UQTk: Uncertainty Quantification Toolkit 3.1.5
gproc.h
Go to the documentation of this file.
1/* =====================================================================================
2
3 The UQ Toolkit (UQTk) version 3.1.5
4 Copyright (2024) NTESS
5 https://www.sandia.gov/UQToolkit/
6 https://github.com/sandialabs/UQTk
7
8 Copyright 2024 National Technology & Engineering Solutions of Sandia, LLC (NTESS).
9 Under the terms of Contract DE-NA0003525 with NTESS, the U.S. Government
10 retains certain rights in this software.
11
12 This file is part of The UQ Toolkit (UQTk)
13
14 UQTk is open source software: you can redistribute it and/or modify
15 it under the terms of BSD 3-Clause License
16
17 UQTk is distributed in the hope that it will be useful,
18 but WITHOUT ANY WARRANTY; without even the implied warranty of
19 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
20 BSD 3 Clause License for more details.
21
22 You should have received a copy of the BSD 3 Clause License
23 along with UQTk. If not, see https://choosealicense.com/licenses/bsd-3-clause/.
24
25 Questions? Contact the UQTk Developers at https://github.com/sandialabs/UQTk/discussions
26 Sandia National Laboratories, Livermore, CA, USA
27===================================================================================== */
31
32#ifndef GPROC_H_SEEN
33#define GPROC_H_SEEN
34
35#include "Array1D.h"
36#include "Array2D.h"
37#include "PCSet.h"
38
41class Gproc {
42public:
43
46 Gproc(const string covtype, PCSet *PCModel, Array1D<double>& param);
48 ~Gproc() {};
49
51 void SetupPrior();
53 void SetupData(Array2D<double>& xdata, Array1D<double>& ydata,Array1D<double>& datavar);
55 void setCorrParam(Array1D<double> param){param_=param; return;}
56
59 void BuildGP();
64 void BuildGP_inv();
67 void EvalGP(Array2D<double>& xgrid, string msc, Array1D<double>& mst);
72 void EvalGP_inv(Array2D<double>& xgrid, string msc, Array1D<double>& mst);
74 int getNpt() const {return npt_;}
76 int getNdim() const {return ndim_;}
78 int getNPC() const {return npc_;}
80 double getAl() const {return al_;}
82 double getBe() const {return be_;}
84 double getSig2hat() const {return sig2hat_;}
86 void getVst(Array2D<double>& vst) {vst=Vst_; return;}
88 void getA(Array2D<double>& acor) {acor=A_; return;}
90 void getParam(Array1D<double>& param) {param=param_; return;}
92 void getCov(Array2D<double>& cov) {cov=cov_;}
94 void getVar(Array1D<double>& var) {var=var_;}
96 void getXYCov(Array2D<double>& xgrid,Array2D<double>& xycov);
98 void getSttPars(Array1D<double>& sttmat);
100 void findBestCorrParam();
101
102
103 private:
104
109
116
118 int npc_;
123 //double sig2f_;
125 double al_;
127 double be_;
129 double sig2hat_;
130
132 int npt_;
134 int ndim_;
136 string covType_;
139
148
150
165
166
167};
168#endif /* GPROC_H_SEEN */
1D Array class for any type T
2D Array class for any type T
Header file for the Multivariate PC class.
Stores data of any type T in a 1D array.
Definition Array1D.h:61
Stores data of any type T in a 2D array.
Definition Array2D.h:60
Class for Gaussian processes.
Definition gproc.h:41
Array1D< double > bhat_
Definition gproc.h:159
Array1D< double > ydata_
ydata array
Definition gproc.h:113
void getVst(Array2D< double > &vst)
Get , an auxiliary matrix.
Definition gproc.h:86
Array2D< double > Vst_
Definition gproc.h:158
void EvalGP(Array2D< double > &xgrid, string msc, Array1D< double > &mst)
Evaluate the Gaussian Process at a given grid msc controls whether only mean will be computed,...
Definition gproc.cpp:205
Array2D< double > Ht_
Definition gproc.h:151
double be_
Prior parameter .
Definition gproc.h:127
void BuildGP()
Build Gaussian Process regressor, i.e. compute internally all necessary matrices and vectors that des...
Definition gproc.cpp:99
void getCov(Array2D< double > &cov)
Get the posterior covariance matrix.
Definition gproc.h:92
void getParam(Array1D< double > &param)
Get the roughness parameters.
Definition gproc.h:90
void setCorrParam(Array1D< double > param)
Set the roughness parameter vector.
Definition gproc.h:55
void findBestCorrParam()
Function to find the best values for roughness parameters.
Definition gproc.cpp:511
Array2D< double > xdata_
xdata array
Definition gproc.h:111
Array2D< double > A_
Definition gproc.h:152
Array1D< double > HtAinvd_
Definition gproc.h:155
void SetupPrior()
Setup the prior.
Definition gproc.cpp:68
Array2D< double > AinvH_
Definition gproc.h:156
double sig2hat_
Posterior variance factor.
Definition gproc.h:129
double getSig2hat() const
Get Sigma-hat-squared, i.e. the posterior variance factor.
Definition gproc.h:84
string covType_
Covariance type, only 'SqExp' implemented so far.
Definition gproc.h:136
void BuildGP_inv()
Build Gaussian Process regressor, i.e. compute internally all necessary matrices and vectors that des...
Definition gproc.cpp:155
Array1D< double > Vinvz_
Definition gproc.h:154
int getNPC() const
Get the number of basis terms in the trend.
Definition gproc.h:78
Gproc(const string covtype, PCSet *PCModel, Array1D< double > &param)
Constructor: initialize with covariance type, trend function basis and roughness parameter vector.
Definition gproc.cpp:52
void getA(Array2D< double > &acor)
Get the correlation matrix .
Definition gproc.h:88
Array2D< double > Vstinv_
Definition gproc.h:163
Array1D< double > Ainvd_
Definition gproc.h:153
void getXYCov(Array2D< double > &xgrid, Array2D< double > &xycov)
Get the covariance in a different format, with the x,x' values.
Definition gproc.cpp:435
void getVar(Array1D< double > &var)
Get the posterior variance vector.
Definition gproc.h:94
Array2D< double > Ainv_
Definition gproc.h:152
int ndim_
Dimensionality.
Definition gproc.h:134
void computeDataCov_(Array2D< double > &xdata, Array1D< double > &param, Array2D< double > &A)
Compute the data covariance .
Definition gproc.cpp:489
Array1D< double > yHbhat_
Definition gproc.h:161
Array2D< double > H_
Auxiliary matrices or vectors, see the UQTk Manual.
Definition gproc.h:151
Array2D< double > Vinv_
Inverse of the mean trend coefficient prior covariance.
Definition gproc.h:120
void getSttPars(Array1D< double > &sttmat)
Get the Student-t parameters.
Definition gproc.cpp:418
int npc_
Number of bases in the mean trend.
Definition gproc.h:118
double covariance(Array1D< double > &x1, Array1D< double > &x2, Array1D< double > &param)
Prior covariance function.
Definition gproc.cpp:460
Array1D< double > var_
Variance of the Student-t posterior.
Definition gproc.h:143
int getNpt() const
Get the number of data points.
Definition gproc.h:74
Array1D< double > param_
Roughness parameter vector.
Definition gproc.h:147
double getBe() const
Get beta parameter.
Definition gproc.h:82
Array1D< double > Hbhat_
Definition gproc.h:160
double al_
Prior parameter .
Definition gproc.h:125
~Gproc()
Destructor: cleans up all memory and destroys object.
Definition gproc.h:48
void EvalGP_inv(Array2D< double > &xgrid, string msc, Array1D< double > &mst)
Evaluate the Gaussian Process at a given grid msc controls whether only mean will be computed,...
Definition gproc.cpp:328
int npt_
Number of data points.
Definition gproc.h:132
Array1D< double > AinvyHbhat_
Definition gproc.h:162
Array2D< double > cov_
Covariance of the Student-t posterior.
Definition gproc.h:145
Array1D< double > z_
Prior mean of the mean trend.
Definition gproc.h:122
Array1D< double > dataVar_
Data noise 'nugget'.
Definition gproc.h:115
int getNdim() const
Get the dimensionality.
Definition gproc.h:76
PCSet * PCModel_
Basis set for the trend function.
Definition gproc.h:138
double getAl() const
Get alpha parameter.
Definition gproc.h:80
Array1D< double > mst_
Mean of the Student-t posterior.
Definition gproc.h:141
void SetupData(Array2D< double > &xdata, Array1D< double > &ydata, Array1D< double > &datavar)
Setup the data.
Definition gproc.cpp:87
Array2D< double > HtAinvH_
Definition gproc.h:157
Defines and initializes PC basis function set and provides functions to manipulate PC expansions defi...
Definition PCSet.h:73
static double x1[]
Definition gkpclib.cpp:36