42 #ifndef STOKHOS_SGMODELEVALUATOR_INTERLACED_HPP 43 #define STOKHOS_SGMODELEVALUATOR_INTERLACED_HPP 47 #include "EpetraExt_ModelEvaluator.h" 48 #include "EpetraExt_MultiComm.h" 49 #include "EpetraExt_BlockVector.h" 51 #include "Teuchos_RCP.hpp" 52 #include "Teuchos_Array.hpp" 53 #include "Teuchos_ParameterList.hpp" 86 const Teuchos::RCP<EpetraExt::ModelEvaluator>&
me,
91 const Teuchos::RCP<Teuchos::ParameterList>&
params,
101 Teuchos::RCP<const Epetra_Map>
get_x_map()
const;
104 Teuchos::RCP<const Epetra_Map>
get_p_map(
int l)
const;
107 Teuchos::RCP<const Teuchos::Array<std::string> >
111 Teuchos::RCP<const Epetra_Vector>
get_x_init()
const;
114 Teuchos::RCP<const Epetra_Vector>
get_p_init(
int l)
const;
120 Teuchos::RCP<const Epetra_Map>
get_f_map()
const;
123 Teuchos::RCP<const Epetra_Map>
get_g_map(
int l)
const;
126 Teuchos::RCP<Epetra_Operator>
create_W()
const;
138 void evalModel(
const InArgs& inArgs,
const OutArgs& outArgs)
const;
146 Teuchos::RCP<const Stokhos::EpetraVectorOrthogPoly>
get_x_sg_init()
const;
152 Teuchos::RCP<const Stokhos::EpetraVectorOrthogPoly>
get_p_sg_init(
int l)
const;
179 Teuchos::RCP<Stokhos::EpetraVectorOrthogPoly>
181 const Epetra_Vector* v = NULL)
const;
184 Teuchos::RCP<Stokhos::EpetraVectorOrthogPoly>
186 const Epetra_Vector* v = NULL)
const;
189 Teuchos::RCP<Stokhos::EpetraMultiVectorOrthogPoly>
191 Epetra_DataAccess CV = Copy,
192 const Epetra_MultiVector* v = NULL)
const;
195 Teuchos::RCP<Stokhos::EpetraMultiVectorOrthogPoly>
197 Epetra_DataAccess CV = Copy,
198 const Epetra_MultiVector* v = NULL)
const;
201 Teuchos::RCP<Stokhos::EpetraVectorOrthogPoly>
203 const Epetra_Vector* v = NULL)
const;
206 Teuchos::RCP<Stokhos::EpetraVectorOrthogPoly>
208 const Epetra_Vector* v = NULL)
const;
211 Teuchos::RCP<Stokhos::EpetraVectorOrthogPoly>
213 const Epetra_Vector* v = NULL)
const;
216 Teuchos::RCP<Stokhos::EpetraMultiVectorOrthogPoly>
218 const Epetra_MultiVector* v = NULL)
const;
221 Teuchos::RCP<Stokhos::EpetraMultiVectorOrthogPoly>
223 const Epetra_MultiVector* v = NULL)
const;
226 Teuchos::RCP<Stokhos::EpetraVectorOrthogPoly>
228 const Epetra_Vector* v = NULL)
const;
231 Teuchos::RCP<Stokhos::EpetraMultiVectorOrthogPoly>
233 const Epetra_MultiVector* v = NULL)
const;
238 static Teuchos::RCP<Epetra_Map>
239 buildInterlaceMap(
const Epetra_BlockMap & determ_map,
const Epetra_BlockMap & stocha_map);
254 Teuchos::RCP<EpetraExt::ModelEvaluator>
me;
257 Teuchos::RCP<const Stokhos::OrthogPolyBasis<int, double> >
sg_basis;
260 Teuchos::RCP<const Stokhos::Quadrature<int,double> >
sg_quad;
263 Teuchos::RCP<Stokhos::OrthogPolyExpansion<int,double> >
sg_exp;
266 Teuchos::RCP<Teuchos::ParameterList>
params;
281 Teuchos::RCP<const Epetra_Map>
x_map;
284 Teuchos::RCP<const Epetra_Map>
f_map;
290 Teuchos::RCP<const EpetraExt::MultiComm>
sg_comm;
335 Teuchos::Array< Teuchos::RCP<const Epetra_Map> >
sg_p_map;
338 Teuchos::Array< Teuchos::RCP< Teuchos::Array<std::string> > >
sg_p_names;
350 Teuchos::Array< Teuchos::RCP<const Epetra_Map> >
sg_g_map;
359 mutable Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly >
f_sg_blocks;
362 mutable Teuchos::RCP< Stokhos::EpetraOperatorOrthogPoly >
W_sg_blocks;
364 mutable Teuchos::Array< Teuchos::RCP< Stokhos::EpetraMultiVectorOrthogPoly > >
dfdp_sg_blocks;
370 mutable Teuchos::Array< Teuchos::RCP< Stokhos::EpetraMultiVectorOrthogPoly > >
dgdx_sg_blocks;
373 Teuchos::RCP<Stokhos::EpetraVectorOrthogPoly>
sg_x_init;
376 Teuchos::Array< Teuchos::RCP<Stokhos::EpetraVectorOrthogPoly> >
sg_p_init;
382 mutable Teuchos::RCP<Stokhos::SGOperator>
my_W;
385 mutable Teuchos::RCP<Epetra_Vector>
my_x;
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_Vector > get_p_init(int l) const
Return initial parameters.
Teuchos::Array< int > get_g_sg_map_indices() const
Get indices of SG responses.
Teuchos::RCP< const Stokhos::EpetraVectorOrthogPoly > get_x_sg_init() const
Return initial SG x.
void evalModel(const InArgs &inArgs, const OutArgs &outArgs) const
Evaluate model on InArgs.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > x_sg_blocks
x stochastic Galerkin components
static Teuchos::RCP< Epetra_Map > buildInterlaceMap(const Epetra_BlockMap &determ_map, const Epetra_BlockMap &stocha_map)
Teuchos::RCP< const Epetra_Map > get_f_map() const
Return residual vector map.
Teuchos::RCP< const Epetra_Map > x_map
Underlying unknown map.
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.
bool supports_x
Whether we support x (and thus f and W)
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< const Stokhos::EpetraSparse3Tensor > serialCijk
Serial Epetra Cijk for dgdx*.
Teuchos::RCP< const Epetra_Map > interlace_x_map
Block SG unknown map.
unsigned int num_sg_blocks
Number of stochastic blocks.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > x_dot_sg_blocks
x_dot stochastic Galerkin components
Teuchos::RCP< const Epetra_Vector > get_x_init() const
Return initial solution.
Teuchos::RCP< const Stokhos::Quadrature< int, double > > sg_quad
Stochastic Galerkin quadrature.
int num_g_sg
Number of stochastic response vectors.
Teuchos::RCP< Epetra_Import > interlace_overlapped_x_importer
Importer from SG to SG-overlapped maps.
int num_g
Number of response vectors of underlying model evaluator.
OutArgs createOutArgs() const
Create OutArgs.
Teuchos::RCP< Stokhos::EpetraOperatorOrthogPoly > W_sg_blocks
W stochastic Galerkin components.
Teuchos::RCP< const Epetra_Map > interlace_overlapped_f_map
Block SG overlapped residual 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.
SGModelEvaluator_Interlaced(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< const Stokhos::EpetraSparse3Tensor > epetraCijk
Epetra Cijk.
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.
Teuchos::RCP< const Stokhos::OrthogPolyBasis< int, double > > sg_basis
Stochastic Galerkin basis.
Teuchos::Array< Teuchos::RCP< Stokhos::EpetraMultiVectorOrthogPoly > > dgdx_sg_blocks
dg/dx stochastic Galerkin components
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.
Abstract base class for orthogonal polynomial-based expansions.
Teuchos::RCP< const Epetra_Import > get_x_sg_importer() const
Return x sg importer.
static void copyToPolyOrthogVector(const Epetra_Vector &x, Stokhos::EpetraVectorOrthogPoly &x_sg)
const IndexType const IndexType const IndexType const IndexType const ValueType const ValueType * x
InArgs createInArgs() const
Create InArgs.
void set_x_sg_init(const Stokhos::EpetraVectorOrthogPoly &x_sg_in)
Set initial solution polynomial.
Teuchos::RCP< const Epetra_BlockMap > overlapped_stoch_row_map
Overlapped map for stochastic blocks (local map)
int num_p
Number of parameter vectors of underlying model evaluator.
Teuchos::RCP< const Epetra_BlockMap > get_overlap_stochastic_map() const
Return overlap stochastic map.
Abstract base class for quadrature methods.
A container class storing an orthogonal polynomial whose coefficients are vectors, operators, or in general any type that would have an expensive copy constructor.
Teuchos::RCP< const Epetra_BlockMap > get_x_sg_overlap_map() const
Return x sg overlap map.
Teuchos::Array< Teuchos::RCP< Teuchos::Array< std::string > > > sg_p_names
SG coefficient parameter names.
Teuchos::RCP< const EpetraExt::MultiComm > sg_comm
Parallel SG communicator.
Top-level namespace for Stokhos classes and functions.
Teuchos::RCP< const Stokhos::ParallelData > sg_parallel_data
Parallel SG data.
Teuchos::RCP< Stokhos::SGOperator > my_W
W pointer for evaluating W with f.
Teuchos::Array< Teuchos::RCP< const Epetra_Map > > get_g_sg_base_maps() const
Get base maps of SG responses.
Teuchos::Array< Teuchos::RCP< const Epetra_Map > > sg_g_map
Block SG response map.
Teuchos::Array< Teuchos::RCP< Stokhos::EpetraMultiVectorOrthogPoly > > dgdx_dot_sg_blocks
dg/dxdot stochastic Galerkin components
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.
Teuchos::RCP< const Epetra_BlockMap > stoch_row_map
Map for stochastic blocks.
Teuchos::RCP< Epetra_Vector > my_x
x pointer for evaluating preconditioner
Teuchos::RCP< Epetra_Operator > create_W() const
Create W = alpha*M + beta*J matrix.
Teuchos::RCP< const Epetra_Map > interlace_overlapped_x_map
Block SG overlapped unknown map.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > sg_x_init
SG initial x.
unsigned int num_p_blocks
Number of p stochastic blocks (may be smaller than num_sg_blocks)
Teuchos::RCP< const Stokhos::EpetraVectorOrthogPoly > get_p_sg_init(int l) const
Return initial SG parameters.
void set_p_sg_init(int i, const Stokhos::EpetraVectorOrthogPoly &p_sg_in)
Set initial parameter polynomial.
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< Epetra_Export > interlace_overlapped_f_exporter
Exporter from SG-overlapped to SG maps.
Teuchos::RCP< const Epetra_Map > interlace_f_map
Block SG residual map.
Teuchos::RCP< const Epetra_Map > get_p_map(int l) const
Return parameter vector map.
bool eval_W_with_f
Whether to always evaluate W with f.
Teuchos::Array< Teuchos::RCP< Stokhos::EpetraMultiVectorOrthogPoly > > dfdp_sg_blocks
Teuchos::RCP< const Epetra_BlockMap > overlapped_stoch_p_map
Overlapped map for p stochastic blocks (local map)
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::RCP< const Epetra_Map > get_g_map(int l) const
Return response map.
unsigned int num_W_blocks
Number of W stochastic blocks (may be smaller than num_sg_blocks)
Teuchos::RCP< Teuchos::ParameterList > params
Algorithmic parameters.
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::Array< Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > > sg_p_init
SG initial p.
Teuchos::Array< int > get_p_sg_map_indices() const
Get indices of SG parameters.
Teuchos::RCP< const Epetra_Map > f_map
Underlying residual map.
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< const Teuchos::Array< std::string > > get_p_names(int l) const
Return array of parameter names.
Teuchos::RCP< EpetraExt::ModelEvaluator > me
Underlying model evaluator.
Teuchos::RCP< const Epetra_Map > get_x_map() const
Return solution vector map.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > f_sg_blocks
f stochastic Galerkin components
int num_p_sg
Number of stochastic parameter vectors.
Teuchos::Array< Teuchos::RCP< const Epetra_Map > > sg_p_map
Block SG parameter map.
Teuchos::Array< int > sg_p_index_map
Index map between block-p and p_sg maps.
Nonlinear, stochastic Galerkin ModelEvaluator that constructs a interlaced Jacobian.
Teuchos::RCP< Stokhos::OrthogPolyExpansion< int, double > > sg_exp
Stochastic Galerkin expansion.
Teuchos::Array< int > sg_g_index_map
Index map between block-g and g_sg maps.
static void copyToInterlacedVector(const Stokhos::EpetraVectorOrthogPoly &x_sg, Epetra_Vector &x)