Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_SGQuadModelEvaluator.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_SGQUADMODELEVALUATOR_HPP
11 #define STOKHOS_SGQUADMODELEVALUATOR_HPP
12 
13 #include "EpetraExt_ModelEvaluator.h"
14 
15 #include "Teuchos_RCP.hpp"
16 #include "Teuchos_Array.hpp"
17 
18 namespace Stokhos {
19 
31  class SGQuadModelEvaluator : public EpetraExt::ModelEvaluator {
32  public:
33 
34  // Constructor
37 
40 
43 
46 
49 
52 
55  get_p_names(int l) const;
56 
59 
62 
65 
67  InArgs createInArgs() const;
68 
70  OutArgs createOutArgs() const;
71 
73  void evalModel(const InArgs& inArgs, const OutArgs& outArgs) const;
74 
76 
77  protected:
78 
81 
83  int num_p;
84 
86  int num_g;
87 
90 
93 
96 
99 
102 
105 
108 
111 
114 
117 
118  };
119 
120 }
121 
122 #endif // STOKHOS_SGQUADMODELEVALUATOR_HPP
int num_p
Number of parameter vectors.
Teuchos::RCP< const Epetra_Map > get_x_map() const
Return solution vector map.
int num_g
Number of response vectors.
Teuchos::RCP< const Epetra_Vector > get_x_init() const
Return initial solution.
ModelEvaluator adaptor that implements the stochastic Galerkin residual and Jacobian computations usi...
Teuchos::Array< EpetraExt::ModelEvaluator::Derivative > dgdx_dot_qp
Response derivative.
Teuchos::Array< EpetraExt::ModelEvaluator::Derivative > dgdx_qp
Response derivative.
Teuchos::RCP< const Teuchos::Array< std::string > > get_p_names(int l) const
Return array of parameter names.
InArgs createInArgs() const
Create InArgs.
Teuchos::Array< Teuchos::RCP< Epetra_Vector > > g_qp
Response vectors.
Teuchos::Array< Teuchos::RCP< Epetra_Vector > > p_qp
Parameter vectors.
Teuchos::RCP< const Epetra_Map > get_p_map(int l) const
Return parameter vector map.
Teuchos::Array< EpetraExt::ModelEvaluator::Derivative > dfdp_qp
Residual derivatives.
Teuchos::RCP< Epetra_Vector > x_dot_qp
Time derivative vector.
Teuchos::RCP< Epetra_Operator > create_W() const
Create W = alpha*M + beta*J matrix.
Teuchos::RCP< Epetra_Vector > f_qp
Residual vector.
Teuchos::RCP< Epetra_Vector > x_qp
Solution vector.
Teuchos::RCP< Epetra_Operator > W_qp
W operator.
Teuchos::Array< Teuchos::Array< EpetraExt::ModelEvaluator::Derivative > > dgdp_qp
Response sensitivities.
Teuchos::RCP< EpetraExt::ModelEvaluator > me
Underlying model evaluator.
Teuchos::RCP< const Epetra_Vector > get_p_init(int l) const
Return initial parameters.
Teuchos::RCP< const Epetra_Map > get_f_map() const
Return residual vector map.
OutArgs createOutArgs() const
Create OutArgs.
SGQuadModelEvaluator(const Teuchos::RCP< EpetraExt::ModelEvaluator > &me)
void evalModel(const InArgs &inArgs, const OutArgs &outArgs) const
Evaluate model on InArgs.
Teuchos::RCP< const Epetra_Map > get_g_map(int l) const
Return observation vector map.