Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_SGModelEvaluator.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Stokhos Package
4 //
5 // Copyright 2009 NTESS and the Stokhos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef STOKHOS_SGMODELEVALUATOR_HPP
11 #define STOKHOS_SGMODELEVALUATOR_HPP
12 
13 #include <vector>
14 
15 #include "EpetraExt_ModelEvaluator.h"
16 #include "EpetraExt_MultiComm.h"
17 #include "EpetraExt_BlockVector.h"
18 
20 #include "Teuchos_RCP.hpp"
21 #include "Teuchos_Array.hpp"
23 #include "Stokhos_ParallelData.hpp"
29 #include "Stokhos_SGOperator.hpp"
32 
33 namespace Stokhos {
34 
36 
59  public:
60 
61  // Constructor
69  bool scaleOP = true);
70 
73 
76 
79 
82 
85 
88  get_p_names(int l) const;
89 
92 
95 
98 
101 
104 
107 
109  Teuchos::RCP<Epetra_Operator> create_DgDp_op(int j, int i) const;
110 
113 
115  InArgs createInArgs() const;
116 
118  OutArgs createOutArgs() const;
119 
121  void evalModel(const InArgs& inArgs, const OutArgs& outArgs) const;
122 
124 
127 
129  void set_x_sg_init(const Stokhos::EpetraVectorOrthogPoly& x_sg_in);
130 
133 
135  void set_p_sg_init(int i, const Stokhos::EpetraVectorOrthogPoly& p_sg_in);
136 
139 
141 
145 
147 
151 
154 
157 
160 
163 
167  const Epetra_Vector* v = NULL) const;
168 
172  const Epetra_Vector* v = NULL) const;
173 
176  create_x_mv_sg(int num_vecs,
177  Epetra_DataAccess CV = Copy,
178  const Epetra_MultiVector* v = NULL) const;
179 
182  create_x_mv_sg_overlap(int num_vecs,
183  Epetra_DataAccess CV = Copy,
184  const Epetra_MultiVector* v = NULL) const;
185 
188  create_p_sg(int l, Epetra_DataAccess CV = Copy,
189  const Epetra_Vector* v = NULL) const;
190 
194  const Epetra_Vector* v = NULL) const;
195 
199  const Epetra_Vector* v = NULL) const;
200 
203  create_f_mv_sg(int num_vecs, Epetra_DataAccess CV = Copy,
204  const Epetra_MultiVector* v = NULL) const;
205 
208  create_f_mv_sg_overlap(int num_vecs, Epetra_DataAccess CV = Copy,
209  const Epetra_MultiVector* v = NULL) const;
210 
213  create_g_sg(int l, Epetra_DataAccess CV = Copy,
214  const Epetra_Vector* v = NULL) const;
215 
218  create_g_mv_sg(int l, int num_vecs, Epetra_DataAccess CV = Copy,
219  const Epetra_MultiVector* v = NULL) const;
220 
222 
225  import_solution(const Epetra_Vector& x) const;
226 
229  import_solution_poly(const Epetra_Vector& x) const;
230 
233  export_solution(const Epetra_Vector& x_overlapped) const;
234 
237  import_residual(const Epetra_Vector& f) const;
238 
241  export_residual(const Epetra_Vector& f_overlapped) const;
242 
243  protected:
244 
247 
250 
253 
256 
259 
261  unsigned int num_sg_blocks;
262 
264  unsigned int num_W_blocks;
265 
267  unsigned int num_p_blocks;
268 
271 
274 
277 
280 
283 
286 
289 
292 
295 
298 
301 
304 
307 
310 
313 
316 
318  int num_p;
319 
321  int num_p_sg;
322 
325 
328 
331 
333  int num_g;
334 
336  int num_g_sg;
337 
340 
343 
346 
349 
352 
355 
358 
361 
364 
367 
370 
371  bool scaleOP;
372 
375 
376  };
377 
378 }
379 
380 #endif // FEAPP_MODELEVALUATOR_HPP
Teuchos::RCP< const Stokhos::EpetraSparse3Tensor > epetraCijk
Epetra Cijk.
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< 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< Stokhos::EpetraOperatorOrthogPoly > W_sg_blocks
W stochastic Galerkin components.
Teuchos::RCP< const Epetra_Map > get_x_map() const
Return solution vector 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.
Teuchos::RCP< const Epetra_Map > sg_overlapped_f_map
Block SG overlapped residual map.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > sg_x_init
SG initial x.
Teuchos::Array< int > sg_g_index_map
Index map between block-g and g_sg maps.
unsigned int num_sg_blocks
Number of stochastic blocks.
Teuchos::RCP< const Epetra_Map > get_g_map(int l) const
Return response map.
Teuchos::RCP< const Epetra_Map > x_map
Underlying unknown map.
Teuchos::RCP< const Epetra_Map > get_f_map() const
Return residual vector map.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > x_sg_blocks
x stochastic Galerkin components
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > f_sg_blocks
f stochastic Galerkin components
int num_p_sg
Number of stochastic parameter vectors.
bool eval_W_with_f
Whether to always evaluate W with f.
void evalModel(const InArgs &inArgs, const OutArgs &outArgs) const
Evaluate model on InArgs.
Teuchos::RCP< const Epetra_BlockMap > get_overlap_stochastic_map() const
Return overlap stochastic map.
Teuchos::RCP< const Epetra_Map > sg_overlapped_x_map
Block SG overlapped unknown map.
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)
Teuchos::Array< Teuchos::RCP< const Epetra_Map > > sg_g_map
Block SG response map.
Teuchos::RCP< const Stokhos::EpetraVectorOrthogPoly > get_p_sg_init(int l) const
Return initial SG parameters.
Teuchos::RCP< Epetra_Operator > create_DgDp_op(int j, int i) const
Create SG operator representing dg/dp.
unsigned int num_p_blocks
Number of p stochastic blocks (may be smaller than num_sg_blocks)
Teuchos::RCP< const Stokhos::ParallelData > sg_parallel_data
Parallel SG data.
Teuchos::RCP< const Epetra_BlockMap > overlapped_stoch_p_map
Overlapped map for p stochastic blocks (local map)
Nonlinear, stochastic Galerkin ModelEvaluator.
Teuchos::Array< Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > > sg_p_init
SG initial p.
Teuchos::RCP< const Stokhos::EpetraSparse3Tensor > serialCijk
Serial Epetra Cijk for dgdx*.
Teuchos::Array< int > get_p_sg_map_indices() const
Get indices of SG parameters.
Teuchos::RCP< const Epetra_Map > sg_f_map
Block SG residual map.
OutArgs createOutArgs() const
Create OutArgs.
Abstract base class for orthogonal polynomial-based expansions.
Teuchos::RCP< EpetraExt::BlockVector > import_solution(const Epetra_Vector &x) const
Import parallel solution vector.
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::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.
int num_g_sg
Number of stochastic response vectors.
Teuchos::RCP< const Epetra_BlockMap > stoch_row_map
Map for stochastic blocks.
Teuchos::RCP< const Epetra_BlockMap > get_x_sg_overlap_map() const
Return x sg overlap map.
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< Stokhos::SGPreconditionerFactory > sg_prec_factory
Preconditioner factory.
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 Stokhos::OrthogPolyBasis< int, double > > sg_basis
Stochastic Galerkin basis.
Teuchos::Array< Teuchos::RCP< const Epetra_Map > > sg_p_map
Block SG parameter map.
bool supports_x
Whether we support x (and thus f and W)
void set_p_sg_init(int i, const Stokhos::EpetraVectorOrthogPoly &p_sg_in)
Set initial parameter polynomial.
Teuchos::RCP< Stokhos::SGOperator > my_W
W pointer for evaluating W with f.
Teuchos::RCP< const Stokhos::EpetraVectorOrthogPoly > get_x_sg_init() const
Return initial SG x.
Teuchos::Array< int > get_g_sg_map_indices() const
Get indices of SG responses.
Teuchos::RCP< Epetra_Export > sg_overlapped_f_exporter
Exporter from SG-overlapped to SG maps.
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_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::Array< Teuchos::RCP< Teuchos::Array< std::string > > > sg_p_names
SG coefficient parameter names.
Teuchos::RCP< EpetraExt::ModelEvaluator > me
Underlying model evaluator.
Teuchos::RCP< const Epetra_Map > get_p_map(int l) const
Return parameter vector 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< 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.
unsigned int num_W_blocks
Number of W stochastic blocks (may be smaller than num_sg_blocks)
InArgs createInArgs() const
Create InArgs.
void set_x_sg_init(const Stokhos::EpetraVectorOrthogPoly &x_sg_in)
Set initial solution polynomial.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > import_solution_poly(const Epetra_Vector &x) const
Import parallel solution vector.
Teuchos::RCP< EpetraExt::BlockVector > import_residual(const Epetra_Vector &f) const
Import parallel residual vector.
Teuchos::RCP< const Epetra_Import > get_x_sg_importer() const
Return x sg importer.
Teuchos::RCP< const Epetra_BlockMap > overlapped_stoch_row_map
Overlapped map for stochastic blocks (local map)
ScalarType f(const Teuchos::Array< ScalarType > &x, double a, double b)
int num_g
Number of response vectors of underlying model evaluator.
Teuchos::RCP< Epetra_Operator > create_DgDx_dot_op(int j) const
Create SG operator representing dg/dxdot.
int num_p
Number of parameter vectors of underlying model evaluator.
Teuchos::RCP< Epetra_Import > sg_overlapped_x_importer
Importer from SG to SG-overlapped maps.
Teuchos::RCP< Epetra_Operator > create_DgDx_op(int j) const
Create SG operator representing dg/dx.
Base class for stochastic Galerkin model evaluators.
Teuchos::RCP< Epetra_Vector > my_x
x pointer for evaluating preconditioner
Teuchos::RCP< Teuchos::ParameterList > params
Algorithmic parameters.
Teuchos::RCP< const Teuchos::Array< std::string > > get_p_names(int l) const
Return array of parameter names.
Teuchos::RCP< EpetraExt::BlockVector > export_residual(const Epetra_Vector &f_overlapped) const
Export parallel residual vector.
Teuchos::RCP< EpetraExt::BlockVector > export_solution(const Epetra_Vector &x_overlapped) const
Export parallel solution vector.
Teuchos::RCP< Epetra_Operator > create_W() const
Create W = alpha*M + beta*J matrix.
Teuchos::Array< int > sg_p_index_map
Index map between block-p and p_sg maps.
Epetra_DataAccess
Teuchos::Array< Teuchos::RCP< const Epetra_Map > > get_g_sg_base_maps() const
Get base maps of SG responses.
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< 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_Vector > get_p_init(int l) const
Return initial parameters.
Teuchos::RCP< Stokhos::OrthogPolyExpansion< int, double > > sg_exp
Stochastic Galerkin expansion.
Teuchos::RCP< EpetraExt::ModelEvaluator::Preconditioner > create_WPrec() const
Create preconditioner operator.
Teuchos::RCP< Epetra_Operator > create_DfDp_op(int i) const
Create SG operator representing df/dp.
Teuchos::RCP< const Stokhos::Quadrature< int, double > > sg_quad
Stochastic Galerkin quadrature.
Teuchos::RCP< const Epetra_Map > sg_x_map
Block SG unknown map.
Teuchos::RCP< const EpetraExt::MultiComm > sg_comm
Parallel SG communicator.