UQTk: Uncertainty Quantification Toolkit 3.1.5
PCSet Class Reference

Defines and initializes PC basis function set and provides functions to manipulate PC expansions defined on this basis set. More...

#include <PCSet.h>

Public Member Functions

 PCSet (const string sp_type, const int order, const int n_dim, const string pc_type, const double alpha=0.0, const double betta=1.0)
 Constructor: initializes the PC basis set for the order, number of dimensions and type that are passed in.
 
 PCSet (const string sp_type, const int order, const int n_dim, const string pc_type, const string pc_seq, const double alpha=0.0, const double betta=1.0)
 Constructor: initializes the PC basis set for the order, number of dimensions and type that are passed in. It also customizes the multiindex sequence ( lexicographical-lex, colexicographical-colex, reverse lexicographical-revlex, and reverse clexicographical-revcolex).
 
 PCSet (const string sp_type, const Array1D< int > &maxOrders, const int n_dim, const string pc_type, const double alpha=0.0, const double betta=1.0)
 Constructor: initializes the PC basis set ordered in an HDMR fashion given order per each HDMR rank (univariate, bivariate, etc...)
 
 PCSet (const string sp_type, const Array2D< int > &customMultiIndex, const string pc_type, const double alpha=0.0, const double betta=1.0)
 Constructor: initializes the PC basis set for a given custom multiIndex.
 
 ~PCSet ()
 Destructor: cleans up all memory and destroys object.
 
void dPhi_alpha (Array1D< double > &x, Array1D< int > &alpha, Array1D< double > &grad)
 Set Gradient and Hessian operators.
 
void dPhi (Array1D< double > &x, Array2D< int > &mindex, Array1D< double > &grad, Array1D< double > &ck)
 Evaluate Gradient at a single d-dim point.
 
void dPhi (Array2D< double > &x, Array2D< int > &mindex, Array2D< double > &grad, Array1D< double > &ck)
 Evaluate Gradient at a multiple d-dim x points.
 
void ddPhi_alpha (Array1D< double > &x, Array1D< int > &alpha, Array2D< double > &grad)
 Evaluate Hessian at a single d-dim point and d-dim basis polynomial.
 
void ddPhi (Array1D< double > &x, Array2D< int > &mindex, Array2D< double > &grad, Array1D< double > &ck)
 Evaluate Gradient at a single d-dim point.
 
void SetQd1d (Array1D< double > &qdpts1d, Array1D< double > &wghts1d, int nqd)
 Set the quadrature rule.
 
void SetQuadRule (const string grid_type, const string fs_type, int param)
 Set the quadrature points by specifying a grid type, a full/sparse indicator, and an integer parameter.
 
void SetQuadRule (Quad &quadRule)
 Set a custom quadrature rule by pointing to the corresponding object.
 
void PrintMultiIndex () const
 Print information on the screen.
 
void PrintMultiIndexNormSquared () const
 For all terms, print their multi-index and norm^2 on the screen.
 
string GetPCType () const
 Get and set variables/arrays inline.
 
double GetAlpha () const
 Get the value of the parameter alpha.
 
double GetBeta () const
 Get the value of the parameter beta.
 
void GetMultiIndex (Array2D< int > &mindex) const
 Get the multiindex (return Array2D)
 
void GetMultiIndex (int *mindex) const
 Get the multiindex (return double *)
 
void GetNormSq (Array1D< double > &normsq) const
 Get the norm-squared.
 
int GetNumberPCTerms () const
 Get the number of terms in a PC expansion of this order and dimension.
 
int GetNDim () const
 Get the PC dimensionality.
 
int GetOrder () const
 Get the PC order.
 
int GetNQuadPoints () const
 Get the number of quadrature points.
 
void GetQuadPoints (Array2D< double > &quad) const
 Get the quadrature points.
 
void GetQuadPointsWeights (Array2D< double > &quad, Array1D< double > &wghts) const
 Get the quadrature points and weights.
 
void GetQuadPoints (double *quad) const
 Get the quadrature points folded into a one-dimensional array quad.
 
void GetQuadWeights (Array1D< double > &wghts) const
 Get the quadrature weights.
 
void GetQuadWeights (double *wghts) const
 Get the quadrature weights folded into a one-dimensional array wghts.
 
void GetPsi (Array2D< double > &psi) const
 Get the values of the basis polynomials evaluated at the quadrature points.
 
void GetPsi (double *psi) const
 Get the polynomials evaluated at the quadrature points folded into a one-dimensional array psi.
 
void GetPsiSq (Array1D< double > &psisq) const
 Get the basis polynomial norms-squared in an array class object psisq.
 
void GetPsiSq (double *psisq) const
 Get the basis polynomial norms-squared in a double* array psisq.
 
double GetTaylorTolerance () const
 Get relative tolerance for Taylor series approximations.
 
void SetTaylorTolerance (const double &rTol)
 Set relative tolerance for Taylor series approximations.
 
int GetTaylorTermsMax () const
 Get maximum number of terms in Taylor series approximations.
 
void SetTaylorTermsMax (const int &maxTerm)
 Set maximum number of terms in Taylor series approximations.
 
void SetLogCompMethod (const LogCompMethod &logMethod)
 Set method of computing the log function.
 
double GetGMRESDivTolerance () const
 Get relative tolerance for GMRES in Div routine.
 
void SetGMRESDivTolerance (const double &rTol)
 Set the relative tolerance for GMRES in Div routine.
 
void InitMeanStDv (const double &m, const double &s, double *p) const
 Intrusive arithmetics.
 
void InitMeanStDv (const double &m, const double &s, Array1D< double > &p) const
 Initializes a PC expansion p in Array1D<double> format to have the same distribution as the underlying PC germ, but with a specified mean m and standard deviation s.
 
void Copy (double *p1, const double *p2) const
 Copy PC expansion p2 into p1 (i.e. p1 = p2).
 
void Copy (Array1D< double > &p1, const Array1D< double > &p2) const
 Copy PC expansion p2 into p1 (i.e. p1 = p2).
 
void Add (const double *p1, const double *p2, double *p3) const
 Add two PC expansions given by double* arguments p1 and p2, and return the result in p3.
 
void Add (const Array1D< double > &p1, const Array1D< double > &p2, Array1D< double > &p3) const
 Add two PC expansions given by Array1D arguments p1 and p2, and return the result in p3.
 
void AddInPlace (double *p1, const double *p2) const
 Add PC expansions given by double* argument p2 to p1 and return the result in p1.
 
void AddInPlace (Array1D< double > &p1, const Array1D< double > &p2) const
 Add PC expansions given by Array1D argument p2 to p1 and return the result in p1.
 
void Multiply (const double *p1, const double &a, double *p2) const
 Multiply PC expansion p1 with scalar a and return the result in p2. All PCEs are in double* format.
 
void Multiply (const Array1D< double > &p1, const double &a, Array1D< double > &p2) const
 Multiply PC expansion p1 with scalar a and return the result in p2. All PCEs are in Array1D format.
 
void MultiplyInPlace (double *p1, const double &a) const
 Multiply PC expansions given by double* argument p1 with scalar a and return the result in p1.
 
void MultiplyInPlace (Array1D< double > &p1, const double &a) const
 Multiply PC expansions given by Array1D argument p1 with scalar a and return the result in p1.
 
void Subtract (const double *p1, const double *p2, double *p3) const
 Subtract PC expansion p2 from p1, and return the result in p3, with all arguments given as double*.
 
void Subtract (const Array1D< double > &p1, const Array1D< double > &p2, Array1D< double > &p3) const
 Subtract PC expansion p2 from p1, and return the result in p3, with all arguments given as Array1D structures.
 
void SubtractInPlace (double *p1, const double *p2) const
 Subtract PC expansion p2 from p1, and return the result in p1, with all arguments given as double*.
 
void SubtractInPlace (Array1D< double > &p1, const Array1D< double > &p2) const
 Subtract PC expansion p2 from p1, and return the result in p1, with all arguments given as Array1D structures.
 
void Prod (const double *p1, const double *p2, double *p3) const
 Multiply two PC expansions given by double* arguments p1 and p2, and return the result in p3.
 
void Prod (const Array1D< double > &p1, const Array1D< double > &p2, Array1D< double > &p3) const
 Multipy two PC expansions given by Array1D arguments p1 and p2, and return the result in p3.
 
void Prod3 (const double *p1, const double *p2, const double *p3, double *p4) const
 Multiply three PC expansions given by double* arguments p1, p2, and p3, and return the result in p4.
 
