Doxygen 1.9.1
Toolkit for Adaptive Stochastic Modeling and Non-Intrusive ApproximatioN: Tasmanian v8.2 (development)
tsgGridCore.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 __TSG_BASE_CLASS_HPP
32 #define __TSG_BASE_CLASS_HPP
33 
46 #include "tsgCudaLoadStructures.hpp"
48 
49 #ifndef __TASMANIAN_DOXYGEN_SKIP
50 namespace TasGrid{
51 
52 class BaseCanonicalGrid{
53 public:
54  BaseCanonicalGrid(AccelerationContext const *acc) : acceleration(acc), num_dimensions(0), num_outputs(0){}
55  BaseCanonicalGrid(AccelerationContext const *acc, int cnum_dimensions, int cnum_outputs, MultiIndexSet &&cpoints, MultiIndexSet &&cneeded, StorageSet &&cvalues)
56  : acceleration(acc), num_dimensions(cnum_dimensions), num_outputs(cnum_outputs),
57  points(std::forward<MultiIndexSet>(cpoints)), needed(std::forward<MultiIndexSet>(cneeded)),
58  values(std::forward<StorageSet>(cvalues)){}
59  BaseCanonicalGrid(AccelerationContext const *acc, BaseCanonicalGrid const &other, int ibegin, int iend)
60  : acceleration(acc), num_dimensions(other.num_dimensions), num_outputs(iend - ibegin),
61  points(other.points), needed(other.needed),
62  values((iend - ibegin == other.num_outputs) ? other.values : other.values.splitValues(ibegin, iend))
63  {}
64  virtual ~BaseCanonicalGrid() = default;
65 
66  virtual bool isGlobal() const{ return false; }
67  virtual bool isSequence() const{ return false; }
68  virtual bool isLocalPolynomial() const{ return false; }
69  virtual bool isWavelet() const{ return false; }
70  virtual bool isFourier() const{ return false; }
71 
72  int getNumDimensions() const{ return num_dimensions; }
73  int getNumOutputs() const{ return num_outputs; }
74  virtual TypeOneDRule getRule() const = 0;
75 
76  int getNumLoaded() const{ return (num_outputs == 0) ? 0 : points.getNumIndexes(); }
77  int getNumNeeded() const{ return needed.getNumIndexes(); }
78  int getNumPoints() const{ return ((points.empty()) ? needed.getNumIndexes() : points.getNumIndexes()); }
79  const int* getPointIndexes() const{ return ((points.empty()) ? needed.getIndex(0) : points.getIndex(0)); }
80 
81  virtual void write(std::ostream&, bool) const = 0;
82 
83  virtual void getLoadedPoints(double *x) const = 0;
84  virtual void getNeededPoints(double *x) const = 0;
85  virtual void getPoints(double *x) const = 0;
86 
87  virtual void getQuadratureWeights(double weights[]) const = 0;
88  virtual void getInterpolationWeights(const double x[], double weights[]) const = 0;
89  virtual void getDifferentiationWeights(const double x[], double weights[]) const = 0;
90 
91  virtual void loadNeededValues(const double *vals) = 0;
92  const double* getLoadedValues() const{ return (points.empty()) ? nullptr : values.getValues(0); }
93 
94  virtual void evaluate(const double x[], double y[]) const = 0;
95  virtual void integrate(double q[], double *conformal_correction) const = 0;
96  virtual void differentiate(const double x[], double jacobian[]) const = 0;
97 
98  virtual void evaluateBatch(const double x[], int num_x, double y[]) const = 0;
99 
100  virtual void evaluateBatchGPU(const double[], int, double[]) const = 0;
101  virtual void evaluateBatchGPU(const float[], int, float[]) const = 0;
102  virtual void evaluateHierarchicalFunctionsGPU(const double[], int, double[]) const = 0;
103  virtual void evaluateHierarchicalFunctionsGPU(const float[], int, float[]) const = 0;
104 
105  virtual void clearRefinement() = 0;
106  virtual void mergeRefinement() = 0;
107 
108  virtual void beginConstruction() = 0;
109  virtual void writeConstructionData(std::ostream&, bool) const = 0;
110  virtual void readConstructionData(std::istream&, bool) = 0;
111  virtual void loadConstructedPoint(const double[], const std::vector<double> &) = 0;
112  virtual void loadConstructedPoint(const double[], int, const double[]) = 0;
113  virtual void finishConstruction() = 0;
114 
115  virtual void evaluateHierarchicalFunctions(const double x[], int num_x, double y[]) const = 0; // add acceleration here
116  virtual std::vector<double> getSupport() const{ return std::vector<double>(Utils::size_mult(getNumPoints(), getNumDimensions()), 2.0); }
117  virtual void setHierarchicalCoefficients(const double c[]) = 0;
118  virtual void integrateHierarchicalFunctions(double integrals[]) const = 0;
119 
120  virtual void updateAccelerationData(AccelerationContext::ChangeType change) const = 0;
121 
122 protected:
123  AccelerationContext const *acceleration;
124  int num_dimensions, num_outputs;
125  MultiIndexSet points;
126  MultiIndexSet needed;
127  StorageSet values;
128 };
129 
130 // For purposes of reading grids with older file formats
131 // Specialize the reader class for each grid type
132 // Befriend the reader class and each grid type
133 template<class GridType>
134 struct GridReaderVersion5{ // 5 refers to the file format version, not the Tasmanian version
135  template<typename iomode> static std::unique_ptr<BaseCanonicalGrid> read(AccelerationContext const*, std::istream &){
136  return std::unique_ptr<BaseCanonicalGrid>();
137  }
138 };
139 
140 // Factory reader method that instantiates the GridReaderVersion5 class
141 template<class GridType, typename iomode> std::unique_ptr<GridType> readGridVersion5(AccelerationContext const *acc, std::istream &is, iomode){
142  return GridReaderVersion5<GridType>::template read<iomode>(acc, is);
143 }
144 
145 }
146 #endif
147 
148 #endif
TypeOneDRule
Used to specify the one dimensional family of rules that induces the sparse grid.
Definition: tsgEnumerates.hpp:285
void loadNeededValues(std::function< void(double const x[], double y[], size_t thread_id)> model, TasmanianSparseGrid &grid, size_t num_threads)
Loads the current grid with model values, does not perform any refinement.
Definition: tsgLoadNeededValues.hpp:104
size_t size_mult(IntA a, IntB b)
Converts two integer-like variables to size_t and returns the product..
Definition: tsgUtils.hpp:82
bool isGlobal(TypeOneDRule rule)
Return True if the rule has basis with Lagrange polynomials with global support, fit for GridGlobal,...
bool isFourier(TypeOneDRule rule)
Return True if the rule has trigonometric basis.
bool isLocalPolynomial(TypeOneDRule rule)
Return True if the rule has polynomial basis with local support, fit for GridLocalPolynomial,...
bool isWavelet(TypeOneDRule rule)
Return True if the rule has wavelet basis with local support.
bool isSequence(TypeOneDRule rule)
Return True if the rule is nested and has single-node growth, fit for GridSequence,...
int getNumPoints(int level, TypeOneDRule rule)
Return the number of points for the rule at the level, includes all global rules.
Encapsulates the Tasmanian Sparse Grid module.
Definition: TasmanianSparseGrid.hpp:68
ChangeType
Defines the types of acceleration context updates so they can be linked to acceleration cache updates...
Definition: tsgAcceleratedDataStructures.hpp:593
Class and data types used for dynamic construction of Global grids.
Algorithms for manipulating multi-indexes defined by hierarchy rules.