44 #ifndef ROL_STDVECTOR_H
45 #define ROL_STDVECTOR_H
49 #include <initializer_list>
59 template <
class Real,
class Element=Real>
73 std_vec_( makePtr<std::vector<Element>>(ilist) ) {}
81 std::invalid_argument,
82 "Error: Vectors must have the same dimension." );
85 const std::vector<Element>& xval = *ex.
getVector();
86 std::copy(xval.begin(),xval.end(),
std_vec_->begin());
92 std::invalid_argument,
93 "Error: Vectors must have the same dimension." );
96 const std::vector<Element>& xval = *ex.
getVector();
99 (*std_vec_)[i] += xval[i];
106 std::invalid_argument,
107 "Error: Vectors must have the same dimension." );
110 const std::vector<Element>& xval = *ex.
getVector();
113 (*std_vec_)[i] += alpha*xval[i];
118 for(
auto& e : *
std_vec_ ) e *= alpha;
128 std::invalid_argument,
129 "Error: Vectors must have the same dimension." );
132 const std::vector<Element>& xval = *ex.
getVector();
139 return std::inner_product(
std_vec_->begin(),
std_vec_->end(), xval.begin(), Real(0));
144 val = std::sqrt(
dot(*
this) );
148 virtual Ptr<Vector<Real> >
clone()
const {
149 return makePtr<StdVector>( makePtr<std::vector<Element>>(
std_vec_->size(),
static_cast<Element
>(0)));
160 Ptr<Vector<Real> >
basis(
const int i )
const {
162 ROL_TEST_FOR_EXCEPTION( i >=
dimension() || i<0,
163 std::invalid_argument,
164 "Error: Basis index must be between 0 and vector dimension." );
166 (*staticPtrCast<StdVector>(e)->
getVector())[i] = 1.0;
171 return static_cast<int>(
std_vec_->size());
174 void applyUnary(
const Elementwise::UnaryFunction<Real> &f ) {
179 for(
auto& e : *
std_vec_ ) e = f.apply(e);
185 std::invalid_argument,
186 "Error: Vectors must have the same dimension." );
189 const std::vector<Element>& xval = *ex.
getVector();
192 (*std_vec_)[i] = f.apply((*
std_vec_)[i],xval[i]);
197 Real
reduce(
const Elementwise::ReductionOp<Real> &r )
const {
198 Real result = r.initialValue();
211 void randomize(
const Real l = 0.0,
const Real u = 1.0 ) {
220 auto get_rand = [a,b]( Real& e ) {
221 auto x =
static_cast<Real
>(rand())/static_cast<Real>(RAND_MAX);
227 virtual void print( std::ostream &outStream )
const {
232 for(
auto e : *
std_vec_ ) outStream << e <<
" ";
233 outStream << std::endl;
typename PV< Real >::size_type size_type
void axpy(const Real alpha, const Vector< Real > &x)
Compute where .
void scale(const Real alpha)
Compute where .
StdVector(const int dim, const Element val=0.0)
virtual int dimension() const
Return dimension of the vector space.
Ptr< Vector< Real > > basis(const int i) const
Return i-th basis vector.
StdVector(const Ptr< std::vector< Element >> &std_vec)
void randomize(const Real l=0.0, const Real u=1.0)
Set vector to be uniform random between [l,u].
Ptr< const std::vector< Element > > getVector() const
typename std::vector< Real >::size_type size_type
void applyBinary(const Elementwise::BinaryFunction< Real > &f, const Vector< Real > &x)
Defines the linear algebra or vector space interface.
virtual void print(std::ostream &outStream) const
Real norm() const
Returns where .
Provides the ROL::Vector interface for scalar values, to be used, for example, with scalar constraint...
void plus(const Vector< Real > &x)
Compute , where .
virtual Ptr< Vector< Real > > clone() const
Clone to make a new (uninitialized) vector.
Real reduce(const Elementwise::ReductionOp< Real > &r) const
int dimension() const
Return dimension of the vector space.
void set(const Vector< Real > &x)
Set where .
void setScalar(const Real C)
Set where .
Ptr< std::vector< Element > > std_vec_
void applyUnary(const Elementwise::UnaryFunction< Real > &f)
virtual Real dot(const Vector< Real > &x) const
Compute where .
Ptr< std::vector< Element > > getVector()
StdVector(std::initializer_list< Element > ilist)