test_single_stratimikos_solver.cpp

The bottom of this page shows a simple test program for how the Thyra::DefaultRealLinearSolverBuilder class is used to take a parameter list (e.g. read from a file) and use it to build Thyra::LinearOpWithSolveFactory and Thyra::PreconditionerFactoryBase objects.

The arguments that the testing program driver that calls this example are:

Teuchos::GlobalMPISession::GlobalMPISession(): started serial run

Echoing the command-line:

./test_single_stratimikos_solver_driver.exe --echo-command-line --help 

Usage: ./test_single_stratimikos_solver_driver.exe [options]
  options:
  --help                         Prints this help message
  --pause-for-debugging          Pauses for user input to allow attaching a debugger
  --echo-command-line            Echo the command-line but continue as normal
  --input-file           string  Input file [Required].
                                 (default: --input-file="")
  --extra-params         string  Extra parameters overriding the parameters read in from --input-file
                                 (default: --extra-params="")

DETAILED DOCUMENTATION:

Testing program for Trilinos (and non-Trilinos) linear solvers access through Thyra.

Sample output looks like:

Teuchos::GlobalMPISession::GlobalMPISession(): started serial run

Echoing the command-line:

./test_single_stratimikos_solver_driver.exe --echo-command-line --input-file=FourByFour.belos.ifpack.xml 


Reading parameters from XML file "FourByFour.belos.ifpack.xml" ...

Echoing input parameters ...
  Matrix File : string = FourByFour.mtx
  Linear Solver Builder -> 
    Linear Solver Type : string = Belos
    Preconditioner Type : string = Ifpack
  LinearOpTester -> 
    [empty list]
  LinearOpWithSolveTester -> 
    [empty list]

Validating top-level input parameters ...

Reading in an epetra matrix A from the file 'FourByFour.mtx' ...

Creating a Thyra::DefaultRealLinearSolverBuilder object ...

