Anasazi  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
AnasaziBlockKrylovSchurSolMgr.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Anasazi: Block Eigensolvers Package
5 // Copyright 2004 Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8 // the U.S. Government retains certain rights in this software.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ***********************************************************************
40 // @HEADER
41 
42 #ifndef ANASAZI_BLOCK_KRYLOV_SCHUR_SOLMGR_HPP
43 #define ANASAZI_BLOCK_KRYLOV_SCHUR_SOLMGR_HPP
44 
48 
49 #include "AnasaziConfigDefs.hpp"
50 #include "AnasaziTypes.hpp"
51 
52 #include "AnasaziEigenproblem.hpp"
53 #include "AnasaziSolverManager.hpp"
54 #include "AnasaziSolverUtils.hpp"
55 
57 #include "AnasaziBasicSort.hpp"
64 #include "AnasaziOutputManager.hpp"
66 
67 #include "Teuchos_BLAS.hpp"
68 #include "Teuchos_LAPACK.hpp"
69 #include "Teuchos_TimeMonitor.hpp"
70 #include "Teuchos_FancyOStream.hpp"
71 
78 
126 namespace Anasazi {
127 
128 
155 template<class ScalarType, class MV, class OP>
156 class BlockKrylovSchurSolMgr : public SolverManager<ScalarType,MV,OP> {
157 
158  private:
162  typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
164 
165  public:
166 
168 
169 
192 
196 
198 
199 
202  return *problem_;
203  }
204 
206  int getNumIters() const {
207  return numIters_;
208  }
209 
212  std::vector<Value<ScalarType> > getRitzValues() const {
213  std::vector<Value<ScalarType> > ret( ritzValues_ );
214  return ret;
215  }
216 
224  return Teuchos::tuple(timerSolve_, timerRestarting_);
225  }
226 
228 
230 
231 
250  ReturnType solve();
251 
254 
257 
260 
263 
265 
266  private:
269 
270  std::string whch_, ortho_;
271  MagnitudeType ortho_kappa_;
272 
273  MagnitudeType convtol_;
274  int maxRestarts_;
275  bool relconvtol_,conjSplit_;
276  int blockSize_, numBlocks_, stepSize_, nevBlocks_, xtra_nevBlocks_;
277  int numIters_;
278  int verbosity_;
279  bool inSituRestart_, dynXtraNev_;
280  int osProc_;
281 
282  std::vector<Value<ScalarType> > ritzValues_;
283 
284  int printNum_;
285  Teuchos::RCP<Teuchos::Time> timerSolve_, timerRestarting_;
286 
288 
291 
292 };
293 
294 
295 // Constructor
296 template<class ScalarType, class MV, class OP>
299  Teuchos::ParameterList &pl ) :
300  problem_(problem),
301  whch_("LM"),
302  ortho_("SVQB"),
303  ortho_kappa_(-1.0),
304  convtol_(0),
305  maxRestarts_(20),
306  relconvtol_(true),
307  conjSplit_(false),
308  blockSize_(0),
309  numBlocks_(0),
310  stepSize_(0),
311  nevBlocks_(0),
312  xtra_nevBlocks_(0),
313  numIters_(0),
314  verbosity_(Anasazi::Errors),
315  inSituRestart_(false),
316  dynXtraNev_(false),
317  osProc_(0),
318  printNum_(-1)
319 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
320  ,timerSolve_(Teuchos::TimeMonitor::getNewTimer("Anasazi: BlockKrylovSchurSolMgr::solve()")),
321  timerRestarting_(Teuchos::TimeMonitor::getNewTimer("Anasazi: BlockKrylovSchurSolMgr restarting"))
322 #endif
323 {
324  TEUCHOS_TEST_FOR_EXCEPTION(problem_ == Teuchos::null, std::invalid_argument, "Problem not given to solver manager.");
325  TEUCHOS_TEST_FOR_EXCEPTION(!problem_->isProblemSet(), std::invalid_argument, "Problem not set.");
326  TEUCHOS_TEST_FOR_EXCEPTION(problem_->getInitVec() == Teuchos::null, std::invalid_argument, "Problem does not contain initial vectors to clone from.");
327 
328  const int nev = problem_->getNEV();
329 
330  // convergence tolerance
331  convtol_ = pl.get("Convergence Tolerance",MT::prec());
332  relconvtol_ = pl.get("Relative Convergence Tolerance",relconvtol_);
333 
334  // maximum number of restarts
335  maxRestarts_ = pl.get("Maximum Restarts",maxRestarts_);
336 
337  // block size: default is 1
338  blockSize_ = pl.get("Block Size",1);
339  TEUCHOS_TEST_FOR_EXCEPTION(blockSize_ <= 0, std::invalid_argument,
340  "Anasazi::BlockKrylovSchurSolMgr: \"Block Size\" must be strictly positive.");
341 
342  // set the number of blocks we need to save to compute the nev eigenvalues of interest.
343  xtra_nevBlocks_ = pl.get("Extra NEV Blocks",0);
344  if (nev%blockSize_) {
345  nevBlocks_ = nev/blockSize_ + 1;
346  } else {
347  nevBlocks_ = nev/blockSize_;
348  }
349 
350  // determine if we should use the dynamic scheme for selecting the current number of retained eigenvalues.
351  // NOTE: This employs a technique similar to ARPACK in that it increases the number of retained eigenvalues
352  // by one for every converged eigenpair to accelerate convergence.
353  if (pl.isParameter("Dynamic Extra NEV")) {
354  if (Teuchos::isParameterType<bool>(pl,"Dynamic Extra NEV")) {
355  dynXtraNev_ = pl.get("Dynamic Extra NEV",dynXtraNev_);
356  } else {
357  dynXtraNev_ = ( Teuchos::getParameter<int>(pl,"Dynamic Extra NEV") != 0 );
358  }
359  }
360 
361  numBlocks_ = pl.get("Num Blocks",3*nevBlocks_);
362  TEUCHOS_TEST_FOR_EXCEPTION(numBlocks_ <= nevBlocks_, std::invalid_argument,
363  "Anasazi::BlockKrylovSchurSolMgr: \"Num Blocks\" must be strictly positive and large enough to compute the requested eigenvalues.");
364 
365  TEUCHOS_TEST_FOR_EXCEPTION(static_cast<ptrdiff_t>(numBlocks_)*blockSize_ > MVT::GetGlobalLength(*problem_->getInitVec()),
366  std::invalid_argument,
367  "Anasazi::BlockKrylovSchurSolMgr: Potentially impossible orthogonality requests. Reduce basis size.");
368 
369  // step size: the default is maxRestarts_*numBlocks_, so that Ritz values are only computed every restart.
370  if (maxRestarts_) {
371  stepSize_ = pl.get("Step Size", (maxRestarts_+1)*(numBlocks_+1));
372  } else {
373  stepSize_ = pl.get("Step Size", numBlocks_+1);
374  }
375  TEUCHOS_TEST_FOR_EXCEPTION(stepSize_ < 1, std::invalid_argument,
376  "Anasazi::BlockKrylovSchurSolMgr: \"Step Size\" must be strictly positive.");
377 
378  // get the sort manager
379  if (pl.isParameter("Sort Manager")) {
380  sort_ = Teuchos::getParameter<Teuchos::RCP<Anasazi::SortManager<MagnitudeType> > >(pl,"Sort Manager");
381  } else {
382  // which values to solve for
383  whch_ = pl.get("Which",whch_);
384  TEUCHOS_TEST_FOR_EXCEPTION(whch_ != "SM" && whch_ != "LM" && whch_ != "SR" && whch_ != "LR" && whch_ != "SI" && whch_ != "LI",
385  std::invalid_argument, "Invalid sorting string.");
386  sort_ = Teuchos::rcp( new BasicSort<MagnitudeType>(whch_) );
387  }
388 
389  // which orthogonalization to use
390  ortho_ = pl.get("Orthogonalization",ortho_);
391  if (ortho_ != "DGKS" && ortho_ != "SVQB") {
392  ortho_ = "SVQB";
393  }
394 
395  // which orthogonalization constant to use
396  ortho_kappa_ = pl.get("Orthogonalization Constant",ortho_kappa_);
397 
398  // Create a formatted output stream to print to.
399  // See if user requests output processor.
400  osProc_ = pl.get("Output Processor", osProc_);
401 
402  // If not passed in by user, it will be chosen based upon operator type.
403  if (pl.isParameter("Output Stream")) {
404  osp_ = Teuchos::getParameter<Teuchos::RCP<Teuchos::FancyOStream> >(pl,"Output Stream");
405  }
406  else {
407  osp_ = OutputStreamTraits<OP>::getOutputStream (*problem_->getOperator(),osProc_);
408  }
409 
410  // verbosity level
411  if (pl.isParameter("Verbosity")) {
412  if (Teuchos::isParameterType<int>(pl,"Verbosity")) {
413  verbosity_ = pl.get("Verbosity", verbosity_);
414  } else {
415  verbosity_ = (int)Teuchos::getParameter<Anasazi::MsgType>(pl,"Verbosity");
416  }
417  }
418 
419  // restarting technique: V*Q or applyHouse(V,H,tau)
420  if (pl.isParameter("In Situ Restarting")) {
421  if (Teuchos::isParameterType<bool>(pl,"In Situ Restarting")) {
422  inSituRestart_ = pl.get("In Situ Restarting",inSituRestart_);
423  } else {
424  inSituRestart_ = ( Teuchos::getParameter<int>(pl,"In Situ Restarting") != 0 );
425  }
426  }
427 
428  printNum_ = pl.get<int>("Print Number of Ritz Values",printNum_);
429 }
430 
431 
432 // solve()
433 template<class ScalarType, class MV, class OP>
436 
437  const int nev = problem_->getNEV();
438  ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
439  ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
440 
443  typedef SolverUtils<ScalarType,MV,OP> msutils;
444 
446  // Output manager
448 
450  // Status tests
451  //
452  // convergence
454  if (globalTest_ == Teuchos::null) {
455  convtest = Teuchos::rcp( new StatusTestResNorm<ScalarType,MV,OP>(convtol_,nev,RITZRES_2NORM,relconvtol_) );
456  }
457  else {
458  convtest = globalTest_;
459  }
461  = Teuchos::rcp( new StatusTestWithOrdering<ScalarType,MV,OP>(convtest,sort_,nev) );
462  // for a non-short-circuited OR test, the order doesn't matter
464  alltests.push_back(ordertest);
465 
466  if (debugTest_ != Teuchos::null) alltests.push_back(debugTest_);
467 
470  // printing StatusTest
472  if ( printer->isVerbosity(Debug) ) {
473  outputtest = Teuchos::rcp( new StatusTestOutput<ScalarType,MV,OP>( printer,combotest,1,Passed+Failed+Undefined ) );
474  }
475  else {
476  outputtest = Teuchos::rcp( new StatusTestOutput<ScalarType,MV,OP>( printer,combotest,1,Passed ) );
477  }
478 
480  // Orthomanager
482  if (ortho_=="SVQB") {
483  ortho = Teuchos::rcp( new SVQBOrthoManager<ScalarType,MV,OP>(problem_->getM()) );
484  } else if (ortho_=="DGKS") {
485  if (ortho_kappa_ <= 0) {
486  ortho = Teuchos::rcp( new BasicOrthoManager<ScalarType,MV,OP>(problem_->getM()) );
487  }
488  else {
489  ortho = Teuchos::rcp( new BasicOrthoManager<ScalarType,MV,OP>(problem_->getM(),ortho_kappa_) );
490  }
491  } else {
492  TEUCHOS_TEST_FOR_EXCEPTION(ortho_!="SVQB"&&ortho_!="DGKS",std::logic_error,"Anasazi::BlockKrylovSchurSolMgr::solve(): Invalid orthogonalization type.");
493  }
494 
496  // Parameter list
498  plist.set("Block Size",blockSize_);
499  plist.set("Num Blocks",numBlocks_);
500  plist.set("Step Size",stepSize_);
501  if (printNum_ == -1) {
502  plist.set("Print Number of Ritz Values",nevBlocks_*blockSize_);
503  }
504  else {
505  plist.set("Print Number of Ritz Values",printNum_);
506  }
507 
509  // BlockKrylovSchur solver
511  = Teuchos::rcp( new BlockKrylovSchur<ScalarType,MV,OP>(problem_,sort_,printer,outputtest,ortho,plist) );
512  // set any auxiliary vectors defined in the problem
513  Teuchos::RCP< const MV > probauxvecs = problem_->getAuxVecs();
514  if (probauxvecs != Teuchos::null) {
515  bks_solver->setAuxVecs( Teuchos::tuple< Teuchos::RCP<const MV> >(probauxvecs) );
516  }
517 
518  // Create workspace for the Krylov basis generated during a restart
519  // Need at most (nevBlocks_*blockSize_+1) for the updated factorization and another block for the current factorization residual block (F).
520  // ---> (nevBlocks_*blockSize_+1) + blockSize_
521  // If Hermitian, this becomes nevBlocks_*blockSize_ + blockSize_
522  // we only need this if there is the possibility of restarting, ex situ
523 
524  // Maximum allowable extra vectors that BKS can keep to accelerate convergence
525  int maxXtraBlocks = 0;
526  if ( dynXtraNev_ ) maxXtraBlocks = ( ( bks_solver->getMaxSubspaceDim() - nev ) / blockSize_ ) / 2;
527 
528  Teuchos::RCP<MV> workMV;
529  if (maxRestarts_ > 0) {
530  if (inSituRestart_==true) {
531  // still need one work vector for applyHouse()
532  workMV = MVT::Clone( *problem_->getInitVec(), 1 );
533  }
534  else { // inSituRestart == false
535  if (problem_->isHermitian()) {
536  workMV = MVT::Clone( *problem_->getInitVec(), (nevBlocks_+maxXtraBlocks)*blockSize_ + blockSize_ );
537  } else {
538  // Increase workspace by 1 because of potential complex conjugate pairs on the Ritz vector boundary
539  workMV = MVT::Clone( *problem_->getInitVec(), (nevBlocks_+maxXtraBlocks)*blockSize_ + 1 + blockSize_ );
540  }
541  }
542  } else {
543  workMV = Teuchos::null;
544  }
545 
546  // go ahead and initialize the solution to nothing in case we throw an exception
548  sol.numVecs = 0;
549  problem_->setSolution(sol);
550 
551  int numRestarts = 0;
552  int cur_nevBlocks = 0;
553 
554  // enter solve() iterations
555  {
556 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
557  Teuchos::TimeMonitor slvtimer(*timerSolve_);
558 #endif
559 
560  // tell bks_solver to iterate
561  while (1) {
562  try {
563  bks_solver->iterate();
564 
566  //
567  // check convergence first
568  //
570  if ( ordertest->getStatus() == Passed ) {
571  // we have convergence
572  // ordertest->whichVecs() tells us which vectors from solver state are the ones we want
573  // ordertest->howMany() will tell us how many
574  break;
575  }
577  //
578  // check for restarting, i.e. the subspace is full
579  //
581  // this is for the Hermitian case, or non-Hermitian conjugate split situation.
582  // --> for the Hermitian case the current subspace dimension needs to match the maximum subspace dimension
583  // --> for the non-Hermitian case:
584  // --> if a conjugate pair was detected in the previous restart then the current subspace dimension needs to match the
585  // maximum subspace dimension (the BKS solver keeps one extra vector if the problem is non-Hermitian).
586  // --> if a conjugate pair was not detected in the previous restart then the current subspace dimension will be one less
587  // than the maximum subspace dimension.
588  else if ( (bks_solver->getCurSubspaceDim() == bks_solver->getMaxSubspaceDim()) ||
589  (!problem_->isHermitian() && !conjSplit_ && (bks_solver->getCurSubspaceDim()+1 == bks_solver->getMaxSubspaceDim())) ) {
590 
591  // Update the Schur form of the projected eigenproblem, then sort it.
592  if (!bks_solver->isSchurCurrent()) {
593  bks_solver->computeSchurForm( true );
594 
595  // Check for convergence, just in case we wait for every restart to check
596  outputtest->checkStatus( &*bks_solver );
597  }
598 
599  // Don't bother to restart if we've converged or reached the maximum number of restarts
600  if ( numRestarts >= maxRestarts_ || ordertest->getStatus() == Passed) {
601  break; // break from while(1){bks_solver->iterate()}
602  }
603 
604  // Start restarting timer and increment counter
605 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
606  Teuchos::TimeMonitor restimer(*timerRestarting_);
607 #endif
608  numRestarts++;
609 
610  int numConv = ordertest->howMany();
611  cur_nevBlocks = nevBlocks_*blockSize_;
612 
613  // Add in extra blocks for restarting if either static or dynamic boundaries are being used.
614  int moreNevBlocks = std::min( maxXtraBlocks, std::max( numConv/blockSize_, xtra_nevBlocks_) );
615  if ( dynXtraNev_ )
616  cur_nevBlocks += moreNevBlocks * blockSize_;
617  else if ( xtra_nevBlocks_ )
618  cur_nevBlocks += xtra_nevBlocks_ * blockSize_;
619 /*
620  int cur_numConv = numConv;
621  while ( (cur_nevBlocks < (_nevBlocks + maxXtraVecs)) && cur_numConv > 0 ) {
622  cur_nevBlocks++;
623  cur_numConv--;
624 */
625 
626  printer->stream(Debug) << " Performing restart number " << numRestarts << " of " << maxRestarts_ << std::endl << std::endl;
627  printer->stream(Debug) << " - Current NEV blocks is " << cur_nevBlocks << ", the minimum is " << nevBlocks_*blockSize_ << std::endl;
628 
629  // Get the most current Ritz values before we continue.
630  ritzValues_ = bks_solver->getRitzValues();
631 
632  // Get the state.
633  BlockKrylovSchurState<ScalarType,MV> oldState = bks_solver->getState();
634 
635  // Get the current dimension of the factorization
636  int curDim = oldState.curDim;
637 
638  // Determine if the storage for the nev eigenvalues of interest splits a complex conjugate pair.
639  std::vector<int> ritzIndex = bks_solver->getRitzIndex();
640  if (ritzIndex[cur_nevBlocks-1]==1) {
641  conjSplit_ = true;
642  cur_nevBlocks++;
643  } else {
644  conjSplit_ = false;
645  }
646 
647  // Print out a warning to the user if complex eigenvalues were found on the boundary of the restart subspace
648  // and the eigenproblem is Hermitian. This solver is not prepared to handle this situation.
649  if (problem_->isHermitian() && conjSplit_)
650  {
651  printer->stream(Warnings)
652  << " Eigenproblem is Hermitian, complex eigenvalues have been detected, and eigenvalues of interest split a conjugate pair!!!"
653  << std::endl
654  << " Block Krylov-Schur eigensolver cannot guarantee correct behavior in this situation, please turn Hermitian flag off!!!"
655  << std::endl;
656  }
657 
658  // Update the Krylov-Schur decomposition
659 
660  // Get a view of the Schur vectors of interest.
661  Teuchos::SerialDenseMatrix<int,ScalarType> Qnev(Teuchos::View, *(oldState.Q), curDim, cur_nevBlocks);
662 
663  // Get a view of the current Krylov basis.
664  std::vector<int> curind( curDim );
665  for (int i=0; i<curDim; i++) { curind[i] = i; }
666  Teuchos::RCP<const MV> basistemp = MVT::CloneView( *(oldState.V), curind );
667 
668  // Compute the new Krylov basis: Vnew = V*Qnev
669  //
670  // this will occur ex situ in workspace allocated for this purpose (tmpMV)
671  // or in situ in the solver's memory space.
672  //
673  // we will also set a pointer for the location that the current factorization residual block (F),
674  // currently located after the current basis in oldstate.V, will be moved to
675  //
676  Teuchos::RCP<MV> newF;
677  if (inSituRestart_) {
678  //
679  // get non-const pointer to solver's basis so we can work in situ
680  Teuchos::RCP<MV> solverbasis = Teuchos::rcp_const_cast<MV>(oldState.V);
682  //
683  // perform Householder QR of copyQnev = Q [D;0], where D is unit diag. We will want D below.
684  std::vector<ScalarType> tau(cur_nevBlocks), work(cur_nevBlocks);
685  int info;
686  lapack.GEQRF(curDim,cur_nevBlocks,copyQnev.values(),copyQnev.stride(),&tau[0],&work[0],work.size(),&info);
687  TEUCHOS_TEST_FOR_EXCEPTION(info != 0,std::logic_error,
688  "Anasazi::BlockKrylovSchurSolMgr::solve(): error calling GEQRF during restarting.");
689  // we need to get the diagonal of D
690  std::vector<ScalarType> d(cur_nevBlocks);
691  for (int j=0; j<copyQnev.numCols(); j++) {
692  d[j] = copyQnev(j,j);
693  }
694  if (printer->isVerbosity(Debug)) {
695  Teuchos::SerialDenseMatrix<int,ScalarType> R(Teuchos::Copy,copyQnev,cur_nevBlocks,cur_nevBlocks);
696  for (int j=0; j<R.numCols(); j++) {
697  R(j,j) = SCT::magnitude(R(j,j)) - 1.0;
698  for (int i=j+1; i<R.numRows(); i++) {
699  R(i,j) = zero;
700  }
701  }
702  printer->stream(Debug) << "||Triangular factor of Su - I||: " << R.normFrobenius() << std::endl;
703  }
704  //
705  // perform implicit V*Qnev
706  // this actually performs V*[Qnev Qtrunc*M] = [newV truncV], for some unitary M
707  // we are interested in only the first cur_nevBlocks vectors of the result
708  curind.resize(curDim);
709  for (int i=0; i<curDim; i++) curind[i] = i;
710  {
711  Teuchos::RCP<MV> oldV = MVT::CloneViewNonConst(*solverbasis,curind);
712  msutils::applyHouse(cur_nevBlocks,*oldV,copyQnev,tau,workMV);
713  }
714  // multiply newV*D
715  // get pointer to new basis
716  curind.resize(cur_nevBlocks);
717  for (int i=0; i<cur_nevBlocks; i++) { curind[i] = i; }
718  {
719  Teuchos::RCP<MV> newV = MVT::CloneViewNonConst( *solverbasis, curind );
720  MVT::MvScale(*newV,d);
721  }
722  // get pointer to new location for F
723  curind.resize(blockSize_);
724  for (int i=0; i<blockSize_; i++) { curind[i] = cur_nevBlocks + i; }
725  newF = MVT::CloneViewNonConst( *solverbasis, curind );
726  }
727  else {
728  // get pointer to first part of work space
729  curind.resize(cur_nevBlocks);
730  for (int i=0; i<cur_nevBlocks; i++) { curind[i] = i; }
731  Teuchos::RCP<MV> tmp_newV = MVT::CloneViewNonConst(*workMV, curind );
732  // perform V*Qnev
733  MVT::MvTimesMatAddMv( one, *basistemp, Qnev, zero, *tmp_newV );
734  tmp_newV = Teuchos::null;
735  // get pointer to new location for F
736  curind.resize(blockSize_);
737  for (int i=0; i<blockSize_; i++) { curind[i] = cur_nevBlocks + i; }
738  newF = MVT::CloneViewNonConst( *workMV, curind );
739  }
740 
741  // Move the current factorization residual block (F) to the last block of newV.
742  curind.resize(blockSize_);
743  for (int i=0; i<blockSize_; i++) { curind[i] = curDim + i; }
744  Teuchos::RCP<const MV> oldF = MVT::CloneView( *(oldState.V), curind );
745  for (int i=0; i<blockSize_; i++) { curind[i] = i; }
746  MVT::SetBlock( *oldF, curind, *newF );
747  newF = Teuchos::null;
748 
749  // Update the Krylov-Schur quasi-triangular matrix.
750  //
751  // Create storage for the new Schur matrix of the Krylov-Schur factorization
752  // Copy over the current quasi-triangular factorization of oldState.H which is stored in oldState.S.
754  Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::Copy, *(oldState.S), cur_nevBlocks+blockSize_, cur_nevBlocks) );
755  //
756  // Get a view of the B block of the current factorization
757  Teuchos::SerialDenseMatrix<int,ScalarType> oldB(Teuchos::View, *(oldState.H), blockSize_, blockSize_, curDim, curDim-blockSize_);
758  //
759  // Get a view of the a block row of the Schur vectors.
760  Teuchos::SerialDenseMatrix<int,ScalarType> subQ(Teuchos::View, *(oldState.Q), blockSize_, cur_nevBlocks, curDim-blockSize_);
761  //
762  // Get a view of the new B block of the updated Krylov-Schur factorization
763  Teuchos::SerialDenseMatrix<int,ScalarType> newB(Teuchos::View, *newH, blockSize_, cur_nevBlocks, cur_nevBlocks);
764  //
765  // Compute the new B block.
766  blas.GEMM( Teuchos::NO_TRANS, Teuchos::NO_TRANS, blockSize_, cur_nevBlocks, blockSize_, one,
767  oldB.values(), oldB.stride(), subQ.values(), subQ.stride(), zero, newB.values(), newB.stride() );
768 
769 
770  //
771  // Set the new state and initialize the solver.
773  if (inSituRestart_) {
774  newstate.V = oldState.V;
775  } else {
776  newstate.V = workMV;
777  }
778  newstate.H = newH;
779  newstate.curDim = cur_nevBlocks;
780  bks_solver->initialize(newstate);
781 
782  } // end of restarting
784  //
785  // we returned from iterate(), but none of our status tests Passed.
786  // something is wrong, and it is probably our fault.
787  //
789  else {
790  TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,"Anasazi::BlockKrylovSchurSolMgr::solve(): Invalid return from bks_solver::iterate().");
791  }
792  }
793  catch (const AnasaziError &err) {
794  printer->stream(Errors)
795  << "Anasazi::BlockKrylovSchurSolMgr::solve() caught unexpected exception from Anasazi::BlockKrylovSchur::iterate() at iteration " << bks_solver->getNumIters() << std::endl
796  << err.what() << std::endl
797  << "Anasazi::BlockKrylovSchurSolMgr::solve() returning Unconverged with no solutions." << std::endl;
798  return Unconverged;
799  }
800  }
801 
802  //
803  // free temporary space
804  workMV = Teuchos::null;
805 
806  // Get the most current Ritz values before we return
807  ritzValues_ = bks_solver->getRitzValues();
808 
809  sol.numVecs = ordertest->howMany();
810  printer->stream(Debug) << "ordertest->howMany() : " << sol.numVecs << std::endl;
811  std::vector<int> whichVecs = ordertest->whichVecs();
812 
813  // Place any converged eigenpairs in the solution container.
814  if (sol.numVecs > 0) {
815 
816  // Next determine if there is a conjugate pair on the boundary and resize.
817  std::vector<int> tmpIndex = bks_solver->getRitzIndex();
818  for (int i=0; i<(int)ritzValues_.size(); ++i) {
819  printer->stream(Debug) << ritzValues_[i].realpart << " + i " << ritzValues_[i].imagpart << ", Index = " << tmpIndex[i] << std::endl;
820  }
821  printer->stream(Debug) << "Number of converged eigenpairs (before) = " << sol.numVecs << std::endl;
822  for (int i=0; i<sol.numVecs; ++i) {
823  printer->stream(Debug) << "whichVecs[" << i << "] = " << whichVecs[i] << ", tmpIndex[" << whichVecs[i] << "] = " << tmpIndex[whichVecs[i]] << std::endl;
824  }
825  if (tmpIndex[whichVecs[sol.numVecs-1]]==1) {
826  printer->stream(Debug) << "There is a conjugate pair on the boundary, resizing sol.numVecs" << std::endl;
827  whichVecs.push_back(whichVecs[sol.numVecs-1]+1);
828  sol.numVecs++;
829  for (int i=0; i<sol.numVecs; ++i) {
830  printer->stream(Debug) << "whichVecs[" << i << "] = " << whichVecs[i] << ", tmpIndex[" << whichVecs[i] << "] = " << tmpIndex[whichVecs[i]] << std::endl;
831  }
832  }
833 
834  bool keepMore = false;
835  int numEvecs = sol.numVecs;
836  printer->stream(Debug) << "Number of converged eigenpairs (after) = " << sol.numVecs << std::endl;
837  printer->stream(Debug) << "whichVecs[sol.numVecs-1] > sol.numVecs-1 : " << whichVecs[sol.numVecs-1] << " > " << sol.numVecs-1 << std::endl;
838  if (whichVecs[sol.numVecs-1] > (sol.numVecs-1)) {
839  keepMore = true;
840  numEvecs = whichVecs[sol.numVecs-1]+1; // Add 1 to fix zero-based indexing
841  printer->stream(Debug) << "keepMore = true; numEvecs = " << numEvecs << std::endl;
842  }
843 
844  // Next set the number of Ritz vectors that the iteration must compute and compute them.
845  bks_solver->setNumRitzVectors(numEvecs);
846  bks_solver->computeRitzVectors();
847 
848  // If the leading Ritz pairs are the converged ones, get the information
849  // from the iteration to the solution container. Otherwise copy the necessary
850  // information using 'whichVecs'.
851  if (!keepMore) {
852  sol.index = bks_solver->getRitzIndex();
853  sol.Evals = bks_solver->getRitzValues();
854  sol.Evecs = MVT::CloneCopy( *(bks_solver->getRitzVectors()) );
855  }
856 
857  // Resize based on the number of solutions being returned and set the number of Ritz
858  // vectors for the iteration to compute.
859  sol.Evals.resize(sol.numVecs);
860  sol.index.resize(sol.numVecs);
861 
862  // If the converged Ritz pairs are not the leading ones, copy over the information directly.
863  if (keepMore) {
864  std::vector<Anasazi::Value<ScalarType> > tmpEvals = bks_solver->getRitzValues();
865  for (int vec_i=0; vec_i<sol.numVecs; ++vec_i) {
866  sol.index[vec_i] = tmpIndex[whichVecs[vec_i]];
867  sol.Evals[vec_i] = tmpEvals[whichVecs[vec_i]];
868  }
869  sol.Evecs = MVT::CloneCopy( *(bks_solver->getRitzVectors()), whichVecs );
870  }
871 
872  // Set the solution space to be the Ritz vectors at this time.
873  sol.Espace = sol.Evecs;
874  }
875  }
876 
877  // print final summary
878  bks_solver->currentStatus(printer->stream(FinalSummary));
879 
880  // print timing information
881 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
882  if ( printer->isVerbosity( TimingDetails ) ) {
883  Teuchos::TimeMonitor::summarize( printer->stream( TimingDetails ) );
884  }
885 #endif
886 
887  problem_->setSolution(sol);
888  printer->stream(Debug) << "Returning " << sol.numVecs << " eigenpairs to eigenproblem." << std::endl;
889 
890  // get the number of iterations performed during this solve.
891  numIters_ = bks_solver->getNumIters();
892 
893  if (sol.numVecs < nev) {
894  return Unconverged; // return from BlockKrylovSchurSolMgr::solve()
895  }
896  return Converged; // return from BlockKrylovSchurSolMgr::solve()
897 }
898 
899 
900 template <class ScalarType, class MV, class OP>
901 void
904 {
905  globalTest_ = global;
906 }
907 
908 template <class ScalarType, class MV, class OP>
911 {
912  return globalTest_;
913 }
914 
915 template <class ScalarType, class MV, class OP>
916 void
919 {
920  debugTest_ = debug;
921 }
922 
923 template <class ScalarType, class MV, class OP>
926 {
927  return debugTest_;
928 }
929 
930 } // end Anasazi namespace
931 
932 #endif /* ANASAZI_BLOCK_KRYLOV_SCHUR_SOLMGR_HPP */
ScalarType * values() const
static ptrdiff_t GetGlobalLength(const MV &mv)
Return the number of rows in the given multivector mv.
Pure virtual base class which describes the basic interface for a solver manager. ...
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > H
The current Hessenberg matrix.
std::vector< Value< ScalarType > > Evals
The computed eigenvalues.
A special StatusTest for printing other status tests.
This class defines the interface required by an eigensolver and status test class to compute solution...
An implementation of the Anasazi::SortManager that performs a collection of common sorting techniques...
T & get(const std::string &name, T def_value)
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Virtual base class which defines basic traits for the operator type.
Teuchos::RCP< MV > Evecs
The computed eigenvectors.
ReturnType solve()
This method performs possibly repeated calls to the underlying eigensolver&#39;s iterate() routine until ...
Status test for forming logical combinations of other status tests.
The Anasazi::SolverManager is a templated virtual base class that defines the basic interface that an...
void setDebugStatusTest(const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &debug)
Set the status test for debugging.
Basic implementation of the Anasazi::SortManager class.
An implementation of the Anasazi::MatOrthoManager that performs orthogonalization using the SVQB iter...
Structure to contain pointers to BlockKrylovSchur state variables.
void GEMM(ETransp transa, ETransp transb, const OrdinalType &m, const OrdinalType &n, const OrdinalType &k, const alpha_type alpha, const A_type *A, const OrdinalType &lda, const B_type *B, const OrdinalType &ldb, const beta_type beta, ScalarType *C, const OrdinalType &ldc) const
An exception class parent to all Anasazi exceptions.
Output managers remove the need for the eigensolver to know any information about the required output...
ScalarTraits< ScalarType >::magnitudeType normFrobenius() const
The Anasazi::BlockKrylovSchurSolMgr provides a flexible solver manager over the BlockKrylovSchur eige...
A status test for testing the norm of the eigenvectors residuals along with a set of auxiliary eigenv...
Anasazi&#39;s templated, static class providing utilities for the solvers.
bool isParameter(const std::string &name) const
int numVecs
The number of computed eigenpairs.
std::vector< Value< ScalarType > > getRitzValues() const
Return the Ritz values from the most recent solve.
void setGlobalStatusTest(const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &global)
Set the status test defining global convergence.
BlockKrylovSchurSolMgr(const Teuchos::RCP< Eigenproblem< ScalarType, MV, OP > > &problem, Teuchos::ParameterList &pl)
Basic constructor for BlockKrylovSchurSolMgr.
Teuchos::RCP< MV > Espace
An orthonormal basis for the computed eigenspace.
Abstract class definition for Anasazi Output Managers.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Abstract base class which defines the interface required by an eigensolver and status test class to c...
static void summarize(Ptr< const Comm< int > > comm, std::ostream &out=std::cout, const bool alwaysWriteLocal=false, const bool writeGlobalStats=true, const bool writeZeroTimers=true, const ECounterSetOp setOp=Intersection, const std::string &filter="", const bool ignoreZeroTimers=false)
static magnitudeType prec()
ReturnType
Enumerated type used to pass back information from a solver manager.
Output managers remove the need for the eigensolver to know any information about the required output...
A status test for testing the norm of the eigenvectors residuals.
Traits class which defines basic operations on multivectors.
std::vector< int > index
An index into Evecs to allow compressed storage of eigenvectors for real, non-Hermitian problems...
const Eigenproblem< ScalarType, MV, OP > & getProblem() const
Return the eigenvalue problem.
Teuchos::Array< Teuchos::RCP< Teuchos::Time > > getTimers() const
Return the timers for this object.
Anasazi header file which uses auto-configuration information to include necessary C++ headers...
Orthogonalization manager based on the SVQB technique described in &quot;A Block Orthogonalization Procedu...
Struct for storing an eigenproblem solution.
const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > & getGlobalStatusTest() const
Get the status test defining global convergence.
OrdinalType numCols() const
void push_back(const value_type &x)
void GEQRF(const OrdinalType &m, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, ScalarType *TAU, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *info) const
Teuchos::RCP< const MulVec > V
The current Krylov basis.
Special StatusTest for printing status tests.
A status test for testing the norm of the eigenvectors residuals along with a set of auxiliary eigenv...
Status test for forming logical combinations of other status tests.
An implementation of the Anasazi::MatOrthoManager that performs orthogonalization using (potentially)...
Types and exceptions used within Anasazi solvers and interfaces.
Abstract class definition for Anasazi output stream.
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > Q
The current Schur vectors of the valid part of H.
Implementation of a block Krylov-Schur eigensolver.
This class implements the block Krylov-Schur iteration, for solving linear eigenvalue problems...
int curDim
The current dimension of the reduction.
Common interface of stopping criteria for Anasazi&#39;s solvers.
A status test for testing the norm of the eigenvectors residuals.
Basic implementation of the Anasazi::OrthoManager class.
const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > & getDebugStatusTest() const
Get the status test for debugging.
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > S
The current Schur form reduction of the valid part of H.
OrdinalType stride() const
int getNumIters() const
Get the iteration count for the most recent call to solve().
OrdinalType numRows() const
Class which provides internal utilities for the Anasazi solvers.