void Prod3 (const Array1D< double > &p1, const Array1D< double > &p2, const Array1D< double > &p3, Array1D< double > &p4) const
 Multipy three PC expansions given by Array1D arguments p1, p2, and p3, and return the result in p4.
 
void Polyn (const double *polycf, int npoly, const double *p1, double *p2) const
 Evaluates a polynomial of PC that is given in double* argument p1. Polynomial coefficients are given in double* argument polycf of size npoly. The output PC is contained in double* argument p2.
 
void Polyn (const Array1D< double > &polycf, const Array1D< double > &p1, Array1D< double > &p2) const
 Evaluates a polynomial of PC that is given by Array1D argument p1. Polynomial coefficients are given in the Array1D argument polycf. The output PC is contained in Array1D argument p2.
 
void PolynMulti (const Array1D< double > &polycf, const Array2D< int > &mindex, const Array2D< double > &p1, Array1D< double > &p2) const
 Evaluates a multivariate polynomial of a set of PC inputs given by Array2D argument p1 (each column of p1 is a PC input). Polynomial coefficients are given in Array1D argument polycf. Multiindex set for the multivariate polynomial is given in Array2D argument mindex. The output PC is contained in Array1D argument p2.
 
void Exp (const double *p1, double *p2) const
 Take the exp() of the PC expansion given by double* argument p1, and return the result in p2.
 
void Exp (const Array1D< double > &p1, Array1D< double > &p2) const
 Take the exp() of the PC expansion given by Array1D argument p1, and return the result in p2.
 
void Log (const double *p1, double *p2) const
 Take the natural logarithm log() of the PC expansion given by double* argument p1, and return the result in p2. The logarithm is evaluated either via Taylor series or via integration depending on the value of parameter logMethod_.
 
void Log (const Array1D< double > &p1, Array1D< double > &p2) const
 Take the natural logarithm, log(), of the PC expansion given by Array1D argument p1, and return the result in Array1D argument p2.
 
void Log10 (const double *p1, double *p2) const
 Take the logarithm to base 10 of the PC expansion given by double* argument p1, and return the result in p2.
 
void Log10 (const Array1D< double > &p1, Array1D< double > &p2) const
 Take the logarithm to base 10 of the PC expansion given by Array1D argument p1, and return the result in Array1D argument p2.
 
void RPow (const double *p1, double *p2, const double &a) const
 Evaluate power a (a real number) of PC expansion given by double* argument p1, and return the result in p2. The power is computed as p1^a = exp(a*log(p1)), where log(p1) is evaluated either via Taylor series or via integration depending on the value of parameter logMethod_.
 
void RPow (const Array1D< double > &p1, Array1D< double > &p2, const double &a) const
 Evaluate power a (a real number) of PC expansion given by Array1D argument p1, and return the result in Array1D argument p2.
 
void IPow (const double *p1, double *p2, const int &ia) const
 Evaluate power ia (an integer number) of PC expansion given by double* argument p1, and return the result in p2.
 
void IPow (const Array1D< double > &p1, Array1D< double > &p2, const int &ia) const
 Evaluate power ia (an integer number) of PC expansion given by Array1D argument p1, and return the result in Array1D argument p2.
 
void Inv (const double *p1, double *p2) const
 Evaluate the inverse of PC expansion given by double* argument p1, and return the result in p2. The inverse is computed using the division function.
 
void Inv (const Array1D< double > &p1, Array1D< double > &p2) const
 Evaluate the inverse of PC expansion given by Array1D argument p1, and return the result in Array1D argument p2.
 
void Div (const double *p1, const double *p2, double *p3) const
 Divide the PC expansion p1 by p2, and return the result in p3 (All arguments in double* format)
 
void Div (const Array1D< double > &p1, const Array1D< double > &p2, Array1D< double > &p3) const
 Divide the PC expansion p1 by p2, and return the result in p3 (All arguments in Array1D<double> format)
 
double StDv (const double *p) const
 Returns the standard deviation of PC expansion p in a double* format.
 
double StDv (const Array1D< double > &p) const
 Returns the standard deviation of PC expansion p (Argument in Array1D<double> format)
 
double GetModesRMS (const double *p) const
 Compute the rms average of the PC coefficients (i.e. the square root of the average of the square of the PC coefficients, not taking into account any basis functions). (Arguments in double* format)
 
double GetModesRMS (const Array1D< double > &p) const
 Compute the rms average of the PC coefficients (i.e. the square root of the average of the square of the PC coefficients, not taking into account any basis functions). (Arguments in Array1D<double> format)
 
void Derivative (const double *p1, double *p2) const
 Computes derivatives of univariate PC given by coefficients p1 returns coefficient vector of the derivative in p2.
 
void Derivative (const Array1D< double > &p1, Array1D< double > &p2) const
 Computes derivatives of univariate PC given by coefficients p1 returns coefficient vector of the derivative in p2.
 
int GetNumTripleProd () const
 Returns number of triple products.
 
void GetTripleProd (int *nTriple, int *iProd, int *jProd, double *Cijk) const
 Returns triple products indices (int*‍/double* version)
 
void GetTripleProd (Array1D< int > &nTriple, Array1D< int > &iProd, Array1D< int > &jProd, Array1D< double > &Cijk) const
 Returns triple products indices (Array version)
 
int GetNumQuadProd () const
 Returns number of quad products.
 
void GetQuadProd (int *nQuad, int *iProd, int *jProd, int *kProd, double *Cijkl) const
 Returns quad products indices (int*‍/double* version)
 
void GetQuadProd (Array1D< int > &nQuad, Array1D< int > &iProd, Array1D< int > &jProd, Array1D< int > &kProd, Array1D< double > &Cijkl) const
 Returns quad products indices (Array version)
 
void SeedBasisRandNumGen (const int &seed) const
 Random sample generator functions.
 
void DrawSampleSet (const Array1D< double > &p, Array1D< double > &samples)
 Draw a set of samples from the PC expansion p, and return the result in the array samples. All arguments are in Array1D<double> format The number of samples requested is assumed to be the size of the samples array.
 
void DrawSampleSet (const double *p, double *samples, const int &nSamples)
 Draw a set of samples from the PC expansion given in double* argument p, and return the result in double* array samples. The number of samples requested is the argument nSamples.
 
void DrawSampleVar (Array2D< double > &samples) const
 Draw a set of samples of the underlying germ random variable.
 
void DrawSampleVar (double *samples, const int &nS, const int &nD) const
 
double EvalPC (const Array1D< double > &p, Array1D< double > &randVarSamples)
 PC evaluation functionalities.
 
double EvalPC (const double *p, const double *randVarSamples)
 Evaluate the given PC expansion p, at the specified values of the random variables, randVarSamples. All arguments in const double* format.
 
void EvalPCAtCustPoints (Array1D< double > &xch, Array2D< double > &custPoints, Array1D< double > &p)
 Evaluate the given PC expansion at given set of points with given coefficient vector and return the values in an 1D Array in the first argument.
 
void EvalBasisAtCustPts (const Array2D< double > &custPoints, Array2D< double > &psi)
 Evaluate Basis Functions at given points custPoints and return in the array psi.
 
void EvalBasisAtCustPts (const double *custPoints, const int npts, double *psi)
 
void GalerkProjection (const Array1D< double > &fcn, Array1D< double > &ck)
 Galerkin projection functionalities.
 
void GalerkProjectionMC (const Array2D< double > &x, const Array1D< double > &fcn, Array1D< double > &ck)
 Galerkin Projection via Monte-Carlo integration.
 
int ComputeOrders (Array1D< int > &orders)
 Multiindex parsing functionalities.
 
int ComputeEffDims (int *effdim)
 Computes the effective dimensionality of each basis term, i.e., the number of dimensions that enter with a non-zero degree. also returns the maximal dimensionality among all basis terms.
 
int ComputeEffDims (Array1D< int > &effdim)
 Computes the effective dimensionality of each basis term, i.e., the number of dimensions that enter with a non-zero degree. also returns the maximal dimensionality among all basis terms.
 
void EncodeMindex (Array1D< Array2D< int > > &sp_mindex)
 Encode multiIndex into a 'sparse' format where the bases are ordered by their effective dimensionality. The i-th element in sp_mindex stores all the bases that have effective dimensionality equal to i. Also, only non-zero components are stored.
 
double ComputeMean (const double *coef)
 Moment/sensitivity extraction given coefficients.
 
double ComputeMean (Array1D< double > &coef)
 Compute the mean of the PC given coefficient array coef(seeking the zero-th order multiindex)
 
