NOX  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Public Types | Public Member Functions | List of all members
NOX::Abstract::Vector Class Referenceabstract

NOX's pure abstract vector interface for vectors that are used by the nonlinear solver. More...

#include <NOX_Abstract_Vector.H>

Inheritance diagram for NOX::Abstract::Vector:
Inheritance graph
[legend]

Public Types

enum  NormType { TwoNorm, OneNorm, MaxNorm }
 Norm types used in norm() calculations. More...
 

Public Member Functions

 Vector ()
 Abstract Vector constructor (does nothing)
 
virtual ~Vector ()
 Abstract Vector destructor (does nothing)
 
virtual NOX::size_type length () const =0
 Return the length of vector. More...
 
virtual void print (std::ostream &stream) const
 Print the vector. To be used for debugging only.
 
virtual NOX::Abstract::Vectorinit (double gamma)=0
 Initialize every element of this vector with gamma. More...
 
virtual NOX::Abstract::Vectorrandom (bool useSeed=false, int seed=1)=0
 Initialize each element of this vector with a random value. More...
 
virtual NOX::Abstract::Vectorabs (const NOX::Abstract::Vector &y)=0
 Put element-wise absolute values of source vector y into this vector. More...
 
virtual NOX::Abstract::Vectoroperator= (const NOX::Abstract::Vector &y)=0
 Copy source vector y into this vector. More...
 
virtual NOX::Abstract::Vectorreciprocal (const NOX::Abstract::Vector &y)=0
 Put element-wise reciprocal of source vector y into this vector. More...
 
virtual NOX::Abstract::Vectorscale (double gamma)=0
 Scale each element of this vector by gamma. More...
 
virtual NOX::Abstract::Vectorscale (const NOX::Abstract::Vector &a)=0
 Scale this vector element-by-element by the vector a. More...
 
virtual NOX::Abstract::Vectorupdate (double alpha, const NOX::Abstract::Vector &a, double gamma=0.0)=0
 Compute x = (alpha * a) + (gamma * x) where x is this vector. More...
 
virtual NOX::Abstract::Vectorupdate (double alpha, const NOX::Abstract::Vector &a, double beta, const NOX::Abstract::Vector &b, double gamma=0.0)=0
 Compute x = (alpha * a) + (beta * b) + (gamma * x) where x is this vector. More...
 
virtual Teuchos::RCP
< NOX::Abstract::Vector
clone (NOX::CopyType type=NOX::DeepCopy) const =0
 Create a new Vector of the same underlying type by cloning "this", and return a pointer to the new vector. More...
 
virtual Teuchos::RCP
< NOX::Abstract::MultiVector
createMultiVector (const NOX::Abstract::Vector *const *vecs, int numVecs, NOX::CopyType type=NOX::DeepCopy) const
 Create a MultiVector with numVecs+1 columns out of an array of Vectors. The vector stored under this will be the first column with the remaining numVecs columns given by vecs. More...
 
virtual Teuchos::RCP
< NOX::Abstract::MultiVector
createMultiVector (int numVecs, NOX::CopyType type=NOX::DeepCopy) const
 Create a MultiVector with numVecs columns. More...
 
virtual double norm (NOX::Abstract::Vector::NormType type=NOX::Abstract::Vector::TwoNorm) const =0
 Norm. More...
 
virtual double norm (const NOX::Abstract::Vector &weights) const =0
 Weighted 2-Norm. More...
 
virtual double innerProduct (const NOX::Abstract::Vector &y) const =0
 Inner product with y. More...
 

Detailed Description

NOX's pure abstract vector interface for vectors that are used by the nonlinear solver.

This class is a member of the namespace NOX::Abstract.

The user should implement their own concrete implementation of this class or use one of the implementations provided by us.

Author
Tammy Kolda (SNL 8950), Roger Pawlowski (SNL 9233)

Member Enumeration Documentation

Norm types used in norm() calculations.

Enumerator
TwoNorm 

