00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043 #include "BelosConfigDefs.hpp"
00044 #include "BelosLinearProblem.hpp"
00045 #include "BelosOutputManager.hpp"
00046 #include "BelosStatusTestMaxIters.hpp"
00047 #include "BelosStatusTestResNorm.hpp"
00048 #include "BelosStatusTestOutputter.hpp"
00049 #include "BelosStatusTestCombo.hpp"
00050 #include "BelosEpetraAdapter.hpp"
00051 #include "BelosBlockCG.hpp"
00052 #include "createEpetraProblem.hpp"
00053 #include "Trilinos_Util.h"
00054 #include "Ifpack_CrsIct.h"
00055 #include "Epetra_CrsMatrix.h"
00056 #include "Epetra_Map.h"
00057 #include "Teuchos_CommandLineProcessor.hpp"
00058
00059 int main(int argc, char *argv[]) {
00060
00061 #ifdef EPETRA_MPI
00062
00063 MPI_Init(&argc,&argv);
00064 Belos::MPIFinalize mpiFinalize;
00065 #endif
00066
00067 using Teuchos::RefCountPtr;
00068 using Teuchos::rcp;
00069
00070
00071
00072 bool verbose = false, proc_verbose = false;
00073 int frequency = -1;
00074 int numrhs = 15;
00075 int blockSize = 10;
00076 std::string filename("bcsstk14.hb");
00077 double tol = 1.0e-5;
00078
00079 Teuchos::CommandLineProcessor cmdp(false,true);
00080 cmdp.setOption("verbose","quiet",&verbose,"Print messages and results.");
00081 cmdp.setOption("frequency",&frequency,"Solvers frequency for printing residuals (#iters).");
00082 cmdp.setOption("filename",&filename,"Filename for Harwell-Boeing test matrix.");
00083 cmdp.setOption("tol",&tol,"Relative residual tolerance used by CG solver.");
00084 cmdp.setOption("num-rhs",&numrhs,"Number of right-hand sides to be solved for.");
00085 cmdp.setOption("block-size",&blockSize,"Block size to be used by CG solver.");
00086 if (cmdp.parse(argc,argv) != Teuchos::CommandLineProcessor::PARSE_SUCCESSFUL) {
00087 return -1;
00088 }
00089 if (!verbose)
00090 frequency = -1;
00091
00092
00093
00094 int MyPID;
00095 RefCountPtr<Epetra_CrsMatrix> A;
00096 int return_val =Belos::createEpetraProblem(filename,NULL,&A,NULL,NULL,&MyPID);
00097 if(return_val != 0) return return_val;
00098 proc_verbose = ( verbose && (MyPID==0) );
00099
00100
00101
00102 if (proc_verbose) cout << endl << endl;
00103 if (proc_verbose) cout << "Constructing ICT preconditioner" << endl;
00104 int Lfill = 0;
00105
00106 if (proc_verbose) cout << "Using Lfill = " << Lfill << endl;
00107 int Overlap = 0;
00108
00109 if (proc_verbose) cout << "Using Level Overlap = " << Overlap << endl;
00110 double Athresh = 0.0;
00111
00112 if (proc_verbose) cout << "Using Absolute Threshold Value of " << Athresh << endl;
00113 double Rthresh = 1.0;
00114
00115 if (proc_verbose) cout << "Using Relative Threshold Value of " << Rthresh << endl;
00116 double dropTol = 1.0e-6;
00117
00118 Teuchos::RefCountPtr<Ifpack_CrsIct> ICT;
00119
00120 if (Lfill > -1) {
00121 ICT = Teuchos::rcp( new Ifpack_CrsIct(*A, dropTol, Lfill) );
00122 ICT->SetAbsoluteThreshold(Athresh);
00123 ICT->SetRelativeThreshold(Rthresh);
00124 int initerr = ICT->InitValues(*A);
00125 if (initerr != 0) cout << "InitValues error = " << initerr;
00126 assert(ICT->Factor() == 0);
00127 }
00128
00129 bool transA = false;
00130 double Cond_Est;
00131 ICT->Condest(transA, Cond_Est);
00132 if (proc_verbose) {
00133 cout << "Condition number estimate for this preconditoner = " << Cond_Est << endl;
00134 cout << endl;
00135 }
00136
00137
00138
00139 typedef double ST;
00140 typedef Belos::Operator<ST> OP;
00141 typedef Belos::MultiVec<ST> MV;
00142 typedef Belos::OperatorTraits<ST,MV,OP> OPT;
00143 typedef Belos::MultiVecTraits<ST,MV> MVT;
00144
00145
00146
00147 Belos::EpetraOp Amat( A );
00148
00149
00150
00151 Belos::EpetraPrecOp Prec( ICT );
00152
00153
00154
00155
00156 const Epetra_Map &Map = A->RowMap();
00157 const int NumGlobalElements = Map.NumGlobalElements();
00158 int maxits = NumGlobalElements/blockSize - 1;
00159
00160
00161
00162 Belos::EpetraMultiVec soln(Map, numrhs);
00163 Belos::EpetraMultiVec rhs(Map, numrhs);
00164 rhs.MvRandom();
00165
00166
00167
00168 Belos::LinearProblem<ST,MV,OP> My_LP( rcp(&Amat, false), rcp(&soln, false), rcp(&rhs,false) );
00169 My_LP.SetLeftPrec( rcp(&Prec,false) );
00170 My_LP.SetBlockSize( blockSize );
00171
00172
00173
00174 Belos::OutputManager<ST> My_OM( MyPID );
00175 if (verbose)
00176 My_OM.SetVerbosity( Belos::Errors + Belos::Warnings
00177 + Belos::TimingDetails + Belos::FinalSummary );
00178
00179 Belos::StatusTestMaxIters<ST,MV,OP> test1( maxits );
00180 Belos::StatusTestResNorm<ST,MV,OP> test2( tol );
00181 Belos::StatusTestOutputter<ST,MV,OP> test3( frequency, false );
00182 test3.set_resNormStatusTest( rcp(&test2,false) );
00183 test3.set_outputManager( rcp(&My_OM,false) );
00184 Belos::StatusTestCombo<ST,MV,OP> My_Test( Belos::StatusTestCombo<ST,MV,OP>::OR, test1, test3 );
00185
00186
00187
00188
00189
00190
00191 Belos::BlockCG<ST,MV,OP> MyBlockCG( rcp(&My_LP, false), rcp(&My_Test,false), rcp(&My_OM,false));
00192
00193
00194
00195 if (proc_verbose) {
00196 cout << endl << endl;
00197 cout << "Dimension of matrix: " << NumGlobalElements << endl;
00198 cout << "Number of right-hand sides: " << numrhs << endl;
00199 cout << "Block size used by solver: " << blockSize << endl;
00200 cout << "Max number of CG iterations: " << maxits << endl;
00201 cout << "Relative residual tolerance: " << tol << endl;
00202 cout << endl;
00203 }
00204
00205 if (proc_verbose) {
00206 cout << endl << endl;
00207 cout << "Running Block CG -- please wait" << endl;
00208 cout << (numrhs+blockSize-1)/blockSize
00209 << " pass(es) through the solver required to solve for " << endl;
00210 cout << numrhs << " right-hand side(s) -- using a block size of " << blockSize
00211 << endl << endl;
00212 }
00213
00214 MyBlockCG.Solve();
00215
00216
00217
00218
00219 std::vector<ST> actual_resids( numrhs );
00220 std::vector<ST> rhs_norm( numrhs );
00221 Belos::EpetraMultiVec resid( Map, numrhs );
00222 OPT::Apply( Amat, soln, resid );
00223 MVT::MvAddMv( -1.0, resid, 1.0, rhs, resid );
00224 MVT::MvNorm( resid, &actual_resids );
00225 MVT::MvNorm( rhs, &rhs_norm );
00226 if (proc_verbose) {
00227 cout<< "---------- Actual Residuals (normalized) ----------"<<endl<<endl;
00228 for (int i=0; i<numrhs; i++) {
00229 cout<<"Problem "<<i<<" : \t"<< actual_resids[i]/rhs_norm[i] <<endl;
00230 }
00231 }
00232
00233 if (My_Test.GetStatus()!=Belos::Converged) {
00234 if (proc_verbose)
00235 cout << "End Result: TEST FAILED" << endl;
00236 return -1;
00237 }
00238
00239
00240
00241 if (proc_verbose)
00242 cout << "End Result: TEST PASSED" << endl;
00243 return 0;
00244
00245 }
00246