double ComputeVarFrac (const double *coef, double *varfrac)
 Compute the variance fractions of each basis term given coefficients in double *coef; returns the variance fractions in the double *varfrac.
 
double ComputeVarFrac (Array1D< double > &coef, Array1D< double > &varfrac)
 Compute the variance fractions of each basis term given coefficient array coef; returns the variance fractions in the array varfrac.
 
void ComputeMainSens (Array1D< double > &coef, Array1D< double > &mainsens)
 Compute main effect sensitivity (Sobol) indices given coefficient array coef; returns the indices in the array mainsens.
 
void ComputeTotSens (Array1D< double > &coef, Array1D< double > &totsens)
 Compute total effect sensitivity (Sobol) indices given coefficient array coeff; returns the indices in the array totsens.
 
void ComputeJointSens (Array1D< double > &coef, Array2D< double > &jointsens)
 Compute joint effect sensitivity (Sobol) indices given coefficient array coeff; returns the indices in the array jointsens.
 
void SetVerbosity (int verbosity)
 Other.
 
void EvalNormSq (Array1D< double > &normsq)
 Evaluate norms-squared of all bases and return in the array normsq.
 
void EvalNormSq (double *normsq, const int npc)
 
void EvalNormSqExact (Array1D< double > &normsq)
 Evaluate norms-squared analytically of all bases and return in the array normsq.
 
bool IsInDomain (double x)
 Check if the point x is in the PC domain.
 

Private Types

typedef std::map< int, PCSet * > OMap_t
 Definition of a map to connect integer indexes with pointers to this class.
 

Private Member Functions

 PCSet ()
 Dummy default constructor, which should not be used as it is not well defined Therefore we make it private so it is not accessible.
 
 PCSet (const PCSet &obj)
 Dummy copy constructor, which should not be used as it is currently not well defined. Therefore we make it private so it is not accessible.
 
void ComputeMaxOrdPerDim ()
 Compute maximal order per dimension and fill in the array maxOrdPerDim_.
 
void Initialize (const string ordertype)
 Initialization of the appropriate variables.
 
void InitISP ()
 Initialize quadrature for computing triple products(ISP) and orthogonal projection(NISP)
 
void InitNISP ()
 Initialize variables that are needed only in non-intrusive computations.
 
void EvalBasisProd3 ()
 Evaluate the expectation of product of three basis functions.
 
void EvalBasisProd4 ()
 Evaluate the expectation of product of four basis functions.
 
void GMRESMatrixVectorProd (const double *x, const double *a, double *y) const
 Actual C++ implementation of the matric vector multiplication for GMRES for the division operation.
 
void LogTaylor (const double *p1, double *p2) const
 Computes natural logarithm using Taylor expansion: N p2 = ln(p1) = ln(p1Mean) + sum d n=1 n.
 
void LogInt (const double *p1, double *p2) const
 Computes natural logarithm by numerical integration: calculate p2=ln(p1) by integrating du=dx/x to get ln(x)
 
int LogIntRhs (sunrealtype t, N_Vector y, N_Vector ydot, void *f_data) const
 Evaluates rhs necessary to compute natural logarithm via integration.
 
int Check_CVflag (void *flagvalue, const char *funcname, int opt) const
 Check cvode return for errors.
 

Static Private Member Functions

static void GMRESMatrixVectorProdWrapper (int *n, double *x, double *y, int *nelt, int *ia, int *ja, double *a, int *obj)
 Wrapper for Matrix-vector multiplication routine to be called by GMRES.
 
static void GMRESPreCondWrapper (int *n, double *r, double *z, int *nelt, int *ia, int *ja, double *a, int *obj, double *rwork, int *iwork)
 Wrapper for preconditioner routine to be called by GMRES.
 
static int LogIntRhsWrapper (sunrealtype t, N_Vector y, N_Vector ydot, void *f_data)
 Wrapper for LogIntRhs. The first component of f_data pointer carries an integer handle identifying the appropriate PC object.
 

Private Attributes

int uqtkverbose_
 Verbosity level.
 
string spType_
 String indicator of ISP or NISP implementation type.
 
string pcType_
 String indicator of PC type.
 
string pcSeq_
 String indicator of multiindex ordering.
 
PCBasisp_basis_
 Pointer to the class that defines the basis type and functions.
 
int order_
 Order of the PC representation.
 
int maxorddim_
 Maximal order within all dimensions.
 
Array1D< int > maxOrders_
 Array of maximum orders requested if custom(HDMR) ordering is requested.
 
Array1D< int > maxOrdPerDim_
 Array of maximum orders per dimension.
 
const int nDim_
 Number of stochastic dimensions (degrees of freedom) in the PC representation.
 
int nQuadPoints_
 Number of quadrature points used.
 
int nPCTerms_
 Total number of terms in the PC expansions.
 
double rTolTaylor_
 Relative tolerance for Taylor series approximations.
 
int maxTermTaylor_
 Max number of terms in Taylor series approximations.
 
double SMALL_
 Tolerance to avoid floating-point errors.
 
double rTolGMRESDiv_
 GMRES tolerance in Div()
 
Array2D< double > psi_
 Array to store basis functions evaluated at quadrature points for each order: psi_(iqp,ipc) contains the value of the polynomial chaos ipc-th basis at the location of quadrature point iqp.
 
Array1D< double > psiSq_
 Array with the norms squared of the basis functions, corresponding to each term in the PC expansion.
 
Array2D< double > quadPoints_
 Array to store quadrature points.
 
Array1D< double > quadWeights_
 Array to store quadrature weights.
 
Array2D< int > quadIndices_
 Array to store quadrature point indexing; useful for nested rules.
 
Array2D< int > multiIndex_
 Array to store multi-index: multiIndex_(ipc,idim) contains the order of the basis function associated with dimension idim, for the ipc-th term in the PC expansion.
 
Array1D< Array1D< int > > iProd2_
 i-indices of <\Psi_i \Psi_j \Psi_k> terms that are not zero, for all k
 
Array1D< Array1D< int > > jProd2_
 j-indices of <\Psi_i \Psi_j \Psi_k> terms that are not zero, for all k
 
Array1D< Array1D< double > > psiIJKProd2_
 <\Psi_i \Psi_j \Psi_k> terms that are not zero, for all k
 
Array1D< Array1D< int > > iProd3_
 i-indices of <\Psi_i \Psi_j \Psi_k \Psi_l> terms that are not zero, for all l
 
Array1D< Array1D< int > > jProd3_
 j-indices of <\Psi_i \Psi_j \Psi_k \Psi_l> terms that are not zero, for all l
 
Array1D< Array1D< int > > kProd3_
 k-indices of <\Psi_i \Psi_j \Psi_k \Psi_l> terms that are not zero, for all l
 
Array1D< Array1D< double > > psiIJKLProd3_
 <\Psi_i \Psi_j \Psi_k \Psi_l> terms that are not zero, for all l
 
LogCompMethod logMethod_
 Flag for method to compute log: TaylorSeries or Integration.
 
int CVmaxord_
 CVODE parameter: maximal order.
 
int CVmaxnumsteps_
 CVODE parameter: maximal number of steps.
 
double CVinitstep_
 CVODE parameter: initial step size.
 
double CVmaxstep_
 CVODE parameter: maximal step size.
 
double CVrelt_
 CVODE parameter: relative tolerance.
 
double CVabst_
 CVODE parameter: absolute tolerance.
 
int my_index_
 Index of this class.
 
int narg_
 Number of free parameters to specify the basis.
 
double alpha_
 Parameter alpha for PCs that require a parameter (LG,SW,JB)
 
double beta_
 Parameter beta for PCs that require two parameters (SW,JB)
 

Static Private Attributes

static int next_index_ = 0
 index of next object in map
 
static OMap_tomap_ = NULL
 Map to connect integer indexes with pointers to this class.
 

Detailed Description

Defines and initializes PC basis function set and provides functions to manipulate PC expansions defined on this basis set.

Member Typedef Documentation

◆ OMap_t

std::map<int, PCSet*> PCSet::OMap_t
private

Definition of a map to connect integer indexes with pointers to this class.

Constructor & Destructor Documentation

◆ PCSet() [1/6]

PCSet::PCSet ( const string sp_type,
const int order,
const int n_dim,
const string pc_type,
const double alpha = 0.0,
const double betta = 1.0 )

Constructor: initializes the PC basis set for the order, number of dimensions and type that are passed in.

Implementation type sp_type has three options "ISP" (intrusive methods), "NISP" (non-intrusive), or "NISPnoq" (non-intrusive without quadrature initialization)

