44#ifndef _fei_Trilinos_Helpers_hpp_
45#define _fei_Trilinos_Helpers_hpp_
47#include "fei_trilinos_macros.hpp"
50#include <fei_Include_Trilinos.hpp>
53#include <fei_SharedPtr.hpp>
55#include <fei_LinearProblemManager.hpp>
56#include <fei_VectorSpace.hpp>
57#include <fei_Reducer.hpp>
58#include <fei_MatrixGraph.hpp>
60namespace Trilinos_Helpers {
69 Epetra_Map create_Epetra_Map(MPI_Comm comm,
70 const std::vector<int>& local_eqns);
78 create_Epetra_BlockMap(
const fei::SharedPtr<fei::VectorSpace>& vecspace);
81 create_Epetra_CrsGraph(
const fei::SharedPtr<fei::MatrixGraph>& matgraph,
83 bool orderRowsWithLocalColsFirst=
false);
85 fei::SharedPtr<fei::Matrix>
86 create_from_Epetra_Matrix(fei::SharedPtr<fei::MatrixGraph> matrixGraph,
87 bool blockEntryMatrix,
88 fei::SharedPtr<fei::Reducer> reducer,
89 bool orderRowsWithLocalColsFirst=
false);
91 fei::SharedPtr<fei::Matrix>
92 create_from_LPM_EpetraBasic(fei::SharedPtr<fei::MatrixGraph> matrixGraph,
93 bool blockEntryMatrix,
94 fei::SharedPtr<fei::Reducer> reducer,
95 fei::SharedPtr<fei::LinearProblemManager>
102 void copy_parameterset(
const fei::ParameterSet& paramset,
103 Teuchos::ParameterList& paramlist);
108 void copy_parameterlist(
const Teuchos::ParameterList& paramlist,
109 fei::ParameterSet& paramset);
111#ifdef HAVE_FEI_EPETRA
116 get_Epetra_MultiVector(fei::Vector* feivec,
bool soln_vec);
121 Epetra_VbrMatrix* get_Epetra_VbrMatrix(fei::Matrix* feimat);
126 Epetra_CrsMatrix* get_Epetra_CrsMatrix(fei::Matrix* feimat);
132 void get_Epetra_pointers(fei::SharedPtr<fei::Matrix> feiA,
133 fei::SharedPtr<fei::Vector> feix,
134 fei::SharedPtr<fei::Vector> feib,
135 Epetra_CrsMatrix*& crsA,
136 Epetra_Operator*& opA,
137 Epetra_MultiVector*& x,
138 Epetra_MultiVector*& b);
141 int zero_Epetra_VbrMatrix(Epetra_VbrMatrix* mat);
147 double* getBeginPointer(
const fei::SharedPtr<fei::Matrix>& feiA)
149 Epetra_CrsMatrix* A = get_Epetra_CrsMatrix(feiA.
get());
158 int getOffsetG(
const fei::SharedPtr<fei::Matrix>& feiA,
159 int global_row_index,
int global_col_index)
161 Epetra_CrsMatrix* A = get_Epetra_CrsMatrix(feiA.
get());
162 const Epetra_Map& erowmap = A->
RowMap();
163 const Epetra_Map& ecolmap = A->
ColMap();
164 int local_row = erowmap.
LID(global_row_index);
165 int local_col = ecolmap.
LID(global_col_index);
172 int* row_ptr = &colIndices[rowOffsets[local_row]];
173 int* end_row = &colIndices[rowOffsets[local_row+1]];
176 for(; row_ptr != end_row; ++row_ptr) {
177 if (*row_ptr == local_col)
break;
181 return rowOffsets[local_row] + col_offset;
189 int getOffsetL(
const fei::SharedPtr<fei::Matrix>& feiA,
190 int local_row_index,
int local_col_index)
192 Epetra_CrsMatrix* A = get_Epetra_CrsMatrix(feiA.
get());
199 int* row_ptr = &colIndices[rowOffsets[local_row_index]];
200 int* end_row = &colIndices[rowOffsets[local_row_index+1]];
203 for(; row_ptr != end_row; ++row_ptr) {
204 if (*row_ptr == local_col_index)
break;
208 return rowOffsets[local_row_index] + col_offset;
const Epetra_Map & RowMap() const
int ExtractCrsDataPointers(int *&IndexOffset, int *&Indices, double *&Values_in) const
const Epetra_Map & ColMap() const