Belos Version of the Day
Loading...
Searching...
No Matches
BelosBlockFGmresIter.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Belos: Block Linear Solvers Package
4//
5// Copyright 2004-2016 NTESS and the Belos contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef BELOS_BLOCK_FGMRES_ITER_HPP
11#define BELOS_BLOCK_FGMRES_ITER_HPP
12
17#include "BelosConfigDefs.hpp"
18#include "BelosTypes.hpp"
20
24#include "BelosStatusTest.hpp"
27
28#include "Teuchos_BLAS.hpp"
29#include "Teuchos_SerialDenseMatrix.hpp"
30#include "Teuchos_SerialDenseVector.hpp"
31#include "Teuchos_ScalarTraits.hpp"
32#include "Teuchos_ParameterList.hpp"
33#include "Teuchos_TimeMonitor.hpp"
34
48namespace Belos {
49
50template<class ScalarType, class MV, class OP>
51class BlockFGmresIter : virtual public GmresIteration<ScalarType,MV,OP> {
52
53 public:
54
55 //
56 // Convenience typedefs
57 //
60 typedef Teuchos::ScalarTraits<ScalarType> SCT;
61 typedef typename SCT::magnitudeType MagnitudeType;
62
64
65
76 const Teuchos::RCP<OutputManager<ScalarType> > &printer,
77 const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > &tester,
78 const Teuchos::RCP<MatOrthoManager<ScalarType,MV,OP> > &ortho,
79 Teuchos::ParameterList &params );
80
82 virtual ~BlockFGmresIter() {};
84
85
87
88
110 void iterate();
111
134
143
153 state.curDim = curDim_;
154 state.V = V_;
155 state.Z = Z_;
156 state.H = H_;
157 state.R = R_;
158 state.z = z_;
159 return state;
160 }
161
163
164
166
167
169 int getNumIters() const { return iter_; }
170
172 void resetNumIters( int iter = 0 ) { iter_ = iter; }
173
176 Teuchos::RCP<const MV> getNativeResiduals( std::vector<MagnitudeType> *norms ) const;
177
179
184 Teuchos::RCP<MV> getCurrentUpdate() const;
185
187
190 void updateLSQR( int dim = -1 );
191
193 int getCurSubspaceDim() const {
194 if (!initialized_) return 0;
195 return curDim_;
196 };
197
199 int getMaxSubspaceDim() const { return blockSize_*numBlocks_; }
200
202
203
205
206
208 const LinearProblem<ScalarType,MV,OP>& getProblem() const { return *lp_; }
209
211 int getBlockSize() const { return blockSize_; }
212
214 void setBlockSize(int blockSize) { setSize( blockSize, numBlocks_ ); }
215
217 int getNumBlocks() const { return numBlocks_; }
218
220 void setNumBlocks(int numBlocks) { setSize( blockSize_, numBlocks ); }
221
228 void setSize(int blockSize, int numBlocks);
229
231 bool isInitialized() { return initialized_; }
232
234
235 private:
236
237 //
238 // Internal methods
239 //
241 void setStateSize();
242
243 //
244 // Classes inputed through constructor that define the linear problem to be solved.
245 //
246 const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > lp_;
247 const Teuchos::RCP<OutputManager<ScalarType> > om_;
248 const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > stest_;
249 const Teuchos::RCP<OrthoManager<ScalarType,MV> > ortho_;
250
251 //
252 // Algorithmic parameters
253 //
254 // blockSize_ is the solver block size.
255 // It controls the number of vectors added to the basis on each iteration.
256 int blockSize_;
257 // numBlocks_ is the size of the allocated space for the Krylov basis, in blocks.
258 int numBlocks_;
259
260 // Storage for QR factorization of the least squares system.
261 Teuchos::SerialDenseVector<int,ScalarType> beta, sn;
262 Teuchos::SerialDenseVector<int,MagnitudeType> cs;
263
264 //
265 // Current solver state
266 //
267 // initialized_ specifies that the basis vectors have been initialized and the iterate() routine
268 // is capable of running; _initialize is controlled by the initialize() member method
269 // For the implications of the state of initialized_, please see documentation for initialize()
270 bool initialized_;
271
272 // stateStorageInitialized_ specified that the state storage has be initialized to the current
273 // blockSize_ and numBlocks_. This initialization may be postponed if the linear problem was
274 // generated without the right-hand side or solution vectors.
275 bool stateStorageInitialized_;
276
277 // keepHessenberg_ specifies that the upper Hessenberg matrix should be stored separately
278 // from the QR-factored least squares system (R_). When false, H_ and R_ point to the
279 // same object and only the QR-rotated form is available via getState().
280 bool keepHessenberg_;
281
282 // Current subspace dimension, and number of iterations performed.
283 int curDim_, iter_;
284
285 //
286 // State Storage
287 //
288 Teuchos::RCP<MV> V_;
289 Teuchos::RCP<MV> Z_;
290 //
291 // Projected matrices
292 // H_ : Projected matrix from the Krylov factorization AV = VH + FE^T
293 //
294 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > H_;
295 //
296 // QR decomposition of Projected matrices for solving the least squares system HY = B.
297 // R_: Upper triangular reduction of H
298 // z_: Q applied to right-hand side of the least squares system
299 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > R_;
300 Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > z_;
301};
302
304 // Constructor.
305 template<class ScalarType, class MV, class OP>
308 const Teuchos::RCP<OutputManager<ScalarType> > &printer,
309 const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > &tester,
310 const Teuchos::RCP<MatOrthoManager<ScalarType,MV,OP> > &ortho,
311 Teuchos::ParameterList &params ):
312 lp_(problem),
313 om_(printer),
314 stest_(tester),
315 ortho_(ortho),
316 blockSize_(0),
317 numBlocks_(0),
318 initialized_(false),
319 stateStorageInitialized_(false),
320 keepHessenberg_(false),
321 curDim_(0),
322 iter_(0)
323 {
324 // Find out whether we are saving the Hessenberg matrix separately from R.
325 if (om_->isVerbosity(Debug))
326 keepHessenberg_ = true;
327 else
328 keepHessenberg_ = params.get("Keep Hessenberg", false);
329
330 // Get the maximum number of blocks allowed for this Krylov subspace
332 ! params.isParameter ("Num Blocks"), std::invalid_argument,
333 "Belos::BlockFGmresIter::constructor: mandatory parameter 'Num Blocks' is not specified.");
334 const int nb = params.get<int> ("Num Blocks");
335
336 // Set the block size and allocate data.
337 const int bs = params.get ("Block Size", 1);
338 setSize (bs, nb);
339 }
340
342 // Set the block size and make necessary adjustments.
343 template <class ScalarType, class MV, class OP>
345 {
346 // This routine only allocates space; it doesn't not perform any computation
347 // any change in size will invalidate the state of the solver.
348
349 TEUCHOS_TEST_FOR_EXCEPTION(numBlocks <= 0 || blockSize <= 0, std::invalid_argument, "Belos::BlockFGmresIter::setSize was passed a non-positive argument.");
350 if (blockSize == blockSize_ && numBlocks == numBlocks_) {
351 // do nothing
352 return;
353 }
354
355 if (blockSize!=blockSize_ || numBlocks!=numBlocks_)
356 stateStorageInitialized_ = false;
357
358 blockSize_ = blockSize;
359 numBlocks_ = numBlocks;
360
361 initialized_ = false;
362 curDim_ = 0;
363
364 // Use the current blockSize_ and numBlocks_ to initialize the state storage.
365 setStateSize();
366
367 }
368
370 // Setup the state storage.
371 template <class ScalarType, class MV, class OP>
373 {
374 using Teuchos::RCP;
375 using Teuchos::rcp;
376 typedef Teuchos::SerialDenseMatrix<int, ScalarType> SDM;
377
378 if (! stateStorageInitialized_) {
379 // Check if there is any multivector to clone from.
380 RCP<const MV> lhsMV = lp_->getLHS();
381 RCP<const MV> rhsMV = lp_->getRHS();
382 if (lhsMV == Teuchos::null && rhsMV == Teuchos::null) {
383 stateStorageInitialized_ = false;
384 return;
385 }
386 else {
388 // blockSize*numBlocks dependent
389 //
390 int newsd = blockSize_*(numBlocks_+1);
391
392 if (blockSize_==1) {
393 cs.resize (newsd);
394 sn.resize (newsd);
395 }
396 else {
397 beta.resize (newsd);
398 }
399
400 // Initialize the state storage
402 blockSize_ * static_cast<ptrdiff_t> (numBlocks_) > MVT::GetGlobalLength (*rhsMV),
403 std::invalid_argument, "Belos::BlockFGmresIter::setStateSize(): "
404 "Cannot generate a Krylov basis with dimension larger the operator!");
405
406 // If the subspace has not be initialized before, generate it using the LHS or RHS from lp_.
407 if (V_ == Teuchos::null) {
408 // Get the multivector that is not null.
409 RCP<const MV> tmp = (rhsMV != Teuchos::null) ? rhsMV : lhsMV;
411 tmp == Teuchos::null, std::invalid_argument,
412 "Belos::BlockFGmresIter::setStateSize(): "
413 "linear problem does not specify multivectors to clone from.");
414 V_ = MVT::Clone (*tmp, newsd);
415 }
416 else {
417 // Generate V_ by cloning itself ONLY if more space is needed.
418 if (MVT::GetNumberVecs (*V_) < newsd) {
419 RCP<const MV> tmp = V_;
420 V_ = MVT::Clone (*tmp, newsd);
421 }
422 }
423
424 if (Z_ == Teuchos::null) {
425 // Get the multivector that is not null.
426 RCP<const MV> tmp = (rhsMV != Teuchos::null) ? rhsMV : lhsMV;
428 tmp == Teuchos::null, std::invalid_argument,
429 "Belos::BlockFGmresIter::setStateSize(): "
430 "linear problem does not specify multivectors to clone from.");
431 Z_ = MVT::Clone (*tmp, newsd);
432 }
433 else {
434 // Generate Z_ by cloning itself ONLY if more space is needed.
435 if (MVT::GetNumberVecs (*Z_) < newsd) {
436 RCP<const MV> tmp = Z_;
437 Z_ = MVT::Clone (*tmp, newsd);
438 }
439 }
440
441 // R_ holds the QR-factored least squares system. Always allocated.
442 if (R_ == Teuchos::null) {
443 R_ = rcp (new SDM (newsd, newsd-blockSize_));
444 }
445 else {
446 R_->shapeUninitialized (newsd, newsd - blockSize_);
447 }
448
449 // H_ holds the raw (pre-QR) upper Hessenberg matrix.
450 // When keepHessenberg_ is false, H_ and R_ point to the same object
451 // (matching BlockGmresIter behavior).
452 if (keepHessenberg_) {
453 if (H_ == Teuchos::null) {
454 H_ = rcp (new SDM (newsd, newsd-blockSize_));
455 }
456 else {
457 H_->shapeUninitialized (newsd, newsd - blockSize_);
458 }
459 }
460 else {
461 H_ = R_;
462 }
463
464 // Generate z_ only if it doesn't exist, otherwise resize it.
465 if (z_ == Teuchos::null) {
466 z_ = rcp (new SDM (newsd, blockSize_));
467 }
468 else {
469 z_->shapeUninitialized (newsd, blockSize_);
470 }
471
472 // State storage has now been initialized.
473 stateStorageInitialized_ = true;
474 }
475 }
476 }
477
478
479 template <class ScalarType, class MV, class OP>
480 Teuchos::RCP<MV>
482 {
483 typedef Teuchos::SerialDenseMatrix<int, ScalarType> SDM;
484
485 Teuchos::RCP<MV> currentUpdate = Teuchos::null;
486 if (curDim_ == 0) {
487 // If this is the first iteration of the Arnoldi factorization,
488 // then there is no update, so return Teuchos::null.
489 return currentUpdate;
490 }
491 else {
492 const ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero ();
493 const ScalarType one = Teuchos::ScalarTraits<ScalarType>::one ();
494 Teuchos::BLAS<int,ScalarType> blas;
495
496 currentUpdate = MVT::Clone (*Z_, blockSize_);
497
498 // Make a view and then copy the RHS of the least squares problem. DON'T OVERWRITE IT!
499 SDM y (Teuchos::Copy, *z_, curDim_, blockSize_);
500
501 // Solve the least squares problem.
502 blas.TRSM (Teuchos::LEFT_SIDE, Teuchos::UPPER_TRI, Teuchos::NO_TRANS,
503 Teuchos::NON_UNIT_DIAG, curDim_, blockSize_, one,
504 R_->values (), R_->stride (), y.values (), y.stride ());
505
506 // Compute the current update.
507 std::vector<int> index (curDim_);
508 for (int i = 0; i < curDim_; ++i) {
509 index[i] = i;
510 }
511 Teuchos::RCP<const MV> Zjp1 = MVT::CloneView (*Z_, index);
512 MVT::MvTimesMatAddMv (one, *Zjp1, y, zero, *currentUpdate);
513 }
514 return currentUpdate;
515 }
516
517
518 template <class ScalarType, class MV, class OP>
519 Teuchos::RCP<const MV>
521 getNativeResiduals (std::vector<MagnitudeType> *norms) const
522 {
523 // NOTE: Make sure the incoming std::vector is the correct size!
524 if (norms != NULL && (int)norms->size() < blockSize_) {
525 norms->resize (blockSize_);
526 }
527
528 if (norms != NULL) {
529 Teuchos::BLAS<int, ScalarType> blas;
530 for (int j = 0; j < blockSize_; ++j) {
531 (*norms)[j] = blas.NRM2 (blockSize_, &(*z_)(curDim_, j), 1);
532 }
533 }
534
535 // FGmres does not return a residual (multi)vector.
536 return Teuchos::null;
537 }
538
539
540 template <class ScalarType, class MV, class OP>
543 {
544 using Teuchos::RCP;
545 using Teuchos::rcp;
546 using std::endl;
547 typedef Teuchos::ScalarTraits<ScalarType> STS;
548 typedef Teuchos::SerialDenseMatrix<int, ScalarType> SDM;
549 const ScalarType ZERO = STS::zero ();
550 const ScalarType ONE = STS::one ();
551
552 // Initialize the state storage if it isn't already.
553 if (! stateStorageInitialized_) {
554 setStateSize ();
555 }
556
558 ! stateStorageInitialized_, std::invalid_argument,
559 "Belos::BlockFGmresIter::initialize(): Cannot initialize state storage!");
560
561 // NOTE: In BlockFGmresIter, V and Z are required!!! Inconsistent
562 // multivectors widths and lengths will not be tolerated, and will
563 // be treated with exceptions.
564 const char errstr[] = "Belos::BlockFGmresIter::initialize(): The given "
565 "multivectors must have a consistent length and width.";
566
567 if (! newstate.V.is_null () && ! newstate.z.is_null ()) {
568
569 // initialize V_,z_, and curDim_
570
572 MVT::GetGlobalLength(*newstate.V) != MVT::GetGlobalLength(*V_),
573 std::invalid_argument, errstr );
575 MVT::GetNumberVecs(*newstate.V) < blockSize_,
576 std::invalid_argument, errstr );
578 newstate.curDim > blockSize_*(numBlocks_+1),
579 std::invalid_argument, errstr );
580
581 curDim_ = newstate.curDim;
582 const int lclDim = MVT::GetNumberVecs(*newstate.V);
583
584 // check size of Z
586 newstate.z->numRows() < curDim_ || newstate.z->numCols() < blockSize_,
587 std::invalid_argument, errstr);
588
589 // copy basis vectors from newstate into V
590 if (newstate.V != V_) {
591 // only copy over the first block and print a warning.
592 if (curDim_ == 0 && lclDim > blockSize_) {
593 std::ostream& warn = om_->stream (Warnings);
594 warn << "Belos::BlockFGmresIter::initialize(): the solver was "
595 << "initialized with a kernel of " << lclDim << endl
596 << "The block size however is only " << blockSize_ << endl
597 << "The last " << lclDim - blockSize_
598 << " vectors will be discarded." << endl;
599 }
600 std::vector<int> nevind (curDim_ + blockSize_);
601 for (int i = 0; i < curDim_ + blockSize_; ++i) {
602 nevind[i] = i;
603 }
604 RCP<const MV> newV = MVT::CloneView (*newstate.V, nevind);
605 RCP<MV> lclV = MVT::CloneViewNonConst (*V_, nevind);
606 MVT::MvAddMv (ONE, *newV, ZERO, *newV, *lclV);
607
608 // done with local pointers
609 lclV = Teuchos::null;
610 }
611
612 // put data into z_, make sure old information is not still hanging around.
613 if (newstate.z != z_) {
614 z_->putScalar();
615 SDM newZ (Teuchos::View, *newstate.z, curDim_ + blockSize_, blockSize_);
617 lclz = rcp (new SDM (Teuchos::View, *z_, curDim_ + blockSize_, blockSize_));
618 lclz->assign (newZ);
619 lclz = Teuchos::null; // done with local pointers
620 }
621 }
622 else {
624 newstate.V == Teuchos::null,std::invalid_argument,
625 "Belos::BlockFGmresIter::initialize(): BlockFGmresStateIterState does not have initial kernel V_0.");
626
628 newstate.z == Teuchos::null,std::invalid_argument,
629 "Belos::BlockFGmresIter::initialize(): BlockFGmresStateIterState does not have initial norms z_0.");
630 }
631
632 // the solver is initialized
633 initialized_ = true;
634 }
635
636
637 template <class ScalarType, class MV, class OP>
639 {
640 using Teuchos::Array;
641 using Teuchos::null;
642 using Teuchos::RCP;
643 using Teuchos::rcp;
644 using Teuchos::View;
645 typedef Teuchos::SerialDenseMatrix<int, ScalarType> SDM;
646
647 // Allocate/initialize data structures
648 if (initialized_ == false) {
649 initialize();
650 }
651
652 // Compute the current search dimension.
653 const int searchDim = blockSize_ * numBlocks_;
654
655 // Iterate until the status test tells us to stop.
656 // Raise an exception if a computed block is not full rank.
657 while (stest_->checkStatus (this) != Passed && curDim_+blockSize_ <= searchDim) {
658 ++iter_;
659
660 // F can be found at the curDim_ block, but the next block is at curDim_ + blockSize_.
661 const int lclDim = curDim_ + blockSize_;
662
663 // Get the current part of the basis.
664 std::vector<int> curind (blockSize_);
665 for (int i = 0; i < blockSize_; ++i) {
666 curind[i] = lclDim + i;
667 }
668 RCP<MV> Vnext = MVT::CloneViewNonConst (*V_, curind);
669
670 // Get a view of the previous vectors.
671 // This is used for orthogonalization and for computing V^H K H.
672 for (int i = 0; i < blockSize_; ++i) {
673 curind[i] = curDim_ + i;
674 }
675 RCP<const MV> Vprev = MVT::CloneView (*V_, curind);
676 RCP<MV> Znext = MVT::CloneViewNonConst (*Z_, curind);
677
678 // Compute the next (multi)vector in the Krylov basis: Znext = M*Vprev
679 lp_->applyRightPrec (*Vprev, *Znext);
680 Vprev = null;
681
682 // Compute the next (multi)vector in the Krylov basis: Vnext = A*Znext
683 lp_->applyOp (*Znext, *Vnext);
684 Znext = null;
685
686 // Remove all previous Krylov basis vectors from Vnext
687 // Get a view of all the previous vectors
688 std::vector<int> prevind (lclDim);
689 for (int i = 0; i < lclDim; ++i) {
690 prevind[i] = i;
691 }
692 Vprev = MVT::CloneView (*V_, prevind);
694
695 // Get a view of the part of the Hessenberg matrix needed to hold the ortho coeffs.
696 // Ortho always writes into H_ (the raw Hessenberg).
697 RCP<SDM> subH = rcp (new SDM (View, *H_, lclDim, blockSize_, 0, curDim_));
699 AsubH.append (subH);
700
701 // Get a view of the part of the Hessenberg matrix needed to hold the norm coeffs.
702 RCP<SDM> subH2 = rcp (new SDM (View, *H_, blockSize_, blockSize_, lclDim, curDim_));
703 const int rank = ortho_->projectAndNormalize (*Vnext, AsubH, subH2, AVprev);
705 rank != blockSize_, GmresIterationOrthoFailure,
706 "Belos::BlockFGmresIter::iterate(): After orthogonalization, the new "
707 "basis block does not have full rank. It contains " << blockSize_
708 << " vector" << (blockSize_ != 1 ? "s" : "")
709 << ", but its rank is " << rank << ".");
710
711 // If keeping the Hessenberg separately, copy the new columns into R_
712 // before updateLSQR() overwrites them with the QR factorization.
713 if (keepHessenberg_) {
714 RCP<SDM> subR = rcp (new SDM (View, *R_, lclDim, blockSize_, 0, curDim_));
715 subR->assign (*subH);
716 RCP<SDM> subR2 = rcp (new SDM (View, *R_, blockSize_, blockSize_, lclDim, curDim_));
717 subR2->assign (*subH2);
718 }
719
720 //
721 // V has been extended, and H has been extended.
722 //
723 // Update the QR factorization of the upper Hessenberg matrix (applied to R_).
724 //
725 updateLSQR ();
726 //
727 // Update basis dim and release all pointers.
728 //
729 Vnext = null;
730 curDim_ += blockSize_;
731 } // end while (statusTest == false)
732 }
733
734
735 template<class ScalarType, class MV, class OP>
737 {
738 typedef Teuchos::ScalarTraits<ScalarType> STS;
739 typedef Teuchos::ScalarTraits<MagnitudeType> STM;
740
741 const ScalarType zero = STS::zero ();
742 const ScalarType two = (STS::one () + STS::one());
744 Teuchos::BLAS<int, ScalarType> blas;
745
746 // Get correct dimension based on input 'dim'. Remember that
747 // orthogonalization failures result in an exit before
748 // updateLSQR() is called. Therefore, it is possible that dim ==
749 // curDim_.
750 int curDim = curDim_;
751 if (dim >= curDim_ && dim < getMaxSubspaceDim ()) {
752 curDim = dim;
753 }
754
755 // Apply previous transformations, and compute new transformation
756 // to reduce upper Hessenberg system to upper triangular form.
757 // The type of transformation we use depends the block size. We
758 // use Givens rotations for a block size of 1, and Householder
759 // reflectors otherwise.
760 if (blockSize_ == 1) {
761 // QR factorization of upper Hessenberg matrix using Givens rotations
762 for (int i = 0; i < curDim; ++i) {
763 // Apply previous Givens rotations to new column of Hessenberg matrix
764 blas.ROT (1, &(*R_)(i, curDim), 1, &(*R_)(i+1, curDim), 1, &cs[i], &sn[i]);
765 }
766 // Calculate new Givens rotation
767 blas.ROTG (&(*R_)(curDim, curDim), &(*R_)(curDim+1, curDim), &cs[curDim], &sn[curDim]);
768 (*R_)(curDim+1, curDim) = zero;
769
770 // Update RHS w/ new transformation
771 blas.ROT (1, &(*z_)(curDim,0), 1, &(*z_)(curDim+1,0), 1, &cs[curDim], &sn[curDim]);
772 }
773 else {
774 // QR factorization of least-squares system using Householder reflectors.
775 for (int j = 0; j < blockSize_; ++j) {
776 // Apply previous Householder reflectors to new block of Hessenberg matrix
777 for (int i = 0; i < curDim + j; ++i) {
778 sigma = blas.DOT (blockSize_, &(*R_)(i+1,i), 1, &(*R_)(i+1,curDim+j), 1);
779 sigma += (*R_)(i,curDim+j);
780 sigma *= beta[i];
781 blas.AXPY (blockSize_, ScalarType(-sigma), &(*R_)(i+1,i), 1, &(*R_)(i+1,curDim+j), 1);
782 (*R_)(i,curDim+j) -= sigma;
783 }
784
785 // Compute new Householder reflector
786 const int maxidx = blas.IAMAX (blockSize_+1, &(*R_)(curDim+j,curDim+j), 1);
787 maxelem = (*R_)(curDim + j + maxidx - 1, curDim + j);
788 for (int i = 0; i < blockSize_ + 1; ++i) {
789 (*R_)(curDim+j+i,curDim+j) /= maxelem;
790 }
791 sigma = blas.DOT (blockSize_, &(*R_)(curDim + j + 1, curDim + j), 1,
792 &(*R_)(curDim + j + 1, curDim + j), 1);
793 if (sigma == zero) {
794 beta[curDim + j] = zero;
795 } else {
796 mu = STS::squareroot ((*R_)(curDim+j,curDim+j)*(*R_)(curDim+j,curDim+j)+sigma);
797 if (STS::real ((*R_)(curDim + j, curDim + j)) < STM::zero ()) {
798 vscale = (*R_)(curDim+j,curDim+j) - mu;
799 } else {
800 vscale = -sigma / ((*R_)(curDim+j, curDim+j) + mu);
801 }
802 beta[curDim+j] = two * vscale * vscale / (sigma + vscale*vscale);
803 (*R_)(curDim+j, curDim+j) = maxelem*mu;
804 for (int i = 0; i < blockSize_; ++i) {
805 (*R_)(curDim+j+1+i,curDim+j) /= vscale;
806 }
807 }
808
809 // Apply new Householder reflector to the right-hand side.
810 for (int i = 0; i < blockSize_; ++i) {
811 sigma = blas.DOT (blockSize_, &(*R_)(curDim+j+1,curDim+j),
812 1, &(*z_)(curDim+j+1,i), 1);
813 sigma += (*z_)(curDim+j,i);
814 sigma *= beta[curDim+j];
815 blas.AXPY (blockSize_, ScalarType(-sigma), &(*R_)(curDim+j+1,curDim+j),
816 1, &(*z_)(curDim+j+1,i), 1);
817 (*z_)(curDim+j,i) -= sigma;
818 }
819 }
820 } // end if (blockSize_ == 1)
821
822 // If the least-squares problem is updated wrt "dim" then update curDim_.
823 if (dim >= curDim_ && dim < getMaxSubspaceDim ()) {
824 curDim_ = dim + blockSize_;
825 }
826 } // end updateLSQR()
827
828} // namespace Belos
829
830#endif /* BELOS_BLOCK_FGMRES_ITER_HPP */
Belos header file which uses auto-configuration information to include necessary C++ headers.
Pure virtual base class which augments the basic interface for a Gmres linear solver iteration.
Class which describes the linear problem to be solved by the iterative solver.
Templated virtual class for providing orthogonalization/orthonormalization methods with matrix-based ...
Declaration of basic traits for the multivector type.
Class which defines basic traits for the operator type.
Class which manages the output and verbosity of the Belos solvers.
Pure virtual base class for defining the status testing capabilities of Belos.
Collection of types and exceptions used within the Belos solvers.
This class implements the block flexible GMRES iteration, where a block Krylov subspace is constructe...
int getNumIters() const
Get the current iteration count.
int getBlockSize() const
Get the blocksize to be used by the iterative solver in solving this linear problem.
OperatorTraits< ScalarType, MV, OP > OPT
Teuchos::RCP< const MV > getNativeResiduals(std::vector< MagnitudeType > *norms) const
Get the norms of the residuals native to the solver.
void iterate()
This method performs block FGmres iterations until the status test indicates the need to stop or an e...
MultiVecTraits< ScalarType, MV > MVT
GmresIterationState< ScalarType, MV > getState() const
Get the current state of the linear solver.
int getNumBlocks() const
Get the maximum number of blocks used by the iterative solver in solving this linear problem.
void initializeGmres(GmresIterationState< ScalarType, MV > &newstate)
Initialize the solver to an iterate, providing a complete state.
SCT::magnitudeType MagnitudeType
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the linear problem.
Teuchos::RCP< MV > getCurrentUpdate() const
Get the current update to the linear system.
void resetNumIters(int iter=0)
Reset the iteration count.
virtual ~BlockFGmresIter()
Destructor.
void updateLSQR(int dim=-1)
Method for updating QR factorization of upper Hessenberg matrix.
void setBlockSize(int blockSize)
Set the blocksize.
void initialize()
Initialize the solver with the initial vectors from the linear problem or random data.
BlockFGmresIter(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem, const Teuchos::RCP< OutputManager< ScalarType > > &printer, const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &tester, const Teuchos::RCP< MatOrthoManager< ScalarType, MV, OP > > &ortho, Teuchos::ParameterList &params)
BlockFGmresIter constructor with linear problem, solver utilities, and parameter list of solver optio...
int getCurSubspaceDim() const
Get the dimension of the search subspace used to generate the current solution to the linear problem.
int getMaxSubspaceDim() const
Get the maximum dimension allocated for the search subspace.
void setSize(int blockSize, int numBlocks)
Set the blocksize and number of blocks to be used by the iterative solver in solving this linear prob...
void setNumBlocks(int numBlocks)
Set the maximum number of blocks used by the iterative solver.
Teuchos::ScalarTraits< ScalarType > SCT
bool isInitialized()
States whether the solver has been initialized or not.
GmresIterationOrthoFailure is thrown when the GmresIteration object is unable to compute independent ...
Alternative run-time polymorphic interface for operators.
Operator()
Default constructor (does nothing).

Generated for Belos by doxygen 1.9.8