44#ifndef GLP_APP_ADV_DIFF_REACT_OPT_MODEL_HPP
45#define GLP_APP_ADV_DIFF_REACT_OPT_MODEL_HPP
53#include "Teuchos_VerboseObject.hpp"
54#include "Teuchos_Array.hpp"
164 ,
public Teuchos::VerboseObject<AdvDiffReactOptModel>
170 const Teuchos::RCP<const Epetra_Comm> &comm
176 ,
const char meshFile[]
180 ,
const double reactionRate
181 ,
const bool normalizeBasis
182 ,
const bool supportDerivatives
186 void set_q( Teuchos::RCP<const Epetra_Vector>
const& q );
189 Teuchos::RCP<GLpApp::GLpYUEpetraDataPool>
getDataPool();
192 Teuchos::RCP<const Epetra_MultiVector>
get_B_bar()
const;
198 Teuchos::RCP<const Epetra_Map>
get_x_map()
const;
200 Teuchos::RCP<const Epetra_Map>
get_f_map()
const;
202 Teuchos::RCP<const Epetra_Map>
get_p_map(
int l)
const;
204 Teuchos::RCP<const Epetra_Map>
get_g_map(
int j)
const;
206 Teuchos::RCP<const Epetra_Vector>
get_x_init()
const;
208 Teuchos::RCP<const Epetra_Vector>
get_p_init(
int l)
const;
218 Teuchos::RCP<Epetra_Operator>
create_W()
const;
249 Teuchos::RCP<GLpApp::GLpYUEpetraDataPool>
dat_;
251 Teuchos::RCP<const Epetra_Vector>
q_;
262 Teuchos::RCP<Epetra_Vector>
x0_;
263 Teuchos::RCP<Epetra_Vector>
xL_;
264 Teuchos::RCP<Epetra_Vector>
xU_;
268 Teuchos::RCP<Epetra_Vector>
gL_;
269 Teuchos::RCP<Epetra_Vector>
gU_;
Base interface for evaluating a stateless "model".
PDE-constrained inverse problem based on a 2D discretization of a diffusion/reaction system.
Teuchos::RCP< Epetra_Operator > create_W() const
Teuchos::RCP< Epetra_Vector > xU_
Teuchos::RCP< const Epetra_Map > get_p_map(int l) const
\breif .
InArgs createInArgs() const
Teuchos::RCP< Epetra_MultiVector > B_bar_
Teuchos::RCP< const Epetra_Map > map_x_
static const int p_rx_idx
RCP_Eptra_Map_Array_t map_p_
Teuchos::RCP< const Epetra_Map > map_p_bar_
Teuchos::RCP< const Epetra_Vector > get_x_init() const
Teuchos::RCP< const Epetra_Vector > get_p_upper_bounds(int l) const
Teuchos::RCP< Epetra_Vector > xL_
RCP_Eptra_Vector_Array_t p0_
Teuchos::RCP< GLpApp::GLpYUEpetraDataPool > dat_
Teuchos::RCP< const Epetra_Map > get_g_map(int j) const
\breif .
RCP_Eptra_Vector_Array_t pU_
Teuchos::RCP< const Epetra_Vector > q_
Teuchos::RCP< const Epetra_Map > map_f_
Teuchos::RCP< Epetra_Vector > gU_
Teuchos::RCP< Epetra_CrsGraph > W_graph_
Teuchos::RCP< GLpApp::GLpYUEpetraDataPool > getDataPool()
Teuchos::RCP< const Epetra_Vector > get_x_upper_bounds() const
Teuchos::RCP< const Epetra_Vector > get_p_lower_bounds(int l) const
Teuchos::RCP< Epetra_Vector > gL_
Teuchos::RCP< Epetra_Operator > create_DfDp_op(int l) const
Teuchos::RCP< const Epetra_Comm > epetra_comm_
Teuchos::RCP< const Epetra_MultiVector > get_B_bar() const
OutArgs createOutArgs() const
Teuchos::RCP< const Epetra_Vector > get_p_init(int l) const
Teuchos::Array< Teuchos::RCP< Epetra_Vector > > RCP_Eptra_Vector_Array_t
Teuchos::Array< Teuchos::RCP< const Epetra_Map > > RCP_Eptra_Map_Array_t
Teuchos::RCP< const Epetra_Map > get_f_map() const
void set_q(Teuchos::RCP< const Epetra_Vector > const &q)
Teuchos::RCP< const Epetra_Map > map_g_
void evalModel(const InArgs &inArgs, const OutArgs &outArgs) const
Teuchos::RCP< const Epetra_Vector > get_x_lower_bounds() const
static const int p_bndy_idx
Teuchos::RCP< const Epetra_Map > get_x_map() const
RCP_Eptra_Vector_Array_t pL_
Teuchos::RCP< Epetra_Vector > x0_