Stokhos Package Browser (Single Doxygen Collection)
Version of the Day
|
Nonlinear, stochastic Galerkin ModelEvaluator. More...
#include <Stokhos_SGModelEvaluator.hpp>
Public Member Functions | |
SGModelEvaluator (const Teuchos::RCP< EpetraExt::ModelEvaluator > &me, const Teuchos::RCP< const Stokhos::OrthogPolyBasis< int, double > > &sg_basis, const Teuchos::RCP< const Stokhos::Quadrature< int, double > > &sg_quad, const Teuchos::RCP< Stokhos::OrthogPolyExpansion< int, double > > &sg_exp, const Teuchos::RCP< const Stokhos::ParallelData > &sg_parallel_data, const Teuchos::RCP< Teuchos::ParameterList > ¶ms, bool scaleOP=true) | |
Teuchos::RCP < EpetraExt::BlockVector > | import_solution (const Epetra_Vector &x) const |
Import parallel solution vector. More... | |
Teuchos::RCP < Stokhos::EpetraVectorOrthogPoly > | import_solution_poly (const Epetra_Vector &x) const |
Import parallel solution vector. More... | |
Teuchos::RCP < EpetraExt::BlockVector > | export_solution (const Epetra_Vector &x_overlapped) const |
Export parallel solution vector. More... | |
Teuchos::RCP < EpetraExt::BlockVector > | import_residual (const Epetra_Vector &f) const |
Import parallel residual vector. More... | |
Teuchos::RCP < EpetraExt::BlockVector > | export_residual (const Epetra_Vector &f_overlapped) const |
Export parallel residual vector. More... | |
Public Member Functions inherited from Stokhos::SGModelEvaluatorBase | |
SGModelEvaluatorBase () | |
virtual | ~SGModelEvaluatorBase () |
Destructor. More... | |
Protected Attributes | |
Teuchos::RCP < EpetraExt::ModelEvaluator > | me |
Underlying model evaluator. More... | |
Teuchos::RCP< const Stokhos::OrthogPolyBasis< int, double > > | sg_basis |
Stochastic Galerkin basis. More... | |
Teuchos::RCP< const Stokhos::Quadrature< int, double > > | sg_quad |
Stochastic Galerkin quadrature. More... | |
Teuchos::RCP < Stokhos::OrthogPolyExpansion < int, double > > | sg_exp |
Stochastic Galerkin expansion. More... | |
Teuchos::RCP < Teuchos::ParameterList > | params |
Algorithmic parameters. More... | |
unsigned int | num_sg_blocks |
Number of stochastic blocks. More... | |
unsigned int | num_W_blocks |
Number of W stochastic blocks (may be smaller than num_sg_blocks) More... | |
unsigned int | num_p_blocks |
Number of p stochastic blocks (may be smaller than num_sg_blocks) More... | |
bool | supports_x |
Whether we support x (and thus f and W) More... | |
Teuchos::RCP< const Epetra_Map > | x_map |
Underlying unknown map. More... | |
Teuchos::RCP< const Epetra_Map > | f_map |
Underlying residual map. More... | |
Teuchos::RCP< const Stokhos::ParallelData > | sg_parallel_data |
Parallel SG data. More... | |
Teuchos::RCP< const EpetraExt::MultiComm > | sg_comm |
Parallel SG communicator. More... | |
Teuchos::RCP< const Stokhos::EpetraSparse3Tensor > | epetraCijk |
Epetra Cijk. More... | |
Teuchos::RCP< const Stokhos::EpetraSparse3Tensor > | serialCijk |
Serial Epetra Cijk for dgdx*. More... | |
Teuchos::RCP< const Epetra_BlockMap > | stoch_row_map |
Map for stochastic blocks. More... | |
Teuchos::RCP< const Epetra_BlockMap > | overlapped_stoch_row_map |
Overlapped map for stochastic blocks (local map) More... | |
Teuchos::RCP< const Epetra_BlockMap > | overlapped_stoch_p_map |
Overlapped map for p stochastic blocks (local map) More... | |
Teuchos::RCP< const Epetra_Map > | sg_x_map |
Block SG unknown map. More... | |
Teuchos::RCP< const Epetra_Map > | sg_overlapped_x_map |
Block SG overlapped unknown map. More... | |
Teuchos::RCP< const Epetra_Map > | sg_f_map |
Block SG residual map. More... | |
Teuchos::RCP< const Epetra_Map > | sg_overlapped_f_map |
Block SG overlapped residual map. More... | |
Teuchos::RCP< Epetra_Import > | sg_overlapped_x_importer |
Importer from SG to SG-overlapped maps. More... | |
Teuchos::RCP< Epetra_Export > | sg_overlapped_f_exporter |
Exporter from SG-overlapped to SG maps. More... | |
int | num_p |
Number of parameter vectors of underlying model evaluator. More... | |
int | num_p_sg |
Number of stochastic parameter vectors. More... | |
Teuchos::Array< int > | sg_p_index_map |
Index map between block-p and p_sg maps. More... | |
Teuchos::Array< Teuchos::RCP < const Epetra_Map > > | sg_p_map |
Block SG parameter map. More... | |
Teuchos::Array< Teuchos::RCP < Teuchos::Array< std::string > > > | sg_p_names |
SG coefficient parameter names. More... | |
int | num_g |
Number of response vectors of underlying model evaluator. More... | |
int | num_g_sg |
Number of stochastic response vectors. More... | |
Teuchos::Array< int > | sg_g_index_map |
Index map between block-g and g_sg maps. More... | |
Teuchos::Array< Teuchos::RCP < const Epetra_Map > > | sg_g_map |
Block SG response map. More... | |
Teuchos::RCP < Stokhos::EpetraVectorOrthogPoly > | x_dot_sg_blocks |
x_dot stochastic Galerkin components More... | |
Teuchos::RCP < Stokhos::EpetraVectorOrthogPoly > | x_sg_blocks |
x stochastic Galerkin components More... | |
Teuchos::RCP < Stokhos::EpetraVectorOrthogPoly > | f_sg_blocks |
f stochastic Galerkin components More... | |
Teuchos::RCP < Stokhos::EpetraOperatorOrthogPoly > | W_sg_blocks |
W stochastic Galerkin components. More... | |
Teuchos::RCP < Stokhos::EpetraVectorOrthogPoly > | sg_x_init |
SG initial x. More... | |
Teuchos::Array< Teuchos::RCP < Stokhos::EpetraVectorOrthogPoly > > | sg_p_init |
SG initial p. More... | |
bool | eval_W_with_f |
Whether to always evaluate W with f. More... | |
Teuchos::RCP< Stokhos::SGOperator > | my_W |
W pointer for evaluating W with f. More... | |
Teuchos::RCP< Epetra_Vector > | my_x |
x pointer for evaluating preconditioner More... | |
bool | scaleOP |
Teuchos::RCP < Stokhos::SGPreconditionerFactory > | sg_prec_factory |
Preconditioner factory. More... | |
Overridden from EpetraExt::ModelEvaluator . | |
Teuchos::RCP< const Epetra_Map > | get_x_map () const |
Return solution vector map. More... | |
Teuchos::RCP< const Epetra_Map > | get_f_map () const |
Return residual vector map. More... | |
Teuchos::RCP< const Epetra_Map > | get_p_map (int l) const |
Return parameter vector map. More... | |
Teuchos::RCP< const Epetra_Map > | get_g_map (int l) const |
Return response map. More... | |
Teuchos::RCP< const Teuchos::Array< std::string > > | get_p_names (int l) const |
Return array of parameter names. More... | |
Teuchos::RCP< const Epetra_Vector > | get_x_init () const |
Return initial solution. More... | |
Teuchos::RCP< const Epetra_Vector > | get_p_init (int l) const |
Return initial parameters. More... | |
Teuchos::RCP< Epetra_Operator > | create_W () const |
Create W = alpha*M + beta*J matrix. More... | |
Teuchos::RCP < EpetraExt::ModelEvaluator::Preconditioner > | create_WPrec () const |
Create preconditioner operator. More... | |
Teuchos::RCP< Epetra_Operator > | create_DgDx_dot_op (int j) const |
Create SG operator representing dg/dxdot. More... | |
Teuchos::RCP< Epetra_Operator > | create_DgDx_op (int j) const |
Create SG operator representing dg/dx. More... | |
Teuchos::RCP< Epetra_Operator > | create_DgDp_op (int j, int i) const |
Create SG operator representing dg/dp. More... | |
Teuchos::RCP< Epetra_Operator > | create_DfDp_op (int i) const |
Create SG operator representing df/dp. More... | |
InArgs | createInArgs () const |
Create InArgs. More... | |
OutArgs | createOutArgs () const |
Create OutArgs. More... | |
void | evalModel (const InArgs &inArgs, const OutArgs &outArgs) const |
Evaluate model on InArgs. More... | |
Nonlinear, stochastic Galerkin ModelEvaluator.
SGModelEvaluator is an implementation of EpetraExt::ModelEvaluator that generates a nonlinear problem from a stochastic Galerkin expansion. It wraps a supplied ModelEvaluator that supports the SG versions of p, x, and possibly x_dot InArgs, and f and W OutArgs, and translates those into a new nonlinear problem. It does so by concatenating all of the SG components of p, x, x_dot, and f into extended block vectors that form the parameters, solution vector, time derivative vector and residual for the new nonlinear problem. For dealing with the W matrix two methods are supported: forming a fully-assembled SG matrix and a "matrix free" method. The choice is selected by setting the "Jacobian Method" parameter of the parameter list supplied to the constructor, which can be either "Fully Assembled" or "Matrix Free". In the first case, the W operator of the underlying model evaluator must be an Epetra_CrsMatrix. In the second case, a preconditioner for the mean block must also be supplied via the "Preconditioner Factory" parameter of this list. This preconditioner factory must implement the Stokhos::PreconditionerFactory interface also supplied in this file. Currently using a preconditioner for the mean is the only option available for preconditioning the SG system when using the matrix-free method.
Definition at line 58 of file Stokhos_SGModelEvaluator.hpp.
Stokhos::SGModelEvaluator::SGModelEvaluator | ( | const Teuchos::RCP< EpetraExt::ModelEvaluator > & | me, |
const Teuchos::RCP< const Stokhos::OrthogPolyBasis< int, double > > & | sg_basis, | ||
const Teuchos::RCP< const Stokhos::Quadrature< int, double > > & | sg_quad, | ||
const Teuchos::RCP< Stokhos::OrthogPolyExpansion< int, double > > & | sg_exp, | ||
const Teuchos::RCP< const Stokhos::ParallelData > & | sg_parallel_data, | ||
const Teuchos::RCP< Teuchos::ParameterList > & | params, | ||
bool | scaleOP = true |
||
) |
Definition at line 24 of file Stokhos_SGModelEvaluator.cpp.
Teuchos::RCP< const Epetra_Map > Stokhos::SGModelEvaluator::get_x_map | ( | ) | const |
Return solution vector map.
Definition at line 232 of file Stokhos_SGModelEvaluator.cpp.
Teuchos::RCP< const Epetra_Map > Stokhos::SGModelEvaluator::get_f_map | ( | ) | const |
Return residual vector map.
Definition at line 238 of file Stokhos_SGModelEvaluator.cpp.
Teuchos::RCP< const Epetra_Map > Stokhos::SGModelEvaluator::get_p_map | ( | int | l | ) | const |
Return parameter vector map.
Definition at line 244 of file Stokhos_SGModelEvaluator.cpp.
Teuchos::RCP< const Epetra_Map > Stokhos::SGModelEvaluator::get_g_map | ( | int | l | ) | const |
Return response map.
Definition at line 257 of file Stokhos_SGModelEvaluator.cpp.
Teuchos::RCP< const Teuchos::Array< std::string > > Stokhos::SGModelEvaluator::get_p_names | ( | int | l | ) | const |
Return array of parameter names.
Definition at line 265 of file Stokhos_SGModelEvaluator.cpp.
Teuchos::RCP< const Epetra_Vector > Stokhos::SGModelEvaluator::get_x_init | ( | ) | const |
Return initial solution.
Definition at line 278 of file Stokhos_SGModelEvaluator.cpp.
Teuchos::RCP< const Epetra_Vector > Stokhos::SGModelEvaluator::get_p_init | ( | int | l | ) | const |
Return initial parameters.
Definition at line 284 of file Stokhos_SGModelEvaluator.cpp.
Teuchos::RCP< Epetra_Operator > Stokhos::SGModelEvaluator::create_W | ( | ) | const |
Create W = alpha*M + beta*J matrix.
Definition at line 297 of file Stokhos_SGModelEvaluator.cpp.
Teuchos::RCP< EpetraExt::ModelEvaluator::Preconditioner > Stokhos::SGModelEvaluator::create_WPrec | ( | ) | const |
Create preconditioner operator.
Definition at line 323 of file Stokhos_SGModelEvaluator.cpp.
Teuchos::RCP< Epetra_Operator > Stokhos::SGModelEvaluator::create_DgDx_dot_op | ( | int | j | ) | const |
Create SG operator representing dg/dxdot.
Definition at line 390 of file Stokhos_SGModelEvaluator.cpp.
Teuchos::RCP< Epetra_Operator > Stokhos::SGModelEvaluator::create_DgDx_op | ( | int | j | ) | const |
Create SG operator representing dg/dx.
Definition at line 337 of file Stokhos_SGModelEvaluator.cpp.
Teuchos::RCP< Epetra_Operator > Stokhos::SGModelEvaluator::create_DgDp_op | ( | int | j, |
int | i | ||
) | const |
Create SG operator representing dg/dp.
Definition at line 443 of file Stokhos_SGModelEvaluator.cpp.
Teuchos::RCP< Epetra_Operator > Stokhos::SGModelEvaluator::create_DfDp_op | ( | int | i | ) | const |
Create SG operator representing df/dp.
Definition at line 521 of file Stokhos_SGModelEvaluator.cpp.
EpetraExt::ModelEvaluator::InArgs Stokhos::SGModelEvaluator::createInArgs | ( | ) | const |
Create InArgs.
Definition at line 596 of file Stokhos_SGModelEvaluator.cpp.
EpetraExt::ModelEvaluator::OutArgs Stokhos::SGModelEvaluator::createOutArgs | ( | ) | const |
Create OutArgs.
Definition at line 618 of file Stokhos_SGModelEvaluator.cpp.
void Stokhos::SGModelEvaluator::evalModel | ( | const InArgs & | inArgs, |
const OutArgs & | outArgs | ||
) | const |
Evaluate model on InArgs.
Definition at line 654 of file Stokhos_SGModelEvaluator.cpp.
|
virtual |
Set initial solution polynomial.
Implements Stokhos::SGModelEvaluatorBase.
Definition at line 1027 of file Stokhos_SGModelEvaluator.cpp.
|
virtual |
Return initial SG x.
Implements Stokhos::SGModelEvaluatorBase.
Definition at line 1034 of file Stokhos_SGModelEvaluator.cpp.
|
virtual |
Set initial parameter polynomial.
Implements Stokhos::SGModelEvaluatorBase.
Definition at line 1040 of file Stokhos_SGModelEvaluator.cpp.
|
virtual |
Return initial SG parameters.
Implements Stokhos::SGModelEvaluatorBase.
Definition at line 1053 of file Stokhos_SGModelEvaluator.cpp.
|
virtual |
Get indices of SG parameters.
These indices determine which parameter vectors support SG
Implements Stokhos::SGModelEvaluatorBase.
Definition at line 1065 of file Stokhos_SGModelEvaluator.cpp.
|
virtual |
Get indices of SG responses.
These indices determine which response vectors support SG
Implements Stokhos::SGModelEvaluatorBase.
Definition at line 1071 of file Stokhos_SGModelEvaluator.cpp.
|
virtual |
Get base maps of SG responses.
Implements Stokhos::SGModelEvaluatorBase.
Definition at line 1077 of file Stokhos_SGModelEvaluator.cpp.
|
virtual |
Return overlap stochastic map.
Implements Stokhos::SGModelEvaluatorBase.
Definition at line 1086 of file Stokhos_SGModelEvaluator.cpp.
|
virtual |
Return x sg overlap map.
Implements Stokhos::SGModelEvaluatorBase.
Definition at line 1092 of file Stokhos_SGModelEvaluator.cpp.
|
virtual |
Return x sg importer.
Implements Stokhos::SGModelEvaluatorBase.
Definition at line 1098 of file Stokhos_SGModelEvaluator.cpp.
|
virtual |
Create vector orthog poly using x map and owned sg map.
Implements Stokhos::SGModelEvaluatorBase.
Definition at line 1104 of file Stokhos_SGModelEvaluator.cpp.
|
virtual |
Create vector orthog poly using x map and overlap sg map.
Implements Stokhos::SGModelEvaluatorBase.
Definition at line 1119 of file Stokhos_SGModelEvaluator.cpp.
|
virtual |
Create vector orthog poly using x map and owned sg map.
Implements Stokhos::SGModelEvaluatorBase.
Definition at line 1135 of file Stokhos_SGModelEvaluator.cpp.
|
virtual |
Create vector orthog poly using x map and overlap sg map.
Implements Stokhos::SGModelEvaluatorBase.
Definition at line 1151 of file Stokhos_SGModelEvaluator.cpp.
|
virtual |
Create vector orthog poly using p map.
Implements Stokhos::SGModelEvaluatorBase.
Definition at line 1169 of file Stokhos_SGModelEvaluator.cpp.
|
virtual |
Create vector orthog poly using f map and owned sg map.
Implements Stokhos::SGModelEvaluatorBase.
Definition at line 1191 of file Stokhos_SGModelEvaluator.cpp.
|
virtual |
Create vector orthog poly using f map and overlap sg map.
Implements Stokhos::SGModelEvaluatorBase.
Definition at line 1206 of file Stokhos_SGModelEvaluator.cpp.
|
virtual |
Create multi-vector orthog poly using f map and owned sg map.
Implements Stokhos::SGModelEvaluatorBase.
Definition at line 1222 of file Stokhos_SGModelEvaluator.cpp.
|
virtual |
Create multi-vector orthog poly using f map and overlap sg map.
Implements Stokhos::SGModelEvaluatorBase.
Definition at line 1240 of file Stokhos_SGModelEvaluator.cpp.
|
virtual |
Create vector orthog poly using g map.
Implements Stokhos::SGModelEvaluatorBase.
Definition at line 1258 of file Stokhos_SGModelEvaluator.cpp.
|
virtual |
Create multi-vector orthog poly using g map.
Implements Stokhos::SGModelEvaluatorBase.
Definition at line 1282 of file Stokhos_SGModelEvaluator.cpp.
Teuchos::RCP< EpetraExt::BlockVector > Stokhos::SGModelEvaluator::import_solution | ( | const Epetra_Vector & | x | ) | const |
Import parallel solution vector.
Definition at line 1307 of file Stokhos_SGModelEvaluator.cpp.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > Stokhos::SGModelEvaluator::import_solution_poly | ( | const Epetra_Vector & | x | ) | const |
Import parallel solution vector.
Definition at line 1316 of file Stokhos_SGModelEvaluator.cpp.
Teuchos::RCP< EpetraExt::BlockVector > Stokhos::SGModelEvaluator::export_solution | ( | const Epetra_Vector & | x_overlapped | ) | const |
Export parallel solution vector.
Definition at line 1328 of file Stokhos_SGModelEvaluator.cpp.
Teuchos::RCP< EpetraExt::BlockVector > Stokhos::SGModelEvaluator::import_residual | ( | const Epetra_Vector & | f | ) | const |
Import parallel residual vector.
Definition at line 1337 of file Stokhos_SGModelEvaluator.cpp.
Teuchos::RCP< EpetraExt::BlockVector > Stokhos::SGModelEvaluator::export_residual | ( | const Epetra_Vector & | f_overlapped | ) | const |
Export parallel residual vector.
Definition at line 1346 of file Stokhos_SGModelEvaluator.cpp.
|
protected |
Underlying model evaluator.
Definition at line 246 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Stochastic Galerkin basis.
Definition at line 249 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Stochastic Galerkin quadrature.
Definition at line 252 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Stochastic Galerkin expansion.
Definition at line 255 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Algorithmic parameters.
Definition at line 258 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Number of stochastic blocks.
Definition at line 261 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Number of W stochastic blocks (may be smaller than num_sg_blocks)
Definition at line 264 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Number of p stochastic blocks (may be smaller than num_sg_blocks)
Definition at line 267 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Whether we support x (and thus f and W)
Definition at line 270 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Underlying unknown map.
Definition at line 273 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Underlying residual map.
Definition at line 276 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Parallel SG data.
Definition at line 279 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Parallel SG communicator.
Definition at line 282 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Epetra Cijk.
Definition at line 285 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Serial Epetra Cijk for dgdx*.
Definition at line 288 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Map for stochastic blocks.
Definition at line 291 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Overlapped map for stochastic blocks (local map)
Definition at line 294 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Overlapped map for p stochastic blocks (local map)
Definition at line 297 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Block SG unknown map.
Definition at line 300 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Block SG overlapped unknown map.
Definition at line 303 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Block SG residual map.
Definition at line 306 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Block SG overlapped residual map.
Definition at line 309 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Importer from SG to SG-overlapped maps.
Definition at line 312 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Exporter from SG-overlapped to SG maps.
Definition at line 315 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Number of parameter vectors of underlying model evaluator.
Definition at line 318 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Number of stochastic parameter vectors.
Definition at line 321 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Index map between block-p and p_sg maps.
Definition at line 324 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Block SG parameter map.
Definition at line 327 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
SG coefficient parameter names.
Definition at line 330 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Number of response vectors of underlying model evaluator.
Definition at line 333 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Number of stochastic response vectors.
Definition at line 336 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Index map between block-g and g_sg maps.
Definition at line 339 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Block SG response map.
Definition at line 342 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
x_dot stochastic Galerkin components
Definition at line 345 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
x stochastic Galerkin components
Definition at line 348 of file Stokhos_SGModelEvaluator.hpp.
|
mutableprotected |
f stochastic Galerkin components
Definition at line 351 of file Stokhos_SGModelEvaluator.hpp.
|
mutableprotected |
W stochastic Galerkin components.
Definition at line 354 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
SG initial x.
Definition at line 357 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
SG initial p.
Definition at line 360 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Whether to always evaluate W with f.
Definition at line 363 of file Stokhos_SGModelEvaluator.hpp.
|
mutableprotected |
W pointer for evaluating W with f.
Definition at line 366 of file Stokhos_SGModelEvaluator.hpp.
|
mutableprotected |
x pointer for evaluating preconditioner
Definition at line 369 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Definition at line 371 of file Stokhos_SGModelEvaluator.hpp.
|
protected |
Preconditioner factory.
Definition at line 374 of file Stokhos_SGModelEvaluator.hpp.