UQTk: Uncertainty Quantification Toolkit  3.1.1
quad.h
Go to the documentation of this file.
1 /* =====================================================================================
2 
3  The UQ Toolkit (UQTk) version 3.1.1
4  Copyright (2021) NTESS
5  https://www.sandia.gov/UQToolkit/
6  https://github.com/sandialabs/UQTk
7 
8  Copyright 2021 National Technology & Engineering Solutions of Sandia, LLC (NTESS).
9  Under the terms of Contract DE-NA0003525 with NTESS, the U.S. Government
10  retains certain rights in this software.
11 
12  This file is part of The UQ Toolkit (UQTk)
13 
14  UQTk is open source software: you can redistribute it and/or modify
15  it under the terms of BSD 3-Clause License
16 
17  UQTk is distributed in the hope that it will be useful,
18  but WITHOUT ANY WARRANTY; without even the implied warranty of
19  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
20  BSD 3 Clause License for more details.
21 
22  You should have received a copy of the BSD 3 Clause License
23  along with UQTk. If not, see https://choosealicense.com/licenses/bsd-3-clause/.
24 
25  Questions? Contact the UQTk Developers at <uqtk-developers@software.sandia.gov>
26  Sandia National Laboratories, Livermore, CA, USA
27 ===================================================================================== */
31 
32 #ifndef QUAD_H_SEEN
33 #define QUAD_H_SEEN
34 
35 #include "Array1D.h"
36 #include "Array2D.h"
37 
38 #include <iostream>
39 #include <string.h>
40 #include <stdio.h>
41 #include <sstream>
42 using namespace std; // needed for python string conversion
43 
46 #define QD_MAX 20000000
47 
54 class Quad {
55 public:
65  Quad(char *grid_type, char *fs_type,int ndim,int param,double alpha=0.0, double betta=1.0);
66 
77  Quad(Array1D<string>& grid_types, char *fs_type, Array1D<int>& param,Array1D<double>& alphas, Array1D<double>& bettas);
78 
80  Quad() {};
82  ~Quad() {};
83 
85  void init();
86 
88  void SetAlpha(double alpha){this->alpha_=alpha; }
90  void SetBeta(double betta){this->beta_=betta; }
91 
93  void SetDomain(Array1D<double>& aa, Array1D<double>& bb);
95  void SetDomain(Array1D<double>& aa);
97  void GetDomain(Array1D<double>& aa, Array1D<double>& bb) const {aa=aa_; bb=bb_;}
99  void GetDomain(Array1D<double>& aa) const {aa=aa_;}
100 
102  void SetRule(Array2D<double>& q, Array1D<double>& w);
107  //void SetRule(Array2D<double>& q, Array1D<double>& w, Array2D<int>& ind, int param);
108 
110  void SetRule();
111 
113  void GetRule(Array2D<double>& q, Array1D<double>& w);
117 
119  void SetQdpts(Array2D<double>& q){ rule_.qdpts=q; return;}
121  void SetWghts(Array1D<double>& w){ rule_.wghts=w; return;}
123  // void SetIndices(Array2D<int>& ind){ rule_.indices=ind; return;}
124 
126  void GetQdpts(Array2D<double>& q){ q=rule_.qdpts; return;}
128  void GetWghts(Array1D<double>& w){ w=rule_.wghts; return;}
130  //void GetIndices(Array2D<int>& ind){ ind=rule_.indices; return;}
131 
133  void SetLevel(int param) {nlevel_=param; return;}
134 
136  void nextLevel();
137 
139  int GetNQ() {return rule_.qdpts.XSize(); }
140 
143  void SetVerbosity(int verbosity) { quadverbose_ = verbosity; }
144 
145 
146 
147  private:
148 
150  Quad(const Quad &) {};
151 
156 
160 
162  double alpha_;
164  double beta_;
169 
171  typedef struct
172  {
177  } QuadRule;
178 
181 
183  int ndim_;
184 
186  int nlevel_;
187 
192 
196 
198  void MultiplyTwoRules(QuadRule *rule1,QuadRule *rule2,QuadRule *rule_prod);
200  void MultiplyManyRules(int nrules, QuadRule *rules, QuadRule *rule_prod);
201  void MultiplyManyRules_(int nrules, QuadRule *rules, QuadRule *rule_prod);
203  void AddTwoRules(QuadRule *rule1,QuadRule *rule2,QuadRule *rule_sum);
205  void SubtractTwoRules(QuadRule *rule1,QuadRule *rule2,QuadRule *rule_sum);
206 
208  void create1DRule(string gridtype,Array1D<double>& qdpts,Array1D<double>& wghts, int ngr, double a, double b);
209 
212  void create1DRule_CC(Array1D<double>& qdpts,Array1D<double>& wghts, int ngr, double a, double b);
214  void create1DRule_LU(Array1D<double>& qdpts,Array1D<double>& wghts, int ngr, double a, double b);
216  void create1DRule_HG(Array1D<double>& qdpts,Array1D<double>& wghts, int ngr);
218  void create1DRule_NC(Array1D<double>& qdpts,Array1D<double>& wghts, int ngr, double a, double b);
220  void create1DRule_NCO(Array1D<double>& qdpts,Array1D<double>& wghts, int ngr, double a, double b);
223  void create1DRule_CCO(Array1D<double>& qdpts,Array1D<double>& wghts, int ngr, double a, double b);
225  void create1DRule_JB(Array1D<double>& qdpts,Array1D<double>& wghts, int ngr, double a, double b);
227  void create1DRule_LG(Array1D<double>& qdpts,Array1D<double>& wghts, int ngr);
229  void create1DRule_SW(Array1D<double>& qdpts,Array1D<double>& wghts,int ngr);
232  void create1DRule_pdf(Array1D<double>& qdpts,Array1D<double>& wghts, int ngr, double a, double b);
235  void create1DRule_GP3(Array1D<double>& qdpts,Array1D<double>& wghts, int ngr, double a, double b);
236 
238  void getMultiIndexLevel(Array2D<int>& multiIndexLevel, int level, int ndim);
239 
244 
246  string grid_type_;
249 
251  string fs_type_;
252 
254  void compressRule(QuadRule *rule);
255 
256 
257 
258 };
259 #endif /* QUAD_H_SEEN */
1D Array class for any type T
2D Array class for any type T
Definition: Array1D.h:472
Definition: Array1D.h:262
Generates quadrature rules.
Definition: quad.h:54
Array1D< string > grid_types_
Vector of grid types: 'CC','CCO','NC','NCO','LU', 'HG', 'JB', 'LG', 'SW', 'pdf', or 'GP3'.
Definition: quad.h:248
string grid_type_
Grid type: 'CC','CCO','NC','NCO','LU', 'HG', 'JB', 'LG', 'SW', 'pdf', or 'GP3'.
Definition: quad.h:246
void SetRule(Array2D< double > &q, Array1D< double > &w, Array2D< int > &ind)
Set the rule externally (quadrature points, weights and indices) Dummy function for backward compatib...
Definition: quad.h:105
Array2D< int > npts_1_all
Definition: quad.h:194
Array1D< int > param_
Definition: quad.h:191
int ndim_
The dimensionality.
Definition: quad.h:183
Array1D< double > alphas_
The first parameter of the rule, if any.
Definition: quad.h:166
void SetQdpts(Array2D< double > &q)
Externally set quadrature points.
Definition: quad.h:119
double alpha_
The first parameter of the rule, if any.
Definition: quad.h:162
void SetBeta(double betta)
Set the parameter beta.
Definition: quad.h:90
~Quad()
Destructor.
Definition: quad.h:82
Array1D< double > bb_
the right endpoints of the domain
Definition: quad.h:155
int maxlevel_
The level for sparse rules, or the number of grid points per dim for full product rules.
Definition: quad.h:190
int GetNQ()
Get the number of quadrature points.
Definition: quad.h:139
int quadverbose_
Verbosity level.
Definition: quad.h:159
void GetDomain(Array1D< double > &aa, Array1D< double > &bb) const
Get the domain endpoints (for compact support domains)
Definition: quad.h:97
void GetQdpts(Array2D< double > &q)
Externally set the indices.
Definition: quad.h:126
void SetLevel(int param)
Get the indices.
Definition: quad.h:133
Quad(const Quad &)
Dummy copy constructor, which should not be used as it is currently not well defined.
Definition: quad.h:150
int nlevel_
The current level, working variable for hierarchical construction.
Definition: quad.h:186
void SetAlpha(double alpha)
Set the parameter alpha.
Definition: quad.h:88
int growth_rule_
Growth rule: exponential(0) or linear(1)
Definition: quad.h:241
Array1D< double > aa_
The left endpoints of the domain.
Definition: quad.h:150
Array2D< QuadRule > qr_1_all
Definition: quad.h:195
QuadRule rule_
The quadrature rule structure.
Definition: quad.h:180
void SetVerbosity(int verbosity)
Set the verbosity level.
Definition: quad.h:143
string fs_type_
Sparseness type (full or sparse)
Definition: quad.h:251
Array1D< double > betas_
The second parameter of the rule, if any.
Definition: quad.h:168
void GetRule(Array2D< double > &q, Array1D< double > &w, Array2D< int > &ind)
Get the quadrature rule with indexing Dummy function for backward compatibility.
Definition: quad.h:116
void GetDomain(Array1D< double > &aa) const
Get the domain endpoint (for semi-infinite domains)
Definition: quad.h:99
double beta_
The second parameter of the rule, if any.
Definition: quad.h:164
void GetWghts(Array1D< double > &w)
Get the weights.
Definition: quad.h:128
Quad()
Constructor: empty.
Definition: quad.h:80
Array1D< int > growth_rules_
Growth rules: exponential(0) or linear(1)
Definition: quad.h:243
void SetWghts(Array1D< double > &w)
Externally set the weights.
Definition: quad.h:121
Rule structure that stores quadrature points, weights and indices.
Definition: quad.h:172
Array2D< double > qdpts
Quadrature points.
Definition: quad.h:174
Array1D< double > wghts
Quadrature weights.
Definition: quad.h:176