Note
alpha and betta are parameters only relevant for LG, JB or SW chaoses

◆ PCSet() [2/6]

PCSet::PCSet ( const string sp_type,
const int order,
const int n_dim,
const string pc_type,
const string pc_seq,
const double alpha = 0.0,
const double betta = 1.0 )

Constructor: initializes the PC basis set for the order, number of dimensions and type that are passed in. It also customizes the multiindex sequence ( lexicographical-lex, colexicographical-colex, reverse lexicographical-revlex, and reverse clexicographical-revcolex).

Implementation type sp_type has three options "ISP" (intrusive methods), "NISP" (non-intrusive), or "NISPnoq" (non-intrusive without quadrature initialization)

Note
alpha and betta are parameters only relevant for LG, JB or SW chaoses

◆ PCSet() [3/6]

PCSet::PCSet ( const string sp_type,
const Array1D< int > & maxOrders,
const int n_dim,
const string pc_type,
const double alpha = 0.0,
const double betta = 1.0 )

Constructor: initializes the PC basis set ordered in an HDMR fashion given order per each HDMR rank (univariate, bivariate, etc...)

Implementation type sp_type has three options "ISP" (intrusive methods), "NISP" (non-intrusive), or "NISPnoq" (non-intrusive without quadrature initialization)

Note
alpha and betta are parameters only relevant for LG, JB or SW chaoses

◆ PCSet() [4/6]

PCSet::PCSet ( const string sp_type,
const Array2D< int > & customMultiIndex,
const string pc_type,
const double alpha = 0.0,
const double betta = 1.0 )

Constructor: initializes the PC basis set for a given custom multiIndex.

Implementation type sp_type has three options "ISP" (intrusive methods), "NISP" (non-intrusive), or "NISPnoq" (non-intrusive without quadrature initialization)

Note
alpha and betta are parameters only relevant for LG, JB or SW chaoses

◆ ~PCSet()

PCSet::~PCSet ( )

Destructor: cleans up all memory and destroys object.

◆ PCSet() [5/6]

PCSet::PCSet ( )
inlineprivate

Dummy default constructor, which should not be used as it is not well defined Therefore we make it private so it is not accessible.

Note
All parameters are intialized to dummy values.

◆ PCSet() [6/6]

PCSet::PCSet ( const PCSet & obj)
inlineprivate

Dummy copy constructor, which should not be used as it is currently not well defined. Therefore we make it private so it is not accessible.

Note
I am not sure actually whether the initialization performed below is legal as it requires access to private data members of the class that is passed in.

Member Function Documentation

◆ Add() [1/2]

void PCSet::Add ( const Array1D< double > & p1,
const Array1D< double > & p2,
Array1D< double > & p3 ) const

Add two PC expansions given by Array1D arguments p1 and p2, and return the result in p3.

Note
Requires the size of the arrays that are passed in to equal the number of PC terms

◆ Add() [2/2]

void PCSet::Add ( const double * p1,
const double * p2,
double * p3 ) const

Add two PC expansions given by double* arguments p1 and p2, and return the result in p3.

◆ AddInPlace() [1/2]

void PCSet::AddInPlace ( Array1D< double > & p1,
const Array1D< double > & p2 ) const

Add PC expansions given by Array1D argument p2 to p1 and return the result in p1.

Note
Requires the size of the arrays that are passed in to equal the number of PC terms

◆ AddInPlace() [2/2]

void PCSet::AddInPlace ( double * p1,
const double * p2 ) const

Add PC expansions given by double* argument p2 to p1 and return the result in p1.

◆ Check_CVflag()

int PCSet::Check_CVflag ( void * flagvalue,
const char * funcname,
int opt ) const
private

Check cvode return for errors.

◆ ComputeEffDims() [1/2]

int PCSet::ComputeEffDims ( Array1D< int > & effdim)

Computes the effective dimensionality of each basis term, i.e., the number of dimensions that enter with a non-zero degree. also returns the maximal dimensionality among all basis terms.

Note
This is not the classical effective dimensionality, since all dimensions can still be involved.

◆ ComputeEffDims() [2/2]

int PCSet::ComputeEffDims ( int * effdim)

Computes the effective dimensionality of each basis term, i.e., the number of dimensions that enter with a non-zero degree. also returns the maximal dimensionality among all basis terms.

Note
This is not the classical effective dimensionality, since all dimensions can still be involved.

◆ ComputeJointSens()

void PCSet::ComputeJointSens ( Array1D< double > & coef,
Array2D< double > & jointsens )

Compute joint effect sensitivity (Sobol) indices given coefficient array coeff; returns the indices in the array jointsens.

Note
jointsens will be populated as a strictly upper-diagonal matrix
Todo
There is no double* version of this function

◆ ComputeMainSens()

void PCSet::ComputeMainSens ( Array1D< double > & coef,
Array1D< double > & mainsens )

Compute main effect sensitivity (Sobol) indices given coefficient array coef; returns the indices in the array mainsens.

Todo
There is no double* version of this function

◆ ComputeMaxOrdPerDim()

void PCSet::ComputeMaxOrdPerDim ( )
private

Compute maximal order per dimension and fill in the array maxOrdPerDim_.

◆ ComputeMean() [1/2]

double PCSet::ComputeMean ( Array1D< double > & coef)

Compute the mean of the PC given coefficient array coef(seeking the zero-th order multiindex)

◆ ComputeMean() [2/2]

double PCSet::ComputeMean ( const double * coef)

Moment/sensitivity extraction given coefficients.

Compute the mean of the PC given coefficients in double *coef (seeking the zero-th order multiindex)

◆ ComputeOrders()

int PCSet::ComputeOrders ( Array1D< int > & orders)

Multiindex parsing functionalities.

Computes the order of each basis term and return it in the array orders, also returns the maximal order

Todo
There is no double* version of this function

◆ ComputeTotSens()

void PCSet::ComputeTotSens ( Array1D< double > & coef,
Array1D< double > & totsens )

Compute total effect sensitivity (Sobol) indices given coefficient array coeff; returns the indices in the array totsens.

Todo
There is no double* version of this function

◆ ComputeVarFrac() [1/2]

double PCSet::ComputeVarFrac ( Array1D< double > & coef,
Array1D< double > & varfrac )

Compute the variance fractions of each basis term given coefficient array coef; returns the variance fractions in the array varfrac.

Note
Also returns the variance
The value for the zeroth order term has a special meaning: it is equal to mean^2/variance or (mean/std)^2.

◆ ComputeVarFrac() [2/2]

double PCSet::ComputeVarFrac ( const double * coef,
double * varfrac )

Compute the variance fractions of each basis term given coefficients in double *coef; returns the variance fractions in the double *varfrac.

Note
Also returns the variance
The value for the zeroth order term has a special meaning: it is equal to mean^2/variance or (mean/std)^2.

◆ Copy() [1/2]

void PCSet::Copy ( Array1D< double > & p1,
const Array1D< double > & p2 ) const

Copy PC expansion p2 into p1 (i.e. p1 = p2).

All arguments in Array format

Note
Requires the size of the arrays that are passed in to equal the number of PC terms

◆ Copy() [2/2]

void PCSet::Copy ( double * p1,
const double * p2 ) const

Copy PC expansion p2 into p1 (i.e. p1 = p2).

All arguments in double* format.

◆ ddPhi()

void PCSet::ddPhi ( Array1D< double > & x,
Array2D< int > & mindex,
Array2D< double > & grad,
Array1D< double > & ck )

Evaluate Gradient at a single d-dim point.

for a PCSet object

◆ ddPhi_alpha()

void PCSet::ddPhi_alpha ( Array1D< double > & x,
Array1D< int > & alpha,
Array2D< double > & grad )

Evaluate Hessian at a single d-dim point and d-dim basis polynomial.

◆ Derivative() [1/2]

void PCSet::Derivative ( const Array1D< double > & p1,
Array1D< double > & p2 ) const

Computes derivatives of univariate PC given by coefficients p1 returns coefficient vector of the derivative in p2.

Note
Makes use of intrusive computations on recursive formulae for derivatives
Todo

Supports LU and HG bases only

Supports only for 1d PCs

◆ Derivative() [2/2]

void PCSet::Derivative ( const double * p1,
double * p2 ) const

Computes derivatives of univariate PC given by coefficients p1 returns coefficient vector of the derivative in p2.

Note
Makes use of intrusive computations on recursive formulae for derivatives
Todo

Supports LU and HG bases only

Supports only for 1d PCs

◆ Div() [1/2]

