10#ifndef OPTIPACK_DIAGONAL_QUADRATIC_RESPONSE_ONLY_MODEL_EVALUATOR_DEF_HPP
11#define OPTIPACK_DIAGONAL_QUADRATIC_RESPONSE_ONLY_MODEL_EVALUATOR_DEF_HPP
14#include "Thyra_DiagonalQuadraticResponseOnlyModelEvaluator_decl.hpp"
15#include "Thyra_DiagonalScalarProd.hpp"
16#include "Thyra_VectorStdOps.hpp"
17#include "Thyra_DefaultSpmdVectorSpace.hpp"
18#include "Thyra_DetachedSpmdVectorView.hpp"
19#include "Teuchos_DefaultComm.hpp"
20#include "Teuchos_CommHelpers.hpp"
21#include "Teuchos_Assert.hpp"
40 :Np_(1), Ng_(1), comm_(comm), localDim_(localDim),
41 nonlinearTermFactor_(0.0), g_offset_(0.0)
67 V_S(diag.
ptr(), ST::one());
73 V_S(s_bar.
ptr(), ST::one());
77 g_offset_ = ST::zero();
86 ps_ = ps.assert_not_null();
107template<
class Scalar>
113 using Teuchos::rcp_dynamic_cast;
118 diag_bar_ = diag_bar.assert_not_null();
128 V_S( s_bar.
ptr(), ST::zero() );
133 rcp_dynamic_cast<Thyra::ScalarProdVectorSpaceBase<Scalar> >(p_space_,
true);
139template<
class Scalar>
141 const Scalar &nonlinearTermFactor)
143 nonlinearTermFactor_ = nonlinearTermFactor;
147template<
class Scalar>
149 const Scalar &g_offset)
151 g_offset_ = g_offset;
158template<
class Scalar>
165template<
class Scalar>
172template<
class Scalar>
185template<
class Scalar>
198template<
class Scalar>
203 MEB::InArgsSetup<Scalar> inArgs;
204 inArgs.setModelEvalDescription(this->
description());
213template<
class Scalar>
215DiagonalQuadraticResponseOnlyModelEvaluator<Scalar>::createOutArgsImpl()
const
218 MEB::OutArgsSetup<Scalar> outArgs;
219 outArgs.setModelEvalDescription(this->description());
220 outArgs.set_Np_Ng(Np_,Ng_);
221 outArgs.setSupports(MEB::OUT_ARG_DgDp, 0 ,0, MEB::DERIV_TRANS_MV_BY_ROW);
226template<
class Scalar>
227void DiagonalQuadraticResponseOnlyModelEvaluator<Scalar>::evalModelImpl(
234 using Teuchos::outArg;
241 const ConstDetachedSpmdVectorView<Scalar> p(inArgs.
get_p(0));
242 const ConstDetachedSpmdVectorView<Scalar> ps(ps_);
243 const ConstDetachedSpmdVectorView<Scalar> diag(diag_);
244 const ConstDetachedSpmdVectorView<Scalar> s_bar(s_bar_);
247 if (!is_null(outArgs.
get_g(0))) {
248 Scalar g_val = ST::zero();
249 for (Ordinal i = 0; i < p.subDim(); ++i) {
250 const Scalar p_ps = p[i] - ps[i];
251 g_val += diag[i] * p_ps*p_ps;
252 if (nonlinearTermFactor_ != ST::zero()) {
253 g_val += nonlinearTermFactor_ * p_ps * p_ps * p_ps;
257 Teuchos::reduceAll<Ordinal, Scalar>(*comm_, Teuchos::REDUCE_SUM,
258 g_val, outArg(global_g_val) );
260 as<Scalar>(0.5) * global_g_val + g_offset_;
264 if (!outArgs.
get_DgDp(0,0).isEmpty()) {
266 get_mv<Scalar>(outArgs.
get_DgDp(0,0),
"DgDp^T", MEB::DERIV_TRANS_MV_BY_ROW);
268 for (Thyra::Ordinal i = 0; i < p.subDim(); ++i) {
269 const Scalar p_ps = p[i] - ps[i];
270 Scalar DgDp_grad_i = diag[i] * p_ps;
271 if (nonlinearTermFactor_ != ST::zero()) {
272 DgDp_grad_i += as<Scalar>(1.5) * nonlinearTermFactor_ * p_ps * p_ps;
274 DgDp_grad[i] = DgDp_grad_i / s_bar[i];
static Teuchos::RCP< const Comm< OrdinalType > > getComm()
virtual std::string description() const
Create an explicit detached non-mutable (const) view of all of the local elements on this process of ...
RCP< DefaultSpmdVectorSpace< Scalar > > defaultSpmdVectorSpace()
Nonmember consturctor that creats an uninitialized vector space.
RCP< DefaultSpmdVectorSpace< Scalar > > locallyReplicatedDefaultSpmdVectorSpace(const RCP< const Teuchos::Comm< Ordinal > > &comm, const Ordinal globalDim)
Nonmember consturctor function that creates a locally-replicated parallel vector space.
Create an explicit detached mutable (non-const) view of all of the local elements on this process of ...
void setDiagonalVector(const RCP< const VectorBase< Scalar > > &diag)
Set the diagonal vector diag.
void setSolutionVector(const RCP< const VectorBase< Scalar > > &ps)
Set the solution vector ps .
ModelEvaluatorBase::InArgs< Scalar > createInArgs() const
void setScalarOffset(const Scalar &g_offset)
Set offset scalar g_offset .
RCP< const VectorSpaceBase< Scalar > > get_g_space(int j) const
void setNonlinearTermFactor(const Scalar &nonlinearTermFactor)
Set nonlinear term factory.
DiagonalQuadraticResponseOnlyModelEvaluator(const int localDim, const RCP< const Teuchos::Comm< Ordinal > > &comm=Teuchos::null)
RCP< const VectorSpaceBase< Scalar > > get_p_space(int l) const
void setDiagonalBarVector(const RCP< const VectorBase< Scalar > > &diag_bar)
Set the diagonal vector diag_bar.
const RCP< const VectorBase< Scalar > > getSolutionVector() const
Get the solution vector ps .
RCP< DiagonalScalarProd< Scalar > > diagonalScalarProd(const RCP< const VectorBase< Scalar > > &s_diag)
Nonmember constructor.
Concrete aggregate class for all input arguments computable by a ModelEvaluator subclass object.
RCP< const VectorBase< Scalar > > get_p(int l) const
Get p(l) where 0 <= l && l < this->Np().
Concrete aggregate class for all output arguments computable by a ModelEvaluator subclass object.
Derivative< Scalar > get_DgDp(int j, int l) const
Precondition: supports(OUT_ARG_DgDp,j,l)==true.
Evaluation< VectorBase< Scalar > > get_g(int j) const
Precondition: supports(OUT_ARG_g)==true..
Base subclass for ModelEvaluator that defines some basic types.
Abstract interface for finite-dimensional dense vectors.
void ele_wise_divide(const Scalar &alpha, const VectorBase< Scalar > &x, const VectorBase< Scalar > &v, const Ptr< VectorBase< Scalar > > &y)
Element-wise division update: y(i) += alpha * x(i) / v(i), i = 0...y->space()->dim()-1.
void V_S(const Ptr< VectorBase< Scalar > > &y, const Scalar &alpha)
y(i) = alpha, i = 0...y->space()->dim()-1.
RCP< VectorBase< Scalar > > createMember(const RCP< const VectorSpaceBase< Scalar > > &vs, const std::string &label="")
Create a vector member from the vector space.
#define TEUCHOS_ASSERT(assertion_test)
#define TEUCHOS_ASSERT_IN_RANGE_UPPER_EXCLUSIVE(index, lower_inclusive, upper_exclusive)
TypeTo as(const TypeFrom &t)