Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_STK_SetupLOWSFactory.cpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Panzer: A partial differential equation assembly
5 // engine for strongly coupled complex multiphysics systems
6 // Copyright (2011) Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Roger P. Pawlowski (rppawlo@sandia.gov) and
39 // Eric C. Cyr (eccyr@sandia.gov)
40 // ***********************************************************************
41 // @HEADER
42 
43 #include "PanzerAdaptersSTK_config.hpp"
47 
48 #include "Teuchos_AbstractFactoryStd.hpp"
49 
50 #include "Stratimikos_DefaultLinearSolverBuilder.hpp"
51 
52 #include "Epetra_MpiComm.h"
53 #include "Epetra_Vector.h"
54 #include "EpetraExt_VectorOut.h"
55 
56 #include "ml_rbm.h"
57 
58 #include "Tpetra_Map.hpp"
59 #include "Tpetra_MultiVector.hpp"
60 
61 #ifdef PANZER_HAVE_TEKO
62 #include "Teko_StratimikosFactory.hpp"
63 #endif
64 
65 #ifdef PANZER_HAVE_MUELU
66 #include <Thyra_MueLuPreconditionerFactory.hpp>
67 #include <Thyra_MueLuRefMaxwellPreconditionerFactory.hpp>
68 #include "Stratimikos_MueLuHelpers.hpp"
69 //#include "MatrixMarket_Tpetra.hpp"
70 #include "Xpetra_MapFactory.hpp"
71 #include "Xpetra_MultiVectorFactory.hpp"
72 #endif
73 
74 #ifdef PANZER_HAVE_IFPACK2
75 #include <Thyra_Ifpack2PreconditionerFactory.hpp>
76 #endif
77 
78 namespace panzer_stk {
79 
80 namespace {
81 
82  bool
83  determineCoordinateField(const panzer::GlobalIndexer & globalIndexer,std::string & fieldName)
84  {
85  std::vector<std::string> elementBlocks;
86  globalIndexer.getElementBlockIds(elementBlocks);
87 
88  // grab fields for first block
89  std::set<int> runningFields;
90  {
91  const std::vector<int> & fields = globalIndexer.getBlockFieldNumbers(elementBlocks[0]);
92  runningFields.insert(fields.begin(),fields.end());
93  }
94 
95  // grab fields for first block
96  for(std::size_t i=1;i<elementBlocks.size();i++) {
97  const std::vector<int> & fields = globalIndexer.getBlockFieldNumbers(elementBlocks[i]);
98 
99  std::set<int> currentFields(runningFields);
100  runningFields.clear();
101  std::set_intersection(fields.begin(),fields.end(),
102  currentFields.begin(),currentFields.end(),
103  std::inserter(runningFields,runningFields.begin()));
104  }
105 
106  if(runningFields.size()<1)
107  return false;
108 
109  fieldName = globalIndexer.getFieldString(*runningFields.begin());
110  return true;
111  }
112 
113  void
114  fillFieldPatternMap(const panzer::DOFManager & globalIndexer,
115  const std::string & fieldName,
116  std::map<std::string,Teuchos::RCP<const panzer::Intrepid2FieldPattern> > & fieldPatterns)
117  {
118  std::vector<std::string> elementBlocks;
119  globalIndexer.getElementBlockIds(elementBlocks);
120 
121  for(std::size_t e=0;e<elementBlocks.size();e++) {
122  std::string blockId = elementBlocks[e];
123 
124  if(globalIndexer.fieldInBlock(fieldName,blockId))
125  fieldPatterns[blockId] =
126  Teuchos::rcp_dynamic_cast<const panzer::Intrepid2FieldPattern>(globalIndexer.getFieldPattern(blockId,fieldName),true);
127  }
128  }
129 
130  void
131  fillFieldPatternMap(const panzer::GlobalIndexer & globalIndexer,
132  const std::string & fieldName,
133  std::map<std::string,Teuchos::RCP<const panzer::Intrepid2FieldPattern> > & fieldPatterns)
134  {
135  using Teuchos::Ptr;
136  using Teuchos::ptrFromRef;
137  using Teuchos::ptr_dynamic_cast;
138  using panzer::DOFManager;
139 
140  // first standard dof manager
141  {
142  Ptr<const DOFManager> dofManager = ptr_dynamic_cast<const DOFManager>(ptrFromRef(globalIndexer));
143 
144  if(dofManager!=Teuchos::null) {
145  fillFieldPatternMap(*dofManager,fieldName,fieldPatterns);
146  return;
147  }
148  }
149  }
150 } // end anonymous namespace
151 
153  buildLOWSFactory(bool blockedAssembly,
154  const Teuchos::RCP<const panzer::GlobalIndexer> & globalIndexer,
155  const Teuchos::RCP<panzer_stk::STKConnManager> & stkConn_manager,
156  int spatialDim,
157  const Teuchos::RCP<const Teuchos::MpiComm<int> > & mpi_comm,
158  const Teuchos::RCP<Teuchos::ParameterList> & strat_params,
159  #ifdef PANZER_HAVE_TEKO
160  const Teuchos::RCP<Teko::RequestHandler> & reqHandler,
161  #endif
162  bool writeCoordinates,
163  bool writeTopo,
164  const Teuchos::RCP<const panzer::GlobalIndexer> & auxGlobalIndexer,
165  bool useCoordinates
166  )
167  {
168  using Teuchos::RCP;
169  using Teuchos::rcp;
170  using Teuchos::rcp_dynamic_cast;
171 
172  Stratimikos::DefaultLinearSolverBuilder linearSolverBuilder;
173 
174  // Note if you want to use new solvers within Teko they have to be added to the solver builer
175  // before teko is added. This is because Teko steals its defaults from the solver its being injected
176  // into!
177 
178  #ifdef PANZER_HAVE_MUELU
179  {
180  Stratimikos::enableMueLu<int,panzer::GlobalOrdinal,panzer::TpetraNodeType>(linearSolverBuilder,"MueLu");
181  Stratimikos::enableMueLuRefMaxwell<int,panzer::GlobalOrdinal,panzer::TpetraNodeType>(linearSolverBuilder,"MueLuRefMaxwell");
182  #ifndef PANZER_HIDE_DEPRECATED_CODE
183  // the next two are only for backwards compatibility
184  Stratimikos::enableMueLu<int,panzer::GlobalOrdinal,panzer::TpetraNodeType>(linearSolverBuilder,"MueLu-Tpetra");
185  Stratimikos::enableMueLuRefMaxwell<int,panzer::GlobalOrdinal,panzer::TpetraNodeType>(linearSolverBuilder,"MueLuRefMaxwell-Tpetra");
186  #endif
187  }
188  #endif // MUELU
189  #ifdef PANZER_HAVE_IFPACK2
190  {
192  typedef Thyra::Ifpack2PreconditionerFactory<Tpetra::CrsMatrix<double, int, panzer::GlobalOrdinal,panzer::TpetraNodeType> > Impl;
193 
194  linearSolverBuilder.setPreconditioningStrategyFactory(Teuchos::abstractFactoryStd<Base, Impl>(), "Ifpack2");
195  }
196  #endif // MUELU
197 
198 
199  #ifdef PANZER_HAVE_TEKO
200  RCP<Teko::RequestHandler> reqHandler_local = reqHandler;
201 
202  if(!blockedAssembly) {
203 
204  std::string fieldName;
205 
206  // try to set request handler from member variable; This is a potential segfault
207  // if its internally stored data (e.g. callback) gets released and all its data
208  // required by ML or whatever gets hosed
209  if(reqHandler_local==Teuchos::null)
210  reqHandler_local = rcp(new Teko::RequestHandler);
211 
212  // add in the coordinate parameter list callback handler
213  if(determineCoordinateField(*globalIndexer,fieldName)) {
214  std::map<std::string,RCP<const panzer::Intrepid2FieldPattern> > fieldPatterns;
215  fillFieldPatternMap(*globalIndexer,fieldName,fieldPatterns);
216 
218  panzer_stk::ParameterListCallback(fieldName,fieldPatterns,stkConn_manager,
219  rcp_dynamic_cast<const panzer::GlobalIndexer>(globalIndexer)));
220  reqHandler_local->addRequestCallback(callback);
221 
222  // determine if you want rigid body null space modes...currently an extremely specialized case!
223  if(strat_params->sublist("Preconditioner Types").isSublist("ML")) {
224 /* COMMENTING THIS OUT FOR NOW, this is causing problems with some of the preconditioners in optimization...not sure why
225 
226  Teuchos::ParameterList & ml_params = strat_params->sublist("Preconditioner Types").sublist("ML").sublist("ML Settings");
227 
228  {
229  // force parameterlistcallback to build coordinates
230  callback->preRequest(Teko::RequestMesg(rcp(new Teuchos::ParameterList())));
231 
232  // extract coordinate vectors
233  std::vector<double> & xcoords = const_cast<std::vector<double> & >(callback->getXCoordsVector());
234  std::vector<double> & ycoords = const_cast<std::vector<double> & >(callback->getYCoordsVector());
235  std::vector<double> & zcoords = const_cast<std::vector<double> & >(callback->getZCoordsVector());
236 
237  ml_params.set<double*>("x-coordinates",&xcoords[0]);
238  ml_params.set<double*>("y-coordinates",&ycoords[0]);
239  ml_params.set<double*>("z-coordinates",&zcoords[0]);
240  }
241 */
242 /*
243  bool useRigidBodyNullSpace = false;
244  if(ml_params.isType<std::string>("null space: type"))
245  useRigidBodyNullSpace = ml_params.get<std::string>("null space: type") == "pre-computed";
246 
247  if(useRigidBodyNullSpace) {
248  // force parameterlistcallback to build coordinates
249  callback->preRequest(Teko::RequestMesg(rcp(new Teuchos::ParameterList())));
250 
251  RCP<std::vector<double> > rbm = rcp(new std::vector<double>);
252  std::vector<double> & rbm_ref = *rbm;
253 
254  // extract coordinate vectors
255  std::vector<double> & xcoords = const_cast<std::vector<double> & >(callback->getXCoordsVector());
256  std::vector<double> & ycoords = const_cast<std::vector<double> & >(callback->getYCoordsVector());
257  std::vector<double> & zcoords = const_cast<std::vector<double> & >(callback->getZCoordsVector());
258 
259  // use ML to build the null space modes for ML
260  int Nnodes = Teuchos::as<int>(xcoords.size());
261  int NscalarDof = 0;
262  int Ndof = spatialDim;
263  int nRBM = spatialDim==3 ? 6 : (spatialDim==2 ? 3 : 1);
264  int rbmSize = Nnodes*(nRBM+NscalarDof)*(Ndof+NscalarDof);
265  rbm_ref.resize(rbmSize);
266 
267  ML_Coord2RBM(Nnodes,&xcoords[0],&ycoords[0],&zcoords[0],&rbm_ref[0],Ndof,NscalarDof);
268 
269  ml_params.set<double*>("null space: vectors",&rbm_ref[0]);
270  ml_params.set<int>("null space: dimension",nRBM);
271 
272  callback->storeExtraVector(rbm);
273  }
274 */
275  }
276 
277  if(writeCoordinates) {
278  // force parameterlistcallback to build coordinates
279  callback->preRequest(Teko::RequestMesg(rcp(new Teuchos::ParameterList())));
280 
281  // extract coordinate vectors
282  const std::vector<double> & xcoords = callback->getXCoordsVector();
283  const std::vector<double> & ycoords = callback->getYCoordsVector();
284  const std::vector<double> & zcoords = callback->getZCoordsVector();
285 
286  // use epetra to write coordinates to matrix market files
287  Epetra_MpiComm ep_comm(*mpi_comm->getRawMpiComm()); // this is OK access to RawMpiComm becase its declared on the stack?
288  // and all users of this object are on the stack (within scope of mpi_comm
289  Epetra_Map map(-1,xcoords.size(),0,ep_comm);
290 
291  RCP<Epetra_Vector> vec;
292  switch(spatialDim) {
293  case 3:
294  vec = rcp(new Epetra_Vector(Copy,map,const_cast<double *>(&zcoords[0])));
295  EpetraExt::VectorToMatrixMarketFile("zcoords.mm",*vec);
296  // Intentional fall-through.
297  case 2:
298  vec = rcp(new Epetra_Vector(Copy,map,const_cast<double *>(&ycoords[0])));
299  EpetraExt::VectorToMatrixMarketFile("ycoords.mm",*vec);
300  // Intentional fall-through.
301  case 1:
302  vec = rcp(new Epetra_Vector(Copy,map,const_cast<double *>(&xcoords[0])));
303  EpetraExt::VectorToMatrixMarketFile("xcoords.mm",*vec);
304  break;
305  default:
306  TEUCHOS_ASSERT(false);
307  }
308  }
309 
310  #ifdef PANZER_HAVE_MUELU
311  if(rcp_dynamic_cast<const panzer::GlobalIndexer>(globalIndexer)!=Teuchos::null
312  && useCoordinates) {
313  if(!writeCoordinates)
314  callback->preRequest(Teko::RequestMesg(rcp(new Teuchos::ParameterList())));
315 
316  typedef Tpetra::Map<int,panzer::GlobalOrdinal,panzer::TpetraNodeType> Map;
317  typedef Tpetra::MultiVector<double,int,panzer::GlobalOrdinal,panzer::TpetraNodeType> MV;
318 
319  // extract coordinate vectors and modify strat_params to include coordinate vectors
320  unsigned dim = Teuchos::as<unsigned>(spatialDim);
321  RCP<MV> coords;
322  for(unsigned d=0;d<dim;d++) {
323  const std::vector<double> & coord = callback->getCoordsVector(d);
324 
325  // no coords vector has been build yet, build one
326  if(coords==Teuchos::null) {
327  if(globalIndexer->getNumFields()==1) {
329  = rcp_dynamic_cast<const panzer::GlobalIndexer>(globalIndexer);
330  std::vector<panzer::GlobalOrdinal> ownedIndices;
331  ugi->getOwnedIndices(ownedIndices);
332  RCP<const Map> coords_map = rcp(new Map(Teuchos::OrdinalTraits<panzer::GlobalOrdinal>::invalid(),ownedIndices,0,mpi_comm));
333  coords = rcp(new MV(coords_map,dim));
334  }
335  else {
336  RCP<const Map> coords_map = rcp(new Map(Teuchos::OrdinalTraits<panzer::GlobalOrdinal>::invalid(),coord.size(),0,mpi_comm));
337  coords = rcp(new MV(coords_map,dim));
338  }
339  }
340 
341  // sanity check the size
342  TEUCHOS_ASSERT(coords->getLocalLength()==coord.size());
343 
344  // fill appropriate coords vector
345  Teuchos::ArrayRCP<double> dest = coords->getDataNonConst(d);
346  for(std::size_t i=0;i<coord.size();i++)
347  dest[i] = coord[i];
348  }
349 
350  // inject coordinates into parameter list
351  Teuchos::ParameterList & muelu_params = strat_params->sublist("Preconditioner Types").sublist("MueLu");
352  muelu_params.set<RCP<MV> >("Coordinates",coords);
353  }
354  #endif
355  }
356  // else write_out_the_mesg("Warning: No unique field determines the coordinates, coordinates unavailable!")
357 
358  Teko::addTekoToStratimikosBuilder(linearSolverBuilder,reqHandler_local);
359  }
360  else {
361  // try to set request handler from member variable
362  if(reqHandler_local==Teuchos::null)
363  reqHandler_local = rcp(new Teko::RequestHandler);
364 
365  std::string fieldName;
366  if(determineCoordinateField(*globalIndexer,fieldName)) {
368  rcp_dynamic_cast<const panzer::BlockedDOFManager>(globalIndexer);
370  rcp_dynamic_cast<const panzer::BlockedDOFManager>(auxGlobalIndexer);
372  rcp(new panzer_stk::ParameterListCallbackBlocked(stkConn_manager,blkDofs,auxBlkDofs));
373  reqHandler_local->addRequestCallback(callback);
374  }
375 
376  Teko::addTekoToStratimikosBuilder(linearSolverBuilder,reqHandler_local);
377 
378  if(writeCoordinates) {
380  rcp_dynamic_cast<const panzer::BlockedDOFManager>(globalIndexer);
381 
382  // loop over blocks
383  const std::vector<RCP<panzer::GlobalIndexer>> & dofVec
384  = blkDofs->getFieldDOFManagers();
385  for(std::size_t i=0;i<dofVec.size();i++) {
386 
387  // add in the coordinate parameter list callback handler
388  TEUCHOS_ASSERT(determineCoordinateField(*dofVec[i],fieldName));
389 
390  std::map<std::string,RCP<const panzer::Intrepid2FieldPattern> > fieldPatterns;
391  fillFieldPatternMap(*dofVec[i],fieldName,fieldPatterns);
392  panzer_stk::ParameterListCallback plCall(fieldName,fieldPatterns,stkConn_manager,dofVec[i]);
393  plCall.buildArrayToVector();
394  plCall.buildCoordinates();
395 
396  // extract coordinate vectors
397  const std::vector<double> & xcoords = plCall.getXCoordsVector();
398  const std::vector<double> & ycoords = plCall.getYCoordsVector();
399  const std::vector<double> & zcoords = plCall.getZCoordsVector();
400 
401  // use epetra to write coordinates to matrix market files
402  Epetra_MpiComm ep_comm(*mpi_comm->getRawMpiComm()); // this is OK access to RawMpiComm becase its declared on the stack?
403  // and all users of this object are on the stack (within scope of mpi_comm
404  Epetra_Map map(-1,xcoords.size(),0,ep_comm);
405 
406  RCP<Epetra_Vector> vec;
407  switch(spatialDim) {
408  case 3:
409  vec = rcp(new Epetra_Vector(Copy,map,const_cast<double *>(&zcoords[0])));
410  EpetraExt::VectorToMatrixMarketFile((fieldName+"_zcoords.mm").c_str(),*vec);
411  // Intentional fall-through.
412  case 2:
413  vec = rcp(new Epetra_Vector(Copy,map,const_cast<double *>(&ycoords[0])));
414  EpetraExt::VectorToMatrixMarketFile((fieldName+"_ycoords.mm").c_str(),*vec);
415  // Intentional fall-through.
416  case 1:
417  vec = rcp(new Epetra_Vector(Copy,map,const_cast<double *>(&xcoords[0])));
418  EpetraExt::VectorToMatrixMarketFile((fieldName+"_xcoords.mm").c_str(),*vec);
419  break;
420  default:
421  TEUCHOS_ASSERT(false);
422  }
423 
424  // TODO add MueLu code...
425  #ifdef PANZER_HAVE_MUELU
426  if(useCoordinates) {
427 
428  typedef Xpetra::Map<int,panzer::GlobalOrdinal> Map;
429  typedef Xpetra::MultiVector<double,int,panzer::GlobalOrdinal> MV;
430 
431  // TODO This is Epetra-specific
432  RCP<const Map> coords_map = Xpetra::MapFactory<int,panzer::GlobalOrdinal>::Build(Xpetra::UseEpetra,
434  //Teuchos::ArrayView<GO>(ownedIndices),
435  xcoords.size(),
436  0,
437  mpi_comm
438  );
439 
440  unsigned dim = Teuchos::as<unsigned>(spatialDim);
441 
442  RCP<MV> coords = Xpetra::MultiVectorFactory<double,int,panzer::GlobalOrdinal>::Build(coords_map,spatialDim);
443 
444  for(unsigned d=0;d<dim;d++) {
445  // sanity check the size
446  TEUCHOS_ASSERT(coords->getLocalLength()==xcoords.size());
447 
448  // fill appropriate coords vector
449  Teuchos::ArrayRCP<double> dest = coords->getDataNonConst(d);
450  for(std::size_t j=0;j<coords->getLocalLength();++j) {
451  if (d == 0) dest[j] = xcoords[j];
452  if (d == 1) dest[j] = ycoords[j];
453  if (d == 2) dest[j] = zcoords[j];
454  }
455  }
456 
457  // TODO This is Epetra-specific
458  // inject coordinates into parameter list
459  Teuchos::ParameterList & muelu_params = strat_params->sublist("Preconditioner Types").sublist("MueLu");
460  muelu_params.set<RCP<MV> >("Coordinates",coords);
461 
462  }
463  #endif
464 
465  } /* end loop over all physical fields */
466  }
467 
468  if(writeTopo) {
469  /*
470  RCP<const panzer::BlockedDOFManager> blkDofs =
471  rcp_dynamic_cast<const panzer::BlockedDOFManager>(globalIndexer);
472 
473  writeTopology(*blkDofs);
474  */
475  TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,
476  "Topology writing is no longer implemented. It needs to be reimplemented for the "
477  "default DOFManager (estimate 2 days with testing)");
478  }
479  }
480  #endif
481 
482  linearSolverBuilder.setParameterList(strat_params);
483  RCP<Thyra::LinearOpWithSolveFactoryBase<double> > lowsFactory = createLinearSolveStrategy(linearSolverBuilder);
484 
485  return lowsFactory;
486  }
487 
488 
490  buildLOWSFactory(bool blockedAssembly,
491  const Teuchos::RCP<const panzer::GlobalIndexer> & globalIndexer,
492  const Teuchos::RCP<panzer::ConnManager> & conn_manager,
493  int spatialDim,
494  const Teuchos::RCP<const Teuchos::MpiComm<int> > & mpi_comm,
495  const Teuchos::RCP<Teuchos::ParameterList> & strat_params,
496  #ifdef PANZER_HAVE_TEKO
497  const Teuchos::RCP<Teko::RequestHandler> & reqHandler,
498  #endif
499  bool writeCoordinates,
500  bool writeTopo,
501  const Teuchos::RCP<const panzer::GlobalIndexer> & auxGlobalIndexer,
502  bool useCoordinates
503  )
504  {
505  #ifdef PANZER_HAVE_TEKO
506  Teuchos::RCP<Teko::RequestHandler> reqHandler_local = reqHandler;
507  if(reqHandler_local==Teuchos::null)
508  reqHandler_local = Teuchos::rcp(new Teko::RequestHandler);
509  #endif
510 
511  auto stk_conn_manager = Teuchos::rcp_dynamic_cast<panzer_stk::STKConnManager>(conn_manager,true);
512 
513  return buildLOWSFactory(blockedAssembly,globalIndexer,stk_conn_manager,spatialDim,mpi_comm,strat_params,
514 #ifdef PANZER_HAVE_TEKO
515  reqHandler_local,
516 #endif
517  writeCoordinates,
518  writeTopo,
519  auxGlobalIndexer,
520  useCoordinates
521  );
522 
523  // should never reach this
524  TEUCHOS_ASSERT(false);
525  return Teuchos::null;
526  }
527 }
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.
void getElementBlockIds(std::vector< std::string > &elementBlockIds) const
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
const std::vector< Teuchos::RCP< GlobalIndexer > > & getFieldDOFManagers() const
virtual void getElementBlockIds(std::vector< std::string > &elementBlockIds) const =0
virtual void getOwnedIndices(std::vector< panzer::GlobalOrdinal > &indices) const =0
Get the set of indices owned by this processor.
virtual int getNumFields() const =0
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
bool isSublist(const std::string &name) const
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)
int VectorToMatrixMarketFile(const char *filename, const Epetra_Vector &A, const char *matrixName=0, const char *matrixDescription=0, bool writeHeader=true)
virtual const std::vector< int > & getBlockFieldNumbers(const std::string &blockId) const =0
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
#define TEUCHOS_ASSERT(assertion_test)
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...