ROL
ROL_Sketch.hpp
Go to the documentation of this file.
1// @HEADER
2// ************************************************************************
3//
4// Rapid Optimization Library (ROL) Package
5// Copyright (2014) Sandia Corporation
6//
7// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8// license for use of this work by or on behalf of the U.S. Government.
9//
10// Redistribution and use in source and binary forms, with or without
11// modification, are permitted provided that the following conditions are
12// met:
13//
14// 1. Redistributions of source code must retain the above copyright
15// notice, this list of conditions and the following disclaimer.
16//
17// 2. Redistributions in binary form must reproduce the above copyright
18// notice, this list of conditions and the following disclaimer in the
19// documentation and/or other materials provided with the distribution.
20//
21// 3. Neither the name of the Corporation nor the names of the
22// contributors may be used to endorse or promote products derived from
23// this software without specific prior written permission.
24//
25// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36//
37// Questions? Contact lead developers:
38// Drew Kouri (dpkouri@sandia.gov) and
39// Denis Ridzal (dridzal@sandia.gov)
40//
41// ************************************************************************
42// @HEADER
43
44#ifndef ROL_SKETCH_H
45#define ROL_SKETCH_H
46
47#include "ROL_Vector.hpp"
48#include "ROL_LinearAlgebra.hpp"
49#include "ROL_LAPACK.hpp"
50#include "ROL_UpdateType.hpp"
51#include "ROL_Types.hpp"
52#include <random>
53#include <chrono>
54
61
62namespace ROL {
63
64template <class Real>
65class Sketch {
66private:
67 std::vector<Ptr<Vector<Real>>> Upsilon_, Phi_, Y_;
68 LA::Matrix<Real> Omega_, Psi_, X_, Z_, C_;
69
71
72 const Real orthTol_;
73 const int orthIt_;
74
75 const bool truncate_;
76
77 LAPACK<int,Real> lapack_;
78
80
81 Ptr<std::ostream> out_;
82
83 Ptr<Elementwise::NormalRandom<Real>> nrand_;
84 Ptr<std::mt19937_64> gen_;
85 Ptr<std::normal_distribution<Real>> dist_;
86
87 int computeP(void) {
88 int INFO(0);
89 if (!flagP_) {
90 // Solve least squares problem using LAPACK
91 int M = ncol_;
92 int N = k_;
93 int K = std::min(M,N);
94 int LDA = M;
95 std::vector<Real> TAU(K);
96 std::vector<Real> WORK(1);
97 int LWORK = -1;
98 // Compute QR factorization of X
99 lapack_.GEQRF(M,N,X_.values(),LDA,&TAU[0],&WORK[0],LWORK,&INFO);
100 LWORK = static_cast<int>(WORK[0]);
101 WORK.resize(LWORK);
102 lapack_.GEQRF(M,N,X_.values(),LDA,&TAU[0],&WORK[0],LWORK,&INFO);
103 // Generate Q
104 LWORK = -1;
105 lapack_.ORGQR(M,N,K,X_.values(),LDA,&TAU[0],&WORK[0],LWORK,&INFO);
106 LWORK = static_cast<int>(WORK[0]);
107 WORK.resize(LWORK);
108 lapack_.ORGQR(M,N,K,X_.values(),LDA,&TAU[0],&WORK[0],LWORK,&INFO);
109 flagP_ = true;
110 }
111 return INFO;
112 }
113
114 void mgs2(std::vector<Ptr<Vector<Real>>> &Y) const {
115 const Real one(1);
116 Real rjj(0), rij(0);
117 std::vector<Real> normQ(k_,0);
118 bool flag(true);
119 for (int j = 0; j < k_; ++j) {
120 rjj = Y[j]->norm();
121 if (rjj > ROL_EPSILON<Real>()) { // Ignore update if Y[i] is zero.
122 for (int k = 0; k < orthIt_; ++k) {
123 for (int i = 0; i < j; ++i) {
124 rij = Y[i]->dot(*Y[j]);
125 Y[j]->axpy(-rij,*Y[i]);
126 }
127 normQ[j] = Y[j]->norm();
128 flag = true;
129 for (int i = 0; i < j; ++i) {
130 rij = std::abs(Y[i]->dot(*Y[j]));
131 if (rij > orthTol_*normQ[j]*normQ[i]) {
132 flag = false;
133 break;
134 }
135 }
136 if (flag) {
137 break;
138 }
139 }
140 }
141 rjj = normQ[j];
142 Y[j]->scale(one/rjj);
143 }
144 }
145
146 int computeQ(void) {
147 if (!flagQ_) {
148 mgs2(Y_);
149 flagQ_ = true;
150 }
151 return 0;
152 }
153
154 int LSsolver(LA::Matrix<Real> &A, LA::Matrix<Real> &B, const bool trans = false) const {
155 int flag(0);
156 char TRANS = (trans ? 'T' : 'N');
157 int M = A.numRows();
158 int N = A.numCols();
159 int NRHS = B.numCols();
160 int LDA = M;
161 int LDB = std::max(M,N);
162 std::vector<Real> WORK(1);
163 int LWORK = -1;
164 int INFO;
165 lapack_.GELS(TRANS,M,N,NRHS,A.values(),LDA,B.values(),LDB,&WORK[0],LWORK,&INFO);
166 flag += INFO;
167 LWORK = static_cast<int>(WORK[0]);
168 WORK.resize(LWORK);
169 lapack_.GELS(TRANS,M,N,NRHS,A.values(),LDA,B.values(),LDB,&WORK[0],LWORK,&INFO);
170 flag += INFO;
171 return flag;
172 }
173
174 int lowRankApprox(LA::Matrix<Real> &A, const int r) const {
175 const Real zero(0);
176 char JOBU = 'S';
177 char JOBVT = 'S';
178 int M = A.numRows();
179 int N = A.numCols();
180 int K = std::min(M,N);
181 int LDA = M;
182 std::vector<Real> S(K);
183 LA::Matrix<Real> U(M,K);
184 int LDU = M;
185 LA::Matrix<Real> VT(K,N);
186 int LDVT = K;
187 std::vector<Real> WORK(1), WORK0(1);
188 int LWORK = -1;
189 int INFO;
190 lapack_.GESVD(JOBU,JOBVT,M,N,A.values(),LDA,&S[0],U.values(),LDU,VT.values(),LDVT,&WORK[0],LWORK,&WORK0[0],&INFO);
191 LWORK = static_cast<int>(WORK[0]);
192 WORK.resize(LWORK);
193 lapack_.GESVD(JOBU,JOBVT,M,N,A.values(),LDA,&S[0],U.values(),LDU,VT.values(),LDVT,&WORK[0],LWORK,&WORK0[0],&INFO);
194 for (int i = 0; i < M; ++i) {
195 for (int j = 0; j < N; ++j) {
196 A(i,j) = zero;
197 for (int k = 0; k < r; ++k) {
198 A(i,j) += S[k] * U(i,k) * VT(k,j);
199 }
200 }
201 }
202 return INFO;
203 }
204
205 int computeC(void) {
206 int infoP(0), infoQ(0), infoLS1(0), infoLS2(0), infoLRA(0);
207 infoP = computeP();
208 infoQ = computeQ();
209 if (!flagC_) {
210 const Real zero(0);
211 LA::Matrix<Real> L(s_,k_), R(s_,k_);
212 for (int i = 0; i < s_; ++i) {
213 for (int j = 0; j < k_; ++j) {
214 L(i,j) = Phi_[i]->dot(*Y_[j]);
215 R(i,j) = zero;
216 for (int k = 0; k < ncol_; ++k) {
217 R(i,j) += Psi_(k,i) * X_(k,j);
218 }
219 }
220 }
221 // Solve least squares problems using LAPACK
222 infoLS1 = LSsolver(L,Z_,false);
223 LA::Matrix<Real> Z(s_,k_);
224 for (int i = 0; i < k_; ++i) {
225 for (int j = 0; j < s_; ++j) {
226 Z(j,i) = Z_(i,j);
227 }
228 }
229 infoLS2 = LSsolver(R,Z,false);
230 for (int i = 0; i < k_; ++i) {
231 for (int j = 0; j < k_; ++j) {
232 C_(j,i) = Z(i,j);
233 }
234 }
235 // Compute best rank r approximation
236 if (truncate_) {
237 infoLRA = lowRankApprox(C_,rank_);
238 }
239 // Set flag
240 flagC_ = true;
241 }
242 return std::abs(infoP)+std::abs(infoQ)+std::abs(infoLS1)
243 +std::abs(infoLS2)+std::abs(infoLRA);
244 }
245
246 void reset(void) {
247 const Real zero(0);
248 // Randomize Upsilon, Omega, Psi and Phi, and zero X, Y and Z
249 X_.scale(zero); Z_.scale(zero); C_.scale(zero);
250 for (int i = 0; i < s_; ++i) {
251 Phi_[i]->applyUnary(*nrand_);
252 for (int j = 0; j < ncol_; ++j) {
253 Psi_(j,i) = (*dist_)(*gen_);
254 }
255 }
256 for (int i = 0; i < k_; ++i) {
257 Y_[i]->zero();
258 Upsilon_[i]->applyUnary(*nrand_);
259 for (int j = 0; j < ncol_; ++j) {
260 Omega_(j,i) = (*dist_)(*gen_);
261 }
262 }
263 }
264
265// void reset(void) {
266// const Real zero(0), one(1), a(2), b(-1);
267// Real x(0);
268// // Randomize Upsilon, Omega, Psi and Phi, and zero X, Y and Z
269// X_.scale(zero); Z_.scale(zero); C_.scale(zero);
270// for (int i = 0; i < s_; ++i) {
271// Phi_[i]->randomize(-one,one);
272// for (int j = 0; j < ncol_; ++j) {
273// x = static_cast<Real>(rand())/static_cast<Real>(RAND_MAX);
274// Psi_(j,i) = a*x + b;
275// }
276// }
277// for (int i = 0; i < k_; ++i) {
278// Y_[i]->zero();
279// Upsilon_[i]->randomize(-one,one);
280// for (int j = 0; j < ncol_; ++j) {
281// x = static_cast<Real>(rand())/static_cast<Real>(RAND_MAX);
282// Omega_(j,i) = a*x + b;
283// }
284// }
285// }
286
287public:
288 virtual ~Sketch(void) {}
289
290 Sketch(const Vector<Real> &x, const int ncol, const int rank,
291 const Real orthTol = 1e-8, const int orthIt = 2,
292 const bool truncate = false,
293 const unsigned dom_seed = 0, const unsigned rng_seed = 0)
294 : ncol_(ncol), orthTol_(orthTol), orthIt_(orthIt),
295 truncate_(truncate), flagP_(false), flagQ_(false), flagC_(false),
296 out_(nullPtr) {
297 Real mu(0), sig(1);
298 nrand_ = makePtr<Elementwise::NormalRandom<Real>>(mu,sig,dom_seed);
299 unsigned seed = rng_seed;
300 if (seed == 0) {
301 seed = std::chrono::system_clock::now().time_since_epoch().count();
302 }
303 gen_ = makePtr<std::mt19937_64>(seed);
304 dist_ = makePtr<std::normal_distribution<Real>>(mu,sig);
305 // Compute reduced dimensions
306 maxRank_ = std::min(ncol_, x.dimension());
307 rank_ = std::min(rank, maxRank_);
308 k_ = std::min(2*rank_+1, maxRank_);
309 s_ = std::min(2*k_ +1, maxRank_);
310 // Initialize matrix storage
311 Upsilon_.clear(); Phi_.clear(); Y_.clear();
312 Omega_.reshape(ncol_,k_); Psi_.reshape(ncol_,s_);
313 X_.reshape(ncol_,k_); Z_.reshape(s_,s_); C_.reshape(k_,k_);
314 for (int i = 0; i < k_; ++i) {
315 Y_.push_back(x.clone());
316 Upsilon_.push_back(x.clone());
317 }
318 for (int i = 0; i < s_; ++i) {
319 Phi_.push_back(x.clone());
320 }
321 // Randomize Psi and Omega and zero W and Y
322 reset();
323 }
324
325 void setStream(Ptr<std::ostream> &out) {
326 out_ = out;
327 }
328
329 void setRank(const int rank) {
330 rank_ = std::min(rank, maxRank_);
331 // Compute reduced dimensions
332 int sold = s_, kold = k_;
333 k_ = std::min(2*rank_+1, maxRank_);
334 s_ = std::min(2*k_ +1, maxRank_);
335 Omega_.reshape(ncol_,k_); Psi_.reshape(ncol_,s_);
336 X_.reshape(ncol_,k_); Z_.reshape(s_,s_); C_.reshape(k_,k_);
337 if (s_ > sold) {
338 for (int i = sold; i < s_; ++i) {
339 Phi_.push_back(Phi_[0]->clone());
340 }
341 }
342 if (k_ > kold) {
343 for (int i = kold; i < k_; ++i) {
344 Y_.push_back(Y_[0]->clone());
345 Upsilon_.push_back(Upsilon_[0]->clone());
346 }
347 }
348 update();
349 if ( out_ != nullPtr ) {
350 *out_ << std::string(80,'=') << std::endl;
351 *out_ << " ROL::Sketch::setRank" << std::endl;
352 *out_ << " **** Rank = " << rank_ << std::endl;
353 *out_ << " **** k = " << k_ << std::endl;
354 *out_ << " **** s = " << s_ << std::endl;
355 *out_ << std::string(80,'=') << std::endl;
356 }
357 }
358
359 void update(void) {
360 flagP_ = false;
361 flagQ_ = false;
362 flagC_ = false;
363 // Randomize Psi and Omega and zero W and Y
364 reset();
365 }
366
367 void update(UpdateType type) {
368 switch(type) {
373 case UpdateType::Temp:
374 flagP_ = false; flagQ_ = false; flagC_ = false;
375 reset();
376 break;
377 default:
378 break;
379 }
380 }
381
382 int advance(const Real nu, Vector<Real> &h, const int col, const Real eta = 1.0) {
383 // Check to see if col is less than ncol_
384 if ( col >= ncol_ || col < 0 ) {
385 // Input column index out of range!
386 return 1;
387 }
388 if (!flagP_ && !flagQ_ && !flagC_) {
389 for (int i = 0; i < k_; ++i) {
390 // Update X
391 for (int j = 0; j < ncol_; ++j) {
392 X_(j,i) *= eta;
393 }
394 X_(col,i) += nu*h.dot(*Upsilon_[i]);
395 // Update Y
396 Y_[i]->scale(eta);
397 Y_[i]->axpy(nu*Omega_(col,i),h);
398 }
399 // Update Z
400 Real hphi(0);
401 for (int i = 0; i < s_; ++i) {
402 hphi = h.dot(*Phi_[i]);
403 for (int j = 0; j < s_; ++j) {
404 Z_(i,j) *= eta;
405 Z_(i,j) += nu*Psi_(col,j)*hphi;
406 }
407 }
408 if ( out_ != nullPtr ) {
409 *out_ << std::string(80,'=') << std::endl;
410 *out_ << " ROL::Sketch::advance" << std::endl;
411 *out_ << " **** col = " << col << std::endl;
412 *out_ << " **** norm(h) = " << h.norm() << std::endl;
413 *out_ << std::string(80,'=') << std::endl;
414 }
415 }
416 else {
417 // Reconstruct has already been called!
418 return 1;
419 }
420 return 0;
421 }
422
423 int reconstruct(Vector<Real> &a, const int col) {
424 // Check to see if col is less than ncol_
425 if ( col >= ncol_ || col < 0 ) {
426 // Input column index out of range!
427 return 2;
428 }
429 const Real zero(0);
430 int flag(0);
431 // Compute QR factorization of X store in X
432 flag = computeP();
433 if (flag > 0 ) {
434 return 3;
435 }
436 // Compute QR factorization of Y store in Y
437 flag = computeQ();
438 if (flag > 0 ) {
439 return 4;
440 }
441 // Compute (Phi Q)\Z/(Psi P)* store in C
442 flag = computeC();
443 if (flag > 0 ) {
444 return 5;
445 }
446 // Recover sketch
447 a.zero();
448 Real coeff(0);
449 for (int i = 0; i < k_; ++i) {
450 coeff = zero;
451 for (int j = 0; j < k_; ++j) {
452 coeff += C_(i,j) * X_(col,j);
453 }
454 a.axpy(coeff,*Y_[i]);
455 }
456 if ( out_ != nullPtr ) {
457 *out_ << std::string(80,'=') << std::endl;
458 *out_ << " ROL::Sketch::reconstruct" << std::endl;
459 *out_ << " **** col = " << col << std::endl;
460 *out_ << " **** norm(a) = " << a.norm() << std::endl;
461 *out_ << std::string(80,'=') << std::endl;
462 }
463 return 0;
464 }
465
466 bool test(const int rank, std::ostream &outStream = std::cout, const int verbosity = 0) {
467 const Real one(1), tol(std::sqrt(ROL_EPSILON<Real>()));
468 // Initialize low rank factors
469 std::vector<Ptr<Vector<Real>>> U(rank);
470 LA::Matrix<Real> V(ncol_,rank);
471 for (int i = 0; i < rank; ++i) {
472 U[i] = Y_[0]->clone();
473 U[i]->randomize(-one,one);
474 for (int j = 0; j < ncol_; ++j) {
475 V(j,i) = static_cast<Real>(rand())/static_cast<Real>(RAND_MAX);
476 }
477 }
478 // Initialize A and build sketch
479 update();
480 std::vector<Ptr<Vector<Real>>> A(ncol_);
481 for (int i = 0; i < ncol_; ++i) {
482 A[i] = Y_[0]->clone(); A[i]->zero();
483 for (int j = 0; j < rank; ++j) {
484 A[i]->axpy(V(i,j),*U[j]);
485 }
486 advance(one,*A[i],i,one);
487 }
488 // Test QR decomposition of X
489 bool flagP = testP(outStream, verbosity);
490 // Test QR decomposition of Y
491 bool flagQ = testQ(outStream, verbosity);
492 // Test reconstruction of A
493 Real nerr(0), maxerr(0);
494 Ptr<Vector<Real>> err = Y_[0]->clone();
495 for (int i = 0; i < ncol_; ++i) {
496 reconstruct(*err,i);
497 err->axpy(-one,*A[i]);
498 nerr = err->norm();
499 maxerr = (nerr > maxerr ? nerr : maxerr);
500 }
501 if (verbosity > 0) {
502 std::ios_base::fmtflags oflags(outStream.flags());
503 outStream << std::scientific << std::setprecision(3) << std::endl;
504 outStream << " TEST RECONSTRUCTION: Max Error = "
505 << std::setw(12) << std::right << maxerr
506 << std::endl << std::endl;
507 outStream.flags(oflags);
508 }
509 return flagP & flagQ & (maxerr < tol ? true : false);
510 }
511
512private:
513
514 // Test functions
515 bool testQ(std::ostream &outStream = std::cout, const int verbosity = 0) {
516 const Real one(1), tol(std::sqrt(ROL_EPSILON<Real>()));
517 computeQ();
518 Real qij(0), err(0), maxerr(0);
519 std::ios_base::fmtflags oflags(outStream.flags());
520 if (verbosity > 0) {
521 outStream << std::scientific << std::setprecision(3);
522 }
523 if (verbosity > 1) {
524 outStream << std::endl
525 << " Printing Q'Q...This should be approximately equal to I"
526 << std::endl << std::endl;
527 }
528 for (int i = 0; i < k_; ++i) {
529 for (int j = 0; j < k_; ++j) {
530 qij = Y_[i]->dot(*Y_[j]);
531 err = (i==j ? std::abs(qij-one) : std::abs(qij));
532 maxerr = (err > maxerr ? err : maxerr);
533 if (verbosity > 1) {
534 outStream << std::setw(12) << std::right << qij;
535 }
536 }
537 if (verbosity > 1) {
538 outStream << std::endl;
539 }
540 if (maxerr > tol) {
541 break;
542 }
543 }
544 if (verbosity > 0) {
545 outStream << std::endl << " TEST ORTHOGONALIZATION: Max Error = "
546 << std::setw(12) << std::right << maxerr
547 << std::endl;
548 outStream.flags(oflags);
549 }
550 return (maxerr < tol ? true : false);
551 }
552
553 bool testP(std::ostream &outStream = std::cout, const int verbosity = 0) {
554 const Real zero(0), one(1), tol(std::sqrt(ROL_EPSILON<Real>()));
555 computeP();
556 Real qij(0), err(0), maxerr(0);
557 std::ios_base::fmtflags oflags(outStream.flags());
558 if (verbosity > 0) {
559 outStream << std::scientific << std::setprecision(3);
560 }
561 if (verbosity > 1) {
562 outStream << std::endl
563 << " Printing P'P...This should be approximately equal to I"
564 << std::endl << std::endl;
565 }
566 for (int i = 0; i < k_; ++i) {
567 for (int j = 0; j < k_; ++j) {
568 qij = zero;
569 for (int k = 0; k < ncol_; ++k) {
570 qij += X_(k,i) * X_(k,j);
571 }
572 err = (i==j ? std::abs(qij-one) : std::abs(qij));
573 maxerr = (err > maxerr ? err : maxerr);
574 if (verbosity > 1) {
575 outStream << std::setw(12) << std::right << qij;
576 }
577 }
578 if (verbosity > 1) {
579 outStream << std::endl;
580 }
581 if (maxerr > tol) {
582 break;
583 }
584 }
585 if (verbosity > 0) {
586 outStream << std::endl << " TEST ORTHOGONALIZATION: Max Error = "
587 << std::setw(12) << std::right << maxerr
588 << std::endl;
589 outStream.flags(oflags);
590 }
591 return (maxerr < tol ? true : false);
592 }
593
594}; // class Sketch
595
596} // namespace ROL
597
598#endif
Vector< Real > V
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0 zero)()
Contains definitions of custom data types in ROL.
int computeC(void)
LA::Matrix< Real > C_
Ptr< std::mt19937_64 > gen_
std::vector< Ptr< Vector< Real > > > Upsilon_
void mgs2(std::vector< Ptr< Vector< Real > > > &Y) const
LA::Matrix< Real > X_
LA::Matrix< Real > Psi_
void update(UpdateType type)
virtual ~Sketch(void)
int reconstruct(Vector< Real > &a, const int col)
bool testP(std::ostream &outStream=std::cout, const int verbosity=0)
LA::Matrix< Real > Omega_
int advance(const Real nu, Vector< Real > &h, const int col, const Real eta=1.0)
const Real orthTol_
LAPACK< int, Real > lapack_
Ptr< std::ostream > out_
void update(void)
std::vector< Ptr< Vector< Real > > > Y_
Ptr< std::normal_distribution< Real > > dist_
LA::Matrix< Real > Z_
int computeP(void)
void setStream(Ptr< std::ostream > &out)
int LSsolver(LA::Matrix< Real > &A, LA::Matrix< Real > &B, const bool trans=false) const
Ptr< Elementwise::NormalRandom< Real > > nrand_
const bool truncate_
std::vector< Ptr< Vector< Real > > > Phi_
bool testQ(std::ostream &outStream=std::cout, const int verbosity=0)
int lowRankApprox(LA::Matrix< Real > &A, const int r) const
void setRank(const int rank)
bool test(const int rank, std::ostream &outStream=std::cout, const int verbosity=0)
const int orthIt_
int computeQ(void)
void reset(void)
Sketch(const Vector< Real > &x, const int ncol, const int rank, const Real orthTol=1e-8, const int orthIt=2, const bool truncate=false, const unsigned dom_seed=0, const unsigned rng_seed=0)
Defines the linear algebra or vector space interface.
virtual Real norm() const =0
Returns where .
virtual void zero()
Set to zero vector.
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
virtual int dimension() const
Return dimension of the vector space.
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
virtual Real dot(const Vector &x) const =0
Compute where .
Real ROL_EPSILON(void)
Platform-dependent machine epsilon.
Definition ROL_Types.hpp:91