UQTk: Uncertainty Quantification Toolkit 3.1.5
|
Contains all basis type specific definitions and operations needed to generate a PCSet. More...
#include <PCBasis.h>
Public Member Functions | |
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. | |
~PCBasis () | |
Destructor. | |
void | Init1dQuadPoints (int qdpts) |
Initialize the quadrature points and weights and store the information in arrays quadPoints_, quadWeights_,quadIndices_. | |
void | Eval1dBasisAtQuadPoints () |
Evaluate polynomial 1d basis functions at quadrature points and store in the private variable psi1d_. | |
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 custPoints Returns the evaluations in the first argument psi, where the number of rows are the number of points, and columns correspond to successive orders. | |
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 functions for all orders in the passed Array1D array (indexed by their order), also returns the highest-order value. | |
double | EvalBasis (const double &xi, const int kord, double *basisEvals) const |
Evaluate 1d basis functions for the given value of random variable xi. Return the value of the basis functions for all orders in the passed double * array (indexed by their order), also returns the highest-order value. | |
void | Eval1dNormSq_Exact (int kord) |
Evaluate the norms (squared) of the basis functions exactly and stores in the private array psi1dSqExact_. | |
void | EvalDerivBasis (const double &xi, Array1D< double > &basisDEvals) |
Evaluate derivative of 1d non-normalized Legendre basis. | |
void | Eval1dDerivBasisAtCustPoints (Array2D< double > &dpsi, int kord, const Array1D< double > &custPoints) |
void | Eval2ndDerivBasis (const double &xi, Array1D< double > &ddP) |
void | Eval2ndDerivCustPoints (Array2D< double > &psi, int kord, Array1D< double > &custPoints) |
void | Get1dNormsSq (Array1D< double > &psi1dSq) const |
Get the norms-squared of the basis functions. Returns the values for each basis function in the passed Array1D array. | |
void | Get1dNormsSqExact (Array1D< double > &psi1dSqExact) const |
Get the analytic norms-squared of the basis functions. Returns the values for each basis function in the passed Array1D array. | |
void | GetRandSample (Array1D< double > &randSamples) |
Get samples of the random variables associated with the current PC basis functions and return them in the 1D array randSamples. Take as many samples as the length of the array randSamples. | |
void | GetRandSample (double *randSamples, const int &nSamp) |
Get nSamp samples of the random variables associated with the current PC basis functions and return them in the double* randSamples. | |
int | GetSeed () const |
Get the random number generator seed. | |
void | SeedRandNumGen (const int &seed) |
Function to (re)seed the random number generator used to sample the Basis functions. | |
void | GetQuadRule (Array2D< double > &qPoints, Array1D< double > &qWeights, Array2D< int > &qIndices) |
Get the quadrature integration information. | |
void | GetQuadPoints (Array2D< double > &quadPoints) const |
Get the quadrature points in the passed Array2D array. | |
void | GetQuadWeights (Array1D< double > &quadWeights) const |
Get the quadrature weights in the passed Array1D array. | |
void | GetQuadIndices (Array2D< int > &quadIndices) const |
Get the quadrature points' indices in the passed Array1D array. | |
void | GetBasisAtQuadPoints (Array2D< double > &psi1d) const |
Get the basis values at quadrature points in the passed Array2D array. | |
string | GetPCType () const |
Get the PC type. | |
double | GetAlpha () const |
Get the value of the parameter alpha. | |
double | GetBeta () const |
Get the value of the parameter beta. | |
Private Member Functions | |
PCBasis (const PCBasis &obj) | |
Dummy default constructor, which should not be used as it is not well defined Therefore we make it private so it is not accessible. | |
void | Eval1dNormSq (int kord) |
Evaluate the norms (squared) of the basis functions and stores in the private array psi1dSq_. | |
double | NormSq_Exact (int kord) |
Evaluate 1d norm of order kord exactly. | |
Private Attributes | |
string | type_ |
String indicator of type of basis functions used. | |
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 only for nested rules. | |
Array2D< double > | psi1d_ |
Array to store basis functions evaluated at quadrature points for each order: psi1d_(iqp,iord) contains the value of the polynomial chaos basis of order iord at the location of quadrature point iqp. | |
Array1D< double > | psi1dSq_ |
Array with the norms squared of the 1D basis functions for each order. | |
Array1D< double > | psi1dSqExact_ |
Array with the exact norms squared of the 1D basis functions for each order. | |
int | maxord_ |
Maximal order of any dimension. | |
int | narg_ |
Number of 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) | |
dsfmt_t | rnstate_ |
Random sequence state for dsfmt. | |
int | rSeed_ |
The seed used for the random number generators that sample the xi's in the basis functions. | |
Contains all basis type specific definitions and operations needed to generate a PCSet.
PCBasis::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.
Currently, the only valid types are Hermite-Gaussian, denoted with "HG", Legendre-Uniform, denoted with "LU", or Laguerre-Gamma, denoted with "LG". (Where the shape parameter for the Gamma distribution is alpha + 1 = 2)
|
inline |
Destructor.
|
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.
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 PCBasis::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 custPoints Returns the evaluations in the first argument psi, where the number of rows are the number of points, and columns correspond to successive orders.
void PCBasis::Eval1dBasisAtQuadPoints | ( | ) |
Evaluate polynomial 1d basis functions at quadrature points and store in the private variable psi1d_.
void PCBasis::Eval1dDerivBasisAtCustPoints | ( | Array2D< double > & | dpsi, |
int | kord, | ||
const Array1D< double > & | custPoints ) |
|
private |
Evaluate the norms (squared) of the basis functions and stores in the private array psi1dSq_.
void PCBasis::Eval1dNormSq_Exact | ( | int | kord | ) |
Evaluate the norms (squared) of the basis functions exactly and stores in the private array psi1dSqExact_.
void PCBasis::Eval2ndDerivBasis | ( | const double & | xi, |
Array1D< double > & | ddP ) |
void PCBasis::Eval2ndDerivCustPoints | ( | Array2D< double > & | psi, |
int | kord, | ||
Array1D< double > & | custPoints ) |
double PCBasis::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 functions for all orders in the passed Array1D array (indexed by their order), also returns the highest-order value.
double PCBasis::EvalBasis | ( | const double & | xi, |
const int | kord, | ||
double * | basisEvals ) const |
Evaluate 1d basis functions for the given value of random variable xi. Return the value of the basis functions for all orders in the passed double * array (indexed by their order), also returns the highest-order value.
void PCBasis::EvalDerivBasis | ( | const double & | xi, |
Array1D< double > & | basisDEvals ) |
Evaluate derivative of 1d non-normalized Legendre basis.
|
inline |
Get the norms-squared of the basis functions. Returns the values for each basis function in the passed Array1D array.
|
inline |
Get the analytic norms-squared of the basis functions. Returns the values for each basis function in the passed Array1D array.
|
inline |
Get the value of the parameter alpha.
|
inline |
Get the basis values at quadrature points in the passed Array2D array.
|
inline |
Get the value of the parameter beta.
|
inline |
Get the PC type.
|
inline |
Get the quadrature points' indices in the passed Array1D array.
|
inline |
Get the quadrature points in the passed Array2D array.
void PCBasis::GetQuadRule | ( | Array2D< double > & | qPoints, |
Array1D< double > & | qWeights, | ||
Array2D< int > & | qIndices ) |
Get the quadrature integration information.
|
inline |
Get the quadrature weights in the passed Array1D array.
void PCBasis::GetRandSample | ( | Array1D< double > & | randSamples | ) |
Get samples of the random variables associated with the current PC basis functions and return them in the 1D array randSamples. Take as many samples as the length of the array randSamples.
void PCBasis::GetRandSample | ( | double * | randSamples, |
const int & | nSamp ) |
Get nSamp samples of the random variables associated with the current PC basis functions and return them in the double* randSamples.
|
inline |
Get the random number generator seed.
void PCBasis::Init1dQuadPoints | ( | int | qdpts | ) |
Initialize the quadrature points and weights and store the information in arrays quadPoints_, quadWeights_,quadIndices_.
|
private |
Evaluate 1d norm of order kord exactly.
void PCBasis::SeedRandNumGen | ( | const int & | seed | ) |
Function to (re)seed the random number generator used to sample the Basis functions.
|
private |
Parameter alpha for PCs that require a parameter (LG,SW,JB)
|
private |
Parameter beta for PCs that require two parameters (SW,JB)
|
private |
Maximal order of any dimension.
|
private |
Number of parameters to specify the basis.
|
private |
Array to store basis functions evaluated at quadrature points for each order: psi1d_(iqp,iord) contains the value of the polynomial chaos basis of order iord at the location of quadrature point iqp.
|
private |
Array with the norms squared of the 1D basis functions for each order.
|
private |
Array with the exact norms squared of the 1D basis functions for each order.
|
private |
Array to store quadrature point indexing; useful only for nested rules.
|
private |
Array to store quadrature points.
|
private |
Array to store quadrature weights.
|
private |
Random sequence state for dsfmt.
|
private |
The seed used for the random number generators that sample the xi's in the basis functions.
This seed is set to 1 during the class construction and can be reset with the SeedRandNumGen function
|
private |
String indicator of type of basis functions used.