$extrastylesheet
Dakota
Version 6.2
|
Derived nondeterministic class that generates N-dimensional Smolyak sparse grids for numerical evaluation of expectation integrals over independent standard random variables. More...
Public Member Functions | |
NonDSparseGrid (Model &model, const UShortArray &ssg_level_seq, const RealVector &dim_pref, short exp_coeffs_soln_approach, short driver_mode, short growth_rate=Pecos::MODERATE_RESTRICTED_GROWTH, short refine_control=Pecos::NO_CONTROL, bool track_uniq_prod_wts=true, bool track_colloc_indices=true) | |
void | increment_grid () |
increment ssgDriver::ssgLevel | |
void | increment_grid_weights (const RealVector &aniso_wts) |
update ssgDriver::ssgAnisoLevelWts and increment ssgDriver::ssgLevel based on specified anisotropic weighting | |
void | increment_specification_sequence () |
advance to next nevel in ssgLevelSeqSpec sequence | |
const std::set< UShortArray > & | active_multi_index () const |
returns SparseGridDriver::active_multi_index() | |
const std::set< UShortArray > & | old_multi_index () const |
returns SparseGridDriver::old_multi_index() | |
void | print_smolyak_multi_index () const |
invokes SparseGridDriver::print_smolyak_multi_index() | |
void | initialize_sets () |
invokes SparseGridDriver::initialize_sets() | |
void | update_reference () |
invokes SparseGridDriver::update_reference() | |
void | increment_set (const UShortArray &set) |
invokes SparseGridDriver::push_trial_set() | |
int | increment_size () const |
invokes SparseGridDriver::unique_trial_points() | |
void | restore_set () |
invokes SparseGridDriver::restore_set() | |
void | evaluate_set () |
invokes SparseGridDriver::compute_trial_grid() | |
void | decrement_set () |
invokes SparseGridDriver::pop_trial_set() | |
void | update_sets (const UShortArray &set_star) |
invokes SparseGridDriver::update_sets() | |
void | print_final_sets (bool converged_within_tol) |
invokes SparseGridDriver::print_final_sets(bool) | |
void | finalize_sets () |
invokes SparseGridDriver::finalize_sets() | |
void | evaluate_grid_increment () |
invokes SparseGridDriver::evaluate_grid_increment() | |
int | num_samples () const |
Protected Member Functions | |
NonDSparseGrid (ProblemDescDB &problem_db, Model &model) | |
constructor | |
~NonDSparseGrid () | |
destructor | |
void | initialize_grid (const std::vector< Pecos::BasisPolynomial > &poly_basis) |
initialize integration grid by drawing from polynomial basis settings | |
void | get_parameter_sets (Model &model) |
Returns one block of samples (ndim * num_samples) | |
void | reset () |
restore initial state for repeated sub-iterator executions | |
void | sampling_reset (int min_samples, bool all_data_flag, bool stats_flag) |
Private Attributes | |
Pecos::SparseGridDriver * | ssgDriver |
convenience pointer to the numIntDriver representation | |
UShortArray | ssgLevelSeqSpec |
the user specification for the Smolyak sparse grid level, defining a sequence of refinement levels. | |
unsigned short | ssgLevelRef |
reference point (e.g., lower bound) for the Smolyak sparse grid level maintained within ssgDriver |
Derived nondeterministic class that generates N-dimensional Smolyak sparse grids for numerical evaluation of expectation integrals over independent standard random variables.
This class is used by NonDPolynomialChaos and NonDStochCollocation, but could also be used for general numerical integration of moments. It employs 1-D Clenshaw-Curtis and Gaussian quadrature rules within Smolyak sparse grids.
NonDSparseGrid | ( | Model & | model, |
const UShortArray & | ssg_level_seq, | ||
const RealVector & | dim_pref, | ||
short | exp_coeffs_soln_approach, | ||
short | driver_mode, | ||
short | growth_rate = Pecos::MODERATE_RESTRICTED_GROWTH , |
||
short | refine_control = Pecos::NO_CONTROL , |
||
bool | track_uniq_prod_wts = true , |
||
bool | track_colloc_indices = true |
||
) |
This alternate constructor is used for on-the-fly generation and evaluation of sparse grids within PCE and SC.
References NonDIntegration::numIntDriver, and NonDSparseGrid::ssgDriver.
NonDSparseGrid | ( | ProblemDescDB & | problem_db, |
Model & | model | ||
) | [protected] |
constructor
This constructor is called for a standard letter-envelope iterator instantiation. In this case, set_db_list_nodes has been called and probDescDB can be queried for settings from the method specification. It is not currently used, as there is not a separate sparse_grid method specification.
References Model::aleatory_distribution_parameters(), NonDIntegration::check_variables(), Iterator::convergenceTol, NonDIntegration::dimPrefSpec, ProblemDescDB::get_bool(), ProblemDescDB::get_short(), ProblemDescDB::get_ushort(), Iterator::iteratedModel, Iterator::maxEvalConcurrency, Iterator::maxIterations, NonD::natafTransform, NonDIntegration::numIntDriver, Iterator::outputLevel, Iterator::probDescDB, NonDSparseGrid::ssgDriver, and NonDSparseGrid::ssgLevelRef.
int num_samples | ( | ) | const [inline, virtual] |
Return current number of evaluation points. Since the calculation of samples, collocation points, etc. might be costly, provide a default implementation here that backs out from the maxEvalConcurrency.
Reimplemented from Analyzer.
References NonDSparseGrid::ssgDriver.
void sampling_reset | ( | int | min_samples, |
bool | all_data_flag, | ||
bool | stats_flag | ||
) | [protected, virtual] |
used by DataFitSurrModel::build_global() to publish the minimum number of points needed from the sparse grid routine in order to build a particular global approximation.
Reimplemented from Iterator.
References NonDSparseGrid::ssgDriver, and NonDSparseGrid::ssgLevelRef.