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_Interlaced.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Stokhos Package
5 // Copyright (2009) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Eric T. Phipps (etphipp@sandia.gov).
38 //
39 // ***********************************************************************
40 // @HEADER
41 
42 #ifndef STOKHOS_SGMODELEVALUATOR_INTERLACED_HPP
43 #define STOKHOS_SGMODELEVALUATOR_INTERLACED_HPP
44 
45 #include <vector>
46 
47 #include "EpetraExt_ModelEvaluator.h"
48 #include "EpetraExt_MultiComm.h"
49 #include "EpetraExt_BlockVector.h"
50 
52 #include "Teuchos_RCP.hpp"
53 #include "Teuchos_Array.hpp"
55 #include "Stokhos_ParallelData.hpp"
62 #include "Stokhos_SGOperator.hpp"
64 
65 namespace Stokhos {
66 
68 
83  public:
84 
85  // Constructor
93  bool scaleOP = true);
94 
97 
98  // inputs
100 
103 
106 
109  get_p_names(int l) const;
110 
113 
116 
117  // outputs
119 
122 
125 
128 
129  // ????
131 
133  InArgs createInArgs() const;
134 
136  OutArgs createOutArgs() const;
137 
139  void evalModel(const InArgs& inArgs, const OutArgs& outArgs) const;
140 
142 
145 
147  void set_x_sg_init(const Stokhos::EpetraVectorOrthogPoly& x_sg_in);
148 
151 
153  void set_p_sg_init(int i, const Stokhos::EpetraVectorOrthogPoly& p_sg_in);
154 
157 
159 
163 
165 
169 
172 
175 
178 
181 
185  const Epetra_Vector* v = NULL) const;
186 
190  const Epetra_Vector* v = NULL) const;
191 
194  create_x_mv_sg(int num_vecs,
195  Epetra_DataAccess CV = Copy,
196  const Epetra_MultiVector* v = NULL) const;
197 
200  create_x_mv_sg_overlap(int num_vecs,
201  Epetra_DataAccess CV = Copy,
202  const Epetra_MultiVector* v = NULL) const;
203 
206  create_p_sg(int l, Epetra_DataAccess CV = Copy,
207  const Epetra_Vector* v = NULL) const;
208 
212  const Epetra_Vector* v = NULL) const;
213 
217  const Epetra_Vector* v = NULL) const;
218 
221  create_f_mv_sg(int num_vecs, Epetra_DataAccess CV = Copy,
222  const Epetra_MultiVector* v = NULL) const;
223 
226  create_f_mv_sg_overlap(int num_vecs, Epetra_DataAccess CV = Copy,
227  const Epetra_MultiVector* v = NULL) const;
228 
231  create_g_sg(int l, Epetra_DataAccess CV = Copy,
232  const Epetra_Vector* v = NULL) const;
233 
236  create_g_mv_sg(int l, int num_vecs, Epetra_DataAccess CV = Copy,
237  const Epetra_MultiVector* v = NULL) const;
238 
240 
245  buildInterlaceMap(const Epetra_BlockMap & determ_map,const Epetra_BlockMap & stocha_map);
246 
250  Epetra_Vector & x);
251 
254  static void copyToPolyOrthogVector(const Epetra_Vector & x,
256 
257  protected:
258 
261 
264 
267 
270 
273 
275  unsigned int num_sg_blocks;
276 
278  unsigned int num_W_blocks;
279 
281  unsigned int num_p_blocks;
282 
285 
288 
291 
294 
297 
300 
303 
306 
309 
312 
315 
318 
321 
324 
327 
330 
332  int num_p;
333 
335  int num_p_sg;
336 
339 
342 
345 
347  int num_g;
348 
350  int num_g_sg;
351 
354 
357 
360 
363 
366 
369 
371 
374 
377 
380 
383 
386 
389 
392 
393  bool scaleOP;
394 
395  };
396 
397 }
398 
399 #endif
Teuchos::RCP< const Epetra_Map > get_x_map() const
Return solution vector map.
Teuchos::RCP< const Epetra_Map > get_g_map(int l) const
Return response map.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > x_sg_blocks
x 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.
static Teuchos::RCP< Epetra_Map > buildInterlaceMap(const Epetra_BlockMap &determ_map, const Epetra_BlockMap &stocha_map)
Teuchos::RCP< const Epetra_Map > x_map
Underlying unknown map.
bool supports_x
Whether we support x (and thus f and W)
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.
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.
int num_g
Number of response vectors of underlying model evaluator.
Teuchos::RCP< Stokhos::EpetraOperatorOrthogPoly > W_sg_blocks
W stochastic Galerkin components.
Teuchos::Array< Teuchos::RCP< const Epetra_Map > > get_g_sg_base_maps() const
Get base maps of SG responses.
Teuchos::RCP< const Epetra_Map > interlace_overlapped_f_map
Block SG overlapped residual map.
Teuchos::RCP< const Stokhos::EpetraVectorOrthogPoly > get_p_sg_init(int l) const
Return initial SG parameters.
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 > &params, bool scaleOP=true)
Teuchos::RCP< const Stokhos::EpetraSparse3Tensor > epetraCijk
Epetra Cijk.
Teuchos::Array< int > get_g_sg_map_indices() const
Get indices of SG responses.
Teuchos::RCP< const Epetra_BlockMap > get_overlap_stochastic_map() const
Return overlap stochastic map.
Teuchos::RCP< const Epetra_Vector > get_p_init(int l) const
Return initial parameters.
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::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< int > get_p_sg_map_indices() const
Get indices of SG parameters.
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< const Epetra_BlockMap > get_x_sg_overlap_map() const
Return x sg overlap map.
Abstract base class for orthogonal polynomial-based expansions.
Teuchos::RCP< const Epetra_Map > get_f_map() const
Return residual vector map.
Teuchos::RCP< const Epetra_Map > get_p_map(int l) const
Return parameter vector map.
static void copyToPolyOrthogVector(const Epetra_Vector &x, Stokhos::EpetraVectorOrthogPoly &x_sg)
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< 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.
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::Array< Teuchos::RCP< Teuchos::Array< std::string > > > sg_p_names
SG coefficient parameter names.
Teuchos::RCP< const EpetraExt::MultiComm > sg_comm
Parallel SG communicator.
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 Stokhos::ParallelData > sg_parallel_data
Parallel SG data.
Teuchos::RCP< Stokhos::SGOperator > my_W
W pointer for evaluating W with f.
Teuchos::RCP< Epetra_Operator > create_W() const
Create W = alpha*M + beta*J matrix.
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< const Epetra_BlockMap > stoch_row_map
Map for stochastic blocks.
Teuchos::RCP< Epetra_Vector > my_x
x pointer for evaluating preconditioner
Teuchos::RCP< const Epetra_Map > interlace_overlapped_x_map
Block SG overlapped unknown map.
Teuchos::RCP< const Epetra_Import > get_x_sg_importer() const
Return x sg importer.
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< 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.
void set_p_sg_init(int i, const Stokhos::EpetraVectorOrthogPoly &p_sg_in)
Set initial parameter polynomial.
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.
void evalModel(const InArgs &inArgs, const OutArgs &outArgs) const
Evaluate model on InArgs.
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)
unsigned int num_W_blocks
Number of W stochastic blocks (may be smaller than num_sg_blocks)
Base class for stochastic Galerkin model evaluators.
Teuchos::RCP< Teuchos::ParameterList > params
Algorithmic parameters.
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::EpetraVectorOrthogPoly > get_x_sg_init() const
Return initial SG x.
Teuchos::Array< Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > > sg_p_init
SG initial p.
Teuchos::RCP< const Epetra_Map > f_map
Underlying residual 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< 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::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.
Epetra_DataAccess
Teuchos::RCP< EpetraExt::ModelEvaluator > me
Underlying model evaluator.
Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > f_sg_blocks
f stochastic Galerkin components
Teuchos::RCP< const Teuchos::Array< std::string > > get_p_names(int l) const
Return array of parameter names.
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)