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

The master-class that represents an instance of a Tasmanian sparse grid. More...

#include <TasmanianSparseGrid.hpp>

Public Types

using EvaluateCallable = std::function< void(std::vector< double > const &, std::vector< double > &)>
 Signature compatible with TasDREAM::DreamPDF, TasDREAM::DreamModel amd TasDREAM::DreamMergedLikelyModel.
 
using DomainInsideSignature = std::function< bool(std::vector< double > const &)>
 Signature of the domain inside lambda, identical to TasDREAM::DreamDomain.
 

Public Member Functions

 TasmanianSparseGrid ()
 Default constructor, creates and empty grid.
 
 TasmanianSparseGrid (const TasmanianSparseGrid &source)
 Copy constructor, note that the selected acceleration mode is not copied, acceleration is reset to the default.
 
 TasmanianSparseGrid (TasmanianSparseGrid &&source)=default
 Move constructor, the selected acceleration mode is also carried over.
 
 ~TasmanianSparseGrid ()=default
 Destructor, releases all resources.
 
TasmanianSparseGridoperator= (TasmanianSparseGrid const &source)
 Copy assignment, note that the selected acceleration mode is not copied, acceleration is reset to the default.
 
TasmanianSparseGridoperator= (TasmanianSparseGrid &&source)=default
 Move assignment, the selected acceleration mode is also carried over.
 
void write (const char *filename, bool binary=mode_binary) const
 Write the grid to the given filename using either binary or ASCII format.
 
void read (const char *filename)
 Read the grid from the given filename, automatically detect the format.
 
void write (std::ostream &ofs, bool binary=mode_binary) const
 Write the grid to the given stream ofs using either binary or ASCII format.
 
void read (std::istream &ifs, bool binary=mode_binary)
 Read the grid from the given stream ifs using either binary or ASCII format.
 
void write (std::string const &fname, bool binary=mode_binary) const
 Overload that works directly with std::string.
 
void read (std::string const &fname)
 Overload that works directly with std::string.
 
void makeGlobalGrid (int dimensions, int outputs, int depth, TypeDepth type, TypeOneDRule rule, std::vector< int > const &anisotropic_weights, double alpha=0.0, double beta=0.0, const char *custom_filename=nullptr, std::vector< int > const &level_limits=std::vector< int >())
 Make a Global Grid using Lagrange polynomials with support over the entire domain. More...
 
void makeGlobalGrid (int dimensions, int outputs, int depth, TypeDepth type, TypeOneDRule rule, const int *anisotropic_weights=nullptr, double alpha=0.0, double beta=0.0, const char *custom_filename=nullptr, const int *level_limits=nullptr)
 Overload using raw-arrays. More...
 
void makeGlobalGrid (int dimensions, int outputs, int depth, TypeDepth type, CustomTabulated &&crule, std::vector< int > const &anisotropic_weights, std::vector< int > const &level_limits=std::vector< int >())
 Overload make a Global Grid using the provided custom rule. More...
 
void makeGlobalGrid (int dimensions, int outputs, int depth, TypeDepth type, CustomTabulated &&crule, const int *anisotropic_weights=nullptr, const int *level_limits=nullptr)
 Overload using raw-arrays. More...
 
void makeSequenceGrid (int dimensions, int outputs, int depth, TypeDepth type, TypeOneDRule rule, std::vector< int > const &anisotropic_weights, std::vector< int > const &level_limits=std::vector< int >())
 Make Sequence Grid using Newton polynomials with support over the entire domain. More...
 
void makeSequenceGrid (int dimensions, int outputs, int depth, TypeDepth type, TypeOneDRule rule, const int *anisotropic_weights=nullptr, const int *level_limits=nullptr)
 Overload using raw-arrays. More...
 
void makeLocalPolynomialGrid (int dimensions, int outputs, int depth, int order, TypeOneDRule rule, std::vector< int > const &level_limits)
 Make Local Polynomial Grid using piece-wise polynomials with decreasing (compact) support. More...
 
void makeLocalPolynomialGrid (int dimensions, int outputs, int depth, int order=1, TypeOneDRule rule=rule_localp, const int *level_limits=nullptr)
 Overload using raw-arrays. More...
 
void makeWaveletGrid (int dimensions, int outputs, int depth, int order, std::vector< int > const &level_limits)
 Make a Wavelet grid using local hierarchical wavelet basis. More...
 
void makeWaveletGrid (int dimensions, int outputs, int depth, int order=1, const int *level_limits=nullptr)
 Overload using raw-arrays. More...
 
void makeFourierGrid (int dimensions, int outputs, int depth, TypeDepth type, std::vector< int > const &anisotropic_weights, std::vector< int > const &level_limits=std::vector< int >())
 Make a Fourier grid using trigonometric basis with support over the entire domain. More...
 
void makeFourierGrid (int dimensions, int outputs, int depth, TypeDepth type, const int *anisotropic_weights=nullptr, const int *level_limits=nullptr)
 Overload using raw-arrays. More...
 
void copyGrid (const TasmanianSparseGrid *source, int outputs_begin=0, int outputs_end=-1)
 Replace the grid with a copy of the source, does not copy the acceleration options. More...
 
void copyGrid (TasmanianSparseGrid const &source, int outputs_begin=0, int outputs_end=-1)
 Overload using pass-by-reference as opposed to a pointer. More...
 
void updateGlobalGrid (int depth, TypeDepth type, std::vector< int > const &anisotropic_weights, std::vector< int > const &level_limits=std::vector< int >())
 Construct a new grid and merge it with the current one. More...
 
void updateGlobalGrid (int depth, TypeDepth type, const int *anisotropic_weights=nullptr, const int *level_limits=nullptr)
 Overload using raw-arrays. More...
 
void updateSequenceGrid (int depth, TypeDepth type, std::vector< int > const &anisotropic_weights, std::vector< int > const &level_limits=std::vector< int >())
 Construct a new grid and merge it with the current one. More...
 
void updateSequenceGrid (int depth, TypeDepth type, const int *anisotropic_weights=nullptr, const int *level_limits=nullptr)
 Overload using raw-arrays. More...
 
void updateFourierGrid (int depth, TypeDepth type, std::vector< int > const &anisotropic_weights, std::vector< int > const &level_limits=std::vector< int >())
 Construct a new grid and merge it with the current one. More...
 
void updateFourierGrid (int depth, TypeDepth type, const int *anisotropic_weights=nullptr, const int *level_limits=nullptr)
 Overload using raw-arrays. More...
 
void updateGrid (int depth, TypeDepth type, std::vector< int > const &anisotropic_weights, std::vector< int > const &level_limits=std::vector< int >())
 Based on the grid type, calls updateGlobalGrid(), updateSequenceGrid() or updateFourierGrid(). More...
 
void updateGrid (int depth, TypeDepth type, const int *anisotropic_weights=nullptr, const int *level_limits=nullptr)
 Overload using raw-arrays. More...
 
double getAlpha () const
 Return the alpha parameter in the call to makeGlobalGrid(), or return 0 if the grid is not Global.
 
double getBeta () const
 Return the beta parameter in the call to makeGlobalGrid(), or return 0 if the grid is not Global.
 
int getOrder () const
 Return the order parameter in the call to makeLocalPolynomialGrid() or makeWaveletGrid(), or return -1 for all other grid types.
 
int getNumDimensions () const
 Return the dimensions of the grid, i.e., number of model inputs.
 
int getNumOutputs () const
 Return the outputs of the grid, i.e., number of model outputs.
 
TypeOneDRule getRule () const
 Return the underlying TasGrid::TypeOneDRule that gives the points and basis functions.
 
const char * getCustomRuleDescription () const
 Return the character string that is the description of the user-provided tabulated rule. More...
 
int getNumLoaded () const
 Return the number of points already associated with model values via loadNeededValues().
 
int getNumNeeded () const
 Return the number of points that should be provided to the next call of loadNeededValues().
 
int getNumPoints () const
 Returns getNumLoaded() if positive, otherwise returns getNumNeeded(), see getPoints().
 
std::vector< double > getLoadedPoints () const
 Return the points already associated with model values. More...
 
void getLoadedPoints (std::vector< double > &x) const
 Overload that accepts the vector as a parameter. More...
 
void getLoadedPoints (double x[]) const
 Overload using raw-array, writes the loaded points to the first getNumLoaded() times getNumDimensions() entries of the array. More...
 
std::vector< double > getNeededPoints () const
 Return the points that require model values. More...
 
void getNeededPoints (std::vector< double > &x) const
 Overload that accepts the vector as a parameter. More...
 
void getNeededPoints (double *x) const
 Overload using raw-array, writes the loaded points to the first getNumNeeded() times getNumDimensions() entries of the array. More...
 
std::vector< double > getPoints () const
 Returns the loaded points if any, otherwise returns the needed points. More...
 
void getPoints (std::vector< double > &x) const
 Overload that accepts the vector as a parameter. More...
 
void getPoints (double x[]) const
 Overload that accepts the raw array as an input. More...
 
std::vector< double > getQuadratureWeights () const
 Returns a vector of size getNumPoints() of the quadrature weights of the grid. More...
 
void getQuadratureWeights (std::vector< double > &weights) const
 Overload that accepts the vector as a parameter. More...
 
void getQuadratureWeights (double weights[]) const
 Overload that accepts the raw array as an input. More...
 
std::vector< double > getInterpolationWeights (std::vector< double > const &x) const
 Returns the weights of the model outputs that combine to construct the approximation value at x. More...
 
std::vector< double > getInterpolationWeights (double const x[]) const
 Overload that uses raw-array, does not check the array size. More...
 
void getInterpolationWeights (const std::vector< double > &x, std::vector< double > &weights) const
 Overload that uses the vector as a parameter. More...
 
void getInterpolationWeights (const double x[], double weights[]) const
 Overload that uses raw-array, does not check the array size. More...
 
std::vector< double > getDifferentiationWeights (std::vector< double > const &x) const
 Returns the weights of the model outputs that combine to construct the approximate Jacobian matrix (derivative) at x. More...
 
std::vector< double > getDifferentiationWeights (double const x[]) const
 Overload that uses raw-array, does not check the array size. More...
 
void getDifferentiationWeights (const std::vector< double > &x, std::vector< double > &weights) const
 Overload that uses the vector as a parameter. More...
 
void getDifferentiationWeights (const double x[], double weights[]) const
 Overload that uses raw-array, does not check the array size. More...
 
void loadNeededValues (std::vector< double > const &vals)
 Provides the values of the model outputs at the needed points, or overwrites the currently loaded points. More...
 
void loadNeededValues (const double *vals)
 Overload that uses a raw-array, does not check the array size. More...
 
void loadNeededPoints (std::vector< double > const &vals)
 Alias of loadNeededValues().
 
void loadNeededPoints (const double *vals)
 Overload that uses a raw-array, does not check the array size. More...
 
const double * getLoadedValues () const
 Returns the model values that have been loaded in the gird. More...
 
void evaluate (std::vector< double > const &x, std::vector< double > &y) const
 Computes the value of the interpolant (or point-wise approximation) at the given point x. More...
 
void evaluate (const double x[], double y[]) const
 Overload that uses raw-arrays, does not check the array size. More...
 
template<typename FloatType >
void evaluateBatch (std::vector< FloatType > const &x, std::vector< FloatType > &y) const
 Computes the value of the interpolant (or point-wise approximation) for a batch of points. More...
 
void evaluateBatch (const double x[], int num_x, double y[]) const
 Overload that uses raw-arrays. More...
 
void evaluateBatch (const float x[], int num_x, float y[]) const
 Overload using single precision and GPU/CUDA acceleration. More...
 
template<typename FloatType >
void evaluateBatchGPU (const FloatType gpu_x[], int cpu_num_x, FloatType gpu_y[]) const
 Overload that uses GPU raw-arrays. More...
 
template<typename FloatType >
void evaluateFast (const FloatType x[], FloatType y[]) const
 Equivalent to evaluate() with enabled acceleration or evaluateBatch() with a batch of one point. More...
 
template<typename FloatType >
void evaluateFast (std::vector< FloatType > const &x, std::vector< FloatType > &y) const
 Alias to evaluateBatch(). More...
 
void integrate (std::vector< double > &q) const
 Computes the integral of each model output over the sparse grid domain. More...
 
void integrate (double q[]) const
 Overload that uses a raw-array. More...
 
std::vector< double > integrate () const
 Overload that returns a vector.
 
void differentiate (std::vector< double > const &x, std::vector< double > &jacobian) const
 Computes the derivative (if available) of the surrogate model at an input point. More...
 
std::vector< double > differentiate (std::vector< double > const &x) const
 Same as TasmanianSparseGrid::differentiate() but returns the jacobian.
 
void differentiate (const double x[], double jacobian[]) const
 Overload that uses a raw-array. More...
 
bool isGlobal () const
 Returns true if the grid is of type global, false otherwise.
 
bool isSequence () const
 Returns true if the grid is of type sequence, false otherwise.
 
bool isLocalPolynomial () const
 Returns true if the grid is of type local polynomial, false otherwise.
 
bool isWavelet () const
 Returns true if the grid is of type wavelet, false otherwise.
 
bool isFourier () const
 Returns true if the grid is of type Fourier, false otherwise.
 