Use the 2-norm.

OneNorm 

Use the 1-norm.

MaxNorm 

Use the max-norm, a.k.a. the $\infty$-norm.

Member Function Documentation

virtual NOX::Abstract::Vector& NOX::Abstract::Vector::abs ( const NOX::Abstract::Vector y)
pure virtual
virtual Teuchos::RCP<NOX::Abstract::Vector> NOX::Abstract::Vector::clone ( NOX::CopyType  type = NOX::DeepCopy) const
pure virtual

Create a new Vector of the same underlying type by cloning "this", and return a pointer to the new vector.

If type is NOX::DeepCopy, then we need to create an exact replica of "this". Otherwise, if type is NOX::ShapeCopy, we need only replicate the shape of "this" (the memory is allocated for the objects, but the current values are not copied into the vector). Note that there is no assumption that a vector created by ShapeCopy is initialized to zeros.

Returns
Pointer to newly created vector or NULL if clone is not supported.

Implemented in NOX::Epetra::Vector, NOX::Thyra::Vector, NOX::Petsc::Vector, NOX::LAPACK::Vector, LOCA::Extended::Vector, LOCA::Hopf::MooreSpence::ExtendedVector, LOCA::Pitchfork::MooreSpence::ExtendedVector, LOCA::TurningPoint::MooreSpence::ExtendedVector, LOCA::PhaseTransition::ExtendedVector, LOCA::Hopf::ComplexVector, and LOCA::MultiContinuation::ExtendedVector.

Referenced by LOCA::Bifurcation::TPBord::StatusTest::NullVectorNormWRMS::checkStatus(), LOCA::Bifurcation::PitchforkBord::StatusTest::NullVectorNormWRMS::checkStatus(), NOX::StatusTest::NormUpdate::checkStatus(), NOX::StatusTest::NormWRMS::checkStatus(), LOCA::SingularJacobianSolve::ItRef::compute(), LOCA::SingularJacobianSolve::Nic::compute(), LOCA::SingularJacobianSolve::NicDay::compute(), NOX::Direction::NonlinearCG::compute(), LOCA::DerivUtils::computeDCeDxa(), LOCA::DerivUtils::computeDJnDxa(), LOCA::DerivUtils::computeDwtCeDp(), LOCA::DerivUtils::computeDwtCeDx(), LOCA::DerivUtils::computeDwtJnDp(), LOCA::DerivUtils::computeDwtJnDx(), NOX::Direction::Utils::InexactNewton::computeForcingTerm(), LOCA::SingularJacobianSolve::ItRef::computeMulti(), LOCA::SingularJacobianSolve::Nic::computeMulti(), LOCA::SingularJacobianSolve::NicDay::computeMulti(), NOX::MeritFunction::SumOfSquares::computeQuadraticMinimizer(), NOX::MeritFunction::SumOfSquares::computeQuadraticModel(), LOCA::Epetra::Group::computeScaledDotProduct(), NOX::LineSearch::Utils::Slope::computeSlope(), NOX::MeritFunction::SumOfSquares::computeSlope(), NOX::LineSearch::Utils::Slope::computeSlopeWithOutJac(), LOCA::MultiContinuation::ConstraintModelEvaluator::ConstraintModelEvaluator(), NOX::Solver::TensorBased::getDirectionalDerivative(), LOCA::TurningPoint::MooreSpence::ExtendedGroup::getLengthVector(), NOX::Solver::TensorBased::getNormModelResidual(), LOCA::Homotopy::Group::Group(), LOCA::Epetra::Group::Group(), LOCA::Hopf::MooreSpence::ExtendedGroup::init(), LOCA::Pitchfork::MooreSpence::ExtendedGroup::init(), LOCA::TurningPoint::MooreSpence::ExtendedGroup::init(), NOX::MultiVector::MultiVector(), LOCA::Epetra::Group::operator=(), LOCA::DerivUtils::perturbXVec(), LOCA::AnasaziOperator::JacobianInverse::rayleighQuotient(), NOX::Direction::Broyden::BroydenMemoryUnit::reset(), NOX::Direction::Newton::resetForcingTerm(), LOCA::Epetra::Group::setScaleVector(), and LOCA::Extended::Vector::setVector().

