43#ifndef __TASMANIAN_SPARSE_GRID_HPP
44#define __TASMANIAN_SPARSE_GRID_HPP
46#include "tsgGridGlobal.hpp"
47#include "tsgGridLocalPolynomial.hpp"
48#include "tsgGridWavelet.hpp"
49#include "tsgGridFourier.hpp"
333 void read(
const char *filename);
343 void read(std::string
const& fname){
read(fname.c_str()); }
392 std::vector<int>
const &anisotropic_weights,
double alpha = 0.0,
double beta = 0.0,
393 const char* custom_filename =
nullptr, std::vector<int>
const &level_limits = std::vector<int>());
402 const int *anisotropic_weights =
nullptr,
double alpha = 0.0,
double beta = 0.0,
403 const char* custom_filename =
nullptr,
const int *level_limits =
nullptr);
410 std::vector<int>
const &anisotropic_weights, std::vector<int>
const &level_limits = std::vector<int>());
417 const int *anisotropic_weights =
nullptr,
const int *level_limits =
nullptr);
437 std::vector<int>
const &anisotropic_weights, std::vector<int>
const &level_limits = std::vector<int>());
446 const int *anisotropic_weights =
nullptr,
const int *level_limits =
nullptr);
502 void makeWaveletGrid(
int dimensions,
int outputs,
int depth,
int order, std::vector<int>
const &level_limits);
510 void makeWaveletGrid(
int dimensions,
int outputs,
int depth,
int order = 1,
const int *level_limits =
nullptr);
522 std::vector<int>
const &anisotropic_weights, std::vector<int>
const &level_limits = std::vector<int>());
530 void makeFourierGrid(
int dimensions,
int outputs,
int depth,
TypeDepth type,
const int* anisotropic_weights =
nullptr,
const int* level_limits =
nullptr);
571 void updateGlobalGrid(
int depth,
TypeDepth type, std::vector<int>
const &anisotropic_weights, std::vector<int>
const &level_limits = std::vector<int>());
592 void updateSequenceGrid(
int depth,
TypeDepth type, std::vector<int>
const &anisotropic_weights, std::vector<int>
const &level_limits = std::vector<int>());
613 void updateFourierGrid(
int depth,
TypeDepth type, std::vector<int>
const &anisotropic_weights, std::vector<int>
const &level_limits = std::vector<int>());
625 void updateGrid(
int depth,
TypeDepth type, std::vector<int>
const &anisotropic_weights, std::vector<int>
const &level_limits = std::vector<int>());
631 void updateGrid(
int depth,
TypeDepth type,
const int *anisotropic_weights =
nullptr,
const int *level_limits =
nullptr);
917 void evaluate(std::vector<double>
const &x, std::vector<double> &y)
const;
957 template<
typename FloatType>
void evaluateBatch(std::vector<FloatType>
const &x, std::vector<FloatType> &y)
const;
1007 template<
typename FloatType>
void evaluateBatchGPU(
const FloatType gpu_x[],
int cpu_num_x, FloatType gpu_y[])
const;
1054 std::vector<double> result;
1064 void differentiate(std::vector<double>
const &x, std::vector<double> &jacobian)
const;
1069 std::vector<double> jacobian;
1081 bool isGlobal()
const{
return base && base->isGlobal(); }
1349 std::vector<int>
const &level_limits, std::vector<double>
const &scale_correction = std::vector<double>());
1362 std::vector<double>
const &scale_correction = std::vector<double>()) {
1436 std::vector<int>
const &level_limits = std::vector<int>());
1471 std::vector<int>
const &level_limits = std::vector<int>(),
1472 std::vector<double>
const &scale_correction = std::vector<double>());
1604 std::vector<double> y;
1685 std::vector<double> integrals;
1871 static int getNumGPUs(){
return AccelerationMeta::getNumGpuDevices(); }
1908 template<
typename FloatType>
1934 template<
typename FloatType>
1938 using EvaluateCallable = std::function<void(std::vector<double>
const&, std::vector<double>&)>;
1947 return [&](std::vector<double>
const &x, std::vector<double> &y)->
void{
evaluateBatch(x, y); };
1962 return [](std::vector<double>
const &)->
bool{
return true; };
1964 if (domain_transform_a.empty()){
1965 return [](std::vector<double>
const &x)->
bool{
1966 for(
auto const &v : x)
if (v < 0.0)
return false;
1971 #if __cplusplus > 201103L
1972 return [dta=domain_transform_a, dims](std::vector<double>
const &x)->
bool{
1973 for(
size_t i=0; i<dims; i++)
if (x[i] < dta[i])
return false;
1977 return [=](std::vector<double>
const &x)->
bool{
1978 for(
size_t i=0; i<dims; i++)
if (x[i] < domain_transform_a[i])
return false;
1984 if (domain_transform_a.empty()){
1986 return [](std::vector<double>
const &x)->
bool{
1987 for(
auto const &v : x)
if ((v < 0.0) || (v > 1.0))
return false;
1991 return [](std::vector<double>
const &x)->
bool{
1992 for(
auto const &v : x)
if ((v < -1.0) || (v > 1.0))
return false;
1998 #if __cplusplus > 201103L
1999 return [dta=domain_transform_a, dtb=domain_transform_b, dims](std::vector<double>
const &x)->
bool{
2000 for(
size_t i=0; i<dims; i++)
2001 if (x[i] < dta[i] or x[i] > dtb[i])
return false;
2005 return [=](std::vector<double>
const &x)->
bool{
2006 for(
size_t i=0; i<dims; i++)
2007 if (x[i] < domain_transform_a[i] or x[i] > domain_transform_b[i])
return false;
2046 #ifndef __TASMANIAN_DOXYGEN_SKIP_INTERNAL
2088 template<
class T>
inline T*
get(){
return dynamic_cast<T*
>(base.get()); }
2096 template<
class T>
inline T
const*
get()
const{
return dynamic_cast<T const*
>(base.get()); }
2107 #ifndef __TASMANIAN_DOXYGEN_SKIP_INTERNAL
2188 template<
typename T>
2246 std::unique_ptr<AccelerationContext> acceleration;
2248 std::unique_ptr<BaseCanonicalGrid> base;
2250 std::vector<double> domain_transform_a, domain_transform_b;
2251 std::vector<int> conformal_asin_power;
2252 std::vector<int> llimits;
2254 bool using_dynamic_construction;
2256 mutable std::unique_ptr<AccelerationDomainTransform> acc_domain;
2286inline TasmanianSparseGrid
2288 std::vector<int>
const &anisotropic_weights = std::vector<int>(),
double alpha = 0.0,
double beta = 0.0,
2289 const char* custom_filename =
nullptr, std::vector<int>
const &level_limits = std::vector<int>()){
2291 grid.
makeGlobalGrid(dimensions, outputs, depth, type, rule, anisotropic_weights, alpha, beta, custom_filename, level_limits);
2301inline TasmanianSparseGrid
2303 std::vector<int>
const &anisotropic_weights = std::vector<int>(), std::vector<int>
const &level_limits = std::vector<int>()){
2305 grid.
makeSequenceGrid(dimensions, outputs, depth, type, rule, anisotropic_weights, level_limits);
2315inline TasmanianSparseGrid
2328inline TasmanianSparseGrid
2329makeWaveletGrid(
int dimensions,
int outputs,
int depth,
int order = 1, std::vector<int>
const &level_limits = std::vector<int>()){
2331 grid.
makeWaveletGrid(dimensions, outputs, depth, order, level_limits);
2341inline TasmanianSparseGrid
2343 std::vector<int>
const &anisotropic_weights = std::vector<int>(), std::vector<int>
const &level_limits = std::vector<int>()){
2345 grid.
makeFourierGrid(dimensions, outputs, depth, type, anisotropic_weights, level_limits);
2359 grid.
read(filename);
2384 grid.
copyGrid(source, outputs_begin, outputs_end);
Class providing manipulation of custom tabulated rules, file I/O and structured access to the points,...
Definition tsgCoreOneDimensional.hpp:66
Generic 2D data structure divided into contiguous strips of fixed length (similar to a matrix).
Definition tsgIndexSets.hpp:104
Template class that wraps around a single GPU array, providing functionality that mimics std::vector.
Definition tsgAcceleratedDataStructures.hpp:95
The master-class that represents an instance of a Tasmanian sparse grid.
Definition TasmanianSparseGrid.hpp:293
bool isGlobal() const
Returns true if the grid is of type global, false otherwise.
Definition TasmanianSparseGrid.hpp:1081
void setAnisotropicRefinement(TypeDepth type, int min_growth, int output, const int *level_limits=nullptr)
Overload using raw-array.
std::vector< int > estimateAnisotropicCoefficients(TypeDepth type, int output) const
Estimate the anisotropic rates of coefficient decay for different direction.
Definition TasmanianSparseGrid.hpp:1272
int getGPUID() const
Returns the currently set CUDA device.
Definition TasmanianSparseGrid.hpp:1860
void getDifferentiationWeights(const std::vector< double > &x, std::vector< double > &weights) const
Overload that uses the vector as a parameter.
T * get()
Return a reference to the internal grid.
Definition TasmanianSparseGrid.hpp:2088
std::vector< double > getNeededPoints() const
Return the points that require model values.
Definition TasmanianSparseGrid.hpp:695
void mapTransformedToCanonical(int num_dimensions, int num_points, TypeOneDRule rule, FloatType x[]) const
Apply the inverse map from transformed to canonical points.
const T * formCanonicalPointsGPU(const T *gpu_x, int num_x, GpuVector< T > &gpu_x_temp) const
Returns a CUDA raw-array with the canonical points, linear transform only.
int getNumPoints() const
Returns getNumLoaded() if positive, otherwise returns getNumNeeded(), see getPoints().
Definition TasmanianSparseGrid.hpp:661
void getDomainTransform(double a[], double b[]) const
Returns the values of the two vectors used to call setDomainTransform().
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.
void updateSequenceGrid(int depth, TypeDepth type, const int *anisotropic_weights=nullptr, const int *level_limits=nullptr)
Overload using raw-arrays.
void getLoadedPoints(std::vector< double > &x) const
Overload that accepts the vector as a parameter.
Definition TasmanianSparseGrid.hpp:678
TypeOneDRule getRule() const
Return the underlying TasGrid::TypeOneDRule that gives the points and basis functions.
int getNumOutputs() const
Return the outputs of the grid, i.e., number of model outputs.
Definition TasmanianSparseGrid.hpp:644
void setConformalTransformASIN(std::vector< int > const &truncation)
Set conformal transformation using truncated Maclaurin series of the arcsin() function.
void formTransformedPoints(int num_points, double x[]) const
Applies both linear and non-linear transformation to the canonical points.
void finishConstruction()
End the procedure, clears flags and unused constructed points, can go back to using regular refinemen...
T const * get() const
Overload using const.
Definition TasmanianSparseGrid.hpp:2096
static int getVersionMinor()
Return the library minor version.
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.
void write(const char *filename, bool binary=mode_binary) const
Write the grid to the given filename using either binary or ASCII format.
void readBinary(std::istream &ifs)
Read the grid from a stream using binary format.
void getLoadedPoints(double x[]) const
Overload using raw-array, writes the loaded points to the first getNumLoaded() times getNumDimensions...
const double * getLoadedValues() const
Returns the model values that have been loaded in the gird.
Definition TasmanianSparseGrid.hpp:899
void getDifferentiationWeights(const double x[], double weights[]) const
Overload that uses raw-array, does not check the array size.
double getBeta() const
Return the beta parameter in the call to makeGlobalGrid(), or return 0 if the grid is not Global.
Definition TasmanianSparseGrid.hpp:636
void getNeededPoints(std::vector< double > &x) const
Overload that accepts the vector as a parameter.
Definition TasmanianSparseGrid.hpp:702
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.
static bool isAccelerationAvailable(TypeAcceleration acc)
Returns whether a specific mode can be enabled.
static bool isHipEnabled()
Returns true if compiled with HIP support, e.g., Tasmanian_ENABLE_HIP=ON.
void evaluateFast(const FloatType x[], FloatType y[]) const
Equivalent to evaluate() with enabled acceleration or evaluateBatch() with a batch of one point.
Definition TasmanianSparseGrid.hpp:1021
std::vector< int > getLevelLimits() const
Return the currently set level limits.
Definition TasmanianSparseGrid.hpp:1207
void differentiate(const double x[], double jacobian[]) const
Overload that uses a raw-array.
void mapConformalTransformedToCanonical(int num_dimensions, int num_points, Data2D< FloatType > &x) const
Applies the inverse non-linear transformation to the points.
static bool isDpcppEnabled()
Returns true if compiled with DPC++ support, e.g., Tasmanian_ENABLE_DPCPP=ON.
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.
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.
void loadNeededValues(const double *vals)
Overload that uses a raw-array, does not check the array size.
void evaluate(const double x[], double y[]) const
Overload that uses raw-arrays, does not check the array size.
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.
std::vector< double > getHierarchicalSupport() const
Returns the support of the hierarchical functions.
void setHierarchicalCoefficients(const double c[])
Overload that uses raw-arrays.
Definition TasmanianSparseGrid.hpp:1569
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.
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.
std::vector< double > getQuadratureWeights() const
Returns a vector of size getNumPoints() of the quadrature weights of the grid.
Definition TasmanianSparseGrid.hpp:747
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).
void readAscii(std::istream &ifs)
Read the grid from a stream using ASCII format.
const FloatType * formCanonicalPoints(const FloatType *x, Data2D< FloatType > &x_temp, int num_x) const
Returns a raw-array with the canonical points, combines both transformations.
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.
TasmanianSparseGrid & operator=(TasmanianSparseGrid const &source)
Copy assignment, note that the selected acceleration mode is not copied, acceleration is reset to the...
void getNeededPoints(double *x) const
Overload using raw-array, writes the loaded points to the first getNumNeeded() times getNumDimensions...
void evaluateBatch(const float x[], int num_x, float y[]) const
Overload using single precision and GPU/CUDA acceleration.
static const char * getLicense()
Return a hard-coded character string with a brief statement of the license.
void makeLocalPolynomialGrid(int dimensions, int outputs, int depth, int order=1, TypeOneDRule rule=rule_localp, const int *level_limits=nullptr)
Overload using raw-arrays.
std::vector< double > getAnisotropicRefinement(TypeDepth type, int min_growth, int output, const std::vector< int > &level_limits)
Call setAnisotropicRefinement() and then getNeededPoints().
Definition TasmanianSparseGrid.hpp:1248
void getQuadratureWeights(std::vector< double > &weights) const
Overload that accepts the vector as a parameter.
Definition TasmanianSparseGrid.hpp:753
static bool isCudaEnabled()
Returns true if compiled with CUDA support, e.g., Tasmanian_ENABLE_CUDA=ON.
void updateGrid(int depth, TypeDepth type, const int *anisotropic_weights=nullptr, const int *level_limits=nullptr)
Overload using raw-arrays.
void integrateHierarchicalFunctions(double integrals[]) const
Overload using raw-arrays.
~TasmanianSparseGrid()=default
Destructor, releases all resources.
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 write(std::string const &fname, bool binary=mode_binary) const
Overload that works directly with std::string.
Definition TasmanianSparseGrid.hpp:341
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.
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 ru...
double getQuadratureScale(int num_dimensions, TypeOneDRule rule) const
Returns the quadrature scale factor associated with the linear transform.
std::vector< double > diffCanonicalTransform() const
Calculates the Jacobian matrix (derivative) of the transform defined by formCanonicalPoints()....
void getPoints(std::vector< double > &x) const
Overload that accepts the vector as a parameter.
Definition TasmanianSparseGrid.hpp:730
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.
bool isSetDomainTransfrom() const
Returns true if a linear domain transformation has been set, false otherwise.
void getDomainTransform(std::vector< double > &a, std::vector< double > &b) const
Returns the two vectors used to call setDomainTransform().
void integrateHierarchicalFunctions(std::vector< double > &integrals) const
Overload where the vector is passes as a parameter.
Definition TasmanianSparseGrid.hpp:1695
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.
void writeAscii(std::ostream &ofs) const
Write the grid to a stream using ASCII format.
std::vector< double > integrate() const
Overload that returns a vector.
Definition TasmanianSparseGrid.hpp:1053
int getNumNeeded() const
Return the number of points that should be provided to the next call of loadNeededValues().
Definition TasmanianSparseGrid.hpp:659
void read(std::string const &fname)
Overload that works directly with std::string.
Definition TasmanianSparseGrid.hpp:343
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().
Definition TasmanianSparseGrid.hpp:1361
void clearDomainTransform()
Removes the domain transformation.
std::vector< double > getPoints() const
Returns the loaded points if any, otherwise returns the needed points.
Definition TasmanianSparseGrid.hpp:724
void evaluateHierarchicalFunctions(const double x[], int num_x, double y[]) const
Array overload, the inputs must have pre-allocated and correct size.
TypeAcceleration getAccelerationType() const
Returns the current effective acceleration mode.
Definition TasmanianSparseGrid.hpp:1834
void clear()
Reset the grid to empty.
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.
void setRocSparseHandle(void *handle)
Takes a user provided cuSparse handle.
void clearLevelLimits()
Removes the currently set level limits.
Definition TasmanianSparseGrid.hpp:1200
int getNumLoaded() const
Return the number of points already associated with model values via loadNeededValues().
Definition TasmanianSparseGrid.hpp:657
void evaluateFast(std::vector< FloatType > const &x, std::vector< FloatType > &y) const
Alias to evaluateBatch().
Definition TasmanianSparseGrid.hpp:1027
std::function< void(std::vector< double > const &, std::vector< double > &)> EvaluateCallable
Signature compatible with TasDREAM::DreamPDF, TasDREAM::DreamModel amd TasDREAM::DreamMergedLikelyMod...
Definition TasmanianSparseGrid.hpp:1938
void estimateAnisotropicCoefficients(TypeDepth type, int output, std::vector< int > &weights) const
Overload that writes the result to a parameter.
std::vector< int > getConformalTransformASIN() const
Fills the array with the values of the set transformation.
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.
void clearConformalTransform()
Removes any currently set transformation.
std::vector< double > differentiate(std::vector< double > const &x) const
Same as TasmanianSparseGrid::differentiate() but returns the jacobian.
Definition TasmanianSparseGrid.hpp:1068
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.
static int getNumGPUs()
Return the number of visible CUDA devices.
Definition TasmanianSparseGrid.hpp:1871
TasmanianSparseGrid(const TasmanianSparseGrid &source)
Copy constructor, note that the selected acceleration mode is not copied, acceleration is reset to th...
void beginConstruction()
Begin a dynamic construction procedure.
static std::string getGPUName(int gpu)
Return the CUDA device name.
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 getNeededPoin...
Definition TasmanianSparseGrid.hpp:1307
bool isSequence() const
Returns true if the grid is of type sequence, false otherwise.
Definition TasmanianSparseGrid.hpp:1083
void getQuadratureWeights(double weights[]) const
Overload that accepts the raw array as an input.
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.
void setCuSolverHandle(void *handle)
Takes a user provided cuSparse handle.
std::vector< int > getGlobalPolynomialSpace(bool interpolation) const
Returns the powers of the polynomial that span the effective basis, Global and Sequence grids only.
void mapConformalCanonicalToTransformed(int num_dimensions, int num_points, double x[]) const
Applies the non-linear transformation to the points.
void integrate(std::vector< double > &q) const
Computes the integral of each model output over the sparse grid domain.
bool isFourier() const
Returns true if the grid is of type Fourier, false otherwise.
Definition TasmanianSparseGrid.hpp:1089
static int getVersionMajor()
Return the library major version.
void read(std::istream &ifs, bool binary=mode_binary)
Read the grid from the given stream ifs using either binary or ASCII format.
void getInterpolationWeights(const std::vector< double > &x, std::vector< double > &weights) const
Overload that uses the vector as a parameter.
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.
const char * getCustomRuleDescription() const
Return the character string that is the description of the user-provided tabulated rule.
static int getGPUMemory(int gpu)
Return the available device memory, in units of MB.
std::vector< double > integrateHierarchicalFunctions() const
Returns the integrals of the hierarchical basis functions.
Definition TasmanianSparseGrid.hpp:1684
void printStats(std::ostream &os=std::cout) const
Prints short human-readable text describing the grid properties.
void clearRefinement()
Remove all needed points from the grid.
void evaluateSparseHierarchicalFunctionsStatic(const double x[], int num_x, int pntr[], int indx[], double vals[]) const
Assumes pre-allocated arrays, otherwise identical to evaluateSparseHierarchicalFunctions()
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.
void mergeRefinement()
Merges the loaded and needed points into a single grid, resets all loaded values to zero.
void loadNeededPoints(const double *vals)
Overload that uses a raw-array, does not check the array size.
Definition TasmanianSparseGrid.hpp:890
TasmanianSparseGrid & operator=(TasmanianSparseGrid &&source)=default
Move assignment, the selected acceleration mode is also carried over.
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.
std::vector< double > getDifferentiationWeights(double const x[]) const
Overload that uses raw-array, does not check the array size.
Definition TasmanianSparseGrid.hpp:835
void copyGrid(TasmanianSparseGrid const &source, int outputs_begin=0, int outputs_end=-1)
Overload using pass-by-reference as opposed to a pointer.
Definition TasmanianSparseGrid.hpp:555
void updateFourierGrid(int depth, TypeDepth type, const int *anisotropic_weights=nullptr, const int *level_limits=nullptr)
Overload using raw-arrays.
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 ...
void getPoints(double x[]) const
Overload that accepts the raw array as an input.
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 (d...
static const char * getGitCommitHash()
Return the git hash string, will use a placeholder if the git command was not found on compile time o...
std::vector< double > getInterpolationWeights(double const x[]) const
Overload that uses raw-array, does not check the array size.
Definition TasmanianSparseGrid.hpp:795
void enableAcceleration(TypeAcceleration acc)
Change the current acceleration mode to the one specified.
static const char * getVersion()
Return a hard-coded character string with the version in format "Major.Minor".
void getHierarchicalCoefficientsStatic(double *coeff) const
Copies the coefficients to the pre-allocated array, intended for internal use. .
Definition TasmanianSparseGrid.hpp:1542
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.
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).
bool isSetConformalTransformASIN() const
Returns true if conformal transformation has been set.
void evaluateBatch(const double x[], int num_x, double y[]) const
Overload that uses raw-arrays.
void read(const char *filename)
Read the grid from the given filename, automatically detect the format.
void makeFourierGrid(int dimensions, int outputs, int depth, TypeDepth type, const int *anisotropic_weights=nullptr, const int *level_limits=nullptr)
Overload using raw-arrays.
void setAnisotropicRefinement(TypeDepth type, int min_growth, int output, const std::vector< int > &level_limits)
Set refinement using anisotropic estimate for the optimal points.
void writeBinary(std::ostream &ofs) const
Write the grid to a stream using binary format.
const int * getPointsIndexes() const
Return a reference to the multi-indexes of the loaded points.
double getAlpha() const
Return the alpha parameter in the call to makeGlobalGrid(), or return 0 if the grid is not Global.
Definition TasmanianSparseGrid.hpp:634
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().
void favorSparseAcceleration(bool favor)
Set the preferred back-end algorithm for Local Polynomial grids.
static bool isOpenMPEnabled()
Returns true if compiled with OpenMP support, e.g., Tasmanian_ENABLE_OPENMP=ON.
void setSurplusRefinement(double tolerance, TypeRefinement criteria, int output=-1, const int *level_limits=nullptr, const double *scale_correction=nullptr)
Overload that uses raw-arrays.
void setSycleQueue(void *queue)
Takes a user provided sycl::queue handle.
void setSurplusRefinement(double tolerance, int output, const int *level_limits=nullptr)
Overload that uses array for the level limits.
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)
static const char * getCmakeCxxFlags()
Return the CMAKE_BUILD_TYPE and CMAKE_CXX_FLAGS used in the configuration.
DomainInsideSignature getDomainInside() const
Returns a lambda object that satisfies the TasDREAM::DreamDomain signature.
Definition TasmanianSparseGrid.hpp:1959
AccelerationContext const * getAccelerationContext() const
Allows the addon methods to use the acceleration context.
Definition TasmanianSparseGrid.hpp:2103
void setDomainTransform(std::vector< double > const &a, std::vector< double > const &b)
Set a linear domain transformation.
bool isUsingConstruction() const
Returns true if the dynamic construction procedure has been initialized, false otherwise.
Definition TasmanianSparseGrid.hpp:1403
void setCuBlasHandle(void *handle)
Takes a user provided cuBlas handle.
const int * getNeededIndexes() const
Return a reference to the multi-indexes of the needed points.
void evaluateHierarchicalFunctions(std::vector< double > const &x, std::vector< double > &y) const
Computes the values of the hierarchical function basis at the specified points.
bool empty() const
Returns true if the grid is empty (no type), false otherwise.
Definition TasmanianSparseGrid.hpp:1093
TasmanianSparseGrid()
Default constructor, creates and empty grid.
void makeWaveletGrid(int dimensions, int outputs, int depth, int order=1, const int *level_limits=nullptr)
Overload using raw-arrays.
void loadNeededPoints(std::vector< double > const &vals)
Alias of loadNeededValues().
Definition TasmanianSparseGrid.hpp:884
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.
void setCuSparseHandle(void *handle)
Takes a user provided cuSparse handle.
void evaluateBatchGPU(const FloatType gpu_x[], int cpu_num_x, FloatType gpu_y[]) const
Overload that uses GPU raw-arrays.
void mapCanonicalToTransformed(int num_dimensions, int num_points, TypeOneDRule rule, double x[]) const
Maps canonical points to the transformed equivalent.
void setHierarchicalCoefficients(const std::vector< double > &c)
Overwrites the current set of coefficients (and loaded values) with the ones provided.
void updateGlobalGrid(int depth, TypeDepth type, const int *anisotropic_weights=nullptr, const int *level_limits=nullptr)
Overload using raw-arrays.
void setDomainTransform(const double a[], const double b[])
Overload using raw-arrays.
const double * getHierarchicalCoefficients() const
Return a reference to the internal data-structure that stores the hierarchical coefficients.
TasmanianSparseGrid(TasmanianSparseGrid &&source)=default
Move constructor, the selected acceleration mode is also carried over.
std::vector< double > getLoadedPoints() const
Return the points already associated with model values.
Definition TasmanianSparseGrid.hpp:671
std::function< bool(std::vector< double > const &)> DomainInsideSignature
Signature of the domain inside lambda, identical to TasDREAM::DreamDomain.
Definition TasmanianSparseGrid.hpp:1951
int getNumDimensions() const
Return the dimensions of the grid, i.e., number of model inputs.
Definition TasmanianSparseGrid.hpp:642
void enableAcceleration(TypeAcceleration acc, int new_gpu_id)
Combines calls the enableAcceleration(), getGPUID() and allows for user provided handles.
bool isLocalPolynomial() const
Returns true if the grid is of type local polynomial, false otherwise.
Definition TasmanianSparseGrid.hpp:1085
int evaluateSparseHierarchicalFunctionsGetNZ(const double x[], int num_x) const
Count the number of non-zeros in a call to evaluateSparseHierarchicalFunctions().
void getInterpolationWeights(const double x[], double weights[]) const
Overload that uses raw-array, does not check the array size.
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.
bool isEmpty() const
Returns true if the grid is empty (no type), false otherwise.
Definition TasmanianSparseGrid.hpp:1091
bool isWavelet() const
Returns true if the grid is of type wavelet, false otherwise.
Definition TasmanianSparseGrid.hpp:1087
void setRocBlasHandle(void *handle)
Takes a user provided cuBlas handle.
void mapConformalWeights(int num_dimensions, int num_points, double weights[]) const
Computes the quadrature weight correction for the conformal map.
void loadNeededValues(std::vector< double > const &vals)
Provides the values of the model outputs at the needed points, or overwrites the currently loaded poi...
void loadConstructedPoints(std::vector< double > const &x, std::vector< double > const &y)
Add pairs of points with associated model values.
void setGPUID(int new_gpu_id)
Select the current CUDA device.
std::vector< double > evaluateHierarchicalFunctions(std::vector< double > const &x) const
Overload that returns the result.
Definition TasmanianSparseGrid.hpp:1603
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.
int getOrder() const
Return the order parameter in the call to makeLocalPolynomialGrid() or makeWaveletGrid(),...
Definition TasmanianSparseGrid.hpp:638
void integrate(double q[]) const
Overload that uses a raw-array.
TypeOneDRule
Used to specify the one dimensional family of rules that induces the sparse grid.
Definition tsgEnumerates.hpp:285
TypeAcceleration
Modes of acceleration.
Definition tsgEnumerates.hpp:551
TypeDepth
Used by Global Sequence and Fourier grids, indicates the selection criteria.
Definition tsgEnumerates.hpp:203
TypeRefinement
Refinement strategy for local polynomial and wavelet grids.
Definition tsgEnumerates.hpp:425
@ rule_gausshermite
Non-nested rule optimized for integral of the form .
Definition tsgEnumerates.hpp:355
@ rule_gausslaguerreodd
Same as rule_gausslaguerre but using only odd levels, partially mitigates the non-nested issues.
Definition tsgEnumerates.hpp:353
@ rule_localp
Nested rule with a hierarchy of uniformly distributed nodes and functions with compact support.
Definition tsgEnumerates.hpp:362
@ rule_gausslaguerre
Non-nested rule optimized for integral of the form .
Definition tsgEnumerates.hpp:351
@ rule_gausshermiteodd
Same as rule_gausshermite but using only odd levels, partially mitigates the non-nested issues.
Definition tsgEnumerates.hpp:357
TasmanianSparseGrid makeEmpty()
Returns an empty sparse grid.
Definition TasmanianSparseGrid.hpp:2270
TasmanianSparseGrid makeWaveletGrid(int dimensions, int outputs, int depth, int order=1, std::vector< int > const &level_limits=std::vector< int >())
Factory method, creates a new grid and calls TasmanianSparseGrid::makeWaveletGrid().
Definition TasmanianSparseGrid.hpp:2329
TasmanianSparseGrid makeSequenceGrid(int dimensions, int outputs, int depth, TypeDepth type, TypeOneDRule rule, std::vector< int > const &anisotropic_weights=std::vector< int >(), std::vector< int > const &level_limits=std::vector< int >())
Factory method, creates a new grid and calls TasmanianSparseGrid::makeSequenceGrid().
Definition TasmanianSparseGrid.hpp:2302
TasmanianSparseGrid makeGlobalGrid(int dimensions, int outputs, int depth, TypeDepth type, TypeOneDRule rule, std::vector< int > const &anisotropic_weights=std::vector< int >(), double alpha=0.0, double beta=0.0, const char *custom_filename=nullptr, std::vector< int > const &level_limits=std::vector< int >())
Factory method, creates a new grid and calls TasmanianSparseGrid::makeGlobalGrid().
Definition TasmanianSparseGrid.hpp:2287
TasmanianSparseGrid readGrid(const char *filename)
Factory method, creates a new grid and calls TasmanianSparseGrid::read().
Definition TasmanianSparseGrid.hpp:2357
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:2342
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:2316
constexpr bool mode_binary
Constant allowing for more expressive selection of ascii and binary mode in IO methods.
Definition tsgIOHelpers.hpp:68
TasmanianSparseGrid copyGrid(TasmanianSparseGrid const &source, int outputs_begin=0, int outputs_end=-1)
Returns a grid that is a copy of the source.
Definition TasmanianSparseGrid.hpp:2382
size_t size_mult(IntA a, IntB b)
Converts two integer-like variables to size_t and returns the product..
Definition tsgUtils.hpp:82
Encapsulates the Tasmanian Sparse Grid module.
Definition TasmanianSparseGrid.hpp:68
Wrapper class around GPU device ID, acceleration type and GpuEngine.
Definition tsgAcceleratedDataStructures.hpp:576