31 #ifndef __TASMANIAN_ADDONS_LOADUNSTRUCTURED_HPP
32 #define __TASMANIAN_ADDONS_LOADUNSTRUCTURED_HPP
93 template<
typename scalar_type>
99 int num_equations = (tolerance > 0.0) ? num_data + num_points : num_data;
104 if (tolerance > 0.0){
105 double correction = std::sqrt(tolerance);
106 constexpr
long long esize = (std::is_same<scalar_type, double>::value) ? 1 : 2;
107 long long num_total =
static_cast<long long>(num_points) *
static_cast<long long>(num_points) * esize;
111 long long stride =
static_cast<long long>(num_points + 1) * esize;
116 reinterpret_cast<double*
>(model_values +
Utils::size_mult(num_data, num_outputs)));
120 basis_matrix.
data(), num_outputs, model_values);
134 template<
typename scalar_type>
139 int num_equations = (tolerance > 0.0) ? num_data + num_points : num_data;
147 int mem_usage = 268435456 / ((std::is_same<double, scalar_type>::value) ? 1 : 2);
149 int num_batch = std::min(num_data, mem_usage / num_points);
150 if (num_batch < 128){
157 for(
int i = 0; i < num_data; i += num_batch){
158 int num_this_batch = std::min(num_batch, num_data - i);
169 if (tolerance > 0.0){
170 double correction = std::sqrt(tolerance);
172 basis_matrix.
getStrip(i + num_data)[i] = correction;
175 auto icoeff = coefficients.
begin();
177 *icoeff++ = model_values[i];
194 template<
typename scalar_type>
198 if (grid.
empty())
throw std::runtime_error(
"Cannot use loadUnstructuredDataL2() with an empty grid.");
204 throw std::runtime_error(
"The loadUnstructuredDataL2() method cannot be used with acceleration mode accel_none.");
209 int num_equations = (tolerance > 0.0) ? num_data + num_points : num_data;
215 GpuVector<double> gpu_points(acceleration, num_dimensions, num_data, data_points);
218 generateCoefficientsGPU<scalar_type>(gpu_points.
data(), num_data, gpu_values.
data(), tolerance, grid);
221 return generateCoefficients<scalar_type>(data_points, num_data, model_values, tolerance, grid);
226 if (std::is_same<scalar_type, std::complex<double>>::value){
228 auto icoeff = coefficients.
begin();
230 real_coeffs[i] = std::real(*icoeff++);
231 icoeff = coefficients.
begin();
233 real_coeffs[i] = std::imag(*icoeff++);
285 loadUnstructuredDataL2tmpl<std::complex<double>>(data_points, num_data, model_values, tolerance, grid);
287 loadUnstructuredDataL2tmpl<double>(data_points, num_data, model_values, tolerance, grid);
301 if (grid.
empty())
throw std::runtime_error(
"Cannot use loadUnstructuredDataL2() with an empty grid.");
302 int num_data =
static_cast<int>(data_points.size() / grid.
getNumDimensions());
304 throw std::runtime_error(
"In loadUnstructuredDataL2(), provided more points than data.");
Generic 2D data structure divided into contiguous strips of fixed length (similar to a matrix).
Definition: tsgIndexSets.hpp:104
std::vector< T >::iterator begin()
Returns an iterator to the beginning of the internal data.
Definition: tsgIndexSets.hpp:167
T * data()
Returns a reference to the internal data.
Definition: tsgIndexSets.hpp:153
T * getStrip(int i)
Returns a reference to the i-th strip.
Definition: tsgIndexSets.hpp:141
Template class that wraps around a single GPU array, providing functionality that mimics std::vector.
Definition: tsgAcceleratedDataStructures.hpp:95
void unload(AccelerationContext const *acc, std::vector< T > &cpu_data) const
Copy the data from the GPU array to cpu_data, the cpu_data will be resized and overwritten.
Definition: tsgAcceleratedDataStructures.hpp:196
T * data()
Get a reference to the GPU array, which an be used as input to GPU libraries and kernels.
Definition: tsgAcceleratedDataStructures.hpp:145
void load(AccelerationContext const *acc, const std::vector< T > &cpu_data)
Copy the content of cpu_data to the GPU device, all pre-existing data is deleted and the vector is re...
Definition: tsgAcceleratedDataStructures.hpp:157
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
int getNumOutputs() const
Return the outputs of the grid, i.e., number of model outputs.
Definition: TasmanianSparseGrid.hpp:644
AccelerationContext const * getAccelerationContext() const
Allows the addon methods to use the acceleration context.
Definition: TasmanianSparseGrid.hpp:2088
int getNumNeeded() const
Return the number of points that should be provided to the next call of loadNeededValues().
Definition: TasmanianSparseGrid.hpp:659
bool isFourier() const
Returns true if the grid is of type Fourier, false otherwise.
Definition: TasmanianSparseGrid.hpp:1089
static int getGPUMemory(int gpu)
Return the available device memory, in units of MB.
void mergeRefinement()
Merges the loaded and needed points into a single grid, resets all loaded values to zero.
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).
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
void setHierarchicalCoefficients(const std::vector< double > &c)
Overwrites the current set of coefficients (and loaded values) with the ones provided.
int getNumDimensions() const
Return the dimensions of the grid, i.e., number of model inputs.
Definition: TasmanianSparseGrid.hpp:642
bool isLocalPolynomial() const
Returns true if the grid is of type local polynomial, false otherwise.
Definition: TasmanianSparseGrid.hpp:1085
bool isWavelet() const
Returns true if the grid is of type wavelet, false otherwise.
Definition: TasmanianSparseGrid.hpp:1087
int getOrder() const
Return the order parameter in the call to makeLocalPolynomialGrid() or makeWaveletGrid(),...
Definition: TasmanianSparseGrid.hpp:638
@ accel_none
Usually the slowest mode, uses only OpenMP multi-threading, but optimized for memory and could be the...
Definition: tsgEnumerates.hpp:553
@ accel_gpu_magma
Same the CUDA option but uses the UTK MAGMA library for the linear algebra operations.
Definition: tsgEnumerates.hpp:563
@ accel_gpu_cuda
Similar to the cuBLAS option but also uses a set of Tasmanian custom GPU kernels.
Definition: tsgEnumerates.hpp:561
void load_n(AccelerationContext const *acc, T const *cpu_data, size_t num_entries, T *gpu_data)
Similar to copy_n, copies the data from the CPU to the GPU.
void fillDataGPU(AccelerationContext const *acc, double value, long long N, long long stride, double data[])
Fills the data with the provided real number at the given stride.
void generateCoefficientsGPU(double const data_points[], int num_data, scalar_type model_values[], double tolerance, TasmanianSparseGrid &grid)
Generates the coefficients and solves the least-squares problem.
Definition: tsgLoadUnstructuredPoints.hpp:94
Data2D< scalar_type > generateCoefficients(double const data_points[], int num_data, double const model_values[], double tolerance, TasmanianSparseGrid &grid)
Generates the coefficients and solves the least-squares problem.
Definition: tsgLoadUnstructuredPoints.hpp:135
void loadUnstructuredDataL2tmpl(double const data_points[], int num_data, double const model_values[], double tolerance, TasmanianSparseGrid &grid)
Template implementation that handles the case of Fourier grids vs. all other types.
Definition: tsgLoadUnstructuredPoints.hpp:195
bool hasGPUBasis(TasmanianSparseGrid const &grid)
Returns true if the grid has a GPU algorithm for computing the hierarchical basis.
Definition: tsgLoadUnstructuredPoints.hpp:76
void loadUnstructuredDataL2(double const data_points[], int num_data, double const model_values[], double tolerance, TasmanianSparseGrid &grid)
Construct a sparse grid surrogate using a least-squares fit.
Definition: tsgLoadUnstructuredPoints.hpp:282
size_t size_mult(IntA a, IntB b)
Converts two integer-like variables to size_t and returns the product..
Definition: tsgUtils.hpp:82
void solvesLeastSquares(AccelerationContext const *acceleration, int n, int m, scalar_type A[], int nrhs, scalar_type B[])
Least squares solver, operates on multiple right-hand sides and row-major matrices.
void solvesLeastSquaresGPU(AccelerationContext const *acceleration, int n, int m, scalar_type A[], int nrhs, scalar_type B[])
Overload that accepts arrays on the GPU device.
Encapsulates the Tasmanian Sparse Grid module.
Definition: TasmanianSparseGrid.hpp:68
Wrapper class around GPU device ID, acceleration type and GpuEngine.
Definition: tsgAcceleratedDataStructures.hpp:576
TypeAcceleration mode
The current active acceleration mode.
Definition: tsgAcceleratedDataStructures.hpp:607
int device
If using a GPU acceleration mode, holds the active device.
Definition: tsgAcceleratedDataStructures.hpp:611
void setDevice() const
Set default device.
Definition: tsgAcceleratedDataStructures.hpp:717
Templates for non-adaptive sampling from lambda models.