void PCSet::Div ( const Array1D< double > & p1,
const Array1D< double > & p2,
Array1D< double > & p3 ) const

Divide the PC expansion p1 by p2, and return the result in p3 (All arguments in Array1D<double> format)

Note
Requires the size of the arrays that are passed in to equal the number of PC terms

◆ Div() [2/2]

void PCSet::Div ( const double * p1,
const double * p2,
double * p3 ) const

Divide the PC expansion p1 by p2, and return the result in p3 (All arguments in double* format)

The "division" p3 = p1/p2 is performed by solving the system of equations p2*p3 = p1 for the unknown p3.

Note
When GMRES is used to solve this system of equations (based on a preprocessor flag in the source code for this routine), a relative tolerance criterium is used that is set by default to 1.e-8, and can be changed with SetGMRESDivTolerance().
Todo
Remove duplication of data and parameters that was required for enforcing imposed "const" constraints on some of the arguments and the class data members when they are being passed to fortran.

◆ dPhi() [1/2]

void PCSet::dPhi ( Array1D< double > & x,
Array2D< int > & mindex,
Array1D< double > & grad,
Array1D< double > & ck )

Evaluate Gradient at a single d-dim point.

for a PCSet object

◆ dPhi() [2/2]

void PCSet::dPhi ( Array2D< double > & x,
Array2D< int > & mindex,
Array2D< double > & grad,
Array1D< double > & ck )

Evaluate Gradient at a multiple d-dim x points.

for a PCSet object

◆ dPhi_alpha()

void PCSet::dPhi_alpha ( Array1D< double > & x,
Array1D< int > & alpha,
Array1D< double > & grad )

Set Gradient and Hessian operators.

Evaluate Gradient at a single d-dim point and d-dim basis polynomial

◆ DrawSampleSet() [1/2]

void PCSet::DrawSampleSet ( const Array1D< double > & p,
Array1D< double > & samples )

Draw a set of samples from the PC expansion p, and return the result in the array samples. All arguments are in Array1D<double> format The number of samples requested is assumed to be the size of the samples array.

Note
The size of the array p that is passed in needs to equal the number of PC terms

◆ DrawSampleSet() [2/2]

void PCSet::DrawSampleSet ( const double * p,
double * samples,
const int & nSamples )

Draw a set of samples from the PC expansion given in double* argument p, and return the result in double* array samples. The number of samples requested is the argument nSamples.

◆ DrawSampleVar() [1/2]

void PCSet::DrawSampleVar ( Array2D< double > & samples) const

Draw a set of samples of the underlying germ random variable.

Todo
There is no double* version of this function

◆ DrawSampleVar() [2/2]

void PCSet::DrawSampleVar ( double * samples,
const int & nS,
const int & nD ) const

◆ EncodeMindex()

void PCSet::EncodeMindex ( Array1D< Array2D< int > > & sp_mindex)

Encode multiIndex into a 'sparse' format where the bases are ordered by their effective dimensionality. The i-th element in sp_mindex stores all the bases that have effective dimensionality equal to i. Also, only non-zero components are stored.

Todo
There is no double* version of this function

◆ EvalBasisAtCustPts() [1/2]

void PCSet::EvalBasisAtCustPts ( const Array2D< double > & custPoints,
Array2D< double > & psi )

Evaluate Basis Functions at given points custPoints and return in the array psi.

Todo
There is no double* version of this function

◆ EvalBasisAtCustPts() [2/2]

void PCSet::EvalBasisAtCustPts ( const double * custPoints,
const int npts,
double * psi )

◆ EvalBasisProd3()

void PCSet::EvalBasisProd3 ( )
private

Evaluate the expectation of product of three basis functions.

◆ EvalBasisProd4()

void PCSet::EvalBasisProd4 ( )
private

Evaluate the expectation of product of four basis functions.

◆ EvalNormSq() [1/2]

void PCSet::EvalNormSq ( Array1D< double > & normsq)

Evaluate norms-squared of all bases and return in the array normsq.

Todo
There is no double* version of this function

◆ EvalNormSq() [2/2]

void PCSet::EvalNormSq ( double * normsq,
const int npc )

◆ EvalNormSqExact()

void PCSet::EvalNormSqExact ( Array1D< double > & normsq)

Evaluate norms-squared analytically of all bases and return in the array normsq.

Todo
There is no double* version of this function
Note
Custom PCs do not have this capability

◆ EvalPC() [1/2]

double PCSet::EvalPC ( const Array1D< double > & p,
Array1D< double > & randVarSamples )

PC evaluation functionalities.

Evaluate the given PC expansion p, at the specified values of the random variables, randVarSamples. All arguments in const Array1D<double> format

Note
The number of elements in p needs to match the number of terms in the PC expansions in this PCSet.
The number of elements in randVarSamples needs to match the number of dimensions in the PC expansion.

◆ EvalPC() [2/2]

double PCSet::EvalPC ( const double * p,
const double * randVarSamples )

Evaluate the given PC expansion p, at the specified values of the random variables, randVarSamples. All arguments in const double* format.

Note
The number of elements in p is assumed to match the number of terms in the PC expansions in this PCSet.
The number of elements in randVarSamples is assumed to match the number of dimensions in the PC expansion.

◆ EvalPCAtCustPoints()

void PCSet::EvalPCAtCustPoints ( Array1D< double > & xch,
Array2D< double > & custPoints,
Array1D< double > & p )

Evaluate the given PC expansion at given set of points with given coefficient vector and return the values in an 1D Array in the first argument.

Todo
There is no double* version of this function

◆ Exp() [1/2]

void PCSet::Exp ( const Array1D< double > & p1,
Array1D< double > & p2 ) const

Take the exp() of the PC expansion given by Array1D argument p1, and return the result in p2.

Note
Requires the size of the arrays that are passed in to equal the number of PC terms

◆ Exp() [2/2]

void PCSet::Exp ( const double * p1,
double * p2 ) const

Take the exp() of the PC expansion given by double* argument p1, and return the result in p2.

Relies on Taylor series expansion: exp(x) = 1 + x + x^2/2! + x^3/3! + ... However, for efficiency and to avoid overflow, the terms are computed as d_i = d_{i-1}*x/i. Also, to reduce the number of terms needed in the series, we subtract the mean out of a random variable u as u = u_0 + (u-u_0) and exp(u) = exp(u_0)*exp(u-u_0), where exp(u_0) can be computed with the regular exp(double& ) function

Note
The Taylor series is truncated after a tolerance criterium is achieved on the relative error defined as the max absolute value of the PC coefficients in the last added term, divided by the mean of exp(p1). The tolerance is set to 1.e-6 by default and can be changed with SetTaylorTolerance().
The maximum number of terms in the Taylor series is set by default to 500 and can be changed with SetTaylorTermsMax()

◆ GalerkProjection()

void PCSet::GalerkProjection ( const Array1D< double > & fcn,
Array1D< double > & ck )

Galerkin projection functionalities.

Performs (NISP) Galerkin projection, given function evaluations at quadrature points Returns in the coefficient vector in the second argument

Note
User should make sure that the function HAS BEEN evaluated at the correct quadrature points by first extracting the quadrature points and evaluating the function externally
Todo

Overload this with forward function pointers

There is no double* version of this function

◆ GalerkProjectionMC()

void PCSet::GalerkProjectionMC ( const Array2D< double > & x,
const Array1D< double > & fcn,
Array1D< double > & ck )

Galerkin Projection via Monte-Carlo integration.

Note
User should make sure that the function HAS BEEN evaluated at the correct sampling points by first sampling the proper PC germ distribution and evaluating the function externally
Todo

Overload this with forward function pointers

There is no double* version of this function

◆ GetAlpha()

double PCSet::GetAlpha ( ) const
inline

Get the value of the parameter alpha.

◆ GetBeta()

double PCSet::GetBeta ( ) const
inline

Get the value of the parameter beta.

◆ GetGMRESDivTolerance()

double PCSet::GetGMRESDivTolerance ( ) const
inline

Get relative tolerance for GMRES in Div routine.

◆ GetModesRMS() [1/2]

double PCSet::GetModesRMS ( const Array1D< double > & p) const

Compute the rms average of the PC coefficients (i.e. the square root of the average of the square of the PC coefficients, not taking into account any basis functions). (Arguments in Array1D<double> format)

Note
Requires the size of the array that is passed in to equal the number of PC terms

◆ GetModesRMS() [2/2]

double PCSet::GetModesRMS ( const double * p) const

Compute the rms average of the PC coefficients (i.e. the square root of the average of the square of the PC coefficients, not taking into account any basis functions). (Arguments in double* format)

