|
EpetraExt
Development
|
A simple serial example model which includes a parameter subvector and a response function that can be used to define an optimization problem. More...
#include <EpetraModelEval4DOpt.hpp>

Public Member Functions | |
| EpetraModelEval4DOpt (const double xt0=1.0, const double xt1=1.0, const double pt0=2.0, const double pt1=0.0, const double d=10.0, const double x00=1.0, const double x01=1.0, const double p00=2.0, const double p01=0.0) | |
| | |
| void | setSupportDerivs (bool supportDerivs) |
| | |
| void | set_p_bounds (double pL0, double pL1, double pU0, double pU1) |
| | |
| void | set_x_bounds (double xL0, double xL1, double xU0, double xU1) |
| | |
Overridden from EpetraExt::ModelEvaluator . | |
| Teuchos::RefCountPtr< const Epetra_Map > | get_x_map () const |
| | |
| Teuchos::RefCountPtr< const Epetra_Map > | get_f_map () const |
| | |
| Teuchos::RefCountPtr< const Epetra_Map > | get_p_map (int l) const |
| . | |
| Teuchos::RefCountPtr< const Epetra_Map > | get_g_map (int j) const |
| . | |
| Teuchos::RefCountPtr< const Epetra_Vector > | get_x_init () const |
| | |
| Teuchos::RefCountPtr< const Epetra_Vector > | get_p_init (int l) const |
| | |
| Teuchos::RefCountPtr< const Epetra_Vector > | get_x_lower_bounds () const |
| | |
| Teuchos::RefCountPtr< const Epetra_Vector > | get_x_upper_bounds () const |
| | |
| Teuchos::RefCountPtr< const Epetra_Vector > | get_p_lower_bounds (int l) const |
| | |
| Teuchos::RefCountPtr< const Epetra_Vector > | get_p_upper_bounds (int l) const |
| | |
| Teuchos::RefCountPtr < Epetra_Operator > | create_W () const |
| | |
| InArgs | createInArgs () const |
| | |
| OutArgs | createOutArgs () const |
| | |
| void | evalModel (const InArgs &inArgs, const OutArgs &outArgs) const |
| | |
A simple serial example model which includes a parameter subvector and a response function that can be used to define an optimization problem.
Represents the model:
f[0] = x[0] + x[1]*x[1] - p[0];
f[1] = d_ * ( x[0]*x[0] - x[1] - p[1] );
g[0] = 0.5 * ( sqr(x[0]-xt0_) + sqr(x[1]-xt1_) + sqr(p[0]-pt0_) + sqr(p[1]-pt1_) );
where there is just one state vector x = [ x[0], x[1] ] and one parameter subvector p = [ p[0], p[1] ].
See the function evalModel() for more details.
Definition at line 72 of file EpetraModelEval4DOpt.hpp.
| EpetraModelEval4DOpt::EpetraModelEval4DOpt | ( | const double | xt0 = 1.0, |
| const double | xt1 = 1.0, |
||
| const double | pt0 = 2.0, |
||
| const double | pt1 = 0.0, |
||
| const double | d = 10.0, |
||
| const double | x00 = 1.0, |
||
| const double | x01 = 1.0, |
||
| const double | p00 = 2.0, |
||
| const double | p01 = 0.0 |
||
| ) |
Definition at line 61 of file EpetraModelEval4DOpt.cpp.
| void EpetraModelEval4DOpt::setSupportDerivs | ( | bool | supportDerivs | ) |
Definition at line 111 of file EpetraModelEval4DOpt.cpp.
| void EpetraModelEval4DOpt::set_p_bounds | ( | double | pL0, |
| double | pL1, | ||
| double | pU0, | ||
| double | pU1 | ||
| ) |
Definition at line 117 of file EpetraModelEval4DOpt.cpp.
| void EpetraModelEval4DOpt::set_x_bounds | ( | double | xL0, |
| double | xL1, | ||
| double | xU0, | ||
| double | xU1 | ||
| ) |
Definition at line 127 of file EpetraModelEval4DOpt.cpp.
| Teuchos::RefCountPtr< const Epetra_Map > EpetraModelEval4DOpt::get_x_map | ( | ) | const [virtual] |
Implements EpetraExt::ModelEvaluator.
Definition at line 140 of file EpetraModelEval4DOpt.cpp.
| Teuchos::RefCountPtr< const Epetra_Map > EpetraModelEval4DOpt::get_f_map | ( | ) | const [virtual] |
Implements EpetraExt::ModelEvaluator.
Definition at line 146 of file EpetraModelEval4DOpt.cpp.
| Teuchos::RefCountPtr< const Epetra_Map > EpetraModelEval4DOpt::get_p_map | ( | int | l | ) | const [virtual] |
.
Reimplemented from EpetraExt::ModelEvaluator.
Definition at line 152 of file EpetraModelEval4DOpt.cpp.
| Teuchos::RefCountPtr< const Epetra_Map > EpetraModelEval4DOpt::get_g_map | ( | int | j | ) | const [virtual] |
.
Reimplemented from EpetraExt::ModelEvaluator.
Definition at line 159 of file EpetraModelEval4DOpt.cpp.
| Teuchos::RefCountPtr< const Epetra_Vector > EpetraModelEval4DOpt::get_x_init | ( | ) | const [virtual] |
Reimplemented from EpetraExt::ModelEvaluator.
Definition at line 166 of file EpetraModelEval4DOpt.cpp.
| Teuchos::RefCountPtr< const Epetra_Vector > EpetraModelEval4DOpt::get_p_init | ( | int | l | ) | const [virtual] |
Reimplemented from EpetraExt::ModelEvaluator.
Definition at line 172 of file EpetraModelEval4DOpt.cpp.
| Teuchos::RefCountPtr< const Epetra_Vector > EpetraModelEval4DOpt::get_x_lower_bounds | ( | ) | const [virtual] |
Reimplemented from EpetraExt::ModelEvaluator.
Definition at line 179 of file EpetraModelEval4DOpt.cpp.
| Teuchos::RefCountPtr< const Epetra_Vector > EpetraModelEval4DOpt::get_x_upper_bounds | ( | ) | const [virtual] |
Reimplemented from EpetraExt::ModelEvaluator.
Definition at line 185 of file EpetraModelEval4DOpt.cpp.
| Teuchos::RefCountPtr< const Epetra_Vector > EpetraModelEval4DOpt::get_p_lower_bounds | ( | int | l | ) | const [virtual] |
Reimplemented from EpetraExt::ModelEvaluator.
Definition at line 191 of file EpetraModelEval4DOpt.cpp.
| Teuchos::RefCountPtr< const Epetra_Vector > EpetraModelEval4DOpt::get_p_upper_bounds | ( | int | l | ) | const [virtual] |
Reimplemented from EpetraExt::ModelEvaluator.
Definition at line 198 of file EpetraModelEval4DOpt.cpp.
| Teuchos::RefCountPtr< Epetra_Operator > EpetraModelEval4DOpt::create_W | ( | ) | const [virtual] |
Reimplemented from EpetraExt::ModelEvaluator.
Definition at line 205 of file EpetraModelEval4DOpt.cpp.
| EpetraExt::ModelEvaluator::InArgs EpetraModelEval4DOpt::createInArgs | ( | ) | const [virtual] |
Implements EpetraExt::ModelEvaluator.
Definition at line 211 of file EpetraModelEval4DOpt.cpp.
| EpetraExt::ModelEvaluator::OutArgs EpetraModelEval4DOpt::createOutArgs | ( | ) | const [virtual] |
Implements EpetraExt::ModelEvaluator.
Definition at line 221 of file EpetraModelEval4DOpt.cpp.
| void EpetraModelEval4DOpt::evalModel | ( | const InArgs & | inArgs, |
| const OutArgs & | outArgs | ||
| ) | const [virtual] |
Implements EpetraExt::ModelEvaluator.
Definition at line 265 of file EpetraModelEval4DOpt.cpp.
1.7.6.1