67 const std::string fwdScalingVecName =
"fwdScalingVec";
75 void assertModelVarScalings(
81 (varScalings.
supports(EME::IN_ARG_x) && varScalings.
supports(EME::IN_ARG_x_dot))
84 "Error, if scaling for x is given and x_dot is supported, then\n"
85 "the scaling for x_dot must also be set and must be the same scaling\n"
89 (varScalings.
supports(EME::IN_ARG_x) && varScalings.
supports(EME::IN_ARG_x_dotdot))
92 "Error, if scaling for x is given and x_dotdot is supported, then\n"
93 "the scaling for x_dotdot must also be set and must be the same scaling\n"
98 #endif // TEUCHOS_DEBUG
103 template<
class InArgsVectorGetterSetter>
105 InArgsVectorGetterSetter vecGetterSetter,
118 using Teuchos::rcp_const_cast;
125 RCP<const Epetra_Vector>
126 orig_vec = vecGetterSetter.getVector(origVars);
128 RCP<const Epetra_Vector>
129 inv_s_vec = vecGetterSetter.getVector(varScalings);
133 vecGetterSetter.getVector(*scaledVars) );
137 Ptr<const RCP<const Epetra_Vector> > fwd_s_vec =
138 Teuchos::getOptionalEmbeddedObj<Epetra_Vector, RCP<const Epetra_Vector> >(
146 scaled_vec->Multiply( 1.0, **fwd_s_vec, *orig_vec, 0.0 );
151 *orig_vec, *inv_s_vec, &*scaled_vec );
153 vecGetterSetter.setVector( scaled_vec, scaledVars );
156 vecGetterSetter.setVector( orig_vec, scaledVars );
160 vecGetterSetter.setVector( null, scaledVars );
168 template<
class InArgsVectorGetterSetter>
169 void scaleModelBound(
170 InArgsVectorGetterSetter vecGetterSetter,
185 using Teuchos::rcp_const_cast;
192 RCP<const Epetra_Vector>
193 orig_lower_vec = vecGetterSetter.getVector(origLowerBounds);
194 if ( !
is_null(orig_lower_vec) ) {
195 RCP<const Epetra_Vector>
196 inv_s_vec = vecGetterSetter.getVector(varScalings);
201 vecGetterSetter.setVector( orig_lower_vec, scaledLowerBounds );
205 vecGetterSetter.setVector( null, scaledLowerBounds );
208 RCP<const Epetra_Vector>
209 orig_upper_vec = vecGetterSetter.getVector(origUpperBounds);
210 if ( !
is_null(orig_upper_vec) ) {
211 RCP<const Epetra_Vector>
212 inv_s_vec = vecGetterSetter.getVector(varScalings);
217 vecGetterSetter.setVector( orig_upper_vec, scaledUpperBounds );
221 vecGetterSetter.setVector( null, scaledUpperBounds );
229 template<
class InArgsVectorGetterSetter>
230 void unscaleModelVar(
231 InArgsVectorGetterSetter vecGetterSetter,
243 using Teuchos::rcp_const_cast;
251 RCP<const Epetra_Vector>
252 scaled_vec = vecGetterSetter.getVector(scaledVars);
254 RCP<const Epetra_Vector>
255 inv_s_vec = vecGetterSetter.getVector(varScalings);
259 vecGetterSetter.getVector(*origVars) );
263 *scaled_vec, *inv_s_vec, &*orig_vec );
265 *out <<
"\nUnscaled vector "<<vecGetterSetter.getName()<<
":\n";
267 orig_vec->Print(*out);
269 vecGetterSetter.setVector( orig_vec, origVars );
272 vecGetterSetter.setVector( scaled_vec, origVars );
276 vecGetterSetter.setVector( null, origVars );
284 template<
class OutArgsVectorGetterSetter>
286 OutArgsVectorGetterSetter vecGetterSetter,
296 func = vecGetterSetter.getVector(origFuncs);
299 funcScaling = vecGetterSetter.getVector(funcScalings);
304 vecGetterSetter.setVector( func, scaledFuncs );
326 if(nominalValues->
supports(EME::IN_ARG_x)) {
330 if(nominalValues->
supports(EME::IN_ARG_x_dot)) {
334 if(nominalValues->
supports(EME::IN_ARG_x_dotdot)) {
338 for(
int l = 0; l < nominalValues->
Np(); ++l ) {
342 if(nominalValues->
supports(EME::IN_ARG_t)) {
367 if(lowerBounds->
supports(EME::IN_ARG_x)) {
372 for(
int l = 0; l < lowerBounds->
Np(); ++l ) {
377 if(lowerBounds->
supports(EME::IN_ARG_t)) {
397 assertModelVarScalings(varScalings);
400 if (origVars.
supports(EME::IN_ARG_x)) {
405 if (origVars.
supports(EME::IN_ARG_x_dot)) {
410 if (origVars.
supports(EME::IN_ARG_x_dotdot)) {
415 const int Np = origVars.
Np();
416 for (
int l = 0; l < Np; ++l ) {
421 if (origVars.
supports(EME::IN_ARG_x_poly)) {
424 "Error, can't hanlde scaling of x_poly yet!"
429 if (origVars.
supports(EME::IN_ARG_x_dot_poly)) {
432 "Error, can't hanlde scaling of x_dot_poly yet!"
437 if (origVars.
supports(EME::IN_ARG_x_dotdot_poly)) {
440 "Error, can't hanlde scaling of x_dotdot_poly yet!"
445 if (origVars.
supports(EME::IN_ARG_t)) {
447 varScalings.
get_t() > 0.0, std::logic_error,
448 "Error, can't hanlde scaling of t yet!"
453 if (origVars.
supports(EME::IN_ARG_alpha)) {
455 varScalings.
get_alpha() > 0.0, std::logic_error,
456 "Error, can't hanlde scaling of alpha yet!"
461 if (origVars.
supports(EME::IN_ARG_beta)) {
463 varScalings.
get_beta() > 0.0, std::logic_error,
464 "Error, can't hanlde scaling of beta yet!"
491 assertModelVarScalings(varScalings);
494 if (origLowerBounds.
supports(EME::IN_ARG_x)) {
497 varScalings, scaledLowerBounds, scaledUpperBounds,
501 if (origLowerBounds.
supports(EME::IN_ARG_x_dot)) {
504 varScalings, scaledLowerBounds, scaledUpperBounds,
508 if (origLowerBounds.
supports(EME::IN_ARG_x_dotdot)) {
511 varScalings, scaledLowerBounds, scaledUpperBounds,
515 const int Np = origLowerBounds.
Np();
516 for (
int l = 0; l < Np; ++l ) {
519 varScalings, scaledLowerBounds, scaledUpperBounds,
543 assertModelVarScalings(varScalings);
549 RCP<const Epetra_Vector> inv_s_x;
550 if ( scaledVars.
supports(EME::IN_ARG_x) &&
553 *out <<
"\nState inverse scaling vector inv_s_x:\n";
555 inv_s_x->Print(*out);
561 if (scaledVars.
supports(EME::IN_ARG_x_dot)) {
566 if (scaledVars.
supports(EME::IN_ARG_x_dotdot)) {
571 if (scaledVars.
supports(EME::IN_ARG_x)) {
576 const int Np = scaledVars.
Np();
577 for (
int l = 0; l < Np; ++l ) {
582 if (scaledVars.
supports(EME::IN_ARG_x_dot_poly)) {
585 "Error, can't hanlde unscaling of x_dot_poly yet!"
590 if (scaledVars.
supports(EME::IN_ARG_x_dotdot_poly)) {
593 "Error, can't hanlde unscaling of x_dotdot_poly yet!"
598 if (scaledVars.
supports(EME::IN_ARG_x_poly)) {
601 "Error, can't hanlde unscaling of x_poly yet!"
606 if (scaledVars.
supports(EME::IN_ARG_t)) {
608 varScalings.
get_t() > 0.0, std::logic_error,
609 "Error, can't hanlde unscaling of t yet!"
614 if (scaledVars.
supports(EME::IN_ARG_alpha)) {
616 varScalings.
get_alpha() > 0.0, std::logic_error,
617 "Error, can't hanlde unscaling of alpha yet!"
622 if (scaledVars.
supports(EME::IN_ARG_beta)) {
624 varScalings.
get_beta() > 0.0, std::logic_error,
625 "Error, can't hanlde unscaling of beta yet!"
638 bool *allFuncsWhereScaled,
649 *allFuncsWhereScaled =
true;
651 const int Np = origFuncs.
Np();
652 const int Ng = origFuncs.
Ng();
662 origFuncs.
supports(EME::OUT_ARG_f_poly)
668 "Error, we can't handle scaling of f_poly yet!"
674 for (
int j = 0; j < Ng; ++j ) {
680 RCP<Epetra_Operator> W;
682 bool didScaling =
false;
688 scaledFuncs->
set_W(W);
690 *allFuncsWhereScaled =
false;
694 for (
int l = 0; l < Np; ++l ) {
695 EME::Derivative orig_DfDp_l;
697 !origFuncs.
supports(EME::OUT_ARG_DfDp,l).none()
698 && !(orig_DfDp_l=origFuncs.
get_DfDp(l)).isEmpty()
701 EME::Derivative scaled_DfDp_l;
702 bool didScaling =
false;
705 &scaled_DfDp_l, &didScaling
708 scaledFuncs->
set_DfDp(l,scaled_DfDp_l);
710 *allFuncsWhereScaled =
false;
716 for (
int j = 0; j < Ng; ++j ) {
718 EME::Derivative orig_DgDx_dot_j;
720 !origFuncs.
supports(EME::OUT_ARG_DgDx_dot,j).none()
721 && !(orig_DgDx_dot_j=origFuncs.
get_DgDx_dot(j)).isEmpty()
724 EME::Derivative scaled_DgDx_dot_j;
725 bool didScaling =
false;
728 &scaled_DgDx_dot_j, &didScaling
733 *allFuncsWhereScaled =
false;
736 EME::Derivative orig_DgDx_dotdot_j;
738 !origFuncs.
supports(EME::OUT_ARG_DgDx_dotdot,j).none()
742 EME::Derivative scaled_DgDx_dotdot_j;
743 bool didScaling =
false;
746 &scaled_DgDx_dotdot_j, &didScaling
751 *allFuncsWhereScaled =
false;
754 EME::Derivative orig_DgDx_j;
756 !origFuncs.
supports(EME::OUT_ARG_DgDx,j).none()
757 && !(orig_DgDx_j=origFuncs.
get_DgDx(j)).isEmpty()
760 EME::Derivative scaled_DgDx_j;
761 bool didScaling =
false;
764 &scaled_DgDx_j, &didScaling
767 scaledFuncs->
set_DgDx(j,scaled_DgDx_j);
769 *allFuncsWhereScaled =
false;
772 for (
int l = 0; l < Np; ++l ) {
773 EME::Derivative orig_DgDp_j_l;
775 !origFuncs.
supports(EME::OUT_ARG_DgDp,j,l).none()
776 && !(orig_DgDp_j_l=origFuncs.
get_DgDp(j,l)).isEmpty()
779 EME::Derivative scaled_DgDp_j_l;
780 bool didScaling =
false;
783 &scaled_DgDp_j_l, &didScaling
786 scaledFuncs->
set_DgDp(j,l,scaled_DgDp_j_l);
788 *allFuncsWhereScaled =
false;
802 Teuchos::rcpWithEmbeddedObj(
806 invScalingVector->Reciprocal(*scalingVector);
807 return invScalingVector;
826 for (
int i = 0; i < localDim; ++i )
827 (*scaledVars)[i] = origVars[i] / invVarScaling[i];
851 scaledVars->Multiply( 1.0, invVarScaling, origVars, 0.0 );
861 funcs->Multiply( 1.0, fwdFuncScaling, *funcs, 0.0 );
879 if (funcDerivRowMatrix) {
881 funcDerivRowMatrix->
LeftScale(*fwdFuncScaling);
883 funcDerivRowMatrix->
RightScale(*invVarScaling);
904 const RCP<Epetra_MultiVector>
906 const EME::EDerivativeMultiVectorOrientation
909 if ( funcDerivMv_orientation == EME::DERIV_MV_BY_COL )
911 if (fwdFuncScaling) {
912 funcDerivMv->Multiply(1.0, *fwdFuncScaling, *funcDerivMv, 0.0);
919 else if ( funcDerivMv_orientation == EME::DERIV_TRANS_MV_BY_ROW )
922 funcDerivMv->Multiply(1.0, *invVarScaling, *funcDerivMv, 0.0);
924 if (fwdFuncScaling) {
932 *scaledFuncDeriv = EME::Derivative(funcDerivMv,funcDerivMv_orientation);
940 &*funcDerivOp, didScaling );
942 *scaledFuncDeriv = EME::Derivative(funcDerivOp);
Class that gets and sets p(l) in an InArgs object.
Teuchos::RCP< const Teuchos::Polynomial< Epetra_Vector > > get_x_dot_poly() const
Get time derivative vector Taylor polynomial.
bool supports(EOutArgsMembers arg) const
virtual InArgs createInArgs() const =0
virtual int RightScale(const Epetra_Vector &x)=0
Evaluation< Epetra_Vector > get_g(int j) const
Get g(j) where 0 <= j && j < this->Ng().
void set_DgDx(int j, const Derivative &DgDx_j)
bool is_null(const boost::shared_ptr< T > &p)
Derivative get_DfDp(int l) const
Teuchos::RCP< const Epetra_Vector > get_x_dotdot() const
void set_x_dotdot_poly(const Teuchos::RCP< const Teuchos::Polynomial< Epetra_Vector > > &x_dotdot_poly)
virtual Teuchos::RCP< const Epetra_Vector > get_p_init(int l) const
bool SameAs(const Epetra_BlockMap &Map) const
void set_beta(double beta)
Teuchos::RCP< Teuchos::Polynomial< Epetra_Vector > > get_f_poly() const
Get residual vector Taylor polynomial.
void set_x(const Teuchos::RCP< const Epetra_Vector > &x)
void set_x_dot(const Teuchos::RCP< const Epetra_Vector > &x_dot)
virtual Teuchos::RCP< const Epetra_Vector > get_p_upper_bounds(int l) const
virtual double get_t_upper_bound() const
Teuchos::RCP< Epetra_Operator > get_W() const
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void set_DfDp(int l, const Derivative &DfDp_l)
void set_DgDx_dotdot(int j, const Derivative &DgDx_dotdot_j)
Derivative get_DgDx_dot(int j) const
Class that gets and sets x_dotdot in an InArgs object.
Teuchos::RCP< const Teuchos::Polynomial< Epetra_Vector > > get_x_poly() const
Get solution vector Taylor polynomial.
Teuchos::RCP< Epetra_Operator > getLinearOp() const
Teuchos::RCP< const Epetra_Vector > get_p(int l) const
virtual int LeftScale(const Epetra_Vector &x)=0
Simple aggregate class that stores a derivative object as a general linear operator or as a multi-vec...
Class that gets and sets x_dot in an InArgs object.
Derivative get_DgDx(int j) const
virtual Teuchos::RCP< const Epetra_Vector > get_x_upper_bounds() const
void set_x_dotdot(const Teuchos::RCP< const Epetra_Vector > &x_dotdot)
void set_DgDx_dot(int j, const Derivative &DgDx_dot_j)
Teuchos::RCP< Epetra_MultiVector > getMultiVector() const
int NumMyElements() const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
virtual Teuchos::RCP< const Epetra_Vector > get_x_lower_bounds() const
Teuchos::RCP< const Teuchos::Polynomial< Epetra_Vector > > get_x_dotdot_poly() const
virtual Teuchos::RCP< const Epetra_Vector > get_x_dotdot_init() const
void set_p(int l, const Teuchos::RCP< const Epetra_Vector > &p_l)
Derivative get_DgDp(int j, int l) const
Class that gets and sets f in an OutArgs object.
virtual const Epetra_BlockMap & Map() const =0
Teuchos::RCP< const Epetra_Vector > get_x_dot() const
bool supports(EInArgsMembers arg) const
void set_W(const Teuchos::RCP< Epetra_Operator > &W)
virtual Teuchos::RCP< const Epetra_Vector > get_x_dot_init() const
Class that gets and sets g(j) in an OutArgs object.
Derivative get_DgDx_dotdot(int j) const
void set_DgDp(int j, int l, const Derivative &DgDp_j_l)
TEUCHOSCORE_LIB_DLL_EXPORT bool includesVerbLevel(const EVerbosityLevel verbLevel, const EVerbosityLevel requestedVerbLevel, const bool isDefaultLevel=false)
void set_f_poly(const Teuchos::RCP< Teuchos::Polynomial< Epetra_Vector > > &f_poly)
Set residual vector Taylor polynomial.
TypeTo implicit_cast(const TypeFrom &t)
EDerivativeMultiVectorOrientation getMultiVectorOrientation() const
void set_x_poly(const Teuchos::RCP< const Teuchos::Polynomial< Epetra_Vector > > &x_poly)
virtual Teuchos::RCP< const Epetra_Vector > get_p_lower_bounds(int l) const
virtual double get_t_lower_bound() const
void set_alpha(double alpha)
Class that gets and sets x in an InArgs object.
virtual Teuchos::RCP< const Epetra_Vector > get_x_init() const
virtual double get_t_init() const
void set_x_dot_poly(const Teuchos::RCP< const Teuchos::Polynomial< Epetra_Vector > > &x_dot_poly)
Set time derivative vector Taylor polynomial.
Evaluation< Epetra_Vector > get_f() const
Teuchos::RCP< const Epetra_Vector > get_x() const
Set solution vector Taylor polynomial.
Base interface for evaluating a stateless "model".
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)