10#ifndef STRATIMIKOS_LINEARSOLVERBUILDER_DEF_HPP 
   11#define STRATIMIKOS_LINEARSOLVERBUILDER_DEF_HPP 
   15#include "Stratimikos_InternalConfig.h" 
   16#include "Thyra_DelayedLinearOpWithSolveFactory.hpp" 
   17#include "Teuchos_AbstractFactoryStd.hpp" 
   18#include "Teuchos_CommandLineProcessor.hpp" 
   19#include "Teuchos_XMLParameterListHelpers.hpp" 
   20#include "Teuchos_GlobalMPISession.hpp" 
   22#ifdef HAVE_STRATIMIKOS_AMESOS 
   23#  include "Thyra_AmesosLinearOpWithSolveFactory.hpp" 
   25#ifdef HAVE_STRATIMIKOS_AMESOS2 
   26#  include "Thyra_Amesos2LinearOpWithSolveFactory.hpp" 
   28#if defined(HAVE_STRATIMIKOS_EPETRAEXT) && defined(HAVE_STRATIMIKOS_AZTECOO) 
   29#  include "Thyra_AztecOOLinearOpWithSolveFactory.hpp" 
   31#ifdef HAVE_STRATIMIKOS_BELOS 
   32#  include "Thyra_BelosLinearOpWithSolveFactory.hpp" 
   34#ifdef HAVE_STRATIMIKOS_IFPACK 
   35#  include "Thyra_IfpackPreconditionerFactory.hpp" 
   37#if defined(HAVE_STRATIMIKOS_IFPACK2) && defined(HAVE_STRATIMIKOS_THYRATPETRAADAPTERS) 
   38#  include "Thyra_Ifpack2PreconditionerFactory.hpp" 
   39#  include "Tpetra_CrsMatrix.hpp" 
   41#ifdef HAVE_STRATIMIKOS_ML 
   42#  include "Thyra_MLPreconditionerFactory.hpp" 
   49const std::string LinearSolverType_name    = 
"Linear Solver Type";
 
   50const std::string LinearSolverTypes_name   = 
"Linear Solver Types";
 
   51const std::string PreconditionerType_name    = 
"Preconditioner Type";
 
   52const std::string PreconditionerTypes_name   = 
"Preconditioner Types";
 
   53const std::string None_name = 
"None";
 
   54const std::string EnableDelayedSolverConstruction_name = 
"Enable Delayed Solver Construction";
 
   55const bool EnableDelayedSolverConstruction_default = 
