45 #include "Teuchos_ScalarTraits.hpp"
46 #include "Epetra_SerialComm.h"
47 #include "Epetra_CrsMatrix.h"
49 #include "Epetra_MpiComm.h"
53 inline double sqr(
const double& s ) {
return s*s; }
58 Teuchos::RCP<Epetra_Comm> epetra_comm
70 :isInitialized_(false), epetra_comm_(epetra_comm),
71 xt0_(xt0),xt1_(xt1),pt0_(pt0),pt1_(pt1),d_(d)
75 const int nx = 2, np = 2, ng = 1, nq = 1;
77 map_x_ = rcp(
new Epetra_Map(nx,0,*epetra_comm_));
78 map_p_ = rcp(
new Epetra_Map(np,0,*epetra_comm_));
79 map_q_ = rcp(
new Epetra_Map(nq,0,*epetra_comm_));
80 map_g_ = rcp(
new Epetra_Map(ng,0,*epetra_comm_));
83 const double inf = 1e+50;
87 x0_ = rcp(
new Epetra_Vector(*map_x_)); (*x0_)[0] = x00; (*x0_)[1] = x01;
90 p0_ = rcp(
new Epetra_Vector(*map_p_)); (*p0_)[0] = p00; (*p0_)[1] = p01;
100 int indices[nx] = { 0, 1 };
101 for(
int i = 0; i < nx; ++i )
102 W_graph_->InsertGlobalIndices(i,nx,indices);
104 W_graph_->FillComplete();
106 isInitialized_ =
true;
111 double pL0,
double pL1,
double pU0,
double pU1
121 double xL0,
double xL1,
double xU0,
double xU1
132 Teuchos::RCP<const Epetra_Map>
138 Teuchos::RCP<const Epetra_Map>
144 Teuchos::RCP<const Epetra_Map>
147 TEUCHOS_TEST_FOR_EXCEPT(l>1);
148 if (l==0)
return map_p_;
152 Teuchos::RCP<const Epetra_Map>
155 TEUCHOS_TEST_FOR_EXCEPT(j!=0);
159 Teuchos::RCP<const Epetra_Vector>
165 Teuchos::RCP<const Epetra_Vector>
168 TEUCHOS_TEST_FOR_EXCEPT(l>1);
169 if (l==0)
return p0_;
173 Teuchos::RCP<const Epetra_Vector>
179 Teuchos::RCP<const Epetra_Vector>
185 Teuchos::RCP<const Epetra_Vector>
188 TEUCHOS_TEST_FOR_EXCEPT(l>1);
189 if (l==0)
return pL_;
193 Teuchos::RCP<const Epetra_Vector>
196 TEUCHOS_TEST_FOR_EXCEPT(l>1);
197 if (l==0)
return pU_;
201 Teuchos::RCP<Epetra_Operator>
261 using Teuchos::dyn_cast;
262 using Teuchos::rcp_dynamic_cast;
266 Teuchos::RCP<const Epetra_Vector> p_in = inArgs.
get_p(0);
267 Teuchos::RCP<const Epetra_Vector> q_in = inArgs.
get_p(1);
287 f[0] = x[0] + x[1]*x[1] - p[0] + q[0];
288 f[1] = d_ * ( x[0]*x[0] - x[1] - p[1] );
292 g[0] = 0.5 * ( sqr(x[0]-xt0_) + sqr(x[1]-xt1_) + sqr(p[0]-pt0_) + sqr(p[1]-pt1_) );
306 values[0] = 1.0; indexes[0] = 0;
307 values[1] = 2.0*x[1]; indexes[1] = 1;
310 values[0] = 2.0*d_*x[0]; indexes[0] = 0;
311 values[1] = -d_; indexes[1] = 1;
322 DgDx_trans[0] = x[0]-xt0_;
323 DgDx_trans[1] = x[1]-xt1_;
327 DgDp_trans[0] = p[0]-pt0_;
328 DgDp_trans[1] = p[1]-pt1_;
InArgs createInArgs() const
void set_W_properties(const DerivativeProperties &properties)
Evaluation< Epetra_Vector > get_g(int j) const
Get g(j) where 0 <= j && j < this->Ng().
void setSupports(EOutArgsMembers arg, bool supports=true)
Teuchos::RCP< const Epetra_Vector > get_p_init(int l) const
virtual int SumIntoGlobalValues(int GlobalRow, int NumEntries, const double *Values, const int *Indices)
Teuchos::RCP< const Epetra_Vector > get_p_upper_bounds(int l) const
Teuchos::RCP< const Epetra_Map > get_x_map() const
Teuchos::RCP< Epetra_MultiVector > get_DgDp_mv(const int j, const int l, const ModelEvaluator::OutArgs &outArgs, const ModelEvaluator::EDerivativeMultiVectorOrientation mvOrientation)
Teuchos::RCP< Epetra_Operator > get_W() const
void setSupports(EInArgsMembers arg, bool supports=true)
Teuchos::RCP< Epetra_Operator > create_W() const
Teuchos::RCP< const Epetra_Vector > get_p(int l) const
Teuchos::RCP< const Epetra_Map > get_p_map(int l) const
.
Teuchos::RCP< const Epetra_Map > get_f_map() const
Teuchos::RCP< const Epetra_Vector > get_p_lower_bounds(int l) const
Teuchos::RCP< const Epetra_Map > get_g_map(int j) const
.
int PutScalar(double ScalarConstant)
void set_DgDp_properties(int j, int l, const DerivativeProperties &properties)
void set_DfDp_properties(int l, const DerivativeProperties &properties)
void setModelEvalDescription(const std::string &modelEvalDescription)
Teuchos::RCP< Epetra_MultiVector > get_DfDp_mv(const int l, const ModelEvaluator::OutArgs &outArgs)
Teuchos::RCP< Epetra_MultiVector > get_DgDx_mv(const int j, const ModelEvaluator::OutArgs &outArgs, const ModelEvaluator::EDerivativeMultiVectorOrientation mvOrientation)
void setModelEvalDescription(const std::string &modelEvalDescription)
void set_p_bounds(double pL0, double pL1, double pU0, double pU1)
Teuchos::RCP< const Epetra_Vector > get_x_upper_bounds() const
EpetraMultiPointModelEval4DOpt(Teuchos::RCP< Epetra_Comm > epetra_comm, 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, const double q0=0.0)
Teuchos::RCP< const Epetra_Vector > get_x_lower_bounds() const
void set_DgDx_properties(int j, const DerivativeProperties &properties)
void set_x_bounds(double xL0, double xL1, double xU0, double xU1)
void evalModel(const InArgs &inArgs, const OutArgs &outArgs) const
Teuchos::RCP< const Epetra_Vector > get_x_init() const
Evaluation< Epetra_Vector > get_f() const
Teuchos::RCP< const Epetra_Vector > get_x() const
Set solution vector Taylor polynomial.
void set_Np_Ng(int Np, int Ng)
OutArgs createOutArgs() const