Doxygen 1.9.8
Toolkit for Adaptive Stochastic Modeling and Non-Intrusive ApproximatioN: Tasmanian v8.2
 
Loading...
Searching...
No Matches
tsgLoadUnstructuredPoints.hpp
Go to the documentation of this file.
1/*
2 * Copyright (c) 2017, Miroslav Stoyanov
3 *
4 * This file is part of
5 * Toolkit for Adaptive Stochastic Modeling And Non-Intrusive ApproximatioN: TASMANIAN
6 *
7 * Redistribution and use in source and binary forms, with or without modification, are permitted provided that the following conditions are met:
8 *
9 * 1. Redistributions of source code must retain the above copyright notice, this list of conditions and the following disclaimer.
10 *
11 * 2. Redistributions in binary form must reproduce the above copyright notice, this list of conditions
12 * and the following disclaimer in the documentation and/or other materials provided with the distribution.
13 *
14 * 3. Neither the name of the copyright holder nor the names of its contributors may be used to endorse
15 * or promote products derived from this software without specific prior written permission.
16 *
17 * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES,
18 * INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED.
19 * IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY,
20 * OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA,
21 * OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY,
22 * OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
23 *
24 * UT-BATTELLE, LLC AND THE UNITED STATES GOVERNMENT MAKE NO REPRESENTATIONS AND DISCLAIM ALL WARRANTIES, BOTH EXPRESSED AND IMPLIED.
25 * THERE ARE NO EXPRESS OR IMPLIED WARRANTIES OF MERCHANTABILITY OR FITNESS FOR A PARTICULAR PURPOSE, OR THAT THE USE OF THE SOFTWARE WILL NOT INFRINGE ANY PATENT,
26 * COPYRIGHT, TRADEMARK, OR OTHER PROPRIETARY RIGHTS, OR THAT THE SOFTWARE WILL ACCOMPLISH THE INTENDED RESULTS OR THAT THE SOFTWARE OR ITS USE WILL NOT RESULT IN INJURY OR DAMAGE.
27 * THE USER ASSUMES RESPONSIBILITY FOR ALL LIABILITIES, PENALTIES, FINES, CLAIMS, CAUSES OF ACTION, AND COSTS AND EXPENSES, CAUSED BY, RESULTING FROM OR ARISING OUT OF,
28 * IN WHOLE OR IN PART THE USE, STORAGE OR DISPOSAL OF THE SOFTWARE.
29 */
30
31#ifndef __TASMANIAN_ADDONS_LOADUNSTRUCTURED_HPP
32#define __TASMANIAN_ADDONS_LOADUNSTRUCTURED_HPP
33
46
67namespace TasGrid{
68
76inline bool hasGPUBasis(TasmanianSparseGrid const &grid){
77 return AccelerationMeta::isAvailable(accel_gpu_cuda)
78 and not (grid.isLocalPolynomial() and ((grid.getOrder() < 0) or (grid.getOrder() > 2)))
79 and not (grid.isWavelet() and grid.getOrder() == 3);
80}
81
93template<typename scalar_type>
94void generateCoefficientsGPU(double const data_points[], int num_data, scalar_type model_values[],
95 double tolerance, TasmanianSparseGrid &grid){
96 AccelerationContext const *acceleration = grid.getAccelerationContext();
97 int num_outputs = grid.getNumOutputs();
98 int num_points = grid.getNumPoints();
99 int num_equations = (tolerance > 0.0) ? num_data + num_points : num_data;
100
101 GpuVector<scalar_type> basis_matrix(acceleration, num_equations, num_points);
102 grid.evaluateHierarchicalFunctionsGPU(data_points, num_data, reinterpret_cast<double*>(basis_matrix.data()));
103
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;
108 TasGpu::fillDataGPU(acceleration, 0.0, num_total, 1,
109 reinterpret_cast<double*>(basis_matrix.data() + Utils::size_mult(num_data, num_points)));
110
111 long long stride = static_cast<long long>(num_points + 1) * esize;
112 TasGpu::fillDataGPU(acceleration, correction, num_points, stride,
113 reinterpret_cast<double*>(basis_matrix.data() + Utils::size_mult(num_data, num_points)));
114
115 TasGpu::fillDataGPU(acceleration, 0.0, esize * num_outputs * num_points, 1,
116 reinterpret_cast<double*>(model_values + Utils::size_mult(num_data, num_outputs)));
117 }
118
119 TasmanianDenseSolver::solvesLeastSquaresGPU(acceleration, num_equations, num_points,
120 basis_matrix.data(), num_outputs, model_values);
121}
122
134template<typename scalar_type>
135Data2D<scalar_type> generateCoefficients(double const data_points[], int num_data, double const model_values[], double tolerance, TasmanianSparseGrid &grid){
136 int num_dimensions = grid.getNumDimensions();
137 int num_outputs = grid.getNumOutputs();
138 int num_points = grid.getNumPoints();
139 int num_equations = (tolerance > 0.0) ? num_data + num_points : num_data;
140
141 AccelerationContext const *acceleration = grid.getAccelerationContext();
142
143 Data2D<scalar_type> basis_matrix(num_points, num_equations, 0.0);
144 Data2D<scalar_type> coefficients(num_outputs, num_equations, 0.0);
145
146 if (acceleration->mode == accel_gpu_magma and hasGPUBasis(grid)){
147 int mem_usage = 268435456 / ((std::is_same<double, scalar_type>::value) ? 1 : 2); // 2GB of double numbers
148 if (grid.getGPUMemory(acceleration->device) < 2048) mem_usage /= 2; // use minimum of 1GB
149 int num_batch = std::min(num_data, mem_usage / num_points);
150 if (num_batch < 128){
151 // cannot use the GPU with less than 128 entries, might as well use the CPU
152 grid.evaluateHierarchicalFunctions(data_points, num_data, reinterpret_cast<double*>(basis_matrix.data()));
153 }else{
154 acceleration->setDevice();
155 GpuVector<double> gpu_points(acceleration, num_dimensions, num_batch);
156 GpuVector<scalar_type> gpu_matrix(acceleration, num_points, num_batch);
157 for(int i = 0; i < num_data; i += num_batch){
158 int num_this_batch = std::min(num_batch, num_data - i);
159 gpu_points.load(acceleration, Utils::size_mult(num_this_batch, num_dimensions),
160 data_points + Utils::size_mult(i, num_dimensions));
161 grid.evaluateHierarchicalFunctionsGPU(gpu_points.data(), num_this_batch, reinterpret_cast<double*>(gpu_matrix.data()));
162 gpu_matrix.unload(acceleration, Utils::size_mult(num_this_batch, num_points), basis_matrix.getStrip(i));
163 }
164 }
165 }else{
166 grid.evaluateHierarchicalFunctions(data_points, num_data, reinterpret_cast<double*>(basis_matrix.data()));
167 }
168
169 if (tolerance > 0.0){
170 double correction = std::sqrt(tolerance);
171 for(int i=0; i<grid.getNumPoints(); i++)
172 basis_matrix.getStrip(i + num_data)[i] = correction;
173 }
174
175 auto icoeff = coefficients.begin();
176 for(size_t i=0; i<Utils::size_mult(num_data, grid.getNumOutputs()); i++)
177 *icoeff++ = model_values[i];
178
179 TasmanianDenseSolver::solvesLeastSquares(acceleration, num_equations, grid.getNumPoints(),
180 basis_matrix.data(), grid.getNumOutputs(), coefficients.data());
181
182 return coefficients;
183}
184
194template<typename scalar_type>
195inline void loadUnstructuredDataL2tmpl(double const data_points[], int num_data, double const model_values[],
196 double tolerance, TasmanianSparseGrid &grid){
197
198 if (grid.empty()) throw std::runtime_error("Cannot use loadUnstructuredDataL2() with an empty grid.");
199 if (grid.getNumNeeded() != 0)
200 grid.mergeRefinement();
201
202 AccelerationContext const *acceleration = grid.getAccelerationContext();
203 if (acceleration->mode == accel_none)
204 throw std::runtime_error("The loadUnstructuredDataL2() method cannot be used with acceleration mode accel_none.");
205
206 int num_dimensions = grid.getNumDimensions();
207 int num_outputs = grid.getNumOutputs();
208 int num_points = grid.getNumPoints();
209 int num_equations = (tolerance > 0.0) ? num_data + num_points : num_data;
210
211 Data2D<scalar_type> coefficients =
212 [&]()->Data2D<scalar_type>{
213 if (acceleration->mode == accel_gpu_cuda and hasGPUBasis(grid)){
214 acceleration->setDevice();
215 GpuVector<double> gpu_points(acceleration, num_dimensions, num_data, data_points);
216 GpuVector<scalar_type> gpu_values(acceleration, num_outputs, num_equations);
217 TasGpu::load_n(acceleration, model_values, Utils::size_mult(num_outputs, num_data), gpu_values.data());
218 generateCoefficientsGPU<scalar_type>(gpu_points.data(), num_data, gpu_values.data(), tolerance, grid);
219 return Data2D<scalar_type>(num_outputs, num_equations, gpu_values.unload(acceleration));
220 }else{
221 return generateCoefficients<scalar_type>(data_points, num_data, model_values, tolerance, grid);
222 }
223 }();
224
225 // do the real-complex split (used in Fourier grids)
226 if (std::is_same<scalar_type, std::complex<double>>::value){
227 std::vector<double> real_coeffs(Utils::size_mult(2 * grid.getNumOutputs(), grid.getNumPoints()));
228 auto icoeff = coefficients.begin();
229 for(size_t i=0; i<Utils::size_mult(grid.getNumOutputs(), grid.getNumPoints()); i++)
230 real_coeffs[i] = std::real(*icoeff++);
231 icoeff = coefficients.begin();
232 for(size_t i=Utils::size_mult(grid.getNumOutputs(), grid.getNumPoints()); i<Utils::size_mult(2 * grid.getNumOutputs(), grid.getNumPoints()); i++)
233 real_coeffs[i] = std::imag(*icoeff++);
234 grid.setHierarchicalCoefficients(real_coeffs.data());
235 }else{
236 grid.setHierarchicalCoefficients(reinterpret_cast<double*>(coefficients.data()));
237 }
238}
239
282inline void loadUnstructuredDataL2(double const data_points[], int num_data, double const model_values[],
283 double tolerance, TasmanianSparseGrid &grid){
284 if (grid.isFourier()){
285 loadUnstructuredDataL2tmpl<std::complex<double>>(data_points, num_data, model_values, tolerance, grid);
286 }else{
287 loadUnstructuredDataL2tmpl<double>(data_points, num_data, model_values, tolerance, grid);
288 }
289}
290
299inline void loadUnstructuredDataL2(std::vector<double> const &data_points, std::vector<double> const &model_values,
300 double tolerance, TasmanianSparseGrid &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());
303 if (model_values.size() < Utils::size_mult(num_data, grid.getNumOutputs()))
304 throw std::runtime_error("In loadUnstructuredDataL2(), provided more points than data.");
305 loadUnstructuredDataL2(data_points.data(), num_data, model_values.data(), tolerance, grid);
306}
307
308}
309
310#endif
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 * getStrip(int i)
Returns a reference to the i-th strip.
Definition tsgIndexSets.hpp:141
T * data()
Returns a reference to the internal data.
Definition tsgIndexSets.hpp:153
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
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).
AccelerationContext const * getAccelerationContext() const
Allows the addon methods to use the acceleration context.
Definition TasmanianSparseGrid.hpp:2103
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 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
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
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
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:727
Templates for non-adaptive sampling from lambda models.