48#include "ROL_LinearAlgebra.hpp"
49#include "ROL_LAPACK.hpp"
83 Ptr<Elementwise::NormalRandom<Real>>
nrand_;
85 Ptr<std::normal_distribution<Real>>
dist_;
93 int K = std::min(M,N);
95 std::vector<Real> TAU(K);
96 std::vector<Real> WORK(1);
99 lapack_.GEQRF(M,N,
X_.values(),LDA,&TAU[0],&WORK[0],LWORK,&INFO);
100 LWORK =
static_cast<int>(WORK[0]);
102 lapack_.GEQRF(M,N,
X_.values(),LDA,&TAU[0],&WORK[0],LWORK,&INFO);
105 lapack_.ORGQR(M,N,K,
X_.values(),LDA,&TAU[0],&WORK[0],LWORK,&INFO);
106 LWORK =
static_cast<int>(WORK[0]);
108 lapack_.ORGQR(M,N,K,
X_.values(),LDA,&TAU[0],&WORK[0],LWORK,&INFO);
117 std::vector<Real> normQ(
k_,0);
119 for (
int j = 0; j <
k_; ++j) {
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]);
127 normQ[j] = Y[j]->norm();
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]) {
142 Y[j]->scale(one/rjj);
154 int LSsolver(LA::Matrix<Real> &A, LA::Matrix<Real> &B,
const bool trans =
false)
const {
156 char TRANS = (trans ?
'T' :
'N');
159 int NRHS = B.numCols();
161 int LDB = std::max(M,N);
162 std::vector<Real> WORK(1);
165 lapack_.GELS(TRANS,M,N,NRHS,A.values(),LDA,B.values(),LDB,&WORK[0],LWORK,&INFO);
167 LWORK =
static_cast<int>(WORK[0]);
169 lapack_.GELS(TRANS,M,N,NRHS,A.values(),LDA,B.values(),LDB,&WORK[0],LWORK,&INFO);
180 int K = std::min(M,N);
182 std::vector<Real> S(K);
183 LA::Matrix<Real> U(M,K);
185 LA::Matrix<Real> VT(K,N);
187 std::vector<Real> WORK(1), WORK0(1);
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]);
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) {
197 for (
int k = 0; k < r; ++k) {
198 A(i,j) += S[k] * U(i,k) * VT(k,j);
206 int infoP(0), infoQ(0), infoLS1(0), infoLS2(0), infoLRA(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]);
216 for (
int k = 0; k <
ncol_; ++k) {
217 R(i,j) +=
Psi_(k,i) *
X_(k,j);
223 LA::Matrix<Real> Z(
s_,
k_);
224 for (
int i = 0; i <
k_; ++i) {
225 for (
int j = 0; j <
s_; ++j) {
230 for (
int i = 0; i <
k_; ++i) {
231 for (
int j = 0; j <
k_; ++j) {
242 return std::abs(infoP)+std::abs(infoQ)+std::abs(infoLS1)
243 +std::abs(infoLS2)+std::abs(infoLRA);
250 for (
int i = 0; i <
s_; ++i) {
252 for (
int j = 0; j <
ncol_; ++j) {
253 Psi_(j,i) = (*dist_)(*gen_);
256 for (
int i = 0; i <
k_; ++i) {
259 for (
int j = 0; j <
ncol_; ++j) {
260 Omega_(j,i) = (*dist_)(*gen_);
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)
298 nrand_ = makePtr<Elementwise::NormalRandom<Real>>(mu,sig,dom_seed);
299 unsigned seed = rng_seed;
301 seed = std::chrono::system_clock::now().time_since_epoch().count();
303 gen_ = makePtr<std::mt19937_64>(seed);
304 dist_ = makePtr<std::normal_distribution<Real>>(mu,sig);
314 for (
int i = 0; i <
k_; ++i) {
318 for (
int i = 0; i <
s_; ++i) {
332 int sold =
s_, kold =
k_;
338 for (
int i = sold; i <
s_; ++i) {
343 for (
int i = kold; i <
k_; ++i) {
344 Y_.push_back(
Y_[0]->clone());
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;
384 if ( col >=
ncol_ || col < 0 ) {
389 for (
int i = 0; i <
k_; ++i) {
391 for (
int j = 0; j <
ncol_; ++j) {
401 for (
int i = 0; i <
s_; ++i) {
403 for (
int j = 0; j <
s_; ++j) {
405 Z_(i,j) += nu*
Psi_(col,j)*hphi;
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;
425 if ( col >=
ncol_ || col < 0 ) {
449 for (
int i = 0; i <
k_; ++i) {
451 for (
int j = 0; j <
k_; ++j) {
452 coeff +=
C_(i,j) *
X_(col,j);
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;
466 bool test(
const int rank, std::ostream &outStream = std::cout,
const int verbosity = 0) {
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);
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]);
489 bool flagP =
testP(outStream, verbosity);
491 bool flagQ =
testQ(outStream, verbosity);
493 Real nerr(0), maxerr(0);
494 Ptr<Vector<Real>> err =
Y_[0]->clone();
495 for (
int i = 0; i <
ncol_; ++i) {
497 err->axpy(-one,*A[i]);
499 maxerr = (nerr > maxerr ? nerr : maxerr);
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);
509 return flagP & flagQ & (maxerr < tol ? true :
false);
515 bool testQ(std::ostream &outStream = std::cout,
const int verbosity = 0) {
518 Real qij(0), err(0), maxerr(0);
519 std::ios_base::fmtflags oflags(outStream.flags());
521 outStream << std::scientific << std::setprecision(3);
524 outStream << std::endl
525 <<
" Printing Q'Q...This should be approximately equal to I"
526 << std::endl << std::endl;
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);
534 outStream << std::setw(12) << std::right << qij;
538 outStream << std::endl;
545 outStream << std::endl <<
" TEST ORTHOGONALIZATION: Max Error = "
546 << std::setw(12) << std::right << maxerr
548 outStream.flags(oflags);
550 return (maxerr < tol ?
true :
false);
553 bool testP(std::ostream &outStream = std::cout,
const int verbosity = 0) {
556 Real qij(0), err(0), maxerr(0);
557 std::ios_base::fmtflags oflags(outStream.flags());
559 outStream << std::scientific << std::setprecision(3);
562 outStream << std::endl
563 <<
" Printing P'P...This should be approximately equal to I"
564 << std::endl << std::endl;
566 for (
int i = 0; i <
k_; ++i) {
567 for (
int j = 0; j <
k_; ++j) {
569 for (
int k = 0; k <
ncol_; ++k) {
570 qij +=
X_(k,i) *
X_(k,j);
572 err = (i==j ? std::abs(qij-one) : std::abs(qij));
573 maxerr = (err > maxerr ? err : maxerr);
575 outStream << std::setw(12) << std::right << qij;
579 outStream << std::endl;
586 outStream << std::endl <<
" TEST ORTHOGONALIZATION: Max Error = "
587 << std::setw(12) << std::right << maxerr
589 outStream.flags(oflags);
591 return (maxerr < tol ?
true :
false);
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0 zero)()
Contains definitions of custom data types in ROL.
Ptr< std::mt19937_64 > gen_
std::vector< Ptr< Vector< Real > > > Upsilon_
void mgs2(std::vector< Ptr< Vector< Real > > > &Y) const
void update(UpdateType type)
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)
LAPACK< int, Real > lapack_
std::vector< Ptr< Vector< Real > > > Y_
Ptr< std::normal_distribution< Real > > dist_
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_
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)
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.