bool isEmpty () const
 Returns true if the grid is empty (no type), false otherwise.
 
bool empty () const
 Returns true if the grid is empty (no type), false otherwise.
 
void setDomainTransform (std::vector< double > const &a, std::vector< double > const &b)
 Set a linear domain transformation. More...
 
void setDomainTransform (const double a[], const double b[])
 Overload using raw-arrays. More...
 
bool isSetDomainTransfrom () const
 Returns true if a linear domain transformation has been set, false otherwise. More...
 
void clearDomainTransform ()
 Removes the domain transformation. More...
 
void getDomainTransform (std::vector< double > &a, std::vector< double > &b) const
 Returns the two vectors used to call setDomainTransform(). More...
 
void getDomainTransform (double a[], double b[]) const
 Returns the values of the two vectors used to call setDomainTransform(). More...
 
void setConformalTransformASIN (std::vector< int > const &truncation)
 Set conformal transformation using truncated Maclaurin series of the arcsin() function. More...
 
bool isSetConformalTransformASIN () const
 Returns true if conformal transformation has been set.
 
void clearConformalTransform ()
 Removes any currently set transformation.
 
std::vector< int > getConformalTransformASIN () const
 Fills the array with the values of the set transformation.
 
void clearLevelLimits ()
 Removes the currently set level limits. More...
 
std::vector< int > getLevelLimits () const
 Return the currently set level limits. More...
 
void setAnisotropicRefinement (TypeDepth type, int min_growth, int output, const std::vector< int > &level_limits)
 Set refinement using anisotropic estimate for the optimal points. More...
 
void setAnisotropicRefinement (TypeDepth type, int min_growth, int output, const int *level_limits=nullptr)
 Overload using raw-array. More...
 
std::vector< double > getAnisotropicRefinement (TypeDepth type, int min_growth, int output, const std::vector< int > &level_limits)
 Call setAnisotropicRefinement() and then getNeededPoints().
 
std::vector< int > estimateAnisotropicCoefficients (TypeDepth type, int output) const
 Estimate the anisotropic rates of coefficient decay for different direction. More...
 
void estimateAnisotropicCoefficients (TypeDepth type, int output, std::vector< int > &weights) const
 Overload that writes the result to a parameter. More...
 
void setSurplusRefinement (double tolerance, int output, std::vector< int > const &level_limits)
 Refine the grid based on the surplus coefficients, Sequence grids and Global grids with a sequence rule. More...
 
void setSurplusRefinement (double tolerance, int output, const int *level_limits=nullptr)
 Overload that uses array for the level limits.
 
std::vector< double > getSurplusRefinement (double tolerance, int output, std::vector< int > const &level_limits)
 Call setSurplusRefinement() for Sequence and Global grids with a sequence rule and then getNeededPoints().
 
void setSurplusRefinement (double tolerance, TypeRefinement criteria, int output, std::vector< int > const &level_limits, std::vector< double > const &scale_correction=std::vector< double >())
 Refine the grid based on the surplus coefficients, Local-Polynomial and Wavelet grids. More...
 
void setSurplusRefinement (double tolerance, TypeRefinement criteria, int output=-1, const int *level_limits=nullptr, const double *scale_correction=nullptr)
 Overload that uses raw-arrays. More...
 
std::vector< double > getSurplusRefinement (double tolerance, TypeRefinement criteria, int output, std::vector< int > const &level_limits, std::vector< double > const &scale_correction=std::vector< double >())
 Call setSurplusRefinement() for Local-Polynomial and Wavelet grids and then getNeededPoints().
 
void clearRefinement ()
 Remove all needed points from the grid. More...
 
void mergeRefinement ()
 Merges the loaded and needed points into a single grid, resets all loaded values to zero. More...
 
void beginConstruction ()
 Begin a dynamic construction procedure. More...
 
bool isUsingConstruction () const
 Returns true if the dynamic construction procedure has been initialized, false otherwise. More...
 
std::vector< double > getCandidateConstructionPoints (TypeDepth type, std::vector< int > const &anisotropic_weights=std::vector< int >(), std::vector< int > const &level_limits=std::vector< int >())
 Generate a sorted list of points weighted by descending importance. More...
 
std::vector< double > getCandidateConstructionPoints (TypeDepth type, int output, std::vector< int > const &level_limits=std::vector< int >())
 Essentially the same as getCandidateConstructionPoints() but the weights are obtained from a call to estimateAnisotropicCoefficients(). More...
 
std::vector< double > getCandidateConstructionPoints (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 >())
 Returns a sorted list of points weighted by descending importance using the hierarchical surpluses. More...
 
void loadConstructedPoints (std::vector< double > const &x, std::vector< double > const &y)
 Add pairs of points with associated model values. More...
 
void loadConstructedPoints (const double x[], int numx, const double y[])
 Same as loadConstructedPoint() but using arrays in place of vectors (array size is not checked) More...
 
void finishConstruction ()
 End the procedure, clears flags and unused constructed points, can go back to using regular refinement. More...
 
const double * getHierarchicalCoefficients () const
 Return a reference to the internal data-structure that stores the hierarchical coefficients. More...
 
void getHierarchicalCoefficientsStatic (double *coeff) const
 
void setHierarchicalCoefficients (const std::vector< double > &c)
 Overwrites the current set of coefficients (and loaded values) with the ones provided. More...
 
void setHierarchicalCoefficients (const double c[])
 Overload that uses raw-arrays. More...
 
void evaluateHierarchicalFunctions (std::vector< double > const &x, std::vector< double > &y) const
 Computes the values of the hierarchical function basis at the specified points. More...
 
std::vector< double > evaluateHierarchicalFunctions (std::vector< double > const &x) const
 Overload that returns the result. More...
 
void evaluateHierarchicalFunctions (const double x[], int num_x, double y[]) const
 Array overload, the inputs must have pre-allocated and correct size. More...
 
void evaluateSparseHierarchicalFunctions (const std::vector< double > &x, std::vector< int > &pntr, std::vector< int > &indx, std::vector< double > &vals) const
 Constructs a sparse matrix with the values of the hierarchical basis functions. More...
 
std::vector< double > getHierarchicalSupport () const
 Returns the support of the hierarchical functions. More...
 
std::vector< double > integrateHierarchicalFunctions () const
 Returns the integrals of the hierarchical basis functions. More...
 
void integrateHierarchicalFunctions (std::vector< double > &integrals) const
 Overload where the vector is passes as a parameter. More...
 
void integrateHierarchicalFunctions (double integrals[]) const
 Overload using raw-arrays. More...
 
std::vector< int > getGlobalPolynomialSpace (bool interpolation) const
 Returns the powers of the polynomial that span the effective basis, Global and Sequence grids only. More...
 
void printStats (std::ostream &os=std::cout) const
 Prints short human-readable text describing the grid properties. More...
 
void enableAcceleration (TypeAcceleration acc)
 Change the current acceleration mode to the one specified. More...
 
void enableAcceleration (TypeAcceleration acc, int new_gpu_id)
 Combines calls the enableAcceleration(), getGPUID() and allows for user provided handles. More...
 
void favorSparseAcceleration (bool favor)
 Set the preferred back-end algorithm for Local Polynomial grids. More...
 
void setCuBlasHandle (void *handle)
 Takes a user provided cuBlas handle. More...
 
void setCuSparseHandle (void *handle)
 Takes a user provided cuSparse handle. More...
 
void setCuSolverHandle (void *handle)
 Takes a user provided cuSparse handle. More...
 
void setRocBlasHandle (void *handle)
 Takes a user provided cuBlas handle. More...
 
void setRocSparseHandle (void *handle)
 Takes a user provided cuSparse handle. More...
 
void setSycleQueue (void *queue)
 Takes a user provided sycl::queue handle. More...
 
TypeAcceleration getAccelerationType () const
 Returns the current effective acceleration mode. More...
 
void setGPUID (int new_gpu_id)
 Select the current CUDA device. More...
 
int getGPUID () const
 Returns the currently set CUDA device. More...
 
template<typename FloatType >
void evaluateHierarchicalFunctionsGPU (const FloatType gpu_x[], int cpu_num_x, FloatType gpu_y[]) const
 Computes the values of the hierarchical function basis at the specified points (CUDA version). More...
 
template<typename FloatType >
void evaluateSparseHierarchicalFunctionsGPU (const FloatType gpu_x[], int cpu_num_x, int *&gpu_pntr, int *&gpu_indx, FloatType *&gpu_vals, int &num_nz) const
 Computes the values of the hierarchical function basis at the specified points (sparse/CUDA version). More...
 
 operator EvaluateCallable () const
 Custom conversion to a callable method using the TasDREAM::DreamPDF signature. More...
 
DomainInsideSignature getDomainInside () const
 Returns a lambda object that satisfies the TasDREAM::DreamDomain signature. More...
 
void removePointsByHierarchicalCoefficient (double tolerance, int output=-1, const double *scale_correction=nullptr)
 Removes all points from the grid that have relative surplus less than the tolerance. More...
 
void removePointsByHierarchicalCoefficient (int num_new_points, int output=-1, const double *scale_correction=nullptr)
 Keeps only the given number of points with largest scaled surpluses. More...
 

Static Public Member Functions

static const char * getVersion ()
 Return a hard-coded character string with the version in format "Major.Minor".
 
static int getVersionMajor ()
 Return the library major version.
 
static int getVersionMinor ()
 Return the library minor version.
 
static const char * getLicense ()
 Return a hard-coded character string with a brief statement of the license.
 
static const char * getGitCommitHash ()
 Return the git hash string, will use a placeholder if the git command was not found on compile time or if building from an official release.
 
static const char * getCmakeCxxFlags ()
 Return the CMAKE_BUILD_TYPE and CMAKE_CXX_FLAGS used in the configuration.
 
static bool isOpenMPEnabled ()
 Returns true if compiled with OpenMP support, e.g., Tasmanian_ENABLE_OPENMP=ON.
 
static bool isCudaEnabled ()
 Returns true if compiled with CUDA support, e.g., Tasmanian_ENABLE_CUDA=ON.
 
static bool isHipEnabled ()
 Returns true if compiled with HIP support, e.g., Tasmanian_ENABLE_HIP=ON.
 
static bool isDpcppEnabled ()
 Returns true if compiled with DPC++ support, e.g., Tasmanian_ENABLE_DPCPP=ON.
 
static bool isAccelerationAvailable (TypeAcceleration acc)
 Returns whether a specific mode can be enabled. More...
 
static int getNumGPUs ()
 Return the number of visible CUDA devices. More...
 
static int getGPUMemory (int gpu)
 Return the available device memory, in units of MB. More...
 
static std::string getGPUName (int gpu)
 Return the CUDA device name. More...
 

Detailed Description

The master-class that represents an instance of a Tasmanian sparse grid.

class TasmanianSparseGrid
An object of this class represents a sparse grid, all aspects of the grid can be accessed through member methods. Each object operates independently to the extend allowed by third-party libraries, e.g., making concurrent calls to OpenBLAS will result in a race condition unless OpenBLAS is compiled with a special flag. The available member methods are split into groups that deal with different aspects of the grid. The object itself wraps around one of five different types of grids, Global, Sequence, Local Polynomial, Wavelet, and Fourier, each using different back-end data structures. The class itself handles the domain transformations, e.g., the translation between the canonical [-1, 1] and a user provided [a, b], as well as mechanisms to select acceleration back-end, common meta I/O, and error checking.
Error Handling
TasmanianSparseGrid throws std::invalid_argument and std::runtime_error exceptions with simple human-readable messages. The std::invalid_argument exception indicates method call arguments of insufficient or incompatible data, e.g., vectors have insufficient size or negative number was given in place of a positive one. The std::runtime_error indicates problems with missing files, incorrect file format, or out-of-order function calls, e.g., adaptive refinement strategies require model values in order to infer the optimal structure of the grid, therefore calling a refinement method before model values are provided will result in a std::runtime_error. Similarly, std::runtime_error will be thrown if a method is called with the wrong grid type, e.g., global anisotropic coefficients cannot be computed for local polynomial grids.

If CUDA acceleration is used, CUDA error may be encountered, e.g., incorrect driver initialization, out-of-memory, etc. Tasmanian will translate CUDA error codes and messages into std::runtime_error exceptions and such exceptions can be raised by any method that uses acceleration.

