ROL
ROL_ProfiledVector.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Rapid Optimization Library (ROL) Package
4//
5// Copyright 2014 NTESS and the ROL contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef ROL_VECTORPROFILER_H
11#define ROL_VECTORPROFILER_H
12
13#include "ROL_Vector.hpp"
14#include <ostream>
15
16namespace ROL {
17
35template<class Ordinal>
37 Ordinal constructor_;
38 Ordinal destructor_;
39 Ordinal plus_;
40 Ordinal scale_;
41 Ordinal dot_;
42 Ordinal norm_;
43 Ordinal clone_;
44 Ordinal axpy_;
45 Ordinal zero_;
46 Ordinal basis_;
47 Ordinal dimension_;
48 Ordinal set_;
49 Ordinal dual_;
50 Ordinal apply_;
51 Ordinal applyUnary_;
52 Ordinal applyBinary_;
53 Ordinal reduce_;
54 Ordinal setScalar_;
55 Ordinal randomize_;
57 constructor_(0), destructor_(0), plus_(0), scale_(0), dot_(0), norm_(0), clone_(0),
58 axpy_(0), zero_(0), basis_(0), dimension_(0), set_(0), dual_(0), apply_(0),
60
61}; // struct VectorFunctionCalls
62
63
64// Forward declaration for friend functions
65template<class Ordinal,class Real>
66class ProfiledVector;
67
68template<class Ordinal,class Real>
72
73template<class Ordinal, class Real>
74void printVectorFunctionCalls( const ProfiledVector<Ordinal,Real> &x, std::ostream &outStream = std::cout ) {
75 outStream << "Total Vector Function Calls" << std::endl;
76 outStream << "---------------------------" << std::endl;
77 outStream << "Constructor : " << x.functionCalls_.constructor_ << std::endl;
78 outStream << "Destructor : " << x.functionCalls_.destructor_ << std::endl;
79 outStream << "set : " << x.functionCalls_.set_ << std::endl;
80 outStream << "plus : " << x.functionCalls_.plus_ << std::endl;
81 outStream << "axpy : " << x.functionCalls_.axpy_ << std::endl;
82 outStream << "scale : " << x.functionCalls_.scale_ << std::endl;
83 outStream << "dot : " << x.functionCalls_.dot_ << std::endl;
84 outStream << "zero : " << x.functionCalls_.zero_ << std::endl;
85 outStream << "norm : " << x.functionCalls_.norm_ << std::endl;
86 outStream << "clone : " << x.functionCalls_.clone_ << std::endl;
87 outStream << "basis : " << x.functionCalls_.basis_ << std::endl;
88 outStream << "dual : " << x.functionCalls_.dual_ << std::endl;
89 outStream << "apply : " << x.functionCalls_.apply_ << std::endl;
90 outStream << "dimension : " << x.functionCalls_.dimension_ << std::endl;
91 outStream << "applyUnary : " << x.functionCalls_.applyUnary_ << std::endl;
92 outStream << "applyBinary : " << x.functionCalls_.applyBinary_ << std::endl;
93 outStream << "reduce : " << x.functionCalls_.reduce_ << std::endl;
94 outStream << "setScalar : " << x.functionCalls_.setScalar_ << std::endl;
95 outStream << "randomize : " << x.functionCalls_.randomize_ << std::endl;
96}
97
98
99
100template<class Ordinal,class Real>
101class ProfiledVector : public Vector<Real> {
102
104
105private:
106 ROL::Ptr<Vector<Real> > v_;
108public:
109
110 ProfiledVector( const ROL::Ptr<Vector<Real> > &v ) {
111 // Make sure that given vector is not itself a ProfiledVector to avoid recursion
112 ROL::Ptr<ProfiledVector> pv = ROL::nullPtr;
113 pv = ROL::dynamicPtrCast<ProfiledVector>(v);
114 ROL_TEST_FOR_EXCEPTION( pv != ROL::nullPtr, std::logic_error, "ProfiledVector class "
115 "cannot encapsulate a ProfiledVector object!");
116
117 v_ = v;
118
119 functionCalls_.constructor_++;
120 }
121
122 virtual ~ProfiledVector() {
123 functionCalls_.destructor_++;
124 }
125
126 void plus( const Vector<Real> &x ) {
127 ROL::Ptr<const V> xp = dynamic_cast<const ProfiledVector&>(x).getVector();
128
129 functionCalls_.plus_++;
130 v_->plus(*xp);
131 }
132
133 void scale( const Real alpha ) {
134 functionCalls_.scale_++;
135 v_->scale(alpha);
136 }
137
138 Real dot( const Vector<Real> &x ) const {
139 ROL::Ptr<const V> xp = dynamic_cast<const ProfiledVector&>(x).getVector();
140 functionCalls_.dot_++;
141 return v_->dot(*xp);
142 }
143
144 Real norm() const {
145 functionCalls_.norm_++;
146 return v_->norm();
147 }
148
149 ROL::Ptr<Vector<Real> > clone() const {
150 functionCalls_.clone_++;
151 return ROL::makePtr<ProfiledVector>( v_->clone() );
152 }
153
154 void axpy( const Real alpha, const Vector<Real> &x ) {
155 ROL::Ptr<const V> xp = dynamic_cast<const ProfiledVector&>(x).getVector();
156 functionCalls_.axpy_++;
157 return v_->axpy(alpha,*xp);
158 }
159
160 void zero() {
161 functionCalls_.zero_++;
162 v_->zero();
163 }
164
165 ROL::Ptr<Vector<Real> > basis( const int i ) const {
166 functionCalls_.basis_++;
167 return ROL::makePtr<ProfiledVector>( v_->basis(i) );
168 }
169
170 int dimension() const {
171 functionCalls_.dimension_++;
172 return v_->dimension();
173 }
174
175 void set( const Vector<Real> &x ) {
176 ROL::Ptr<const V> xp = dynamic_cast<const ProfiledVector&>(x).getVector();
177 functionCalls_.set_++;
178 v_->set(*xp);
179 }
180
181 // TODO: determine the correct way to handle dual when v_ is a generic ROL::Ptr<ROL::Vector>
182 const Vector<Real> & dual() const {
183 functionCalls_.dual_++;
184 return *this;
185 }
186
187 Real apply(const Vector<Real> &x) const {
188 functionCalls_.apply_++;
189 return v_->apply(x);
190 }
191
192 ROL::Ptr<Vector<Real> > getVector() {
193 return v_;
194 }
195
196 ROL::Ptr<const Vector<Real> > getVector() const {
197 return v_;
198 }
199
200 void applyUnary( const Elementwise::UnaryFunction<Real> &f ) {
201 functionCalls_.applyUnary_++;
202 v_->applyUnary(f);
203 }
204
205 void applyBinary( const Elementwise::BinaryFunction<Real> &f, const Vector<Real> &x ) {
206 functionCalls_.applyBinary_++;
207 v_->applyBinary(f,x);
208 }
209
210 Real reduce( const Elementwise::ReductionOp<Real> &r ) const {
211 functionCalls_.reduce_++;
212 return v_->reduce(r);
213 }
214
215 void setScalar( const Real C ) {
216 functionCalls_.setScalar_++;
217 v_->setScalar(C);
218 }
219
220 void randomize( const Real l=0.0, const Real u=1.0) {
221 functionCalls_.randomize_++;
222 v_->randomize(l,u);
223 }
224
225 void print( std::ostream &outStream ) const {
226 v_->print(outStream);
227 }
228
229 friend VectorFunctionCalls<Ordinal> getVectorFunctionCalls<>( const ProfiledVector<Ordinal,Real> & );
230 friend void printVectorFunctionCalls<>( const ProfiledVector<Ordinal,Real> &, std::ostream & );
231
232};
233
234
235} // namespace ROL
236
237#endif // ROL_RANDOMVECTOR_H
By keeping a pointer to this in a derived Vector class, a tally of all methods is kept for profiling ...
void set(const Vector< Real > &x)
Set where .
ROL::Ptr< const Vector< Real > > getVector() const
Real apply(const Vector< Real > &x) const
Apply to a dual vector. This is equivalent to the call .
Real reduce(const Elementwise::ReductionOp< Real > &r) const
int dimension() const
Return dimension of the vector space.
ROL::Ptr< Vector< Real > > basis(const int i) const
Return i-th basis vector.
void scale(const Real alpha)
Compute where .
ROL::Ptr< Vector< Real > > clone() const
Clone to make a new (uninitialized) vector.
void randomize(const Real l=0.0, const Real u=1.0)
Set vector to be uniform random between [l,u].
void plus(const Vector< Real > &x)
Compute , where .
Real dot(const Vector< Real > &x) const
Compute where .
ProfiledVector(const ROL::Ptr< Vector< Real > > &v)
const Vector< Real > & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis,...
void setScalar(const Real C)
Set where .
void axpy(const Real alpha, const Vector< Real > &x)
Compute where .
void print(std::ostream &outStream) const
Real norm() const
Returns where .
ROL::Ptr< Vector< Real > > getVector()
ROL::Ptr< Vector< Real > > v_
void applyUnary(const Elementwise::UnaryFunction< Real > &f)
static VectorFunctionCalls< Ordinal > functionCalls_
void applyBinary(const Elementwise::BinaryFunction< Real > &f, const Vector< Real > &x)
void zero()
Set to zero vector.
Defines the linear algebra or vector space interface.
void printVectorFunctionCalls(const ProfiledVector< Ordinal, Real > &x, std::ostream &outStream=std::cout)
VectorFunctionCalls< Ordinal > getVectorFunctionCalls(const ProfiledVector< Ordinal, Real > &x)