FEI Version of the Day
fei_Trilinos_Helpers.hpp
00001 /*
00002 // @HEADER
00003 // ************************************************************************
00004 //             FEI: Finite Element Interface to Linear Solvers
00005 //                  Copyright (2005) Sandia Corporation.
00006 //
00007 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation, the
00008 // U.S. Government retains certain rights in this software.
00009 //
00010 // Redistribution and use in source and binary forms, with or without
00011 // modification, are permitted provided that the following conditions are
00012 // met:
00013 //
00014 // 1. Redistributions of source code must retain the above copyright
00015 // notice, this list of conditions and the following disclaimer.
00016 //
00017 // 2. Redistributions in binary form must reproduce the above copyright
00018 // notice, this list of conditions and the following disclaimer in the
00019 // documentation and/or other materials provided with the distribution.
00020 //
00021 // 3. Neither the name of the Corporation nor the names of the
00022 // contributors may be used to endorse or promote products derived from
00023 // this software without specific prior written permission.
00024 //
00025 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
00026 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00027 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
00028 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
00029 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
00030 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
00031 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
00032 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
00033 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
00034 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
00035 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00036 //
00037 // Questions? Contact Alan Williams (william@sandia.gov) 
00038 //
00039 // ************************************************************************
00040 // @HEADER
00041 */
00042 
00043 
00044 #ifndef _fei_Trilinos_Helpers_hpp_
00045 #define _fei_Trilinos_Helpers_hpp_
00046 
00047 #include "fei_trilinos_macros.hpp"
00048 #include "fei_fwd.hpp"
00049 
00050 #include <fei_Include_Trilinos.hpp>
00051 
00052 #include <fei_mpi.h>
00053 #include <fei_SharedPtr.hpp>
00054 
00055 #include <fei_LinearProblemManager.hpp>
00056 #include <fei_VectorSpace.hpp>
00057 #include <fei_Reducer.hpp>
00058 #include <fei_MatrixGraph.hpp>
00059 
00060 namespace Trilinos_Helpers {
00061 
00062 #ifdef HAVE_FEI_EPETRA
00063 
00069   Epetra_Map create_Epetra_Map(MPI_Comm comm,
00070                                const std::vector<int>& local_eqns);
00071 
00077   Epetra_BlockMap
00078     create_Epetra_BlockMap(const fei::SharedPtr<fei::VectorSpace>& vecspace);
00079 
00080   Epetra_CrsGraph
00081     create_Epetra_CrsGraph(const fei::SharedPtr<fei::MatrixGraph>& matgraph,
00082                            bool blockEntries,
00083                            bool orderRowsWithLocalColsFirst=false);
00084 
00085   fei::SharedPtr<fei::Matrix>
00086     create_from_Epetra_Matrix(fei::SharedPtr<fei::MatrixGraph> matrixGraph,
00087                               bool blockEntryMatrix,
00088                               fei::SharedPtr<fei::Reducer> reducer,
00089                               bool orderRowsWithLocalColsFirst=false);
00090 
00091   fei::SharedPtr<fei::Matrix>
00092     create_from_LPM_EpetraBasic(fei::SharedPtr<fei::MatrixGraph> matrixGraph,
00093                                  bool blockEntryMatrix,
00094                                  fei::SharedPtr<fei::Reducer> reducer,
00095                                  fei::SharedPtr<fei::LinearProblemManager>
00096                                    lpm_epetrabasic);
00097 #endif
00098 
00102   void copy_parameterset(const fei::ParameterSet& paramset,
00103                          Teuchos::ParameterList& paramlist);
00104 
00108   void copy_parameterlist(const Teuchos::ParameterList& paramlist,
00109                           fei::ParameterSet& paramset);
00110 
00111 #ifdef HAVE_FEI_EPETRA
00112 
00115   Epetra_MultiVector*
00116     get_Epetra_MultiVector(fei::Vector* feivec, bool soln_vec);
00117 
00121   Epetra_VbrMatrix* get_Epetra_VbrMatrix(fei::Matrix* feimat);
00122 
00126   Epetra_CrsMatrix* get_Epetra_CrsMatrix(fei::Matrix* feimat);
00127 
00132   void get_Epetra_pointers(fei::SharedPtr<fei::Matrix> feiA,
00133                            fei::SharedPtr<fei::Vector> feix,
00134                            fei::SharedPtr<fei::Vector> feib,
00135                            Epetra_CrsMatrix*& crsA,
00136                            Epetra_Operator*& opA,
00137                            Epetra_MultiVector*& x,
00138                            Epetra_MultiVector*& b);
00139 
00141   int zero_Epetra_VbrMatrix(Epetra_VbrMatrix* mat);
00142 
00146   inline
00147   double* getBeginPointer(const fei::SharedPtr<fei::Matrix>& feiA)
00148   {
00149     Epetra_CrsMatrix* A = get_Epetra_CrsMatrix(feiA.get());
00150     return (*A)[0];
00151   }
00152 
00157   inline
00158   int getOffsetG(const fei::SharedPtr<fei::Matrix>& feiA,
00159                  int global_row_index, int global_col_index)
00160   {
00161     Epetra_CrsMatrix* A = get_Epetra_CrsMatrix(feiA.get());
00162     const Epetra_Map& erowmap = A->RowMap();
00163     const Epetra_Map& ecolmap = A->ColMap();
00164     int local_row = erowmap.LID(global_row_index);
00165     int local_col = ecolmap.LID(global_col_index);
00166   
00167     int* rowOffsets;
00168     int* colIndices;
00169     double* coefs;
00170     A->ExtractCrsDataPointers(rowOffsets, colIndices, coefs);
00171   
00172     int* row_ptr = &colIndices[rowOffsets[local_row]];
00173     int* end_row = &colIndices[rowOffsets[local_row+1]];
00174   
00175     int col_offset = 0;
00176     for(; row_ptr != end_row; ++row_ptr) {
00177       if (*row_ptr == local_col) break;
00178       ++col_offset;
00179     }
00180   
00181     return rowOffsets[local_row] + col_offset;
00182   }
00183 
00188   inline
00189   int getOffsetL(const fei::SharedPtr<fei::Matrix>& feiA,
00190                  int local_row_index, int local_col_index)
00191   {
00192     Epetra_CrsMatrix* A = get_Epetra_CrsMatrix(feiA.get());
00193   
00194     int* rowOffsets;
00195     int* colIndices;
00196     double* coefs;
00197     A->ExtractCrsDataPointers(rowOffsets, colIndices, coefs);
00198   
00199     int* row_ptr = &colIndices[rowOffsets[local_row_index]];
00200     int* end_row = &colIndices[rowOffsets[local_row_index+1]];
00201   
00202     int col_offset = 0;
00203     for(; row_ptr != end_row; ++row_ptr) {
00204       if (*row_ptr == local_col_index) break;
00205       ++col_offset;
00206     }
00207   
00208     return rowOffsets[local_row_index] + col_offset;
00209   }
00210 
00211 #endif // HAVE_FEI_EPETRA
00212 
00213 }//namespace Trilinos_Helpers
00214 
00215 #endif // _Trilinos_Helpers_hpp_
00216 
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Friends