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 //
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_SGQUADMODELEVALUATOR_HPP
43 #define STOKHOS_SGQUADMODELEVALUATOR_HPP
44 
45 #include "EpetraExt_ModelEvaluator.h"
46 
47 #include "Teuchos_RCP.hpp"
48 #include "Teuchos_Array.hpp"
49 
50 namespace Stokhos {
51 
63  class SGQuadModelEvaluator : public EpetraExt::ModelEvaluator {
64  public:
65 
66  // Constructor
69 
72 
75 
78 
81 
84 
87  get_p_names(int l) const;
88 
91 
94 
97 
99  InArgs createInArgs() const;
100 
102  OutArgs createOutArgs() const;
103 
105  void evalModel(const InArgs& inArgs, const OutArgs& outArgs) const;
106 
108 
109  protected:
110 
113 
115  int num_p;
116 
118  int num_g;
119 
122 
125 
128 
131 
134 
137 
140 
143 
146 
149 
150  };
151 
152 }
153 
154 #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.