Constructors and Copy/Move assignment
TasGrid::TasmanianSparseGrid includes move and copy constructors and assignment operator= overloads. There is also a default constructor that creates an empty grid. Note that move will preserve the enabled acceleration mode and the internal cache data-structures, while copy will reset the acceleration to the default.
Make Grid Methods
Each of the five grid types comes with a separate make command. The first three parameters are always the number of inputs (dimensions), the number of outputs, and the initial depth or level of the grid. The remaining parameters are specific to the type of grid.
Domain Transformations
In most cases, sparse grids are constructed on a canonical domain [-1, 1]. The exceptions are the Fourier grids with domain [0, 1], and the Gauss-Hermite and Gauss-Laguerre rules with unbounded domain, see TasGrid::rule_gausshermite and TasGrid::rule_gausslaguerre for details. Tasmanian can automatically apply a linear domain transformations shifting the canonical interval to an arbitrary [a, b] or a non-linear (conformal) transformation that may accelerate convergence for some models. Each dimension of the grid can have a different linear and/or non-linear transform.
Stream and File I/O
Sparse grid objects can be written to and from files and std::ostream objects. Tasmanian supports both binary and ASCII formats; the binary format has lower overhead in both file size and operations required to read/write, the ASCII format is portable without considerations of endians and for small grids the files are human-readable which helps debugging.
Get Points and Weights
The points on the grid are labeled as "loaded" or "needed" based on whether the associated model values are already provided by the user. The grid points are also associated with quadrature, interpolation, and differentiation weights, where the computed weights correspond to the loaded points unless all points are labeled as needed. If the grid is set to work with zero model outputs, then the distinction is meaningless and the points will be accessed through the generic getPoints() methods.
Load Model Values or Hierarchical Coefficients
In order to construct an interpolant, Tasmanian needed to compute the coefficients of the basis functions. The user can directly provide the coefficients, e.g., computed externally with a regression method using an unstructured set of samples, or the user can provide the model values at the needed points and let Tasmanian do the computations.
Update and Adaptive Refinement
The most efficient approximation schemes adapt the grid (and the associated basis functions) to the target model. Tasmanian provides several adaptive procedures tailored to different types of grids. The refinement requires two-way communication between the model and Tasmanian, i.e., Tasmanian provides an initial set of points, the model provides the values, then Tasmanian provides an updated set of points, and so on. This can be done either in batches of point or in dynamic construction setup where values can be given one at a time in an arbitrary order (see next paragraph). See also the papers referenced in TasGrid::TypeDepth and TasGrid::TypeRefinement.
Dynamic Construction
The standard refinement procedure works in batches, where all model values for the batch have to be computed externally and loaded together in the correct order. The dynamic construction alleviates those restrictions, but comes at an increase cost since Tasmanian has to do more work and store extra data. The dynamic construction (and the internal data-structures) are initiated with the beginConstruction() command and batch refinement cannot be used until the procedure is concluded. See also TasGrid::constructSurrogate() that can take a user defined model and fully automate the construction process.
Using Unstructured Data
The standard sparse grid methods assume that model data is available at the very specific grid points that are chosen according to optimal estimates. However, in some cases, the model inputs cannot be controlled precisely and only randomly samples model data is available. Sparse grids can still produce accurate surrogates using such unstructured data by effectively removing the points and working only with the underlying basis (the basis is still optimal for the corresponding class of functions). The hierarchical basis method allow direct access to the values of the basis functions and the associated coefficients.
Evaluate Methods
The evaluate methods compute the sparse grid approximation to the model for arbitrary point within the domain, i.e., for points that don't necessarily align with the original grid. The batch and fast evaluate method can leverage accelerated linear algebra libraries such as BLAS, cuBLAS and MAGAMA. See the documentation for TasGrid::TypeAcceleration for details.
Acceleration Back-end Selection
Allows specifying the acceleration used for evaluation methods and (in some cases) for computing the basis coefficients during loadNeededValues(). For example, methods are provided to check the number of available CUDA devices and to select an active device on a multi-GPU system.
Set User Provided Handles for Accelerated Linear Algebra
Acceleration frameworks such as CUDA ROCm and oneAPI use handles and queues that are either opaque pointers to internal data-structures (CUDA and ROCm) or queues that must be used for all types of acceleration and data operations. A Tasmanian object will allocate internal handles as needed, but the user can specify the handles manually. The handles must be set after a GPU capable acceleration has been enabled via enableAcceleration(), Tasmanian takes a non-owning reference and the user is responsible for deleting the handle after the Tasmanian object has been destroyed or a non-GPU acceleration is selected, i.e., accel_none or accel_cpu_blas. The handles are accepted as void-pointers (the type is stripped) to keep consistent API since only one GPU framework can be enabled in a single Tasmanian build and the headers for the other frameworks will not be present.
Get Grid Meta-data
Various method that read the number of points, grid type, specifics about the rule and domain transformations and others. Some of the output is in human-readable format for debugging and sanity check.
Get Library Meta-data
Runtime API is provided for the library version and various compile time options. Those methods are static and are included in the class API solely for consistency reasons.

Member Function Documentation

◆ makeGlobalGrid() [1/4]

void TasGrid::TasmanianSparseGrid::makeGlobalGrid ( int  dimensions,
int  outputs,
int  depth,
TypeDepth  type,
TypeOneDRule  rule,
std::vector< int > const &  anisotropic_weights,
double  alpha = 0.0,
double  beta = 0.0,
const char *  custom_filename = nullptr,
std::vector< int > const &  level_limits = std::vector< int >() 
)

Make a Global Grid using Lagrange polynomials with support over the entire domain.

Construct a sparse grid with the given set of parameters.

Parameters
dimensionsis a positive integer specifying the dimension or number of model inputs of the grid. There is no hard restriction on how big the dimension can be; however, for large dimensions, the number of points of the sparse grid grows fast and hence the grid may require a prohibitive amount of memory.
outputsis a non-negative integer specifying the number of outputs for the model. If outputs is zero, then the grid can only generate quadrature and interpolation weights. There is no hard restriction on how many outputs can be handled; however, Tasmanian requires at least outputs times number of points in storage and the computational complexity of evaluate and I/O methods increases linearly with the number of outputs.
depthis a non-negative integer that controls the density of grid points. This is commonly referred to the "level" of the grid and corresponds to the L parameter in the formulas for TasGrid::TypeDepth. There is no hard restriction on how big depth can be; however, the number of points have a strong influence on performance and memory requirements.
typeindicates the tensor selection strategy, see TasGrid::TypeDepth.
ruleis one of the global rules, see TasGrid::TypeOneDRule.
anisotropic_weightsindicates the relative "importance" of the inputs. If an empty vector is provided, the isotropic selection is used, i.e., assuming all inputs are equally important. If non-empty vector is given, then the size must match the number of required weights based on the selection type, i.e., the curved types use 2 x dimensions number of weights while the non-curved types use only dimensions. Integer values are used, but only the relative scale affects the selection, e.g., {2, 1}, {4, 2}, {20, 10} and {200, 100} are equivalent. The first dimensions entries correspond to the $ \xi_1, \xi_2, \cdots, \xi_d $ weights and the second set corrections to the $ \eta_1, \eta_2, \cdots, \eta_d $ weights in the formulas in TasGrid::TypeDepth.
alphaspecifies the $ \alpha $ parameter for the integration weight $ \rho(x) $, ignored when rule is doesn't have such parameter. See TasGrid::TypeOneDRule.
betaspecifies the $ \beta $ parameter for the integration weight $ \rho(x) $, ignored when rule is doesn't have such parameter. See TasGrid::TypeOneDRule.
custom_filenamespecifies the file containing the custom rule description, used only when rule is TasGrid::rule_customtabulated.
level_limitsis either empty or has size dimensions indicating the restriction of the levels for each direction which in turn restricts the number of points. For example, limit 1 when using TasGrid::rule_clenshawcurtis will restricts the grid to the mid-point and end-points of the domain. Each dimension can have a different restriction and negative numbers indicate no restriction for that direction.
Exceptions
std::invalid_argumentwith human readable messages when integers are out of range, the rule is not a global rule, or the vectors have incorrect size.
std::runtime_errorif the custom_filename if missing, or it cannot be opened, or the format is incorrect.

◆ makeGlobalGrid() [2/4]

void TasGrid::TasmanianSparseGrid::makeGlobalGrid ( int  dimensions,
int  outputs,
int  depth,
TypeDepth  type,
TypeOneDRule  rule,
const int *  anisotropic_weights = nullptr,
double  alpha = 0.0,
double  beta = 0.0,
const char *  custom_filename = nullptr,
const int *  level_limits = nullptr 
)

Overload using raw-arrays.

The inputs are the same as makeGlobalGrid() except the vectors are replaced by arrays. Passing nullptr for an array is equivalent to passing an empty vector. Throws the same exceptions, but it does not check if the arrays have the correct size.

◆ makeGlobalGrid() [3/4]

void TasGrid::TasmanianSparseGrid::makeGlobalGrid ( int  dimensions,
int  outputs,
int  depth,
TypeDepth  type,
CustomTabulated &&  crule,
std::vector< int > const &  anisotropic_weights,
std::vector< int > const &  level_limits = std::vector< int >() 
)

Overload make a Global Grid using the provided custom rule.

Compared to makeGlobalGrid(), this uses the provided CustomTabulated rule with rule_customtabulated.

◆ makeGlobalGrid() [4/4]

void TasGrid::TasmanianSparseGrid::makeGlobalGrid ( int  dimensions,
int  outputs,
int  depth,
TypeDepth  type,
CustomTabulated &&  crule,
const int *  anisotropic_weights = nullptr,
const int *  level_limits = nullptr 
)

Overload using raw-arrays.

Same as the other CustomTabulated overload but uses raw-arrays.

◆ makeSequenceGrid() [1/2]

void TasGrid::TasmanianSparseGrid::makeSequenceGrid ( int  dimensions,
int  outputs,
int  depth,
TypeDepth  type,
TypeOneDRule  rule,
std::vector< int > const &  anisotropic_weights,
std::vector< int > const &  level_limits = std::vector< int >() 
)

Make Sequence Grid using Newton polynomials with support over the entire domain.

Mathematically the Sequence and Global grids do not differ in properties; however, the implementation of the Sequence grids uses optimized internal data structures which leads to massive increase in speed when calling evaluate methods at the expense of doubled memory size and increased cost of loadNeededValues(). The inputs are identical to makeGlobalGrid() with the restriction that rule must be one of:

