Doxygen 1.9.1
Toolkit for Adaptive Stochastic Modeling and Non-Intrusive ApproximatioN: Tasmanian v8.1
TasGrid Namespace Reference

Encapsulates the Tasmanian Sparse Grid module. More...

Classes

class  TasmanianSparseGrid
 The master-class that represents an instance of a Tasmanian sparse grid. More...
 
class  Data2D
 
class  MultiIndexSet
 
class  StorageSet
 
struct  HandleDeleter
 
class  GpuVector
 Template class that wraps around a single GPU array, providing functionality that mimics std::vector.
 
struct  GpuEngine
 Wrapper class around calls GPU accelerated linear algebra libraries.
 
class  AccelerationDomainTransform
 
struct  AccelerationContext
 
struct  InternalSyclQueue
 
class  CacheLagrange
 
class  CacheLagrangeDerivative
 
class  CustomTabulated
 Class providing manipulation of custom tabulated rules, file I/O and structured access to the points, weights, and meta-data.
 
class  OneDimensionalWrapper
 A class to cache one dimensional rules, nodes, weight, meta-data, etc.
 
class  TableGaussPatterson
 
struct  NodeData
 
struct  TensorData
 
class  DynamicConstructorDataGlobal
 
struct  SimpleConstructData
 
class  CandidateManager
 
class  CompleteStorage
 
class  VectorToStreamBuffer
 

Typedefs

using int_gpu_lapack = std::int64_t
 
using ModelSignature = std::function< void(std::vector< double > const &x, std::vector< double > &y, size_t thread_id)>
 Signature of a model function to be used in the construction procedures. More...
 
using ModelSignatureMPI = std::function< void(std::vector< double > const &x, std::vector< double > &y)>
 Signature of a model function to be used in the MPI construction procedures. More...
 

Enumerations

enum  TypeIndexRelation { type_abeforeb , type_bbeforea , type_asameb }
 
enum  TypeDepth {
  type_none , type_level , type_curved , type_hyperbolic ,
  type_iptotal , type_qptotal , type_ipcurved , type_qpcurved ,
  type_iphyperbolic , type_qphyperbolic , type_tensor , type_iptensor ,
  type_qptensor
}
 Used by Global Sequence and Fourier grids, indicates the selection criteria. More...
 
enum  TypeOneDRule {
  rule_none , rule_clenshawcurtis , rule_clenshawcurtis0 , rule_fejer2 ,
  rule_chebyshev , rule_chebyshevodd , rule_leja , rule_lejaodd ,
  rule_rleja , rule_rlejadouble2 , rule_rlejadouble4 , rule_rlejaodd ,
  rule_rlejashifted , rule_rlejashiftedeven , rule_rlejashifteddouble , rule_maxlebesgue ,
  rule_maxlebesgueodd , rule_minlebesgue , rule_minlebesgueodd , rule_mindelta ,
  rule_mindeltaodd , rule_gausslegendre , rule_gausslegendreodd , rule_gausspatterson ,
  rule_gausschebyshev1 , rule_gausschebyshev1odd , rule_gausschebyshev2 , rule_gausschebyshev2odd ,
  rule_gaussgegenbauer , rule_gaussgegenbauerodd , rule_gaussjacobi , rule_gaussjacobiodd ,
  rule_gausslaguerre , rule_gausslaguerreodd , rule_gausshermite , rule_gausshermiteodd ,
  rule_customtabulated , rule_localp , rule_localp0 , rule_semilocalp ,
  rule_localpb , rule_wavelet , rule_fourier
}
 Used to specify the one dimensional family of rules that induces the sparse grid. More...
 
enum  TypeRefinement {
  refine_classic , refine_parents_first , refine_direction_selective , refine_fds ,
  refine_stable , refine_none
}
 Refinement strategy for local polynomial and wavelet grids. More...
 
enum  TypeAcceleration {
  accel_none , accel_cpu_blas , accel_gpu_default , accel_gpu_cublas ,
  accel_gpu_cuda , accel_gpu_magma
}
 Modes of acceleration. More...
 

Functions

TasmanianSparseGrid makeEmpty ()
 Returns an empty sparse grid. More...
 
TasmanianSparseGrid makeGlobalGrid (int dimensions, int outputs, int depth, TypeDepth type, TypeOneDRule rule, std::vector< int > const &anisotropic_weights=std::vector< int >(), double alpha=0.0, double beta=0.0, const char *custom_filename=nullptr, std::vector< int > const &level_limits=std::vector< int >())
 Factory method, creates a new grid and calls TasmanianSparseGrid::makeGlobalGrid(). More...
 