Valid parameters for DefaultRealLinearSolverBuilder ...
  Enable Delayed Solver Construction : bool = 0
  Linear Solver Type : string = Amesos
  Preconditioner Type : string = Ifpack
  Linear Solver Types -> 
    Amesos -> 
      Refactorization Policy : string = RepivotOnRefactorization
      Solver Type : string = Klu
      Throw on Preconditioner Input : bool = 1
      Amesos Settings -> 
        AddToDiag : double = 0
        AddZeroToDiag : bool = 0
        ComputeTrueResidual : bool = 0
        ComputeVectorNorms : bool = 0
        DebugLevel : int = 0
        MatrixProperty : string = general
        MaxProcs : int = -1
        NoDestroy : bool = 0
        OutputLevel : int = 1
        PrintTiming : bool = 0
        RcondThreshold : double = 1e-12
        Redistribute : bool = 0
        Refactorize : bool = 0
        Reindex : int = 0
        ScaleMethod : int = 0
        TrustMe : bool = 0
        Lapack -> 
          Equilibrate : bool = 1
        Mumps -> 
          ColScaling : double* = 0
          Equilibrate : bool = 1
          RowScaling : double* = 0
        Pardiso -> 
          IPARM(1) : int = 0
          IPARM(10) : int = 0
          IPARM(11) : int = 0
          IPARM(18) : int = 0
          IPARM(19) : int = 0
          IPARM(2) : int = 0
          IPARM(21) : int = 0
          IPARM(3) : int = 0
          IPARM(4) : int = 0
          IPARM(8) : int = 0
          MSGLVL : int = 0
        Scalapack -> 
          2D distribution : bool = 1
          grid_nb : int = 32
        Superludist -> 
          ColPerm : string = NOT SET
          Equil : bool = 0
          Fact : string = SamePattern
          IterRefine : string = NOT SET
          PrintNonzeros : bool = 0
          ReplaceTinyPivot : bool = 1
          ReuseSymbolic : bool = 0
          RowPerm : string = NOT SET
          perm_c : int* = 0
          perm_r : int* = 0
      VerboseObject -> 
        Output File : string = none
        Verbosity Level : string = default
    AztecOO -> 
      Output Every RHS : bool = 0
      Adjoint Solve -> 
        Max Iterations : int = 400
        Tolerance : double = 1e-06
        AztecOO Settings -> 
          Aztec Preconditioner : string = ilu
          Aztec Solver : string = GMRES
          Convergence Test : string = r0
          Drop Tolerance : double = 0
          Fill Factor : double = 1
          Graph Fill : int = 0
          Ill-Conditioning Threshold : double = 1e+11
          Orthogonalization : string = Classical
          Output Frequency : int = 0
          Overlap : int = 0
          Polynomial Order : int = 3
          RCM Reordering : string = Disabled
          Size of Krylov Subspace : int = 300
          Steps : int = 3
      Forward Solve -> 
        Max Iterations : int = 400
        Tolerance : double = 1e-06
        AztecOO Settings -> 
          Aztec Preconditioner : string = ilu
          Aztec Solver : string = GMRES
          Convergence Test : string = r0
          Drop Tolerance : double = 0
          Fill Factor : double = 1
          Graph Fill : int = 0
          Ill-Conditioning Threshold : double = 1e+11
          Orthogonalization : string = Classical
          Output Frequency : int = 0
          Overlap : int = 0
          Polynomial Order : int = 3
          RCM Reordering : string = Disabled
          Size of Krylov Subspace : int = 300
          Steps : int = 3
      VerboseObject -> 
        Output File : string = none
        Verbosity Level : string = default
    Belos -> 
      Solver Type : string = Block GMRES
      Solver Types -> 
        Block CG -> 
          Adaptive Block Size : bool = 1
          Block Size : int = 1
          Convergence Tolerance : double = 1e-08
          Maximum Iterations : int = 1000
          Orthogonalization : string = DGKS
          Orthogonalization Constant : double = -1
          Output Frequency : int = -1
          Output Stream : RCP<std::ostream> = RCP<std::ostream>{ptr=0xb85280,node=0xb945a0,count=4}
          Show Maximum Residual Norm Only : bool = 0
          Timer Label : string = Belos
          Verbosity : int = 0
        Block GMRES -> 
          Adaptive Block Size : bool = 1
          Block Size : int = 1
          Convergence Tolerance : double = 1e-08
          Explicit Residual Scaling : string = Norm of Initial Residual
          Flexible Gmres : bool = 0
          Implicit Residual Scaling : string = Norm of Preconditioned Initial Residual
          Maximum Iterations : int = 1000
          Maximum Restarts : int = 20
          Num Blocks : int = 300
          Orthogonalization : string = DGKS
          Orthogonalization Constant : double = -1
          Output Frequency : int = -1
          Output Stream : RCP<std::ostream> = RCP<std::ostream>{ptr=0xb85280,node=0xb94360,count=4}
          Show Maximum Residual Norm Only : bool = 0
          Timer Label : string = Belos
          Verbosity : int = 0
        Pseudo Block GMRES -> 
          Adaptive Block Size : bool = 1
          Block Size : int = 1
          Convergence Tolerance : double = 1e-08
          Deflation Quorum : int = 1
          Explicit Residual Scaling : string = Norm of Initial Residual
          Implicit Residual Scaling : string = Norm of Preconditioned Initial Residual
          Maximum Iterations : int = 1000
          Maximum Restarts : int = 20
          Num Blocks : int = 300
          Orthogonalization : string = DGKS
          Orthogonalization Constant : double = -1
          Output Frequency : int = -1
          Output Stream : RCP<std::ostream> = RCP<std::ostream>{ptr=0xb85280,node=0xb94480,count=4}
          Show Maximum Residual Norm Only : bool = 0
          Timer Label : string = Belos
          Verbosity : int = 0
      VerboseObject -> 
        Output File : string = none
        Verbosity Level : string = default
  Preconditioner Types -> 
    Ifpack -> 
      Overlap : int = 0
      Prec Type : string = ILU
      Ifpack Settings -> 
        amesos: solver type : string = Amesos_Klu
        fact: absolute threshold : double = 0
        fact: drop tolerance : double = 0
        fact: ict level-of-fill : double = 1
        fact: ilut level-of-fill : double = 1
        fact: level-of-fill : int = 0
        fact: relative threshold : double = 1
        fact: relax value : double = 0
        fact: sparskit: alph : double = 0
        fact: sparskit: droptol : double = 0
        fact: sparskit: lfil : int = 0
        fact: sparskit: mbloc : int = -1
        fact: sparskit: permtol : double = 0.1
        fact: sparskit: tol : double = 0
        fact: sparskit: type : string = ILUT
        partitioner: local parts : int = 1
        partitioner: overlap : int = 0
        partitioner: print level : int = 0
        partitioner: type : string = greedy
        partitioner: use symmetric graph : bool = 1
        relaxation: damping factor : double = 1
        relaxation: min diagonal value : double = 1
        relaxation: sweeps : int = 1
        relaxation: type : string = Jacobi
        relaxation: zero starting solution : bool = 1
        schwarz: combine mode : string = Zero
        schwarz: compute condest : bool = 1
        schwarz: filter singletons : bool = 0
        schwarz: reordering type : string = none
      VerboseObject -> 
        Output File : string = none
        Verbosity Level : string = default
    ML -> 
      Base Method Defaults : string = DD
      ML Settings -> 
        aggregation: damping factor : double = 1.333
        aggregation: edge prolongator drop threshold : double = 0
        aggregation: local aggregates : int = 1
        aggregation: next-level aggregates per process : int = 128
        aggregation: nodes per aggregate : int = 512
        aggregation: type : string = Uncoupled-MIS
        coarse: max size : int = 128
        coarse: type : string = Amesos-KLU
        default values : string = maxwell
        eigen-analysis: iterations : int = 10
        eigen-analysis: type : string = cg
        increasing or decreasing : string = decreasing
        max levels : int = 10
        prec type : string = MGV
        smoother: Aztec as solver : bool = 0
        smoother: Aztec options : RCP<vector<int>> = RCP<vector<int>>{ptr=0xbac9a0,node=0xbad610,count=2}
        smoother: Aztec params : RCP<vector<double>> = RCP<vector<double>>{ptr=0xbad830,node=0xbacfe0,count=2}
        smoother: Hiptmair efficient symmetric : bool = 1
        smoother: damping factor : double = 1
        smoother: pre or post : string = both
        smoother: sweeps : int = 1
        smoother: type : string = Hiptmair
        subsmoother: Chebyshev alpha : double = 20
        subsmoother: edge sweeps : int = 4
        subsmoother: node sweeps : int = 4
        subsmoother: type : string = Chebyshev