◆ GetMultiIndex() [1/2]

void PCSet::GetMultiIndex ( Array2D< int > & mindex) const
inline

Get the multiindex (return Array2D)

◆ GetMultiIndex() [2/2]

void PCSet::GetMultiIndex ( int * mindex) const

Get the multiindex (return double *)

◆ GetNDim()

int PCSet::GetNDim ( ) const
inline

Get the PC dimensionality.

◆ GetNormSq()

void PCSet::GetNormSq ( Array1D< double > & normsq) const
inline

Get the norm-squared.

Todo
this seems like a duplication, see below GetPsiSq()

◆ GetNQuadPoints()

int PCSet::GetNQuadPoints ( ) const
inline

Get the number of quadrature points.

◆ GetNumberPCTerms()

int PCSet::GetNumberPCTerms ( ) const
inline

Get the number of terms in a PC expansion of this order and dimension.

◆ GetNumQuadProd()

int PCSet::GetNumQuadProd ( ) const

Returns number of quad products.

◆ GetNumTripleProd()

int PCSet::GetNumTripleProd ( ) const

Returns number of triple products.

◆ GetOrder()

int PCSet::GetOrder ( ) const
inline

Get the PC order.

◆ GetPCType()

string PCSet::GetPCType ( ) const
inline

Get and set variables/arrays inline.

Get the PC type

◆ GetPsi() [1/2]

void PCSet::GetPsi ( Array2D< double > & psi) const
inline

Get the values of the basis polynomials evaluated at the quadrature points.

◆ GetPsi() [2/2]

void PCSet::GetPsi ( double * psi) const
inline

Get the polynomials evaluated at the quadrature points folded into a one-dimensional array psi.

◆ GetPsiSq() [1/2]

void PCSet::GetPsiSq ( Array1D< double > & psisq) const
inline

Get the basis polynomial norms-squared in an array class object psisq.

◆ GetPsiSq() [2/2]

void PCSet::GetPsiSq ( double * psisq) const
inline

Get the basis polynomial norms-squared in a double* array psisq.

◆ GetQuadPoints() [1/2]

void PCSet::GetQuadPoints ( Array2D< double > & quad) const
inline

Get the quadrature points.

◆ GetQuadPoints() [2/2]

void PCSet::GetQuadPoints ( double * quad) const
inline

Get the quadrature points folded into a one-dimensional array quad.

◆ GetQuadPointsWeights()

void PCSet::GetQuadPointsWeights ( Array2D< double > & quad,
Array1D< double > & wghts ) const
inline

Get the quadrature points and weights.

◆ GetQuadProd() [1/2]

void PCSet::GetQuadProd ( Array1D< int > & nQuad,
Array1D< int > & iProd,
Array1D< int > & jProd,
Array1D< int > & kProd,
Array1D< double > & Cijkl ) const

Returns quad products indices (Array version)

◆ GetQuadProd() [2/2]

void PCSet::GetQuadProd ( int * nQuad,
int * iProd,
int * jProd,
int * kProd,
double * Cijkl ) const

Returns quad products indices (int*‍/double* version)

◆ GetQuadWeights() [1/2]

void PCSet::GetQuadWeights ( Array1D< double > & wghts) const
inline

Get the quadrature weights.

◆ GetQuadWeights() [2/2]

void PCSet::GetQuadWeights ( double * wghts) const
inline

Get the quadrature weights folded into a one-dimensional array wghts.

◆ GetTaylorTermsMax()

int PCSet::GetTaylorTermsMax ( ) const
inline

Get maximum number of terms in Taylor series approximations.

◆ GetTaylorTolerance()

double PCSet::GetTaylorTolerance ( ) const
inline

Get relative tolerance for Taylor series approximations.

◆ GetTripleProd() [1/2]

void PCSet::GetTripleProd ( Array1D< int > & nTriple,
Array1D< int > & iProd,
Array1D< int > & jProd,
Array1D< double > & Cijk ) const

Returns triple products indices (Array version)

◆ GetTripleProd() [2/2]

void PCSet::GetTripleProd ( int * nTriple,
int * iProd,
int * jProd,
double * Cijk ) const

Returns triple products indices (int*‍/double* version)

◆ GMRESMatrixVectorProd()

void PCSet::GMRESMatrixVectorProd ( const double * x,
const double * a,
double * y ) const
private

Actual C++ implementation of the matric vector multiplication for GMRES for the division operation.

Given the structure of the problem, this boils down to the product between two PC variables.

◆ GMRESMatrixVectorProdWrapper()

static void PCSet::GMRESMatrixVectorProdWrapper ( int * n,
double * x,
double * y,
int * nelt,
int * ia,
int * ja,
double * a,
int * obj )
inlinestaticprivate

Wrapper for Matrix-vector multiplication routine to be called by GMRES.

As GMRES is a Fortran77 routine, this routine is defined as a static function. One of the function arguments (obj) was originally isym, a flag for matrix symmetry, but has been repurposed to carry an integer handle to identify this object.

Note
The matrix vector product here comes down to a product between two PC expansions.

◆ GMRESPreCondWrapper()

static void PCSet::GMRESPreCondWrapper ( int * n,
double * r,
double * z,
int * nelt,
int * ia,
int * ja,
double * a,
int * obj,
double * rwork,
int * iwork )
inlinestaticprivate

Wrapper for preconditioner routine to be called by GMRES.

As GMRES is a Fortran77 routine, this routine is defined as a static function. One of the function arguments (obj) was originally isym, a flag for matrix symmetry, but has been repurposed to carry an integer handle to identify this object.

Note
Since we currently do not use preconditioning, this routine does nothing. It is a place holder for future use.

◆ Initialize()

void PCSet::Initialize ( const string ordertype)
private

Initialization of the appropriate variables.

Note
Intrusive implementation only works with TotalOrder multiindes
Todo
Test and allow intrusive implementation with customized multiindices

◆ InitISP()

void PCSet::InitISP ( )
private

Initialize quadrature for computing triple products(ISP) and orthogonal projection(NISP)

Initialize variables that are needed only in intrusive computations

◆ InitMeanStDv() [1/2]

void PCSet::InitMeanStDv ( const double & m,
const double & s,
Array1D< double > & p ) const

Initializes a PC expansion p in Array1D<double> format to have the same distribution as the underlying PC germ, but with a specified mean m and standard deviation s.

Note
This assumes that the zeroth order term is the first one in the multi-index - this assumption does not hold in general
This function only holds for expansions with one stochastic dimension
All existing coefficient values in p will be overwritten
Todo
Make this function work for general multi-indices, and for any number of stochastic dimensions

◆ InitMeanStDv() [2/2]

void PCSet::InitMeanStDv ( const double & m,
const double & s,
double * p ) const

Intrusive arithmetics.

Initializes a PC expansion p in a double* format to have the same distribution as the underlying PC germ, but with a specified mean m and standard deviation s

Note
This assumes that the zeroth order term is the first one in the multi-index - this assumption does not hold in general
This function only holds for expansions with one stochastic dimension
All existing coefficient values in p will be overwritten
Todo
Make this function work for general multi-indices, and for any number of stochastic dimensions

◆ InitNISP()

void PCSet::InitNISP ( )
private

Initialize variables that are needed only in non-intrusive computations.

◆ Inv() [1/2]

void PCSet::Inv ( const Array1D< double > & p1,
Array1D< double > & p2 ) const

Evaluate the inverse of PC expansion given by Array1D argument p1, and return the result in Array1D argument p2.

Note
Requires the size of the arrays that are passed in to equal the number of PC terms

◆ Inv() [2/2]

void PCSet::Inv ( const double * p1,
double * p2 ) const

Evaluate the inverse of PC expansion given by double* argument p1, and return the result in p2. The inverse is computed using the division function.

◆ IPow() [1/2]

void PCSet::IPow ( const Array1D< double > & p1,
Array1D< double > & p2,
const int & ia ) const

Evaluate power ia (an integer number) of PC expansion given by Array1D argument p1, and return the result in Array1D argument p2.

Note
Requires the size of the arrays that are passed in to equal the number of PC terms

◆ IPow() [2/2]

void PCSet::IPow ( const double * p1,
double * p2,
const int & ia ) const

Evaluate power ia (an integer number) of PC expansion given by double* argument p1, and return the result in p2.

◆ IsInDomain()

bool PCSet::IsInDomain ( double x)

Check if the point x is in the PC domain.

◆ Log() [1/2]

