Doxygen 1.9.8
Toolkit for Adaptive Stochastic Modeling and Non-Intrusive ApproximatioN: Tasmanian v8.2
 
Loading...
Searching...
No Matches
tsgIOHelpers.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_IOHELPERS_HPP
32#define __TASMANIAN_IOHELPERS_HPP
33
34#include "tsgEnumerates.hpp"
35#include <functional> // put here because it causes a conflict with some versinos of CUDA (keep out of tsgAcceleratedStructures.hpp)
36
56namespace TasGrid{
57
62constexpr bool mode_ascii = false;
63
68constexpr bool mode_binary = true;
69
70
76namespace IO{
77
83
89
111
116inline std::map<std::string, TypeOneDRule> getStringRuleMap(){
117 return std::initializer_list<std::pair<std::string const, TypeOneDRule>>{
118 {"none", rule_none},
119 {"clenshaw-curtis", rule_clenshawcurtis},
120 {"clenshaw-curtis-zero", rule_clenshawcurtis0},
121 {"chebyshev", rule_chebyshev},
122 {"chebyshev-odd", rule_chebyshevodd},
123 {"gauss-legendre", rule_gausslegendre},
124 {"gauss-legendre-odd", rule_gausslegendreodd},
125 {"gauss-patterson", rule_gausspatterson},
126 {"leja", rule_leja},
127 {"leja-odd", rule_lejaodd},
128 {"rleja", rule_rleja},
129 {"rleja-double2", rule_rlejadouble2},
130 {"rleja-double4", rule_rlejadouble4},
131 {"rleja-odd", rule_rlejaodd},
132 {"rleja-shifted", rule_rlejashifted},
133 {"rleja-shifted-even", rule_rlejashiftedeven},
134 {"rleja-shifted-double", rule_rlejashifteddouble},
135 {"max-lebesgue", rule_maxlebesgue},
136 {"max-lebesgue-odd", rule_maxlebesgueodd},
137 {"min-lebesgue", rule_minlebesgue},
138 {"min-lebesgue-odd", rule_minlebesgueodd},
139 {"min-delta", rule_mindelta},
140 {"min-delta-odd", rule_mindeltaodd},
141 {"gauss-chebyshev1", rule_gausschebyshev1},
142 {"gauss-chebyshev1-odd", rule_gausschebyshev1odd},
143 {"gauss-chebyshev2", rule_gausschebyshev2},
144 {"gauss-chebyshev2-odd", rule_gausschebyshev2odd},
145 {"fejer2", rule_fejer2},
146 {"gauss-gegenbauer", rule_gaussgegenbauer},
147 {"gauss-gegenbauer-odd", rule_gaussgegenbauerodd},
148 {"gauss-jacobi", rule_gaussjacobi},
149 {"gauss-jacobi-odd", rule_gaussjacobiodd},
150 {"gauss-laguerre", rule_gausslaguerre},
151 {"gauss-laguerre-odd", rule_gausslaguerreodd},
152 {"gauss-hermite", rule_gausshermite},
153 {"gauss-hermite-odd", rule_gausshermiteodd},
154 {"custom-tabulated", rule_customtabulated},
155 {"localp", rule_localp},
156 {"localp-zero", rule_localp0},
157 {"localp-boundary", rule_localpb},
158 {"semi-localp", rule_semilocalp},
159 {"wavelet", rule_wavelet},
160 {"fourier", rule_fourier}};
161}
162
167inline TypeOneDRule getRuleString(std::string const &name){
168 try{
169 return getStringRuleMap().at(name);
170 }catch(std::out_of_range &){
171 return rule_none;
172 }
173}
174
179inline std::string getRuleString(TypeOneDRule rule){
180 auto smap = getStringRuleMap();
181 return std::find_if(smap.begin(), smap.end(),
182 [&](std::pair<std::string, TypeOneDRule> r)->bool{ return (r.second == rule); })->first;
183}
184
204
210 auto rmap = getIntRuleMap();
211 return ((size_t) r < rmap.size()) ? rmap[(size_t) r] : rule_none;
212}
213
218inline int getRuleInt(TypeOneDRule rule){
219 auto rmap = getIntRuleMap();
220 return (int) std::distance(rmap.begin(), std::find_if(rmap.begin(), rmap.end(),
221 [&](TypeOneDRule r)->bool{ return (r == rule); }));
222}
223
228inline std::map<std::string, TypeDepth> getStringToDepthMap(){
229 return std::initializer_list<std::pair<std::string const, TypeDepth>>{
230 {"level", type_level},
231 {"curved", type_curved},
232 {"iptotal", type_iptotal},
233 {"ipcurved", type_ipcurved},
234 {"qptotal", type_qptotal},
235 {"qpcurved", type_qpcurved},
236 {"hyperbolic", type_hyperbolic},
237 {"iphyperbolic", type_iphyperbolic},
238 {"qphyperbolic", type_qphyperbolic},
239 {"tensor", type_tensor},
240 {"iptensor", type_iptensor},
241 {"qptensor", type_qptensor}};
242}
247inline TypeDepth getDepthTypeString(std::string const &name){
248 try{
249 return getStringToDepthMap().at(name);
250 }catch(std::out_of_range &){
251 return type_none;
252 }
253}
254
265
270inline std::map<std::string, TypeRefinement> getStringToRefinementMap(){
271 return std::initializer_list<std::pair<std::string const, TypeRefinement>>{
272 {"classic", refine_classic},
273 {"parents", refine_parents_first},
274 {"direction", refine_direction_selective},
275 {"fds", refine_fds},
276 {"stable", refine_stable}};
277}
282inline TypeRefinement getTypeRefinementString(std::string const &name){
283 try{
284 return getStringToRefinementMap().at(name);
285 }catch(std::out_of_range &){
286 return refine_none;
287 }
288}
299
304template<bool iomode, IOPad pad>
305void writeFlag(bool flag, std::ostream &os){
306 if (iomode == mode_ascii){
307 os << ((flag) ? "1" : "0");
308 if ((pad == pad_rspace) || ((pad == pad_auto) && flag)) os << " ";
309 if ((pad == pad_line) || ((pad == pad_auto) && !flag)) os << std::endl;
310 }else{
311 char cflag = ((flag) ? 'y' : 'n');
312 os.write(&cflag, sizeof(char));
313 }
314}
315
320template<typename iomode>
321bool readFlag(std::istream &os){
322 if (std::is_same<iomode, mode_ascii_type>::value){
323 int flag;
324 os >> flag;
325 return (flag != 0);
326 }else{
327 char cflag;
328 os.read(&cflag, sizeof(char));
329 return (cflag == 'y');
330 }
331}
332
337template<bool iomode, IOPad pad, typename VecType>
338void writeVector(const std::vector<VecType> &x, std::ostream &os){
339 if (iomode == mode_ascii){
340 if (pad == pad_lspace)
341 for(auto i : x) os << " " << i;
342 if (pad == pad_rspace)
343 for(auto i : x) os << i << " ";
344 if ((pad == pad_none) || (pad == pad_line)){
345 os << x[0];
346 for(size_t i = 1; i < x.size(); i++) os << " " << x[i];
347 if (pad == pad_line) os << std::endl;
348 }
349 }else{
350 os.write((char*) x.data(), x.size() * sizeof(VecType));
351 }
352}
353
358template<typename iomode, typename VecType>
359void readVector(std::istream &is, std::vector<VecType> &x){
360 if (std::is_same<iomode, mode_ascii_type>::value){
361 for(auto &i : x) is >> i;
362 }else{
363 is.read((char*) x.data(), x.size() * sizeof(VecType));
364 }
365}
366
371template<typename iomode, typename VecType, typename SizeType>
372std::vector<VecType> readVector(std::istream &is, SizeType num_entries){
373 std::vector<VecType> x((size_t) num_entries);
375 return x;
376}
377
382template<bool iomode, IOPad pad, typename... Vals>
383void writeNumbers(std::ostream &os, Vals... vals){
384 std::vector<typename std::tuple_element<0, std::tuple<Vals...>>::type> values = {vals...};
386}
387
392template<typename iomode, typename Val>
393Val readNumber(std::istream &is){
394 Val v;
395 if (std::is_same<iomode, mode_ascii_type>::value){
396 is >> v;
397 }else{
398 is.read((char*) &v, sizeof(Val));
399 }
400 return v;
401}
402
407template<bool iomode>
408void writeRule(TypeOneDRule rule, std::ostream &os){
409 if (iomode == mode_ascii){
410 os << getRuleString(rule) << std::endl;
411 }else{
412 int r = getRuleInt(rule);
413 os.write((char*) &r, sizeof(int));
414 }
415}
416
421template<typename iomode>
422TypeOneDRule readRule(std::istream &is){
423 if (std::is_same<iomode, mode_ascii_type>::value){
424 std::string T;
425 is >> T;
426 return getRuleString(T);
427 }else{
429 }
430}
431
432}
433
434}
435
436#endif
TypeOneDRule
Used to specify the one dimensional family of rules that induces the sparse grid.
Definition tsgEnumerates.hpp:285
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_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_localp
Nested rule with a hierarchy of uniformly distributed nodes and functions with compact support.
Definition tsgEnumerates.hpp:362
@ 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_wavelet
Wavelet basis with uniformly distributed nodes (primarily for internal use).
Definition tsgEnumerates.hpp:370
@ 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_localpb
Variation of rule_localp focusing nodes on the boundary instead of the interior.
Definition tsgEnumerates.hpp:368
@ 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_localp0
Variation of rule_localp assuming the model is zero at the domain boundary.
Definition tsgEnumerates.hpp:364
@ rule_fourier
Trigonometric basis with uniformly distributed nodes (primarily for internal use).
Definition tsgEnumerates.hpp:372
@ 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_semilocalp
Variation of rule_localp using increased support in exchange for higher order basis (better for smoot...
Definition tsgEnumerates.hpp:366
@ rule_gaussjacobi
Non-nested rule optimized for integral of the form .
Definition tsgEnumerates.hpp:347
@ type_curved
Ignoring the polynomial space, use rules with index .
Definition tsgEnumerates.hpp:213
@ type_iptensor
Make a dense tensor grid with interpolation range that includes .
Definition tsgEnumerates.hpp:245
@ type_qptotal
Total degree polynomial space for quadrature/integration.
Definition tsgEnumerates.hpp:224
@ type_qpcurved
Curved polynomial space for quadrature/integration.
Definition tsgEnumerates.hpp:231
@ type_level
Ignoring the polynomial space, use rules with index .
Definition tsgEnumerates.hpp:209
@ type_qphyperbolic
Hyperbolic cross section polynomial space for quadrature/integration.
Definition tsgEnumerates.hpp:238
@ type_none
Null type, should never be used for input, indicates an error of some sort.
Definition tsgEnumerates.hpp:205
@ type_iptotal
Total degree polynomial space for interpolation, i.e., the span of .
Definition tsgEnumerates.hpp:221
@ type_qptensor
Make a dense tensor grid with quadrature/integration range that includes .
Definition tsgEnumerates.hpp:249
@ type_hyperbolic
Ignoring the polynomial space, use rules with index .
Definition tsgEnumerates.hpp:217
@ type_iphyperbolic
Hyperbolic cross section polynomial space for interpolation, i.e., the span of .
Definition tsgEnumerates.hpp:235
@ type_tensor
Make a dense tensor grid with rules indexed by .
Definition tsgEnumerates.hpp:241
@ type_ipcurved
Curved polynomial space for interpolation, i.e., the span of .
Definition tsgEnumerates.hpp:228
@ refine_parents_first
Isotropic refinement adding children only if the parents are already included.
Definition tsgEnumerates.hpp:429
@ refine_classic
Isotropic refinement using only the children and disregarding missing parents.
Definition tsgEnumerates.hpp:427
@ refine_direction_selective
Anisotropic refinement using only the children and disregarding missing parents.
Definition tsgEnumerates.hpp:431
@ refine_fds
Anisotropic refinement adding children only if the parents are already included.
Definition tsgEnumerates.hpp:433
@ refine_none
Null method, should never be used as input.
Definition tsgEnumerates.hpp:437
@ refine_stable
Isotropic refinement that ensures the points maintain lower-complete structures.
Definition tsgEnumerates.hpp:435
Val readNumber(std::istream &is)
Read a single number, used to read ints (and potentially cast to size_t) or read a double.
Definition tsgIOHelpers.hpp:393
TypeRefinement getTypeRefinementString(std::string const &name)
Map the string to the enumerate hierarchical refinement strategy, used in command line and Python.
Definition tsgIOHelpers.hpp:282
std::vector< TypeOneDRule > getIntRuleMap()
Creates a map with int (used by Fortran and binary I/O) mapped to TypeOneDRule enums.
Definition tsgIOHelpers.hpp:189
void writeVector(const std::vector< VecType > &x, std::ostream &os)
Write the vector to the stream, the vector cannot be empty.
Definition tsgIOHelpers.hpp:338
bool readFlag(std::istream &os)
Read a flag, ascii uses 0 and 1, binary uses characters y and n (counter intuitive,...
Definition tsgIOHelpers.hpp:321
std::map< std::string, TypeOneDRule > getStringRuleMap()
Creates a map with std::string rule names (used by C/Python/CLI) mapped to TypeOneDRule enums.
Definition tsgIOHelpers.hpp:116
TypeDepth getDepthTypeString(std::string const &name)
Map the string to the enumerate multi-index selection strategy, used in command line and Python.
Definition tsgIOHelpers.hpp:247
TypeDepth getDepthTypeInt(int t)
Map the integer to the enumerate multi-index selection strategy, used in Fortran.
Definition tsgIOHelpers.hpp:259
void writeFlag(bool flag, std::ostream &os)
Write the flag to file, ascii uses 0 and 1, binary uses characters y and n (counter intuitive,...
Definition tsgIOHelpers.hpp:305
TypeOneDRule getRuleString(std::string const &name)
Map the string rule name to the enumerate, used in ASCII I/O, command line and Python.
Definition tsgIOHelpers.hpp:167
void writeRule(TypeOneDRule rule, std::ostream &os)
Write a rule.
Definition tsgIOHelpers.hpp:408
std::map< std::string, TypeRefinement > getStringToRefinementMap()
Creates a map with std::string rule names (used by C/Python/CLI) mapped to TypeRefinement enums.
Definition tsgIOHelpers.hpp:270
std::map< std::string, TypeDepth > getStringToDepthMap()
Creates a map with std::string rule names (used by C/Python/CLI) mapped to TypeDepth enums.
Definition tsgIOHelpers.hpp:228
void readVector(std::istream &is, std::vector< VecType > &x)
Read the vector from the stream, the size must already be set.
Definition tsgIOHelpers.hpp:359
void writeNumbers(std::ostream &os, Vals... vals)
Write a bunch of numbers with the same type.
Definition tsgIOHelpers.hpp:383
IOPad
Indicate the type of padding to use, none, space, new-line, etc.
Definition tsgIOHelpers.hpp:99
TypeOneDRule readRule(std::istream &is)
Read a rule.
Definition tsgIOHelpers.hpp:422
TypeRefinement getTypeRefinementInt(int refinement)
Map the integer to the enumerate hierarchical refinement strategy, used by Fortran.
Definition tsgIOHelpers.hpp:293
TypeOneDRule getRuleInt(int r)
Map the int rule index to the enumerate, used in Fortran and binary IO.
Definition tsgIOHelpers.hpp:209
@ pad_lspace
Pad with space at the beginning.
Definition tsgIOHelpers.hpp:105
@ pad_auto
Pad with space if the flag is true, newline if false.
Definition tsgIOHelpers.hpp:109
@ pad_line
Pad with new line.
Definition tsgIOHelpers.hpp:107
@ pad_rspace
Pad with space at the end.
Definition tsgIOHelpers.hpp:103
@ pad_none
Do not add padding.
Definition tsgIOHelpers.hpp:101
constexpr bool mode_binary
Constant allowing for more expressive selection of ascii and binary mode in IO methods.
Definition tsgIOHelpers.hpp:68
constexpr bool mode_ascii
Constant allowing for more expressive selection of ascii and binary mode in IO methods.
Definition tsgIOHelpers.hpp:62
Encapsulates the Tasmanian Sparse Grid module.
Definition TasmanianSparseGrid.hpp:68
Type indicating ascii I/O mode.
Definition tsgIOHelpers.hpp:82
Type indicating binary I/O mode.
Definition tsgIOHelpers.hpp:88
Omnipresent enumerate types.