Creating the LinearOpWithSolveFactoryBase object lowsFactory ...

lowsFactory described as:
  Thyra::BelosLinearOpWithSolveFactory


Running example use cases ...
  
  Running example use cases for a LinearOpWithSolveFactoryBase object ...
    
    Performing a single linear solve ...
    
    Solve status:
      solveStatus = SOLVE_STATUS_CONVERGED
      achievedTol = 1e-08
      message:
        The Belos solver of type "Belos::BlockGmresSolMgr<...,double>{Ortho Type='DGKS', Block Size=1, Num Blocks=300, Max Restarts=20}" returned a solve status of "SOLVE_STATUS_CONVERGED with total CPU time of 0.003241 sec
      extraParameters: NONE
    
    Performing a solve, changing the operator, then performing another solve ...
    
    Performing a solve, changing the operator in a very small way, then performing another solve ...
    
    Performing a solve, changing the operator in a major way, then performing another solve ...

Printing the parameter list (showing what was used) ...
  Matrix File : string = FourByFour.mtx
  Linear Solver Builder -> 
    Enable Delayed Solver Construction : bool = 0   [default]
    Linear Solver Type : string = Belos
    Preconditioner Type : string = Ifpack
    Linear Solver Types -> 
      Belos -> 
        Solver Type : string = Block GMRES   [default]
        Solver Types -> 
          Block CG -> 
            Adaptive Block Size : bool = 1   [unused]
            Block Size : int = 1   [unused]
            Convergence Tolerance : double = 1e-08   [unused]
            Maximum Iterations : int = 1000   [unused]
            Orthogonalization : string = DGKS   [unused]
            Orthogonalization Constant : double = -1   [unused]
            Output Frequency : int = -1   [unused]
            Output Stream : RCP<std::ostream> = RCP<std::ostream>{ptr=0xb85280,node=0xb945a0,count=6}   [unused]
            Show Maximum Residual Norm Only : bool = 0   [unused]
            Timer Label : string = Belos   [unused]
            Verbosity : int = 0   [unused]
          Block GMRES -> 
            Adaptive Block Size : bool = 1   [unused]
            Block Size : int = 1
            Convergence Tolerance : double = 1e-08
            Explicit Residual Scaling : string = Norm of Initial Residual
            Flexible Gmres : bool = 0
            Implicit Residual Scaling : string = Norm of Preconditioned Initial Residual
            Maximum Iterations : int = 1000
            Maximum Restarts : int = 20
            Num Blocks : int = 300
            Orthogonalization : string = DGKS
            Orthogonalization Constant : double = -1
            Output Frequency : int = -1   [unused]
            Output Stream : RCP<std::ostream> = RCP<std::ostream>{ptr=0xb85280,node=0xb94360,count=6}
            Show Maximum Residual Norm Only : bool = 0
            Timer Label : string = Belos
            Verbosity : int = 0
          Pseudo Block GMRES -> 
            Adaptive Block Size : bool = 1   [unused]
            Block Size : int = 1   [unused]
            Convergence Tolerance : double = 1e-08   [unused]
            Deflation Quorum : int = 1   [unused]
            Explicit Residual Scaling : string = Norm of Initial Residual   [unused]
            Implicit Residual Scaling : string = Norm of Preconditioned Initial Residual   [unused]
            Maximum Iterations : int = 1000   [unused]
            Maximum Restarts : int = 20   [unused]
            Num Blocks : int = 300   [unused]
            Orthogonalization : string = DGKS   [unused]
            Orthogonalization Constant : double = -1   [unused]
            Output Frequency : int = -1   [unused]
            Output Stream : RCP<std::ostream> = RCP<std::ostream>{ptr=0xb85280,node=0xb94480,count=6}   [unused]
            Show Maximum Residual Norm Only : bool = 0   [unused]
            Timer Label : string = Belos   [unused]
            Verbosity : int = 0   [unused]
        VerboseObject -> 
          Output File : string = none
          Verbosity Level : string = default
    Preconditioner Types -> 
      Ifpack -> 
        Overlap : int = 0   [default]
        Prec Type : string = ILU   [default]
        Ifpack Settings -> 
          schwarz: combine mode : string = Zero   [default]
          schwarz: compute condest : bool = 1   [default]
          schwarz: filter singletons : bool = 0   [default]
          schwarz: reordering type : string = none   [default]
        VerboseObject -> 
          Output File : string = none   [default]
          Verbosity Level : string = default   [default]
  LinearOpTester -> 
    [empty list]
  LinearOpWithSolveTester -> 
    [empty list]

