|
Stratimikos Package Browser (Single Doxygen Collection)
Version of the Day
|
00001 00002 // @HEADER 00003 // *********************************************************************** 00004 // 00005 // Stratimikos: Thyra-based strategies for linear solvers 00006 // Copyright (2006) Sandia Corporation 00007 // 00008 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive 00009 // license for use of this work by or on behalf of the U.S. Government. 00010 // 00011 // Redistribution and use in source and binary forms, with or without 00012 // modification, are permitted provided that the following conditions are 00013 // met: 00014 // 00015 // 1. Redistributions of source code must retain the above copyright 00016 // notice, this list of conditions and the following disclaimer. 00017 // 00018 // 2. Redistributions in binary form must reproduce the above copyright 00019 // notice, this list of conditions and the following disclaimer in the 00020 // documentation and/or other materials provided with the distribution. 00021 // 00022 // 3. Neither the name of the Corporation nor the names of the 00023 // contributors may be used to endorse or promote products derived from 00024 // this software without specific prior written permission. 00025 // 00026 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY 00027 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE 00028 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR 00029 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE 00030 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, 00031 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, 00032 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR 00033 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF 00034 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING 00035 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS 00036 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. 00037 // 00038 // Questions? Contact Roscoe A. Bartlett (rabartl@sandia.gov) 00039 // 00040 // *********************************************************************** 00041 // @HEADER 00042 // 00043 // This driver reads a problem from a Harwell-Boeing (HB) file into an 00044 // Epetra_CrsMatrix. This matrix is then converted into a Thyra linear operator 00045 // through the Thyra-Epetra adapters. 00046 // The right-hand-side from the HB file is used instead of random vectors. 00047 // The initial guesses are all set to zero. 00048 // 00049 // 00050 00051 // Thyra includes 00052 #include "Thyra_BelosLinearOpWithSolveFactory.hpp" 00053 #include "Thyra_LinearOpWithSolveFactoryHelpers.hpp" 00054 #include "Thyra_EpetraLinearOp.hpp" 00055 00056 // Epetra includes 00057 #include "Epetra_SerialComm.h" 00058 #include "EpetraExt_readEpetraLinearSystem.h" 00059 00060 // Ifpack includes 00061 #ifdef HAVE_BELOS_IFPACK 00062 #include "Thyra_IfpackPreconditionerFactory.hpp" 00063 #endif 00064 00065 // Teuchos includes 00066 #include "Teuchos_ParameterList.hpp" 00067 #include "Teuchos_VerboseObject.hpp" 00068 00069 int main(int argc, char* argv[]) 00070 { 00071 // 00072 // Get a default output stream from the Teuchos::VerboseObjectBase 00073 // 00074 Teuchos::RCP<Teuchos::FancyOStream> 00075 out = Teuchos::VerboseObjectBase::getDefaultOStream(); 00076 // 00077 // Set the parameters for the Belos LOWS Factory and create a parameter list. 00078 // 00079 int blockSize = 2; 00080 int maxIterations = 400; 00081 int maxRestarts = 25; 00082 int gmresKrylovLength = 25; 00083 int outputFrequency = 1; 00084 bool outputMaxResOnly = true; 00085 double maxResid = 1e-6; 00086 00087 Teuchos::RCP<Teuchos::ParameterList> 00088 belosLOWSFPL = Teuchos::rcp( new Teuchos::ParameterList() ); 00089 00090 belosLOWSFPL->set("Solver Type","Block GMRES"); 00091 00092 Teuchos::ParameterList& belosLOWSFPL_solver = 00093 belosLOWSFPL->sublist("Solver Types"); 00094 00095 Teuchos::ParameterList& belosLOWSFPL_gmres = 00096 belosLOWSFPL_solver.sublist("Block GMRES"); 00097 00098 belosLOWSFPL_gmres.set("Maximum Iterations",int(maxIterations)); 00099 belosLOWSFPL_gmres.set("Convergence Tolerance",double(maxResid)); 00100 belosLOWSFPL_gmres.set("Maximum Restarts",int(maxRestarts)); 00101 belosLOWSFPL_gmres.set("Block Size",int(blockSize)); 00102 belosLOWSFPL_gmres.set("Num Blocks",int(gmresKrylovLength)); 00103 belosLOWSFPL_gmres.set("Output Frequency",int(outputFrequency)); 00104 belosLOWSFPL_gmres.set("Show Maximum Residual Norm Only",bool(outputMaxResOnly)); 00105 00106 #ifdef HAVE_BELOS_IFPACK 00107 // 00108 // Set the parameters for the Ifpack Preconditioner Factory and create parameter list 00109 // 00110 Teuchos::ParameterList &ifpackPFSL = belosLOWSFPL->sublist("IfpackPreconditionerFactory"); 00111 00112 ifpackPFSL.set("Overlap",int(2)); 00113 ifpackPFSL.set("Prec Type","ILUT"); 00114 #endif 00115 00116 // Whether the linear solver succeeded. 00117 // (this will be set during the residual check at the end) 00118 bool success = true; 00119 00120 // Number of random right-hand sides we will be solving for. 00121 int numRhs = 5; 00122 00123 // Name of input matrix file 00124 std::string matrixFile = "orsirr1.hb"; 00125 00126 // Read in the matrix file (can be *.mtx, *.hb, etc.) 00127 Epetra_SerialComm comm; 00128 Teuchos::RCP<Epetra_CrsMatrix> epetra_A; 00129 EpetraExt::readEpetraLinearSystem( matrixFile, comm, &epetra_A ); 00130 00131 // Create a Thyra linear operator (A) using the Epetra_CrsMatrix (epetra_A). 00132 Teuchos::RCP<const Thyra::LinearOpBase<double> > 00133 A = Thyra::epetraLinearOp(epetra_A); 00134 00135 // Get the domain space for the Thyra linear operator 00136 Teuchos::RCP<const Thyra::VectorSpaceBase<double> > domain = A->domain(); 00137 00138 // Create the Belos LOWS factory. 00139 Teuchos::RCP<Thyra::LinearOpWithSolveFactoryBase<double> > 00140 belosLOWSFactory = Teuchos::rcp(new Thyra::BelosLinearOpWithSolveFactory<double>()); 00141 00142 #ifdef HAVE_BELOS_IFPACK 00143 00144 // Set the preconditioner factory for the LOWS factory. 00145 belosLOWSFactory->setPreconditionerFactory( 00146 Teuchos::rcp(new Thyra::IfpackPreconditionerFactory()) 00147 ,"IfpackPreconditionerFactory" 00148 ); 00149 #endif 00150 00151 // Set the parameter list to specify the behavior of the factory. 00152 belosLOWSFactory->setParameterList( belosLOWSFPL ); 00153 00154 // Set the output stream and the verbosity level (prints to std::cout by defualt) 00155 belosLOWSFactory->setVerbLevel(Teuchos::VERB_LOW); 00156 00157 // Create a BelosLinearOpWithSolve object from the Belos LOWS factory. 00158 Teuchos::RCP<Thyra::LinearOpWithSolveBase<double> > 00159 nsA = belosLOWSFactory->createOp(); 00160 00161 // Initialize the BelosLinearOpWithSolve object with the Thyra linear operator. 00162 Thyra::initializeOp<double>( *belosLOWSFactory, A, &*nsA ); 00163 00164 // Create a right-hand side with numRhs vectors in it and randomize it. 00165 Teuchos::RCP< Thyra::MultiVectorBase<double> > 00166 b = Thyra::createMembers(domain, numRhs); 00167 Thyra::seed_randomize<double>(0); 00168 Thyra::randomize(-1.0, 1.0, &*b); 00169 00170 // Create an initial std::vector with numRhs vectors in it and initialize it to zero. 00171 Teuchos::RCP< Thyra::MultiVectorBase<double> > 00172 x = Thyra::createMembers(domain, numRhs); 00173 Thyra::assign(&*x, 0.0); 00174 00175 // Perform solve using the linear operator to get the approximate solution of Ax=b, 00176 // where b is the right-hand side and x is the left-hand side. 00177 Thyra::SolveStatus<double> solveStatus; 00178 solveStatus = Thyra::solve( *nsA, Thyra::NONCONJ_ELE, *b, &*x ); 00179 00180 // Print out status of solve. 00181 *out << "\nBelos LOWS Status: "<< solveStatus << std::endl; 00182 00183 // 00184 // Compute residual and double check convergence. 00185 // 00186 std::vector<double> norm_b(numRhs), norm_res(numRhs); 00187 Teuchos::RCP< Thyra::MultiVectorBase<double> > 00188 y = Thyra::createMembers(domain, numRhs); 00189 00190 // Compute the column norms of the right-hand side b. 00191 Thyra::norms_2( *b, &norm_b[0] ); 00192 00193 // Compute y=A*x, where x is the solution from the linear solver. 00194 A->apply( Thyra::NONCONJ_ELE, *x, &*y ); 00195 00196 // Compute A*x-b = y-b 00197 Thyra::update( -1.0, *b, &*y ); 00198 00199 // Compute the column norms of A*x-b. 00200 Thyra::norms_2( *y, &norm_res[0] ); 00201 00202 // Print out the final relative residual norms. 00203 double rel_res = 0.0; 00204 *out << "Final relative residual norms" << std::endl; 00205 for (int i=0; i<numRhs; ++i) { 00206 rel_res = norm_res[i]/norm_b[i]; 00207 if (rel_res > maxResid) 00208 success = false; 00209 *out << "RHS " << i+1 << " : " 00210 << std::setw(16) << std::right << rel_res << std::endl; 00211 } 00212 00213 return ( success ? 0 : 1 ); 00214 }
1.7.6.1