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 #include "Thyra_BelosLinearOpWithSolveFactory.hpp"
00041 #include "Thyra_LinearOpWithSolveFactoryHelpers.hpp"
00042 #include "Thyra_TpetraLinearOp.hpp"
00043
00044
00045 #include "Tpetra_ElementSpace.hpp"
00046 #include "Tpetra_VectorSpace.hpp"
00047 #include "Tpetra_CisMatrix.hpp"
00048
00049
00050 #include "Teuchos_ParameterList.hpp"
00051 #include "Teuchos_VerboseObject.hpp"
00052 #include "Teuchos_GlobalMPISession.hpp"
00053
00054 #ifdef HAVE_MPI
00055 # include "Tpetra_MpiPlatform.hpp"
00056 #else
00057 # include "Tpetra_SerialPlatform.hpp"
00058 #endif
00059
00060 int main(int argc, char* argv[])
00061 {
00062
00063
00064
00065 Teuchos::RefCountPtr<Teuchos::FancyOStream>
00066 out = Teuchos::VerboseObjectBase::getDefaultOStream();
00067
00068 Teuchos::GlobalMPISession mpiSession(&argc,&argv);
00069
00070
00071
00072
00073
00074
00075
00076 #ifdef HAVE_COMPLEX
00077 typedef std::complex<double> ST;
00078 #elif HAVE_COMPLEX_H
00079 typedef ::complex<double> ST;
00080 #else
00081 typedef double ST;
00082 #endif
00083 typedef Teuchos::ScalarTraits<ST>::magnitudeType MT;
00084 typedef int OT;
00085
00086
00087
00088
00089
00090
00091
00092 #ifdef HAVE_MPI
00093 MPI_Comm mpiComm = MPI_COMM_WORLD;
00094 const Tpetra::MpiPlatform<OT,OT> ordinalPlatform(mpiComm);
00095 const Tpetra::MpiPlatform<OT,ST> scalarPlatform(mpiComm);
00096 #else
00097 const Tpetra::SerialPlatform<OT,OT> ordinalPlatform;
00098 const Tpetra::SerialPlatform<OT,ST> scalarPlatform;
00099 #endif
00100
00101
00102 OT globalDim = 500;
00103
00104
00105 const Tpetra::ElementSpace<OT> elementSpace(globalDim,0,ordinalPlatform);
00106 const Tpetra::VectorSpace<OT,ST> vectorSpace(elementSpace,scalarPlatform);
00107
00108
00109 RefCountPtr<Tpetra::CisMatrix<OT,ST> >
00110 tpetra_A = rcp(new Tpetra::CisMatrix<OT,ST>(vectorSpace));
00111
00112
00113 const int numMyElements = vectorSpace.getNumMyEntries();
00114 const std::vector<int> &myGlobalElements = vectorSpace.elementSpace().getMyGlobalElements();
00115
00116
00117 const ST offDiag = -1.0, diag = 2.0;
00118 int numEntries; ST values[3]; int indexes[3];
00119 for( int k = 0; k < numMyElements; ++k ) {
00120 const int rowIndex = myGlobalElements[k];
00121 if( rowIndex == 0 ) {
00122 numEntries = 2;
00123 values[0] = diag; values[1] = offDiag;
00124 indexes[0] = 0; indexes[1] = 1;
00125 }
00126 else if( rowIndex == globalDim - 1 ) {
00127 numEntries = 2;
00128 values[0] = offDiag; values[1] = diag;
00129 indexes[0] = globalDim-2; indexes[1] = globalDim-1;
00130 }
00131 else {
00132 numEntries = 3;
00133 values[0] = offDiag; values[1] = diag; values[2] = offDiag;
00134 indexes[0] = rowIndex-1; indexes[1] = rowIndex; indexes[2] = rowIndex+1;
00135 }
00136 tpetra_A->submitEntries(Tpetra::Insert,rowIndex,numEntries,values,indexes);
00137 }
00138
00139
00140 tpetra_A->fillComplete();
00141
00142
00143 RefCountPtr<Thyra::LinearOpBase<ST> >
00144 A = Teuchos::rcp( new Thyra::TpetraLinearOp<OT,ST>(
00145 Teuchos::rcp_implicit_cast<Tpetra::Operator<OT,ST> >(tpetra_A) ) );
00146
00147
00148
00149
00150 int blockSize = 1;
00151 int maxIterations = globalDim;
00152 int maxRestarts = 0;
00153 int gmresKrylovLength = globalDim;
00154 int outputFrequency = 100;
00155 bool outputMaxResOnly = true;
00156 MT maxResid = 1e-3;
00157
00158 ST one = Teuchos::ScalarTraits<ST>::one();
00159 ST zero = Teuchos::ScalarTraits<ST>::zero();
00160
00161 Teuchos::RefCountPtr<Teuchos::ParameterList>
00162 belosLOWSFPL = Teuchos::rcp( new Teuchos::ParameterList() );
00163
00164 belosLOWSFPL->set("Solver Type","GMRES");
00165 belosLOWSFPL->set("Max Iters",int(maxIterations));
00166 belosLOWSFPL->set("Default Rel Res Norm",double(maxResid));
00167 belosLOWSFPL->set("Max Restarts",int(maxRestarts));
00168 belosLOWSFPL->set("Block Size",int(blockSize));
00169 belosLOWSFPL->sublist("GMRES").set("Max Number of Krylov Vectors",int(gmresKrylovLength*blockSize));
00170 belosLOWSFPL->sublist("GMRES").set("Variant","Standard");
00171 Teuchos::ParameterList &outputterSL = belosLOWSFPL->sublist("Outputter");
00172 outputterSL.set("Output Frequency",int(outputFrequency));
00173 outputterSL.set("Output Max Res Only",bool(outputMaxResOnly));
00174
00175
00176
00177 bool success = true;
00178
00179
00180
00181
00182
00183
00184
00185
00186 int numRhs = 5;
00187
00188
00189 Teuchos::RefCountPtr<Thyra::MultiVectorBase<ST> > x, b;
00190
00191 if (numRhs==1) {
00192
00193
00194
00195
00196
00197 Teuchos::RefCountPtr<Tpetra::Vector<OT,ST> > tpetra_b =
00198 Teuchos::rcp( new Tpetra::Vector<OT,ST>(vectorSpace) );
00199
00200
00201 tpetra_b->setAllToRandom();
00202
00203
00204 b = Thyra::create_Vector(tpetra_b);
00205
00206
00207 Teuchos::RefCountPtr<Tpetra::Vector<OT,ST> > tpetra_x =
00208 Teuchos::rcp( new Tpetra::Vector<OT,ST>(vectorSpace) );
00209
00210
00211 tpetra_x->setAllToScalar( zero );
00212
00213
00214 x = Thyra::create_Vector(tpetra_x);
00215
00216 }
00217 else {
00218
00219
00220
00221
00222
00223
00224
00225 Teuchos::RefCountPtr<const Thyra::VectorSpaceBase<ST> > domain = A->domain();
00226
00227
00228 x = Thyra::createMembers(domain, numRhs);
00229 b = Thyra::createMembers(domain, numRhs);
00230
00231
00232
00233
00234 for ( int j=0; j<numRhs; ++j ) {
00235
00236
00237
00238 Teuchos::RefCountPtr<Tpetra::Vector<OT,ST> >
00239 tpetra_b_j = Thyra::get_Tpetra_Vector(vectorSpace,b->col(j));
00240 tpetra_b_j->setAllToRandom();
00241
00242
00243
00244 Teuchos::RefCountPtr<Tpetra::Vector<OT,ST> >
00245 tpetra_x_j = Thyra::get_Tpetra_Vector(vectorSpace,x->col(j));
00246 tpetra_x_j->setAllToScalar( zero );
00247
00248
00249
00250
00251 for ( int i=0; i<numMyElements; ++i ) {
00252
00253
00254
00255 const int rowIndex = myGlobalElements[ i ];
00256
00257
00258 (*tpetra_x_j)[ rowIndex ] = zero;
00259 }
00260 }
00261 }
00262
00263
00264
00265
00266
00267
00268
00269
00270
00271
00272
00273 Teuchos::RefCountPtr<Thyra::LinearOpWithSolveFactoryBase<ST> >
00274 belosLOWSFactory = Teuchos::rcp(new Thyra::BelosLinearOpWithSolveFactory<ST>());
00275
00276
00277
00278 belosLOWSFactory->setVerbLevel(Teuchos::VERB_LOW);
00279
00280
00281 belosLOWSFactory->setParameterList( belosLOWSFPL );
00282
00283
00284 Teuchos::RefCountPtr<Thyra::LinearOpWithSolveBase<ST> >
00285 nsA = belosLOWSFactory->createOp();
00286
00287
00288 Thyra::initializeOp<ST>( *belosLOWSFactory, A, &*nsA );
00289
00290
00291
00292 Thyra::SolveStatus<ST> solveStatus;
00293 solveStatus = Thyra::solve( *nsA, Thyra::NONCONJ_ELE, *b, &*x );
00294
00295
00296 *out << "\nBelos LOWS Status: "<< solveStatus << endl;
00297
00298
00299
00300
00301
00302
00303
00304
00305 std::vector<ST> norm_b(numRhs), norm_res(numRhs);
00306
00307 for ( int j=0; j<numRhs; ++j ) {
00308
00309
00310
00311 Teuchos::RefCountPtr<Tpetra::Vector<OT,ST> >
00312 tpetra_x_j = Thyra::get_Tpetra_Vector(vectorSpace,x->col(j));
00313
00314 Teuchos::RefCountPtr<Tpetra::Vector<OT,ST> >
00315 tpetra_b_j = Thyra::get_Tpetra_Vector(vectorSpace,b->col(j));
00316
00317
00318 norm_b[j] = tpetra_b_j->norm2();
00319
00320
00321 Tpetra::Vector<OT,ST> y(vectorSpace);
00322 tpetra_A->apply( *tpetra_x_j, y );
00323
00324
00325 y.update( one, *tpetra_b_j, -one );
00326
00327
00328 norm_res[j] = y.norm2();
00329 }
00330
00331
00332 MT rel_res = 0.0;
00333 *out << "Final relative residual norms" << endl;
00334 for (int i=0; i<numRhs; ++i) {
00335 rel_res = Teuchos::ScalarTraits<ST>::real(norm_res[i]/norm_b[i]);
00336 if (rel_res > maxResid)
00337 success = false;
00338 *out << "RHS " << i+1 << " : "
00339 << std::setw(16) << std::right << rel_res << endl;
00340 }
00341
00342 return ( success ? 0 : 1 );
00343 }