Stratimikos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Thyra_AztecOOLinearOpWithSolveFactory.cpp
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Stratimikos: Thyra-based strategies for linear solvers
5 // Copyright (2006) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Roscoe A. Bartlett (rabartl@sandia.gov)
38 //
39 // ***********************************************************************
40 // @HEADER
41 
42 
43 #ifndef SUN_CXX
44 
45 
46 #include "Thyra_AztecOOLinearOpWithSolveFactory.hpp"
47 #include "Thyra_AztecOOLinearOpWithSolve.hpp"
48 #include "Thyra_PreconditionerFactoryHelpers.hpp"
49 #include "Thyra_EpetraOperatorViewExtractorStd.hpp"
50 #include "Thyra_ScaledAdjointLinearOpBase.hpp"
51 #include "Thyra_EpetraLinearOpBase.hpp"
52 #include "Thyra_EpetraOperatorWrapper.hpp"
54 #include "Teuchos_VerboseObjectParameterListHelpers.hpp"
55 #include "Teuchos_ParameterList.hpp"
56 #include "Teuchos_dyn_cast.hpp"
57 #include "AztecOOParameterList.hpp"
58 
59 
60 namespace {
61 
62 
63 const std::string AOOLOWSF_epetraPrecOp_str
64 = "AOOLOWSF::epetraPrecOp";
65 const std::string AOOLOWSF_aztec_epetra_epetraFwdOp_str
66 = "AOOLOWSF::aztec_epetra_epetraFwdOp";
67 const std::string AOOLOWSF_aztec_epetra_epetraAdjOp_str
68 = "AOOLOWSF::aztec_epetra_epetraAdjOp";
69 const std::string AOOLOWSF_rowmatrix_epetraFwdOp_str
70 = "AOOLOWSF::rowmatrix_epetraFwdOp";
71 const std::string AOOLOWSF_rowmatrix_epetraPrecOp_str
72 = "AOOLOWSF::rowmatrix_epetraPrecOp";
73 const std::string AOOLOWSF_aztec_fwd_epetra_epetraPrecOp_str
74 = "AOOLOWSF::aztec_fwd_epetra_epetraPrecOp";
75 const std::string AOOLOWSF_aztec_adj_epetra_epetraPrecOp_str
76 = "AOOLOWSF::aztec_adj_epetra_epetraPrecOp";
77 const std::string AOOLOWSF_setPrecondtionerOperator_str
78 = "AOOLOWSF::setPrecondtionerOperator";
79 const std::string AOOLOWSF_constructedAztecPreconditoner_str
80 = "AOOLOWSF::constructedAztecPreconditoner";
81 
82 
83 const std::string ForwardSolve_name = "Forward Solve";
84 const std::string AdjointSolve_name = "Adjoint Solve";
85 const std::string MaxIterations_name = "Max Iterations";
86 const int MaxIterations_default = 400;
87 const std::string Tolerance_name = "Tolerance";
88 const double Tolerance_default = 1e-6;
89 const std::string OutputEveryRhs_name = "Output Every RHS";
90 const bool OutputEveryRhs_default = false;
91 const std::string AztecOO_Settings_name = "AztecOO Settings";
92 
93 
94 } // namespace
95 
96 
97 namespace Thyra {
98 
99 
100 // Constructors/initializers/accessors
101 
102 
104  Teuchos::RCP<Teuchos::ParameterList> const& paramList
105  )
106  :epetraFwdOpViewExtractor_(Teuchos::rcp(new EpetraOperatorViewExtractorStd()))
107  ,defaultFwdMaxIterations_(MaxIterations_default)
108  ,defaultFwdTolerance_(Tolerance_default)
109  ,defaultAdjMaxIterations_(MaxIterations_default)
110  ,defaultAdjTolerance_(Tolerance_default)
111  ,outputEveryRhs_(OutputEveryRhs_default)
112  ,useAztecPrec_(false)
113 {
114  updateThisValidParamList();
115  if(paramList.get())
116  setParameterList(paramList);
117 }
118 
119 
120 // Overridden from LinearOpWithSolveFactoryBase
121 
122 
124 {
125  return true;
126 }
127 
128 
130  const Teuchos::RCP<PreconditionerFactoryBase<double> > &precFactory,
131  const std::string &precFactoryName
132  )
133 {
134  TEUCHOS_TEST_FOR_EXCEPT(!precFactory.get());
136  precFactoryValidPL = precFactory->getValidParameters();
137  const std::string _precFactoryName =
138  ( precFactoryName != ""
139  ? precFactoryName
140  : ( precFactoryValidPL.get()
141  ? precFactoryValidPL->name()
142  : "GENERIC PRECONDITIONER FACTORY"
143  )
144  );
145  precFactory_ = precFactory;
146  precFactoryName_ = _precFactoryName;
147  updateThisValidParamList();
148 }
149 
150 
153 {
154  return precFactory_;
155 }
156 
157 
160  std::string *precFactoryName
161  )
162 {
163  if(precFactory) *precFactory = precFactory_;
164  if(precFactoryName) *precFactoryName = precFactoryName_;
165  precFactory_ = Teuchos::null;
166  precFactoryName_ = "";
167  updateThisValidParamList();
168 }
169 
170 
172  const LinearOpSourceBase<double> &fwdOpSrc
173  ) const
174 {
175  return epetraFwdOpViewExtractor_->isCompatible(*fwdOpSrc.getOp());
176 }
177 
178 
181 {
183 }
184 
185 
187  const Teuchos::RCP<const LinearOpSourceBase<double> > &fwdOpSrc,
189  const ESupportSolveUse /* supportSolveUse */
190  ) const
191 {
192  this->initializeOp_impl(fwdOpSrc,Teuchos::null,Teuchos::null,false,Op);
193 }
194 
195 
197  const Teuchos::RCP<const LinearOpSourceBase<double> > &fwdOpSrc,
199  ) const
200 {
201  this->initializeOp_impl(fwdOpSrc,Teuchos::null,Teuchos::null,true,Op);
202 }
203 
204 
206  const EPreconditionerInputType precOpType
207  ) const
208 {
209  const_cast<bool&>(useAztecPrec_) = (
210  paramList_.get()
211  &&
212  paramList_->sublist(ForwardSolve_name).sublist(AztecOO_Settings_name).get(
213  "Aztec Preconditioner","none"
214  )!="none"
215  );
216  switch(precOpType) {
218  return true;
220  return useAztecPrec_;
221  default:
223  }
225 }
226 
227 
229  const Teuchos::RCP<const LinearOpSourceBase<double> > &fwdOpSrc,
230  const Teuchos::RCP<const PreconditionerBase<double> > &prec,
232  const ESupportSolveUse /* supportSolveUse */
233  ) const
234 {
235  TEUCHOS_TEST_FOR_EXCEPT(prec.get()==NULL);
236  this->initializeOp_impl(fwdOpSrc,prec,Teuchos::null,false,Op);
237 }
238 
239 
241  const Teuchos::RCP<const LinearOpSourceBase<double> > &fwdOpSrc,
242  const Teuchos::RCP<const LinearOpSourceBase<double> > &approxFwdOpSrc,
244  const ESupportSolveUse /* supportSolveUse */
245  ) const
246 {
247  TEUCHOS_TEST_FOR_EXCEPT(approxFwdOpSrc.get()==NULL);
248  TEUCHOS_TEST_FOR_EXCEPT(approxFwdOpSrc->getOp().get()==NULL);
249  this->initializeOp_impl(fwdOpSrc,Teuchos::null,approxFwdOpSrc,false,Op);
250 }
251 
252 
255  Teuchos::RCP<const LinearOpSourceBase<double> > *fwdOpSrc,
257  Teuchos::RCP<const LinearOpSourceBase<double> > *approxFwdOpSrc,
258  ESupportSolveUse * /* supportSolveUse */
259  ) const
260 {
261 #ifdef TEUCHOS_DEBUG
262  TEUCHOS_TEST_FOR_EXCEPT(Op==NULL);
263 #endif
266  // Extract and unset the fwdOP and approxFwdOp objects
268  _fwdOpSrc = aztecOp->extract_fwdOpSrc(), // Will be null if not initialized!
269  _approxFwdOpSrc = aztecOp->extract_approxFwdOpSrc(); // Will be null if not set
270  if(fwdOpSrc) *fwdOpSrc = _fwdOpSrc;
271  if(approxFwdOpSrc) *approxFwdOpSrc = _approxFwdOpSrc;
272  // Only extract and uset the prec object if it is external. If it is
273  // internal, then we need to hold on to this so that we can reinitialize it
274  // later.
275  if(aztecOp->isExternalPrec()) {
277  _prec = aztecOp->extract_prec(); // Will be null if not external preconditioner was set
278  if(prec) *prec = _prec;
279  }
280  // ToDo: Extract the Epetra_Operator views what where used to initialize the
281  // forward and adjoint solvers! This is needed to make this totally
282  // stateless.
283 }
284 
285 
286 // Overridden from ParameterListAcceptor
287 
288 
290  Teuchos::RCP<Teuchos::ParameterList> const& paramList
291  )
292 {
293  TEUCHOS_TEST_FOR_EXCEPT(paramList.get()==NULL);
294  paramList->validateParameters(*this->getValidParameters());
295  paramList_ = paramList;
296  //
297  outputEveryRhs_ = paramList_->get(OutputEveryRhs_name,OutputEveryRhs_default);
298  // Foward Solve parameters
300  &fwdSolvePL = paramList_->sublist(ForwardSolve_name);
301  defaultFwdMaxIterations_ = fwdSolvePL.get(MaxIterations_name,defaultFwdMaxIterations_);
302  defaultFwdTolerance_ = fwdSolvePL.get(Tolerance_name,defaultFwdTolerance_);
303  // Adjoint Solve parameters
304  if( !paramList_->getPtr<Teuchos::ParameterList>(AdjointSolve_name) ) {
305  // If adjoint solve sublist is not set, then use the forward solve parameters
306  paramList_->sublist(AdjointSolve_name).setParameters(fwdSolvePL);
307  }
309  &adjSolvePL = paramList_->sublist(AdjointSolve_name);
310  defaultAdjMaxIterations_ = adjSolvePL.get(MaxIterations_name,defaultAdjMaxIterations_);
311  defaultAdjTolerance_ = adjSolvePL.get(Tolerance_name,defaultAdjTolerance_);
312  //
313  if(precFactory_.get()) {
314  // Only reset the PF's PL if the sublist exists or the PF does not already
315  // have a PL. We don't want to overwrite an externally set PL for the PF
316  // if we don't have a nested sublist defined here!
317  const bool nestedPFSublistExists = paramList_->isSublist(precFactoryName_);
318  const bool alreadyHasSublist = !is_null(precFactory_->getParameterList());
319  if( nestedPFSublistExists || !alreadyHasSublist ) {
320  precFactory_->setParameterList(Teuchos::sublist(paramList_,precFactoryName_));
321  }
322  }
323  Teuchos::readVerboseObjectSublist(&*paramList_,this);
324 }
325 
326 
329 {
330  return paramList_;
331 }
332 
333 
336 {
337  Teuchos::RCP<Teuchos::ParameterList> _paramList = paramList_;
338  paramList_ = Teuchos::null;
339  return _paramList;
340 }
341 
342 
345 {
346  return paramList_;
347 }
348 
349 
352 {
353  return thisValidParamList_;
354 }
355 
356 
357 // Public functions overridden from Teuchos::Describable
358 
359 
361 {
362  std::ostringstream oss;
363  oss << "Thyra::AztecOOLinearOpWithSolveFactory{";
364  oss << "precFactory=";
365  if(!is_null(precFactory_))
366  oss << precFactory_->description();
367  else
368  oss << "NULL";
369  oss << "}";
370  return oss.str();
371 }
372 
373 
374 // private
375 
376 
378 AztecOOLinearOpWithSolveFactory::generateAndGetValidParameters()
379 {
380  static Teuchos::RCP<Teuchos::ParameterList> validParamList;
381  if(validParamList.get()==NULL) {
382  validParamList = Teuchos::rcp(
383  new Teuchos::ParameterList("AztecOOLinearOpWithSolveFactory"));
384  validParamList->set(
385  OutputEveryRhs_name,OutputEveryRhs_default
386  ,"Determines if output is created for each individual RHS (true or 1) or if output\n"
387  "is just created for an entire set of RHSs (false or 0)."
388  );
390  aztecParamList = getValidAztecOOParameters();
392  &fwdSolvePL = validParamList->sublist(
393  ForwardSolve_name, false
394  ,"Gives the options for the forward solve."
395  );
396  fwdSolvePL.set(
397  Tolerance_name,Tolerance_default
398  ,"The tolerence used in the convergence check (see the convergence test\n"
399  "in the sublist \"" + AztecOO_Settings_name + "\")"
400  );
401  fwdSolvePL.set(
402  MaxIterations_name,MaxIterations_default
403  ,"The maximum number of iterations the AztecOO solver is allowed to perform."
404  );
405  fwdSolvePL.sublist(
406  AztecOO_Settings_name,false
407  ,"Sets the parameters on the AztecOO object itself."
408  ).setParameters(*aztecParamList);
410  &adjSolvePL = validParamList->sublist(
411  AdjointSolve_name, false
412  ,"The options for the adjoint solve.\n"
413  "If this sublist is missing then the parameters from the\n"
414  "\""+ForwardSolve_name+"\" sublist are used instead."
415  );
416  // Make the adjoint solve have same defaults as forward solve
417  adjSolvePL.setParameters(fwdSolvePL);
418  }
419  return validParamList;
420 }
421 
422 
423 void AztecOOLinearOpWithSolveFactory::updateThisValidParamList()
424 {
425  thisValidParamList_ = Teuchos::rcp(
426  new Teuchos::ParameterList(*generateAndGetValidParameters())
427  );
428  if(precFactory_.get()) {
430  precFactoryValidParamList = precFactory_->getValidParameters();
431  if(precFactoryValidParamList.get()) {
432  thisValidParamList_->sublist(precFactoryName_).setParameters(
433  *precFactoryValidParamList);
434  }
435  }
436  Teuchos::setupVerboseObjectSublist(&*thisValidParamList_);
437 }
438 
439 
440 void AztecOOLinearOpWithSolveFactory::initializeOp_impl(
441  const Teuchos::RCP<const LinearOpSourceBase<double> > &fwdOpSrc,
442  const Teuchos::RCP<const PreconditionerBase<double> > &prec,
443  const Teuchos::RCP<const LinearOpSourceBase<double> > &approxFwdOpSrc,
444  const bool reusePrec,
445  LinearOpWithSolveBase<double> *Op
446  ) const
447 {
448  using Teuchos::RCP;
449  using Teuchos::null;
450  using Teuchos::rcp;
451  using Teuchos::rcp_dynamic_cast;
452  using Teuchos::rcp_const_cast;
453  using Teuchos::set_extra_data;
454  using Teuchos::get_optional_extra_data;
455  using Teuchos::get_optional_nonconst_extra_data;
456  using Teuchos::outArg;
457  typedef EpetraExt::ProductOperator PO;
458 
459  const Teuchos::RCP<Teuchos::FancyOStream> out = this->getOStream();
460  const Teuchos::EVerbosityLevel verbLevel = this->getVerbLevel();
461  Teuchos::OSTab tab(out);
462  if(out.get() && static_cast<int>(verbLevel) > static_cast<int>(Teuchos::VERB_LOW))
463  *out << "\nEntering Thyra::AztecOOLinearOpWithSolveFactory::initializeOp_impl(...) ...\n";
464 
466  VOTSPF precFactoryOutputTempState(precFactory_,out,verbLevel);
467 
468 #ifdef TEUCHOS_DEBUG
469  TEUCHOS_TEST_FOR_EXCEPT(Op==NULL);
470  TEUCHOS_TEST_FOR_EXCEPT(fwdOpSrc.get()==NULL);
471  TEUCHOS_TEST_FOR_EXCEPT(fwdOpSrc->getOp().get()==NULL);
472 #endif
473 
474  //
475  // Determine whether the operators are EpetraLinearOp objects. If so, we're
476  // good to go. If not, we need to wrap it as an Epetra_Operator with some
477  // invasive code.
478  //
480  tmpFwdOp = fwdOpSrc->getOp(),
481  tmpApproxFwdOp = ( approxFwdOpSrc.get() ? approxFwdOpSrc->getOp() : Teuchos::null );
484  if ( dynamic_cast<const EpetraLinearOpBase*>(tmpFwdOp.get())!=0 )
485  {
486  fwdOp = tmpFwdOp;
487  approxFwdOp = tmpApproxFwdOp;
488  }
489  else
490  {
491  fwdOp = makeEpetraWrapper(tmpFwdOp);
492  if (
493  tmpApproxFwdOp.get()
494  &&
495  dynamic_cast<const EpetraLinearOpBase*>(&*tmpApproxFwdOp.get())
496  )
497  {
498  approxFwdOp = makeEpetraWrapper(tmpApproxFwdOp);
499  }
500  }
501 
502  //
503  // Get the AztecOOLinearOpWithSolve object
504  //
505  AztecOOLinearOpWithSolve
506  *aztecOp = &Teuchos::dyn_cast<AztecOOLinearOpWithSolve>(*Op);
507 
508  //
509  // Unwrap and get the forward operator or matrix
510  //
511  Teuchos::RCP<const Epetra_Operator> epetra_epetraFwdOp;
512  EOpTransp epetra_epetraFwdOpTransp;
513  EApplyEpetraOpAs epetra_epetraFwdOpApplyAs;
514  EAdjointEpetraOp epetra_epetraFwdOpAdjointSupport;
515  double epetra_epetraFwdOpScalar;
516  epetraFwdOpViewExtractor_->getEpetraOpView(
517  fwdOp,
518  outArg(epetra_epetraFwdOp), outArg(epetra_epetraFwdOpTransp),
519  outArg(epetra_epetraFwdOpApplyAs), outArg(epetra_epetraFwdOpAdjointSupport),
520  outArg(epetra_epetraFwdOpScalar)
521  );
523  epetra_epetraFwdOp.get()==NULL, std::logic_error
524  ,"Error, The input fwdOp object must be fully initialized "
525  "before calling this function!"
526  );
527 
528  //
529  // Get the preconditioner object to use
530  //
533  if (prec.get()) {
534  // We will be used the passed in external preconditioner
535  precUsed = prec;
536  }
537  else if (precFactory_.get() ) {
538  // We will be creating our own preconditioner using an externally set
539  // preconditioner factory
540  myPrec =
541  ( !aztecOp->isExternalPrec()
542  ? Teuchos::rcp_const_cast<PreconditionerBase<double> >(
543  aztecOp->extract_prec())
544  : Teuchos::null
545  );
546  if(myPrec.get()) {
547  // ToDo: Get the forward operator and validate that it is the same
548  // operator that is used here!
549  }
550  else {
551  myPrec = precFactory_->createPrec();
552  }
553  precFactory_->initializePrec(fwdOpSrc,&*myPrec);
554  precUsed = myPrec;
555  }
556 
557  //
558  // Unwrap and get the preconditioner operator
559  //
560  RCP<const LinearOpBase<double> > rightPrecOp;
561  if (precUsed.get()) {
562  RCP<const LinearOpBase<double> > unspecified = precUsed->getUnspecifiedPrecOp();
563  RCP<const LinearOpBase<double> > left = precUsed->getLeftPrecOp();
564  RCP<const LinearOpBase<double> > right = precUsed->getRightPrecOp();
566  !( left.get() || right.get() || unspecified.get() ), std::logic_error
567  ,"Error, at least one preconditoner linear operator objects must be set!"
568  );
569  if(unspecified.get()) {
570  rightPrecOp = unspecified;
571  }
572  else {
573  // Set a left, right or split preconditioner
575  left.get(),std::logic_error
576  ,"Error, we can not currently handle a left"
577  " preconditioner with the AztecOO/Thyra adapters!"
578  );
579  rightPrecOp = right;
580  }
581  }
582  double wrappedPrecOpScalar = 0.0;
583  EOpTransp wrappedPrecOpTransp = NOTRANS;
584  RCP<const LinearOpBase<double> > wrappedPrecOp = null;
585  RCP<const EpetraLinearOpBase> epetraPrecOp;
586  Teuchos::RCP<const Epetra_Operator> epetra_epetraPrecOp;
587  EOpTransp epetra_epetraPrecOpTransp;
588  EApplyEpetraOpAs epetra_epetraPrecOpApplyAs;
589  EAdjointEpetraOp epetra_epetraPrecOpAdjointSupport;
590  EOpTransp overall_epetra_epetraPrecOpTransp=NOTRANS;
591  if(rightPrecOp.get()) {
592  RCP<const LinearOpBase<double> > tmpWrappedPrecOp;
593  unwrap(
594  rightPrecOp,&wrappedPrecOpScalar,&wrappedPrecOpTransp,&tmpWrappedPrecOp);
595  if( dynamic_cast<const EpetraLinearOpBase*>(&*tmpWrappedPrecOp) ) {
596  wrappedPrecOp = tmpWrappedPrecOp;
597  }
598  else {
599  wrappedPrecOp = makeEpetraWrapper(tmpWrappedPrecOp);
600  }
601  epetraPrecOp = rcp_dynamic_cast<const EpetraLinearOpBase>(
602  wrappedPrecOp,true);
603  epetraPrecOp->getEpetraOpView(
604  outArg(epetra_epetraPrecOp), outArg(epetra_epetraPrecOpTransp),
605  outArg(epetra_epetraPrecOpApplyAs), outArg(epetra_epetraPrecOpAdjointSupport));
607  epetra_epetraPrecOp.get()==NULL,std::logic_error
608  ,"Error, The input prec object and its embedded preconditioner"
609  " operator must be fully initialized before calling this function!"
610  );
611  // 2007/08/10: rabartl: This next set_extra_data(...) call is likely to be
612  // setting up a circular reference! Since epetra_epetraPrecOp was
613  // gotten from epetraPrecOp, if you set epetraPrecOp as extra data
614  // on the RCP epetra_epetraPrecOp then you have a circular reference!
615  //set_extra_data(
616  // epetraPrecOp, AOOLOWSF_epetraPrecOp_str, &epetra_epetraPrecOp,
617  // Teuchos::POST_DESTROY, false );
618  overall_epetra_epetraPrecOpTransp
619  = trans_trans(
620  real_trans(wrappedPrecOpTransp),
621  real_trans(epetra_epetraPrecOpTransp)
622  );
623  }
624 
625  //
626  // Unwrap and get the approximate forward operator to be used to generate a
627  // preconditioner
628  //
629  if(approxFwdOp.get()) {
630  // Note, here we just use the same members data that would be set for an
631  // extenral preconditioner operator since it is not getting used.
632  unwrap(approxFwdOp,&wrappedPrecOpScalar,&wrappedPrecOpTransp,&wrappedPrecOp);
633  epetraPrecOp = rcp_dynamic_cast<const EpetraLinearOpBase>(
634  wrappedPrecOp,true);
635  epetraPrecOp->getEpetraOpView(
636  outArg(epetra_epetraPrecOp), outArg(epetra_epetraPrecOpTransp),
637  outArg(epetra_epetraPrecOpApplyAs), outArg(epetra_epetraPrecOpAdjointSupport)
638  );
640  epetra_epetraPrecOp.get()==NULL,std::logic_error
641  ,"Error, The input approxFwdOp object must be fully initialized"
642  " before calling this function!"
643  );
644  // 2007/08/10: rabartl: This next set_extra_data(...) call is likely to be
645  // setting up a circular reference! Since epetra_epetraPrecOp was
646  // gotten from epetraPrecOp, if you set epetraPrecOp as extra data
647  // on the RCP epetra_epetraPrecOp then you have a circular reference!
648  //set_extra_data(
649  // epetraPrecOp, AOOLOWSF_epetraPrecOp_str, &epetra_epetraPrecOp,
650  // Teuchos::POST_DESTROY, false
651  // );
652  overall_epetra_epetraPrecOpTransp
653  = trans_trans(
654  real_trans(wrappedPrecOpTransp),
655  real_trans(epetra_epetraPrecOpTransp)
656  );
657  }
658 
659  //
660  // Determine if the forward and preconditioner operators are a row matrices
661  // or not
662  //
663  RCP<const Epetra_RowMatrix>
664  rowmatrix_epetraFwdOp = rcp_dynamic_cast<const Epetra_RowMatrix>(
665  epetra_epetraFwdOp),
666  rowmatrix_epetraPrecOp = rcp_dynamic_cast<const Epetra_RowMatrix>(
667  epetra_epetraPrecOp);
668  //
669  // Determine the type of preconditoner
670  //
671  // Update useAztecPrec_, input value does not matter
673  enum ELocalPrecType {
674  PT_NONE, PT_AZTEC_FROM_OP, PT_AZTEC_FROM_APPROX_FWD_MATRIX,
675  PT_FROM_PREC_OP, PT_UPPER_BOUND
676  };
677  ELocalPrecType localPrecType = PT_UPPER_BOUND;
678  if( precUsed.get()==NULL && approxFwdOp.get()==NULL && !useAztecPrec_ ) {
679  // No preconditioning at all!
680  localPrecType = PT_NONE;
681  }
682  else if( precUsed.get()==NULL && approxFwdOp.get()==NULL && useAztecPrec_ ) {
683  // We are using the forward matrix for the preconditioner using aztec
684  // preconditioners
685  localPrecType = PT_AZTEC_FROM_OP;
686  }
687  else if( approxFwdOp.get() && useAztecPrec_ ) {
688  // The preconditioner comes from the input as a matrix and we are using
689  // aztec preconditioners
690  localPrecType = PT_AZTEC_FROM_APPROX_FWD_MATRIX;
691  }
692  else if( precUsed.get() ) {
693  // The preconditioner comes as an external operator so let's use it as
694  // such
695  localPrecType = PT_FROM_PREC_OP;
696  }
698  (localPrecType == PT_UPPER_BOUND, std::logic_error,
699  "AztecOOLinearOpWithSolveFactory::initializeOp_impl(...): "
700  "localPrecType == PT_UPPER_BOUND. This means that previously, "
701  "this value might have been used uninitialized. "
702  "Please report this bug to the Stratimikos developers.");
703 
704  //
705  // Determine if aztecOp already contains solvers and if we need to
706  // reinitialize or not
707  //
708  RCP<AztecOO> aztecFwdSolver, aztecAdjSolver;
709  bool startingOver;
710  {
711  // Let's assume that fwdOp, prec and/or approxFwdOp are compatible with
712  // the already created AztecOO objects. If they are not, then the client
713  // should have created a new LOWSB object from scratch!
717  bool old_isExternalPrec;
719  Teuchos::RCP<AztecOO> old_aztecFwdSolver;
720  Teuchos::RCP<AztecOO> old_aztecAdjSolver;
721  double old_aztecSolverScalar;
722  aztecOp->uninitialize(
723  &old_fwdOp
724  ,&old_fwdOpSrc
725  ,&old_prec
726  ,&old_isExternalPrec
727  ,&old_approxFwdOpSrc
728  ,&old_aztecFwdSolver
729  ,NULL
730  ,&old_aztecAdjSolver
731  ,NULL
732  ,&old_aztecSolverScalar
733  );
734  if( old_aztecFwdSolver.get()==NULL ) {
735  // This has never been initialized before
736  startingOver = true;
737  }
738  else {
739  // Let's assume that fwdOp, prec and/or approxFwdOp are compatible with
740  // the already created AztecOO objects. If they are not, then the
741  // client should have created a new LOWSB object from scratch!
742  aztecFwdSolver = old_aztecFwdSolver;
743  aztecAdjSolver = old_aztecAdjSolver;
744  startingOver = false;
745  // We must wipe out the old preconditoner if we are not reusing the
746  // preconditioner
747  Ptr<bool> constructedAztecPreconditioner;
748  if(
749  !reusePrec
750  &&
751  !is_null(constructedAztecPreconditioner = get_optional_nonconst_extra_data<bool>(
752  aztecFwdSolver, "AOOLOWSF::constructedAztecPreconditoner") )
753  &&
754  *constructedAztecPreconditioner
755  )
756  {
757  aztecFwdSolver->DestroyPreconditioner();
758  *constructedAztecPreconditioner = false;
759  }
760  // We must see if we set an external preconditioner but will not do so
761  // again in which case we must blow away AztecOO and start over again!
762  Ptr<bool> setPreconditionerOperator;
763  if(
764  localPrecType != PT_FROM_PREC_OP
765  && !is_null( setPreconditionerOperator = get_optional_nonconst_extra_data<bool>(
766  aztecFwdSolver,"AOOLOWSF::setPreconditonerOperator") )
767  && *setPreconditionerOperator
768  )
769  {
770  // We must start over again since there is no way to unset an external
771  // preconditioner!
772  startingOver = true;
773  }
774  }
775  }
776 
777  //
778  // Create the AztecOO solvers if we are starting over
779  //
780  startingOver = true; // ToDo: Remove this and figure out why this is not working!
781  if(startingOver) {
782  // Forward solver
783  aztecFwdSolver = rcp(new AztecOO());
784  aztecFwdSolver->SetAztecOption(AZ_diagnostics,AZ_none); // This was turned off in NOX?
785  aztecFwdSolver->SetAztecOption(AZ_keep_info,1);
786  // Adjoint solver (if supported)
787  if(
788  epetra_epetraFwdOpAdjointSupport==EPETRA_OP_ADJOINT_SUPPORTED
789  &&
790  localPrecType!=PT_AZTEC_FROM_OP && localPrecType!=PT_AZTEC_FROM_APPROX_FWD_MATRIX
791  )
792  {
793  aztecAdjSolver = rcp(new AztecOO());
794  aztecAdjSolver->SetAztecOption(AZ_diagnostics,AZ_none);
795  //aztecAdjSolver->SetAztecOption(AZ_keep_info,1);
796  }
797  }
798 
799  //
800  // Set the options on the AztecOO solvers
801  //
802  if( startingOver ) {
803  if(paramList_.get())
804  setAztecOOParameters(
805  &paramList_->sublist(ForwardSolve_name).sublist(AztecOO_Settings_name),
806  &*aztecFwdSolver
807  );
808  if(aztecAdjSolver.get() && paramList_.get())
809  setAztecOOParameters(
810  &paramList_->sublist(AdjointSolve_name).sublist(AztecOO_Settings_name),
811  &*aztecAdjSolver
812  );
813  }
814 
815  //
816  // Process the forward operator
817  //
818  RCP<const Epetra_Operator>
819  aztec_epetra_epetraFwdOp,
820  aztec_epetra_epetraAdjOp;
821  // Forward solve
822  RCP<const Epetra_Operator>
823  epetraOps[]
824  = { epetra_epetraFwdOp };
826  epetraOpsTransp[]
827  = { epetra_epetraFwdOpTransp==NOTRANS ? Teuchos::NO_TRANS : Teuchos::TRANS };
828  PO::EApplyMode
829  epetraOpsApplyMode[]
830  = { epetra_epetraFwdOpApplyAs==EPETRA_OP_APPLY_APPLY
831  ? PO::APPLY_MODE_APPLY
832  : PO::APPLY_MODE_APPLY_INVERSE };
833  if(
834  epetraOpsTransp[0] == Teuchos::NO_TRANS
835  &&
836  epetraOpsApplyMode[0] == PO::APPLY_MODE_APPLY
837  )
838  {
839  aztec_epetra_epetraFwdOp = epetra_epetraFwdOp;
840  }
841  else
842  {
843  aztec_epetra_epetraFwdOp = rcp(
844  new PO(1,epetraOps,epetraOpsTransp,epetraOpsApplyMode));
845  }
846  if(
847  startingOver
848  ||
849  aztec_epetra_epetraFwdOp.get() != aztecFwdSolver->GetUserOperator()
850  )
851  {
852  // Here we will be careful not to reset the forward operator in fears that
853  // it will blow out the internally created stuff.
854  aztecFwdSolver->SetUserOperator(
855  const_cast<Epetra_Operator*>(&*aztec_epetra_epetraFwdOp));
856  set_extra_data(
857  aztec_epetra_epetraFwdOp, AOOLOWSF_aztec_epetra_epetraFwdOp_str,
858  Teuchos::inOutArg(aztecFwdSolver), Teuchos::POST_DESTROY, false
859  );
860  }
861  // Adjoint solve
862  if( aztecAdjSolver.get() ) {
863  epetraOpsTransp[0] = (
864  epetra_epetraFwdOpTransp==NOTRANS
867  );
868  if(
869  epetraOpsTransp[0] == Teuchos::NO_TRANS
870  &&
871  epetraOpsApplyMode[0] == PO::APPLY_MODE_APPLY
872  )
873  {
874  aztec_epetra_epetraAdjOp = epetra_epetraFwdOp;
875  }
876  else {
877  aztec_epetra_epetraAdjOp = rcp(
878  new PO(1,epetraOps,epetraOpsTransp,epetraOpsApplyMode));
879  }
880  aztecAdjSolver->SetUserOperator(
881  const_cast<Epetra_Operator*>(&*aztec_epetra_epetraAdjOp));
882  set_extra_data(
883  aztec_epetra_epetraAdjOp, AOOLOWSF_aztec_epetra_epetraAdjOp_str,
884  Teuchos::inOutArg(aztecAdjSolver), Teuchos::POST_DESTROY, false
885  );
886  }
887 
888  //
889  // Process the preconditioner
890  //
891  RCP<const Epetra_Operator>
892  aztec_fwd_epetra_epetraPrecOp,
893  aztec_adj_epetra_epetraPrecOp;
894  bool setAztecPreconditioner = false;
895  switch(localPrecType) {
896  case PT_NONE: {
897  //
898  // No preconditioning at all!
899  //
900  break;
901  }
902  case PT_AZTEC_FROM_OP: {
903  //
904  // We are using the forward matrix for the preconditioner using aztec
905  // preconditioners
906  //
907  if( startingOver || !reusePrec ) {
909  rowmatrix_epetraFwdOp.get()==NULL, std::logic_error,
910  "AztecOOLinearOpWithSolveFactory::initializeOp_impl(...): "
911  "Error, There is no preconditioner given by client, but the client "
912  "passed in an Epetra_Operator for the forward operator of type \'"
913  <<typeName(*epetra_epetraFwdOp)<<"\' that does not "
914  "support the Epetra_RowMatrix interface!"
915  );
917  epetra_epetraFwdOpTransp!=NOTRANS, std::logic_error,
918  "AztecOOLinearOpWithSolveFactory::initializeOp_impl(...):"
919  " Error, There is no preconditioner given by client and the client "
920  "passed in an Epetra_RowMatrix for the forward operator but the "
921  "overall transpose is not NOTRANS and therefore we can can just "
922  "hand this over to aztec without making a copy which is not supported here!"
923  );
924  aztecFwdSolver->SetPrecMatrix(
925  const_cast<Epetra_RowMatrix*>(&*rowmatrix_epetraFwdOp));
926  set_extra_data(
927  rowmatrix_epetraFwdOp, AOOLOWSF_rowmatrix_epetraFwdOp_str,
928  Teuchos::inOutArg(aztecFwdSolver), Teuchos::POST_DESTROY, false
929  );
930  }
931  setAztecPreconditioner = true;
932  break;
933  }
934  case PT_AZTEC_FROM_APPROX_FWD_MATRIX: {
935  //
936  // The preconditioner comes from the input as a matrix and we are using
937  // aztec preconditioners
938  //
939  if( startingOver || !reusePrec ) {
941  rowmatrix_epetraPrecOp.get()==NULL, std::logic_error
942  ,"AztecOOLinearOpWithSolveFactor::initializeOp_impl(...): The client "
943  "passed in an Epetra_Operator for the preconditioner matrix of type \'"
944  <<typeName(*epetra_epetraPrecOp)<<"\' that does not "
945  "support the Epetra_RowMatrix interface!"
946  );
948  overall_epetra_epetraPrecOpTransp!=NOTRANS, std::logic_error
949  ,"AztecOOLinearOpWithSolveFactor::initializeOp_impl(...): Error, The client "
950  "passed in an Epetra_RowMatrix for the preconditoner matrix but the overall "
951  "transpose is not NOTRANS and therefore we can can just "
952  "hand this over to aztec without making a copy which is not supported here!"
953  );
954  aztecFwdSolver->SetPrecMatrix(
955  const_cast<Epetra_RowMatrix*>(&*rowmatrix_epetraPrecOp));
956  set_extra_data(
957  rowmatrix_epetraPrecOp, AOOLOWSF_rowmatrix_epetraPrecOp_str,
958  Teuchos::inOutArg(aztecFwdSolver), Teuchos::POST_DESTROY, false
959  );
960  }
961  setAztecPreconditioner = true;
962  break;
963  }
964  case PT_FROM_PREC_OP: {
965  //
966  // The preconditioner comes as an operator so let's use it as such
967  //
968  // Forward solve
969  RCP<const Epetra_Operator>
970  theEpetraOps[]
971  = { epetra_epetraPrecOp };
973  theEpetraOpsTransp[]
974  = { overall_epetra_epetraPrecOpTransp==NOTRANS
976  : Teuchos::TRANS };
977  // Here we must toggle the apply mode since aztecoo applies the
978  // preconditioner using ApplyInverse(...)
979  PO::EApplyMode
980  theEpetraOpsApplyMode[]
981  = { epetra_epetraPrecOpApplyAs==EPETRA_OP_APPLY_APPLY
982  ? PO::APPLY_MODE_APPLY_INVERSE
983  : PO::APPLY_MODE_APPLY };
984  if(
985  theEpetraOpsTransp[0] == Teuchos::NO_TRANS
986  &&
987  epetra_epetraPrecOpApplyAs==EPETRA_OP_APPLY_APPLY_INVERSE
988  )
989  {
990  aztec_fwd_epetra_epetraPrecOp = epetra_epetraPrecOp;
991  }
992  else {
993  aztec_fwd_epetra_epetraPrecOp = rcp(new PO(1,theEpetraOps,theEpetraOpsTransp,theEpetraOpsApplyMode));
994  }
995  aztecFwdSolver->SetPrecOperator(
996  const_cast<Epetra_Operator*>(&*aztec_fwd_epetra_epetraPrecOp));
997  set_extra_data(
998  aztec_fwd_epetra_epetraPrecOp, AOOLOWSF_aztec_fwd_epetra_epetraPrecOp_str,
999  Teuchos::inOutArg(aztecFwdSolver), Teuchos::POST_DESTROY, false
1000  );
1001  // Adjoint solve
1002  if(
1003  aztecAdjSolver.get()
1004  &&
1005  epetra_epetraPrecOpAdjointSupport == EPETRA_OP_ADJOINT_SUPPORTED
1006  )
1007  {
1008  theEpetraOpsTransp[0] = (
1009  overall_epetra_epetraPrecOpTransp==NOTRANS
1010  ? Teuchos::TRANS
1012  );
1013  if(
1014  theEpetraOpsTransp[0] == Teuchos::NO_TRANS
1015  &&
1016  epetra_epetraPrecOpApplyAs==EPETRA_OP_APPLY_APPLY_INVERSE
1017  )
1018  {
1019  aztec_adj_epetra_epetraPrecOp = epetra_epetraPrecOp;
1020  }
1021  else {
1022  aztec_adj_epetra_epetraPrecOp = rcp(
1023  new PO(1,theEpetraOps,theEpetraOpsTransp,theEpetraOpsApplyMode));
1024  }
1025  aztecAdjSolver->SetPrecOperator(
1026  const_cast<Epetra_Operator*>(&*aztec_adj_epetra_epetraPrecOp));
1027  set_extra_data(
1028  aztec_adj_epetra_epetraPrecOp, AOOLOWSF_aztec_adj_epetra_epetraPrecOp_str,
1029  Teuchos::inOutArg(aztecAdjSolver), Teuchos::POST_DESTROY, false
1030  );
1031  set_extra_data<bool>(
1032  true, AOOLOWSF_setPrecondtionerOperator_str,
1033  Teuchos::inOutArg(aztecFwdSolver), Teuchos::POST_DESTROY, false
1034  );
1035  }
1036  break;
1037  }
1038  default:
1040  }
1041 
1042  //
1043  // Initialize the interal aztec preconditioner
1044  //
1045  if(setAztecPreconditioner) {
1046  if( startingOver || !reusePrec ) {
1047  double condNumEst = -1.0;
1048  TEUCHOS_TEST_FOR_EXCEPT(0!=aztecFwdSolver->ConstructPreconditioner(condNumEst));
1049  //aztecFwdSolver->SetAztecOption(AZ_pre_calc, AZ_calc);
1050  set_extra_data<bool>(
1051  true, AOOLOWSF_constructedAztecPreconditoner_str,
1052  Teuchos::inOutArg(aztecFwdSolver), Teuchos::POST_DESTROY, false
1053  );
1054  }
1055  else {
1056  //aztecFwdSolver->SetAztecOption(AZ_pre_calc, AZ_reuse);
1057  }
1058  }
1059 
1060  //
1061  // Initialize the AztecOOLinearOpWithSolve object and set its options
1062  //
1063  if(aztecAdjSolver.get() && aztecAdjSolver->GetPrecOperator()) {
1064  aztecOp->initialize(
1065  fwdOp, fwdOpSrc,precUsed, prec.get()!=NULL, approxFwdOpSrc,
1066  aztecFwdSolver, true, aztecAdjSolver, true, epetra_epetraFwdOpScalar
1067  );
1068  }
1069  else {
1070  aztecOp->initialize(
1071  fwdOp, fwdOpSrc, precUsed, prec.get()!=NULL, approxFwdOpSrc,
1072  aztecFwdSolver, true, null, false, epetra_epetraFwdOpScalar
1073  );
1074  }
1075  aztecOp->fwdDefaultMaxIterations(defaultFwdMaxIterations_);
1076  aztecOp->fwdDefaultTol(defaultFwdTolerance_);
1077  aztecOp->adjDefaultMaxIterations(defaultAdjMaxIterations_);
1078  aztecOp->adjDefaultTol(defaultAdjTolerance_);
1079  aztecOp->outputEveryRhs(outputEveryRhs_);
1080  aztecOp->setOStream(this->getOStream());
1081  if(!is_null(this->getOverridingOStream()))
1082  aztecOp->setOverridingOStream(this->getOverridingOStream());
1083  aztecOp->setVerbLevel(this->getVerbLevel());
1084 
1085 #ifdef TEUCHOS_DEBUG
1086  if(paramList_.get())
1087  paramList_->validateParameters(*this->getValidParameters());
1088 #endif
1089 
1090  if(out.get() && static_cast<int>(verbLevel) > static_cast<int>(Teuchos::VERB_LOW))
1091  *out << "\nLeaving Thyra::AztecOOLinearOpWithSolveFactory::initializeOp_impl(...) ...\n";
1092 
1093 }
1094 
1095 
1096 } // namespace Thyra
1097 
1098 
1099 #endif // SUN_CXX
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
void initializeAndReuseOp(const Teuchos::RCP< const LinearOpSourceBase< double > > &fwdOpSrc, LinearOpWithSolveBase< double > *Op) const
bool is_null(const boost::shared_ptr< T > &p)
EOpTransp
void unwrap(const LinearOpBase< Scalar > &Op, Scalar *scalar, EOpTransp *transp, const LinearOpBase< Scalar > **origOp)
Teuchos::RCP< const Teuchos::ParameterList > getParameterList() const
T & get(const std::string &name, T def_value)
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
RCP< const LinearOpSourceBase< double > > extract_fwdOpSrc()
Extract the forward LinearOpBase&lt;double&gt; object so that it can be modified.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
bool isExternalPrec() const
Determine if the preconditioner was external or not.
T_To & dyn_cast(T_From &from)
AztecOOLinearOpWithSolveFactory(Teuchos::RCP< Teuchos::ParameterList > const &paramList=Teuchos::null)
Construct uninitialized.
void initializeApproxPreconditionedOp(const Teuchos::RCP< const LinearOpSourceBase< double > > &fwdOpSrc, const Teuchos::RCP< const LinearOpSourceBase< double > > &approxFwdOpSrc, LinearOpWithSolveBase< double > *Op, const ESupportSolveUse supportSolveUse) const
EOpTransp real_trans(EOpTransp transp)
T * get() const
bool supportsPreconditionerInputType(const EPreconditionerInputType precOpType) const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
bool isCompatible(const LinearOpSourceBase< double > &fwdOpSrc) const
bool isSublist(const std::string &name) const
EAdjointEpetraOp
Teuchos::RCP< Teuchos::ParameterList > unsetParameterList()
EOpTransp trans_trans(EOpTransp trans1, EOpTransp trans2)
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
int DestroyPreconditioner()
void unsetPreconditionerFactory(Teuchos::RCP< PreconditionerFactoryBase< double > > *precFactory, std::string *precFactoryName)
void initializeOp(const Teuchos::RCP< const LinearOpSourceBase< double > > &fwdOpSrc, LinearOpWithSolveBase< double > *Op, const ESupportSolveUse supportSolveUse) const
T * getPtr(const std::string &name)
void setParameterList(Teuchos::RCP< Teuchos::ParameterList > const &paramList)
Teuchos::RCP< PreconditionerFactoryBase< double > > getPreconditionerFactory() const
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
void initializePreconditionedOp(const Teuchos::RCP< const LinearOpSourceBase< double > > &fwdOpSrc, const Teuchos::RCP< const PreconditionerBase< double > > &prec, LinearOpWithSolveBase< double > *Op, const ESupportSolveUse supportSolveUse) const
#define TEUCHOS_UNREACHABLE_RETURN(dummyReturnVal)
RCP< const PreconditionerBase< double > > extract_prec()
Extract the preconditioner.
EApplyEpetraOpAs
EPreconditionerInputType
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
virtual Teuchos::RCP< const LinearOpBase< Scalar > > getOp() const =0
void setPreconditionerFactory(const Teuchos::RCP< PreconditionerFactoryBase< double > > &precFactory, const std::string &precFactoryName)
Teuchos::RCP< Teuchos::ParameterList > getNonconstParameterList()
Teuchos::RCP< LinearOpWithSolveBase< double > > createOp() const
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)
std::string typeName(const T &t)
Concrete LinearOpWithSolveBase subclass implemented using AztecOO.
RCP< const LinearOpSourceBase< double > > extract_approxFwdOpSrc()
Extract the approximate forward LinearOpBase&lt;double&gt; object used to build the preconditioner.

Generated on Fri Jun 5 2020 10:13:28 for Stratimikos by doxygen 1.8.5