42#ifndef STOKHOS_SGMODELEVALUATOR_HPP
43#define STOKHOS_SGMODELEVALUATOR_HPP
47#include "EpetraExt_ModelEvaluator.h"
48#include "EpetraExt_MultiComm.h"
49#include "EpetraExt_BlockVector.h"
52#include "Teuchos_RCP.hpp"
53#include "Teuchos_Array.hpp"
54#include "Teuchos_ParameterList.hpp"
95 const Teuchos::RCP<EpetraExt::ModelEvaluator>&
me,
100 const Teuchos::RCP<Teuchos::ParameterList>&
params,
107 Teuchos::RCP<const Epetra_Map>
get_x_map()
const;
110 Teuchos::RCP<const Epetra_Map>
get_f_map()
const;
113 Teuchos::RCP<const Epetra_Map>
get_p_map(
int l)
const;
116 Teuchos::RCP<const Epetra_Map>
get_g_map(
int l)
const;
119 Teuchos::RCP<const Teuchos::Array<std::string> >
123 Teuchos::RCP<const Epetra_Vector>
get_x_init()
const;
126 Teuchos::RCP<const Epetra_Vector>
get_p_init(
int l)
const;
129 Teuchos::RCP<Epetra_Operator>
create_W()
const;
132 Teuchos::RCP<EpetraExt::ModelEvaluator::Preconditioner>
create_WPrec()
const;
153 void evalModel(
const InArgs& inArgs,
const OutArgs& outArgs)
const;
164 Teuchos::RCP<const Stokhos::EpetraVectorOrthogPoly>
get_x_sg_init()
const;
170 Teuchos::RCP<const Stokhos::EpetraVectorOrthogPoly>
get_p_sg_init(
int l)
const;
197 Teuchos::RCP<Stokhos::EpetraVectorOrthogPoly>
202 Teuchos::RCP<Stokhos::EpetraVectorOrthogPoly>
207 Teuchos::RCP<Stokhos::EpetraMultiVectorOrthogPoly>
213 Teuchos::RCP<Stokhos::EpetraMultiVectorOrthogPoly>
219 Teuchos::RCP<Stokhos::EpetraVectorOrthogPoly>
224 Teuchos::RCP<Stokhos::EpetraVectorOrthogPoly>
229 Teuchos::RCP<Stokhos::EpetraVectorOrthogPoly>
234 Teuchos::RCP<Stokhos::EpetraMultiVectorOrthogPoly>
239 Teuchos::RCP<Stokhos::EpetraMultiVectorOrthogPoly>
244 Teuchos::RCP<Stokhos::EpetraVectorOrthogPoly>
249 Teuchos::RCP<Stokhos::EpetraMultiVectorOrthogPoly>
256 Teuchos::RCP<EpetraExt::BlockVector>
260 Teuchos::RCP<Stokhos::EpetraVectorOrthogPoly>
264 Teuchos::RCP<EpetraExt::BlockVector>
268 Teuchos::RCP<EpetraExt::BlockVector>
272 Teuchos::RCP<EpetraExt::BlockVector>
278 Teuchos::RCP<EpetraExt::ModelEvaluator>
me;
281 Teuchos::RCP<const Stokhos::OrthogPolyBasis<int, double> >
sg_basis;
284 Teuchos::RCP<const Stokhos::Quadrature<int,double> >
sg_quad;
287 Teuchos::RCP<Stokhos::OrthogPolyExpansion<int,double> >
sg_exp;
290 Teuchos::RCP<Teuchos::ParameterList>
params;
305 Teuchos::RCP<const Epetra_Map>
x_map;
308 Teuchos::RCP<const Epetra_Map>
f_map;
314 Teuchos::RCP<const EpetraExt::MultiComm>
sg_comm;
359 Teuchos::Array< Teuchos::RCP<const Epetra_Map> >
sg_p_map;
362 Teuchos::Array< Teuchos::RCP< Teuchos::Array<std::string> > >
sg_p_names;
374 Teuchos::Array< Teuchos::RCP<const Epetra_Map> >
sg_g_map;
383 mutable Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly >
f_sg_blocks;
386 mutable Teuchos::RCP< Stokhos::EpetraOperatorOrthogPoly >
W_sg_blocks;
389 Teuchos::RCP<Stokhos::EpetraVectorOrthogPoly>
sg_x_init;
392 Teuchos::Array< Teuchos::RCP<Stokhos::EpetraVectorOrthogPoly> >
sg_p_init;
398 mutable Teuchos::RCP<Stokhos::SGOperator>
my_W;
401 mutable Teuchos::RCP<Epetra_Vector>
my_x;
A container class storing an orthogonal polynomial whose coefficients are vectors,...
Abstract base class for multivariate orthogonal polynomials.
Abstract base class for orthogonal polynomial-based expansions.
Abstract base class for quadrature methods.
Base class for stochastic Galerkin model evaluators.
Nonlinear, stochastic Galerkin ModelEvaluator.
Teuchos::RCP< const Stokhos::OrthogPolyBasis< int, double > > sg_basis
Stochastic Galerkin basis.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > create_x_sg_overlap(Epetra_DataAccess CV=Copy, const Epetra_Vector *v=NULL) const
Create vector orthog poly using x map and overlap sg map.
Teuchos::RCP< Stokhos::EpetraMultiVectorOrthogPoly > create_f_mv_sg_overlap(int num_vecs, Epetra_DataAccess CV=Copy, const Epetra_MultiVector *v=NULL) const
Create multi-vector orthog poly using f map and overlap sg map.
OutArgs createOutArgs() const
Create OutArgs.
Teuchos::RCP< EpetraExt::BlockVector > export_residual(const Epetra_Vector &f_overlapped) const
Export parallel residual vector.
Teuchos::RCP< const Stokhos::ParallelData > sg_parallel_data
Parallel SG data.
Teuchos::RCP< const Epetra_Map > sg_f_map
Block SG residual map.
int num_p_sg
Number of stochastic parameter vectors.
Teuchos::RCP< Epetra_Operator > create_DgDx_dot_op(int j) const
Create SG operator representing dg/dxdot.
Teuchos::RCP< EpetraExt::ModelEvaluator > me
Underlying model evaluator.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > sg_x_init
SG initial x.
int num_g_sg
Number of stochastic response vectors.
Teuchos::RCP< const Teuchos::Array< std::string > > get_p_names(int l) const
Return array of parameter names.
Teuchos::RCP< Stokhos::SGPreconditionerFactory > sg_prec_factory
Preconditioner factory.
Teuchos::RCP< const Stokhos::EpetraSparse3Tensor > serialCijk
Serial Epetra Cijk for dgdx*.
Teuchos::RCP< const Epetra_Map > get_g_map(int l) const
Return response map.
Teuchos::Array< int > get_g_sg_map_indices() const
Get indices of SG responses.
Teuchos::RCP< const Stokhos::Quadrature< int, double > > sg_quad
Stochastic Galerkin quadrature.
Teuchos::RCP< const Epetra_Map > f_map
Underlying residual map.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > create_f_sg_overlap(Epetra_DataAccess CV=Copy, const Epetra_Vector *v=NULL) const
Create vector orthog poly using f map and overlap sg map.
Teuchos::RCP< const Epetra_Map > get_p_map(int l) const
Return parameter vector map.
Teuchos::RCP< const Epetra_Vector > get_x_init() const
Return initial solution.
unsigned int num_sg_blocks
Number of stochastic blocks.
Teuchos::RCP< const Epetra_Map > sg_overlapped_x_map
Block SG overlapped unknown map.
Teuchos::RCP< Stokhos::EpetraMultiVectorOrthogPoly > create_f_mv_sg(int num_vecs, Epetra_DataAccess CV=Copy, const Epetra_MultiVector *v=NULL) const
Create multi-vector orthog poly using f map and owned sg map.
Teuchos::RCP< const Epetra_BlockMap > overlapped_stoch_row_map
Overlapped map for stochastic blocks (local map)
InArgs createInArgs() const
Create InArgs.
Teuchos::RCP< const Epetra_Map > sg_x_map
Block SG unknown map.
Teuchos::RCP< EpetraExt::BlockVector > import_solution(const Epetra_Vector &x) const
Import parallel solution vector.
Teuchos::RCP< Stokhos::EpetraMultiVectorOrthogPoly > create_x_mv_sg_overlap(int num_vecs, Epetra_DataAccess CV=Copy, const Epetra_MultiVector *v=NULL) const
Create vector orthog poly using x map and overlap sg map.
Teuchos::RCP< Stokhos::EpetraMultiVectorOrthogPoly > create_x_mv_sg(int num_vecs, Epetra_DataAccess CV=Copy, const Epetra_MultiVector *v=NULL) const
Create vector orthog poly using x map and owned sg map.
Teuchos::RCP< Epetra_Operator > create_W() const
Create W = alpha*M + beta*J matrix.
void set_x_sg_init(const Stokhos::EpetraVectorOrthogPoly &x_sg_in)
Set initial solution polynomial.
Teuchos::Array< Teuchos::RCP< const Epetra_Map > > sg_g_map
Block SG response map.
Teuchos::RCP< Epetra_Import > sg_overlapped_x_importer
Importer from SG to SG-overlapped maps.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > create_x_sg(Epetra_DataAccess CV=Copy, const Epetra_Vector *v=NULL) const
Create vector orthog poly using x map and owned sg map.
Teuchos::RCP< const Epetra_BlockMap > overlapped_stoch_p_map
Overlapped map for p stochastic blocks (local map)
bool eval_W_with_f
Whether to always evaluate W with f.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > create_g_sg(int l, Epetra_DataAccess CV=Copy, const Epetra_Vector *v=NULL) const
Create vector orthog poly using g map.
Teuchos::RCP< Stokhos::OrthogPolyExpansion< int, double > > sg_exp
Stochastic Galerkin expansion.
Teuchos::RCP< Teuchos::ParameterList > params
Algorithmic parameters.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > x_sg_blocks
x stochastic Galerkin components
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > create_p_sg(int l, Epetra_DataAccess CV=Copy, const Epetra_Vector *v=NULL) const
Create vector orthog poly using p map.
Teuchos::Array< int > sg_g_index_map
Index map between block-g and g_sg maps.
Teuchos::Array< Teuchos::RCP< const Epetra_Map > > sg_p_map
Block SG parameter map.
Teuchos::RCP< EpetraExt::BlockVector > export_solution(const Epetra_Vector &x_overlapped) const
Export parallel solution vector.
Teuchos::RCP< const Epetra_Vector > get_p_init(int l) const
Return initial parameters.
Teuchos::Array< Teuchos::RCP< Teuchos::Array< std::string > > > sg_p_names
SG coefficient parameter names.
Teuchos::RCP< Epetra_Vector > my_x
x pointer for evaluating preconditioner
Teuchos::RCP< Epetra_Export > sg_overlapped_f_exporter
Exporter from SG-overlapped to SG maps.
Teuchos::Array< int > sg_p_index_map
Index map between block-p and p_sg maps.
Teuchos::RCP< const Epetra_Import > get_x_sg_importer() const
Return x sg importer.
Teuchos::RCP< Epetra_Operator > create_DgDx_op(int j) const
Create SG operator representing dg/dx.
Teuchos::RCP< const EpetraExt::MultiComm > sg_comm
Parallel SG communicator.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > import_solution_poly(const Epetra_Vector &x) const
Import parallel solution vector.
Teuchos::RCP< const Epetra_Map > get_x_map() const
Return solution vector map.
Teuchos::RCP< Epetra_Operator > create_DgDp_op(int j, int i) const
Create SG operator representing dg/dp.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > x_dot_sg_blocks
x_dot stochastic Galerkin components
Teuchos::RCP< EpetraExt::ModelEvaluator::Preconditioner > create_WPrec() const
Create preconditioner operator.
Teuchos::RCP< const Epetra_BlockMap > get_overlap_stochastic_map() const
Return overlap stochastic map.
Teuchos::RCP< Stokhos::EpetraOperatorOrthogPoly > W_sg_blocks
W stochastic Galerkin components.
Teuchos::Array< int > get_p_sg_map_indices() const
Get indices of SG parameters.
Teuchos::RCP< const Epetra_BlockMap > stoch_row_map
Map for stochastic blocks.
Teuchos::Array< Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > > sg_p_init
SG initial p.
Teuchos::RCP< const Stokhos::EpetraVectorOrthogPoly > get_x_sg_init() const
Return initial SG x.
Teuchos::RCP< Stokhos::SGOperator > my_W
W pointer for evaluating W with f.
Teuchos::RCP< const Epetra_Map > get_f_map() const
Return residual vector map.
void set_p_sg_init(int i, const Stokhos::EpetraVectorOrthogPoly &p_sg_in)
Set initial parameter polynomial.
Teuchos::RCP< Stokhos::EpetraMultiVectorOrthogPoly > create_g_mv_sg(int l, int num_vecs, Epetra_DataAccess CV=Copy, const Epetra_MultiVector *v=NULL) const
Create multi-vector orthog poly using g map.
Teuchos::RCP< const Epetra_BlockMap > get_x_sg_overlap_map() const
Return x sg overlap map.
int num_p
Number of parameter vectors of underlying model evaluator.
void evalModel(const InArgs &inArgs, const OutArgs &outArgs) const
Evaluate model on InArgs.
Teuchos::RCP< const Stokhos::EpetraSparse3Tensor > epetraCijk
Epetra Cijk.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > create_f_sg(Epetra_DataAccess CV=Copy, const Epetra_Vector *v=NULL) const
Create vector orthog poly using f map and owned sg map.
bool supports_x
Whether we support x (and thus f and W)
Teuchos::Array< Teuchos::RCP< const Epetra_Map > > get_g_sg_base_maps() const
Get base maps of SG responses.
unsigned int num_W_blocks
Number of W stochastic blocks (may be smaller than num_sg_blocks)
Teuchos::RCP< EpetraExt::BlockVector > import_residual(const Epetra_Vector &f) const
Import parallel residual vector.
Teuchos::RCP< const Epetra_Map > x_map
Underlying unknown map.
Teuchos::RCP< Epetra_Operator > create_DfDp_op(int i) const
Create SG operator representing df/dp.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > f_sg_blocks
f stochastic Galerkin components
Teuchos::RCP< const Stokhos::EpetraVectorOrthogPoly > get_p_sg_init(int l) const
Return initial SG parameters.
unsigned int num_p_blocks
Number of p stochastic blocks (may be smaller than num_sg_blocks)
int num_g
Number of response vectors of underlying model evaluator.
Teuchos::RCP< const Epetra_Map > sg_overlapped_f_map
Block SG overlapped residual map.
ScalarType f(const Teuchos::Array< ScalarType > &x, double a, double b)
Top-level namespace for Stokhos classes and functions.