TasmanianSparseGrid makeSequenceGrid (int dimensions, int outputs, int depth, TypeDepth type, TypeOneDRule rule, std::vector< int > const &anisotropic_weights=std::vector< int >(), std::vector< int > const &level_limits=std::vector< int >())
 Factory method, creates a new grid and calls TasmanianSparseGrid::makeSequenceGrid(). More...
 
TasmanianSparseGrid makeLocalPolynomialGrid (int dimensions, int outputs, int depth, int order=1, TypeOneDRule rule=rule_localp, std::vector< int > const &level_limits=std::vector< int >())
 Factory method, creates a new grid and calls TasmanianSparseGrid::makeLocalPolynomialGrid(). More...
 
TasmanianSparseGrid makeWaveletGrid (int dimensions, int outputs, int depth, int order=1, std::vector< int > const &level_limits=std::vector< int >())
 Factory method, creates a new grid and calls TasmanianSparseGrid::makeWaveletGrid(). More...
 
TasmanianSparseGrid makeFourierGrid (int dimensions, int outputs, int depth, TypeDepth type, std::vector< int > const &anisotropic_weights=std::vector< int >(), std::vector< int > const &level_limits=std::vector< int >())
 Factory method, creates a new grid and calls TasmanianSparseGrid::makeFourierGrid(). More...
 
TasmanianSparseGrid readGrid (const char *filename)
 Factory method, creates a new grid and calls TasmanianSparseGrid::read(). More...
 
TasmanianSparseGrid readGrid (std::string const &filename)
 Overload using std::string. More...
 
TasmanianSparseGrid copyGrid (TasmanianSparseGrid const &source, int outputs_begin=0, int outputs_end=-1)
 Returns a grid that is a copy of the source. More...
 
template<typename T >
std::vector< T > spltVector2D (std::vector< T > const &x, size_t stride, int ibegin, int iend)
 
template<typename >
void deleteHandle (int *)
 
CustomTabulated getSubrules (CustomTabulated &ct, int start_index, int stride, std::string description)
 
template<class T >
std::vector< const T * > makeReverseReferenceVector (const std::forward_list< T > &list)
 
template<bool use_ascii>
void writeNodeDataList (const std::forward_list< NodeData > &data, std::ostream &os)
 
template<typename iomode >
std::forward_list< NodeData > readNodeDataList (std::istream &is, size_t num_dimensions, size_t num_outputs)
 
template<typename iomode >
std::forward_list< TensorData > readTensorDataList (std::istream &is, size_t num_dimensions)
 
template<class RuleLike >
std::vector< double > listToNodes (std::forward_list< NodeData > const &node_list, size_t num_dimensions, RuleLike const &rule)
 
template<typename callable_method >
std::vector< double > listToLocalNodes (std::forward_list< NodeData > const &node_list, size_t num_dimensions, callable_method rule)
 
double poly_eval (const std::vector< double > &roots, double x)
 
double lagrange_eval (size_t idx, const std::vector< double > &roots, double x)
 
template<bool is_symmetric>
void getGaussNodesAndWeights (const int n, const std::vector< double > &ref_points, const std::vector< double > &ref_weights, std::vector< std::vector< double >> &points_cache, std::vector< std::vector< double >> &weights_cache)
 
TasGrid::CustomTabulated getShiftedExoticQuadrature (const int n, const double shift, const std::vector< double > &shifted_weights, const std::vector< double > &ref_points, const char *description, const bool is_symmetric=false)
 
void shiftReferenceWeights (std::vector< double > const &vals, double shift, std::vector< double > &ref_weights)
 
TasGrid::CustomTabulated getExoticQuadrature (const int num_levels, const double shift, const TasGrid::TasmanianSparseGrid &grid, const char *description, const bool is_symmetric=false)
 Constructs an exotic quadrature from a sparse grid object. More...
 
TasGrid::CustomTabulated getExoticQuadrature (const int num_levels, const double shift, std::function< double(double)> weight_fn, const int num_ref_points, const char *description, const bool is_symmetric=false)
 Similar to getExoticQuadrature() but the weight function is defined by a lambda expression. More...
 
template<bool parallel_construction, bool use_initial_guess>
void constructCommon (ModelSignature model, size_t max_num_points, size_t num_parallel_jobs, size_t max_samples_per_job, TasmanianSparseGrid &grid, std::function< std::vector< double >(TasmanianSparseGrid &)> candidates, std::string const &checkpoint_filename)
 