Teuchos::RCP< NOX::Abstract::MultiVector > NOX::Abstract::Vector::createMultiVector ( const NOX::Abstract::Vector *const *  vecs,
int  numVecs,
NOX::CopyType  type = NOX::DeepCopy 
) const
virtual

Create a MultiVector with numVecs+1 columns out of an array of Vectors. The vector stored under this will be the first column with the remaining numVecs columns given by vecs.

The default implementation creates a generic NOX::MultiVector with either Shape or Deep copies of the supplied vectors.

Reimplemented in NOX::Epetra::Vector, NOX::Thyra::Vector, and LOCA::Extended::Vector.

References Teuchos::rcp().

Referenced by LOCA::MultiContinuation::ConstrainedGroup::applyJacobian(), LOCA::Homotopy::DeflatedGroup::applyJacobian(), LOCA::Hopf::MooreSpence::ExtendedGroup::applyJacobian(), LOCA::Pitchfork::MooreSpence::ExtendedGroup::applyJacobian(), LOCA::TurningPoint::MooreSpence::ExtendedGroup::applyJacobian(), LOCA::Pitchfork::MinimallyAugmented::ExtendedGroup::applyJacobian(), LOCA::Hopf::MinimallyAugmented::ExtendedGroup::applyJacobian(), LOCA::MultiContinuation::ConstrainedGroup::applyJacobianInverse(), LOCA::Homotopy::DeflatedGroup::applyJacobianInverse(), LOCA::Hopf::MooreSpence::ExtendedGroup::applyJacobianInverse(), LOCA::Pitchfork::MooreSpence::ExtendedGroup::applyJacobianInverse(), LOCA::Pitchfork::MinimallyAugmented::ExtendedGroup::applyJacobianInverse(), LOCA::TurningPoint::MooreSpence::ExtendedGroup::applyJacobianInverse(), LOCA::Hopf::MinimallyAugmented::ExtendedGroup::applyJacobianInverse(), LOCA::MultiContinuation::ConstrainedGroup::applyJacobianTranspose(), LOCA::Homotopy::DeflatedGroup::applyJacobianTranspose(), LOCA::Hopf::MooreSpence::ExtendedGroup::applyJacobianTranspose(), LOCA::Pitchfork::MooreSpence::ExtendedGroup::applyJacobianTranspose(), LOCA::TurningPoint::MooreSpence::ExtendedGroup::applyJacobianTranspose(), LOCA::Pitchfork::MinimallyAugmented::ExtendedGroup::applyJacobianTranspose(), LOCA::Hopf::MinimallyAugmented::ExtendedGroup::applyJacobianTranspose(), LOCA::MultiContinuation::ConstrainedGroup::applyJacobianTransposeInverse(), LOCA::TurningPoint::MooreSpence::ExtendedGroup::applyJacobianTransposeInverse(), LOCA::Hopf::ComplexMultiVector::ComplexMultiVector(), LOCA::MultiPredictor::Tangent::compute(), LOCA::Eigensolver::DGGEVStrategy::computeEigenvalues(), LOCA::Eigensolver::AnasaziStrategy::computeEigenvalues(), LOCA::Hopf::MinimallyAugmented::Constraint::Constraint(), LOCA::MultiContinuation::ConstraintModelEvaluator::ConstraintModelEvaluator(), LOCA::Pitchfork::MooreSpence::ExtendedGroup::ExtendedGroup(), LOCA::Hopf::MooreSpence::ExtendedGroup::ExtendedGroup(), LOCA::TurningPoint::MooreSpence::ExtendedGroup::ExtendedGroup(), LOCA::PhaseTransition::ExtendedMultiVector::ExtendedMultiVector(), LOCA::MultiContinuation::ExtendedMultiVector::ExtendedMultiVector(), LOCA::Pitchfork::MooreSpence::ExtendedMultiVector::ExtendedMultiVector(), LOCA::TurningPoint::MooreSpence::ExtendedMultiVector::ExtendedMultiVector(), LOCA::Hopf::MooreSpence::ExtendedMultiVector::ExtendedMultiVector(), LOCA::AnasaziOperator::ShiftInvert::rayleighQuotient(), LOCA::Epetra::AnasaziOperator::Floquet::rayleighQuotient(), LOCA::AnasaziOperator::ShiftInvert2Matrix::rayleighQuotient(), LOCA::AnasaziOperator::Cayley::rayleighQuotient(), LOCA::AnasaziOperator::Cayley2Matrix::rayleighQuotient(), and LOCA::Pitchfork::MooreSpence::PhippsBordering::setBlocks().