@ rule_minlebesgue
A greedy sequence rule with nodes added to minimize the Lebesgue constant.
Definition: tsgEnumerates.hpp:321
@ rule_leja
Classic sequence rule, moderate Lebesgue constant growth (empirical result only).
Definition: tsgEnumerates.hpp:299
@ rule_mindelta
A greedy sequence rule with nodes added to minimize the norm of the surplus operator.
Definition: tsgEnumerates.hpp:325
@ rule_maxlebesgue
A greedy sequence rule with nodes placed at the maximum of the Lebesgue function.
Definition: tsgEnumerates.hpp:317
@ rule_rlejashifted
Similar sequence to rule_rleja but with nodes strictly in the interior.
Definition: tsgEnumerates.hpp:311
@ rule_rleja
Classic sequence rule based on complex analysis, moderate Lebesgue constant growth (theoretically pro...
Definition: tsgEnumerates.hpp:303
Exceptions
std::invalid_argumentwith human readable messages when integers are out of range, the rule is not a sequence rule, or the vectors have incorrect size.

◆ makeSequenceGrid() [2/2]

void TasGrid::TasmanianSparseGrid::makeSequenceGrid ( int  dimensions,
int  outputs,
int  depth,
TypeDepth  type,
TypeOneDRule  rule,
const int *  anisotropic_weights = nullptr,
const int *  level_limits = nullptr 
)

Overload using raw-arrays.

The inputs are the same as makeSequenceGrid() except the vectors are replaced by arrays. Passing nullptr for an array is equivalent to passing an empty vector. Throws the same exceptions, but it does not check if the arrays have the correct size.

◆ makeLocalPolynomialGrid() [1/2]

void TasGrid::TasmanianSparseGrid::makeLocalPolynomialGrid ( int  dimensions,
int  outputs,
int  depth,
int  order,
TypeOneDRule  rule,
std::vector< int > const &  level_limits 
)

Make Local Polynomial Grid using piece-wise polynomials with decreasing (compact) support.

Creates a grid based on one of the local hierarchical piece-wise polynomial rules. Local grids can be used for integration, but in many cases the quadrature weights will be zero for some points.

Parameters
dimensionsis a positive integer specifying the dimension or number of model inputs of the grid. See also makeGlobalGrid().
outputsis a non-negative integer specifying the number of outputs for the model. Unlike makeGlobalGrid() the storage requirement is at least twice the number of points times the number of outputs.
depthis a non-negative integer that controls the density of grid points, i.e., the "level". The initial construction of the local grids uses tensor selection equivalent to TasGrid::type_level and depth is the L parameter in the formula.
orderis an integer no smaller than -1 indicating the largest polynomial order of the basis functions, i.e., 1 indicates the use of constant and linear functions only, while 2 would allow quadratics (if enough points are present). Using -1 indicates using the largest possible order for each point. See the papers referenced in the TasGrid::TypeRefinement description.
ruleis one of the local polynomial rules, e.g.,
@ rule_localp
Nested rule with a hierarchy of uniformly distributed nodes and functions with compact support.
Definition: tsgEnumerates.hpp:362
@ rule_localpb
Variation of rule_localp focusing nodes on the boundary instead of the interior.
Definition: tsgEnumerates.hpp:368
@ rule_localp0
Variation of rule_localp assuming the model is zero at the domain boundary.
Definition: tsgEnumerates.hpp:364
@ rule_semilocalp
Variation of rule_localp using increased support in exchange for higher order basis (better for smoot...
Definition: tsgEnumerates.hpp:366
Note that using order 0 with any rule will automatically switch to a special piece-wise constant hierarchy.
level_limitsis identical to that for makeGlobalGrid().
Exceptions
std::invalid_argumentwith human readable messages when integers are out of range, the rule is not a local polynomial rule, or the vector has incorrect size.

◆ makeLocalPolynomialGrid() [2/2]

void TasGrid::TasmanianSparseGrid::makeLocalPolynomialGrid ( int  dimensions,
int  outputs,
int  depth,
int  order = 1,
TypeOneDRule  rule = rule_localp,
const int *  level_limits = nullptr 
)

Overload using raw-arrays.

The inputs are the same as makeLocalPolynomialGrid() except the vector is replaced by an array. Passing nullptr for an array is equivalent to passing an empty vector. Throws the same exceptions, but this does not check if the array has the correct size.

◆ makeWaveletGrid() [1/2]

void TasGrid::TasmanianSparseGrid::makeWaveletGrid ( int  dimensions,
int  outputs,
int  depth,
int  order,
std::vector< int > const &  level_limits 
)

Make a Wavelet grid using local hierarchical wavelet basis.

Wavelets are specialized functions that form a Riesz basis and can offer better convergence in an adaptive refinement process. See the pdf manual for more details.

Parameters
dimensionsis a positive integer specifying the dimension or number of model inputs of the grid. See also makeGlobalGrid().
outputsis a non-negative integer specifying the number of outputs for the model. Unlike makeGlobalGrid() the storage requirement is at least twice the number of points times the number of outputs.
depthis a non-negative integer that controls the density of grid points, i.e., the "level". The initial construction of the local grids uses tensor selection equivalent to TasGrid::type_level and depth is the L parameter in the formula, but unlike the local polynomial hierarchies the zero-th level has 3 or 5 points and hence the grid density is overall higher.
orderis the wavelet approximation order, implemented orders are only 1 and 3 (wavelet order cannot be an even number).
level_limitsis identical to that for makeGlobalGrid().
Exceptions
std::invalid_argumentwith human readable messages when integers are out of range, or the vector has incorrect size.

◆ makeWaveletGrid() [2/2]

void TasGrid::TasmanianSparseGrid::makeWaveletGrid ( int  dimensions,
int  outputs,
int  depth,
int  order = 1,
const int *  level_limits = nullptr 
)

Overload using raw-arrays.

The inputs are the same as makeWaveletGrid() except the vector is replaced by an array. Passing nullptr for an array is equivalent to passing an empty vector. Throws the same exceptions, but this does not check if the array has the correct size.

◆ makeFourierGrid() [1/2]

void TasGrid::TasmanianSparseGrid::makeFourierGrid ( int  dimensions,
int  outputs,
int  depth,
TypeDepth  type,
std::vector< int > const &  anisotropic_weights,
std::vector< int > const &  level_limits = std::vector< int >() 
)

Make a Fourier grid using trigonometric basis with support over the entire domain.

The trigonometric basis guarantees that the interpolant is periodic (in all derivatives) across the domain boundary. This is the only grid that defaults to a canonical interval [0, 1] (as opposed to [-1, 1]). The parameters and error messages are identical to those used in makeGlobalGrid(), but this always uses TasGrid::rule_fourier and the polynomial order in the selection is replaced by the frequency of the trigonometric basis.

◆ makeFourierGrid() [2/2]

void TasGrid::TasmanianSparseGrid::makeFourierGrid ( int  dimensions,
int  outputs,
int  depth,
TypeDepth  type,
const int *  anisotropic_weights = nullptr,
const int *  level_limits = nullptr 
)

Overload using raw-arrays.

The inputs are the same as makeFourierGrid() except the vectors are replaced by arrays. Passing nullptr for an array is equivalent to passing an empty vector. Throws the same exceptions, but it does not check if the arrays have the correct size.

◆ copyGrid() [1/2]

void TasGrid::TasmanianSparseGrid::copyGrid ( const TasmanianSparseGrid source,
int  outputs_begin = 0,
int  outputs_end = -1 
)

Replace the grid with a copy of the source, does not copy the acceleration options.

Using the default inputs is equivalent to the assignment operator, otherwise this allows to copy all points and inputs of the grid but only a sub-range of the outputs.

Parameters
sourcethe grid to copy from.
outputs_beginthe first output of the new grid.
outputs_endis the first not-included output, i.e., the logic is similar to std::copy(). If the last parameter is outside of the range, then include all outputs from output_begin till the end of the outputs.

◆ copyGrid() [2/2]

void TasGrid::TasmanianSparseGrid::copyGrid ( TasmanianSparseGrid const &  source,
int  outputs_begin = 0,
int  outputs_end = -1 
)
inline

Overload using pass-by-reference as opposed to a pointer.

The pointer version is supported for backwards compatibility, this is the preferred way to use the copy command.

◆ updateGlobalGrid() [1/2]

void TasGrid::TasmanianSparseGrid::updateGlobalGrid ( int  depth,
TypeDepth  type,
std::vector< int > const &  anisotropic_weights,
std::vector< int > const &  level_limits = std::vector< int >() 
)

Construct a new grid and merge it with the current one.

This method is used for refinement with user specified anisotropic weights. As such, it can be called only for global grids with a nested rule. If there are no loaded points (or if the number of outputs is zero) then this is equivalent to calling makeGlobalGrid() (i.e., completely replacing the grid) using the current dimensions, outputs, and rule parameters and the new values for depth, type, anisotropic coefficients and level limits. If there are loaded model values, then the new grid will be added to the existing one without removing any existing points.

The parameters used and the thrown exceptions are identical to those in the call to makeGlobalGrid(). In addition, std::runtime_error is thrown if the current grid is not Global.

◆ updateGlobalGrid() [2/2]

void TasGrid::TasmanianSparseGrid::updateGlobalGrid ( int  depth,
TypeDepth  type,
const int *  anisotropic_weights = nullptr,
const int *  level_limits = nullptr 
)

Overload using raw-arrays.

Array dimensions are not checked, otherwise identical to updateGlobalGrid().

◆ updateSequenceGrid() [1/2]

void TasGrid::TasmanianSparseGrid::updateSequenceGrid ( int  depth,
TypeDepth  type,
std::vector< int > const &  anisotropic_weights,
std::vector< int > const &  level_limits = std::vector< int >() 
)

Construct a new grid and merge it with the current one.

This method is used for refinement with user specified anisotropic weights. If there are no loaded points (or if the number of outputs is zero) then this is equivalent to calling makeSequenceGrid() (i.e., completely replacing the grid) using the current dimensions, outputs, and rule parameters and the new values for depth, type, anisotropic coefficients and level limits. If there are loaded model values, then the new grid will be added to the existing one without removing any existing points.

The parameters used and the thrown exceptions are identical to those in the call to makeSequenceGrid(). In addition, std::runtime_error is thrown if the current grid is not Sequence.

◆ updateSequenceGrid() [2/2]

void TasGrid::TasmanianSparseGrid::updateSequenceGrid ( int  depth,
TypeDepth  type,
const int *  anisotropic_weights = nullptr,
const int *  level_limits = nullptr 
)

Overload using raw-arrays.

Array dimensions are not checked, otherwise identical to updateSequenceGrid().

◆ updateFourierGrid() [1/2]

void TasGrid::TasmanianSparseGrid::updateFourierGrid ( int  depth,
TypeDepth  type,
std::vector< int > const &  anisotropic_weights,
std::vector< int > const &  level_limits = std::vector< int >() 
)

Construct a new grid and merge it with the current one.

This method is used for refinement with user specified anisotropic weights. If there are no loaded points (or if the number of outputs is zero) then this is equivalent to calling makeFourierGrid() (i.e., completely replacing the grid) using the current dimensions, outputs, and rule parameters and the new values for depth, type, anisotropic coefficients and level limits. If there are loaded model values, then the new grid will be added to the existing one without removing any existing points.

The parameters used and the thrown exceptions are identical to those in the call to makeFourierGrid(). In addition, std::runtime_error is thrown if the current grid is not Fourier.

◆ updateFourierGrid() [2/2]

void TasGrid::TasmanianSparseGrid::updateFourierGrid ( int  depth,
TypeDepth  type,
const int *  anisotropic_weights = nullptr,
const int *  level_limits = nullptr 
)

Overload using raw-arrays.

Array dimensions are not checked, otherwise identical to updateFourierGrid().

◆ updateGrid() [1/2]

void TasGrid::TasmanianSparseGrid::updateGrid ( int  depth,
TypeDepth  type,
std::vector< int > const &  anisotropic_weights,
std::vector< int > const &  level_limits = std::vector< int >() 
)

Based on the grid type, calls updateGlobalGrid(), updateSequenceGrid() or updateFourierGrid().

Exceptions
std::runtime_errorif the grid is not Global, Sequence or Fourier.

◆ updateGrid() [2/2]

void TasGrid::TasmanianSparseGrid::updateGrid ( int  depth,
TypeDepth  type,
const int *  anisotropic_weights = nullptr,
const int *  level_limits = nullptr 
)

Overload using raw-arrays.

Array dimensions are not checked, otherwise identical to updateFourierGrid().

◆ getCustomRuleDescription()

const char* TasGrid::TasmanianSparseGrid::getCustomRuleDescription ( ) const

Return the character string that is the description of the user-provided tabulated rule.

User-provided rules, i.e., Global grids with TasGrid::rule_customtabulated have a description string that is returned by this method. All other grids and rules will return an empty string.

◆ getLoadedPoints() [1/3]

std::vector<double> TasGrid::TasmanianSparseGrid::getLoadedPoints ( ) const
inline

Return the points already associated with model values.

Returns a vector of size getNumLoaded() times getNumDimensions() with the coordinates of the loaded points. The vector is logically divided into strips of length getNumDimensions(), each strip corresponds to a single point.

◆ getLoadedPoints() [2/3]

void TasGrid::TasmanianSparseGrid::getLoadedPoints ( std::vector< double > &  x) const
inline

Overload that accepts the vector as a parameter.

The vector is resized to getNumLoaded() times getNumDimensions() and overwritten with the loaded points. See getLoadedPoints().

◆ getLoadedPoints() [3/3]

void TasGrid::TasmanianSparseGrid::getLoadedPoints ( double  x[]) const

Overload using raw-array, writes the loaded points to the first getNumLoaded() times getNumDimensions() entries of the array.

Assumes the array size if at least getNumLoaded() times getNumDimensions(), overwrites the entries with the loaded points in the same format as in getLoadedPoints().

◆ getNeededPoints() [1/3]

std::vector<double> TasGrid::TasmanianSparseGrid::getNeededPoints ( ) const
inline

Return the points that require model values.

Returns a vector of size getNumNeeded() times getNumDimensions() with the coordinates of the needed points. The vector is logically divided into strips of length getNumDimensions(), each strip corresponds to a single point.

◆ getNeededPoints() [2/3]

void TasGrid::TasmanianSparseGrid::getNeededPoints ( std::vector< double > &  x) const
inline

Overload that accepts the vector as a parameter.

The vector is resized to getNumNeeded() times getNumDimensions() and overwritten with the needed points. See getNeededPoints().

◆ getNeededPoints() [3/3]

void TasGrid::TasmanianSparseGrid::getNeededPoints ( double *  x) const

Overload using raw-array, writes the loaded points to the first getNumNeeded() times getNumDimensions() entries of the array.

Assumes the array size if at least getNumNeeded() times getNumDimensions(), overwrites the entries with the loaded points in the same format as in getNeededPoints().

◆ getPoints() [1/3]

std::vector<double> TasGrid::TasmanianSparseGrid::getPoints ( ) const
inline

Returns the loaded points if any, otherwise returns the needed points.

Grid operations that do not require model values, e.g., getQuadratureWeights(), operate with the loaded points, but in some cases all points may be needed, e.g., right after a make method is called. Thus, points is an alias to the loaded points unless there are no loaded points, in which case points aliases to the needed points. The overloads of this method have the same behavior as the corresponding overload of the other get points methods.

Note that if getNumPoints() is zero then the grid is empty().

◆ getPoints() [2/3]

void TasGrid::TasmanianSparseGrid::getPoints ( std::vector< double > &  x) const
inline

Overload that accepts the vector as a parameter.

See getNeededPoints().

◆ getPoints() [3/3]

void TasGrid::TasmanianSparseGrid::getPoints ( double  x[]) const

Overload that accepts the raw array as an input.

See getNeededPoints().

◆ getQuadratureWeights() [1/3]

std::vector<double> TasGrid::TasmanianSparseGrid::getQuadratureWeights ( ) const
inline

Returns a vector of size getNumPoints() of the quadrature weights of the grid.

The quadrature is designed to work with weight of constant 1 unless the grid rule is associated with a special weight. See TasGrid::TypeOneDRule for details.

Returns
A vector of size getNumPoints() holding the quadrature weights; the order of the weights matches the getPoints().

◆ getQuadratureWeights() [2/3]

void TasGrid::TasmanianSparseGrid::getQuadratureWeights ( std::vector< double > &  weights) const
inline

Overload that accepts the vector as a parameter.

See getQuadratureWeights().

◆ getQuadratureWeights() [3/3]

void TasGrid::TasmanianSparseGrid::getQuadratureWeights ( double  weights[]) const

Overload that accepts the raw array as an input.

See getQuadratureWeights().

◆ getInterpolationWeights() [1/4]

std::vector<double> TasGrid::TasmanianSparseGrid::getInterpolationWeights ( std::vector< double > const &  x) const

Returns the weights of the model outputs that combine to construct the approximation value at x.

The sum of the model values times the weights will produce the approximation at x. For problems where the model outputs can be represented by a vector, it is recommended to use loadNeededValues() and evaluate() methods which have much better performance. However, not all models can be easily represented as vector valued functions, e.g., the discretization of the operators in a parametrized partial differential equation can result in sparse matrices with very different fill. Therefore, Tasmanian offers the option to compute these weights and leave to the user to compute the corresponding weighted sum, e.g., the matrix for each grid point is stored independently and the action of the parametrized operator onto the vector is approximated as a weighed linear combination of the individual matrix vector results.

Parameters
xis a vector of size getNumDimensions() with the coordinates of the point of interest in the transformed domain.
Returns
A vector of size getNumPoints() with the interpolation weights, the order of the weights matches the order of the getPoints().

Note that using a vector x outside of the domain will result in undefined behavior, but will not throw an exception.

Exceptions
std::runtime_errorif x has an incorrect size.

◆ getInterpolationWeights() [2/4]

std::vector<double> TasGrid::TasmanianSparseGrid::getInterpolationWeights ( double const  x[]) const
inline

Overload that uses raw-array, does not check the array size.

Identical to getInterpolationWeights() but does not throw if x is larger than getNumDimensions(); however, using shorter x is undefined behavior and will likely segfault.

◆ getInterpolationWeights() [3/4]

void TasGrid::TasmanianSparseGrid::getInterpolationWeights ( const std::vector< double > &  x,
std::vector< double > &  weights 
) const

Overload that uses the vector as a parameter.

Identical to getInterpolationWeights() but the weights vector will be resized to size getNumPoints() and the entries will be overwritten with the output of getInterpolationWeights().

◆ getInterpolationWeights() [4/4]

void TasGrid::TasmanianSparseGrid::getInterpolationWeights ( const double  x[],
double  weights[] 
) const

Overload that uses raw-array, does not check the array size.

Identical to getInterpolationWeights() but does not throw if x is larger than getNumDimensions(); the length of x must be at least getNumDimensions() and the length of weighs must be at least getNumPoints().

◆ getDifferentiationWeights() [1/4]

std::vector<double> TasGrid::TasmanianSparseGrid::getDifferentiationWeights ( std::vector< double > const &  x) const

Returns the weights of the model outputs that combine to construct the approximate Jacobian matrix (derivative) at x.

The Jacobian of the k-th output at x is the sum of the model values of the k-th output times the Jacobian of the nodal functions at x.

Parameters
xis a vector of size getNumDimensions() with the coordinates of the point of interest in the transformed domain.
Returns
A vector of size getNumPoints() * getNumDimensions() with the differentiation weights, the order of the weights matches the order of the getPoints().

Note that using a vector x outside of the domain will result in undefined behavior, but will not throw an exception.

Exceptions
std::runtime_errorif x has an incorrect size.

◆ getDifferentiationWeights() [2/4]

std::vector<double> TasGrid::TasmanianSparseGrid::getDifferentiationWeights ( double const  x[]) const
inline

Overload that uses raw-array, does not check the array size.

Identical to getDifferentiationWeights() but does not throw if x is larger than getNumDimensions(); however, using shorter x is undefined behavior and will likely segfault.

◆ getDifferentiationWeights() [3/4]

void TasGrid::TasmanianSparseGrid::getDifferentiationWeights ( const std::vector< double > &  x,
std::vector< double > &  weights 
) const

Overload that uses the vector as a parameter.

Identical to getDifferentiationWeights() but the weights vector will be resized to size getNumPoints() * getNumDimensions() and the entries will be overwritten with the output of getDifferentiationWeights().

◆ getDifferentiationWeights() [4/4]

void TasGrid::TasmanianSparseGrid::getDifferentiationWeights ( const double  x[],
double  weights[] 
) const

Overload that uses raw-array, does not check the array size.

Identical to getDifferentiationWeights() but does not throw if x is larger than getNumDimensions(); the length of x must be at least getNumDimensions() and the length of weighs must be at least getNumPoints() * getNumDimensions().

◆ loadNeededValues() [1/2]

void TasGrid::TasmanianSparseGrid::loadNeededValues ( std::vector< double > const &  vals)

Provides the values of the model outputs at the needed points, or overwrites the currently loaded points.

In order to construct an interpolant, Tasmanian needs the values of the model outputs at each grid point. If there are needed points (i.e., getNumNeeded() is not zero), then vals must correspond to the model outputs at the needed points, otherwise, it must correspond to the loaded points.

Parameters
valsA vector that is logically divided into strips of size getNumOutputs() each strip corresponding to a single point. The order of the outputs must match the order of the points from either getNeededPoints() or getLoadedPoints(). If getNumNeeded() is zero, the total size must be getNumOutputs() times getNumLoaded(), otherwise, it must be getNumOutputs() times getNumNeeded().
Exceptions
std::runtime_errorif vals has an incorrect size.

Note: The needed points can always be cleared with clearRefinement() and new needed points can be assigned with setAnisotropicRefinement() or setSurplusRefinement().

◆ loadNeededValues() [2/2]

void TasGrid::TasmanianSparseGrid::loadNeededValues ( const double *  vals)

Overload that uses a raw-array, does not check the array size.

Identical to loadNeededPoints() but does not throw if vals has an incorrect size (but will segfault).

◆ loadNeededPoints()

void TasGrid::TasmanianSparseGrid::loadNeededPoints ( const double *  vals)
inline

Overload that uses a raw-array, does not check the array size.

Identical to loadNeededPoints() but does not throw if vals has an incorrect size (but will segfault).

◆ getLoadedValues()

const double* TasGrid::TasmanianSparseGrid::getLoadedValues ( ) const
inline

Returns the model values that have been loaded in the gird.

Returns a pointer to the internal data-structures, which must not be modified and will be invalidated by any operation that affects the loaded points, e.g., mergeRefinement() or loadNeededValues(). The model values will follow the internal Tasmanian order, identical to getLoadedPoints().

◆ evaluate() [1/2]

void TasGrid::TasmanianSparseGrid::evaluate ( std::vector< double > const &  x,
std::vector< double > &  y 
) const

Computes the value of the interpolant (or point-wise approximation) at the given point x.

This is the reference implementation that does not use any acceleration mode even if one is set. As a result, the calls to evaluate() are thread-safe but potentially slow.

Parameters
[in]xindicates the coordinate entries of a point within the domain of the sparse grid, the size must be equal to getNumDimensions().
[out]ywill contain the approximated model outputs corresponding to x, the vector will be resized to getNumOutputs().
Exceptions
std::runtime_errorif x has an incorrect size.

Note: calling evaluate() for a point outside of the domain is Mathematically incorrect, even though no exception will be generated.

◆ evaluate() [2/2]

void TasGrid::TasmanianSparseGrid::evaluate ( const double  x[],
double  y[] 
) const

Overload that uses raw-arrays, does not check the array size.

Identical to evaluate() but does not throw, assumes x has size at least getNumDimensions() and y is at least getNumOutputs(), will segfault if either is too short.

◆ evaluateBatch() [1/3]

template<typename FloatType >
void TasGrid::TasmanianSparseGrid::evaluateBatch ( std::vector< FloatType > const &  x,
std::vector< FloatType > &  y 
) const

Computes the value of the interpolant (or point-wise approximation) for a batch of points.

Identical to calling evaluate() for each point defined by x, but uses the specified acceleration mode and is potentially much faster, see TasGrid::TypeAcceleration for details.

Template Parameters
FloatTypeis either float or double indicating the precision to use.
Parameters
[in]xis logically divided into strips of size getNumDimensions() defining the coordinates of points within the sparse grid domain, see also evaluate().
[out]yis logically divided into x.size() / getNumDimensions() number of strips each with length getNumOutputs(), provides the approximated model outputs for each point defined by x. The vector will be resized, if the original size is incorrect.
Exceptions
std::runtime_errorif instantiated with float and the current acceleration mode is neither CUDA nor MAGMA, see TasGrid::TypeAcceleration.

Notes: this does not check if x.size() divides evenly.

The batch call:

grid.evaluateBatch(x, y);

is Mathematically equivalent to:

for(int i=0; i<x.size() / grid.getNumDimensions(); i++)
grid.evaluate(&x[i * grid.getNumDimensions()],
&y[i * grid.getNumOutputs()]);

However, depending on the acceleration mode, the performance can be significantly different.

◆ evaluateBatch() [2/3]

void TasGrid::TasmanianSparseGrid::evaluateBatch ( const double  x[],
int  num_x,
double  y[] 
) const

Overload that uses raw-arrays.

Raw-arrays do not provide size, thus the user must specify the number of points and the arrays are assumed to have sufficient size (otherwise the call will segfault). See also evaluate() and evaluateBatch().

Parameters
[in]xis logically organized into num_x strips of length getNumDimensions(), each strip will hold the coordinates of a point within the sparse grid domain.
[in]num_xis the total number of points stored in x.
[out]yis logically organized into num_x strips of length getNumOutputs(), each strip will be overwritten with the corresponding approximated model outputs.

The following two calls are equivalent:

grid.evaluateBatch(x, y); // using containers
grid.evaluateBatch(x.data(), x.size() / grid.getNumDimensions(), y.data()); // using raw arrays

◆ evaluateBatch() [3/3]

void TasGrid::TasmanianSparseGrid::evaluateBatch ( const float  x[],
int  num_x,
float  y[] 
) const

Overload using single precision and GPU/CUDA acceleration.

Works identical to the other raw-array overload but works only with CUDA and MAGMA acceleration modes.

Parameters
[in]xsee the double precision array overload
[in]num_xsee the double precision array overload
[in]ysee the double precision array overload
Exceptions
std::runtime_errorif the acceleration mode is not CUDA or MAGMA, see TasGrid::TypeAcceleration.
std::runtime_errorfrom failed calls to evaluateBatchGPU().

◆ evaluateBatchGPU()

template<typename FloatType >
void TasGrid::TasmanianSparseGrid::evaluateBatchGPU ( const FloatType  gpu_x[],
int  cpu_num_x,
FloatType  gpu_y[] 
) const

Overload that uses GPU raw-arrays.

Identical to the raw-array version of evaluateBatch(), but gpu_x and gpu_y must point to memory allocated on the CUDA device matching getGPUID(). Requires that Tasmanian was compiled with CUDA support and CUDA (or MAGAMA) acceleration mode has been enabled.

Exceptions
std::runtime_errorif Tasmanian was not build with -DTasmanian_ENABLE_CUDA=ON.
std::runtime_errorif calling for a grid type that doesn't have appropriate CUDA kernels, e.g., local polynomial or wavelet grids with order more than 2.

◆ evaluateFast() [1/2]

template<typename FloatType >
void TasGrid::TasmanianSparseGrid::evaluateFast ( const FloatType  x[],
FloatType  y[] 
) const
inline

Equivalent to evaluate() with enabled acceleration or evaluateBatch() with a batch of one point.

Some use cases still require performance but cannot batch multiple points together. The evaluateFast() method will work the same as evaluate(), but will use the specified acceleration mode.

Note: in older versions of Tasmanian, this function used to call a different set of algorithms optimized for single point evaluations; currently, evaluateBatch() will automatically switch to the appropriate mode depending on x.size() and num_x. Thus, this method is now an alias to evaluateBatch().

◆ evaluateFast() [2/2]

template<typename FloatType >
void TasGrid::TasmanianSparseGrid::evaluateFast ( std::vector< FloatType > const &  x,
std::vector< FloatType > &  y 
) const
inline

Alias to evaluateBatch().

Provided for consistency and backwards compatibility.

◆ integrate() [1/2]

void TasGrid::TasmanianSparseGrid::integrate ( std::vector< double > &  q) const

Computes the integral of each model output over the sparse grid domain.

The integration weight is assumed 1 unless another weight is associated with the underlying one dimensional rule, see TasGrid::TypeOneDRule for details.

Parameters
qwill be resized to getNumOutputs() and overwritten with the approximated integral of each model output.

The output of integrate() is Mathematically equivalent to calling getQuadratureWeights() and computing the sum of model values times weights. However, when the model values have been loaded integrate() is faster and more convenient.

◆ integrate() [2/2]

void TasGrid::TasmanianSparseGrid::integrate ( double  q[]) const

Overload that uses a raw-array.

Equivalent to integrate() but q must have sufficient size to write the result.

Parameters
qmust have size of at least getNumOutputs().

◆ differentiate() [1/2]

void TasGrid::TasmanianSparseGrid::differentiate ( std::vector< double > const &  x,
std::vector< double > &  jacobian 
) const

Computes the derivative (if available) of the surrogate model at an input point.

Parameters
xis the point of interest where the Jacobian should be evaluated.
jacobianwill be resized to getNumOutputs() * getNumDimensions() and overwritten with the Jacobian matrix.

◆ differentiate() [2/2]

void TasGrid::TasmanianSparseGrid::differentiate ( const double  x[],
double  jacobian[] 
) const

Overload that uses a raw-array.

Equivalent to differentiate() but jacobian must have sufficient size to write the result.

◆ setDomainTransform() [1/2]

void TasGrid::TasmanianSparseGrid::setDomainTransform ( std::vector< double > const &  a,
std::vector< double > const &  b 
)

Set a linear domain transformation.

By default integration and interpolation are performed on a canonical interval [-1, 1], with the exception of Fourier grids using [0, 1] and some Gauss rules, see TasGrid::TypeOneDRule. The linear transformation will shift the interval to an arbitrary [a, b] (or shift and scale for the unbounded case). Different values can be specified for each dimension and the transformation will be automatically applied to every operation that uses points, e.g., getPoints() and evaluate() will return/accept only transformed points.

Setting or changing the transformation will change the points and weights, therefore, the transformations should be set immediately after calling a make command and before calling get-points or computing model values. Changing the transformation will not throw, but will likely invalidate the loaded data and should be used with extreme caution (the validity of the underlying Mathematics is left to the user to analyze).

Parameters
awith size getNumDimensions() specifies the a transformation parameter for each input
bwith size getNumDimensions() specifies the b transformation parameter for each input
Exceptions
std::runtime_errorif the grid is empty.
std::invalid_argumentif either input has incorrect size.

◆ setDomainTransform() [2/2]

void TasGrid::TasmanianSparseGrid::setDomainTransform ( const double  a[],
const double  b[] 
)

Overload using raw-arrays.

Identical to setDomainTransform() but does not check for the size of a and b, although it still checks if the grid is empty.

◆ isSetDomainTransfrom()

bool TasGrid::TasmanianSparseGrid::isSetDomainTransfrom ( ) const

Returns true if a linear domain transformation has been set, false otherwise.

Allows to check if there is a set transformation or if working on the canonical domain.

◆ clearDomainTransform()

void TasGrid::TasmanianSparseGrid::clearDomainTransform ( )

Removes the domain transformation.

Use with extreme caution, see setDomainTransform().

◆ getDomainTransform() [1/2]

void TasGrid::TasmanianSparseGrid::getDomainTransform ( std::vector< double > &  a,
std::vector< double > &  b 
) const

Returns the two vectors used to call setDomainTransform().

If the grid is empty or if no transformation has been set, this will return empty vectors.

Parameters
[out]awill be resized to getNumDimensions() and overwritten with the a vector of setDomainTransform().
[out]bwill be resized to getNumDimensions() and overwritten with the b vector of setDomainTransform().

Note: this works the same regardless of which setDomainTransform() overload has been used.

◆ getDomainTransform() [2/2]

void TasGrid::TasmanianSparseGrid::getDomainTransform ( double  a[],
double  b[] 
) const

Returns the values of the two vectors used to call setDomainTransform().

Assuming that the inputs have sufficient size, will overwrite the getNumDimensions() entries with the domain transformation.

Exceptions
std::runtime_errorif the grid is empty or if the domain transformation has not been set.

Note: this works the same regardless of which setDomainTransform() overload has been used.

◆ setConformalTransformASIN()

void TasGrid::TasmanianSparseGrid::setConformalTransformASIN ( std::vector< int > const &  truncation)

Set conformal transformation using truncated Maclaurin series of the arcsin() function.

Conformal transformations apply a non-linear map to the points, weights and basis functions of a grid, in an attempt to accelerate convergence for some types of functions. The objective is to expand the area of analytic extension of the target function and thus accelerate convergence; however, if applied to the wrong function, conformal transformations can deteriorate the accuracy. Characterizing the functions that would most benefit from conformal maps is an active area of research.

The truncated Maclaurin series of the arcsin() work well with functions that have a pole close to the edges of the domain. See:
P. Jantsch, C. G. Webster, Sparse Grid Quadrature Rules Based on Conformal Mappings, Sparse Grids and Applications - Miami 2016 pp 117–134.

Parameters
truncationis a vector of size getNumDimensions() that indicates the number of terms to keep in each direction. The more terms, the more pronounced the transformation becomes; however, too many terms can in fact produce a pole even worse than the original.

◆ clearLevelLimits()

void TasGrid::TasmanianSparseGrid::clearLevelLimits ( )
inline

Removes the currently set level limits.

Once level limits have been set (by either the make or set-refinement commands), the limits will be used for all follow-on commands unless either overwritten with a new set of limits of cleared with this command.

◆ getLevelLimits()

std::vector<int> TasGrid::TasmanianSparseGrid::getLevelLimits ( ) const
inline

Return the currently set level limits.

Returns the limits that have been set with the last command. If no limits have been set, an empty vector will be returned.

◆ setAnisotropicRefinement() [1/2]

void TasGrid::TasmanianSparseGrid::setAnisotropicRefinement ( TypeDepth  type,
int  min_growth,
int  output,
const std::vector< int > &  level_limits 
)

Set refinement using anisotropic estimate for the optimal points.

Computes the anisotropic coefficients based on the current set of loaded points, then adds more points to the grid by selecting the points with largest coefficients. The new points are set to needed. See also estimateAnisotropicCoefficients()

Parameters
typeindicates the type of estimate to use, e.g., total degree or curved; regardless of the specified rule the interpolation estimate is computed (not quadrature).
min_growthis the minimum number of points to add to the grid, e.g., if the model can be executed in parallel then minimum number of points is needed to ensure occupancy for all computing resources. The value of min_growth can never be less than 1, but the actual number of points can exceed the min_growth depending on the weights and growth of the one dimensional rule.
outputindicates which output to use to compute the estimate, using -1 indicates to use all outputs (possible with Sequence and Fourier grids only, Global grids require a specific output).
level_limits(if not empty) will be used to overwrite the currently set limits. The limits must be either empty or have size getNumDimensions(); if empty, the current set of limits will be used.
Exceptions
std::runtime_errorif called during dynamic construction, or there are no loaded points, or there are no outputs.
std::invalid_argumentif min_growth is not positive, output is out of range, or level_limits has the wrong size.

◆ setAnisotropicRefinement() [2/2]

void TasGrid::TasmanianSparseGrid::setAnisotropicRefinement ( TypeDepth  type,
int  min_growth,
int  output,
const int *  level_limits = nullptr 
)

Overload using raw-array.

Identical to setAnisotropicRefinement() with the exception that level_limits is a raw array, empty is equivalent to nullptr and the size is not checked.

◆ estimateAnisotropicCoefficients() [1/2]

std::vector<int> TasGrid::TasmanianSparseGrid::estimateAnisotropicCoefficients ( TypeDepth  type,
int  output 
) const
inline

Estimate the anisotropic rates of coefficient decay for different direction.

Available for Global, Sequence and Fourier grids, the anisotropic coefficients describe the space of the dominant basis functions needed to construct optimal approximation to the model. See the documentation of TasGrid::TypeDepth of the different formulas. The estimate requires that values have been loaded in the grid, i.e., getNumOutputs() > 0 and getNumLoaded() > 0.

Parameters
typeis the assumed type of optimal space, i.e., total degree, total degree with log-correction, or hyperbolic cross section.
outputdetermines the output to use for the inference of the coefficients, for Global grid it must be between 0 and getNumOutputs() -1, in other cases, it can also be set to -1 to indicate "all outputs" or for each basis use the larges coefficient among all outputs.
Returns
the xi and eta (if used) parameters of the estimate indicated by type.
Exceptions
std::runtime_errorif the grid is empty, has no outputs or no loaded points.
std::invalid_argumentif the output is out of range.

◆ estimateAnisotropicCoefficients() [2/2]

void TasGrid::TasmanianSparseGrid::estimateAnisotropicCoefficients ( TypeDepth  type,
int  output,
std::vector< int > &  weights 
) const

Overload that writes the result to a parameter.

The inputs are equivalent to estimateAnisotropicCoefficients() but the result is returned into weights.

◆ setSurplusRefinement() [1/3]

void TasGrid::TasmanianSparseGrid::setSurplusRefinement ( double  tolerance,
int  output,
std::vector< int > const &  level_limits 
)

Refine the grid based on the surplus coefficients, Sequence grids and Global grids with a sequence rule.

Using the relative magnitude of the surplus coefficients, add more points to the grid and set them to needed. The approach differs from the local (polynomial or wavelet) case in the interpretation of the surpluses, instead of local-spacial estimate the surpluses are interpreted as indicators of needed polynomial basis functions.

Parameters
toleranceindicates the cutoff threshold for refinement, points will not be included once the magnitude of the relative surplus drops below the tolerance.
outputindicates the output to use for the surpluses, the Sequence grid accepts -1 to indicate the use of all outputs.
level_limitsindicates a new set of limits, if empty the currently set limits will be used.
Exceptions
std::runtime_errorif the called during the construction process, if the grid is empty, of if the grid has no outputs or values.
std::invalid_argumentif the output is out of range, tolerance is negative, or if the level_limits has the wrong size.

◆ setSurplusRefinement() [2/3]

void TasGrid::TasmanianSparseGrid::setSurplusRefinement ( double  tolerance,
TypeRefinement  criteria,
int  output,
std::vector< int > const &  level_limits,
std::vector< double > const &  scale_correction = std::vector< double >() 
)

Refine the grid based on the surplus coefficients, Local-Polynomial and Wavelet grids.

Using the relative magnitude of the surplus coefficients, add more points to the grid and set them to needed. This method uses the hierarchical and local structure of the one dimensional rule and interprets the surplus as a local error indicator. The refinement can be combined with tests for completeness of the hierarchy (prioritizing parents or even the entire ancestry) or local anisotropy that can manifest even in globally isotropic cases. See TasGrid::TypeRefinement for details.

If this method is called on a Global or a Sequence grid, the criteria will be ignored and the method will use the Global/Sequence variant.

Parameters
toleranceindicates the cutoff threshold for refinement, points will not be included once the magnitude of the relative surplus drops below the tolerance.
criteriaindicates how to prioritize the hierarchy and/or local anisotropy.
outputindicates a specific output to use for the refinement, by default (when -1) all model outputs will be considered together.
level_limitsindicates a new set of limits, if empty the currently set limits will be used.
scale_correctionis a set of weights that would multiply the surpluses before the tolerance test; the correction can be used to modify the threshold test (e.g., multiply the surplus by the integral of the basis functions) or to guide the refinement towards certain regions of the domain. The values of the correction terms are organized in order that matches the order of getNumLoaded() and there is one weight per active output (either 1 or getNumOutputs() depending whether output is -1 or positive). If the scale correction is empty, then no correction is used (i.e., using correction of 1.0).
Exceptions
std::runtime_errorif called during construction, or if the grid is empty or has no outputs or loaded points, of if the grid has incompatible type (not local polynomial or wavelet).
std::invalid_argumentif output is out of range, of if the level limits and/or scale correction have incorrect size.

◆ setSurplusRefinement() [3/3]

void TasGrid::TasmanianSparseGrid::setSurplusRefinement ( double  tolerance,
TypeRefinement  criteria,
int  output = -1,
const int *  level_limits = nullptr,
const double *  scale_correction = nullptr 
)

Overload that uses raw-arrays.

The scale_correction is a potentially large vector and using a raw array avoids a large data-copy when calling from C/Python/Fortran interfaces. Otherwise the method behaves the same, but does not throw if the arrays have incorrect size (will probably segfault).

◆ clearRefinement()

void TasGrid::TasmanianSparseGrid::clearRefinement ( )

Remove all needed points from the grid.

Once a refinement is set, but before the new values have been loaded, the refinement can be canceled with this command. After this call, all needed points will be erased.

◆ mergeRefinement()

void TasGrid::TasmanianSparseGrid::mergeRefinement ( )

Merges the loaded and needed points into a single grid, resets all loaded values to zero.

This method allows refinement to be used in cases where the model values cannot be computed at the grid points, e.g., when working with random or unstructured data. Once a refinement is set, the new points can be merged without model values which will result in a larger grid with all (previous and new) loaded values reset to 0.0. Afterwards, the setHierarchicalCoefficients() method can be used to set a new set of coefficients, e.g., inferred from the data and the hierarchical basis values.

◆ beginConstruction()

void TasGrid::TasmanianSparseGrid::beginConstruction ( )

Begin a dynamic construction procedure.

Initializes the internal data-structures needed for the construction procedure (which is pretty cheap) and makes a call to clearRefinement().

Note: after this call, setSurplusRefinement() and setAnisotropicRefinement() cannot be called until finishConstruction() is issued.

Note: the construction process can be initiated before any model values have been loaded, in such case, the initial set of points will always come first in a call to getCandidateConstructionPoints().

◆ isUsingConstruction()

bool TasGrid::TasmanianSparseGrid::isUsingConstruction ( ) const
inline

Returns true if the dynamic construction procedure has been initialized, false otherwise.

Simply returns the internal flag.

◆ getCandidateConstructionPoints() [1/3]

std::vector<double> TasGrid::TasmanianSparseGrid::getCandidateConstructionPoints ( TypeDepth  type,
std::vector< int > const &  anisotropic_weights = std::vector< int >(),
std::vector< int > const &  level_limits = std::vector< int >() 
)

Generate a sorted list of points weighted by descending importance.

The importance is inferred from the user provided anisotropic weighs using the formula specified by the type, see TasGrid::TypeDepth. The full tensor types will fallback to type_level.

Unlike the batch procedures in setAnisotropicRefinement(), there is no expectation that the entire batch is processed, in fact only a small subset of the most important points should be computed and loaded, then the getCandidateConstructionPoints() should be called again and it will return an updated list of points. Enough points should be used from the top of the list to ensure occupancy across computing resources, e.g., CPU cores or MPI ranks, but the bottom set of points will contribute less and less to the overall accuracy.

Parameters
typesets the formula to use when weighting the potential points, see TasGrid::TypeDepth.
anisotropic_weightsare the xi and eta parameters for the formula, the vector must have the correct size, either getNumDimensions() or twice as much to handle the curved weights.
level_limits(if not empty) will be used to overwrite the currently set limits. The limits must be either empty or have size getNumDimensions(); if empty, the current set of limits will be used.
Returns
A vector organized in strips of length getNumDimensions() that indicate the coordinates of the points to use as model inputs, unlike the getPoints() command, the points here are arranged in decreasing importance.
Exceptions
std::runtime_errorif the grid is empty, Local Polynomial or this is called before beginConstruction().
std::invalid_argumentif the anisotropic_weights or level_limits have incorrect size.

◆ getCandidateConstructionPoints() [2/3]

std::vector<double> TasGrid::TasmanianSparseGrid::getCandidateConstructionPoints ( TypeDepth  type,
int  output,
std::vector< int > const &  level_limits = std::vector< int >() 
)

Essentially the same as getCandidateConstructionPoints() but the weights are obtained from a call to estimateAnisotropicCoefficients().

This method is the construction equivalent to setAnisotropicRefinement(). One notable difference is that this function will not throw if there are no loaded points, instead isotropic coefficient will be used until enough points are loaded so that estimates for the coefficients can be computed.

Parameters
typesets the formula to use when weighting the potential points, see TasGrid::TypeDepth.
outputindicates which coefficients will be used to estimate the anisotropic decay rate, when working with Sequence and Fourier grids this can be set to -1 to use all outputs.
level_limits(if not empty) will be used to overwrite the currently set limits. The limits must be either empty or have size getNumDimensions(); if empty, the current set of limits will be used.
Returns
a vector organized in strips of length getNumDimensions() that indicate the coordinates of the points to use as model inputs, unlike the getPoints() command, the points here are arranged in decreasing importance.
Exceptions
std::runtime_errorif the grid is empty, Local Polynomial or this is called before beginConstruction().
std::invalid_argumentif the level_limits have incorrect size or output is out of range.

◆ getCandidateConstructionPoints() [3/3]

std::vector<double> TasGrid::TasmanianSparseGrid::getCandidateConstructionPoints ( 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 >() 
)

Returns a sorted list of points weighted by descending importance using the hierarchical surpluses.

This is the construction equivalent to the Local Polynomial setSurplusRefinement(). The inputs are the same as setSurplusRefinement() except the returned points will be sorted by decreasing surpluses. Similar to the other getCandidateConstructionPoints() overloads, this can be called before any model values have been loaded.

◆ loadConstructedPoints() [1/2]

void TasGrid::TasmanianSparseGrid::loadConstructedPoints ( std::vector< double > const &  x,
std::vector< double > const &  y 
)

Add pairs of points with associated model values.

This is the construction equivalent to loadNeededValues(), the main difference is that any number of points can be loaded here and the points can be in any arbitrary order (they have to correspond to the model values in this call only).

Parameters
xis a vector with strips of size getNumDimensions() indicating the points where the model values were computed. The points do not have to be in any order; however, every point has to match a potential grid points.
yis a vector with strips of size getNumOutputs() indicating the model outputs corresponding each of the provided points.
Exceptions
std::runtime_errorif the number of strips in y are less than those in x.

Note: regardless of the grid type and rule, as the depth/level increases the points become dense in the domain; thus every point is theoretically a potential grid points. However, if a random point is chosen than the level may easily overflow the range of the int type, the safest option is to make sure the points match the ones returned by getCandidateConstructionPoints().

◆ loadConstructedPoints() [2/2]

void TasGrid::TasmanianSparseGrid::loadConstructedPoints ( const double  x[],
int  numx,
const double  y[] 
)

Same as loadConstructedPoint() but using arrays in place of vectors (array size is not checked)

Does not throw on incorrect array size, but will likely segfault.

◆ finishConstruction()

void TasGrid::TasmanianSparseGrid::finishConstruction ( )

End the procedure, clears flags and unused constructed points, can go back to using regular refinement.

After this call, the construction methods getCandidateConstructionPoints() and loadConstructedPoint() cannot be used until a new call to beginConstruction().

Note: finalizing the construction process can potentially lead to loss of data. The constructed points can be loaded in any order, but the points need to satisfy constraints before being incorporated within a grid. Some grids require that points form a lower complete set, or the point must form a connected graph, or there must be enough points to complete a tensor (when the rules grow by more than one point). In such scenarios, the points and values are stored in a temporary structure until enough data is present to add them to the grid. Finalizing the construction will delete all data for points that have not been incorporated yet.

◆ getHierarchicalCoefficients()

const double* TasGrid::TasmanianSparseGrid::getHierarchicalCoefficients ( ) const

Return a reference to the internal data-structure that stores the hierarchical coefficients.

All types of grids (except Global grids), use a hierarchical basis representation where the interpolant is expressed as a set of coefficients times the basis functions. For Local Polynomial and Wavelet grids, the coefficients are commonly known as hierarchical surpluses; for Sequence grid, the coefficients correspond to the Newton polynomials; for Fourier grids, the coefficients are the linear combination of discrete Fourier transform coefficients.

The returned array has getNumLoaded() strips of size getNumOutputs() where each strip corresponds to one basis function. In the case of Fourier grids, the coefficients are complex numbers and (for performance reasons) the real and complex parts are stored separately, i.e., the first getNumLoaded() strips hold the real parts and there is a second set of getNumLoaded() strips that hold the complex components. If the grid is empty or there are no outputs or no loaded points, then this returns nullptr.

Note: modifying the coefficients through this pointer leads to undefined behavior, use setHierarchicalCoefficients() instead.

◆ setHierarchicalCoefficients() [1/2]

void TasGrid::TasmanianSparseGrid::setHierarchicalCoefficients ( const std::vector< double > &  c)

Overwrites the current set of coefficients (and loaded values) with the ones provided.

Discards the current set of loaded values and the associated hierarchical coefficients, and replaces both with the data provided here. The coefficients are overwritten, while the values are inferred, i.e., the opposed from the use case of loadNeededValues() where the model values are provided and the coefficients are computed.

Parameters
cis a vector of getNumLoaded() strips of size getNumOutputs(), each strip corresponds to the coefficients of one basis function. Fourier coefficients are an exceptions as they require twice as many strips where the first set corresponds to the real components and the second set corresponds to the complex components of the coefficients.
Exceptions
std::runtime_erroris the number of coefficients is incorrect.

◆ setHierarchicalCoefficients() [2/2]

void TasGrid::TasmanianSparseGrid::setHierarchicalCoefficients ( const double  c[])
inline

Overload that uses raw-arrays.

Identical to setHierarchicalCoefficients() but does not check for the size of the array.

◆ evaluateHierarchicalFunctions() [1/3]

void TasGrid::TasmanianSparseGrid::evaluateHierarchicalFunctions ( std::vector< double > const &  x,
std::vector< double > &  y 
) const

Computes the values of the hierarchical function basis at the specified points.

The method for getInterpolationWeights() computes the values of the nodal basis function, e.g., Lagrange polynomials, while this computes the values of the hierarchical functions, e.g., Newton polynomials. The call is optimized for multiple points at a time. The output consists of strips of size getNumPoints() and one strip corresponds to one point provided in x. Effectively this constructs a matrix in column major format where each column corresponds to a single input point x and each row corresponds to a single grid point.

When working with Fourier grids, the real and complex part of each basis is interlaced, i.e., the strip size is twice as large; see the array overload.

Parameters
[in]xhas strips of size getNumDimensions() indicating a point in the domain where the basis is to be computed.
[out]ywill have the same number of strips with size getNumPoints() (or twice that to accommodate the real and complex parts of the basis).

Note: if the the output from getHierarchicalCoefficients() is interpreted as a matrix in column major format with leading dimension of getNumOutputs(), and if the output of this call is interpreted as a matrix in column major format with leading dimension getNumPoints(), then the result of the product of the coefficients times the basis will be equal to the result of evaluateBatch().

◆ evaluateHierarchicalFunctions() [2/3]

std::vector<double> TasGrid::TasmanianSparseGrid::evaluateHierarchicalFunctions ( std::vector< double > const &  x) const
inline

Overload that returns the result.

Useful for direct initialization of vectors.

◆ evaluateHierarchicalFunctions() [3/3]

void TasGrid::TasmanianSparseGrid::evaluateHierarchicalFunctions ( const double  x[],
int  num_x,
double  y[] 
) const

Array overload, the inputs must have pre-allocated and correct size.

The size of x must be num_x times getNumDimensions() and the size of y must be num_x times getNumPoints() (or twice that for the Fourier grids).

Example of returning the result in a vector of complex numbers:

int num_x = ...; // set the number of points
std::vector<double> x(grid.getNumDimensions() * num_x);
// initialize x
std::vector<std::complex<double>> y(grid.getNumPoints() * num_x);
grid.evaluateHierarchicalFunctions(x.data(), num_x, reinterpret_cast<double*>(y.data()));
// at this point y is loaded with the complex numbers
The master-class that represents an instance of a Tasmanian sparse grid.
Definition: TasmanianSparseGrid.hpp:293
int getNumPoints() const
Returns getNumLoaded() if positive, otherwise returns getNumNeeded(), see getPoints().
Definition: TasmanianSparseGrid.hpp:661
void evaluateHierarchicalFunctions(std::vector< double > const &x, std::vector< double > &y) const
Computes the values of the hierarchical function basis at the specified points.
int getNumDimensions() const
Return the dimensions of the grid, i.e., number of model inputs.
Definition: TasmanianSparseGrid.hpp:642
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().
Definition: TasmanianSparseGrid.hpp:2327

◆ evaluateSparseHierarchicalFunctions()

void TasGrid::TasmanianSparseGrid::evaluateSparseHierarchicalFunctions ( const std::vector< double > &  x,
std::vector< int > &  pntr,
std::vector< int > &  indx,
std::vector< double > &  vals 
) const

Constructs a sparse matrix with the values of the hierarchical basis functions.

The local basis functions associated with Local Polynomial and Wavelet grids means that only a small set of the basis functions will be non-zero at any given point in the domain. Effectively, the matrix build by evaluateHierarchicalFunctions() will be sparse. This method generates the sparse matrix in compressed format in non-zeroes per entry in x.

Parameters
[in]xhas strips of size getNumDimensions() indicating a point in the domain where the basis is to be computed.
[out]pntrwill have size one more than the number of strips in x, pntr[i] indicate the starting offsets of the entries for the i-th point, (the first entry will always be zero) and pntr.back() is the total number of non-zeros.
[out]indxholds the indexes of the supported basis functions, i.e., indx[pntr[i]] ... indx[pntr[i+1] - 1] are the indexes of the basis functions supported over the i-th point.
[out]valsare the numeric values of the basis functions corresponding to the indexes in indx.

◆ getHierarchicalSupport()

std::vector<double> TasGrid::TasmanianSparseGrid::getHierarchicalSupport ( ) const

Returns the support of the hierarchical functions.

The Global, Sequence and Fourier grids rely on basis functions with support over the entire domain. However, local-polynomial and wavelet grids have basis that is restricted to a hypercube centered at the corresponding point and with equal support in each direction. This method returns the length for each basis function and direction.

Returns
either an empty vector (for an empty grid), or a vector of size getNumDimensions() times getNumPoints() organized in strips of length getNumDimensions(). Each strip corresponds to a basis function and if (in a some direction for some basis function) an x-point falls away from the corresponding grid point at a distance larger than the support, then the basis will evaluate to zero.

Note: the support of all basis functions is logically restricted to the grid domain, i.e., the effective support equals the intersection of the hypercube and the domain.

◆ integrateHierarchicalFunctions() [1/3]

std::vector<double> TasGrid::TasmanianSparseGrid::integrateHierarchicalFunctions ( ) const
inline

Returns the integrals of the hierarchical basis functions.

Returns a vector of size getNumPoints() that holds the integral of each basis function. If the output of getHierarchicalCoefficients() is interpreted as a matrix, the product of that matrix times this vector is the same as the output of TasmanianSparseGrid::integrate().

One use case of the integrals is to add scale correction to the surplus refinement, e.g., rescale the coefficients by the integral of the basis function:

// ... load values ... //
auto correction = grid.integrateHierarchicalFunctions();
grid.setSurplusRefinement(1.E-4, TasGrid::refine_classic, 0, {}, correction);
@ refine_classic
Isotropic refinement using only the children and disregarding missing parents.
Definition: tsgEnumerates.hpp:427
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().
Definition: TasmanianSparseGrid.hpp:2301

◆ integrateHierarchicalFunctions() [2/3]

void TasGrid::TasmanianSparseGrid::integrateHierarchicalFunctions ( std::vector< double > &  integrals) const
inline

Overload where the vector is passes as a parameter.

Parameters
integralswill be resized to getNumPoints() and overwritten with the integrals of the hierarchical basis functions.

◆ integrateHierarchicalFunctions() [3/3]

void TasGrid::TasmanianSparseGrid::integrateHierarchicalFunctions ( double  integrals[]) const

Overload using raw-arrays.

Useful mostly for the C/Python/Fortran interfaces.

Parameters
integralsmust have size getNumPoints() and will be overwritten with the integrals of the hierarchical basis functions.

◆ getGlobalPolynomialSpace()

std::vector<int> TasGrid::TasmanianSparseGrid::getGlobalPolynomialSpace ( bool  interpolation) const

Returns the powers of the polynomial that span the effective basis, Global and Sequence grids only.

Different rules have different growth in number of points and different exactness with respect to integration and interpolation. This method returns the actual polynomial space spanned by the basis used in the grid.

Parameters
interpolationdetermines whether to look for the space of exactly interpolated polynomials (true), or the space of exactly integrated polynomials (false).
Returns
a vector with getNumPoints() number of strips each having getNumDimensions() entries, each integer corresponds to the polynomial power in that direction, e.g., for a two dimensional grid (0, 0, 1, 0) indicates the constant in all directions and linear in the first input constant in the second.
Exceptions
std::runtime_errorif the grid is not Global or Sequence.

◆ printStats()

void TasGrid::TasmanianSparseGrid::printStats ( std::ostream &  os = std::cout) const

Prints short human-readable text describing the grid properties.

Among the outputs are type of the grid, number of inputs/outputs/points, one dimensional rule, and currently selected acceleration type.

Parameters
osindicates the stream to use for the text output, defaults to std::cout.

◆ enableAcceleration() [1/2]

void TasGrid::TasmanianSparseGrid::enableAcceleration ( TypeAcceleration  acc)

Change the current acceleration mode to the one specified.

Sets a new acceleration mode and releases the cached data-structures for the old mode.

Parameters
accis the new acceleration mode, see TasGrid::TypeAcceleration.

Note: any acceleration method can be set regardless whether it is available. If the selected mode has not been enabled in CMake, this will select the "next-best" mode, thus enableAcceleration() is a suggestion rather than a hard request. See getAccelerationType().

◆ enableAcceleration() [2/2]

void TasGrid::TasmanianSparseGrid::enableAcceleration ( TypeAcceleration  acc,
int  new_gpu_id 
)

Combines calls the enableAcceleration(), getGPUID() and allows for user provided handles.

The purpose of this method is to allow for one-show setup of the acceleration mode and gpu_id.

Parameters
accis a new acceleration mode, see TasGrid::TypeAcceleration.
new_gpu_idis the new device id to use for acceleration, the number must be between 0 and getNumGPUs() - 1.

◆ favorSparseAcceleration()

void TasGrid::TasmanianSparseGrid::favorSparseAcceleration ( bool  favor)

Set the preferred back-end algorithm for Local Polynomial grids.

Usually the Local Polynomial grids use sparse data-structures and sparse linear algebra, but if the fill not sufficiently small then dense methods can yield better performance at the expense of higher memory footprint. In TasGrid::accel_cpu_blas mode, the dense algorithm will be used when the fill increases above the 10% threshold. In CUDA mode, even if the fill is computed, switching modes incurs a large overhead and thus the sparse algorithm is always favored. Sparse computations work better on newer GPU architectures (e.g., Pascal and Volta) and consumer (gaming) cards with reduced double-precision capabilities, but older devices coupled with small grids may work better with the dense algorithm. Hence, Tasmanian includes the manual option to select the desired algorithm.

Parameters
favorif set to true will force the sparse back-end, if set to false will force the dense back-end. If the mode has been forced already, calling the method with the converse favor will reset mode to the default (automatic) selection.

◆ setCuBlasHandle()

void TasGrid::TasmanianSparseGrid::setCuBlasHandle ( void *  handle)

Takes a user provided cuBlas handle.

Parameters
handlemust be a valid and initialized cublasHandle_t
Exceptions
std::runtime_errorif CUDA is not enabled in CMake and enableAcceleration()

◆ setCuSparseHandle()

void TasGrid::TasmanianSparseGrid::setCuSparseHandle ( void *  handle)

Takes a user provided cuSparse handle.

Parameters
handlemust be a valid and initialized cusparseHandle_t
Exceptions
std::runtime_errorif CUDA is not enabled in CMake and enableAcceleration()

◆ setCuSolverHandle()

void TasGrid::TasmanianSparseGrid::setCuSolverHandle ( void *  handle)

Takes a user provided cuSparse handle.

Parameters
handlemust be a valid and initialized cusolverDnHandle_t
Exceptions
std::runtime_errorif CUDA is not enabled in CMake and enableAcceleration()

◆ setRocBlasHandle()

void TasGrid::TasmanianSparseGrid::setRocBlasHandle ( void *  handle)

Takes a user provided cuBlas handle.

Parameters
handlemust be a valid and initialized rocblas_handle
Exceptions
std::runtime_errorif HIP is not enabled in CMake and enableAcceleration()

◆ setRocSparseHandle()

void TasGrid::TasmanianSparseGrid::setRocSparseHandle ( void *  handle)

Takes a user provided cuSparse handle.

Parameters
handlemust be a valid and initialized rocsparse_handle
Exceptions
std::runtime_errorif HIP is not enabled in CMake and enableAcceleration()

◆ setSycleQueue()

void TasGrid::TasmanianSparseGrid::setSycleQueue ( void *  queue)

Takes a user provided sycl::queue handle.

Parameters
queuemust be a valid and initialized sycl::queue
Exceptions
std::runtime_errorif DPC++ is not enabled in CMake and enableAcceleration()

◆ getAccelerationType()

TypeAcceleration TasGrid::TasmanianSparseGrid::getAccelerationType ( ) const
inline

Returns the current effective acceleration mode.

Returns the acceleration mode that will be used, i.e., the one selected internally based on the request made in enableAcceleration().

◆ isAccelerationAvailable()

static bool TasGrid::TasmanianSparseGrid::isAccelerationAvailable ( TypeAcceleration  acc)
static

Returns whether a specific mode can be enabled.

Based on the CMake compile time options, this method returns true for all acceleration modes that will be enabled, and false for all modes that will be replaced by a fallback alternative.

◆ setGPUID()

void TasGrid::TasmanianSparseGrid::setGPUID ( int  new_gpu_id)

Select the current CUDA device.

Select the CUDA device to be used for the CUDA acceleration types, default device is 0.

Parameters
new_gpu_idis the CUDA device ID of the new device, the number is between 0 and getNumGPUs() - 1.
Exceptions
std::runtime_errorif the new_gpu_id is out of range.

◆ getGPUID()

int TasGrid::TasmanianSparseGrid::getGPUID ( ) const
inline

Returns the currently set CUDA device.

Does not throw if CUDA is not enabled in CMake, instead the default 0 is returned.

◆ getNumGPUs()

static int TasGrid::TasmanianSparseGrid::getNumGPUs ( )
inlinestatic

Return the number of visible CUDA devices.

Simple wrapper to cudaGetDeviceCount().

Use the tasgrid command line tool to see all available devices:

tasgrid -v

◆ getGPUMemory()

static int TasGrid::TasmanianSparseGrid::getGPUMemory ( int  gpu)
static

Return the available device memory, in units of MB.

Simple wrapper to cudaGetDeviceProperties() returning totalGlobalMem divided by 2^20.

Parameters
gpuis the CUDA device ID to be queried, if the device is out of range then 0 will be returned (i.e., Tasmanian will not throw).

◆ getGPUName()

static std::string TasGrid::TasmanianSparseGrid::getGPUName ( int  gpu)
static

Return the CUDA device name.

Simple wrapper to cudaGetDeviceProperties() returning "name".

Parameters
gpuis the CUDA device ID to be queried, if the device is out of range then empty string will be returned (i.e., Tasmanian will not throw).
Returns
the CUDA device name.

◆ evaluateHierarchicalFunctionsGPU()

template<typename FloatType >
void TasGrid::TasmanianSparseGrid::evaluateHierarchicalFunctionsGPU ( const FloatType  gpu_x[],
int  cpu_num_x,
FloatType  gpu_y[] 
) const

Computes the values of the hierarchical function basis at the specified points (CUDA version).

Equivalent to evaluateHierarchicalFunctions() but using arrays allocated on the CUDA device.

Template Parameters
FloatTypemust be either float or double to indicate the precision used by the CUDA kernels.
Parameters
gpu_xmust have size getNumDimensions() times cpu_num_x and must be allocated on the currently set CUDA device.
cpu_num_xis an integer (located on the CPU memory) indicating the number of points.
gpu_ymust have size getNumPoints() times cpu_num_x and must be allocated on the currently set CUDA device.
Exceptions
std::runtime_errorif the grid is Global or Wavelet or if the currently set acceleration mode is not compatible, i.e., TasGrid::accel_none or TasGrid::accel_cpu_blas. Also, if CUDA has not been enabled at compile time.

Note: will not work for LocalPolynomial grids with order bigger than 2 or Wavelets with order 3.

◆ evaluateSparseHierarchicalFunctionsGPU()

template<typename FloatType >
void TasGrid::TasmanianSparseGrid::evaluateSparseHierarchicalFunctionsGPU ( const FloatType  gpu_x[],
int  cpu_num_x,
int *&  gpu_pntr,
int *&  gpu_indx,
FloatType *&  gpu_vals,
int &  num_nz 
) const

Computes the values of the hierarchical function basis at the specified points (sparse/CUDA version).

Equivalent to evaluateSparseHierarchicalFunctions() but using arrays allocated on the CUDA device.

Template Parameters
FloatTypemust be either float or double to indicate the precision used by the CUDA kernels.
Parameters
[in]gpu_xmust have size getNumDimensions() times cpu_num_x and must be allocated on the currently set CUDA device.
[in]cpu_num_xis an integer (located on the CPU memory) indicating the number of points.
[out]gpu_pntrwill be allocated to size cpu_num_x + 1 and will hold the offsets of indexes for each point, the current memory will not be freed.
[out]gpu_indxwill be allocated to the number of non-zeros and will hold the indexes of the non-zeros for each point, the current memory will not be freed.
[out]gpu_valswill be allocated to the number of non-zeros and will hold the non-zero values of the basis function, the current memory will not be freed.
[out]num_nzis an integer located on the CPU, will be overwritten to the total number of non-zeros.
Exceptions
std::runtime_errorif the grid is not Local Polynomial or if the currently set acceleration mode is not compatible, i.e., TasGrid::accel_none or TasGrid::accel_cpu_blas. Also, if CUDA has not been enabled at compile time.

Note: will not work for LocalPolynomial grids with order bigger than 2.

◆ operator EvaluateCallable()

TasGrid::TasmanianSparseGrid::operator EvaluateCallable ( ) const
inline

Custom conversion to a callable method using the TasDREAM::DreamPDF signature.

This conversion allows an instance of TasmanianSparseGrid to be passed as input to any method that expects TasDREAM::DreamPDF, TasDREAM::DreamModel or TasDREAM::DreamMergedLikelyModel.

◆ getDomainInside()

DomainInsideSignature TasGrid::TasmanianSparseGrid::getDomainInside ( ) const
inline

Returns a lambda object that satisfies the TasDREAM::DreamDomain signature.

This method allows for the domain currently set in the grid to be passed to the Tasmanian DREAM sampling templates in place of the inside() object.

◆ removePointsByHierarchicalCoefficient() [1/2]

void TasGrid::TasmanianSparseGrid::removePointsByHierarchicalCoefficient ( double  tolerance,
int  output = -1,
const double *  scale_correction = nullptr 
)

Removes all points from the grid that have relative surplus less than the tolerance.

The output and scale_correction variables have the same effects as in the call to setSurplusRefinement(). The purpose of this call is to reduce the number of points and thus the memory footprint of the grid. As such, points will be removed with no regard of preserving lower completeness or connectivity of the hierarchical graphs; therefore, it is possible that the grid no longer has a valid state with respect to the update and refinement algorithms. Calling loadNeededValues() or any refinement or construction method after the removal of points may lead to undefined behavior; get, evaluate and file I/O methods are safe to call.

Parameters
tolerancethe cut-off tolerance for the point removal.
outputis the output to use for the tolerance test, can be set to -1 to use all outputs.
scale_correctionis the same as in the call to setSurplusRefinement().
Exceptions
std::runtime_errorif the grid is not Local Polynomial.

◆ removePointsByHierarchicalCoefficient() [2/2]

void TasGrid::TasmanianSparseGrid::removePointsByHierarchicalCoefficient ( int  num_new_points,
int  output = -1,
const double *  scale_correction = nullptr 
)

Keeps only the given number of points with largest scaled surpluses.

Similar to removePointsByHierarchicalCoefficient(), but the points are not removed based on a comparison to a tolerance. Instead, only the given number of points is kept so that the remaining points have the largest scaled surplus coefficients.

Parameters
num_new_pointsthe number of points to keep in the grid.
outputis the output to use for the tolerance test, can be set to -1 to use all outputs.
scale_correctionis the same as in the call to setSurplusRefinement().
Exceptions
std::runtime_errorif the grid is not Local Polynomial.

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