Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
epetra/SimpleME.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 SIMPLE_ME_HPP
11 #define SIMPLE_ME_HPP
12 
13 #include "Teuchos_RCP.hpp"
14 #include "Teuchos_Array.hpp"
15 
16 #include "EpetraExt_ModelEvaluator.h"
17 
18 #include "Epetra_Map.h"
19 #include "Epetra_LocalMap.h"
20 #include "Epetra_Import.h"
21 #include "Epetra_CrsGraph.h"
22 #include "Epetra_CrsMatrix.h"
23 
24 /* Simple model evaluator demonstrating how to use the Stokhos
25  * ModelEvaluator to solve a nonlinear stochastic Galerkin problem.
26  * It represents the simple function
27  *
28  * f(x,a) = | a^2 - x_0 |
29  * | x_1^2 - x_0 |
30  *
31  * where x = [x_0 x_1]^T and a is a parameter. It has a root at x = [ a^2, a].
32  * The parameter a may be represented by a given polynomial chaos expansion,
33  * and the corresponding expansion for x computed by Stokhos.
34  */
35 class SimpleME : public EpetraExt::ModelEvaluator {
36 public:
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 
84 
87 
90 
93 
96 
99 
102 
105 
106 };
107 
108 #endif // SIMPLE_ME_HPP
Teuchos::RCP< const Epetra_Vector > get_x_init() const
Return initial solution.
Teuchos::RCP< Epetra_Vector > x_init
Initial guess.
Teuchos::RCP< const Epetra_Map > get_f_map() const
Return residual vector map.
SimpleME(const Teuchos::RCP< const Epetra_Comm > &comm)
Constructor.
Teuchos::RCP< const Epetra_Map > get_x_map() const
Return solution vector map.
Teuchos::RCP< Epetra_Vector > p_init
Initial parameters.
InArgs createInArgs() const
Create InArgs.
Teuchos::RCP< Epetra_Vector > x_overlapped
Overlapped solution vector.
Teuchos::RCP< Epetra_Map > x_map
Solution vector map.
void evalModel(const InArgs &inArgs, const OutArgs &outArgs) const
Evaluate model on InArgs.
Teuchos::RCP< Epetra_Map > p_map
Parameter vector map.
Teuchos::RCP< Teuchos::Array< std::string > > p_names
Parameter names.
Teuchos::RCP< Epetra_CrsGraph > graph
Jacobian graph.
Teuchos::RCP< const Epetra_Map > get_p_map(int l) const
Return parameter vector map.
Teuchos::RCP< const Teuchos::Array< std::string > > get_p_names(int l) const
Return array of parameter names.
Teuchos::RCP< Epetra_Import > importer
Importer to overlapped distribution.
OutArgs createOutArgs() const
Create OutArgs.
Teuchos::RCP< Epetra_Operator > create_W() const
Create W = alpha*M + beta*J matrix.
Teuchos::RCP< Epetra_Map > x_overlapped_map
Overlapped solution vector map.
Teuchos::RCP< const Epetra_Vector > get_p_init(int l) const
Return initial parameters.