Amesos2 - Direct Sparse Solver Interfaces Version of the Day
Amesos2_EpetraCrsMatrix_MatrixAdapter_def.hpp
1// @HEADER
2// *****************************************************************************
3// Amesos2: Templated Direct Sparse Solver Package
4//
5// Copyright 2011 NTESS and the Amesos2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10
11#ifndef AMESOS2_EPETRACRSMATRIX_MATRIXADAPTER_DEF_HPP
12#define AMESOS2_EPETRACRSMATRIX_MATRIXADAPTER_DEF_HPP
13
14#include <Epetra_LocalMap.h>
15#include <Epetra_Import.h>
16
19#include "Amesos2_MatrixAdapter_def.hpp"
20
21
22namespace Amesos2 {
23
24 ConcreteMatrixAdapter<Epetra_CrsMatrix>::ConcreteMatrixAdapter(RCP<Epetra_CrsMatrix> m)
25 : AbstractConcreteMatrixAdapter<Epetra_RowMatrix,Epetra_CrsMatrix>(m) // CrsMatrix inherits from RowMatrix virtually, so a dynamic cast is necessary
26 {}
27
28 Teuchos::RCP<const MatrixAdapter<Epetra_CrsMatrix> >
29 ConcreteMatrixAdapter<Epetra_CrsMatrix>::get_impl(const Teuchos::Ptr<const map_t> map, EDistribution distribution) const
30 {
31 using Teuchos::as;
32 using Teuchos::rcp;
33 using Teuchos::rcpFromRef;
34
35 RCP<const Epetra_Map> o_map, t_map;
36 o_map = rcpFromRef(this->mat_->RowMap());
37 t_map = Util::tpetra_map_to_epetra_map<local_ordinal_t,global_ordinal_t,global_size_t,node_t>(*map);
38
39 const int maxRowNNZ = 0;
40 RCP<Epetra_CrsMatrix> t_mat = rcp(new Epetra_CrsMatrix(Copy, *t_map, maxRowNNZ));
41
42 Epetra_Import importer(*t_map, *o_map);
43 t_mat->Import(*(this->mat_), importer, Insert);
44 t_mat->FillComplete();
45
46 // Case for non-contiguous GIDs
47 if ( distribution == CONTIGUOUS_AND_ROOTED ) {
48
49 auto myRank = map->getComm()->getRank();
50
51 const int global_num_contiguous_entries = t_mat->NumGlobalRows();
52 const int local_num_contiguous_entries = (myRank == 0) ? t_mat->NumGlobalRows() : 0;
53
54 RCP<const Epetra_Map> contiguousRowMap = rcp( new Epetra_Map(global_num_contiguous_entries, local_num_contiguous_entries, 0, (t_mat->Comm() ) ) );
55 RCP<const Epetra_Map> contiguousColMap = rcp( new Epetra_Map(global_num_contiguous_entries, local_num_contiguous_entries, 0, (t_mat->Comm() ) ) );
56 RCP<const Epetra_Map> contiguousDomainMap = rcp( new Epetra_Map(global_num_contiguous_entries, local_num_contiguous_entries, 0, (t_mat->Comm() ) ) );
57 RCP<const Epetra_Map> contiguousRangeMap = rcp( new Epetra_Map(global_num_contiguous_entries, local_num_contiguous_entries, 0, (t_mat->Comm() ) ) );
58
59 RCP<Epetra_CrsMatrix> contiguous_t_mat = rcp( new Epetra_CrsMatrix(Epetra_DataAccess::Copy, *contiguousRowMap, *contiguousColMap, t_mat->MaxNumEntries()) );
60
61 // fill local sparse matrix on rank zero
62 if(myRank == 0) {
63 int num_entries;
64 int *indices;
65 double *values;
66 for (int row = 0; row < t_mat->NumMyRows(); row++) {
67 t_mat->ExtractMyRowView(row, num_entries, values, indices);
68 contiguous_t_mat->InsertMyValues(row, num_entries, values, indices);
69 }
70 }
71
72 contiguous_t_mat->FillComplete(*contiguousDomainMap, *contiguousRangeMap);
73
74 return rcp (new ConcreteMatrixAdapter<Epetra_CrsMatrix> (contiguous_t_mat));
75 }
76
77 return( rcp(new ConcreteMatrixAdapter<Epetra_CrsMatrix>(t_mat)) );
78 }
79
80 Teuchos::RCP<const MatrixAdapter<Epetra_CrsMatrix> >
81 ConcreteMatrixAdapter<Epetra_CrsMatrix>::reindex_impl(Teuchos::RCP<const map_t> &contigRowMap,
82 Teuchos::RCP<const map_t> &contigColMap,
83 const EPhase /*current_phase*/) const
84 {
85 #if defined(HAVE_AMESOS2_EPETRAEXT)
86 using Teuchos::RCP;
87 using Teuchos::rcp;
88 using Teuchos::rcpFromRef;
89 auto CrsMatrix = const_cast<Epetra_CrsMatrix *>(this->mat_.getRawPtr());
90 if(!CrsMatrix) {
91 TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Amesos2_EpetraCrsMatrix_MatrixAdapter requires CsrMatrix to reindex matrices.");
92 }
93
94 // Map
95 RCP<const Epetra_Map> OriginalMap = rcpFromRef(CrsMatrix->RowMap());
96 int NumGlobalElements = OriginalMap->NumGlobalElements();
97 int NumMyElements = OriginalMap->NumMyElements();
98 auto ReindexMap = rcp( new Epetra_Map( NumGlobalElements, NumMyElements, 0, OriginalMap->Comm() ) );
99
100 // Matrix
101 StdIndex_ = rcp( new EpetraExt::CrsMatrix_Reindex( *ReindexMap ) );
102 ContigMat_ = rcpFromRef((*StdIndex_)( *CrsMatrix ));
103 if(!ContigMat_) {
104 TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "Amesos2_EpetraCrsMatrix_MatrixAdapter reindexing failed.");
105 }
106
107 auto reindexMat = rcp( new ConcreteMatrixAdapter<Epetra_CrsMatrix>(ContigMat_));
108 contigRowMap = reindexMat->getRowMap();
109 contigColMap = reindexMat->getColMap();
110
111 return reindexMat;
112 #else
113 TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error, "ConcreteMatrixAdapter<Epetra_CrsMatrix> requires EpetraExt to reindex matrices.");
114 #endif
115 }
116
117
118 void
119 ConcreteMatrixAdapter<Epetra_CrsMatrix>::describe (Teuchos::FancyOStream& os,
120 const Teuchos::EVerbosityLevel verbLevel) const
121 {
122 this->mat_->Print(*(os.getOStream()));
123 }
124} // end namespace Amesos2
125
126#endif // AMESOS2_EPETRACRSMATRIX_MATRIXADAPTER_DEF_HPP
Specialization of the ConcreteMatrixAdapter for Epetra_CrsMatrix. Inherits all its functionality from...
Definitions for the Epetra_RowMatrix abstract adapter.