NOX
Development
|
Bordered system solver strategy based on Householder transformations. More...
#include <LOCA_BorderedSolver_EpetraHouseholder.H>
Protected Types | |
enum | PRECONDITIONER_METHOD { JACOBIAN, SMW } |
Enumerated type indicating preconditioner method. | |
Protected Attributes | |
Teuchos::RCP< LOCA::GlobalData > | globalData |
Global data object. | |
Teuchos::RCP < Teuchos::ParameterList > | solverParams |
Solver parameters. | |
Teuchos::RCP< LOCA::Epetra::Group > | grp |
Pointer to group storing J. | |
Teuchos::RCP< const LOCA::BorderedSolver::AbstractOperator > | op |
Teuchos::RCP< const NOX::Abstract::MultiVector > | A |
Pointer to A block. | |
Teuchos::RCP< const NOX::Abstract::MultiVector > | B |
Pointer to B block. | |
Teuchos::RCP< const NOX::Abstract::MultiVector::DenseMatrix > | C |
Pointer to C block. | |
Teuchos::RCP< const LOCA::MultiContinuation::ConstraintInterfaceMVDX > | constraints |
Pointer to constraint interface. | |
LOCA::BorderedSolver::HouseholderQR | qrFact |
QR Factorization object. | |
Teuchos::RCP < NOX::Abstract::MultiVector > | house_x |
Solution component of Householder multivec. | |
NOX::Abstract::MultiVector::DenseMatrix | house_p |
Parameter component of Householder multivec. | |
NOX::Abstract::MultiVector::DenseMatrix | T |
T matrix in compact WY representation. | |
NOX::Abstract::MultiVector::DenseMatrix | R |
R matrix in QR factorization. | |
Teuchos::RCP < NOX::Abstract::MultiVector > | U |
U matrix in low-rank update form P = J + U*V^T. | |
Teuchos::RCP < NOX::Abstract::MultiVector > | V |
V matrix in low-rank update form P = J + U*V^T. | |
Teuchos::RCP < NOX::Abstract::MultiVector > | house_x_trans |
Solution component of Householder multivec for transposed system. | |
NOX::Abstract::MultiVector::DenseMatrix | house_p_trans |
Parameter component of Householder multivec for transposed system. | |
NOX::Abstract::MultiVector::DenseMatrix | T_trans |
T matrix in compact WY representation for transposed system. | |
NOX::Abstract::MultiVector::DenseMatrix | R_trans |
R matrix in QR factorization for transposed system. | |
Teuchos::RCP < NOX::Abstract::MultiVector > | U_trans |
U matrix in low-rank update form P = J + U*V^T for transposed system. | |
Teuchos::RCP < NOX::Abstract::MultiVector > | V_trans |
V matrix in low-rank update form P = J + U*V^T for transposed system. | |
Teuchos::RCP< const NOX::Abstract::MultiVector > | Ablock |
Pointer to A block as an Epetra multivector. | |
Teuchos::RCP< const NOX::Abstract::MultiVector > | Bblock |
Pointer to B block as an Epetra multivector. | |
Teuchos::RCP < NOX::Abstract::MultiVector > | Ascaled |
Pointer to scaled A block. | |
Teuchos::RCP < NOX::Abstract::MultiVector > | Bscaled |
Pointer to scaled B block. | |
Teuchos::RCP < NOX::Abstract::MultiVector::DenseMatrix > | Cscaled |
Pointer to scaled C block. | |
Teuchos::RCP < NOX::Epetra::LinearSystem > | linSys |
Pointer to linear system. | |
Teuchos::RCP< Epetra_Operator > | epetraOp |
Pointer to Epetra operator. | |
Teuchos::RCP< const Epetra_BlockMap > | baseMap |
Pointer to base map for block vectors. | |
Teuchos::RCP< const Epetra_BlockMap > | globalMap |
Pointer to global map for block vectors. | |
int | numConstraints |
Number of constraint equations. | |
bool | isZeroA |
flag indicating whether A block is zero | |
bool | isZeroB |
flag indicating whether B block is zero | |
bool | isZeroC |
flag indicating whether C block is zero | |
bool | isValidForSolve |
Flag indicating whether constraint factorization for solve has been computed. | |
bool | isValidForTransposeSolve |
Flag indicating whether constraint factorization for transpoe solve has been computed. | |
Teuchos::BLAS< int, double > | dblas |
BLAS Wrappers. | |
bool | scale_rows |
Whether we should scale augmented rows to have unit 2-norm. | |
std::vector< double > | scale_vals |
Scale values for each row. | |
PRECONDITIONER_METHOD | precMethod |
Preconditioner method. | |
bool | includeUV |
Flag indicating whether to include U*V^T terms in preconditioner. | |
bool | use_P_For_Prec |
Flag indicating whether to use P = J + U*V^T in preconditioner. | |
bool | isComplex |
Flag indicating whether we are doing a complex solve. | |
double | omega |
Frequency for complex systems. | |
Bordered system solver strategy based on Householder transformations.
This class solves the extended system of equations
using Householder tranformations. The algorithm works as follows: First consider a slightly rearranged version of the extended system of equations:
Let
be the QR decomposition of the constraints matrix where and . Define
then the extended system of equations is equivalent to
and hence
This last equation equation can be written
where is given by
and
We then recover and by
It can be further shown that the operator above can be written
where , and . The equation is solved using an iterative solver using the definition of above, in this case AztecOO. The system is preconditioned using the preconditioner for . The operator is generated using the standard Householder QR algorithm (Algorithm 5.2.1, G. Golub and C. Van Loan, "Matrix Computations," 3rd Edition, Johns Hopkins, Baltimore, 1996) and is stored using the compact WY representation: (see R. Schreiber and C. Van Loan, "A Storage-Efficient WY Representation for Products of Householder Transformations," SIAM J. Sci. Stat. Comput., Vol. 10, No. 1, pp. 53-57, January 1989).
The operator representing is encapsulated in the class LOCA::Epetra::LowRankUpdateRowMatrix if is an Epetra_RowMatrix and LOCA::Epetra::LowRankUpdateOp otherwise. If the row matrix version is available can be scaled and also used to construct a preconditioner. If "Include UV In Preconditioner" is true as discussed below, the and terms will be included when computing this preconditioner, which can help stability when is nearly singular.
The class is intialized via the solverParams
parameter list argument to the constructor. The parameters this class recognizes are:
LOCA::BorderedSolver::EpetraHouseholder::EpetraHouseholder | ( | const Teuchos::RCP< LOCA::GlobalData > & | global_data, |
const Teuchos::RCP< LOCA::Parameter::SublistParser > & | topParams, | ||
const Teuchos::RCP< Teuchos::ParameterList > & | solverParams | ||
) |
Constructor.
global_data | [in] Global data object |
topParams | [in] Parsed top-level parameter list |
solverParams | [in] Bordered solver parameters as described above |
References Teuchos::ParameterList::get(), globalData, includeUV, LOCA::GlobalData::locaErrorCheck, precMethod, scale_rows, solverParams, and use_P_For_Prec.
|
virtual |
Computed extended matrix-multivector product.
Computes
Implements LOCA::BorderedSolver::AbstractStrategy.
References NOX::Abstract::Group::Failed, Teuchos::SerialDenseMatrix< OrdinalType, ScalarType >::multiply(), Teuchos::NO_TRANS, and NOX::Abstract::MultiVector::update().
|
virtual |
Solves the extended system using the technique described above.
The params argument is the linear solver parameters. If isZeroF or isZeroG is true, than the corresponding F or G pointers may be NULL.
Note that if either the A or B blocks are zero, the system is solved using a simple block elimination scheme instead of the Householder scheme.
Implements LOCA::BorderedSolver::AbstractStrategy.
References Teuchos::SerialDenseMatrix< OrdinalType, ScalarType >::assign(), Teuchos::ParameterList::get(), Teuchos::RCP< T >::get(), NOX::Abstract::MultiVector::init(), Teuchos::SerialDenseMatrix< OrdinalType, ScalarType >::numCols(), Teuchos::SerialDenseMatrix< OrdinalType, ScalarType >::numRows(), NOX::Abstract::Group::Ok, Teuchos::SerialDenseMatrix< OrdinalType, ScalarType >::putScalar(), Teuchos::rcp(), LOCA::BorderedSolver::UpperTriangularBlockElimination::solve(), and LOCA::BorderedSolver::LowerTriangularBlockElimination::solve().
|
virtual |
Solves the transpose of the extended system as defined above.
The params argument is the linear solver parameters.
Implements LOCA::BorderedSolver::AbstractStrategy.
References Teuchos::SerialDenseMatrix< OrdinalType, ScalarType >::assign(), Teuchos::ParameterList::get(), Teuchos::RCP< T >::get(), NOX::Abstract::MultiVector::init(), Teuchos::SerialDenseMatrix< OrdinalType, ScalarType >::numCols(), Teuchos::SerialDenseMatrix< OrdinalType, ScalarType >::numRows(), NOX::Abstract::Group::Ok, Teuchos::SerialDenseMatrix< OrdinalType, ScalarType >::putScalar(), Teuchos::rcp(), LOCA::BorderedSolver::UpperTriangularBlockElimination::solveTranspose(), and LOCA::BorderedSolver::LowerTriangularBlockElimination::solveTranspose().
|
virtual |
Computed extended matrix transpose-multivector product.
Computes
Implements LOCA::BorderedSolver::AbstractStrategy.
References NOX::Abstract::Group::Failed, Teuchos::SerialDenseMatrix< OrdinalType, ScalarType >::multiply(), NOX::Abstract::MultiVector::multiply(), Teuchos::NO_TRANS, and Teuchos::TRANS.
|
virtual |
Intialize solver for a solve.
This should be called after setMatrixBlocks(), but before applyInverse().
Implements LOCA::BorderedSolver::AbstractStrategy.
References NOX::Abstract::Group::Ok, and NOX::ShapeCopy.
|
virtual |
Intialize solver for a transpose solve.
This should be called after setMatrixBlocks(), but before applyInverseTranspose().
Implements LOCA::BorderedSolver::AbstractStrategy.
References NOX::Abstract::Group::Ok, and NOX::ShapeCopy.
|
virtual |
Set blocks.
The blockA
or blockC
pointer may be null if either is zero. Whether block B is zero will be determined by querying blockB
via ConstraintInterface::isConstraintDerivativesXZero.
Implements LOCA::BorderedSolver::AbstractStrategy.
References Teuchos::RCP< T >::get(), LOCA::Epetra::Group::getComplexLinearSystem(), LOCA::BorderedSolver::JacobianOperator::getGroup(), NOX::Epetra::Group::getLinearSystem(), Teuchos::SerialDenseMatrix< OrdinalType, ScalarType >::putScalar(), and Teuchos::rcp().