Teuchos::RCP< NOX::Abstract::MultiVector > NOX::Abstract::Vector::createMultiVector ( int  numVecs,
NOX::CopyType  type = NOX::DeepCopy 
) const
virtual

Create a MultiVector with numVecs columns.

The default implementation creates a generic NOX::MultiVector with either Shape or Deep copies of the supplied vector.

Reimplemented in NOX::Epetra::Vector, NOX::Thyra::Vector, and LOCA::Extended::Vector.

References Teuchos::rcp().

virtual NOX::Abstract::Vector& NOX::Abstract::Vector::init ( double  gamma)
pure virtual
virtual double NOX::Abstract::Vector::innerProduct ( const NOX::Abstract::Vector y) const
pure virtual

Inner product with y.

Here x represents this vector, and we compute its inner product with y as follows:

\[ \langle x,y \rangle = \sum_{i=1}^n x_i y_i \]

Returns
$\langle x,y \rangle$

Implemented in NOX::Epetra::Vector, LOCA::Extended::Vector, NOX::Thyra::Vector, NOX::Petsc::Vector, and NOX::LAPACK::Vector.

Referenced by LOCA::SingularJacobianSolve::Nic::compute(), LOCA::SingularJacobianSolve::NicDay::compute(), NOX::LineSearch::NonlinearCG::compute(), NOX::Direction::NonlinearCG::compute(), NOX::Direction::Broyden::compute(), NOX::Solver::TensorBased::computeCurvilinearStep(), LOCA::DerivUtils::computeDwtCeDp(), LOCA::DerivUtils::computeDwtJnDp(), LOCA::SingularJacobianSolve::Nic::computeMulti(), LOCA::SingularJacobianSolve::NicDay::computeMulti(), NOX::MeritFunction::SumOfSquares::computeQuadraticMinimizer(), LOCA::MultiContinuation::AbstractGroup::computeScaledDotProduct(), LOCA::Thyra::Group::computeScaledDotProduct(), LOCA::LAPACK::Group::computeScaledDotProduct(), LOCA::Epetra::Group::computeScaledDotProduct(), NOX::LineSearch::Utils::Slope::computeSlope(), NOX::MeritFunction::SumOfSquares::computeSlope(), NOX::Solver::TensorBased::getDirectionalDerivative(), NOX::Solver::TensorBased::getNormModelResidual(), LOCA::Pitchfork::MooreSpence::AbstractGroup::innerProduct(), NOX::Solver::InexactTrustRegionBased::iterateInexact(), NOX::Solver::InexactTrustRegionBased::iterateStandard(), NOX::Solver::AndersonAcceleration::qrAdd(), LOCA::AnasaziOperator::JacobianInverse::rayleighQuotient(), LOCA::AnasaziOperator::ShiftInvert::rayleighQuotient(), LOCA::Epetra::AnasaziOperator::Floquet::rayleighQuotient(), LOCA::AnasaziOperator::ShiftInvert2Matrix::rayleighQuotient(), LOCA::AnasaziOperator::Cayley::rayleighQuotient(), LOCA::AnasaziOperator::Cayley2Matrix::rayleighQuotient(), and NOX::Solver::TrustRegionBased::step().