void PCSet::Log ( const Array1D< double > & p1,
Array1D< double > & p2 ) const

Take the natural logarithm, log(), of the PC expansion given by Array1D argument p1, and return the result in Array1D argument p2.

Note
Requires the size of the arrays that are passed in to equal the number of PC terms

◆ Log() [2/2]

void PCSet::Log ( const double * p1,
double * p2 ) const

Take the natural logarithm log() of the PC expansion given by double* argument p1, and return the result in p2. The logarithm is evaluated either via Taylor series or via integration depending on the value of parameter logMethod_.

◆ Log10() [1/2]

void PCSet::Log10 ( const Array1D< double > & p1,
Array1D< double > & p2 ) const

Take the logarithm to base 10 of the PC expansion given by Array1D argument p1, and return the result in Array1D argument p2.

Note
Requires the size of the arrays that are passed in to equal the number of PC terms

◆ Log10() [2/2]

void PCSet::Log10 ( const double * p1,
double * p2 ) const

Take the logarithm to base 10 of the PC expansion given by double* argument p1, and return the result in p2.

First use Log() to compute the natural logarithm and then divide it by log(10)

◆ LogInt()

void PCSet::LogInt ( const double * p1,
double * p2 ) const
private

Computes natural logarithm by numerical integration: calculate p2=ln(p1) by integrating du=dx/x to get ln(x)

◆ LogIntRhs()

int PCSet::LogIntRhs ( sunrealtype t,
N_Vector y,
N_Vector ydot,
void * f_data ) const
private

Evaluates rhs necessary to compute natural logarithm via integration.

◆ LogIntRhsWrapper()

static int PCSet::LogIntRhsWrapper ( sunrealtype t,
N_Vector y,
N_Vector ydot,
void * f_data )
inlinestaticprivate

Wrapper for LogIntRhs. The first component of f_data pointer carries an integer handle identifying the appropriate PC object.

Todo
Why is this function a static int instead of static void? Should there be a return statement at the end?

◆ LogTaylor()

void PCSet::LogTaylor ( const double * p1,
double * p2 ) const
private

Computes natural logarithm using Taylor expansion: N p2 = ln(p1) = ln(p1Mean) + sum d n=1 n.

(n+1) (-1) n p1 where d = -— *x , and x = ---— - 1 n n p1Mean

Note
See Exp notes for info related to tolerance and maximum number of terms criteria for truncating the Taylor series

◆ Multiply() [1/2]

void PCSet::Multiply ( const Array1D< double > & p1,
const double & a,
Array1D< double > & p2 ) const

Multiply PC expansion p1 with scalar a and return the result in p2. All PCEs are in Array1D format.

Note
Requires the size of the arrays that are passed in to equal the number of PC terms

◆ Multiply() [2/2]

void PCSet::Multiply ( const double * p1,
const double & a,
double * p2 ) const

Multiply PC expansion p1 with scalar a and return the result in p2. All PCEs are in double* format.

◆ MultiplyInPlace() [1/2]

void PCSet::MultiplyInPlace ( Array1D< double > & p1,
const double & a ) const

Multiply PC expansions given by Array1D argument p1 with scalar a and return the result in p1.

Note
Requires the size of the arrays that are passed in to equal the number of PC terms

◆ MultiplyInPlace() [2/2]

void PCSet::MultiplyInPlace ( double * p1,
const double & a ) const

Multiply PC expansions given by double* argument p1 with scalar a and return the result in p1.

◆ Polyn() [1/2]

void PCSet::Polyn ( const Array1D< double > & polycf,
const Array1D< double > & p1,
Array1D< double > & p2 ) const

Evaluates a polynomial of PC that is given by Array1D argument p1. Polynomial coefficients are given in the Array1D argument polycf. The output PC is contained in Array1D argument p2.

Note
Requires the size of array p1 to equal the number of PC terms

◆ Polyn() [2/2]

void PCSet::Polyn ( const double * polycf,
int npoly,
const double * p1,
double * p2 ) const

Evaluates a polynomial of PC that is given in double* argument p1. Polynomial coefficients are given in double* argument polycf of size npoly. The output PC is contained in double* argument p2.

Note
Recursive algorithm is implemented.

◆ PolynMulti()

void PCSet::PolynMulti ( const Array1D< double > & polycf,
const Array2D< int > & mindex,
const Array2D< double > & p1,
Array1D< double > & p2 ) const

Evaluates a multivariate polynomial of a set of PC inputs given by Array2D argument p1 (each column of p1 is a PC input). Polynomial coefficients are given in Array1D argument polycf. Multiindex set for the multivariate polynomial is given in Array2D argument mindex. The output PC is contained in Array1D argument p2.

Note
Requires the size of the array polycf to equal the first dimension of argument mindex
Requires the size of the array p1 to equal (the number of PC terms) X (second dimension of argument mindex)
Uses a recursive algorithm
Out of convenience, this function so far is implemented for Array classes, not double* arrays.
Todo
A double* version should be added.

◆ PrintMultiIndex()

void PCSet::PrintMultiIndex ( ) const

Print information on the screen.

Print the multi-indices for all terms on the screen

◆ PrintMultiIndexNormSquared()

void PCSet::PrintMultiIndexNormSquared ( ) const

For all terms, print their multi-index and norm^2 on the screen.

◆ Prod() [1/2]

void PCSet::Prod ( const Array1D< double > & p1,
const Array1D< double > & p2,
Array1D< double > & p3 ) const

Multipy two PC expansions given by Array1D arguments p1 and p2, and return the result in p3.

Note
Requires the size of the input arrays to equal the number of PC terms

◆ Prod() [2/2]

void PCSet::Prod ( const double * p1,
const double * p2,
double * p3 ) const

Multiply two PC expansions given by double* arguments p1 and p2, and return the result in p3.

◆ Prod3() [1/2]

void PCSet::Prod3 ( const Array1D< double > & p1,
const Array1D< double > & p2,
const Array1D< double > & p3,
Array1D< double > & p4 ) const

Multipy three PC expansions given by Array1D arguments p1, p2, and p3, and return the result in p4.

Note
Requires the size of the input arrays to equal the number of PC terms

◆ Prod3() [2/2]

void PCSet::Prod3 ( const double * p1,
const double * p2,
const double * p3,
double * p4 ) const

Multiply three PC expansions given by double* arguments p1, p2, and p3, and return the result in p4.

◆ RPow() [1/2]

void PCSet::RPow ( const Array1D< double > & p1,
Array1D< double > & p2,
const double & a ) const

Evaluate power a (a real number) of PC expansion given by Array1D argument p1, and return the result in Array1D argument p2.

Note
Requires the size of the arrays that are passed in to equal the number of PC terms

◆ RPow() [2/2]

void PCSet::RPow ( const double * p1,
double * p2,
const double & a ) const

Evaluate power a (a real number) of PC expansion given by double* argument p1, and return the result in p2. The power is computed as p1^a = exp(a*log(p1)), where log(p1) is evaluated either via Taylor series or via integration depending on the value of parameter logMethod_.

◆ SeedBasisRandNumGen()

void PCSet::SeedBasisRandNumGen ( const int & seed) const

Random sample generator functions.

Reseed the random number generator used for the sampling of the PC variables and expansions

◆ SetGMRESDivTolerance()

void PCSet::SetGMRESDivTolerance ( const double & rTol)
inline

Set the relative tolerance for GMRES in Div routine.

◆ SetLogCompMethod()

void PCSet::SetLogCompMethod ( const LogCompMethod & logMethod)
inline

Set method of computing the log function.

Use the argument TaylorSeries to select the Taylor series approach or Integration to select the integration method.

◆ SetQd1d()

void PCSet::SetQd1d ( Array1D< double > & qdpts1d,
Array1D< double > & wghts1d,
int nqd )

Set the quadrature rule.

Obtain 1d quadrature points and weights

Note
This is used in triple or quadruple product computation for which the default quadrature is not enough

◆ SetQuadRule() [1/2]

void PCSet::SetQuadRule ( const string grid_type,
const string fs_type,
int param )

Set the quadrature points by specifying a grid type, a full/sparse indicator, and an integer parameter.

Full/sparse switch fs_type can be either 'full' or 'sparse' The parameter param is the number of points per dimension for full quadrature, and the level for sparse quadrature Options for grid_type are, besides the standard PC types, 'CC' (Clenshaw-Curtis), 'CCO' (Clenshaw-Curtis open), 'NC' (Newton-Cotes), 'NCO' (Newton-Cotes open), where open means that endpoints are expluded

