Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_SGModelEvaluatorBase.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_SGMODELEVALUATORBASE_HPP
11 #define STOKHOS_SGMODELEVALUATORBASE_HPP
12 
13 #include "EpetraExt_ModelEvaluator.h"
14 #include "EpetraExt_BlockVector.h"
15 
16 #include "Teuchos_RCP.hpp"
17 #include "Teuchos_Array.hpp"
19 
20 namespace Stokhos {
21 
23  class SGModelEvaluatorBase : public virtual EpetraExt::ModelEvaluator {
24  public:
25 
26  // Constructor
28 
30  virtual ~SGModelEvaluatorBase() {}
31 
33  virtual void set_x_sg_init(const Stokhos::EpetraVectorOrthogPoly& x_sg_in) = 0;
34 
37 
39  virtual void set_p_sg_init(int i, const Stokhos::EpetraVectorOrthogPoly& p_sg_in) = 0;
40 
43 
45 
48  virtual Teuchos::Array<int> get_p_sg_map_indices() const = 0;
49 
51 
54  virtual Teuchos::Array<int> get_g_sg_map_indices() const = 0;
55 
58 
61 
64 
67 
71  const Epetra_Vector* v = NULL) const = 0;
72 
76  const Epetra_Vector* v = NULL) const = 0;
77 
80  create_x_mv_sg(int num_vecs,
81  Epetra_DataAccess CV = Copy,
82  const Epetra_MultiVector* v = NULL) const = 0;
83 
86  create_x_mv_sg_overlap(int num_vecs,
87  Epetra_DataAccess CV = Copy,
88  const Epetra_MultiVector* v = NULL) const = 0;
89 
92  create_p_sg(int l, Epetra_DataAccess CV = Copy,
93  const Epetra_Vector* v = NULL) const = 0;
94 
98  const Epetra_Vector* v = NULL) const = 0;
99 
103  const Epetra_Vector* v = NULL) const = 0;
104 
107  create_f_mv_sg(int num_vecs, Epetra_DataAccess CV = Copy,
108  const Epetra_MultiVector* v = NULL) const = 0;
109 
112  create_f_mv_sg_overlap(int num_vecs, Epetra_DataAccess CV = Copy,
113  const Epetra_MultiVector* v = NULL) const = 0;
114 
117  create_g_sg(int l, Epetra_DataAccess CV = Copy,
118  const Epetra_Vector* v = NULL) const = 0;
119 
122  create_g_mv_sg(int l, int num_vecs, Epetra_DataAccess CV = Copy,
123  const Epetra_MultiVector* v = NULL) const = 0;
124 
125  };
126 
127 }
128 
129 #endif // STOKHOS_SGMODELEVALUATORBASE_HPP
virtual Teuchos::Array< Teuchos::RCP< const Epetra_Map > > get_g_sg_base_maps() const =0
Get base maps of SG responses.
virtual Teuchos::RCP< const Epetra_BlockMap > get_overlap_stochastic_map() const =0
Return overlap stochastic map.
virtual Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > create_p_sg(int l, Epetra_DataAccess CV=Copy, const Epetra_Vector *v=NULL) const =0
Create vector orthog poly using p map.
virtual Teuchos::RCP< Stokhos::EpetraMultiVectorOrthogPoly > create_g_mv_sg(int l, int num_vecs, Epetra_DataAccess CV=Copy, const Epetra_MultiVector *v=NULL) const =0
Create multi-vector orthog poly using g map.
virtual Teuchos::Array< int > get_g_sg_map_indices() const =0
Get indices of SG responses.
virtual Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > create_g_sg(int l, Epetra_DataAccess CV=Copy, const Epetra_Vector *v=NULL) const =0
Create vector orthog poly using g map.
virtual Teuchos::Array< int > get_p_sg_map_indices() const =0
Get indices of SG parameters.
virtual Teuchos::RCP< const Epetra_BlockMap > get_x_sg_overlap_map() const =0
Return x sg overlap 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.
virtual Teuchos::RCP< Stokhos::EpetraMultiVectorOrthogPoly > create_f_mv_sg(int num_vecs, Epetra_DataAccess CV=Copy, const Epetra_MultiVector *v=NULL) const =0
Create multi-vector orthog poly using f map and owned sg map.
virtual Teuchos::RCP< Stokhos::EpetraMultiVectorOrthogPoly > create_x_mv_sg_overlap(int num_vecs, Epetra_DataAccess CV=Copy, const Epetra_MultiVector *v=NULL) const =0
Create vector orthog poly using x map and overlap sg map.
virtual Teuchos::RCP< Stokhos::EpetraMultiVectorOrthogPoly > create_x_mv_sg(int num_vecs, Epetra_DataAccess CV=Copy, const Epetra_MultiVector *v=NULL) const =0
Create vector orthog poly using x map and owned sg map.
virtual void set_x_sg_init(const Stokhos::EpetraVectorOrthogPoly &x_sg_in)=0
Set initial solution polynomial.
virtual Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > create_f_sg(Epetra_DataAccess CV=Copy, const Epetra_Vector *v=NULL) const =0
Create vector orthog poly using f map and owned sg map.
virtual void set_p_sg_init(int i, const Stokhos::EpetraVectorOrthogPoly &p_sg_in)=0
Set initial parameter polynomial.
virtual Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > create_f_sg_overlap(Epetra_DataAccess CV=Copy, const Epetra_Vector *v=NULL) const =0
Create vector orthog poly using f map and overlap sg map.
virtual Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > create_x_sg_overlap(Epetra_DataAccess CV=Copy, const Epetra_Vector *v=NULL) const =0
Create vector orthog poly using x map and overlap sg map.
virtual Teuchos::RCP< const Stokhos::EpetraVectorOrthogPoly > get_x_sg_init() const =0
Return initial SG x.
virtual Teuchos::RCP< Stokhos::EpetraMultiVectorOrthogPoly > create_f_mv_sg_overlap(int num_vecs, Epetra_DataAccess CV=Copy, const Epetra_MultiVector *v=NULL) const =0
Create multi-vector orthog poly using f map and overlap sg map.
Base class for stochastic Galerkin model evaluators.
virtual Teuchos::RCP< const Epetra_Import > get_x_sg_importer() const =0
Return x sg importer.
Epetra_DataAccess
virtual Teuchos::RCP< Stokhos::EpetraVectorOrthogPoly > create_x_sg(Epetra_DataAccess CV=Copy, const Epetra_Vector *v=NULL) const =0
Create vector orthog poly using x map and owned sg map.
virtual Teuchos::RCP< const Stokhos::EpetraVectorOrthogPoly > get_p_sg_init(int l) const =0
Return initial SG parameters.