virtual NOX::size_type NOX::Abstract::Vector::length ( ) const
pure virtual
virtual double NOX::Abstract::Vector::norm ( NOX::Abstract::Vector::NormType  type = NOX::Abstract::Vector::TwoNorm) const
pure virtual

Norm.

Here x represents this vector, and we compute its norm as follows: for each NOX::Abstract::Vector::NormType:

Returns
$\|x\|$

Implemented in NOX::Epetra::Vector, LOCA::Extended::Vector, NOX::Thyra::Vector, NOX::Petsc::Vector, and NOX::LAPACK::Vector.

Referenced by NOX::StatusTest::RelativeNormF::checkStatus(), NOX::StatusTest::FiniteValue::checkStatus(), LOCA::Bifurcation::TPBord::StatusTest::NullVectorNormWRMS::checkStatus(), LOCA::Bifurcation::PitchforkBord::StatusTest::NullVectorNormWRMS::checkStatus(), NOX::StatusTest::NormUpdate::checkStatus(), NOX::LineSearch::SafeguardedDirection::compute(), NOX::LineSearch::SafeguardedStep::compute(), NOX::Direction::SteepestDescent::compute(), NOX::Solver::InexactTrustRegionBased::computeNorm(), NOX::Thyra::WeightedMeritFunction::computeSlope(), NOX::LineSearch::Utils::Slope::computeSlopeWithOutJac(), NOX::Solver::TensorBased::computeTensorDirection(), LOCA::DerivUtils::epsVector(), NOX::Solver::TensorBased::getNormModelResidual(), NOX::Solver::TensorBased::printDirectionInfo(), NOX::Solver::TrustRegionBased::printUpdate(), LOCA::MultiContinuation::AbstractGroup::projectToDraw(), NOX::Solver::AndersonAcceleration::qrAdd(), NOX::Direction::Broyden::BroydenMemoryUnit::reset(), LOCA::TurningPoint::MooreSpence::ExtendedGroup::scaleNullVector(), LOCA::TurningPoint::MinimallyAugmented::Constraint::scaleNullVectors(), and NOX::Solver::TrustRegionBased::step().

virtual double NOX::Abstract::Vector::norm ( const NOX::Abstract::Vector weights) const
pure virtual

Weighted 2-Norm.

Here x represents this vector, and we compute its weighted norm as follows:

\[ \|x\|_w = \sqrt{\sum_{i=1}^{n} w_i \; x_i^2} \]

Returns
$ \|x\|_w $

Implemented in NOX::Epetra::Vector, LOCA::Extended::Vector, NOX::Thyra::Vector, NOX::Petsc::Vector, and NOX::LAPACK::Vector.

virtual NOX::Abstract::Vector& NOX::Abstract::Vector::operator= ( const NOX::Abstract::Vector y)
pure virtual
virtual NOX::Abstract::Vector& NOX::Abstract::Vector::random ( bool  useSeed = false,
int  seed = 1 
)
pure virtual

Initialize each element of this vector with a random value.

If useSeed is true, uses the value of seed to seed the random number generator before filling the entries of this vector. So, if two calls are made where useSeed is true and seed is the same, then the vectors returned should be the same.

Default implementation throw an error. Only referenced by LOCA methods.

Returns
Reference to this object

Implemented in LOCA::Extended::Vector, NOX::Epetra::Vector, NOX::Thyra::Vector, NOX::Petsc::Vector, and NOX::LAPACK::Vector.

Referenced by LOCA::Homotopy::Group::Group(), LOCA::Hopf::MooreSpence::ExtendedGroup::init(), LOCA::Pitchfork::MooreSpence::ExtendedGroup::init(), and LOCA::TurningPoint::MooreSpence::ExtendedGroup::init().

