ROL
Public Member Functions | List of all members
ROL::Vector< Real > Class Template Referenceabstract

Defines the linear algebra or vector space interface. More...

#include <ROL_Vector.hpp>

+ Inheritance diagram for ROL::Vector< Real >:

Public Member Functions

virtual ~Vector ()
 
virtual void plus (const Vector &x)=0
 Compute \(y \leftarrow y + x\), where \(y = \mathtt{*this}\). More...
 
virtual void scale (const Real alpha)=0
 Compute \(y \leftarrow \alpha y\) where \(y = \mathtt{*this}\). More...
 
virtual Real dot (const Vector &x) const =0
 Compute \( \langle y,x \rangle \) where \(y = \mathtt{*this}\). More...
 
virtual Real norm () const =0
 Returns \( \| y \| \) where \(y = \mathtt{*this}\). More...
 
virtual Teuchos::RCP< Vectorclone () const =0
 Clone to make a new (uninitialized) vector. More...
 
virtual void axpy (const Real alpha, const Vector &x)
 Compute \(y \leftarrow \alpha x + y\) where \(y = \mathtt{*this}\). More...
 
virtual void zero ()
 Set to zero vector. More...
 
virtual Teuchos::RCP< Vectorbasis (const int i) const
 Return i-th basis vector. More...
 
virtual int dimension () const
 Return dimension of the vector space. More...
 
virtual void set (const Vector &x)
 Set \(y \leftarrow x\) where \(y = \mathtt{*this}\). More...
 
virtual const Vectordual () const
 Return dual representation of \(\mathtt{*this}\), for example, the result of applying a Riesz map, or change of basis, or change of memory layout. More...
 
virtual std::vector< Real > checkVector (const Vector< Real > &x, const Vector< Real > &y, const bool printToStream=true, std::ostream &outStream=std::cout) const
 Verify vector-space methods. More...
 

Detailed Description

template<class Real>
class ROL::Vector< Real >

Defines the linear algebra or vector space interface.

The basic linear algebra interface, to be implemented by the user, includes:

The dot product can represent an inner product (in Hilbert space) or a duality pairing (in general Banach space).

There are additional virtual member functions that can be overloaded for computational efficiency.

Definition at line 72 of file ROL_Vector.hpp.

Constructor & Destructor Documentation

template<class Real>
virtual ROL::Vector< Real >::~Vector ( )
inlinevirtual

Definition at line 75 of file ROL_Vector.hpp.

Member Function Documentation

template<class Real>
virtual void ROL::Vector< Real >::plus ( const Vector< Real > &  x)
pure virtual

Compute \(y \leftarrow y + x\), where \(y = \mathtt{*this}\).

Parameters
[in]xis the vector to be added to \(\mathtt{*this}\).

On return \(\mathtt{*this} = \mathtt{*this} + x\).


Implemented in ConDualStdVector< Real, Element >, ConDualStdVector< Real, Element >, ConDualStdVector< Real, Element >, ConStdVector< Real, Element >, ConStdVector< Real, Element >, ConStdVector< Real, Element >, OptDualStdVector< Real, Element >, OptDualStdVector< Real, Element >, OptDualStdVector< Real, Element >, OptDualStdVector< Real, Element >, OptStdVector< Real, Element >, OptStdVector< Real, Element >, OptStdVector< Real, Element >, OptStdVector< Real, Element >, ROL::CArrayVector< Real, Element >, ROL::StdVector< Real, Element >, ROL::Vector_SimOpt< Real >, and ROL::CVaRVector< Real >.

Referenced by ROL::CompositeStepSQP< Real >::accept(), ROL::Bundle< Real >::aggregate(), ROL::ProjectedHessian< Real >::apply(), ROL::ProjectedPreconditioner< Real >::apply(), ROL::ProjectedHessian< Real >::applyInverse(), ROL::EqualityConstraint_SimOpt< Real >::applyJacobian(), ROL::Vector< Real >::axpy(), ROL::PrimalDualActiveSetStep< Real >::compute(), ROL::LineSearchStep< Real >::compute(), ROL::CompositeStepSQP< Real >::computeLagrangeMultiplier(), ROL::BoundConstraint< Real >::computeProjectedStep(), ROL::Constraints< Real >::computeProjectedStep(), ROL::MeanVariance< Real >::getGradient(), ROL::ExpUtility< Real >::getHessVec(), ROL::MeanVariance< Real >::getHessVec(), ROL::Reduced_ParametrizedObjective_SimOpt< Real >::gradient(), ROL::Reduced_Objective_SimOpt< Real >::gradient(), ROL::ZOO::Objective_PoissonInversion< Real >::gradient(), ROL::Reduced_ParametrizedObjective_SimOpt< Real >::hessVec(), ROL::Reduced_Objective_SimOpt< Real >::hessVec(), ROL::ProjectedObjective< Real >::reducedHessVec(), ROL::ProjectedObjective< Real >::reducedInvHessVec(), ROL::ProjectedObjective< Real >::reducedPrecond(), ROL::Vector< Real >::set(), ROL::EqualityConstraint< Real >::solveAugmentedSystem(), Normalization_Constraint< Real >::solveAugmentedSystem(), ROL::CompositeStepSQP< Real >::solveTangentialSubproblem(), ROL::CompositeStepSQP< Real >::update(), ROL::BundleStep< Real >::update(), ROL::PrimalDualActiveSetStep< Real >::update(), and BoundaryValueProblem< Real >::value().