template<bool parallel_construction = TasGrid::mode_parallel, bool initial_guess = no_initial_guess>
void constructSurrogate (ModelSignature model, size_t max_num_points, size_t num_parallel_jobs, size_t max_samples_per_job, TasmanianSparseGrid &grid, double tolerance, TypeRefinement criteria, int output=-1, std::vector< int > const &level_limits=std::vector< int >(), std::string const &checkpoint_filename=std::string())
 Construct a sparse grid surrogate to the model defined by the lambda. More...
 
template<bool parallel_construction = TasGrid::mode_parallel, bool initial_guess = no_initial_guess>
void constructSurrogate (ModelSignature model, size_t max_num_points, size_t num_parallel_jobs, size_t max_samples_per_job, TasmanianSparseGrid &grid, TypeDepth type, std::vector< int > const &anisotropic_weights=std::vector< int >(), std::vector< int > const &level_limits=std::vector< int >(), std::string const &checkpoint_filename=std::string())
 Construct a sparse grid surrogate to the model defined by the lambda. More...
 
template<bool parallel_construction = TasGrid::mode_parallel, bool initial_guess = no_initial_guess>
void constructSurrogate (ModelSignature model, size_t max_num_points, size_t num_parallel_jobs, size_t max_samples_per_job, TasmanianSparseGrid &grid, TypeDepth type, int output, std::vector< int > const &level_limits=std::vector< int >(), std::string const &checkpoint_filename=std::string())
 Construct a sparse grid surrogate to the model defined by the lambda. More...
 
bool hasGPUBasis (TasmanianSparseGrid const &grid)
 
template<typename scalar_type >
void generateCoefficientsGPU (double const data_points[], int num_data, scalar_type model_values[], double tolerance, TasmanianSparseGrid &grid)
 
template<typename scalar_type >
Data2D< scalar_type > generateCoefficients (double const data_points[], int num_data, double const model_values[], double tolerance, TasmanianSparseGrid &grid)
 
template<typename scalar_type >
void loadUnstructuredDataL2tmpl (double const data_points[], int num_data, double const model_values[], double tolerance, TasmanianSparseGrid &grid)
 
void loadUnstructuredDataL2 (double const data_points[], int num_data, double const model_values[], double tolerance, TasmanianSparseGrid &grid)
 Construct a sparse grid surrogate using a least-squares fit. More...
 
void loadUnstructuredDataL2 (std::vector< double > const &data_points, std::vector< double > const &model_values, double tolerance, TasmanianSparseGrid &grid)
 Overload that used vectors and infers the number of data points from the size of the vectors. More...
 
template<bool use_initial_guess>
void mpiConstructCommon (ModelSignatureMPI model, int num_dimensions, int num_outputs, size_t max_num_points, size_t max_samples_per_job, size_t max_num_ranks, int tagx, int tagy, int root, MPI_Comm comm, TasmanianSparseGrid &grid, std::function< std::vector< double >(TasmanianSparseGrid &)> candidates, std::string const &checkpoint_filename)
 
template<bool use_initial_guess = no_initial_guess>
void mpiConstructSurrogate (ModelSignatureMPI model, int num_dimensions, int num_outputs, size_t max_num_points, size_t max_samples_per_job, size_t max_num_ranks, int tagx, int tagy, int root, MPI_Comm comm, TasmanianSparseGrid &grid, double tolerance, TypeRefinement criteria, int output=-1, std::vector< int > const &level_limits=std::vector< int >(), std::vector< double > const &scale_correction=std::vector< double >(), std::string const &checkpoint_filename=std::string())
 Construct a sparse grid surrogate to the model defined by the lambda, MPI version. More...
 
template<bool use_initial_guess = no_initial_guess>
void mpiConstructSurrogate (ModelSignatureMPI model, int num_dimensions, int num_outputs, size_t max_num_points, size_t max_samples_per_job, size_t max_num_ranks, int tagx, int tagy, int root, MPI_Comm comm, TasmanianSparseGrid &grid, TypeDepth type, std::vector< int > const &anisotropic_weights=std::vector< int >(), std::vector< int > const &level_limits=std::vector< int >(), std::string const &checkpoint_filename=std::string())
 Construct a sparse grid surrogate to the model defined by the lambda, MPI version. More...
 
template<bool use_initial_guess = no_initial_guess>
void mpiConstructSurrogate (ModelSignatureMPI model, int num_dimensions, int num_outputs, size_t max_num_points, size_t max_samples_per_job, size_t max_num_ranks, int tagx, int tagy, int root, MPI_Comm comm, TasmanianSparseGrid &grid, TypeDepth type, int output, std::vector< int > const &level_limits=std::vector< int >(), std::string const &checkpoint_filename=std::string())
 Construct a sparse grid surrogate to the model defined by the lambda, MPI version. More...
 
