Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
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 
51 #include "Teuchos_RCP.hpp"
52 #include "Teuchos_Array.hpp"
53 #include "Teuchos_ParameterList.hpp"
54 #include "Stokhos_ParallelData.hpp"
61 #include "Stokhos_SGOperator.hpp"
63 
64 namespace Stokhos {
65 
67 
81  class SGModelEvaluator_Interlaced : public EpetraExt::ModelEvaluator {
82  public:
83 
84  // Constructor
86  const Teuchos::RCP<EpetraExt::ModelEvaluator>& me,
87  const Teuchos::RCP<const Stokhos::OrthogPolyBasis<int,double> >& sg_basis,
88  const Teuchos::RCP<const Stokhos::Quadrature<int,double> >& sg_quad,
90  const Teuchos::RCP<const Stokhos::ParallelData>& sg_parallel_data,
91  const Teuchos::RCP<Teuchos::ParameterList>& params,
92  bool scaleOP = true);
93 
96 
97  // inputs
99 
101  Teuchos::RCP<const Epetra_Map> get_x_map() const;
102 
104  Teuchos::RCP<const Epetra_Map> get_p_map(int l) const;
105 
107  Teuchos::RCP<const Teuchos::Array<std::string> >
108  get_p_names(int l) const;
109 
111  Teuchos::RCP<const Epetra_Vector> get_x_init() const;
112 
114  Teuchos::RCP<const Epetra_Vector> get_p_init(int l) const;
115 
116  // outputs
118 
120  Teuchos::RCP<const Epetra_Map> get_f_map() const;
121 
123  Teuchos::RCP<const Epetra_Map> get_g_map(int l) const;
124 
126  Teuchos::RCP<Epetra_Operator> create_W() const;
127 
128  // ????
130 
132  InArgs createInArgs() const;
133 
135  OutArgs createOutArgs() const;
136 
138  void evalModel(const InArgs& inArgs, const OutArgs& outArgs) const;
139 
141 
143  void set_x_sg_init(const Stokhos::EpetraVectorOrthogPoly& x_sg_in);
144 
146  Teuchos::RCP<const Stokhos::EpetraVectorOrthogPoly> get_x_sg_init() const;
147 
149  void set_p_sg_init(int i, const Stokhos::EpetraVectorOrthogPoly& p_sg_in);
150 
152  Teuchos::RCP<const Stokhos::EpetraVectorOrthogPoly> get_p_sg_init(int l) const;
153 
155 
158  Teuchos::Array<int> get_p_sg_map_indices() const;
159 
161 
164  Teuchos::Array<int> get_g_sg_map_indices() const;
165 
167  Teuchos::Array< Teuchos::RCP<const Epetra_Map> > get_g_sg_base_maps() const;
168 
170  Teuchos::RCP<const Epetra_BlockMap> get_overlap_stochastic_map() const;
171 
173  Teuchos::RCP<const Epetra_BlockMap> get_x_sg_overlap_map() const;
174 
176  Teuchos::RCP<const Epetra_Import> get_x_sg_importer() const;
177 
179  Teuchos::RCP<Stokhos::EpetraVectorOrthogPoly>
180  create_x_sg(Epetra_DataAccess CV = Copy,
181  const Epetra_Vector* v = NULL) const;
182 
184  Teuchos::RCP<Stokhos::EpetraVectorOrthogPoly>
185  create_x_sg_overlap(Epetra_DataAccess CV = Copy,
186  const Epetra_Vector* v = NULL) const;
187 
189  Teuchos::RCP<Stokhos::EpetraMultiVectorOrthogPoly>
190  create_x_mv_sg(int num_vecs,
191  Epetra_DataAccess CV = Copy,
192  const Epetra_MultiVector* v = NULL) const;
193 
195  Teuchos::RCP<Stokhos::EpetraMultiVectorOrthogPoly>
196  create_x_mv_sg_overlap(int num_vecs,
197  Epetra_DataAccess CV = Copy,
198  const Epetra_MultiVector* v = NULL) const;
199 
201  Teuchos::RCP<Stokhos::EpetraVectorOrthogPoly>
202  create_p_sg(int l, Epetra_DataAccess CV = Copy,
203  const Epetra_Vector* v = NULL) const;
204 
206  Teuchos::RCP<Stokhos::EpetraVectorOrthogPoly>
207  create_f_sg(Epetra_DataAccess CV = Copy,
208  const Epetra_Vector* v = NULL) const;
209 
211  Teuchos::RCP<Stokhos::EpetraVectorOrthogPoly>
212  create_f_sg_overlap(Epetra_DataAccess CV = Copy,
213  const Epetra_Vector* v = NULL) const;
214 
216  Teuchos::RCP<Stokhos::EpetraMultiVectorOrthogPoly>
217  create_f_mv_sg(int num_vecs, Epetra_DataAccess CV = Copy,
218  const Epetra_MultiVector* v = NULL) const;
219 
221  Teuchos::RCP<Stokhos::EpetraMultiVectorOrthogPoly>
222  create_f_mv_sg_overlap(int num_vecs, Epetra_DataAccess CV = Copy,
223  const Epetra_MultiVector* v = NULL) const;
224 
226  Teuchos::RCP<Stokhos::EpetraVectorOrthogPoly>
227  create_g_sg(int l, Epetra_DataAccess CV = Copy,
228  const Epetra_Vector* v = NULL) const;
229 
231  Teuchos::RCP<Stokhos::EpetraMultiVectorOrthogPoly>
232  create_g_mv_sg(int l, int num_vecs, Epetra_DataAccess CV = Copy,
233  const Epetra_MultiVector* v = NULL) const;
234 
238  static Teuchos::RCP<Epetra_Map>
239  buildInterlaceMap(const Epetra_BlockMap & determ_map,const Epetra_BlockMap & stocha_map);
240 
244  Epetra_Vector & x);
245 
248  static void copyToPolyOrthogVector(const Epetra_Vector & x,
250 
251  protected:
252 
254  Teuchos::RCP<EpetraExt::ModelEvaluator> me;
255 
257  Teuchos::RCP<const Stokhos::OrthogPolyBasis<int, double> > sg_basis;
258 
260  Teuchos::RCP<const Stokhos::Quadrature<int,double> > sg_quad;
261 
263  Teuchos::RCP<Stokhos::OrthogPolyExpansion<int,double> > sg_exp;
264 
266  Teuchos::RCP<Teuchos::ParameterList> params;
267 
269  unsigned int num_sg_blocks;
270 
272  unsigned int num_W_blocks;
273 
275  unsigned int num_p_blocks;
276 
279 
281  Teuchos::RCP<const Epetra_Map> x_map;
282 
284  Teuchos::RCP<const Epetra_Map> f_map;
285 
287  Teuchos::RCP<const Stokhos::ParallelData> sg_parallel_data;
288 
290  Teuchos::RCP<const EpetraExt::MultiComm> sg_comm;
291 
293  Teuchos::RCP<const Stokhos::EpetraSparse3Tensor> epetraCijk;
294 
296  Teuchos::RCP<const Stokhos::EpetraSparse3Tensor> serialCijk;
297 
299  Teuchos::RCP<const Epetra_BlockMap> stoch_row_map;
300 
302  Teuchos::RCP<const Epetra_BlockMap> overlapped_stoch_row_map;
303 
305  Teuchos::RCP<const Epetra_BlockMap> overlapped_stoch_p_map;
306 
308  Teuchos::RCP<const Epetra_Map> interlace_x_map;
309 
311  Teuchos::RCP<const Epetra_Map> interlace_overlapped_x_map;
312 
314  Teuchos::RCP<const Epetra_Map> interlace_f_map;
315 
317  Teuchos::RCP<const Epetra_Map> interlace_overlapped_f_map;
318 
320  Teuchos::RCP<Epetra_Import> interlace_overlapped_x_importer;
321 
323  Teuchos::RCP<Epetra_Export> interlace_overlapped_f_exporter;
324 
326  int num_p;
327 
329  int num_p_sg;
330 
332  Teuchos::Array<int> sg_p_index_map;
333 
335  Teuchos::Array< Teuchos::RCP<const Epetra_Map> > sg_p_map;
336 
338  Teuchos::Array< Teuchos::RCP< Teuchos::Array<std::string> > > sg_p_names;
339 
341  int num_g;
342 
344  int num_g_sg;
345 
347  Teuchos::Array<int> sg_g_index_map;
348 
350  Teuchos::Array< Teuchos::RCP<const Epetra_Map> > sg_g_map;
351 
353  Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > x_dot_sg_blocks;
354 
356  Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > x_sg_blocks;
357 
359  mutable Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > f_sg_blocks;
360 
362  mutable Teuchos::RCP< Stokhos::EpetraOperatorOrthogPoly > W_sg_blocks;
363 
364  mutable Teuchos::Array< Teuchos::RCP< Stokhos::EpetraMultiVectorOrthogPoly > > dfdp_sg_blocks;
365 
367  mutable Teuchos::Array< Teuchos::RCP< Stokhos::EpetraMultiVectorOrthogPoly > > dgdx_dot_sg_blocks;
368 
370  mutable Teuchos::Array< Teuchos::RCP< Stokhos::EpetraMultiVectorOrthogPoly > > dgdx_sg_blocks;
371 
373  Teuchos::RCP<Stokhos::EpetraVectorOrthogPoly> sg_x_init;
374 
376  Teuchos::Array< Teuchos::RCP<Stokhos::EpetraVectorOrthogPoly> > sg_p_init;
377 
380 
382  mutable Teuchos::RCP<Stokhos::SGOperator> my_W;
383 
385  mutable Teuchos::RCP<Epetra_Vector> my_x;
386 
387  bool scaleOP;
388 
389  };
390 
391 }
392 
393 #endif
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.
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 > &params, 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
Definition: csr_vector.h:260
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)