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_PreconditionerFactoryHelpers.hpp"
00043 #include "Thyra_TpetraLinearOp.hpp"
00044 #include "Thyra_DefaultPreconditioner.hpp"
00045
00046
00047 #include "Tpetra_ElementSpace.hpp"
00048 #include "Tpetra_VectorSpace.hpp"
00049 #include "Tpetra_CisMatrix.hpp"
00050
00051
00052 #include "Teuchos_ParameterList.hpp"
00053 #include "Teuchos_VerboseObject.hpp"
00054 #include "Teuchos_GlobalMPISession.hpp"
00055
00056 #ifdef HAVE_MPI
00057 # include "Tpetra_MpiPlatform.hpp"
00058 #else
00059 # include "Tpetra_SerialPlatform.hpp"
00060 #endif
00061
00062 int main(int argc, char* argv[])
00063 {
00064
00065
00066
00067 Teuchos::RefCountPtr<Teuchos::FancyOStream>
00068 out = Teuchos::VerboseObjectBase::getDefaultOStream();
00069
00070 Teuchos::GlobalMPISession mpiSession(&argc,&argv);
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 ST one = Teuchos::ScalarTraits<ST>::one();
00087 ST zero = Teuchos::ScalarTraits<ST>::zero();
00088
00089
00090
00091
00092
00093 #ifdef HAVE_MPI
00094 MPI_Comm mpiComm = MPI_COMM_WORLD;
00095 const Tpetra::MpiPlatform<OT,OT> ordinalPlatform(mpiComm);
00096 const Tpetra::MpiPlatform<OT,ST> scalarPlatform(mpiComm);
00097 #else
00098 const Tpetra::SerialPlatform<OT,OT> ordinalPlatform;
00099 const Tpetra::SerialPlatform<OT,ST> scalarPlatform;
00100 #endif
00101
00102
00103 OT globalDim = 500;
00104
00105
00106 const Tpetra::ElementSpace<OT> elementSpace(globalDim,0,ordinalPlatform);
00107 const Tpetra::VectorSpace<OT,ST> vectorSpace(elementSpace,scalarPlatform);
00108
00109
00110 RefCountPtr<Tpetra::CisMatrix<OT,ST> >
00111 tpetra_A = rcp(new Tpetra::CisMatrix<OT,ST>(vectorSpace));
00112
00113
00114 RefCountPtr<Tpetra::CisMatrix<OT,ST> >
00115 tpetra_Prec = rcp(new Tpetra::CisMatrix<OT,ST>(vectorSpace));
00116
00117
00118 const int numMyElements = vectorSpace.getNumMyEntries();
00119 const std::vector<int> &myGlobalElements = vectorSpace.elementSpace().getMyGlobalElements();
00120
00121
00122 const ST offDiag = -1.0, diag = 2.0;
00123 int numEntries; ST values[3]; int indexes[3];
00124 ST prec_value[1]; int prec_index[1];
00125 prec_value[0] = one/diag;
00126 for( int k = 0; k < numMyElements; ++k ) {
00127 const int rowIndex = myGlobalElements[k];
00128 prec_index[0] = rowIndex;
00129 if( rowIndex == 0 ) {
00130 numEntries = 2;
00131 values[0] = diag; values[1] = offDiag;
00132 indexes[0] = 0; indexes[1] = 1;
00133 }
00134 else if( rowIndex == globalDim - 1 ) {
00135 numEntries = 2;
00136 values[0] = offDiag; values[1] = diag;
00137 indexes[0] = globalDim-2; indexes[1] = globalDim-1;
00138 }
00139 else {
00140 numEntries = 3;
00141 values[0] = offDiag; values[1] = diag; values[2] = offDiag;
00142 indexes[0] = rowIndex-1; indexes[1] = rowIndex; indexes[2] = rowIndex+1;
00143 }
00144 tpetra_A->submitEntries(Tpetra::Insert,rowIndex,numEntries,values,indexes);
00145 tpetra_Prec->submitEntries(Tpetra::Insert,rowIndex,1,prec_value,prec_index);
00146 }
00147
00148
00149 tpetra_A->fillComplete();
00150 tpetra_Prec->fillComplete();
00151
00152
00153 RefCountPtr<Thyra::LinearOpBase<ST> >
00154 A = Teuchos::rcp( new Thyra::TpetraLinearOp<OT,ST>(
00155 Teuchos::rcp_implicit_cast<Tpetra::Operator<OT,ST> >(tpetra_A) ) );
00156
00157
00158 RefCountPtr<Thyra::LinearOpBase<ST> >
00159 Prec = Teuchos::rcp( new Thyra::TpetraLinearOp<OT,ST>(
00160 Teuchos::rcp_implicit_cast<Tpetra::Operator<OT,ST> >(tpetra_Prec) ) );
00161
00162
00163 RefCountPtr<const Thyra::DefaultPreconditioner<ST> >
00164 DefPrec = Thyra::unspecifiedPrec<ST>(Prec);
00165
00166
00167
00168
00169
00170 int blockSize = 1;
00171 int maxIterations = globalDim;
00172 int maxRestarts = 0;
00173 int gmresKrylovLength = globalDim;
00174 int outputFrequency = 100;
00175 bool outputMaxResOnly = true;
00176 MT maxResid = 1e-3;
00177
00178 Teuchos::RefCountPtr<Teuchos::ParameterList>
00179 belosLOWSFPL = Teuchos::rcp( new Teuchos::ParameterList() );
00180
00181 belosLOWSFPL->set("Solver Type","GMRES");
00182 belosLOWSFPL->set("Max Iters",int(maxIterations));
00183 belosLOWSFPL->set("Default Rel Res Norm",double(maxResid));
00184 belosLOWSFPL->set("Max Restarts",int(maxRestarts));
00185 belosLOWSFPL->set("Block Size",int(blockSize));
00186 belosLOWSFPL->sublist("GMRES").set("Max Number of Krylov Vectors",int(gmresKrylovLength*blockSize));
00187 belosLOWSFPL->sublist("GMRES").set("Variant","Standard");
00188 Teuchos::ParameterList &outputterSL = belosLOWSFPL->sublist("Outputter");
00189 outputterSL.set("Output Frequency",int(outputFrequency));
00190 outputterSL.set("Output Max Res Only",bool(outputMaxResOnly));
00191
00192
00193
00194 bool success = true;
00195
00196
00197 int numRhs = 5;
00198
00199
00200 Teuchos::RefCountPtr<const Thyra::VectorSpaceBase<ST> > domain = A->domain();
00201
00202
00203 Teuchos::RefCountPtr<Thyra::LinearOpWithSolveFactoryBase<ST> >
00204 belosLOWSFactory = Teuchos::rcp(new Thyra::BelosLinearOpWithSolveFactory<ST>());
00205
00206
00207 belosLOWSFactory->setParameterList( belosLOWSFPL );
00208
00209
00210
00211 belosLOWSFactory->setVerbLevel(Teuchos::VERB_LOW);
00212
00213
00214 Teuchos::RefCountPtr<Thyra::LinearOpWithSolveBase<ST> >
00215 nsA = belosLOWSFactory->createOp();
00216
00217
00218
00219 Thyra::initializeOp<ST>( *belosLOWSFactory, A, &*nsA );
00220
00221
00222 Teuchos::RefCountPtr< Thyra::MultiVectorBase<ST> >
00223 b = Thyra::createMembers(domain, numRhs);
00224 Thyra::seed_randomize<ST>(0);
00225 Thyra::randomize(-one, one, &*b);
00226
00227
00228 Teuchos::RefCountPtr< Thyra::MultiVectorBase<ST> >
00229 x = Thyra::createMembers(domain, numRhs);
00230 Thyra::assign(&*x, zero);
00231
00232
00233
00234 Thyra::SolveStatus<ST> solveStatus;
00235 solveStatus = Thyra::solve( *nsA, Thyra::NONCONJ_ELE, *b, &*x );
00236
00237
00238 *out << "\nBelos LOWS Status: "<< solveStatus << endl;
00239
00240
00241
00242
00243 std::vector<MT> norm_b(numRhs), norm_res(numRhs);
00244 Teuchos::RefCountPtr< Thyra::MultiVectorBase<ST> >
00245 y = Thyra::createMembers(domain, numRhs);
00246
00247
00248 Thyra::norms_2( *b, &norm_b[0] );
00249
00250
00251 A->apply( Thyra::NONCONJ_ELE, *x, &*y );
00252
00253
00254 Thyra::update( -one, *b, &*y );
00255
00256
00257 Thyra::norms_2( *y, &norm_res[0] );
00258
00259
00260 MT rel_res = 0.0;
00261 *out << "Final relative residual norms" << endl;
00262 for (int i=0; i<numRhs; ++i) {
00263 rel_res = norm_res[i]/norm_b[i];
00264 if (rel_res > maxResid)
00265 success = false;
00266 *out << "RHS " << i+1 << " : "
00267 << std::setw(16) << std::right << rel_res << endl;
00268 }
00269
00270 return ( success ? 0 : 1 );
00271 }