Congratulations! All of the tests checked out!

Here is the main driver itself:

// @HEADER
// ***********************************************************************
// 
//         Stratimikos: Thyra-based strategies for linear solvers
//                Copyright (2006) Sandia Corporation
// 
// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
// license for use of this work by or on behalf of the U.S. Government.
// 
// This library is free software; you can redistribute it and/or modify
// it under the terms of the GNU Lesser General Public License as
// published by the Free Software Foundation; either version 2.1 of the
// License, or (at your option) any later version.
//  
// This library is distributed in the hope that it will be useful, but
// WITHOUT ANY WARRANTY; without even the implied warranty of
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
// Lesser General Public License for more details.
//  
// You should have received a copy of the GNU Lesser General Public
// License along with this library; if not, write to the Free Software
// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307
// USA
// Questions? Contact Roscoe A. Bartlett (rabartl@sandia.gov) 
// 
// ***********************************************************************
// @HEADER

#include "test_single_stratimikos_solver.hpp"
#include "Thyra_DefaultRealLinearSolverBuilder.hpp"
#include "Thyra_EpetraLinearOp.hpp"
#include "Thyra_LinearOpWithSolveFactoryExamples.hpp"
#include "EpetraExt_readEpetraLinearSystem.h"
#include "Teuchos_ParameterList.hpp"

