Teko  Version of the Day
 All Classes Files Functions Variables Pages
Teko_StridedMappingStrategy.cpp
1 // @HEADER
2 // *****************************************************************************
3 // Teko: A package for block and physics based preconditioning
4 //
5 // Copyright 2010 NTESS and the Teko contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #include "Epetra/Teko_StridedMappingStrategy.hpp"
11 #include "Epetra/Teko_InterlacedEpetra.hpp"
12 #include "Epetra/Teko_EpetraHelpers.hpp"
13 
14 #include "Thyra_EpetraThyraWrappers.hpp"
15 #include "Thyra_EpetraLinearOp.hpp"
16 #include "Thyra_DefaultProductMultiVector.hpp"
17 #include "Thyra_DefaultProductVectorSpace.hpp"
18 #include "Thyra_DefaultSpmdMultiVector.hpp"
19 #include "Thyra_DefaultBlockedLinearOp.hpp"
20 #include "Thyra_get_Epetra_Operator.hpp"
21 
22 using Teuchos::RCP;
23 using Teuchos::rcp;
24 using Teuchos::rcp_dynamic_cast;
25 
26 namespace Teko {
27 namespace Epetra {
28 
29 // Creates a strided mapping strategy. This class is useful
30 // for breaking up nodally ordered matrices (i.e. the unknowns
31 // in a FEM problem are ordered [u0,v0,p0,u1,v1,p1,...]). Current
32 // implimentation only supports a fixed number of variables
33 //
34 // arguments:
35 // vars - Number of different variables
36 // map - original Epetra_Map to be broken up
37 // comm - Epetra_Comm object related to the map
38 //
39 StridedMappingStrategy::StridedMappingStrategy(const std::vector<int>& vars,
40  const RCP<const Epetra_Map>& map,
41  const Epetra_Comm& comm) {
42  rangeMap_ = map;
43  domainMap_ = map;
44  buildBlockTransferData(vars, rangeMap_, comm);
45 }
46 
47 // Virtual function defined in MappingStrategy. This copies
48 // an Epetra_MultiVector into a Thyra::MultiVectorBase with
49 // blocking handled by the strides defined in the constructor.
50 //
51 // arguments:
52 // X - source Epetra_MultiVector
53 // thyra_X - destination Thyra::MultiVectorBase
54 //
55 void StridedMappingStrategy::copyEpetraIntoThyra(
56  const Epetra_MultiVector& X,
57  const Teuchos::Ptr<Thyra::MultiVectorBase<double> >& thyra_X) const {
58  int count = X.NumVectors();
59 
60  std::vector<RCP<Epetra_MultiVector> > subX;
61 
62  // allocate vectors to copy into
63  Strided::buildSubVectors(blockMaps_, subX, count);
64 
65  // copy source vector to X vector
66  Strided::one2many(subX, X, blockImport_);
67 
68  // convert subX to an array of multi vectors
69  Teuchos::Array<RCP<Thyra::MultiVectorBase<double> > > thyra_subX;
70  Teuchos::Ptr<Thyra::ProductMultiVectorBase<double> > prod_X =
71  Teuchos::ptr_dynamic_cast<Thyra::ProductMultiVectorBase<double> >(thyra_X);
72  for (unsigned int i = 0; i < blockMaps_.size(); i++) {
73  RCP<Thyra::DefaultSpmdMultiVector<double> > vec =
74  rcp_dynamic_cast<Thyra::DefaultSpmdMultiVector<double> >(
75  prod_X->getNonconstMultiVectorBlock(i));
76  fillDefaultSpmdMultiVector(vec, subX[i]);
77  }
78 }
79 
80 // Virtual function defined in MappingStrategy. This copies
81 // an Epetra_MultiVector into a Thyra::MultiVectorBase with
82 // blocking handled by the strides defined in the constructor.
83 //
84 // arguments:
85 // thyra_Y - source Thyra::MultiVectorBase
86 // Y - destination Epetra_MultiVector
87 //
88 void StridedMappingStrategy::copyThyraIntoEpetra(
89  const RCP<const Thyra::MultiVectorBase<double> >& thyra_Y, Epetra_MultiVector& Y) const {
90  std::vector<RCP<const Epetra_MultiVector> > subY;
91  RCP<const Thyra::DefaultProductMultiVector<double> > prod_Y =
92  rcp_dynamic_cast<const Thyra::DefaultProductMultiVector<double> >(thyra_Y);
93 
94  // convert thyra product vector to subY
95  for (unsigned int i = 0; i < blockMaps_.size(); i++)
96  subY.push_back(
97  Thyra::get_Epetra_MultiVector(*blockMaps_[i].second, prod_Y->getMultiVectorBlock(i)));
98 
99  // endow the subVectors with required information about the maps
100  Strided::associateSubVectors(blockMaps_, subY);
101 
102  // copy solution vectors to Y vector
103  Strided::many2one(Y, subY, blockExport_);
104 }
105 
106 // this is the core routine that builds the maps
107 // and importers/exporters neccessary for all the
108 // transfers. Currently it simply calls out to the
109 // interlaced epetra functions. (Comment: this
110 // routine should probably be private or protected
111 // ... it is basically the meat of the constructor)
112 //
113 // arguments:
114 // vars - Vector describing the blocking of variables
115 // baseMap - basic map to use in the transfers
116 // comm - Epetra_Comm object
117 //
118 void StridedMappingStrategy::buildBlockTransferData(const std::vector<int>& vars,
119  const Teuchos::RCP<const Epetra_Map>& baseMap,
120  const Epetra_Comm& comm) {
121  // build maps and exporters/importers
122  Strided::buildSubMaps(*baseMap, vars, comm, blockMaps_);
123  Strided::buildExportImport(*baseMap, blockMaps_, blockExport_, blockImport_);
124 }
125 
126 // Builds a blocked Thyra operator that uses the strided
127 // mapping strategy to define sub blocks.
128 //
129 // arguments:
130 // mat - Epetra_CrsMatrix with FillComplete called, this
131 // matrix is assumed to be square, with the same
132 // range and domain maps
133 // returns: Blocked Thyra linear operator with sub blocks
134 // defined by this mapping strategy
135 //
136 const Teuchos::RCP<Thyra::BlockedLinearOpBase<double> > StridedMappingStrategy::buildBlockedThyraOp(
137  const RCP<const Epetra_CrsMatrix>& crsContent, const std::string& label) const {
138  int dim = blockMaps_.size();
139 
140  RCP<Thyra::DefaultBlockedLinearOp<double> > A = Thyra::defaultBlockedLinearOp<double>();
141 
142  A->beginBlockFill(dim, dim);
143  for (int i = 0; i < dim; i++) {
144  for (int j = 0; j < dim; j++) {
145  // label block correctly
146  std::stringstream ss;
147  ss << label << "_" << i << "," << j;
148 
149  // build the blocks and place it the right location
150  A->setNonconstBlock(i, j,
151  Thyra::nonconstEpetraLinearOp(
152  Strided::buildSubBlock(i, j, *crsContent, blockMaps_), ss.str()));
153  }
154  } // end for i
155  A->endBlockFill();
156 
157  return A;
158 }
159 
160 // Rebuilds a blocked Thyra operator that uses the strided
161 // mapping strategy to define sub blocks.
162 //
163 // arguments:
164 // crsContent - Epetra_CrsMatrix with FillComplete called, this
165 // matrix is assumed to be square, with the same
166 // range and domain maps
167 // A - Destination block linear op composed of blocks of
168 // Epetra_CrsMatrix at all relevant locations
169 //
170 void StridedMappingStrategy::rebuildBlockedThyraOp(
171  const RCP<const Epetra_CrsMatrix>& crsContent,
172  const RCP<Thyra::BlockedLinearOpBase<double> >& A) const {
173  int dim = blockMaps_.size();
174 
175  for (int i = 0; i < dim; i++) {
176  for (int j = 0; j < dim; j++) {
177  // get Epetra version of desired block
178  RCP<Thyra::LinearOpBase<double> > Aij = A->getNonconstBlock(i, j);
179  RCP<Epetra_CrsMatrix> eAij =
180  rcp_dynamic_cast<Epetra_CrsMatrix>(Thyra::get_Epetra_Operator(*Aij), true);
181 
182  // rebuild the blocks and place it the right location
183  Strided::rebuildSubBlock(i, j, *crsContent, blockMaps_, *eAij);
184  }
185  } // end for i
186 }
187 
188 } // end namespace Epetra
189 } // end namespace Teko