Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_STK_SetupLOWSFactory.cpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Panzer: A partial differential equation assembly
4// engine for strongly coupled complex multiphysics systems
5//
6// Copyright 2011 NTESS and the Panzer contributors.
7// SPDX-License-Identifier: BSD-3-Clause
8// *****************************************************************************
9// @HEADER
10
11#include "PanzerAdaptersSTK_config.hpp"
12
16
17#include "Teuchos_AbstractFactoryStd.hpp"
18
19#include "Stratimikos_DefaultLinearSolverBuilder.hpp"
20
21#ifdef PANZER_HAVE_EPETRA_STACK
22#include "Epetra_MpiComm.h"
23#include "Epetra_Vector.h"
24#include "EpetraExt_VectorOut.h"
25#endif // PANZER_HAVE_EPETRA_STACK
26
27#include "Tpetra_Map.hpp"
28#include "Tpetra_MultiVector.hpp"
29
30#ifdef PANZER_HAVE_TEKO
31#include "Teko_StratimikosFactory.hpp"
32#endif
33
34#ifdef PANZER_HAVE_MUELU
35#include "Stratimikos_MueLuHelpers.hpp"
36//#include "MatrixMarket_Tpetra.hpp"
37#include "Xpetra_MapFactory.hpp"
38#include "Xpetra_MultiVectorFactory.hpp"
39#endif
40
41namespace panzer_stk {
42
43namespace {
44
45 bool
46 determineCoordinateField(const panzer::GlobalIndexer & globalIndexer,std::string & fieldName)
47 {
48 std::vector<std::string> elementBlocks;
49 globalIndexer.getElementBlockIds(elementBlocks);
50
51 // grab fields for first block
52 std::set<int> runningFields;
53 {
54 const std::vector<int> & fields = globalIndexer.getBlockFieldNumbers(elementBlocks[0]);
55 runningFields.insert(fields.begin(),fields.end());
56 }
57
58 // grab fields for first block
59 for(std::size_t i=1;i<elementBlocks.size();i++) {
60 const std::vector<int> & fields = globalIndexer.getBlockFieldNumbers(elementBlocks[i]);
61
62 std::set<int> currentFields(runningFields);
63 runningFields.clear();
64 std::set_intersection(fields.begin(),fields.end(),
65 currentFields.begin(),currentFields.end(),
66 std::inserter(runningFields,runningFields.begin()));
67 }
68
69 if(runningFields.size()<1)
70 return false;
71
72 fieldName = globalIndexer.getFieldString(*runningFields.begin());
73 return true;
74 }
75
76 void
77 fillFieldPatternMap(const panzer::DOFManager & globalIndexer,
78 const std::string & fieldName,
79 std::map<std::string,Teuchos::RCP<const panzer::Intrepid2FieldPattern> > & fieldPatterns)
80 {
81 std::vector<std::string> elementBlocks;
82 globalIndexer.getElementBlockIds(elementBlocks);
83
84 for(std::size_t e=0;e<elementBlocks.size();e++) {
85 std::string blockId = elementBlocks[e];
86
87 if(globalIndexer.fieldInBlock(fieldName,blockId))
88 fieldPatterns[blockId] =
89 Teuchos::rcp_dynamic_cast<const panzer::Intrepid2FieldPattern>(globalIndexer.getFieldPattern(blockId,fieldName),true);
90 }
91 }
92
93 void
94 fillFieldPatternMap(const panzer::GlobalIndexer & globalIndexer,
95 const std::string & fieldName,
96 std::map<std::string,Teuchos::RCP<const panzer::Intrepid2FieldPattern> > & fieldPatterns)
97 {
98 using Teuchos::Ptr;
99 using Teuchos::ptrFromRef;
100 using Teuchos::ptr_dynamic_cast;
101 using panzer::DOFManager;
102
103 // first standard dof manager
104 {
105 Ptr<const DOFManager> dofManager = ptr_dynamic_cast<const DOFManager>(ptrFromRef(globalIndexer));
106
107 if(dofManager!=Teuchos::null) {
108 fillFieldPatternMap(*dofManager,fieldName,fieldPatterns);
109 return;
110 }
111 }
112 }
113} // end anonymous namespace
114
115 Teuchos::RCP<Thyra::LinearOpWithSolveFactoryBase<double> >
116 buildLOWSFactory(bool blockedAssembly,
117 const Teuchos::RCP<const panzer::GlobalIndexer> & globalIndexer,
118 const Teuchos::RCP<panzer_stk::STKConnManager> & stkConn_manager,
119 int spatialDim,
120 const Teuchos::RCP<const Teuchos::MpiComm<int> > & mpi_comm,
121 const Teuchos::RCP<Teuchos::ParameterList> & strat_params,
122 #ifdef PANZER_HAVE_TEKO
123 const Teuchos::RCP<Teko::RequestHandler> & reqHandler,
124 #endif
125 bool writeCoordinates,
126 bool writeTopo,
127 const Teuchos::RCP<const panzer::GlobalIndexer> & auxGlobalIndexer,
128 bool useCoordinates
129 )
130 {
131 using Teuchos::RCP;
132 using Teuchos::rcp;
133 using Teuchos::rcp_dynamic_cast;
134
135 Stratimikos::DefaultLinearSolverBuilder linearSolverBuilder;
136
137 // Note if you want to use new solvers within Teko they have to be added to the solver builer
138 // before teko is added. This is because Teko steals its defaults from the solver its being injected
139 // into!
140
141 #ifdef PANZER_HAVE_MUELU
142 {
143 Stratimikos::enableMueLu<double,int,panzer::GlobalOrdinal,panzer::TpetraNodeType>(linearSolverBuilder,"MueLu");
144 Stratimikos::enableMueLuRefMaxwell<double,int,panzer::GlobalOrdinal,panzer::TpetraNodeType>(linearSolverBuilder,"MueLuRefMaxwell");
145 Stratimikos::enableMueLuMaxwell1<double,int,panzer::GlobalOrdinal,panzer::TpetraNodeType>(linearSolverBuilder,"MueLuMaxwell1");
146 #ifndef PANZER_HIDE_DEPRECATED_CODE
147 // the next two are only for backwards compatibility
148 Stratimikos::enableMueLu<double,int,panzer::GlobalOrdinal,panzer::TpetraNodeType>(linearSolverBuilder,"MueLu-Tpetra");
149 Stratimikos::enableMueLuRefMaxwell<double,int,panzer::GlobalOrdinal,panzer::TpetraNodeType>(linearSolverBuilder,"MueLuRefMaxwell-Tpetra");
150 #endif
151 }
152 #endif // MUELU
153
154 #ifdef PANZER_HAVE_TEKO
155 RCP<Teko::RequestHandler> reqHandler_local = reqHandler;
156
157 if(!blockedAssembly) {
158
159 std::string fieldName;
160
161 // try to set request handler from member variable; This is a potential segfault
162 // if its internally stored data (e.g. callback) gets released and all its data
163 // required by ML or whatever gets hosed
164 if(reqHandler_local==Teuchos::null)
165 reqHandler_local = rcp(new Teko::RequestHandler);
166
167 // add in the coordinate parameter list callback handler
168 if(determineCoordinateField(*globalIndexer,fieldName)) {
169 std::map<std::string,RCP<const panzer::Intrepid2FieldPattern> > fieldPatterns;
170 fillFieldPatternMap(*globalIndexer,fieldName,fieldPatterns);
171
172 RCP<panzer_stk::ParameterListCallback> callback = rcp(new
173 panzer_stk::ParameterListCallback(fieldName,fieldPatterns,stkConn_manager,
174 rcp_dynamic_cast<const panzer::GlobalIndexer>(globalIndexer)));
175 reqHandler_local->addRequestCallback(callback);
176
177 // determine if you want rigid body null space modes...currently an extremely specialized case!
178 if(strat_params->sublist("Preconditioner Types").isSublist("ML")) {
179/* COMMENTING THIS OUT FOR NOW, this is causing problems with some of the preconditioners in optimization...not sure why
180
181 Teuchos::ParameterList & ml_params = strat_params->sublist("Preconditioner Types").sublist("ML").sublist("ML Settings");
182
183 {
184 // force parameterlistcallback to build coordinates
185 callback->preRequest(Teko::RequestMesg(rcp(new Teuchos::ParameterList())));
186
187 // extract coordinate vectors
188 std::vector<double> & xcoords = const_cast<std::vector<double> & >(callback->getXCoordsVector());
189 std::vector<double> & ycoords = const_cast<std::vector<double> & >(callback->getYCoordsVector());
190 std::vector<double> & zcoords = const_cast<std::vector<double> & >(callback->getZCoordsVector());
191
192 ml_params.set<double*>("x-coordinates",&xcoords[0]);
193 ml_params.set<double*>("y-coordinates",&ycoords[0]);
194 ml_params.set<double*>("z-coordinates",&zcoords[0]);
195 }
196*/
197/*
198 bool useRigidBodyNullSpace = false;
199 if(ml_params.isType<std::string>("null space: type"))
200 useRigidBodyNullSpace = ml_params.get<std::string>("null space: type") == "pre-computed";
201
202 if(useRigidBodyNullSpace) {
203 // force parameterlistcallback to build coordinates
204 callback->preRequest(Teko::RequestMesg(rcp(new Teuchos::ParameterList())));
205
206 RCP<std::vector<double> > rbm = rcp(new std::vector<double>);
207 std::vector<double> & rbm_ref = *rbm;
208
209 // extract coordinate vectors
210 std::vector<double> & xcoords = const_cast<std::vector<double> & >(callback->getXCoordsVector());
211 std::vector<double> & ycoords = const_cast<std::vector<double> & >(callback->getYCoordsVector());
212 std::vector<double> & zcoords = const_cast<std::vector<double> & >(callback->getZCoordsVector());
213
214 // use ML to build the null space modes for ML
215 int Nnodes = Teuchos::as<int>(xcoords.size());
216 int NscalarDof = 0;
217 int Ndof = spatialDim;
218 int nRBM = spatialDim==3 ? 6 : (spatialDim==2 ? 3 : 1);
219 int rbmSize = Nnodes*(nRBM+NscalarDof)*(Ndof+NscalarDof);
220 rbm_ref.resize(rbmSize);
221
222 ML_Coord2RBM(Nnodes,&xcoords[0],&ycoords[0],&zcoords[0],&rbm_ref[0],Ndof,NscalarDof);
223
224 ml_params.set<double*>("null space: vectors",&rbm_ref[0]);
225 ml_params.set<int>("null space: dimension",nRBM);
226
227 callback->storeExtraVector(rbm);
228 }
229*/
230 }
231
232 if(writeCoordinates) {
233#ifdef PANZER_HAVE_EPETRA_STACK
234 // force parameterlistcallback to build coordinates
235 callback->preRequest(Teko::RequestMesg(rcp(new Teuchos::ParameterList())));
236
237 // extract coordinate vectors
238 const std::vector<double> & xcoords = callback->getXCoordsVector();
239 const std::vector<double> & ycoords = callback->getYCoordsVector();
240 const std::vector<double> & zcoords = callback->getZCoordsVector();
241
242 // use epetra to write coordinates to matrix market files
243 Epetra_MpiComm ep_comm(*mpi_comm->getRawMpiComm()); // this is OK access to RawMpiComm becase its declared on the stack?
244 // and all users of this object are on the stack (within scope of mpi_comm
245 Epetra_Map map(-1,xcoords.size(),0,ep_comm);
246
247 RCP<Epetra_Vector> vec;
248 switch(spatialDim) {
249 case 3:
250 vec = rcp(new Epetra_Vector(Copy,map,const_cast<double *>(&zcoords[0])));
251 EpetraExt::VectorToMatrixMarketFile("zcoords.mm",*vec);
252 // Intentional fall-through.
253 case 2:
254 vec = rcp(new Epetra_Vector(Copy,map,const_cast<double *>(&ycoords[0])));
255 EpetraExt::VectorToMatrixMarketFile("ycoords.mm",*vec);
256 // Intentional fall-through.
257 case 1:
258 vec = rcp(new Epetra_Vector(Copy,map,const_cast<double *>(&xcoords[0])));
259 EpetraExt::VectorToMatrixMarketFile("xcoords.mm",*vec);
260 break;
261 default:
262 TEUCHOS_ASSERT(false);
263 }
264#else
265 TEUCHOS_TEST_FOR_EXCEPTION(true,std::runtime_error,"ERROR: Panzer_STK_SetupLOWSFactory.cpp - writeCoordinates not implemented for Tpetra yet!");
266#endif
267 }
268
269#ifdef PANZER_HAVE_MUELU
270 if(rcp_dynamic_cast<const panzer::GlobalIndexer>(globalIndexer)!=Teuchos::null
271 && useCoordinates) {
272 if(!writeCoordinates)
273 callback->preRequest(Teko::RequestMesg(rcp(new Teuchos::ParameterList())));
274
275 typedef Tpetra::Map<int,panzer::GlobalOrdinal,panzer::TpetraNodeType> Map;
276 typedef Tpetra::MultiVector<double,int,panzer::GlobalOrdinal,panzer::TpetraNodeType> MV;
277
278 // extract coordinate vectors and modify strat_params to include coordinate vectors
279 unsigned dim = Teuchos::as<unsigned>(spatialDim);
280 RCP<MV> coords;
281 for(unsigned d=0;d<dim;d++) {
282 const std::vector<double> & coord = callback->getCoordsVector(d);
283
284 // no coords vector has been build yet, build one
285 if(coords==Teuchos::null) {
286 if(globalIndexer->getNumFields()==1) {
287 RCP<const panzer::GlobalIndexer> ugi
288 = rcp_dynamic_cast<const panzer::GlobalIndexer>(globalIndexer);
289 std::vector<panzer::GlobalOrdinal> ownedIndices;
290 ugi->getOwnedIndices(ownedIndices);
291 RCP<const Map> coords_map = rcp(new Map(Teuchos::OrdinalTraits<panzer::GlobalOrdinal>::invalid(),ownedIndices,0,mpi_comm));
292 coords = rcp(new MV(coords_map,dim));
293 }
294 else {
295 RCP<const Map> coords_map = rcp(new Map(Teuchos::OrdinalTraits<panzer::GlobalOrdinal>::invalid(),coord.size(),0,mpi_comm));
296 coords = rcp(new MV(coords_map,dim));
297 }
298 }
299
300 // sanity check the size
301 TEUCHOS_ASSERT(coords->getLocalLength()==coord.size());
302
303 // fill appropriate coords vector
304 Teuchos::ArrayRCP<double> dest = coords->getDataNonConst(d);
305 for(std::size_t i=0;i<coord.size();i++)
306 dest[i] = coord[i];
307 }
308
309 // inject coordinates into parameter list
310 Teuchos::ParameterList & muelu_params = strat_params->sublist("Preconditioner Types").sublist("MueLu");
311 muelu_params.set<RCP<MV> >("Coordinates",coords);
312 }
313 #endif
314 }
315 // else write_out_the_mesg("Warning: No unique field determines the coordinates, coordinates unavailable!")
316
317 Teko::addTekoToStratimikosBuilder(linearSolverBuilder,reqHandler_local);
318 }
319 else {
320 // try to set request handler from member variable
321 if(reqHandler_local==Teuchos::null)
322 reqHandler_local = rcp(new Teko::RequestHandler);
323
324 std::string fieldName;
325 if(determineCoordinateField(*globalIndexer,fieldName)) {
326 RCP<const panzer::BlockedDOFManager> blkDofs =
327 rcp_dynamic_cast<const panzer::BlockedDOFManager>(globalIndexer);
328 RCP<const panzer::BlockedDOFManager> auxBlkDofs =
329 rcp_dynamic_cast<const panzer::BlockedDOFManager>(auxGlobalIndexer);
330 RCP<panzer_stk::ParameterListCallbackBlocked> callback =
331 rcp(new panzer_stk::ParameterListCallbackBlocked(stkConn_manager,blkDofs,auxBlkDofs));
332 reqHandler_local->addRequestCallback(callback);
333 }
334
335 Teko::addTekoToStratimikosBuilder(linearSolverBuilder,reqHandler_local);
336
337 if(writeCoordinates) {
338#ifdef PANZER_HAVE_EPETRA_STACK
339 RCP<const panzer::BlockedDOFManager> blkDofs =
340 rcp_dynamic_cast<const panzer::BlockedDOFManager>(globalIndexer);
341
342 // loop over blocks
343 const std::vector<RCP<panzer::GlobalIndexer>> & dofVec
344 = blkDofs->getFieldDOFManagers();
345 for(std::size_t i=0;i<dofVec.size();i++) {
346
347 // add in the coordinate parameter list callback handler
348 TEUCHOS_ASSERT(determineCoordinateField(*dofVec[i],fieldName));
349
350 std::map<std::string,RCP<const panzer::Intrepid2FieldPattern> > fieldPatterns;
351 fillFieldPatternMap(*dofVec[i],fieldName,fieldPatterns);
352 panzer_stk::ParameterListCallback plCall(fieldName,fieldPatterns,stkConn_manager,dofVec[i]);
353 plCall.buildArrayToVector();
354 plCall.buildCoordinates();
355
356 // extract coordinate vectors
357 const std::vector<double> & xcoords = plCall.getXCoordsVector();
358 const std::vector<double> & ycoords = plCall.getYCoordsVector();
359 const std::vector<double> & zcoords = plCall.getZCoordsVector();
360
361 // use epetra to write coordinates to matrix market files
362 Epetra_MpiComm ep_comm(*mpi_comm->getRawMpiComm()); // this is OK access to RawMpiComm becase its declared on the stack?
363 // and all users of this object are on the stack (within scope of mpi_comm
364 Epetra_Map map(-1,xcoords.size(),0,ep_comm);
365
366 RCP<Epetra_Vector> vec;
367 switch(spatialDim) {
368 case 3:
369 vec = rcp(new Epetra_Vector(Copy,map,const_cast<double *>(&zcoords[0])));
370 EpetraExt::VectorToMatrixMarketFile((fieldName+"_zcoords.mm").c_str(),*vec);
371 // Intentional fall-through.
372 case 2:
373 vec = rcp(new Epetra_Vector(Copy,map,const_cast<double *>(&ycoords[0])));
374 EpetraExt::VectorToMatrixMarketFile((fieldName+"_ycoords.mm").c_str(),*vec);
375 // Intentional fall-through.
376 case 1:
377 vec = rcp(new Epetra_Vector(Copy,map,const_cast<double *>(&xcoords[0])));
378 EpetraExt::VectorToMatrixMarketFile((fieldName+"_xcoords.mm").c_str(),*vec);
379 break;
380 default:
381 TEUCHOS_ASSERT(false);
382 }
383
384 // TODO add MueLu code...
385 #ifdef PANZER_HAVE_MUELU
386 if(useCoordinates) {
387
388 typedef Xpetra::Map<int,panzer::GlobalOrdinal> Map;
389 typedef Xpetra::MultiVector<double,int,panzer::GlobalOrdinal> MV;
390
391 // TODO This is Epetra-specific
392 RCP<const Map> coords_map = Xpetra::MapFactory<int,panzer::GlobalOrdinal>::Build(Xpetra::UseEpetra,
393 Teuchos::OrdinalTraits<panzer::GlobalOrdinal>::invalid(),
394 //Teuchos::ArrayView<GO>(ownedIndices),
395 xcoords.size(),
396 0,
397 mpi_comm
398 );
399
400 unsigned dim = Teuchos::as<unsigned>(spatialDim);
401
402 RCP<MV> coords = Xpetra::MultiVectorFactory<double,int,panzer::GlobalOrdinal>::Build(coords_map,spatialDim);
403
404 for(unsigned d=0;d<dim;d++) {
405 // sanity check the size
406 TEUCHOS_ASSERT(coords->getLocalLength()==xcoords.size());
407
408 // fill appropriate coords vector
409 Teuchos::ArrayRCP<double> dest = coords->getDataNonConst(d);
410 for(std::size_t j=0;j<coords->getLocalLength();++j) {
411 if (d == 0) dest[j] = xcoords[j];
412 if (d == 1) dest[j] = ycoords[j];
413 if (d == 2) dest[j] = zcoords[j];
414 }
415 }
416
417 // TODO This is Epetra-specific
418 // inject coordinates into parameter list
419 Teuchos::ParameterList & muelu_params = strat_params->sublist("Preconditioner Types").sublist("MueLu");
420 muelu_params.set<RCP<MV> >("Coordinates",coords);
421
422 }
423 #endif
424
425 } /* end loop over all physical fields */
426#else // PANZER_HAVE EPETRA
427 TEUCHOS_TEST_FOR_EXCEPTION(true,std::runtime_error,"ERROR: Panzer_STK_SetupLOWSFactory - writeCoordinates not implemented for Tpetra yet!")
428#endif
429 }
430
431 if(writeTopo) {
432 /*
433 RCP<const panzer::BlockedDOFManager> blkDofs =
434 rcp_dynamic_cast<const panzer::BlockedDOFManager>(globalIndexer);
435
436 writeTopology(*blkDofs);
437 */
438 TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,
439 "Topology writing is no longer implemented. It needs to be reimplemented for the "
440 "default DOFManager (estimate 2 days with testing)");
441 }
442 }
443 #endif
444
445 linearSolverBuilder.setParameterList(strat_params);
446 RCP<Thyra::LinearOpWithSolveFactoryBase<double> > lowsFactory = createLinearSolveStrategy(linearSolverBuilder);
447
448 return lowsFactory;
449 }
450
451
452 Teuchos::RCP<Thyra::LinearOpWithSolveFactoryBase<double> >
453 buildLOWSFactory(bool blockedAssembly,
454 const Teuchos::RCP<const panzer::GlobalIndexer> & globalIndexer,
455 const Teuchos::RCP<panzer::ConnManager> & conn_manager,
456 int spatialDim,
457 const Teuchos::RCP<const Teuchos::MpiComm<int> > & mpi_comm,
458 const Teuchos::RCP<Teuchos::ParameterList> & strat_params,
459 #ifdef PANZER_HAVE_TEKO
460 const Teuchos::RCP<Teko::RequestHandler> & reqHandler,
461 #endif
462 bool writeCoordinates,
463 bool writeTopo,
464 const Teuchos::RCP<const panzer::GlobalIndexer> & auxGlobalIndexer,
465 bool useCoordinates
466 )
467 {
468 #ifdef PANZER_HAVE_TEKO
469 Teuchos::RCP<Teko::RequestHandler> reqHandler_local = reqHandler;
470 if(reqHandler_local==Teuchos::null)
471 reqHandler_local = Teuchos::rcp(new Teko::RequestHandler);
472 #endif
473
474 auto stk_conn_manager = Teuchos::rcp_dynamic_cast<panzer_stk::STKConnManager>(conn_manager,true);
475
476 return buildLOWSFactory(blockedAssembly,globalIndexer,stk_conn_manager,spatialDim,mpi_comm,strat_params,
477#ifdef PANZER_HAVE_TEKO
478 reqHandler_local,
479#endif
480 writeCoordinates,
481 writeTopo,
482 auxGlobalIndexer,
483 useCoordinates
484 );
485
486 // should never reach this
487 TEUCHOS_ASSERT(false);
488 return Teuchos::null;
489 }
490}
Copy
Teuchos::RCP< const FieldPattern > getFieldPattern(const std::string &name) const
Find a field pattern stored for a particular block and field number. This will retrive the pattern ad...
void getElementBlockIds(std::vector< std::string > &elementBlockIds) const
bool fieldInBlock(const std::string &field, const std::string &block) const
virtual const std::string & getFieldString(int num) const =0
Reverse lookup of the field string from a field number.
virtual const std::vector< int > & getBlockFieldNumbers(const std::string &blockId) const =0
virtual void getElementBlockIds(std::vector< std::string > &elementBlockIds) const =0
int VectorToMatrixMarketFile(const char *filename, const Epetra_Vector &A, const char *matrixName=0, const char *matrixDescription=0, bool writeHeader=true)
Teuchos::RCP< Thyra::LinearOpWithSolveFactoryBase< double > > buildLOWSFactory(bool blockedAssembly, const Teuchos::RCP< const panzer::GlobalIndexer > &globalIndexer, const Teuchos::RCP< panzer_stk::STKConnManager > &stkConn_manager, int spatialDim, const Teuchos::RCP< const Teuchos::MpiComm< int > > &mpi_comm, const Teuchos::RCP< Teuchos::ParameterList > &strat_params, bool writeCoordinates, bool writeTopo, const Teuchos::RCP< const panzer::GlobalIndexer > &auxGlobalIndexer, bool useCoordinates)