Note
'NC', 'NCO' quadratures are the same as uniformly spaced grids
Todo
Need to improve it

◆ SetQuadRule() [2/2]

void PCSet::SetQuadRule ( Quad & quadRule)

Set a custom quadrature rule by pointing to the corresponding object.

◆ SetTaylorTermsMax()

void PCSet::SetTaylorTermsMax ( const int & maxTerm)
inline

Set maximum number of terms in Taylor series approximations.

◆ SetTaylorTolerance()

void PCSet::SetTaylorTolerance ( const double & rTol)
inline

Set relative tolerance for Taylor series approximations.

◆ SetVerbosity()

void PCSet::SetVerbosity ( int verbosity)
inline

Other.

Set the verbosity level

Note
Currently, the values of 0, 1 and 2 are implemented

◆ StDv() [1/2]

double PCSet::StDv ( const Array1D< double > & p) const

Returns the standard deviation of PC expansion p (Argument in Array1D<double> format)

Note
This assumes that the zeroth order term is the first one in the multi-index - this assumption does not hold in general
Todo
Lift the assumption by looking for the constant term in the multiindex
Note
For a more general implementation, see ComputeVarFrac()

◆ StDv() [2/2]

double PCSet::StDv ( const double * p) const

Returns the standard deviation of PC expansion p in a double* format.

Note
This assumes that the zeroth order term is the first one in the multi-index - this assumption does not hold in general
Todo
Lift the assumption by looking for the constant term in the multiindex

◆ Subtract() [1/2]

void PCSet::Subtract ( const Array1D< double > & p1,
const Array1D< double > & p2,
Array1D< double > & p3 ) const

Subtract PC expansion p2 from p1, and return the result in p3, with all arguments given as Array1D structures.

Note
Requires the size of the arrays that are passed in to equal the number of PC terms

◆ Subtract() [2/2]

void PCSet::Subtract ( const double * p1,
const double * p2,
double * p3 ) const

Subtract PC expansion p2 from p1, and return the result in p3, with all arguments given as double*.

◆ SubtractInPlace() [1/2]

void PCSet::SubtractInPlace ( Array1D< double > & p1,
const Array1D< double > & p2 ) const

Subtract PC expansion p2 from p1, and return the result in p1, with all arguments given as Array1D structures.

Note
Requires the size of the arrays that are passed in to equal the number of PC terms

◆ SubtractInPlace() [2/2]

void PCSet::SubtractInPlace ( double * p1,
const double * p2 ) const

Subtract PC expansion p2 from p1, and return the result in p1, with all arguments given as double*.

Member Data Documentation

◆ alpha_

double PCSet::alpha_
private

Parameter alpha for PCs that require a parameter (LG,SW,JB)

◆ beta_

double PCSet::beta_
private

Parameter beta for PCs that require two parameters (SW,JB)

◆ CVabst_

double PCSet::CVabst_
private

CVODE parameter: absolute tolerance.

◆ CVinitstep_

double PCSet::CVinitstep_
private

CVODE parameter: initial step size.

◆ CVmaxnumsteps_

int PCSet::CVmaxnumsteps_
private

CVODE parameter: maximal number of steps.

◆ CVmaxord_

int PCSet::CVmaxord_
private

CVODE parameter: maximal order.

◆ CVmaxstep_

double PCSet::CVmaxstep_
private

CVODE parameter: maximal step size.

◆ CVrelt_

double PCSet::CVrelt_
private

CVODE parameter: relative tolerance.

◆ iProd2_

Array1D<Array1D<int> > PCSet::iProd2_
private

i-indices of <\Psi_i \Psi_j \Psi_k> terms that are not zero, for all k

Note
Stored as a vector over k, with each element being a vector of i-indices itself

◆ iProd3_

Array1D<Array1D<int> > PCSet::iProd3_
private

i-indices of <\Psi_i \Psi_j \Psi_k \Psi_l> terms that are not zero, for all l

Note
Stored as a vector over l, with each element being a vector of i-indices itself

◆ jProd2_

Array1D<Array1D<int> > PCSet::jProd2_
private

j-indices of <\Psi_i \Psi_j \Psi_k> terms that are not zero, for all k

Note
Stored as a vector over k, with each element being a vector of j-indices itself

◆ jProd3_

Array1D<Array1D<int> > PCSet::jProd3_
private

j-indices of <\Psi_i \Psi_j \Psi_k \Psi_l> terms that are not zero, for all l

Note
Stored as a vector over l, with each element being a vector of j-indices itself

◆ kProd3_

Array1D<Array1D<int> > PCSet::kProd3_
private

k-indices of <\Psi_i \Psi_j \Psi_k \Psi_l> terms that are not zero, for all l

Note
Stored as a vector over l, with each element being a vector of k-indices itself

◆ logMethod_

LogCompMethod PCSet::logMethod_
private

Flag for method to compute log: TaylorSeries or Integration.

◆ maxorddim_

int PCSet::maxorddim_
private

Maximal order within all dimensions.

◆ maxOrders_

Array1D<int> PCSet::maxOrders_
private

Array of maximum orders requested if custom(HDMR) ordering is requested.

◆ maxOrdPerDim_

Array1D<int> PCSet::maxOrdPerDim_
private

Array of maximum orders per dimension.

◆ maxTermTaylor_

int PCSet::maxTermTaylor_
private

Max number of terms in Taylor series approximations.

◆ multiIndex_

Array2D<int> PCSet::multiIndex_
private

Array to store multi-index: multiIndex_(ipc,idim) contains the order of the basis function associated with dimension idim, for the ipc-th term in the PC expansion.

◆ my_index_

int PCSet::my_index_
private

Index of this class.

◆ narg_

int PCSet::narg_
private

Number of free parameters to specify the basis.

◆ nDim_

const int PCSet::nDim_
private

Number of stochastic dimensions (degrees of freedom) in the PC representation.

◆ next_index_

int PCSet::next_index_ = 0
staticprivate

index of next object in map

◆ nPCTerms_

int PCSet::nPCTerms_
private

Total number of terms in the PC expansions.

◆ nQuadPoints_

int PCSet::nQuadPoints_
private

Number of quadrature points used.

◆ omap_

PCSet::OMap_t * PCSet::omap_ = NULL
staticprivate

Map to connect integer indexes with pointers to this class.

◆ order_

int PCSet::order_
private

Order of the PC representation.

◆ p_basis_

PCBasis* PCSet::p_basis_
private

Pointer to the class that defines the basis type and functions.

◆ pcSeq_

string PCSet::pcSeq_
private

String indicator of multiindex ordering.

◆ pcType_

string PCSet::pcType_
private

String indicator of PC type.

◆ psi_

Array2D<double> PCSet::psi_
private

Array to store basis functions evaluated at quadrature points for each order: psi_(iqp,ipc) contains the value of the polynomial chaos ipc-th basis at the location of quadrature point iqp.

◆ psiIJKLProd3_

Array1D<Array1D<double> > PCSet::psiIJKLProd3_
private

<\Psi_i \Psi_j \Psi_k \Psi_l> terms that are not zero, for all l

Note
Stored as a vector over l, with each element being a vector of <\Psi_i \Psi_j \Psi_k \Psi_l> values

◆ psiIJKProd2_

Array1D<Array1D<double> > PCSet::psiIJKProd2_
private

<\Psi_i \Psi_j \Psi_k> terms that are not zero, for all k

Note
Stored as a vector over k, with each element being a vector of <\Psi_i \Psi_j \Psi_k> values

◆ psiSq_

Array1D<double> PCSet::psiSq_
private

Array with the norms squared of the basis functions, corresponding to each term in the PC expansion.

◆ quadIndices_

Array2D<int> PCSet::quadIndices_
private

Array to store quadrature point indexing; useful for nested rules.

◆ quadPoints_

Array2D<double> PCSet::quadPoints_
private

Array to store quadrature points.

◆ quadWeights_

Array1D<double> PCSet::quadWeights_
private

Array to store quadrature weights.

◆ rTolGMRESDiv_

double PCSet::rTolGMRESDiv_
private

GMRES tolerance in Div()

◆ rTolTaylor_

double PCSet::rTolTaylor_
private

Relative tolerance for Taylor series approximations.

◆ SMALL_

double PCSet::SMALL_
private

Tolerance to avoid floating-point errors.

◆ spType_

string PCSet::spType_
private

String indicator of ISP or NISP implementation type.

◆ uqtkverbose_

int PCSet::uqtkverbose_
private

Verbosity level.

Note
Currently the values of 0, 1 or 2 are implemented.

The documentation for this class was generated from the following files: