60 PCBasis(
const string type=
"LU",
const double alpha=0.0,
const double betta=1.0,
const int maxord=10);
98 double EvalBasis(
const double &xi,
const int kord,
double *basisEvals)
const;
1D Array class for any type T
2D Array class for any type T
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
Contains all basis type specific definitions and operations needed to generate a PCSet.
Definition PCBasis.h:47
PCBasis(const string type="LU", const double alpha=0.0, const double betta=1.0, const int maxord=10)
Constructor: initializes the univariate basis type and order.
Definition PCBasis.cpp:48
Array2D< double > quadPoints_
Array to store quadrature points.
Definition PCBasis.h:194
PCBasis(const PCBasis &obj)
Dummy default constructor, which should not be used as it is not well defined Therefore we make it pr...
Definition PCBasis.h:179
double alpha_
Parameter alpha for PCs that require a parameter (LG,SW,JB)
Definition PCBasis.h:221
double EvalBasis(const double &xi, Array1D< double > &basisEvals) const
Evaluate 1d basis functions for the given value of random variable xi. Return the value of the basis ...
Definition PCBasis.cpp:147
void GetRandSample(Array1D< double > &randSamples)
Get samples of the random variables associated with the current PC basis functions and return them in...
Definition PCBasis.cpp:398
Array2D< int > quadIndices_
Array to store quadrature point indexing; useful only for nested rules.
Definition PCBasis.h:200
double beta_
Parameter beta for PCs that require two parameters (SW,JB)
Definition PCBasis.h:224
void GetBasisAtQuadPoints(Array2D< double > &psi1d) const
Get the basis values at quadrature points in the passed Array2D array.
Definition PCBasis.h:157
int GetSeed() const
Get the random number generator seed.
Definition PCBasis.h:137
void Get1dNormsSq(Array1D< double > &psi1dSq) const
Get the norms-squared of the basis functions. Returns the values for each basis function in the passe...
Definition PCBasis.h:119
double GetAlpha() const
Get the value of the parameter alpha.
Definition PCBasis.h:163
double GetBeta() const
Get the value of the parameter beta.
Definition PCBasis.h:166
Array1D< double > psi1dSq_
Array with the norms squared of the 1D basis functions for each order.
Definition PCBasis.h:209
void GetQuadRule(Array2D< double > &qPoints, Array1D< double > &qWeights, Array2D< int > &qIndices)
Get the quadrature integration information.
Definition PCBasis.cpp:458
int maxord_
Maximal order of any dimension.
Definition PCBasis.h:215
~PCBasis()
Destructor.
Definition PCBasis.h:64
void Init1dQuadPoints(int qdpts)
Initialize the quadrature points and weights and store the information in arrays quadPoints_,...
Definition PCBasis.cpp:95
void Eval1dNormSq(int kord)
Evaluate the norms (squared) of the basis functions and stores in the private array psi1dSq_.
Definition PCBasis.cpp:480
void EvalDerivBasis(const double &xi, Array1D< double > &basisDEvals)
Evaluate derivative of 1d non-normalized Legendre basis.
Definition PCBasis.cpp:277
void Eval2ndDerivBasis(const double &xi, Array1D< double > &ddP)
Definition PCBasis.cpp:338
void GetQuadIndices(Array2D< int > &quadIndices) const
Get the quadrature points' indices in the passed Array1D array.
Definition PCBasis.h:154
Array2D< double > psi1d_
Array to store basis functions evaluated at quadrature points for each order: psi1d_(iqp,...
Definition PCBasis.h:206
void GetQuadPoints(Array2D< double > &quadPoints) const
Get the quadrature points in the passed Array2D array.
Definition PCBasis.h:148
void Eval1dNormSq_Exact(int kord)
Evaluate the norms (squared) of the basis functions exactly and stores in the private array psi1dSqEx...
Definition PCBasis.cpp:497
double NormSq_Exact(int kord)
Evaluate 1d norm of order kord exactly.
Definition PCBasis.cpp:510
void Eval1dBasisAtCustPoints(Array2D< double > &psi, int kord, const Array1D< double > &custPoints)
Evaluate polynomial 1d basis functions up to the order kord at custom points given by an array custPo...
Definition PCBasis.cpp:121
void Eval1dDerivBasisAtCustPoints(Array2D< double > &dpsi, int kord, const Array1D< double > &custPoints)
Definition PCBasis.cpp:308
void Get1dNormsSqExact(Array1D< double > &psi1dSqExact) const
Get the analytic norms-squared of the basis functions. Returns the values for each basis function in ...
Definition PCBasis.h:123
string GetPCType() const
Get the PC type.
Definition PCBasis.h:160
string type_
String indicator of type of basis functions used.
Definition PCBasis.h:191
int narg_
Number of parameters to specify the basis.
Definition PCBasis.h:218
void SeedRandNumGen(const int &seed)
Function to (re)seed the random number generator used to sample the Basis functions.
Definition PCBasis.cpp:468
dsfmt_t rnstate_
Random sequence state for dsfmt.
Definition PCBasis.h:228
void Eval2ndDerivCustPoints(Array2D< double > &psi, int kord, Array1D< double > &custPoints)
Definition PCBasis.cpp:374
void Eval1dBasisAtQuadPoints()
Evaluate polynomial 1d basis functions at quadrature points and store in the private variable psi1d_.
Definition PCBasis.cpp:109
void GetQuadWeights(Array1D< double > &quadWeights) const
Get the quadrature weights in the passed Array1D array.
Definition PCBasis.h:151
Array1D< double > psi1dSqExact_
Array with the exact norms squared of the 1D basis functions for each order.
Definition PCBasis.h:212
Array1D< double > quadWeights_
Array to store quadrature weights.
Definition PCBasis.h:197
int rSeed_
The seed used for the random number generators that sample the xi's in the basis functions.
Definition PCBasis.h:236