template<bool use_initial_guess = no_initial_guess>
void mpiConstructWorker (ModelSignatureMPI model, int num_dimensions, int num_outputs, size_t max_samples_per_job, size_t max_num_ranks, int tagx, int tagy, int root, MPI_Comm comm)
 Executes the worker (i.e., non-root) portion of the MPI sampling. More...
 
int getMPIRank (MPI_Comm comm)
 
template<bool binary = TasGrid::mode_binary>
int MPIGridSend (TasmanianSparseGrid const &grid, int destination, int tag, MPI_Comm comm)
 Send a grid to another process in the MPI comm. More...
 
template<bool binary = TasGrid::mode_binary>
int MPIGridRecv (TasmanianSparseGrid &grid, int source, int tag, MPI_Comm comm, MPI_Status *status=MPI_STATUS_IGNORE)
 Receive a grid from another process in the MPI comm. More...
 
template<bool binary = TasGrid::mode_binary>
int MPIGridBcast (TasmanianSparseGrid &grid, int root, MPI_Comm comm)
 Broadcast a grid to all processes in an MPI comm. More...
 
template<bool binary = TasGrid::mode_binary>
int MPIGridScatterOutputs (TasmanianSparseGrid const &source, TasmanianSparseGrid &destination, int root, int tag, MPI_Comm comm)
 Split the grid across the comm where each rank receives an equal portion of the total outputs. More...
 
template<bool parallel_construction = true, bool overwrite_loaded = false>
void loadNeededValues (std::function< void(double const x[], double y[], size_t thread_id)> model, TasmanianSparseGrid &grid, size_t num_threads)
 Loads the current grid with model values, does not perform any refinement. More...
 
template<bool parallel_construction = true, bool overwrite_loaded = false>
void loadNeededValues (std::function< void(std::vector< double > const &x, std::vector< double > &y, size_t thread_id)> model, TasmanianSparseGrid &grid, size_t num_threads)
 Overload that uses vectors for the model inputs and outputs. More...
 
template<bool parallel_construction = true, bool overwrite_loaded = false>
void loadNeededPoints (std::function< void(double const x[], double y[], size_t thread_id)> model, TasmanianSparseGrid &grid, size_t num_threads)
 Alias to loadNeededValues(), array variant.
 
template<bool parallel_construction = true, bool overwrite_loaded = false>
void loadNeededPoints (std::function< void(std::vector< double > const &x, std::vector< double > &y, size_t thread_id)> model, TasmanianSparseGrid &grid, size_t num_threads)
 Alias to loadNeededValues(), vector variant.
 

Variables

constexpr TypeAcceleration accel_gpu_hip = accel_gpu_cuda
 At the front API, the HIP and CUDA options are equivalent, see TasGrid::TypeAcceleration.
 
constexpr TypeAcceleration accel_gpu_rocblas = accel_gpu_cublas
 At the front API, the HIP and CUDA options are equivalent, see TasGrid::TypeAcceleration.
 
constexpr bool mode_ascii = false
 Constant allowing for more expressive selection of ascii and binary mode in IO methods.
 
constexpr bool mode_binary = true
 Constant allowing for more expressive selection of ascii and binary mode in IO methods.
 
InternalSyclQueue test_queue
 
constexpr bool mode_parallel = true
 Allows for expressive calls to TasGrid::constructSurrogate().
 
constexpr bool mode_sequential = false
 Allows for expressive calls to TasGrid::constructSurrogate().
 
constexpr bool with_initial_guess = true
 Allows for expressive calls to TasGrid::constructSurrogate().
 
constexpr bool no_initial_guess = false
 Allows for expressive calls to TasGrid::constructSurrogate().
 

Detailed Description

Encapsulates the Tasmanian Sparse Grid module.

Enumeration Type Documentation

◆ TypeIndexRelation

The standard C++11 algorithms std::merge and std::binary_search use bools as return types when comparing entries, thus in the case of many repeated entries, two comparisons have to be performed. Repeated entries are encountered in many sparse grids algorithms, especially when dealing with nested one dimensional rules. Thus, for performance reasons, Tasmanian implements merging and searching methods that compare multi-indexes and produce three outcomes in a single comparison.

Specifically, see SetManipulations::push_merge_map() and MultiIndexSet::getSlot() and MultiIndexSet::addSortedInsexes()

Enumerator
type_abeforeb 

Indicates that the first multi-index in the call to compare comes before the second.

type_bbeforea 

Indicates that the second multi-index in the call to compare after the second.

type_asameb 

Indicates that the two multi-indexes are the same.