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 #include "AnasaziBlockKrylovSchurSolMgr.hpp"
00035 #include "AnasaziBasicEigenproblem.hpp"
00036 #include "AnasaziConfigDefs.hpp"
00037 #include "AnasaziEpetraAdapter.hpp"
00038 #include "AnasaziBasicSort.hpp"
00039 #include "Epetra_CrsMatrix.h"
00040
00041 #include "Teuchos_LAPACK.hpp"
00042 #include "Teuchos_CommandLineProcessor.hpp"
00043
00044 #ifdef EPETRA_MPI
00045 #include "Epetra_MpiComm.h"
00046 #else
00047 #include "Epetra_SerialComm.h"
00048 #endif
00049 #include "Epetra_Map.h"
00050
00051 int main(int argc, char *argv[]) {
00052
00053 using std::cout;
00054
00055 #ifdef EPETRA_MPI
00056
00057 MPI_Init(&argc,&argv);
00058 Epetra_MpiComm Comm(MPI_COMM_WORLD);
00059 #else
00060 Epetra_SerialComm Comm;
00061 #endif
00062
00063 bool boolret;
00064 int MyPID = Comm.MyPID();
00065
00066 bool verbose = true;
00067 bool debug = false;
00068 std::string which("SM");
00069
00070 Teuchos::CommandLineProcessor cmdp(false,true);
00071 cmdp.setOption("verbose","quiet",&verbose,"Print messages and results.");
00072 cmdp.setOption("debug","nodebug",&debug,"Print debugging information.");
00073 cmdp.setOption("sort",&which,"Targetted eigenvalues (SM,LM,SR,LR,SI,or LI).");
00074 if (cmdp.parse(argc,argv) != Teuchos::CommandLineProcessor::PARSE_SUCCESSFUL) {
00075 #ifdef HAVE_MPI
00076 MPI_Finalize();
00077 #endif
00078 return -1;
00079 }
00080
00081 typedef double ScalarType;
00082 typedef Teuchos::ScalarTraits<ScalarType> SCT;
00083 typedef SCT::magnitudeType MagnitudeType;
00084 typedef Epetra_MultiVector MV;
00085 typedef Epetra_Operator OP;
00086 typedef Anasazi::MultiVecTraits<ScalarType,MV> MVT;
00087 typedef Anasazi::OperatorTraits<ScalarType,MV,OP> OPT;
00088
00089
00090 int nx = 10;
00091 int NumGlobalElements = nx*nx;
00092
00093
00094
00095
00096 Epetra_Map Map(NumGlobalElements, 0, Comm);
00097
00098
00099
00100 int NumMyElements = Map.NumMyElements();
00101
00102 std::vector<int> MyGlobalElements(NumMyElements);
00103 Map.MyGlobalElements(&MyGlobalElements[0]);
00104
00105
00106
00107
00108 std::vector<int> NumNz(NumMyElements);
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118
00119
00120
00121
00122 for (int i=0; i<NumMyElements; i++) {
00123 if (MyGlobalElements[i] == 0 || MyGlobalElements[i] == NumGlobalElements-1 ||
00124 MyGlobalElements[i] == nx-1 || MyGlobalElements[i] == nx*(nx-1) ) {
00125 NumNz[i] = 3;
00126 }
00127 else if (MyGlobalElements[i] < nx || MyGlobalElements[i] > nx*(nx-1) ||
00128 MyGlobalElements[i]%nx == 0 || (MyGlobalElements[i]+1)%nx == 0) {
00129 NumNz[i] = 4;
00130 }
00131 else {
00132 NumNz[i] = 5;
00133 }
00134 }
00135
00136
00137
00138 Teuchos::RCP<Epetra_CrsMatrix> A = Teuchos::rcp( new Epetra_CrsMatrix(Copy, Map, &NumNz[0]) );
00139
00140
00141
00142
00143 double rho = 2*nx+1;
00144
00145
00146 const double one = 1.0;
00147 std::vector<double> Values(4);
00148 std::vector<int> Indices(4);
00149 double h = one /(nx+1);
00150 double h2 = h*h;
00151 double c = 5.0e-01*rho/ h;
00152 Values[0] = -one/h2 - c; Values[1] = -one/h2 + c; Values[2] = -one/h2; Values[3]= -one/h2;
00153 double diag = 4.0 / h2;
00154 int NumEntries, info;
00155
00156 for (int i=0; i<NumMyElements; i++)
00157 {
00158 if (MyGlobalElements[i]==0)
00159 {
00160 Indices[0] = 1;
00161 Indices[1] = nx;
00162 NumEntries = 2;
00163 info = A->InsertGlobalValues(MyGlobalElements[i], NumEntries, &Values[1], &Indices[0]);
00164 assert( info==0 );
00165 }
00166 else if (MyGlobalElements[i] == nx*(nx-1))
00167 {
00168 Indices[0] = nx*(nx-1)+1;
00169 Indices[1] = nx*(nx-2);
00170 NumEntries = 2;
00171 info = A->InsertGlobalValues(MyGlobalElements[i], NumEntries, &Values[1], &Indices[0]);
00172 assert( info==0 );
00173 }
00174 else if (MyGlobalElements[i] == nx-1)
00175 {
00176 Indices[0] = nx-2;
00177 NumEntries = 1;
00178 info = A->InsertGlobalValues(MyGlobalElements[i], NumEntries, &Values[0], &Indices[0]);
00179 assert( info==0 );
00180 Indices[0] = 2*nx-1;
00181 info = A->InsertGlobalValues(MyGlobalElements[i], NumEntries, &Values[2], &Indices[0]);
00182 assert( info==0 );
00183 }
00184 else if (MyGlobalElements[i] == NumGlobalElements-1)
00185 {
00186 Indices[0] = NumGlobalElements-2;
00187 NumEntries = 1;
00188 info = A->InsertGlobalValues(MyGlobalElements[i], NumEntries, &Values[0], &Indices[0]);
00189 assert( info==0 );
00190 Indices[0] = nx*(nx-1)-1;
00191 info = A->InsertGlobalValues(MyGlobalElements[i], NumEntries, &Values[2], &Indices[0]);
00192 assert( info==0 );
00193 }
00194 else if (MyGlobalElements[i] < nx)
00195 {
00196 Indices[0] = MyGlobalElements[i]-1;
00197 Indices[1] = MyGlobalElements[i]+1;
00198 Indices[2] = MyGlobalElements[i]+nx;
00199 NumEntries = 3;
00200 info = A->InsertGlobalValues(MyGlobalElements[i], NumEntries, &Values[0], &Indices[0]);
00201 assert( info==0 );
00202 }
00203 else if (MyGlobalElements[i] > nx*(nx-1))
00204 {
00205 Indices[0] = MyGlobalElements[i]-1;
00206 Indices[1] = MyGlobalElements[i]+1;
00207 Indices[2] = MyGlobalElements[i]-nx;
00208 NumEntries = 3;
00209 info = A->InsertGlobalValues(MyGlobalElements[i], NumEntries, &Values[0], &Indices[0]);
00210 assert( info==0 );
00211 }
00212 else if (MyGlobalElements[i]%nx == 0)
00213 {
00214 Indices[0] = MyGlobalElements[i]+1;
00215 Indices[1] = MyGlobalElements[i]-nx;
00216 Indices[2] = MyGlobalElements[i]+nx;
00217 NumEntries = 3;
00218 info = A->InsertGlobalValues(MyGlobalElements[i], NumEntries, &Values[1], &Indices[0]);
00219 assert( info==0 );
00220 }
00221 else if ((MyGlobalElements[i]+1)%nx == 0)
00222 {
00223 Indices[0] = MyGlobalElements[i]-nx;
00224 Indices[1] = MyGlobalElements[i]+nx;
00225 NumEntries = 2;
00226 info = A->InsertGlobalValues(MyGlobalElements[i], NumEntries, &Values[2], &Indices[0]);
00227 assert( info==0 );
00228 Indices[0] = MyGlobalElements[i]-1;
00229 NumEntries = 1;
00230 info = A->InsertGlobalValues(MyGlobalElements[i], NumEntries, &Values[0], &Indices[0]);
00231 assert( info==0 );
00232 }
00233 else
00234 {
00235 Indices[0] = MyGlobalElements[i]-1;
00236 Indices[1] = MyGlobalElements[i]+1;
00237 Indices[2] = MyGlobalElements[i]-nx;
00238 Indices[3] = MyGlobalElements[i]+nx;
00239 NumEntries = 4;
00240 info = A->InsertGlobalValues(MyGlobalElements[i], NumEntries, &Values[0], &Indices[0]);
00241 assert( info==0 );
00242 }
00243
00244 info = A->InsertGlobalValues(MyGlobalElements[i], 1, &diag, &MyGlobalElements[i]);
00245 assert( info==0 );
00246 }
00247
00248
00249 info = A->FillComplete();
00250 assert( info==0 );
00251 A->SetTracebackMode(1);
00252
00253
00254
00255
00256
00257
00258
00259 int nev = 4;
00260 int blockSize = 1;
00261 int numBlocks = 20;
00262 int maxRestarts = 500;
00263
00264 double tol = 1e-8;
00265
00266
00267
00268
00269
00270 Teuchos::RCP<Anasazi::SortManager<ScalarType,MV,OP> > MySort =
00271 Teuchos::rcp( new Anasazi::BasicSort<ScalarType,MV,OP>( which ) );
00272
00273
00274 int verbosity = Anasazi::Errors + Anasazi::Warnings;
00275 if (verbose) {
00276 verbosity += Anasazi::FinalSummary + Anasazi::TimingDetails;
00277 }
00278 if (debug) {
00279 verbosity += Anasazi::Debug;
00280 }
00281
00282
00283
00284 Teuchos::ParameterList MyPL;
00285 MyPL.set( "Verbosity", verbosity );
00286 MyPL.set( "Sort Manager", MySort );
00287
00288 MyPL.set( "Block Size", blockSize );
00289 MyPL.set( "Num Blocks", numBlocks );
00290 MyPL.set( "Maximum Restarts", maxRestarts );
00291
00292 MyPL.set( "Convergence Tolerance", tol );
00293
00294
00295
00296 Teuchos::RCP<Epetra_MultiVector> ivec = Teuchos::rcp( new Epetra_MultiVector(Map, blockSize) );
00297 ivec->Random();
00298
00299
00300 Teuchos::RCP<Anasazi::BasicEigenproblem<double, MV, OP> > MyProblem =
00301 Teuchos::rcp( new Anasazi::BasicEigenproblem<double, MV, OP>(A, ivec) );
00302
00303
00304 MyProblem->setHermitian(rho==0.0);
00305
00306
00307 MyProblem->setNEV( nev );
00308
00309
00310 boolret = MyProblem->setProblem();
00311 if (boolret != true) {
00312 if (verbose && MyPID == 0) {
00313 cout << "Anasazi::BasicEigenproblem::setProblem() returned with error." << endl;
00314 }
00315 #ifdef HAVE_MPI
00316 MPI_Finalize() ;
00317 #endif
00318 return -1;
00319 }
00320
00321
00322 Anasazi::BlockKrylovSchurSolMgr<double, MV, OP> MySolverMgr(MyProblem, MyPL);
00323
00324
00325 Anasazi::ReturnType returnCode = MySolverMgr.solve();
00326 if (returnCode != Anasazi::Converged && MyPID==0 && verbose) {
00327 cout << "Anasazi::EigensolverMgr::solve() returned unconverged." << endl;
00328 }
00329
00330
00331 std::vector<Anasazi::Value<double> > ritzValues = MySolverMgr.getRitzValues();
00332
00333
00334 if (verbose && MyPID==0) {
00335 int numritz = (int)ritzValues.size();
00336 cout.setf(std::ios_base::right, std::ios_base::adjustfield);
00337 cout<<endl<< "Computed Ritz Values"<< endl;
00338 if (MyProblem->isHermitian()) {
00339 cout<< std::setw(16) << "Real Part"
00340 << endl;
00341 cout<<"-----------------------------------------------------------"<<endl;
00342 for (int i=0; i<numritz; i++) {
00343 cout<< std::setw(16) << ritzValues[i].realpart
00344 << endl;
00345 }
00346 cout<<"-----------------------------------------------------------"<<endl;
00347 }
00348 else {
00349 cout<< std::setw(16) << "Real Part"
00350 << std::setw(16) << "Imag Part"
00351 << endl;
00352 cout<<"-----------------------------------------------------------"<<endl;
00353 for (int i=0; i<numritz; i++) {
00354 cout<< std::setw(16) << ritzValues[i].realpart
00355 << std::setw(16) << ritzValues[i].imagpart
00356 << endl;
00357 }
00358 cout<<"-----------------------------------------------------------"<<endl;
00359 }
00360 }
00361
00362
00363 Anasazi::Eigensolution<ScalarType,MV> sol = MyProblem->getSolution();
00364 std::vector<Anasazi::Value<ScalarType> > evals = sol.Evals;
00365 Teuchos::RCP<MV> evecs = sol.Evecs;
00366 std::vector<int> index = sol.index;
00367 int numev = sol.numVecs;
00368
00369 if (numev > 0) {
00370
00371 Teuchos::LAPACK<int,double> lapack;
00372 std::vector<double> normA(numev);
00373
00374 if (MyProblem->isHermitian()) {
00375
00376 Epetra_MultiVector Aevecs(Map,numev);
00377 Teuchos::SerialDenseMatrix<int,double> B(numev,numev);
00378 B.putScalar(0.0);
00379 for (int i=0; i<numev; i++) {B(i,i) = evals[i].realpart;}
00380
00381
00382 OPT::Apply( *A, *evecs, Aevecs );
00383
00384
00385 MVT::MvTimesMatAddMv( -1.0, *evecs, B, 1.0, Aevecs );
00386 MVT::MvNorm( Aevecs, &normA );
00387
00388
00389 for (int i=0; i<numev; i++) {
00390 normA[i] /= Teuchos::ScalarTraits<double>::magnitude( evals[i].realpart );
00391 }
00392 } else {
00393
00394 int i=0;
00395 std::vector<int> curind(1);
00396 std::vector<double> resnorm(1), tempnrm(1);
00397 Teuchos::RCP<MV> evecr, eveci, tempAevec;
00398 Epetra_MultiVector Aevec(Map,numev);
00399
00400
00401 OPT::Apply( *A, *evecs, Aevec );
00402
00403 Teuchos::SerialDenseMatrix<int,double> Breal(1,1), Bimag(1,1);
00404 while (i<numev) {
00405 if (index[i]==0) {
00406
00407 curind[0] = i;
00408 evecr = MVT::CloneView( *evecs, curind );
00409
00410
00411 tempAevec = MVT::CloneCopy( Aevec, curind );
00412
00413
00414 Breal(0,0) = evals[i].realpart;
00415 MVT::MvTimesMatAddMv( -1.0, *evecr, Breal, 1.0, *tempAevec );
00416
00417
00418 MVT::MvNorm( *tempAevec, &resnorm );
00419 normA[i] = resnorm[0]/Teuchos::ScalarTraits<MagnitudeType>::magnitude( evals[i].realpart );
00420 i++;
00421 } else {
00422
00423 curind[0] = i;
00424 evecr = MVT::CloneView( *evecs, curind );
00425
00426
00427 tempAevec = MVT::CloneCopy( Aevec, curind );
00428
00429
00430 curind[0] = i+1;
00431 eveci = MVT::CloneView( *evecs, curind );
00432
00433
00434 Breal(0,0) = evals[i].realpart;
00435 Bimag(0,0) = evals[i].imagpart;
00436
00437
00438 MVT::MvTimesMatAddMv( -1.0, *evecr, Breal, 1.0, *tempAevec );
00439 MVT::MvTimesMatAddMv( 1.0, *eveci, Bimag, 1.0, *tempAevec );
00440 MVT::MvNorm( *tempAevec, &tempnrm );
00441
00442
00443 tempAevec = MVT::CloneCopy( Aevec, curind );
00444
00445
00446 MVT::MvTimesMatAddMv( -1.0, *evecr, Bimag, 1.0, *tempAevec );
00447 MVT::MvTimesMatAddMv( -1.0, *eveci, Breal, 1.0, *tempAevec );
00448 MVT::MvNorm( *tempAevec, &resnorm );
00449
00450
00451 normA[i] = lapack.LAPY2( tempnrm[i], resnorm[i] ) /
00452 lapack.LAPY2( evals[i].realpart, evals[i].imagpart );
00453 normA[i+1] = normA[i];
00454
00455 i=i+2;
00456 }
00457 }
00458 }
00459
00460
00461 if (verbose && MyPID==0) {
00462 cout.setf(std::ios_base::right, std::ios_base::adjustfield);
00463 cout<<endl<< "Actual Residuals"<<endl;
00464 if (MyProblem->isHermitian()) {
00465 cout<< std::setw(16) << "Real Part"
00466 << std::setw(20) << "Direct Residual"<< endl;
00467 cout<<"-----------------------------------------------------------"<<endl;
00468 for (int i=0; i<numev; i++) {
00469 cout<< std::setw(16) << evals[i].realpart
00470 << std::setw(20) << normA[i] << endl;
00471 }
00472 cout<<"-----------------------------------------------------------"<<endl;
00473 }
00474 else {
00475 cout<< std::setw(16) << "Real Part"
00476 << std::setw(16) << "Imag Part"
00477 << std::setw(20) << "Direct Residual"<< endl;
00478 cout<<"-----------------------------------------------------------"<<endl;
00479 for (int i=0; i<numev; i++) {
00480 cout<< std::setw(16) << evals[i].realpart
00481 << std::setw(16) << evals[i].imagpart
00482 << std::setw(20) << normA[i] << endl;
00483 }
00484 cout<<"-----------------------------------------------------------"<<endl;
00485 }
00486 }
00487 }
00488
00489 #ifdef EPETRA_MPI
00490 MPI_Finalize();
00491 #endif
00492
00493 return 0;
00494 }