virtual NOX::Abstract::Vector& NOX::Abstract::Vector::reciprocal ( const NOX::Abstract::Vector y)
pure virtual

Put element-wise reciprocal of source vector y into this vector.

Here x represents this vector, and we update it as

\[ x_i = \frac{1}{y_i} \quad \mbox{for } i=1,\dots,n \]

Returns
Reference to this object

Implemented in LOCA::Extended::Vector, NOX::Epetra::Vector, NOX::Thyra::Vector, NOX::Petsc::Vector, and NOX::LAPACK::Vector.

Referenced by LOCA::Bifurcation::TPBord::StatusTest::NullVectorNormWRMS::checkStatus(), LOCA::Bifurcation::PitchforkBord::StatusTest::NullVectorNormWRMS::checkStatus(), and NOX::StatusTest::NormWRMS::checkStatus().

virtual NOX::Abstract::Vector& NOX::Abstract::Vector::scale ( double  gamma)
pure virtual

Scale each element of this vector by gamma.

Here x represents this vector, and we update it as

\[ x_i = \gamma x_i \quad \mbox{for } i=1,\dots,n \]

Returns
Reference to this object

Implemented in LOCA::Extended::Vector, NOX::Epetra::Vector, NOX::Thyra::Vector, NOX::Petsc::Vector, and NOX::LAPACK::Vector.

Referenced by LOCA::Bifurcation::TPBord::StatusTest::NullVectorNormWRMS::checkStatus(), LOCA::Bifurcation::PitchforkBord::StatusTest::NullVectorNormWRMS::checkStatus(), NOX::StatusTest::NormWRMS::checkStatus(), NOX::LineSearch::SafeguardedStep::compute(), LOCA::MultiPredictor::Secant::compute(), NOX::Direction::SteepestDescent::compute(), NOX::Direction::NonlinearCG::compute(), NOX::Direction::Broyden::compute(), LOCA::DerivUtils::computeDCeDp(), LOCA::DerivUtils::computeDfDp(), LOCA::DerivUtils::computeDJnDp(), LOCA::DerivUtils::computeDwtCeDx(), LOCA::DerivUtils::computeDwtJDp(), LOCA::DerivUtils::computeDwtJnDx(), NOX::MeritFunction::SumOfSquares::computeQuadraticMinimizer(), LOCA::Epetra::Group::computeScaledDotProduct(), NOX::Thyra::WeightedMeritFunction::computeSlope(), LOCA::TurningPoint::MooreSpence::ExtendedGroup::getLengthVector(), LOCA::Hopf::MooreSpence::ExtendedGroup::init(), LOCA::Pitchfork::MooreSpence::ExtendedGroup::init(), LOCA::TurningPoint::MooreSpence::ExtendedGroup::init(), NOX::Solver::InexactTrustRegionBased::iterateInexact(), NOX::Solver::AndersonAcceleration::qrAdd(), LOCA::AnasaziOperator::ShiftInvert::rayleighQuotient(), LOCA::AnasaziOperator::ShiftInvert2Matrix::rayleighQuotient(), LOCA::AnasaziOperator::Cayley::rayleighQuotient(), LOCA::AnasaziOperator::Cayley2Matrix::rayleighQuotient(), LOCA::Extended::Vector::scale(), LOCA::TurningPoint::MooreSpence::ExtendedGroup::scaleNullVector(), LOCA::TurningPoint::MinimallyAugmented::Constraint::scaleNullVectors(), LOCA::Thyra::Group::scaleVector(), LOCA::LAPACK::Group::scaleVector(), and LOCA::Epetra::Group::scaleVector().

virtual NOX::Abstract::Vector& NOX::Abstract::Vector::scale ( const NOX::Abstract::Vector a)
pure virtual

Scale this vector element-by-element by the vector a.

Here x represents this vector, and we update it as

\[ x_i = x_i \cdot a_i \quad \mbox{for } i=1,\dots,n \]

Returns
Reference to this object