false;
 
   61namespace Stratimikos {
 
   69  const std::string ¶msXmlFileName_in,
 
   70  const std::string &extraParamsXmlString_in,
 
   71  const std::string ¶msUsedXmlOutFileName_in,
 
   72  const std::string ¶msXmlFileNameOption_in,
 
   73  const std::string &extraParamsXmlStringOption_in,
 
   74  const std::string ¶msUsedXmlOutFileNameOption_in,
 
   75  const bool &replaceDuplicateFactories_in
 
   77  :paramsXmlFileName_(paramsXmlFileName_in)
 
   78  ,extraParamsXmlString_(extraParamsXmlString_in)
 
   79  ,paramsUsedXmlOutFileName_(paramsUsedXmlOutFileName_in)
 
   80  ,paramsXmlFileNameOption_(paramsXmlFileNameOption_in)
 
   81  ,extraParamsXmlStringOption_(extraParamsXmlStringOption_in)
 
   82  ,paramsUsedXmlOutFileNameOption_(paramsUsedXmlOutFileNameOption_in)
 
   83  ,replaceDuplicateFactories_(replaceDuplicateFactories_in)
 
   84  ,enableDelayedSolverConstruction_(EnableDelayedSolverConstruction_default)
 
   86  this->initializeDefaults();
 
 
   95  if (nonnull(paramList_)) {
 
   96    paramList_->validateParameters(*this->getValidParameters());
 
 
  102template<
class Scalar>
 
  104  const RCP<
const AbstractFactory<Thyra::LinearOpWithSolveFactoryBase<Scalar> > >
 
  105  &solveStrategyFactory,
 
  106  const std::string &solveStrategyName,
 
  107  const bool makeDefault
 
  110  const int existingNameIdx =
 
  111    this->getAndAssertExistingFactoryNameIdx(
"setLinearSolveStrategyFactory()",
 
  112      validLowsfNames_(), solveStrategyName);
 
  113  if (existingNameIdx >= 0) {
 
  114    validLowsfNames_[existingNameIdx] = solveStrategyName;
 
  115    lowsfArray_[existingNameIdx] = solveStrategyFactory;
 
  118    validLowsfNames_.push_back(solveStrategyName);
 
  119    lowsfArray_.push_back(solveStrategyFactory);
 
  121  validParamList_ = Teuchos::null;
 
  123    setDefaultLinearSolveStrategyFactoryName(solveStrategyName);
 
 
  128template<
class Scalar>
 
  130  const std::string &solveStrategyName)
 
  132  defaultLOWSF_ = solveStrategyName;
 
 
  136template<
class Scalar>
 
  138  const RCP<
const AbstractFactory<Thyra::PreconditionerFactoryBase<Scalar>>>
 
  139    &precStrategyFactory,
 
  140  const std::string &precStrategyName,
 
  141  const bool makeDefault
 
  144  const int existingNameIdx =
 
  145    this->getAndAssertExistingFactoryNameIdx(
"setPreconditioningStrategyFactory()",
 
  146      validPfNames_(), precStrategyName);
 
  147  if (existingNameIdx >= 0) {
 
  148    validPfNames_[existingNameIdx] = precStrategyName;
 
  149    pfArray_[existingNameIdx-1] = precStrategyFactory; 
 
  152    validPfNames_.push_back(precStrategyName);
 
  153    pfArray_.push_back(precStrategyFactory);
 
  155  validParamList_ = Teuchos::null;
 
  157    setDefaultPreconditioningStrategyFactoryName(precStrategyName);
 
 
  162template<
class Scalar>
 
  164  const std::string &precStrategyName)
 
  166  defaultPF_ = precStrategyName;
 
 
  170template<
class Scalar>
 
  173  TEUCHOS_TEST_FOR_EXCEPT(clp==NULL);
 
  175    paramsXmlFileNameOption().c_str(),¶msXmlFileName_
 
  176    ,
"Name of an XML file containing parameters for linear solver " 
  177    "options to be appended first." 
  180    extraParamsXmlStringOption().c_str(),&extraParamsXmlString_
 
  181    ,
"An XML string containing linear solver parameters to be appended second." 
  184    paramsUsedXmlOutFileNameOption().c_str(),¶msUsedXmlOutFileName_
 
  185    ,
"Name of an XML file that can be written with the parameter list after it " 
  186    "has been used on completion of this program." 
 
  191template<
class Scalar>
 
  194  using Teuchos::parameterList;
 
  196  using Teuchos::updateParametersFromXmlFile;
 
  197  using Teuchos::updateParametersFromXmlString;
 
  200  if (!paramList_.get()) {
 
  201    paramList_ = parameterList(
"LinearSolverBuilder");
 
  203  if (paramsXmlFileName().length()) {
 
  205      *out << endl << 
"Reading parameters from XML file \""  
  206           << paramsXmlFileName() << 
"\" ..." << endl;
 
  208    updateParametersFromXmlFile (paramsXmlFileName (), paramList_.ptr());
 
  210  if (extraParamsXmlString().length()) {
 
  212      *out << endl << 
"Appending extra parameters from the XML string \"" 
  213           << extraParamsXmlString() << 
"\" ..." << endl;
 
  215    updateParametersFromXmlString (extraParamsXmlString (), paramList_.ptr());
 
  217  setParameterList(paramList_);
 
 
  221template<
class Scalar>
 
  223  const Thyra::LinearOpWithSolveFactoryBase<Scalar> &,
 
  224  const std::string &outputXmlFileName
 
  227  justInTimeInitialize();
 
  228  const std::string xmlOutputFile =
 
  229    ( outputXmlFileName.length() ? outputXmlFileName : paramsUsedXmlOutFileName() );
 
  230  if (xmlOutputFile.length()) {
 
  231    Teuchos::writeParameterListToXmlFile(*paramList_, xmlOutputFile);
 
 
  236template<
class Scalar>
 
  240  justInTimeInitialize();
 
  241  return lowsfValidator_->getStringValue(*paramList_, LinearSolverType_name,
 
 
  246template<
class Scalar>
 
  250  justInTimeInitialize();
 
  251  return pfValidator_->getStringValue(*paramList_, PreconditionerType_name,
 
 
  259template<
class Scalar>
 
  261  RCP<Teuchos::ParameterList> 
const& paramList
 
  264  TEUCHOS_TEST_FOR_EXCEPT(is_null(paramList));
 
  265  paramList->validateParameters(*this->getValidParameters());
 
  266  paramList_ = paramList;
 
  267  enableDelayedSolverConstruction_ = paramList_->get(
 
  268    EnableDelayedSolverConstruction_name, EnableDelayedSolverConstruction_default );
 
 
  272template<
class Scalar>
 
  273RCP<Teuchos::ParameterList>
 
  280template<
class Scalar>
 
  281RCP<Teuchos::ParameterList>
 
  284  RCP<Teuchos::ParameterList> _paramList = paramList_;
 
  285  paramList_ = Teuchos::null;
 
 
  290template<
class Scalar>
 
  291RCP<const Teuchos::ParameterList>
 
  298template<
class Scalar>
 
  299RCP<const Teuchos::ParameterList>
 
  302  using Teuchos::rcp_implicit_cast;
 
  303  typedef Teuchos::ParameterEntryValidator PEV;
 
  304  if (is_null(validParamList_)) {
 
  305    RCP<Teuchos::ParameterList>
 
  306      validParamList = Teuchos::rcp(
new Teuchos::ParameterList);
 
  308    lowsfValidator_ = Teuchos::rcp(
 
  309      new Teuchos::StringToIntegralParameterEntryValidator<int>(
 
  310        validLowsfNames_,LinearSolverType_name
 
  314      LinearSolverType_name, defaultLOWSF_,
 
  315      (std::string(
"Determines the type of linear solver that will be used.\n")
 
  316        + 
"The parameters for each solver type are specified in the sublist \"" 
  317        + LinearSolverTypes_name + 
"\"").c_str(),
 
  318      rcp_implicit_cast<const PEV>(lowsfValidator_)
 
  320    Teuchos::ParameterList &linearSolverTypesSL = validParamList->sublist(
 
  321      LinearSolverTypes_name,
false,
 
  322      "Sublists for each of the linear solver types set using the parameter\n" 
  323      "\"" + LinearSolverType_name + 
"\".  Note that the options for each\n" 
  324      "linear solver type given below will only be used if linear solvers\n" 
  325      "of that type are created.  It is fine to list parameter sublists for\n" 
  326      "linear solver types that are not used." 
  328    for( 
int i = 0; i < static_cast<int>(lowsfArray_.size()); ++i ) {
 
  330        &lsname = validLowsfNames_[i];
 
  331      const RCP<Thyra::LinearOpWithSolveFactoryBase<Scalar> >
 
  332        lowsf = lowsfArray_[i]->create();
 
  333      linearSolverTypesSL.sublist(lsname).setParameters(*lowsf->getValidParameters()
 
  334        ).disableRecursiveValidation();
 
  337    pfValidator_ = Teuchos::rcp(
 
  338      new Teuchos::StringToIntegralParameterEntryValidator<int>(
 
  339        validPfNames_, PreconditionerType_name ) );
 
  341      PreconditionerType_name, defaultPF_,
 
  342      (std::string(
"Determines the type of preconditioner that will be used.\n")
 
  343        + 
"This option is only meaningful for linear solvers that accept preconditioner" 
  344        + 
" factory objects!\n" 
  345        + 
"The parameters for each preconditioner are specified in the sublist \"" 
  346        + PreconditionerTypes_name + 
"\"").c_str(),
 
  347      rcp_implicit_cast<const PEV>(pfValidator_)
 
  349    Teuchos::ParameterList &precTypesSL = validParamList->sublist(
 
  350        PreconditionerTypes_name,
false,
 
  351        "Sublists for each of the preconditioner types set using the parameter\n" 
  352        "\"" + PreconditionerType_name + 
"\".  Note that the options for each\n" 
  353        "preconditioner type given below will only be used if preconditioners\n" 
  354        "of that type are created.  It is fine to list parameter sublists for\n" 
  355        "preconditioner types that are not used." 
  357    for( 
int i = 0; i < static_cast<int>(pfArray_.size()); ++i ) {
 
  359        &pfname = validPfNames_[i+1]; 
 
  360      const RCP<Thyra::PreconditionerFactoryBase<Scalar> >
 
  361        pf = pfArray_[i]->create();
 
  362      precTypesSL.sublist(pfname).setParameters(*pf->getValidParameters()
 
  363        ).disableRecursiveValidation();
 
  367      EnableDelayedSolverConstruction_name, EnableDelayedSolverConstruction_default,
 
  368      "When this option is set to true, the linear solver factory will be wrapped\n" 
  369      "in a delayed evaluation Decorator factory object.  This results in a delay\n" 
  370      "in the creation of a linear solver (and the associated preconditioner) until\n" 
  371      "the first solve is actually performed.  This helps in cases where it is not\n" 
  372      "known a-priori if a linear solve will be needed on a given linear operator and\n" 
  373      "therefore can significantly improve performance for some types of algorithms\n" 
  374      "such as NOX and LOCA." 
  377    validParamList_ = validParamList;
 
  379  return validParamList_;
 
 
  386template<
class Scalar>
 
  387RCP<Thyra::LinearOpWithSolveFactoryBase<Scalar> >
 
  389  const std::string &linearSolveStrategyName
 
  392  justInTimeInitialize();
 
  395#ifdef THYRA_DEFAULT_REAL_LINEAR_SOLVER_BUILDER_DUMP 
  396  std::cout << 
"\nEntering LinearSolverBuilder" 
  397            << 
"::createLinearSolveStrategy(...) ...\n";
 
  398  std::cout << 
"\nlinearSolveStrategyName = \"" 
  399            << linearSolveStrategyName << 
"\"\n";
 
  400  std::cout << 
"\nlinearSolveStrategyName.length() = " 
  401            << linearSolveStrategyName.length() << 
"\n";
 
  402  std::cout << 
"\ndefaultLOWSF_ = \"" << defaultLOWSF_ << 
"\"\n";
 
  403  std::cout << 
"\nthis->getLinearSolveStrategyName() = \"" 
  404            << this->getLinearSolveStrategyName() << 
"\"\n";
 
  407    lsname = ( linearSolveStrategyName.length()
 
  408             ? linearSolveStrategyName
 
  409             : this->getLinearSolveStrategyName() );
 
  410#ifdef THYRA_DEFAULT_REAL_LINEAR_SOLVER_BUILDER_DUMP 
  411  std::cout << 
"\nlsname = \"" << lsname << 
"\"\n";
 
  416    ls_idx = lowsfValidator_->getIntegralValue(lsname, LinearSolverType_name);
 
  419  RCP<Thyra::LinearOpWithSolveFactoryBase<Scalar> >
 
  420    lowsf = lowsfArray_[ls_idx]->create();
 
  423  if(lowsf->acceptsPreconditionerFactory()) {
 
  424    const std::string &pfName = this->getPreconditionerStrategyName();
 
  425    RCP<Thyra::PreconditionerFactoryBase<Scalar> >
 
  426      pf = this->createPreconditioningStrategy(pfName);
 
  428      lowsf->setPreconditionerFactory(pf,pfName);
 
  433  lowsf->setParameterList(
 
  434    sublist(sublist(paramList_, LinearSolverTypes_name), lsname));
 
  436  if (enableDelayedSolverConstruction_) {
 
  438      new Thyra::DelayedLinearOpWithSolveFactory<Scalar>(lowsf)
 
 
  447template<
class Scalar>
 
  448RCP<Thyra::PreconditionerFactoryBase<Scalar> >
 
  450  const std::string &preconditioningStrategyName
 
  453  justInTimeInitialize();
 
  457    pfname = ( preconditioningStrategyName.length()
 
  458             ? preconditioningStrategyName
 
  459             : this->getPreconditionerStrategyName() );
 
  460  RCP<Thyra::PreconditionerFactoryBase<Scalar> >
 
  465    pf_idx = pfValidator_->getIntegralValue(pfname, PreconditionerType_name);
 
  467    pf = pfArray_[pf_idx-1]->create(); 
 
  468    pf->setParameterList(
 
  469      sublist(sublist(paramList_, PreconditionerTypes_name), pfname));
 
 
  479template<
class Scalar>
 
  484  using Teuchos::abstractFactoryStd;
 
  487  defaultPF_ = None_name;
 
  488  validLowsfNames_.resize(0);
 
  489  validPfNames_.resize(0);
 
  490  validPfNames_.push_back(None_name); 
 
  496#ifdef HAVE_STRATIMIKOS_AMESOS2 
  497  setLinearSolveStrategyFactory(
 
  498    abstractFactoryStd<Thyra::LinearOpWithSolveFactoryBase<Scalar>,
 
  499    Thyra::Amesos2LinearOpWithSolveFactory<Scalar>>(),
 
  504#ifdef HAVE_STRATIMIKOS_BELOS 
  505  setLinearSolveStrategyFactory(
 
  506    abstractFactoryStd<Thyra::LinearOpWithSolveFactoryBase<Scalar>,
 
  512#if defined(HAVE_STRATIMIKOS_IFPACK2) && defined(HAVE_STRATIMIKOS_THYRATPETRAADAPTERS) 
  513  setPreconditioningStrategyFactory(
 
  514    abstractFactoryStd<Thyra::PreconditionerFactoryBase<Scalar>,
 
  515    Thyra::Ifpack2PreconditionerFactory<Tpetra::CrsMatrix<Scalar>>>(),
 
  526void LinearSolverBuilder<double>::initializeDefaults()
 
  528  using Scalar = double;
 
  530  using Teuchos::abstractFactoryStd;
 
  533  defaultPF_ = None_name;
 
  534  validLowsfNames_.resize(0);
 
  535  validPfNames_.resize(0);
 
  536  validPfNames_.push_back(None_name); 
 
  542#ifdef HAVE_STRATIMIKOS_AMESOS2 
  543  setLinearSolveStrategyFactory(
 
  544    abstractFactoryStd<Thyra::LinearOpWithSolveFactoryBase<Scalar>,
 
  545    Thyra::Amesos2LinearOpWithSolveFactory<Scalar>>(),
 
  550#ifdef HAVE_STRATIMIKOS_BELOS 
  551  setLinearSolveStrategyFactory(
 
  552    abstractFactoryStd<Thyra::LinearOpWithSolveFactoryBase<Scalar>,
 
  558#ifdef HAVE_STRATIMIKOS_AMESOS 
  559  setLinearSolveStrategyFactory(
 
  560    abstractFactoryStd<Thyra::LinearOpWithSolveFactoryBase<Scalar>,
 
  566#if defined(HAVE_STRATIMIKOS_EPETRAEXT) && defined(HAVE_STRATIMIKOS_AZTECOO) 
  567  setLinearSolveStrategyFactory(
 
  568    abstractFactoryStd<Thyra::LinearOpWithSolveFactoryBase<Scalar>,
 
  577#ifdef HAVE_STRATIMIKOS_AMESOS 
  578  if (Teuchos::GlobalMPISession::getNProc() == 1) {
 
  579    setDefaultLinearSolveStrategyFactoryName(
"Amesos");
 
  587#ifdef HAVE_STRATIMIKOS_ML 
  588  setPreconditioningStrategyFactory(
 
  589    abstractFactoryStd<Thyra::PreconditionerFactoryBase<Scalar>,
 
  595#if defined(HAVE_STRATIMIKOS_IFPACK2) && defined(HAVE_STRATIMIKOS_THYRATPETRAADAPTERS) 
  596  setPreconditioningStrategyFactory(
 
  597    abstractFactoryStd<Thyra::PreconditionerFactoryBase<Scalar>,
 
  598    Thyra::Ifpack2PreconditionerFactory<Tpetra::CrsMatrix<Scalar>>>(),
 
  603#ifdef HAVE_STRATIMIKOS_IFPACK 
  604  setPreconditioningStrategyFactory(
 
  605    abstractFactoryStd<Thyra::PreconditionerFactoryBase<Scalar>,
 
  616template<
class Scalar>
 
  617void LinearSolverBuilder<Scalar>::justInTimeInitialize()
 const 
  619  paramList_.assert_not_null();
 
  620  if (is_null(validParamList_)) {
 
  622    this->getValidParameters();
 
  627template<
class Scalar>
 
  628int LinearSolverBuilder<Scalar>::getAndAssertExistingFactoryNameIdx(
 
  629  const std::string &setFunctionName, 
const Teuchos::ArrayView<std::string> namesArray,
 
  630  const std::string &name)
 const 
  632  const int existingNameIdx =
 
  633    LinearSolverBuilderHelpers::existingNameIndex(namesArray, name);
 
  634  TEUCHOS_TEST_FOR_EXCEPTION(
 
  635    (!replaceDuplicateFactories_ && existingNameIdx >= 0), std::logic_error,
 
  636    "ERROR: "<<setFunctionName<<
": the name='"<<name<<
"' already exists" 
  637    << 
" at index="<<existingNameIdx<<
"!\n" 
  639    << 
"TIP: To allow duplicates, change the property replaceDuplicateFactories from" 
  640    << 
" false to true by calling <thisObject>.replaceDuplicateFactories(true)" 
  641    << 
" where <thisObject> is of type Stratimikos::LinearSolverBuilder<Scalar>!");
 
  642  return existingNameIdx;
 
  653#define STRATIMIKOS_LINEARSOLVERBUILDER_INSTANT(SCALAR) \ 
  655  template class LinearSolverBuilder<SCALAR >; 
Concrete subclass of Thyra::LinearSolverBuilderBase for creating Thyra::LinearOpWithSolveFactoryBase ...
void setLinearSolveStrategyFactory(const RCP< const AbstractFactory< Thyra::LinearOpWithSolveFactoryBase< Scalar > > > &solveStrategyFactory, const std::string &solveStrategyName, const bool makeDefault=false)
Set a new linear solver strategy factory object.
RCP< Thyra::PreconditionerFactoryBase< Scalar > > createPreconditioningStrategy(const std::string &preconditioningStrategyName) const
RCP< const ParameterList > getParameterList() const
std::string getLinearSolveStrategyName() const
Get the name of the linear solver strategy that will be created on the next call to this->createLinea...
void setDefaultPreconditioningStrategyFactoryName(const std::string &precStrategyName)
Set the default linear solver factory name.
RCP< const ParameterList > getValidParameters() const
LinearSolverBuilder(const std::string ¶msXmlFileName="", const std::string &extraParamsXmlString="", const std::string ¶msUsedXmlOutFileName="", const std::string ¶msXmlFileNameOption="linear-solver-params-file", const std::string &extraParamsXmlStringOption="extra-linear-solver-params", const std::string ¶msUsedXmlOutFileNameOption="linear-solver-params-used-file", const bool &replaceDuplicateFactories=true)
Construct with default parameters.
RCP< Thyra::LinearOpWithSolveFactoryBase< Scalar > > createLinearSolveStrategy(const std::string &linearSolveStrategyName) const
std::string getPreconditionerStrategyName() const
Get the name of the preconditioner strategy that will be created on the next call to this->createPrec...
void setParameterList(RCP< ParameterList > const ¶mList)
void setDefaultLinearSolveStrategyFactoryName(const std::string &solveStrategyName)
Set the default linear solver factory name.
void setupCLP(Teuchos::CommandLineProcessor *clp)
Setup the command-line processor to read in the needed data to extra the parameters from.
RCP< ParameterList > getNonconstParameterList()
void readParameters(std::ostream *out)
Force the parameters to be read from a file and/or an extra XML string.
void setPreconditioningStrategyFactory(const RCP< const AbstractFactory< Thyra::PreconditionerFactoryBase< Scalar > > > &precStrategyFactory, const std::string &precStrategyName, const bool makeDefault=false)
Set a new preconditioner strategy factory object.
RCP< ParameterList > unsetParameterList()
void writeParamsFile(const Thyra::LinearOpWithSolveFactoryBase< Scalar > &lowsFactory, const std::string &outputXmlFileName="") const
Write the parameters list for a LinearOpWithSolveFactoryBase object to a file after the parameters ar...
Concrete LinearOpWithSolveFactoryBase adapter subclass that uses Amesos direct solvers.
LinearOpWithSolveFactoryBase subclass implemented in terms of AztecOO.
LinearOpWithSolveFactoryBase subclass implemented in terms of Belos.
Concrete preconditioner factory subclass based on Ifpack.
Concrete preconditioner factory subclass based on ML.