|
Anasazi
Version of the Day
|
00001 // @HEADER 00002 // *********************************************************************** 00003 // 00004 // Anasazi: Block Eigensolvers Package 00005 // Copyright (2004) Sandia Corporation 00006 // 00007 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive 00008 // license for use of this work by or on behalf of the U.S. Government. 00009 // 00010 // This library is free software; you can redistribute it and/or modify 00011 // it under the terms of the GNU Lesser General Public License as 00012 // published by the Free Software Foundation; either version 2.1 of the 00013 // License, or (at your option) any later version. 00014 // 00015 // This library is distributed in the hope that it will be useful, but 00016 // WITHOUT ANY WARRANTY; without even the implied warranty of 00017 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 00018 // Lesser General Public License for more details. 00019 // 00020 // You should have received a copy of the GNU Lesser General Public 00021 // License along with this library; if not, write to the Free Software 00022 // Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 00023 // USA 00024 // Questions? Contact Michael A. Heroux (maherou@sandia.gov) 00025 // 00026 // *********************************************************************** 00027 // @HEADER 00028 00033 #ifndef ANASAZI_MINRES_HPP 00034 #define ANASAZI_MINRES_HPP 00035 00036 #include "AnasaziConfigDefs.hpp" 00037 #include "Teuchos_TimeMonitor.hpp" 00038 00039 namespace Anasazi { 00040 namespace Experimental { 00041 00042 template <class Scalar, class MV, class OP> 00043 class PseudoBlockMinres 00044 { 00045 typedef Anasazi::MultiVecTraits<Scalar,MV> MVT; 00046 typedef Anasazi::OperatorTraits<Scalar,MV,OP> OPT; 00047 const Scalar ONE; 00048 const Scalar ZERO; 00049 00050 public: 00051 // Constructor 00052 PseudoBlockMinres(Teuchos::RCP<OP> A, Teuchos::RCP<OP> Prec = Teuchos::null); 00053 00054 // Set tolerance and maximum iterations 00055 void setTol(const std::vector<Scalar>& tolerances) { tolerances_ = tolerances; }; 00056 void setMaxIter(const int maxIt) { maxIt_ = maxIt; }; 00057 00058 // Set solution and RHS 00059 void setSol(Teuchos::RCP<MV> X) { X_ = X; }; 00060 void setRHS(Teuchos::RCP<const MV> B) { B_ = B; }; 00061 00062 // Solve the linear system 00063 void solve(); 00064 00065 private: 00066 Teuchos::RCP<OP> A_; 00067 Teuchos::RCP<OP> Prec_; 00068 Teuchos::RCP<MV> X_; 00069 Teuchos::RCP<const MV> B_; 00070 std::vector<Scalar> tolerances_; 00071 int maxIt_; 00072 00073 void symOrtho(Scalar a, Scalar b, Scalar *c, Scalar *s, Scalar *r); 00074 00075 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00076 Teuchos::RCP<Teuchos::Time> AddTime_, ApplyOpTime_, ApplyPrecTime_, AssignTime_, DotTime_, LockTime_, NormTime_, ScaleTime_, TotalTime_; 00077 #endif 00078 }; 00079 00080 00081 00082 template <class Scalar, class MV, class OP> 00083 PseudoBlockMinres<Scalar,MV,OP>::PseudoBlockMinres(Teuchos::RCP<OP> A, Teuchos::RCP<OP> Prec) : 00084 ONE(Teuchos::ScalarTraits<Scalar>::one()), 00085 ZERO(Teuchos::ScalarTraits<Scalar>::zero()) 00086 { 00087 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00088 AddTime_ = Teuchos::TimeMonitor::getNewTimer("Anasazi: PseudoBlockMinres::Add Multivectors"); 00089 ApplyOpTime_ = Teuchos::TimeMonitor::getNewTimer("Anasazi: PseudoBlockMinres::Apply Operator"); 00090 ApplyPrecTime_ = Teuchos::TimeMonitor::getNewTimer("Anasazi: PseudoBlockMinres::Apply Preconditioner"); 00091 AssignTime_ = Teuchos::TimeMonitor::getNewTimer("Anasazi: PseudoBlockMinres::Assignment (no locking)"); 00092 DotTime_ = Teuchos::TimeMonitor::getNewTimer("Anasazi: PseudoBlockMinres::Compute Dot Product"); 00093 LockTime_ = Teuchos::TimeMonitor::getNewTimer("Anasazi: PseudoBlockMinres::Lock Converged Vectors"); 00094 NormTime_ = Teuchos::TimeMonitor::getNewTimer("Anasazi: PseudoBlockMinres::Compute Norm"); 00095 ScaleTime_ = Teuchos::TimeMonitor::getNewTimer("Anasazi: PseudoBlockMinres::Scale Multivector"); 00096 TotalTime_ = Teuchos::TimeMonitor::getNewTimer("Anasazi: PseudoBlockMinres::Total Time"); 00097 #endif 00098 00099 A_ = A; 00100 Prec_ = Prec; 00101 } 00102 00103 00104 00105 template <class Scalar, class MV, class OP> 00106 void PseudoBlockMinres<Scalar,MV,OP>::solve() 00107 { 00108 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00109 Teuchos::TimeMonitor outertimer( *TotalTime_ ); 00110 #endif 00111 00112 // Get number of vectors 00113 int ncols = MVT::GetNumberVecs(*B_); 00114 int newNumConverged; 00115 00116 // Declare some variables 00117 std::vector<Scalar> alpha(ncols), beta, beta1(ncols), phibar, oldBeta(ncols,ZERO), epsln(ncols,ZERO), cs(ncols,-ONE), sn(ncols,ZERO), dbar(ncols,ZERO), oldeps(ncols), delta(ncols), gbar(ncols), phi(ncols), gamma(ncols), tmpvec(ncols); 00118 std::vector<int> indicesToRemove, newlyConverged, unconvergedIndices(ncols); 00119 Teuchos::RCP<MV> V, Y, R1, R2, W, W1, W2, locX, scaleHelper, swapHelper; 00120 00121 // Allocate space for multivectors 00122 V = MVT::Clone(*B_, ncols); 00123 Y = MVT::Clone(*B_, ncols); 00124 R1 = MVT::Clone(*B_, ncols); 00125 R2 = MVT::Clone(*B_, ncols); 00126 W = MVT::Clone(*B_, ncols); 00127 W1 = MVT::Clone(*B_, ncols); 00128 W2 = MVT::Clone(*B_, ncols); 00129 scaleHelper = MVT::Clone(*B_, ncols); 00130 00131 // Reserve space for arrays 00132 indicesToRemove.reserve(ncols); 00133 newlyConverged.reserve(ncols); 00134 00135 // Initialize array of unconverged indices 00136 for(int i=0; i<ncols; i++) 00137 unconvergedIndices[i] = i; 00138 00139 // Get a local copy of X 00140 // We want the vectors to remain contiguous even as things converge 00141 locX = MVT::CloneCopy(*X_); 00142 00143 // Initialize residuals 00144 // R1 = B - AX 00145 { 00146 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00147 Teuchos::TimeMonitor lcltimer( *ApplyOpTime_ ); 00148 #endif 00149 OPT::Apply(*A_,*locX,*R1); 00150 } 00151 { 00152 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00153 Teuchos::TimeMonitor lcltimer( *AddTime_ ); 00154 #endif 00155 MVT::MvAddMv(ONE, *B_, -ONE, *R1, *R1); 00156 } 00157 00158 // R2 = R1 00159 { 00160 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00161 Teuchos::TimeMonitor lcltimer( *AssignTime_ ); 00162 #endif 00163 MVT::Assign(*R1,*R2); 00164 } 00165 00166 // Initialize the W's to 0. 00167 MVT::MvInit (*W); 00168 MVT::MvInit (*W2); 00169 00170 // Y = M\R1 (preconditioned residual) 00171 if(Prec_ != Teuchos::null) 00172 { 00173 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00174 Teuchos::TimeMonitor lcltimer( *ApplyPrecTime_ ); 00175 #endif 00176 00177 OPT::Apply(*Prec_,*R1,*Y); 00178 } 00179 else 00180 { 00181 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00182 Teuchos::TimeMonitor lcltimer( *AssignTime_ ); 00183 #endif 00184 MVT::Assign(*R1,*Y); 00185 } 00186 00187 // beta1 = sqrt(Y'*R1) 00188 { 00189 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00190 Teuchos::TimeMonitor lcltimer( *DotTime_ ); 00191 #endif 00192 MVT::MvDot(*Y,*R1, beta1); 00193 00194 for(size_t i=0; i<beta1.size(); i++) 00195 beta1[i] = sqrt(beta1[i]); 00196 } 00197 00198 // beta = beta1 00199 beta = beta1; 00200 00201 // phibar = beta1 00202 phibar = beta1; 00203 00205 // Begin Lanczos iterations 00206 for(int iter=1; iter <= maxIt_; iter++) 00207 { 00208 // Test convergence 00209 indicesToRemove.clear(); 00210 for(int i=0; i<ncols; i++) 00211 { 00212 Scalar relres = phibar[i]/beta1[i]; 00213 // std::cout << "relative residual[" << unconvergedIndices[i] << "] at iteration " << iter << " = " << relres << std::endl; 00214 00215 // If the vector converged, mark it for termination 00216 // Make sure to add it to 00217 if(relres < tolerances_[i]) 00218 { 00219 indicesToRemove.push_back(i); 00220 } 00221 } 00222 00223 // Check whether anything converged 00224 newNumConverged = indicesToRemove.size(); 00225 if(newNumConverged > 0) 00226 { 00227 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00228 Teuchos::TimeMonitor lcltimer( *LockTime_ ); 00229 #endif 00230 00231 // If something converged, stick the converged vectors in X 00232 newlyConverged.resize(newNumConverged); 00233 for(int i=0; i<newNumConverged; i++) 00234 newlyConverged[i] = unconvergedIndices[indicesToRemove[i]]; 00235 00236 Teuchos::RCP<MV> helperLocX = MVT::CloneCopy(*locX,indicesToRemove); 00237 00238 MVT::SetBlock(*helperLocX,newlyConverged,*X_); 00239 00240 // If everything has converged, we are done 00241 if(newNumConverged == ncols) 00242 return; 00243 00244 // Update unconverged indices 00245 std::vector<int> helperVec(ncols - newNumConverged); 00246 00247 std::set_difference(unconvergedIndices.begin(), unconvergedIndices.end(), newlyConverged.begin(), newlyConverged.end(), helperVec.begin()); 00248 unconvergedIndices = helperVec; 00249 00250 // Get indices of things we want to keep 00251 std::vector<int> thingsToKeep(ncols - newNumConverged); 00252 helperVec.resize(ncols); 00253 for(int i=0; i<ncols; i++) 00254 helperVec[i] = i; 00255 ncols = ncols - newNumConverged; 00256 00257 std::set_difference(helperVec.begin(), helperVec.end(), indicesToRemove.begin(), indicesToRemove.end(), thingsToKeep.begin()); 00258 00259 // Shrink the multivectors 00260 Teuchos::RCP<MV> helperMV; 00261 helperMV = MVT::CloneCopy(*V,thingsToKeep); 00262 V = helperMV; 00263 helperMV = MVT::CloneCopy(*Y,thingsToKeep); 00264 Y = helperMV; 00265 helperMV = MVT::CloneCopy(*R1,thingsToKeep); 00266 R1 = helperMV; 00267 helperMV = MVT::CloneCopy(*R2,thingsToKeep); 00268 R2 = helperMV; 00269 helperMV = MVT::CloneCopy(*W,thingsToKeep); 00270 W = helperMV; 00271 helperMV = MVT::CloneCopy(*W1,thingsToKeep); 00272 W1 = helperMV; 00273 helperMV = MVT::CloneCopy(*W2,thingsToKeep); 00274 W2 = helperMV; 00275 helperMV = MVT::CloneCopy(*locX,thingsToKeep); 00276 locX = helperMV; 00277 scaleHelper = MVT::Clone(*V,ncols); 00278 00279 // Shrink the arrays 00280 alpha.resize(ncols); 00281 oldeps.resize(ncols); 00282 delta.resize(ncols); 00283 gbar.resize(ncols); 00284 phi.resize(ncols); 00285 gamma.resize(ncols); 00286 tmpvec.resize(ncols); 00287 std::vector<Scalar> helperVecS(ncols); 00288 for(int i=0; i<ncols; i++) 00289 helperVecS[i] = beta[thingsToKeep[i]]; 00290 beta = helperVecS; 00291 for(int i=0; i<ncols; i++) 00292 helperVecS[i] = beta1[thingsToKeep[i]]; 00293 beta1 = helperVecS; 00294 for(int i=0; i<ncols; i++) 00295 helperVecS[i] = phibar[thingsToKeep[i]]; 00296 phibar = helperVecS; 00297 for(int i=0; i<ncols; i++) 00298 helperVecS[i] = oldBeta[thingsToKeep[i]]; 00299 oldBeta = helperVecS; 00300 for(int i=0; i<ncols; i++) 00301 helperVecS[i] = epsln[thingsToKeep[i]]; 00302 epsln = helperVecS; 00303 for(int i=0; i<ncols; i++) 00304 helperVecS[i] = cs[thingsToKeep[i]]; 00305 cs = helperVecS; 00306 for(int i=0; i<ncols; i++) 00307 helperVecS[i] = sn[thingsToKeep[i]]; 00308 sn = helperVecS; 00309 for(int i=0; i<ncols; i++) 00310 helperVecS[i] = dbar[thingsToKeep[i]]; 00311 dbar = helperVecS; 00312 00313 // Tell operator about the removed indices 00314 A_->removeIndices(indicesToRemove); 00315 } 00316 00317 // Normalize previous vector 00318 // V = Y / beta 00319 { 00320 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00321 Teuchos::TimeMonitor lcltimer( *AssignTime_ ); 00322 #endif 00323 MVT::Assign(*Y,*V); 00324 } 00325 for(int i=0; i<ncols; i++) 00326 tmpvec[i] = ONE/beta[i]; 00327 { 00328 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00329 Teuchos::TimeMonitor lcltimer( *ScaleTime_ ); 00330 #endif 00331 MVT::MvScale (*V, tmpvec); 00332 } 00333 00334 // Apply operator 00335 // Y = AV 00336 { 00337 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00338 Teuchos::TimeMonitor lcltimer( *ApplyOpTime_ ); 00339 #endif 00340 OPT::Apply(*A_, *V, *Y); 00341 } 00342 00343 if(iter > 1) 00344 { 00345 // Y = Y - beta/oldBeta R1 00346 for(int i=0; i<ncols; i++) 00347 tmpvec[i] = beta[i]/oldBeta[i]; 00348 { 00349 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00350 Teuchos::TimeMonitor lcltimer( *AssignTime_ ); 00351 #endif 00352 MVT::Assign(*R1,*scaleHelper); 00353 } 00354 { 00355 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00356 Teuchos::TimeMonitor lcltimer( *ScaleTime_ ); 00357 #endif 00358 MVT::MvScale(*scaleHelper,tmpvec); 00359 } 00360 { 00361 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00362 Teuchos::TimeMonitor lcltimer( *AddTime_ ); 00363 #endif 00364 MVT::MvAddMv(ONE, *Y, -ONE, *scaleHelper, *Y); 00365 } 00366 } 00367 00368 // alpha = V'*Y 00369 { 00370 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00371 Teuchos::TimeMonitor lcltimer( *DotTime_ ); 00372 #endif 00373 MVT::MvDot(*V, *Y, alpha); 00374 } 00375 00376 // Y = Y - alpha/beta R2 00377 for(int i=0; i<ncols; i++) 00378 tmpvec[i] = alpha[i]/beta[i]; 00379 { 00380 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00381 Teuchos::TimeMonitor lcltimer( *AssignTime_ ); 00382 #endif 00383 MVT::Assign(*R2,*scaleHelper); 00384 } 00385 { 00386 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00387 Teuchos::TimeMonitor lcltimer( *ScaleTime_ ); 00388 #endif 00389 MVT::MvScale(*scaleHelper,tmpvec); 00390 } 00391 { 00392 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00393 Teuchos::TimeMonitor lcltimer( *AddTime_ ); 00394 #endif 00395 MVT::MvAddMv(ONE, *Y, -ONE, *scaleHelper, *Y); 00396 } 00397 00398 // R1 = R2 00399 // R2 = Y 00400 swapHelper = R1; 00401 R1 = R2; 00402 R2 = Y; 00403 Y = swapHelper; 00404 00405 // Y = M\R2 00406 if(Prec_ != Teuchos::null) 00407 { 00408 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00409 Teuchos::TimeMonitor lcltimer( *ApplyPrecTime_ ); 00410 #endif 00411 00412 OPT::Apply(*Prec_,*R2,*Y); 00413 } 00414 else 00415 { 00416 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00417 Teuchos::TimeMonitor lcltimer( *AssignTime_ ); 00418 #endif 00419 MVT::Assign(*R2,*Y); 00420 } 00421 00422 // Get new beta 00423 // beta = sqrt(R2'*Y) 00424 oldBeta = beta; 00425 { 00426 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00427 Teuchos::TimeMonitor lcltimer( *DotTime_ ); 00428 #endif 00429 MVT::MvDot(*Y,*R2, beta); 00430 00431 for(int i=0; i<ncols; i++) 00432 beta[i] = sqrt(beta[i]); 00433 } 00434 00435 // Apply previous rotation 00436 oldeps = epsln; 00437 for(int i=0; i<ncols; i++) 00438 { 00439 delta[i] = cs[i]*dbar[i] + sn[i]*alpha[i]; 00440 gbar[i] = sn[i]*dbar[i] - cs[i]*alpha[i]; 00441 epsln[i] = sn[i]*beta[i]; 00442 dbar[i] = - cs[i]*beta[i]; 00443 } 00444 00445 // Compute the next plane rotation 00446 for(int i=0; i<ncols; i++) 00447 { 00448 symOrtho(gbar[i], beta[i], &cs[i], &sn[i], &gamma[i]); 00449 00450 phi[i] = cs[i]*phibar[i]; 00451 phibar[i] = sn[i]*phibar[i]; 00452 } 00453 00454 // w1 = w2 00455 // w2 = w 00456 swapHelper = W1; 00457 W1 = W2; 00458 W2 = W; 00459 W = swapHelper; 00460 00461 // W = (V - oldeps*W1 - delta*W2) / gamma 00462 { 00463 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00464 Teuchos::TimeMonitor lcltimer( *AssignTime_ ); 00465 #endif 00466 MVT::Assign(*W1,*scaleHelper); 00467 } 00468 { 00469 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00470 Teuchos::TimeMonitor lcltimer( *ScaleTime_ ); 00471 #endif 00472 MVT::MvScale(*scaleHelper,oldeps); 00473 } 00474 { 00475 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00476 Teuchos::TimeMonitor lcltimer( *AddTime_ ); 00477 #endif 00478 MVT::MvAddMv(ONE, *V, -ONE, *scaleHelper, *W); 00479 } 00480 { 00481 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00482 Teuchos::TimeMonitor lcltimer( *AssignTime_ ); 00483 #endif 00484 MVT::Assign(*W2,*scaleHelper); 00485 } 00486 { 00487 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00488 Teuchos::TimeMonitor lcltimer( *ScaleTime_ ); 00489 #endif 00490 MVT::MvScale(*scaleHelper,delta); 00491 } 00492 { 00493 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00494 Teuchos::TimeMonitor lcltimer( *AddTime_ ); 00495 #endif 00496 MVT::MvAddMv(ONE, *W, -ONE, *scaleHelper, *W); 00497 } 00498 for(int i=0; i<ncols; i++) 00499 tmpvec[i] = ONE/gamma[i]; 00500 { 00501 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00502 Teuchos::TimeMonitor lcltimer( *ScaleTime_ ); 00503 #endif 00504 MVT::MvScale(*W,tmpvec); 00505 } 00506 00507 // Update X 00508 // X = X + phi*W 00509 { 00510 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00511 Teuchos::TimeMonitor lcltimer( *AssignTime_ ); 00512 #endif 00513 MVT::Assign(*W,*scaleHelper); 00514 } 00515 { 00516 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00517 Teuchos::TimeMonitor lcltimer( *ScaleTime_ ); 00518 #endif 00519 MVT::MvScale(*scaleHelper,phi); 00520 } 00521 { 00522 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00523 Teuchos::TimeMonitor lcltimer( *AddTime_ ); 00524 #endif 00525 MVT::MvAddMv(ONE, *locX, ONE, *scaleHelper, *locX); 00526 } 00527 } 00528 00529 // Stick unconverged vectors in X 00530 { 00531 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR 00532 Teuchos::TimeMonitor lcltimer( *AssignTime_ ); 00533 #endif 00534 MVT::SetBlock(*locX,unconvergedIndices,*X_); 00535 } 00536 } 00537 00538 template <class Scalar, class MV, class OP> 00539 void PseudoBlockMinres<Scalar,MV,OP>::symOrtho(Scalar a, Scalar b, Scalar *c, Scalar *s, Scalar *r) 00540 { 00541 const Scalar absA = std::abs(a); 00542 const Scalar absB = std::abs(b); 00543 if ( absB == ZERO ) { 00544 *s = ZERO; 00545 *r = absA; 00546 if ( absA == ZERO ) 00547 *c = ONE; 00548 else 00549 *c = a / absA; 00550 } else if ( absA == ZERO ) { 00551 *c = ZERO; 00552 *s = b / absB; 00553 *r = absB; 00554 } else if ( absB >= absA ) { // && a!=0 && b!=0 00555 Scalar tau = a / b; 00556 if ( b < ZERO ) 00557 *s = -ONE / sqrt( ONE+tau*tau ); 00558 else 00559 *s = ONE / sqrt( ONE+tau*tau ); 00560 *c = *s * tau; 00561 *r = b / *s; 00562 } else { // (absA > absB) && a!=0 && b!=0 00563 Scalar tau = b / a; 00564 if ( a < ZERO ) 00565 *c = -ONE / sqrt( ONE+tau*tau ); 00566 else 00567 *c = ONE / sqrt( ONE+tau*tau ); 00568 *s = *c * tau; 00569 *r = a / *c; 00570 } 00571 } 00572 00573 }} // End of namespace 00574 00575 #endif
1.7.6.1