Implemented in LOCA::Extended::Vector, NOX::Epetra::Vector, NOX::Thyra::Vector, NOX::Petsc::Vector, and NOX::LAPACK::Vector.

virtual NOX::Abstract::Vector& NOX::Abstract::Vector::update ( double  alpha,
const NOX::Abstract::Vector a,
double  gamma = 0.0 
)
pure virtual

Compute x = (alpha * a) + (gamma * x) where x is this vector.

Here x represents this vector, and we update it as

\[ x_i = \alpha \; a_i + \gamma \; x_i \quad \mbox{for } i=1,\dots,n \]

Returns
Reference to this object

Implemented in LOCA::Extended::Vector, NOX::Epetra::Vector, NOX::Thyra::Vector, NOX::Petsc::Vector, and NOX::LAPACK::Vector.

Referenced by LOCA::Homotopy::Group::applyJacobian(), LOCA::Homotopy::Group::applyJacobianTranspose(), LOCA::Bifurcation::TPBord::StatusTest::NullVectorNormWRMS::checkStatus(), LOCA::Bifurcation::PitchforkBord::StatusTest::NullVectorNormWRMS::checkStatus(), NOX::StatusTest::NormUpdate::checkStatus(), NOX::StatusTest::NormWRMS::checkStatus(), LOCA::SingularJacobianSolve::ItRef::compute(), LOCA::SingularJacobianSolve::Nic::compute(), LOCA::SingularJacobianSolve::NicDay::compute(), NOX::Direction::NonlinearCG::compute(), NOX::Direction::Broyden::compute(), NOX::Solver::TensorBased::computeCurvilinearStep(), LOCA::DerivUtils::computeDCeDp(), LOCA::DerivUtils::computeDfDp(), LOCA::DerivUtils::computeDJnDp(), LOCA::DerivUtils::computeDwtCeDx(), LOCA::DerivUtils::computeDwtJDp(), LOCA::DerivUtils::computeDwtJnDx(), NOX::Direction::Utils::InexactNewton::computeForcingTerm(), LOCA::SingularJacobianSolve::ItRef::computeMulti(), LOCA::SingularJacobianSolve::Nic::computeMulti(), LOCA::SingularJacobianSolve::NicDay::computeMulti(), NOX::Thyra::WeightedMeritFunction::computeSlope(), NOX::Solver::TensorBased::computeTensorDirection(), NOX::Solver::TensorBased::getNormModelResidual(), LOCA::Homotopy::Group::Group(), LOCA::Hopf::MooreSpence::ExtendedGroup::init(), LOCA::Pitchfork::MooreSpence::ExtendedGroup::init(), LOCA::TurningPoint::MooreSpence::ExtendedGroup::init(), NOX::Solver::InexactTrustRegionBased::iterateInexact(), NOX::Solver::InexactTrustRegionBased::iterateStandard(), LOCA::DerivUtils::perturbXVec(), NOX::Solver::AndersonAcceleration::qrAdd(), NOX::Direction::Newton::resetForcingTerm(), LOCA::TurningPoint::MooreSpence::SalingerBordering::solveTranspose(), LOCA::TurningPoint::MooreSpence::PhippsBordering::solveTranspose(), and NOX::Solver::TrustRegionBased::step().

virtual NOX::Abstract::Vector& NOX::Abstract::Vector::update ( double  alpha,
const NOX::Abstract::Vector a,
double  beta,
const NOX::Abstract::Vector b,
double  gamma = 0.0 
)
pure virtual

Compute x = (alpha * a) + (beta * b) + (gamma * x) where x is this vector.

Here x represents this vector, and we update it as

\[ x_i = \alpha \; a_i + \beta \; b_i + \gamma \; x_i \quad \mbox{for } i=1,\dots,n \]

Returns
Reference to this object

Implemented in LOCA::Extended::Vector, NOX::Epetra::Vector, NOX::Thyra::Vector, NOX::Petsc::Vector, and NOX::LAPACK::Vector.


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