Doxygen 1.9.8
Toolkit for Adaptive Stochastic Modeling and Non-Intrusive ApproximatioN: Tasmanian v8.2
 
Loading...
Searching...
No Matches
tsgOneDimensionalWrapper.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_ONE_DIMENSIONAL_WRAPPER_HPP
32#define __TSG_ONE_DIMENSIONAL_WRAPPER_HPP
33
35
47namespace TasGrid{
48
54public:
56 OneDimensionalWrapper() : isNonNested(false), num_levels(0), rule(rule_none){}
59
68 OneDimensionalWrapper(const CustomTabulated &custom, int max_level, TypeOneDRule crule, double alpha, double beta) :
69 isNonNested(OneDimensionalMeta::isNonNested(crule)),
70 num_levels(max_level + 1),
71 rule(crule),
72 num_points(num_levels),
73 pntr(num_levels + 1, 0), // the first entry must be set to zero
74 weights(num_levels),
75 nodes(num_levels),
76 coeff(num_levels)
77 {
78 if (crule == rule_customtabulated){
79 if (max_level + 1 > custom.getNumLevels()){
80 std::string message = "ERROR: custom-tabulated rule needed with levels ";
81 message += std::to_string(num_levels);
82 message += ", but only ";
83 message += std::to_string(custom.getNumLevels());
84 message += " are provided.";
85 throw std::runtime_error(message);
86 }
87 }
88 // find the points per level and the cumulative pointers
89 if (rule != rule_customtabulated){
90 for(int l=0; l<num_levels; l++){
91 num_points[l] = OneDimensionalMeta::getNumPoints(l, rule);
92 pntr[l+1] = pntr[l] + num_points[l];
93 }
94 }else{
95 for(int l=0; l<num_levels; l++){
96 num_points[l] = custom.getNumPoints(l);
97 pntr[l+1] = pntr[l] + num_points[l];
98 }
99 }
100 int num_total = pntr[num_levels];
101
102 if (isNonNested){
103 indx.reserve(num_total);
104 unique.reserve(num_total);
105
106 for(int l=0; l<num_levels; l++){
107 int n = num_points[l];
108 if ((rule == rule_chebyshev) || (rule == rule_chebyshevodd)){
109 OneDimensionalNodes::getChebyshev(n, weights[l], nodes[l]);
110 }else if ((rule == rule_gausslegendre) || (rule == rule_gausslegendreodd)){
111 OneDimensionalNodes::getGaussLegendre(n, weights[l], nodes[l]);
112 }else if ((rule == rule_gausschebyshev1) || (rule == rule_gausschebyshev1odd)){
113 OneDimensionalNodes::getGaussChebyshev1(n, weights[l], nodes[l]);
114 }else if ((rule == rule_gausschebyshev2) || (rule == rule_gausschebyshev2odd)){
115 OneDimensionalNodes::getGaussChebyshev2(n, weights[l], nodes[l]);
116 }else if ((rule == rule_gaussgegenbauer) || (rule == rule_gaussgegenbauerodd)){
117 OneDimensionalNodes::getGaussJacobi(n, weights[l], nodes[l], alpha, alpha);
118 }else if ((rule == rule_gausshermite) || (rule == rule_gausshermiteodd)){
119 OneDimensionalNodes::getGaussHermite(n, weights[l], nodes[l], alpha);
120 }else if ((rule == rule_gaussjacobi) || (rule == rule_gaussjacobiodd)){
121 OneDimensionalNodes::getGaussJacobi(n, weights[l], nodes[l], alpha, beta);
122 }else if ((rule == rule_gausslaguerre) || (rule == rule_gausslaguerreodd)){
123 OneDimensionalNodes::getGaussLaguerre(n, weights[l], nodes[l], alpha);
124 }else if (rule == rule_customtabulated){
125 custom.getWeightsNodes(l, weights[l], nodes[l]);
126 }
127
128 for(auto x : nodes[l]){
129 int point = -1;
130 for(int j=0; j<(int) unique.size(); j++){
131 if (std::abs(x - unique[j]) < Maths::num_tol){
132 point = j;
133 break;
134 }
135 }
136 if (point == - 1){ // new point found
137 unique.push_back(x);
138 point = (int) (unique.size() - 1);
139 }
140 indx.push_back(point);
141 }
142 }
143
144 unique.shrink_to_fit();
145
146 // make coefficients
147 for(int l=0; l<num_levels; l++){
148 int n = num_points[l];
149 coeff[l].resize(n);
150 double *x = nodes[l].data();
151 for(int i=0; i<n; i++){
152 double c = 1.0;
153 for(int j=0; j<i; j++){
154 c /= (x[i] - x[j]);
155 }
156 for(int j=i+1; j<n; j++){
157 c /= (x[i] - x[j]);
158 }
159 coeff[l][i] = c;
160 }
161 }
162 }else{
163 if (rule == rule_clenshawcurtis){
164 unique = OneDimensionalNodes::getClenshawCurtisNodes(max_level);
165 }else if (rule == rule_clenshawcurtis0){
166 unique = OneDimensionalNodes::getClenshawCurtisNodesZero(max_level);
167 }else if (rule == rule_fejer2){
168 unique = OneDimensionalNodes::getFejer2Nodes(max_level);
169 }else if (rule == rule_gausspatterson){
171 if (num_levels > gp.getNumLevels()){
172 std::string message = "ERROR: gauss-patterson rule needed with level ";
173 message += std::to_string(max_level);
174 message += ", but only ";
175 message += std::to_string(gp.getNumLevels());
176 message += " are hardcoded.";
177 throw std::runtime_error(message);
178 }
179 unique = gp.getNodes(max_level);
180
181 // move here to avoid a warning
182 for(int l=0; l<num_levels; l++){
183 weights[l].resize(num_points[l]);
184 for(int i=0; i<num_points[l]; i++){
185 weights[l][i] = gp.getWeight(l, i);
186 }
187 }
188 }else if (rule == rule_rleja){
189 unique = OneDimensionalNodes::getRLeja(OneDimensionalMeta::getNumPoints(max_level,rule));
190 }else if ((rule == rule_rlejaodd) || (rule == rule_rlejadouble2) || (rule == rule_rlejadouble4)){
191 unique = OneDimensionalNodes::getRLejaCentered(OneDimensionalMeta::getNumPoints(max_level,rule));
192 }else if ((rule == rule_rlejashifted) || (rule == rule_rlejashiftedeven) || (rule == rule_rlejashifteddouble)){
193 unique = OneDimensionalNodes::getRLejaShifted(OneDimensionalMeta::getNumPoints(max_level,rule));
194 }else if ((rule == rule_leja) || (rule == rule_lejaodd)){
195 unique = Optimizer::getGreedyNodes<rule_leja>(OneDimensionalMeta::getNumPoints(max_level, rule));
196 }else if ((rule == rule_maxlebesgue) || (rule == rule_maxlebesgueodd)){
197 unique = Optimizer::getGreedyNodes<rule_maxlebesgue>(OneDimensionalMeta::getNumPoints(max_level, rule));
198 }else if ((rule == rule_minlebesgue) || (rule == rule_minlebesgueodd)){
199 unique = Optimizer::getGreedyNodes<rule_minlebesgue>(OneDimensionalMeta::getNumPoints(max_level, rule));
200 }else if ((rule == rule_mindelta) || (rule == rule_mindeltaodd)){
201 unique = Optimizer::getGreedyNodes<rule_mindelta>(OneDimensionalMeta::getNumPoints(max_level, rule));
202 }else{ // if (rule==rule_fourier)
203 unique = OneDimensionalNodes::getFourierNodes(max_level);
204 }
205
206 for(int l=0; l<num_levels; l++){
207 int n = num_points[l];
208 weights[l].resize(n);
209 coeff[l].resize(n);
210 for(int i=0; i<n; i++){
211 double c = 1.0;
212 for(int j=0; j<i; j++){
213 c /= (unique[i] - unique[j]);
214 }
215 for(int j=i+1; j<n; j++){
216 c /= (unique[i] - unique[j]);
217 }
218 if (rule == rule_clenshawcurtis0){
219 c /= (unique[i] - 1.0) * (unique[i] + 1.0);
220 }
221 coeff[l][i] = c;
222 }
223 }
224
225 if (rule == rule_clenshawcurtis){
226 for(int l=0; l<num_levels; l++){
227 for(int i=0; i<num_points[l]; i++){
228 weights[l][i] = OneDimensionalNodes::getClenshawCurtisWeight(l, i);
229 }
230 }
231 }else if (rule == rule_clenshawcurtis0){
232 for(int l=0; l<num_levels; l++){
233 for(int i=0; i<num_points[l]; i++){
234 weights[l][i] = OneDimensionalNodes::getClenshawCurtisWeightZero(l, i);
235 }
236 }
237 }else if (rule == rule_fejer2){
238 for(int l=0; l<num_levels; l++){
239 for(int i=0; i<num_points[l]; i++){
240 weights[l][i] = OneDimensionalNodes::getFejer2Weight(l, i);
241 }
242 }
243 }else if ((rule == rule_rleja) || (rule == rule_rlejaodd) || (rule == rule_rlejadouble2) || (rule == rule_rlejadouble4) || (rule == rule_leja) || (rule == rule_lejaodd) ||
244 (rule == rule_rlejashifted) || (rule == rule_rlejashiftedeven) || (rule == rule_rlejashifteddouble) ||
245 (rule == rule_maxlebesgue) || (rule == rule_maxlebesgueodd) || (rule == rule_minlebesgue) || (rule == rule_minlebesgueodd) || (rule == rule_mindelta) || (rule == rule_mindeltaodd)){
246 // if missing the quadrature weights
247 int n = 1 + num_points[max_level] / 2; // number of Gauss-Legendre points needed to integrate the basis functions
248 std::vector<double> lag_x, lag_w;
249 OneDimensionalNodes::getGaussLegendre(n, lag_w, lag_x);
250 for(int l=0; l<num_levels; l++){
251 std::fill(weights[l].begin(), weights[l].end(), 0.0);
252 int npl = num_points[l];
253 std::vector<double> v(npl);
254 for(int i=0; i<n; i++){
255 v[0] = 1.0;
256 for(int j=0; j<npl-1; j++){
257 v[j+1] = (lag_x[i] - unique[j]) * v[j];
258 }
259 v[npl-1] *= coeff[l][npl-1];
260 double s = 1.0;
261 for(int j=npl-2; j>=0; j--){
262 s *= (lag_x[i] - unique[j+1]);
263 v[j] *= s * coeff[l][j];
264 }
265 for(int j=0; j<npl; j++){
266 weights[l][j] += lag_w[i] * v[j];
267 }
268 }
269 }
270 }
271 }
272 }
273
275 OneDimensionalWrapper(int max_level, TypeOneDRule crule, double alpha, double beta) :
276 OneDimensionalWrapper(CustomTabulated(), max_level, crule, alpha, beta){}
277
279 int getNumPoints(int level) const{ return num_points[level]; }
281 int getPointIndex(int level, int j) const{ return indx[pntr[level] + j]; }
283 int getLevel(int point) {
284 // assumes we have the correct number of levels loaded
285 int l = 0;
286 while(point >= num_points[l]) l++;
287 return l;
288 }
289 std::vector<int> getLevels(std::vector<int> const &point){
290 std::vector<int> result(point.size());
291 for(size_t i=0; i<point.size(); i++) result[i] = getLevel(point[i]);
292 return result;
293 }
294
296 int getNumNodes() const{ return (int) unique.size(); }
298 double getNode(int j) const{ return unique[j]; }
300 double getWeight(int level, int j) const{ return weights[level][j]; }
301
303 const double* getNodes(int level) const{ return (isNonNested) ? nodes[level].data() : unique.data(); }
305 const double* getCoefficients(int level) const{ return coeff[level].data(); }
306
308 const std::vector<int>& getPointsCount() const{ return pntr; }
309
311 TypeOneDRule getRule() const{ return rule; }
313 int getNumLevels() const{ return num_levels; }
314
316 const std::vector<double>& getUnique() const{ return unique; }
318 std::vector<double> getAllNodes() const{ return Utils::mergeVectors(nodes); }
320 std::vector<double> getAllCoeff() const{ return Utils::mergeVectors(coeff); }
322 const std::vector<int>& getNumNodesPerLevel() const{ return num_points; }
324 std::vector<int> getOffsetNodesPerLevel() const{
325 std::vector<int> offsets(num_points.size());
326 offsets[0] = 0;
327 for(size_t i=1; i<num_points.size(); i++)
328 offsets[i] = offsets[i-1] + num_points[i-1];
329 return offsets;
330 }
331
332private:
333 bool isNonNested;
334 int num_levels;
335 TypeOneDRule rule;
336
337 std::vector<int> num_points; // give the number of points per level
338 std::vector<int> pntr; // gives the cumulative offset of each level
339 std::vector<int> indx; // gives a list of the points associated with each level
340
341 std::vector<std::vector<double>> weights; // contains the weight associated with each level
342 std::vector<std::vector<double>> nodes; // contains all nodes for each level
343 std::vector<double> unique; // contains the x-coordinate of each sample point
344
345 std::vector<std::vector<double>> coeff; // the coefficients of the Lagrange
346};
347
348}
349
350#endif
Class providing manipulation of custom tabulated rules, file I/O and structured access to the points,...
Definition tsgCoreOneDimensional.hpp:66
int getNumPoints(int level) const
Returns the number of points associated with the selected level.
Definition tsgCoreOneDimensional.hpp:94
int getNumLevels() const
Returns the number of loaded levels.
Definition tsgCoreOneDimensional.hpp:92
void getWeightsNodes(int level, std::vector< double > &w, std::vector< double > &x) const
Get the points x and quadrature weights w associated with the rule at the level.
A class to cache one dimensional rules, nodes, weight, meta-data, etc.
Definition tsgOneDimensionalWrapper.hpp:53
int getLevel(int point)
Returns the first level associated with this point (nested rules only!).
Definition tsgOneDimensionalWrapper.hpp:283
const std::vector< int > & getPointsCount() const
Get a reference to the offsets of all point counts used by the CacheLagrange class.
Definition tsgOneDimensionalWrapper.hpp:308
const double * getCoefficients(int level) const
Get an array of the Lagrange coefficients associated with the level, the order matches getNodes().
Definition tsgOneDimensionalWrapper.hpp:305
const std::vector< int > & getNumNodesPerLevel() const
Return reference to the number of nodes per level.
Definition tsgOneDimensionalWrapper.hpp:322
OneDimensionalWrapper(int max_level, TypeOneDRule crule, double alpha, double beta)
Overload that skips the custom rule altogether, more convenient in Fourier grids.
Definition tsgOneDimensionalWrapper.hpp:275
const std::vector< double > & getUnique() const
Return reference to all unique nodes.
Definition tsgOneDimensionalWrapper.hpp:316
double getWeight(int level, int j) const
Get the quadrature weight of the j-th node on the level (for non-nested rules, using index local to t...
Definition tsgOneDimensionalWrapper.hpp:300
std::vector< int > getOffsetNodesPerLevel() const
Return cumulative points per level.
Definition tsgOneDimensionalWrapper.hpp:324
OneDimensionalWrapper(const CustomTabulated &custom, int max_level, TypeOneDRule crule, double alpha, double beta)
Load the cache.
Definition tsgOneDimensionalWrapper.hpp:68
int getNumLevels() const
Get the number of cached levels.
Definition tsgOneDimensionalWrapper.hpp:313
int getNumNodes() const
Get number of loaded nodes.
Definition tsgOneDimensionalWrapper.hpp:296
TypeOneDRule getRule() const
Get the loaded rule.
Definition tsgOneDimensionalWrapper.hpp:311
std::vector< double > getAllNodes() const
Return all nodes concatenated per level.
Definition tsgOneDimensionalWrapper.hpp:318
std::vector< double > getAllCoeff() const
Return all coefficients concatenated per level.
Definition tsgOneDimensionalWrapper.hpp:320
int getPointIndex(int level, int j) const
Get the global index of point j on level (used only by non-nested rules).
Definition tsgOneDimensionalWrapper.hpp:281
const double * getNodes(int level) const
Get an array of all nodes associated with the level (the array is ordered by local index).
Definition tsgOneDimensionalWrapper.hpp:303
double getNode(int j) const
Get the canonical coordinate of the node with global index j.
Definition tsgOneDimensionalWrapper.hpp:298
~OneDimensionalWrapper()=default
Default destructor.
int getNumPoints(int level) const
Get the number of points for the level, can only query loaded levels.
Definition tsgOneDimensionalWrapper.hpp:279
OneDimensionalWrapper()
Default constructor, create an emptry cache.
Definition tsgOneDimensionalWrapper.hpp:56
Rule with hard-corded tabulated points and weights.
Definition tsgHardCodedTabulatedRules.hpp:64
double getWeight(int level, int point) const
Return the quadrature weight for level and given point.
std::vector< double > getNodes(int level) const
Returns the nodes for the level, note that the nodes are nested.
static int getNumLevels()
Return the number of hard-coded levels.
Definition tsgHardCodedTabulatedRules.hpp:72
TypeOneDRule
Used to specify the one dimensional family of rules that induces the sparse grid.
Definition tsgEnumerates.hpp:285
@ rule_minlebesgue
A greedy sequence rule with nodes added to minimize the Lebesgue constant.
Definition tsgEnumerates.hpp:321
@ rule_gaussgegenbauer
Non-nested rule optimized for integral of the form .
Definition tsgEnumerates.hpp:343
@ rule_gausshermite
Non-nested rule optimized for integral of the form .
Definition tsgEnumerates.hpp:355
@ rule_gausslegendre
Non-nested rule but optimized for integration.
Definition tsgEnumerates.hpp:329
@ rule_fejer2
Similar to rule_clenshawcurtis but with nodes strictly in the interior.
Definition tsgEnumerates.hpp:293
@ rule_leja
Classic sequence rule, moderate Lebesgue constant growth (empirical result only).
Definition tsgEnumerates.hpp:299
@ rule_rlejashifteddouble
Same as rule_rlejashifted but doubling the number of nodes per level, which reduced the Lebesgue cons...
Definition tsgEnumerates.hpp:315
@ rule_gausslaguerreodd
Same as rule_gausslaguerre but using only odd levels, partially mitigates the non-nested issues.
Definition tsgEnumerates.hpp:353
@ rule_chebyshevodd
Same as rule_chebyshev but using only odd levels, partially mitigates the non-nested issues.
Definition tsgEnumerates.hpp:297
@ rule_minlebesgueodd
Same as rule_minlebesgue but using only odd levels, quadrature is more stable.
Definition tsgEnumerates.hpp:323
@ rule_gausslegendreodd
Same as rule_gausslegendre but using only odd levels, partially mitigates the non-nested issues.
Definition tsgEnumerates.hpp:331
@ rule_gausschebyshev1
Non-nested rule optimized for integral of the form .
Definition tsgEnumerates.hpp:335
@ rule_gausslaguerre
Non-nested rule optimized for integral of the form .
Definition tsgEnumerates.hpp:351
@ rule_mindelta
A greedy sequence rule with nodes added to minimize the norm of the surplus operator.
Definition tsgEnumerates.hpp:325
@ rule_gaussgegenbauerodd
Same as rule_gaussgegenbauer but using only odd levels, partially mitigates the non-nested issues.
Definition tsgEnumerates.hpp:345
@ rule_customtabulated
User provided rule, nodes and weights must be provided with a separate file.
Definition tsgEnumerates.hpp:359
@ rule_gausshermiteodd
Same as rule_gausshermite but using only odd levels, partially mitigates the non-nested issues.
Definition tsgEnumerates.hpp:357
@ rule_clenshawcurtis
Classic nested rule using Chebyshev nodes with very low Lebesgue constant.
Definition tsgEnumerates.hpp:289
@ rule_gausschebyshev2odd
Same as rule_gausschebyshev2 but using only odd levels, partially mitigates the non-nested issues.
Definition tsgEnumerates.hpp:341
@ rule_maxlebesgueodd
Same as rule_maxlebesgue but using only odd levels, quadrature is more stable.
Definition tsgEnumerates.hpp:319
@ rule_maxlebesgue
A greedy sequence rule with nodes placed at the maximum of the Lebesgue function.
Definition tsgEnumerates.hpp:317
@ rule_clenshawcurtis0
Same as rule_clenshawcurtis but with modified basis that assumes the model is zero at the boundary.
Definition tsgEnumerates.hpp:291
@ rule_none
Null rule, should never be used as input (default rule for an empty grid).
Definition tsgEnumerates.hpp:287
@ rule_chebyshev
Using Chebyshev nodes with very low Lebesgue constant and slow node growth, but non-nested.
Definition tsgEnumerates.hpp:295
@ rule_gausschebyshev1odd
Same as rule_gausschebyshev1 but using only odd levels, partially mitigates the non-nested issues.
Definition tsgEnumerates.hpp:337
@ rule_rlejashifted
Similar sequence to rule_rleja but with nodes strictly in the interior.
Definition tsgEnumerates.hpp:311
@ rule_gausschebyshev2
Non-nested rule optimized for integral of the form .
Definition tsgEnumerates.hpp:339
@ rule_mindeltaodd
Same as rule_mindelta but using only odd levels, quadrature is more stable.
Definition tsgEnumerates.hpp:327
@ rule_rlejadouble4
Using rule_rleja nodes but doubling the nodes every 4 levels, reduces the Lebesgue constant.
Definition tsgEnumerates.hpp:307
@ rule_gaussjacobiodd
Same as rule_gaussjacobi but using only odd levels, partially mitigates the non-nested issues.
Definition tsgEnumerates.hpp:349
@ rule_lejaodd
Same as rule_leja but using only odd levels, quadrature is more stable.
Definition tsgEnumerates.hpp:301
@ rule_rlejaodd
Same as rule_rleja but using only odd levels, quadrature is more stable.
Definition tsgEnumerates.hpp:309
@ rule_rlejadouble2
Using rule_rleja nodes but doubling the nodes every 2 levels, reduces the Lebesgue constant.
Definition tsgEnumerates.hpp:305
@ rule_rlejashiftedeven
Same as rule_rlejashifted but using only even levels, quadrature is more stable.
Definition tsgEnumerates.hpp:313
@ rule_gausspatterson
Nested rule that is optimized for integration, probably the best integration rule in more than 2 dime...
Definition tsgEnumerates.hpp:333
@ rule_rleja
Classic sequence rule based on complex analysis, moderate Lebesgue constant growth (theoretically pro...
Definition tsgEnumerates.hpp:303
@ rule_gaussjacobi
Non-nested rule optimized for integral of the form .
Definition tsgEnumerates.hpp:347
constexpr double num_tol
Numerical tolerance for various algorithms.
Definition tsgMathUtils.hpp:132
std::vector< T > mergeVectors(std::vector< std::vector< T > > const &vec)
Takes a vector of vectors and returns a single contiguous vector.
Definition tsgUtils.hpp:105
Encapsulates the Tasmanian Sparse Grid module.
Definition TasmanianSparseGrid.hpp:68
Hard-coded nodes and weights.