#ifdef HAVE_MPI
#  include "Epetra_MpiComm.h"
#else
#  include "Epetra_SerialComm.h"
#endif

bool Thyra::test_single_stratimikos_solver(
  Teuchos::ParameterList                  *paramList_inout
  ,const bool                             dumpAll
  ,Teuchos::FancyOStream                  *out
  )
{

  using Teuchos::rcp;
  using Teuchos::RCP;
  using Teuchos::OSTab;
  using Teuchos::ParameterList;
  using Teuchos::getParameter;
  bool success = true;

  try {

    TEST_FOR_EXCEPT(!paramList_inout);

    RCP<ParameterList>
      paramList = rcp(paramList_inout,false);

    if(out) {
      *out << "\nEchoing input parameters ...\n";
      paramList->print(*out,1,true,false);
    }

    // Create list of valid parameter sublists
    Teuchos::ParameterList validParamList("test_single_stratimikos_solver");
    validParamList.set("Matrix File","fileName");
    validParamList.sublist("Linear Solver Builder");
    validParamList.sublist("LinearOpTester");
    validParamList.sublist("LinearOpWithSolveTester");
    
    if(out) *out << "\nValidating top-level input parameters ...\n";
    paramList->validateParameters(validParamList,0);

    const std::string
      &matrixFile = getParameter<std::string>(*paramList,"Matrix File");
    RCP<ParameterList>
      solverBuilderSL  = sublist(paramList,"Linear Solver Builder",true),
      loTesterSL       = sublist(paramList,"LinearOpTester",true),
      lowsTesterSL     = sublist(paramList,"LinearOpWithSolveTester",true);

    if(out) *out << "\nReading in an epetra matrix A from the file \'"<<matrixFile<<"\' ...\n";
  
#ifdef HAVE_MPI
    Epetra_MpiComm comm(MPI_COMM_WORLD);
#else
    Epetra_SerialComm comm;
#endif
    Teuchos::RCP<Epetra_CrsMatrix> epetra_A;
    EpetraExt::readEpetraLinearSystem( matrixFile, comm, &epetra_A );

    Teuchos::RCP<LinearOpBase<double> >
      A = Teuchos::rcp(new EpetraLinearOp(epetra_A));

    if(out) *out << "\nCreating a Thyra::DefaultRealLinearSolverBuilder object ...\n";
    
    RCP<Thyra::LinearSolverBuilderBase<double> >
      linearSolverBuilder = rcp(new Thyra::DefaultRealLinearSolverBuilder);

    if(out) {
      *out << "\nValid parameters for DefaultRealLinearSolverBuilder ...\n";
      linearSolverBuilder->getValidParameters()->print(*out,1,true,false);
    }

    linearSolverBuilder->setParameterList(solverBuilderSL);

    if(out) *out << "\nCreating the LinearOpWithSolveFactoryBase object lowsFactory ...\n";
    RCP<LinearOpWithSolveFactoryBase<double> >
      lowsFactory = createLinearSolveStrategy(*linearSolverBuilder);
    if(out) *out << "\nlowsFactory described as:\n" << describe(*lowsFactory,Teuchos::VERB_MEDIUM) << std::endl;

    if(out) *out << "\nRunning example use cases ...\n";
    
    nonExternallyPreconditionedLinearSolveUseCases(
      *A,*lowsFactory,false,*out
      );
    
    // ToDo: Finish tests!

    if(out) {
      *out << "\nPrinting the parameter list (showing what was used) ...\n";
      paramList->print(*out,1,true,true);
    }
    
  }
  catch( const std::exception &excpt ) {
    std::cerr << "*** Caught standard exception : " << excpt.what() << std::endl;
    success = false;
  }
  
  return success;
  
}

Generated on Wed Jul 22 13:20:34 2009 for Stratimikos by doxygen 1.5.8