template<class Real>
virtual void ROL::Vector< Real >::scale ( const Real  alpha)
pure virtual

Compute \(y \leftarrow \alpha y\) where \(y = \mathtt{*this}\).

Parameters
[in]alphais the scaling of \(\mathtt{*this}\).

On return \(\mathtt{*this} = \alpha (\mathtt{*this}) \).


Implemented in ConDualStdVector< Real, Element >, ConDualStdVector< Real, Element >, ConDualStdVector< Real, Element >, ConStdVector< Real, Element >, ConStdVector< Real, Element >, ConStdVector< Real, Element >, OptDualStdVector< Real, Element >, OptDualStdVector< Real, Element >, OptDualStdVector< Real, Element >, OptDualStdVector< Real, Element >, OptStdVector< Real, Element >, OptStdVector< Real, Element >, ROL::CArrayVector< Real, Element >, OptStdVector< Real, Element >, OptStdVector< Real, Element >, ROL::StdVector< Real, Element >, ROL::Vector_SimOpt< Real >, and ROL::CVaRVector< Real >.

Referenced by ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian_11(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian_12(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian_21(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian_22(), ROL::BarzilaiBorwein< Real >::applyB(), ROL::Secant< Real >::applyB0(), ROL::lDFP< Real >::applyB0(), ROL::BarzilaiBorwein< Real >::applyH(), ROL::lDFP< Real >::applyH0(), ROL::Secant< Real >::applyH0(), ROL::EqualityConstraint< Real >::applyJacobian(), ROL::EqualityConstraint_SimOpt< Real >::applyJacobian_1(), ROL::EqualityConstraint_SimOpt< Real >::applyJacobian_2(), ROL::CauchyPoint< Real >::cauchypoint_M(), ROL::CauchyPoint< Real >::cauchypoint_unc(), ROL::BundleStep< Real >::compute(), ROL::LineSearchStep< Real >::compute(), ROL::CompositeStepSQP< Real >::computeQuasinormalStep(), ROL::ExpUtility< Real >::getGradient(), ROL::MeanVariance< Real >::getGradient(), ROL::ExpUtility< Real >::getHessVec(), ROL::MeanVariance< Real >::getHessVec(), ROL::ZOO::Objective_Zakharov< Real >::gradient(), ROL::Objective< Real >::hessVec(), ROL::ZOO::Objective_DiodeCircuit< Real >::hessVec(), ROL::Objective_SimOpt< Real >::hessVec_11(), ROL::Objective_SimOpt< Real >::hessVec_12(), ROL::Objective_SimOpt< Real >::hessVec_21(), ROL::Objective_SimOpt< Real >::hessVec_22(), ROL::ZOO::Objective_Zakharov< Real >::invHessVec(), ROL::ZOO::Objective_PoissonInversion< Real >::reg_gradient(), ROL::ZOO::Objective_PoissonInversion< Real >::reg_hessVec(), ROL::DogLeg< Real >::run(), ROL::DoubleDogLeg< Real >::run(), Normalization_Constraint< Real >::solveAugmentedSystem(), and ROL::Vector< Real >::zero().

template<class Real>
virtual Real ROL::Vector< Real >::dot ( const Vector< Real > &  x) const
pure virtual

Compute \( \langle y,x \rangle \) where \(y = \mathtt{*this}\).

Parameters
[in]xis the vector that forms the dot product with \(\mathtt{*this}\).
Returns
The number equal to \(\langle \mathtt{*this}, x \rangle\).

Implemented in ConDualStdVector< Real, Element >, ConDualStdVector< Real, Element >, ConDualStdVector< Real, Element >, ConStdVector< Real, Element >, ConStdVector< Real, Element >, ConStdVector< Real, Element >, OptDualStdVector< Real, Element >, OptDualStdVector< Real, Element >, OptDualStdVector< Real, Element >, OptDualStdVector< Real, Element >, OptStdVector< Real, Element >, OptStdVector< Real, Element >, OptStdVector< Real, Element >, OptStdVector< Real, Element >, ROL::CArrayVector< Real, Element >, ROL::StdVector< Real, Element >, ROL::Vector_SimOpt< Real >, and ROL::CVaRVector< Real >.

Referenced by ROL::CompositeStepSQP< Real >::accept(), ROL::lBFGS< Real >::applyB(), ROL::lBFGS< Real >::applyH(), ROL::CauchyPoint< Real >::cauchypoint_CGT(), ROL::CauchyPoint< Real >::cauchypoint_M(), ROL::EqualityConstraint< Real >::checkAdjointConsistencyJacobian(), ROL::EqualityConstraint_SimOpt< Real >::checkAdjointConsistencyJacobian_1(), ROL::EqualityConstraint_SimOpt< Real >::checkAdjointConsistencyJacobian_2(), ROL::Objective< Real >::checkGradient(), ROL::Objective< Real >::checkHessSym(), ROL::BundleStep< Real >::compute(), ROL::LineSearchStep< Real >::compute(), ROL::Objective< Real >::gradient(), ROL::ZOO::Objective_Zakharov< Real >::gradient(), ROL::Objective_SimOpt< Real >::gradient_1(), ROL::Objective_SimOpt< Real >::gradient_2(), ROL::ZOO::Objective_Zakharov< Real >::invHessVec(), ROL::DogLeg< Real >::run(), ROL::DoubleDogLeg< Real >::run(), ROL::EqualityConstraint< Real >::solveAugmentedSystem(), ROL::CompositeStepSQP< Real >::solveTangentialSubproblem(), ROL::LineSearch< Real >::status(), ROL::lSR1< Real >::update(), ROL::Secant< Real >::update(), ROL::TrustRegion< Real >::update(), and ROL::ZOO::Objective_Zakharov< Real >::value().

template<class Real>
virtual Real ROL::Vector< Real >::norm ( ) const
pure virtual

Returns \( \| y \| \) where \(y = \mathtt{*this}\).

Returns
A nonnegative number equal to the norm of \(\mathtt{*this}\).

Implemented in ConDualStdVector< Real, Element >, ConDualStdVector< Real, Element >, ConDualStdVector< Real, Element >, ConStdVector< Real, Element >, ConStdVector< Real, Element >, ConStdVector< Real, Element >, OptDualStdVector< Real, Element >, OptDualStdVector< Real, Element >, OptDualStdVector< Real, Element >, OptDualStdVector< Real, Element >, OptStdVector< Real, Element >, OptStdVector< Real, Element >, OptStdVector< Real, Element >, OptStdVector< Real, Element >, ROL::CArrayVector< Real, Element >, ROL::StdVector< Real, Element >, ROL::Vector_SimOpt< Real >, and ROL::CVaRVector< Real >.

Referenced by ROL::CompositeStepSQP< Real >::accept(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian_11(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian_12(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian_21(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian_22(), ROL::EqualityConstraint< Real >::applyAdjointJacobian(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointJacobian_1(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointJacobian_2(), ROL::EqualityConstraint< Real >::applyJacobian(), ROL::EqualityConstraint_SimOpt< Real >::applyJacobian_1(), ROL::EqualityConstraint_SimOpt< Real >::applyJacobian_2(), ROL::CauchyPoint< Real >::cauchypoint_CGT(), ROL::CauchyPoint< Real >::cauchypoint_M(), ROL::EqualityConstraint_SimOpt< Real >::checkInverseAdjointJacobian_1(), ROL::EqualityConstraint_SimOpt< Real >::checkInverseJacobian_1(), ROL::PrimalDualActiveSetStep< Real >::compute(), ROL::LineSearchStep< Real >::compute(), ROL::TrustRegionStep< Real >::computeCriticalityMeasure(), ROL::Objective< Real >::hessVec(), ROL::ZOO::Objective_DiodeCircuit< Real >::hessVec(), ROL::Objective_SimOpt< Real >::hessVec_11(), ROL::Objective_SimOpt< Real >::hessVec_12(), ROL::Objective_SimOpt< Real >::hessVec_21(), ROL::Objective_SimOpt< Real >::hessVec_22(), ROL::ConjugateGradients< Real >::run(), ROL::ConjugateResiduals< Real >::run(), ROL::TruncatedCG< Real >::run(), EqualityConstraint_BurgersControl< Real >::solve(), ROL::CompositeStepSQP< Real >::solveTangentialSubproblem(), and ROL::PrimalDualActiveSetStep< Real >::update().

template<class Real>
virtual Teuchos::RCP<Vector> ROL::Vector< Real >::clone ( ) const
pure virtual

Clone to make a new (uninitialized) vector.

Returns
A reference-counted pointer to the cloned vector.

Provides the means of allocating temporary memory in ROL.


Implemented in ConDualStdVector< Real, Element >, ConDualStdVector< Real, Element >, ConDualStdVector< Real, Element >, ConStdVector< Real, Element >, ConStdVector< Real, Element >, ConStdVector< Real, Element >, OptDualStdVector< Real, Element >, OptDualStdVector< Real, Element >, OptDualStdVector< Real, Element >, OptDualStdVector< Real, Element >, OptStdVector< Real, Element >, OptStdVector< Real, Element >, OptStdVector< Real, Element >, ROL::CArrayVector< Real, Element >, OptStdVector< Real, Element >, ROL::StdVector< Real, Element >, ROL::Vector_SimOpt< Real >, and ROL::CVaRVector< Real >.

Referenced by ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian_11(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian_12(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian_21(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian_22(), ROL::EqualityConstraint< Real >::applyAdjointJacobian(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointJacobian_1(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointJacobian_2(), ROL::lBFGS< Real >::applyB(), ROL::lDFP< Real >::applyB(), ROL::lSR1< Real >::applyB(), ROL::lBFGS< Real >::applyH(), ROL::lDFP< Real >::applyH(), ROL::lSR1< Real >::applyH(), ROL::EqualityConstraint< Real >::applyJacobian(), ROL::EqualityConstraint_SimOpt< Real >::applyJacobian(), ROL::EqualityConstraint_SimOpt< Real >::applyJacobian_1(), ROL::EqualityConstraint_SimOpt< Real >::applyJacobian_2(), ROL::Vector< Real >::axpy(), ROL::EqualityConstraint< Real >::checkAdjointConsistencyJacobian(), ROL::EqualityConstraint_SimOpt< Real >::checkAdjointConsistencyJacobian_1(), ROL::EqualityConstraint_SimOpt< Real >::checkAdjointConsistencyJacobian_2(), ROL::EqualityConstraint< Real >::checkApplyAdjointHessian(), ROL::EqualityConstraint< Real >::checkApplyAdjointJacobian(), ROL::EqualityConstraint< Real >::checkApplyJacobian(), ROL::Objective< Real >::checkGradient(), ROL::Objective< Real >::checkHessSym(), ROL::Objective< Real >::checkHessVec(), ROL::EqualityConstraint_SimOpt< Real >::checkInverseAdjointJacobian_1(), ROL::EqualityConstraint_SimOpt< Real >::checkInverseJacobian_1(), ROL::EqualityConstraint_SimOpt< Real >::checkSolve(), ROL::Vector< Real >::checkVector(), ROL::computeDenseHessian(), ROL::computeDotMatrix(), ROL::SparseGridGenerator< Real >::computeError(), ROL::BoundConstraint< Real >::computeProjectedGradient(), ROL::Constraints< Real >::computeProjectedGradient(), ROL::Objective< Real >::dirDeriv(), ROL::MeanVariance< Real >::getGradient(), ROL::MeanDeviation< Real >::getGradient(), ROL::ExpUtility< Real >::getHessVec(), ROL::MeanVariance< Real >::getHessVec(), ROL::MeanDeviation< Real >::getHessVec(), ROL::RiskNeutralObjective< Real >::gradient(), ROL::Reduced_ParametrizedObjective_SimOpt< Real >::gradient(), ROL::Reduced_Objective_SimOpt< Real >::gradient(), ROL::Objective_SimOpt< Real >::gradient_1(), ROL::Objective_SimOpt< Real >::gradient_2(), ROL::Objective< Real >::hessVec(), ROL::RiskNeutralObjective< Real >::hessVec(), ROL::Reduced_ParametrizedObjective_SimOpt< Real >::hessVec(), ROL::Reduced_Objective_SimOpt< Real >::hessVec(), ROL::ZOO::Objective_DiodeCircuit< Real >::hessVec(), ROL::Objective_SimOpt< Real >::hessVec_11(), ROL::Objective_SimOpt< Real >::hessVec_12(), ROL::Objective_SimOpt< Real >::hessVec_21(), ROL::Objective_SimOpt< Real >::hessVec_22(), ROL::BackTracking< Real >::initialize(), ROL::IterationScaling< Real >::initialize(), ROL::CubicInterp< Real >::initialize(), ROL::Brents< Real >::initialize(), ROL::Bisection< Real >::initialize(), ROL::GoldenSection< Real >::initialize(), ROL::DogLeg< Real >::initialize(), ROL::DoubleDogLeg< Real >::initialize(), ROL::PathBasedTargetLevel< Real >::initialize(), ROL::TruncatedCG< Real >::initialize(), ROL::CauchyPoint< Real >::initialize(), ROL::Step< Real >::initialize(), ROL::LineSearch< Real >::initialize(), ROL::TrustRegion< Real >::initialize(), ROL::Bundle< Real >::initialize(), ROL::BundleStep< Real >::initialize(), ROL::CompositeStepSQP< Real >::initialize(), ROL::PrimalDualActiveSetStep< Real >::initialize(), ROL::TrustRegionStep< Real >::initialize(), ROL::LineSearchStep< Real >::initialize(), ROL::BoundConstraint< Real >::isFeasible(), ROL::BoundConstraint< Real >::pruneInactive(), ROL::Constraints< Real >::pruneInactive(), ROL::ProjectedObjective< Real >::reducedHessVec(), ROL::ProjectedObjective< Real >::reducedInvHessVec(), ROL::ProjectedObjective< Real >::reducedPrecond(), ROL::ExpUtility< Real >::reset(), ROL::RiskMeasure< Real >::reset(), ROL::MeanVariance< Real >::reset(), ROL::MeanDeviation< Real >::reset(), ROL::MeanVarianceFromTarget< Real >::reset(), ROL::MeanDeviationFromTarget< Real >::reset(), ROL::ConjugateGradients< Real >::run(), ROL::ConjugateResiduals< Real >::run(), ROL::DefaultAlgorithm< Real >::run(), ROL::EqualityConstraint< Real >::solveAugmentedSystem(), ROL::Secant< Real >::test(), ROL::lSR1< Real >::update(), ROL::Secant< Real >::update(), ROL::MeanVariance< Real >::update(), ROL::MeanDeviation< Real >::update(), and Objective_PoissonInversion< Real >::update().

template<class Real>
virtual void ROL::Vector< Real >::axpy ( const Real  alpha,
const Vector< Real > &  x 
)
inlinevirtual

Compute \(y \leftarrow \alpha x + y\) where \(y = \mathtt{*this}\).

Parameters
[in]alphais the scaling of x.
[in]xis a vector.

On return \(\mathtt{*this} = \mathtt{*this} + \alpha x \). Uses clone, set, scale and plus for the computation. Please overload if a more efficient implementation is needed.


Reimplemented in ROL::Vector_SimOpt< Real >, and ROL::CVaRVector< Real >.

Definition at line 141 of file ROL_Vector.hpp.

References ROL::Vector< Real >::clone(), and ROL::Vector< Real >::plus().

Referenced by ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian_11(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian_12(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian_21(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointHessian_22(), ROL::EqualityConstraint< Real >::applyAdjointJacobian(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointJacobian_1(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointJacobian_2(), ROL::lBFGS< Real >::applyB(), ROL::lDFP< Real >::applyB(), ROL::lSR1< Real >::applyB(), ROL::lBFGS< Real >::applyH(), ROL::lDFP< Real >::applyH(), ROL::lSR1< Real >::applyH(), ROL::EqualityConstraint< Real >::applyJacobian(), ROL::EqualityConstraint_SimOpt< Real >::applyJacobian_1(), ROL::EqualityConstraint_SimOpt< Real >::applyJacobian_2(), ROL::CauchyPoint< Real >::cauchypoint_CGT(), ROL::TrustRegionStep< Real >::compute(), ROL::LineSearchStep< Real >::compute(), ROL::BoundConstraint< Real >::computeProjectedStep(), ROL::Constraints< Real >::computeProjectedStep(), ROL::CompositeStepSQP< Real >::computeQuasinormalStep(), ROL::MeanDeviationFromTarget< Real >::getGradient(), ROL::MeanVariance< Real >::getHessVec(), ROL::MeanDeviationFromTarget< Real >::getHessVec(), ROL::Objective< Real >::gradient(), ROL::ZOO::Objective_Zakharov< Real >::gradient(), ROL::Objective_SimOpt< Real >::gradient_1(), ROL::Objective_SimOpt< Real >::gradient_2(), ROL::Objective< Real >::hessVec(), ROL::ZOO::Objective_DiodeCircuit< Real >::hessVec(), ROL::Objective_SimOpt< Real >::hessVec_11(), ROL::Objective_SimOpt< Real >::hessVec_12(), ROL::Objective_SimOpt< Real >::hessVec_21(), ROL::Objective_SimOpt< Real >::hessVec_22(), ROL::ZOO::Objective_Zakharov< Real >::invHessVec(), main(), ROL::BoundConstraint< Real >::pruneInactive(), ROL::Constraints< Real >::pruneInactive(), ROL::ConjugateGradients< Real >::run(), ROL::ConjugateResiduals< Real >::run(), ROL::DogLeg< Real >::run(), ROL::DoubleDogLeg< Real >::run(), ROL::TruncatedCG< Real >::run(), ROL::ZOO::Objective_PoissonInversion< Real >::solve_adjoint_sensitivity_equation(), ROL::TrustRegion< Real >::update(), ROL::LineSearchStep< Real >::update(), and ROL::LineSearch< Real >::updateIterate().

template<class Real>
virtual void ROL::Vector< Real >::zero ( )
inlinevirtual

Set to zero vector.

Uses scale by zero for the computation. Please overload if a more efficient implementation is needed.


Definition at line 155 of file ROL_Vector.hpp.

References ROL::Vector< Real >::scale().

Referenced by ROL::Bundle< Real >::aggregate(), EqualityConstraint_BurgersControl< Real >::applyAdjointHessian_12(), EqualityConstraint_BurgersControl< Real >::applyAdjointHessian_21(), EqualityConstraint_BurgersControl< Real >::applyAdjointHessian_22(), ROL::EqualityConstraint< Real >::applyAdjointJacobian(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointJacobian_1(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointJacobian_2(), ROL::EqualityConstraint< Real >::applyJacobian(), ROL::Vector< Real >::checkVector(), ROL::BundleStep< Real >::compute(), ROL::PrimalDualActiveSetStep< Real >::compute(), ROL::MeanVariance< Real >::getGradient(), ROL::MeanDeviation< Real >::getGradient(), ROL::MeanVariance< Real >::getHessVec(), ROL::MeanDeviation< Real >::getHessVec(), ROL::ZOO::Objective_BVP< Real >::gradient(), ROL::ZOO::Objective_HS25< Real >::gradient(), ROL::Objective< Real >::gradient(), ROL::RiskNeutralObjective< Real >::gradient(), ROL::RiskAverseObjective< Real >::gradient(), ROL::Objective_SimOpt< Real >::gradient_1(), ROL::Objective_SimOpt< Real >::gradient_2(), ROL::Objective< Real >::hessVec(), ROL::RiskNeutralObjective< Real >::hessVec(), ROL::RiskAverseObjective< Real >::hessVec(), ROL::ZOO::Objective_DiodeCircuit< Real >::hessVec(), ROL::Objective_SimOpt< Real >::hessVec_11(), ROL::Objective_SimOpt< Real >::hessVec_12(), Objective_BurgersControl< Real >::hessVec_12(), ROL::Objective_SimOpt< Real >::hessVec_21(), Objective_BurgersControl< Real >::hessVec_21(), ROL::Objective_SimOpt< Real >::hessVec_22(), ROL::ZOO::Objective_HS25< Real >::invHessVec(), main(), ROL::ConjugateGradients< Real >::run(), ROL::ConjugateResiduals< Real >::run(), ROL::TruncatedCG< Real >::run(), ROL::Vector< Real >::set(), ROL::ParametrizedEqualityConstraint_SimOpt< Real >::solve(), ROL::EqualityConstraint< Real >::solveAugmentedSystem(), ROL::CompositeStepSQP< Real >::solveTangentialSubproblem(), and BoundaryValueProblem< Real >::value().

template<class Real>
virtual Teuchos::RCP<Vector> ROL::Vector< Real >::basis ( const int  i) const
inlinevirtual
template<class Real>
virtual int ROL::Vector< Real >::dimension ( ) const
inlinevirtual
template<class Real>
virtual void ROL::Vector< Real >::set ( const Vector< Real > &  x)
inlinevirtual

Set \(y \leftarrow x\) where \(y = \mathtt{*this}\).

Parameters
[in]xis a vector.

On return \(\mathtt{*this} = x\). Uses zero and plus methods for the computation. Please overload if a more efficient implementation is needed.


Reimplemented in ROL::StdVector< Real, Element >.

Definition at line 194 of file ROL_Vector.hpp.

References ROL::Vector< Real >::plus(), and ROL::Vector< Real >::zero().

Referenced by ROL::CompositeStepSQP< Real >::accept(), ROL::Bundle< Real >::aggregate(), ROL::BarzilaiBorwein< Real >::applyB(), ROL::lDFP< Real >::applyB(), ROL::Secant< Real >::applyB0(), ROL::lDFP< Real >::applyB0(), ROL::lSR1< Real >::applyB0(), ROL::lBFGS< Real >::applyH(), ROL::BarzilaiBorwein< Real >::applyH(), ROL::lDFP< Real >::applyH0(), ROL::lSR1< Real >::applyH0(), ROL::Secant< Real >::applyH0(), ROL::LinearOperator< Real >::applyInverse(), ROL::EqualityConstraint< Real >::applyPreconditioner(), ROL::EqualityConstraint_SimOpt< Real >::applyPreconditioner(), ROL::CauchyPoint< Real >::cauchypoint_CGT(), ROL::CauchyPoint< Real >::cauchypoint_M(), ROL::CauchyPoint< Real >::cauchypoint_unc(), ROL::BundleStep< Real >::compute(), ROL::TrustRegionStep< Real >::compute(), ROL::LineSearchStep< Real >::compute(), ROL::CompositeStepSQP< Real >::computeQuasinormalStep(), ROL::ZOO::Objective_Zakharov< Real >::gradient(), ROL::RiskNeutralObjective< Real >::gradient(), ROL::ZOO::Objective_Zakharov< Real >::invHessVec(), main(), ROL::Objective< Real >::precond(), ROL::RiskNeutralObjective< Real >::precond(), ROL::RiskAverseObjective< Real >::precond(), ROL::Reduced_ParametrizedObjective_SimOpt< Real >::precond(), ROL::Reduced_Objective_SimOpt< Real >::precond(), ROL::ZOO::Objective_PoissonInversion< Real >::reg_gradient(), ROL::ZOO::Objective_PoissonInversion< Real >::reg_hessVec(), ROL::DogLeg< Real >::run(), ROL::DoubleDogLeg< Real >::run(), ROL::TruncatedCG< Real >::run(), BoundConstraint_PoissonControl< Real >::setVectorToLowerBound(), ROL::StdBoundConstraint< Real >::setVectorToLowerBound(), BoundConstraint_PoissonInversion< Real >::setVectorToLowerBound(), BoundConstraint_BurgersControl< Real >::setVectorToLowerBound(), ROL::ZOO::BoundConstraint_DiodeCircuit< Real >::setVectorToLowerBound(), BoundConstraint_PoissonControl< Real >::setVectorToUpperBound(), ROL::StdBoundConstraint< Real >::setVectorToUpperBound(), BoundConstraint_PoissonInversion< Real >::setVectorToUpperBound(), BoundConstraint_BurgersControl< Real >::setVectorToUpperBound(), ROL::ZOO::BoundConstraint_DiodeCircuit< Real >::setVectorToUpperBound(), ROL::ZOO::Objective_PoissonInversion< Real >::solve_poisson(), Normalization_Constraint< Real >::solveAugmentedSystem(), ROL::CompositeStepSQP< Real >::solveTangentialSubproblem(), ROL::BatchManager< Real >::sumAll(), ROL::TrustRegionStep< Real >::update(), and ROL::LineSearch< Real >::updateIterate().

template<class Real>
virtual const Vector& ROL::Vector< Real >::dual ( void  ) const
inlinevirtual

Return dual representation of \(\mathtt{*this}\), for example, the result of applying a Riesz map, or change of basis, or change of memory layout.

Returns
A const reference to dual representation.

By default, returns the current object. Please overload if you need a dual representation.


Reimplemented in ConDualStdVector< Real, Element >, ConDualStdVector< Real, Element >, ConDualStdVector< Real, Element >, ConStdVector< Real, Element >, ConStdVector< Real, Element >, ConStdVector< Real, Element >, OptDualStdVector< Real, Element >, OptDualStdVector< Real, Element >, OptDualStdVector< Real, Element >, OptDualStdVector< Real, Element >, OptStdVector< Real, Element >, OptStdVector< Real, Element >, OptStdVector< Real, Element >, OptStdVector< Real, Element >, and ROL::Vector_SimOpt< Real >.

Definition at line 211 of file ROL_Vector.hpp.

Referenced by ROL::CompositeStepSQP< Real >::accept(), ROL::EqualityConstraint< Real >::applyAdjointJacobian(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointJacobian_1(), ROL::EqualityConstraint_SimOpt< Real >::applyAdjointJacobian_2(), ROL::BarzilaiBorwein< Real >::applyB(), ROL::lDFP< Real >::applyB(), ROL::lSR1< Real >::applyB(), ROL::Secant< Real >::applyB0(), ROL::lDFP< Real >::applyB0(), ROL::lSR1< Real >::applyB0(), ROL::lBFGS< Real >::applyH(), ROL::lDFP< Real >::applyH(), ROL::BarzilaiBorwein< Real >::applyH(), ROL::lSR1< Real >::applyH(), ROL::lDFP< Real >::applyH0(), ROL::lSR1< Real >::applyH0(), ROL::Secant< Real >::applyH0(), ROL::EqualityConstraint< Real >::applyPreconditioner(), ROL::CauchyPoint< Real >::cauchypoint_CGT(), ROL::CauchyPoint< Real >::cauchypoint_M(), ROL::CauchyPoint< Real >::cauchypoint_unc(), ROL::EqualityConstraint< Real >::checkAdjointConsistencyJacobian(), ROL::EqualityConstraint_SimOpt< Real >::checkAdjointConsistencyJacobian_1(), ROL::EqualityConstraint_SimOpt< Real >::checkAdjointConsistencyJacobian_2(), ROL::EqualityConstraint< Real >::checkApplyAdjointJacobian(), ROL::Objective< Real >::checkGradient(), ROL::TrustRegionStep< Real >::computeCriticalityMeasure(), ROL::CompositeStepSQP< Real >::computeQuasinormalStep(), ROL::Objective< Real >::precond(), ROL::RiskNeutralObjective< Real >::precond(), ROL::RiskAverseObjective< Real >::precond(), ROL::Reduced_Objective_SimOpt< Real >::precond(), ROL::DogLeg< Real >::run(), ROL::DoubleDogLeg< Real >::run(), ROL::DefaultAlgorithm< Real >::run(), ROL::EqualityConstraint< Real >::solveAugmentedSystem(), ROL::CompositeStepSQP< Real >::solveTangentialSubproblem(), ROL::RiskNeutralObjective< Real >::update(), and ROL::TrustRegion< Real >::update().

template<class Real>
virtual std::vector<Real> ROL::Vector< Real >::checkVector ( const Vector< Real > &  x,
const Vector< Real > &  y,
const bool  printToStream = true,
std::ostream &  outStream = std::cout 
) const
inlinevirtual

Verify vector-space methods.

Parameters
[in]xis a vector.
[in]yis a vector.

Returns a vector of Reals, all of which should be close to zero. They represent consistency errors in the vector space properties, as follows:

  • Commutativity of addition: \(\mathtt{*this} + x = x + \mathtt{*this}\).
  • Associativity of addition: \(\mathtt{*this} + (x + y) = (\mathtt{*this} + x) + y\).
  • Identity element of addition: \(0 + x = x\).
  • Inverse elements of addition: \(\mathtt{*this} + (- \mathtt{*this}) = 0\).
  • Identity element of scalar multiplication: \( 1 \cdot \mathtt{*this} = \mathtt{*this} \).
  • Consistency of scalar multiplication with field multiplication: \(a (b \cdot \mathtt{*this}) = (a b) \cdot \mathtt{*this}\).
  • Distributivity of scalar multiplication with respect to field addition: \((a+b) \cdot \mathtt{*this} = a \cdot \mathtt{*this} + b \cdot \mathtt{*this}\).
  • Distributivity of scalar multiplication with respect to vector addition: \(a \cdot (\mathtt{*this} + x) = a \cdot \mathtt{*this} + a \cdot x\).
  • Commutativity of dot (inner) product over the field of reals: \((\mathtt{*this}, x) = (x, \mathtt{*this})\).
  • Additivity of dot (inner) product: \((\mathtt{*this}, x+y) = (\mathtt{*this}, x) + (\mathtt{*this}, y)\).
  • Consistency of scalar multiplication and norm: \(\|{\mathtt{*this}} / {\|\mathtt{*this}\|} \| = 1\).
  • Reflexivity: \(\mbox{dual}(\mbox{dual}(\mathtt{*this})) = \mathtt{*this}\) .

The consistency errors are defined as the norms or absolute values of the differences between the left-hand side and the right-hand side terms in the above equalities.


Definition at line 242 of file ROL_Vector.hpp.

References ROL::Vector< Real >::clone(), and ROL::Vector< Real >::zero().

Referenced by main().


The documentation for this class was generated from the following file: