Stratimikos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Thyra_AmesosLinearOpWithSolveFactory.cpp
1 // @HEADER
2 // *****************************************************************************
3 // Stratimikos: Thyra-based strategies for linear solvers
4 //
5 // Copyright 2006 NTESS and the Stratimikos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef SUN_CXX
11 
12 #include "Thyra_AmesosLinearOpWithSolveFactory.hpp"
13 
14 #include "Thyra_AmesosLinearOpWithSolve.hpp"
15 #include "Thyra_EpetraOperatorViewExtractorStd.hpp"
16 #include "Amesos.h"
17 #include "Teuchos_dyn_cast.hpp"
18 #include "Teuchos_TimeMonitor.hpp"
19 #include "Teuchos_TypeNameTraits.hpp"
20 #include "Teuchos_VerboseObjectParameterListHelpers.hpp"
21 
22 #ifdef HAVE_AMESOS_KLU
23 #include "Amesos_Klu.h"
24 #endif
25 #ifdef HAVE_AMESOS_PASTIX
26 #include "Amesos_Pastix.h"
27 #endif
28 #ifdef HAVE_AMESOS_LAPACK
29 #include "Amesos_Lapack.h"
30 #endif
31 #ifdef HAVE_AMESOS_MUMPS
32 #include "Amesos_Mumps.h"
33 #endif
34 #ifdef HAVE_AMESOS_SCALAPACK
35 #include "Amesos_Scalapack.h"
36 #endif
37 #ifdef HAVE_AMESOS_UMFPACK
38 #include "Amesos_Umfpack.h"
39 #endif
40 #ifdef HAVE_AMESOS_SUPERLUDIST
41 #include "Amesos_Superludist.h"
42 #endif
43 #ifdef HAVE_AMESOS_SUPERLU
44 #include "Amesos_Superlu.h"
45 #endif
46 #ifdef HAVE_AMESOS_DSCPACK
47 #include "Amesos_Dscpack.h"
48 #endif
49 #ifdef HAVE_AMESOS_PARDISO
50 #include "Amesos_Pardiso.h"
51 #endif
52 #ifdef HAVE_AMESOS_TAUCS
53 #include "Amesos_Taucs.h"
54 #endif
55 #ifdef HAVE_AMESOS_PARAKLETE
56 #include "Amesos_Paraklete.h"
57 #endif
58 
59 namespace {
60 
61 const std::string epetraFwdOp_str = "epetraFwdOp";
62 
63 } // namespace
64 
65 namespace Thyra {
66 
67 
68 // Parameter names for Paramter List
69 
70 const std::string AmesosLinearOpWithSolveFactory::SolverType_name = "Solver Type";
71 
72 const std::string AmesosLinearOpWithSolveFactory::RefactorizationPolicy_name = "Refactorization Policy";
73 
74 const std::string AmesosLinearOpWithSolveFactory::ThrowOnPreconditionerInput_name = "Throw on Preconditioner Input";
75 
76 const std::string AmesosLinearOpWithSolveFactory::Amesos_Settings_name = "Amesos Settings";
77 
78 // Constructors/initializers/accessors
79 
81 {
82 #ifdef TEUCHOS_DEBUG
83  if(paramList_.get())
84  paramList_->validateParameters(
85  *this->getValidParameters(),0 // Only validate this level for now!
86  );
87 #endif
88 }
89 
91  const Amesos::ESolverType solverType
92  ,const Amesos::ERefactorizationPolicy refactorizationPolicy
93  ,const bool throwOnPrecInput
94  )
95  :epetraFwdOpViewExtractor_(Teuchos::rcp(new EpetraOperatorViewExtractorStd()))
96  ,solverType_(solverType)
97  ,refactorizationPolicy_(refactorizationPolicy)
98  ,throwOnPrecInput_(throwOnPrecInput)
99 {}
100 
101 // Overridden from LinearOpWithSolveFactoryBase
102 
104  const LinearOpSourceBase<double> &fwdOpSrc
105  ) const
106 {
107  using Teuchos::outArg;
109  fwdOp = fwdOpSrc.getOp();
110  RCP<const Epetra_Operator> epetraFwdOp;
111  EOpTransp epetraFwdOpTransp;
112  EApplyEpetraOpAs epetraFwdOpApplyAs;
113  EAdjointEpetraOp epetraFwdOpAdjointSupport;
114  double epetraFwdOpScalar;
115  epetraFwdOpViewExtractor_->getEpetraOpView(
116  fwdOp,
117  outArg(epetraFwdOp), outArg(epetraFwdOpTransp),
118  outArg(epetraFwdOpApplyAs), outArg(epetraFwdOpAdjointSupport),
119  outArg(epetraFwdOpScalar)
120  );
121  if( !dynamic_cast<const Epetra_RowMatrix*>(&*epetraFwdOp) )
122  return false;
123  return true;
124 }
125 
128 {
130 }
131 
133  const RCP<const LinearOpSourceBase<double> > &fwdOpSrc
135  ,const ESupportSolveUse /* supportSolveUse */
136  ) const
137 {
138  using Teuchos::outArg;
139  THYRA_FUNC_TIME_MONITOR("Stratimikos: AmesosLOWSF");
140 #ifdef TEUCHOS_DEBUG
141  TEUCHOS_TEST_FOR_EXCEPT(Op==NULL);
142 #endif
144  fwdOp = fwdOpSrc->getOp();
145  //
146  // Unwrap and get the forward Epetra_Operator object
147  //
148  RCP<const Epetra_Operator> epetraFwdOp;
149  EOpTransp epetraFwdOpTransp;
150  EApplyEpetraOpAs epetraFwdOpApplyAs;
151  EAdjointEpetraOp epetraFwdOpAdjointSupport;
152  double epetraFwdOpScalar;
153  epetraFwdOpViewExtractor_->getEpetraOpView(
154  fwdOp,
155  outArg(epetraFwdOp), outArg(epetraFwdOpTransp),
156  outArg(epetraFwdOpApplyAs), outArg(epetraFwdOpAdjointSupport),
157  outArg(epetraFwdOpScalar)
158  );
159  // Get the AmesosLinearOpWithSolve object
161  *amesosOp = &Teuchos::dyn_cast<AmesosLinearOpWithSolve>(*Op);
162  //
163  // Determine if we must start over or not
164  //
165  bool startOver = ( amesosOp->get_amesosSolver()==Teuchos::null );
166  if(!startOver) {
167  startOver =
168  (
169  epetraFwdOpTransp != amesosOp->get_amesosSolverTransp() ||
170  epetraFwdOp.get() != amesosOp->get_epetraLP()->GetOperator()
171  // We must start over if the matrix object changes. This is a
172  // weakness of Amesos but there is nothing I can do about this right
173  // now!
174  );
175  }
176  //
177  // Update the amesos solver
178  //
179  if(startOver) {
180  //
181  // This LOWS object has not be initialized yet or is not compatible with the existing
182  //
183  // so this is where we setup everything from the ground up.
184  //
185  // Create the linear problem and set the operator with memory of RCP to Epetra_Operator view!
187  epetraLP = Teuchos::rcp(new Epetra_LinearProblem());
188  epetraLP->SetOperator(const_cast<Epetra_Operator*>(&*epetraFwdOp));
189  Teuchos::set_extra_data< RCP<const Epetra_Operator> >( epetraFwdOp, epetraFwdOp_str,
190  Teuchos::inOutArg(epetraLP) );
191  // Create the concrete solver
193  amesosSolver;
194  {
195  THYRA_FUNC_TIME_MONITOR_DIFF("Stratimikos: AmesosLOWSF:InitConstruct",
196  InitConstruct);
197  switch(solverType_) {
198  case Thyra::Amesos::LAPACK :
199  amesosSolver = Teuchos::rcp(new Amesos_Lapack(*epetraLP));
200  break;
201 #ifdef HAVE_AMESOS_KLU
202  case Thyra::Amesos::KLU :
203  amesosSolver = Teuchos::rcp(new Amesos_Klu(*epetraLP));
204  break;
205 #endif
206 #ifdef HAVE_AMESOS_PASTIX
207  case Thyra::Amesos::PASTIX :
208  amesosSolver = Teuchos::rcp(new Amesos_Pastix(*epetraLP));
209  break;
210 #endif
211 #ifdef HAVE_AMESOS_MUMPS
212  case Thyra::Amesos::MUMPS :
213  amesosSolver = Teuchos::rcp(new Amesos_Mumps(*epetraLP));
214  break;
215 #endif
216 #ifdef HAVE_AMESOS_SCALAPACK
217  case Thyra::Amesos::SCALAPACK :
218  amesosSolver = Teuchos::rcp(new Amesos_Scalapack(*epetraLP));
219  break;
220 #endif
221 #ifdef HAVE_AMESOS_UMFPACK
222  case Thyra::Amesos::UMFPACK :
223  amesosSolver = Teuchos::rcp(new Amesos_Umfpack(*epetraLP));
224  break;
225 #endif
226 #ifdef HAVE_AMESOS_SUPERLUDIST
227  case Thyra::Amesos::SUPERLUDIST :
228  amesosSolver = Teuchos::rcp(new Amesos_Superludist(*epetraLP));
229  break;
230 #endif
231 #ifdef HAVE_AMESOS_SUPERLU
232  case Thyra::Amesos::SUPERLU :
233  amesosSolver = Teuchos::rcp(new Amesos_Superlu(*epetraLP));
234  break;
235 #endif
236 #ifdef HAVE_AMESOS_DSCPACK
237  case Thyra::Amesos::DSCPACK :
238  amesosSolver = Teuchos::rcp(new Amesos_Dscpack(*epetraLP));
239  break;
240 #endif
241 #ifdef HAVE_AMESOS_PARDISO
242  case Thyra::Amesos::PARDISO :
243  amesosSolver = Teuchos::rcp(new Amesos_Pardiso(*epetraLP));
244  break;
245 #endif
246 #ifdef HAVE_AMESOS_TAUCS
247  case Thyra::Amesos::TAUCS :
248  amesosSolver = Teuchos::rcp(new Amesos_Taucs(*epetraLP));
249  break;
250 #endif
251 #ifdef HAVE_AMESOS_PARAKLETE
252  case Thyra::Amesos::PARAKLETE :
253  amesosSolver = Teuchos::rcp(new Amesos_Paraklete(*epetraLP));
254  break;
255 #endif
256  default:
258  true, std::logic_error
259  ,"Error, the solver type ID = " << solverType_ << " is invalid!"
260  );
261  }
262  }
263  // Set the parameters
264  if(paramList_.get()) amesosSolver->SetParameters(paramList_->sublist("Amesos Settings"));
265  // Do the initial factorization
266  {
267  THYRA_FUNC_TIME_MONITOR_DIFF("Stratimikos: AmesosLOWSF:Symbolic", Symbolic);
268  const int err = amesosSolver->SymbolicFactorization();
270  "Error, SymbolicFactorization() on amesos solver of type \'"<<Teuchos::typeName(*amesosSolver)<<"\'\n"
271  "returned error code "<<err<<"!" );
272  }
273  {
274  THYRA_FUNC_TIME_MONITOR_DIFF("Stratimikos: AmesosLOWSF:Factor", Factor);
275  const int err = amesosSolver->NumericFactorization();
277  "Error, NumericFactorization() on amesos solver of type \'"<<Teuchos::typeName(*amesosSolver)<<"\'\n"
278  "returned error code "<<err<<"!" );
279  }
280  // Initialize the LOWS object and we are done!
281  amesosOp->initialize(fwdOp,fwdOpSrc,epetraLP,amesosSolver,epetraFwdOpTransp,epetraFwdOpScalar);
282  }
283  else {
284  //
285  // This LOWS object has already be initialized once so we must just reset
286  // the matrix and refactor it.
287  //
288  // Get non-const pointers to the linear problem and the amesos solver.
289  // These const-casts are just fine since the amesosOp in non-const.
291  epetraLP = Teuchos::rcp_const_cast<Epetra_LinearProblem>(amesosOp->get_epetraLP());
293  amesosSolver = amesosOp->get_amesosSolver();
294  // Reset the forward operator with memory of RCP to Epetra_Operator view!
295  epetraLP->SetOperator(const_cast<Epetra_Operator*>(&*epetraFwdOp));
296  Teuchos::get_nonconst_extra_data<RCP<const Epetra_Operator> >(epetraLP,epetraFwdOp_str) = epetraFwdOp;
297  // Reset the parameters
298  if(paramList_.get()) amesosSolver->SetParameters(paramList_->sublist(Amesos_Settings_name));
299  // Repivot if asked
300  if(refactorizationPolicy_==Amesos::REPIVOT_ON_REFACTORIZATION) {
301  THYRA_FUNC_TIME_MONITOR_DIFF("Stratimikos: AmesosLOWSF:Symbolic", Symbolic);
302  const int err = amesosSolver->SymbolicFactorization();
304  "Error, SymbolicFactorization() on amesos solver of type \'"<<Teuchos::typeName(*amesosSolver)<<"\'\n"
305  "returned error code "<<err<<"!" );
306  }
307  {
308  THYRA_FUNC_TIME_MONITOR_DIFF("Stratimikos: AmesosLOWSF::Factor", Factor);
309  const int err = amesosSolver->NumericFactorization();
311  "Error, NumericFactorization() on amesos solver of type \'"<<Teuchos::typeName(*amesosSolver)<<"\'\n"
312  "returned error code "<<err<<"!" );
313  }
314  /* ToDo: Put this back in once PrintStatus accepts an std::ostream!
315  OsTab tab(out);
316  amesosSolver->PrintStatus()
317  */
318  // Reinitialize the LOWS object and we are done! (we must do this to get the
319  // possibly new transpose and scaling factors back in)
320  amesosOp->initialize(fwdOp,fwdOpSrc,epetraLP,amesosSolver,epetraFwdOpTransp,epetraFwdOpScalar);
321  }
322  amesosOp->setOStream(this->getOStream());
323  amesosOp->setVerbLevel(this->getVerbLevel());
324 }
325 
327 {
328  return false;
329 }
330 
332  const RCP<const LinearOpSourceBase<double> > &fwdOpSrc
333  ,const RCP<const PreconditionerBase<double> > &/* prec */
335  ,const ESupportSolveUse supportSolveUse
336  ) const
337 {
339  this->throwOnPrecInput_, std::logic_error
340  ,"Error, the concrete implementation described as \'"<<this->description()<<"\' does not support preconditioners "
341  "and has been configured to throw this exception when the initializePreconditionedOp(...) function is called!"
342  );
343  this->initializeOp(fwdOpSrc,Op,supportSolveUse); // Ignore the preconditioner!
344 }
345 
347  const RCP<const LinearOpSourceBase<double> > &fwdOpSrc
348  ,const RCP<const LinearOpSourceBase<double> > &/* approxFwdOpSrc */
350  ,const ESupportSolveUse supportSolveUse
351  ) const
352 {
354  this->throwOnPrecInput_, std::logic_error
355  ,"Error, the concrete implementation described as \'"<<this->description()<<"\' does not support preconditioners "
356  "and has been configured to throw this exception when the initializePreconditionedOp(...) function is called!"
357  );
358  this->initializeOp(fwdOpSrc,Op,supportSolveUse); // Ignore the preconditioner!
359 }
360 
363  ,RCP<const LinearOpSourceBase<double> > *fwdOpSrc
364  ,RCP<const PreconditionerBase<double> > *prec
365  ,RCP<const LinearOpSourceBase<double> > *approxFwdOpSrc
366  ,ESupportSolveUse * /* supportSolveUse */
367  ) const
368 {
369 #ifdef TEUCHOS_DEBUG
370  TEUCHOS_TEST_FOR_EXCEPT(Op==NULL);
371 #endif
373  *amesosOp = &Teuchos::dyn_cast<AmesosLinearOpWithSolve>(*Op);
375  _fwdOpSrc = amesosOp->extract_fwdOpSrc(); // Will be null if uninitialized!
376  if(_fwdOpSrc.get()) {
377  // Erase the Epetra_Operator view of the forward operator!
378  RCP<Epetra_LinearProblem> epetraLP = amesosOp->get_epetraLP();
379  Teuchos::get_nonconst_extra_data< RCP<const Epetra_Operator> >(
380  epetraLP,epetraFwdOp_str
381  )
382  = Teuchos::null;
383  // Note, we did not erase the address of the operator in
384  // epetraLP->GetOperator() since it seems that the amesos solvers do not
385  // recheck the value of GetProblem()->GetOperator() so you had better not
386  // rest this!
387  }
388  if(fwdOpSrc) *fwdOpSrc = _fwdOpSrc; // It is fine if the client does not want this object back!
389  if(prec) *prec = Teuchos::null; // We never keep a preconditioner!
390  if(approxFwdOpSrc) *approxFwdOpSrc = Teuchos::null; // We never keep an approximate fwd operator!
391 }
392 
393 // Overridden from ParameterListAcceptor
394 
396  RCP<Teuchos::ParameterList> const& paramList
397  )
398 {
399  TEUCHOS_TEST_FOR_EXCEPT(paramList.get()==NULL);
400  paramList->validateParameters(*this->getValidParameters(),0); // Only validate this level for now!
401  paramList_ = paramList;
402  solverType_ =
403  Amesos::solverTypeNameToEnumMap.get<Amesos::ESolverType>(
404  paramList_->get(
406  ,Amesos::toString(solverType_)
407  )
408  ,paramList_->name()+"->"+SolverType_name
409  );
410  refactorizationPolicy_ =
411  Amesos::refactorizationPolicyNameToEnumMap.get<Amesos::ERefactorizationPolicy>(
412  paramList_->get(
414  ,Amesos::toString(refactorizationPolicy_)
415  )
416  ,paramList_->name()+"->"+RefactorizationPolicy_name
417  );
418  throwOnPrecInput_ = paramList_->get(ThrowOnPreconditionerInput_name,throwOnPrecInput_);
419  Teuchos::readVerboseObjectSublist(&*paramList_,this);
420 }
421 
424 {
425  return paramList_;
426 }
427 
430 {
431  RCP<Teuchos::ParameterList> _paramList = paramList_;
432  paramList_ = Teuchos::null;
433  return _paramList;
434 }
435 
438 {
439  return paramList_;
440 }
441 
444 {
445  return generateAndGetValidParameters();
446 }
447 
448 // Public functions overridden from Teuchos::Describable
449 
451 {
452  std::ostringstream oss;
453  oss << "Thyra::AmesosLinearOpWithSolveFactory{";
454  oss << "solverType=" << toString(solverType_);
455  oss << "}";
456  return oss.str();
457 }
458 
459 // private
460 
462 AmesosLinearOpWithSolveFactory::generateAndGetValidParameters()
463 {
464  static RCP<Teuchos::ParameterList> validParamList;
465  if(validParamList.get()==NULL) {
466  validParamList = Teuchos::rcp(new Teuchos::ParameterList("Amesos"));
467  validParamList->set(
469 #ifdef HAVE_AMESOS_KLU
470  ,Amesos::toString(Amesos::KLU)
471 #else
472  ,Amesos::toString(Amesos::LAPACK)
473 #endif
474  );
475  validParamList->set(RefactorizationPolicy_name,Amesos::toString(Amesos::REPIVOT_ON_REFACTORIZATION));
476  validParamList->set(ThrowOnPreconditionerInput_name,bool(true));
478  Teuchos::setupVerboseObjectSublist(&*validParamList);
479  }
480  return validParamList;
481 }
482 
483 } // namespace Thyra
484 
485 #endif // SUN_CXX
const std::string & name() const
Teuchos::RCP< const LinearOpSourceBase< double > > extract_fwdOpSrc()
Extract the LinearOpSourceBase&lt;double&gt; object so that it can be modified.
void SetOperator(Epetra_RowMatrix *A)
EOpTransp
bool supportsPreconditionerInputType(const EPreconditionerInputType precOpType) const
Returns false .
T & get(const std::string &name, T def_value)
bool isCompatible(const LinearOpSourceBase< double > &fwdOpSrc) const
Returns true if dynamic_cast&lt;const EpetraLinearOpBase*&gt;(fwdOpSrc)!=NULL .
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Teuchos::RCP< Amesos_BaseSolver > get_amesosSolver() const
T_To & dyn_cast(T_From &from)
static Teuchos::ParameterList GetValidParameters()
T * get() const
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
Teuchos::RCP< const Teuchos::ParameterList > getParameterList() const
void setParameterList(Teuchos::RCP< Teuchos::ParameterList > const &paramList)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< Teuchos::ParameterList > getNonconstParameterList()
EAdjointEpetraOp
void uninitializeOp(LinearOpWithSolveBase< double > *Op, Teuchos::RCP< const LinearOpSourceBase< double > > *fwdOpSrc, Teuchos::RCP< const PreconditionerBase< double > > *prec, Teuchos::RCP< const LinearOpSourceBase< double > > *approxFwdOpSrc, ESupportSolveUse *supportSolveUse) const
TEUCHOSCORE_LIB_DLL_EXPORT std::string toString(const EVerbosityLevel verbLevel)
void initialize(const Teuchos::RCP< const LinearOpBase< double > > &fwdOp, const Teuchos::RCP< const LinearOpSourceBase< double > > &fwdOpSrc, const Teuchos::RCP< Epetra_LinearProblem > &epetraLP, const Teuchos::RCP< Amesos_BaseSolver > &amesosSolver, const EOpTransp amesosSolverTransp, const double amesosSolverScalar)
First initialization.
Concrete LinearOpWithSolveBase subclass that adapts any Amesos_BaseSolver object. ...
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
ParameterList & setParameters(const ParameterList &source)
ESupportSolveUse
Teuchos::RCP< LinearOpWithSolveBase< double > > createOp() const
void initializeOp(const Teuchos::RCP< const LinearOpSourceBase< double > > &fwdOpSrc, LinearOpWithSolveBase< double > *Op, const ESupportSolveUse supportSolveUse) const
void initializePreconditionedOp(const Teuchos::RCP< const LinearOpSourceBase< double > > &fwdOpSrc, const Teuchos::RCP< const PreconditionerBase< double > > &prec, LinearOpWithSolveBase< double > *Op, const ESupportSolveUse supportSolveUse) const
Throws exception if this-&gt;throwOnPrecInput()==true and calls this-&gt;initializeOp(fwdOpSrc,Op) otherwise.
AmesosLinearOpWithSolveFactory(const Amesos::ESolverType solverType=Amesos::LAPACK, const Amesos::ERefactorizationPolicy refactorizationPolicy=Amesos::REPIVOT_ON_REFACTORIZATION, const bool throwOnPrecInput=true)
Constructor which sets the defaults.
EApplyEpetraOpAs
EPreconditionerInputType
Teuchos::RCP< Epetra_LinearProblem > get_epetraLP() const
Epetra_Operator * GetOperator() const
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
virtual Teuchos::RCP< const LinearOpBase< Scalar > > getOp() const =0
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)
std::string typeName(const T &t)
Teuchos::RCP< Teuchos::ParameterList > unsetParameterList()

Generated on Thu Nov 21 2024 09:22:16 for Stratimikos by doxygen 1.8.5