Intrepid
Intrepid_CubatureTensorSorted.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ************************************************************************
3 //
4 // Intrepid Package
5 // Copyright (2007) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Pavel Bochev (pbboche@sandia.gov)
38 // Denis Ridzal (dridzal@sandia.gov), or
39 // Kara Peterson (kjpeter@sandia.gov)
40 //
41 // ************************************************************************
42 // @HEADER
43 
49 #ifndef INTREPID_CUBATURE_TENSORSORTED_HPP
50 #define INTREPID_CUBATURE_TENSORSORTED_HPP
51 
53 
54 namespace Intrepid {
55 
73 /* =========================================================================
74  Structure Definition - quadRule
75  ========================================================================= */
76 //struct quadInfo {
77  /*
78  The quadRule structure contains information about a multidimensional
79  tensor product quadrature rule.
80  */
81  // int dim; // Number of Spatial Dimensions
82  // int maxlevel; // Maximum Order of 1D rules
83  // int * rule; // 1D rules defined in each dimension
84  // int * growth; // 1D growth rules defined in each dimension
85  // int * np; // Number of parameters for each 1D rule
86  // double * param; // Parameters for each 1D rule
87  //};
88 
89 template<class Scalar, class ArrayPoint = FieldContainer<Scalar>, class ArrayWeight = ArrayPoint>
90 class CubatureTensorSorted : public Intrepid::Cubature<Scalar,ArrayPoint,ArrayWeight> {
91 
92 private:
95  typename std::map<std::vector<Scalar>,int> points_; // keys = nodes, values = location of weights
96 
99  std::vector<Scalar> weights_;
100 
104 
108  std::vector<int> degree_;
109 
113 
114 public:
115 
117 
118  // CubatureTensorSorted() {}
119 
120  CubatureTensorSorted(int numPoints = 0, int dimension = 1);
121 
127 
135  CubatureTensorSorted(int dimension, std::vector<int> numPoints1D, std::vector<EIntrepidBurkardt> rule1D, bool isNormalized);
136 
144  CubatureTensorSorted(int dimension, std::vector<int> numPoints1D, std::vector<EIntrepidBurkardt> rule1D, std::vector<EIntrepidGrowth> growth1D, bool isNormalized);
145 
154  CubatureTensorSorted(int dimension, int maxNumPoints, std::vector<EIntrepidBurkardt> rule1D, std::vector<EIntrepidGrowth> growth1D, bool isNormalized);
155 
156  // Access Operator
163  void getCubature(ArrayPoint & cubPoints,
164  ArrayWeight & cubWeights) const;
165 
173  void getCubature(ArrayPoint& cubPoints,
174  ArrayWeight& cubWeights,
175  ArrayPoint& cellCoords) const;
176 
179  int getNumPoints() const;
180 
184  void getAccuracy(std::vector<int> & accuracy) const;
185 
188  int getDimension() const;
189 
192  typename std::map<std::vector<Scalar>,int>::iterator begin();
193 
196  typename std::map<std::vector<Scalar>,int>::iterator end();
197 
200  void insert(typename std::map<std::vector<Scalar>,int>::iterator it,
201  std::vector<Scalar> point, Scalar weight);
202 
205  std::vector<Scalar> getNode(typename std::map<std::vector<Scalar>,int>::iterator it);
206 
209  Scalar getWeight(int node);
210 
213  Scalar getWeight(std::vector<Scalar> point);
214 
218  void update(Scalar alpha2, CubatureTensorSorted<Scalar> & cubRule2, Scalar alpha1);
219 
222  void normalize();
223 
224 };
225 
226 template<class Scalar>
227 CubatureTensorSorted<Scalar> kron_prod(CubatureTensorSorted<Scalar> & rule1,
228  CubatureLineSorted<Scalar> & rule2 );
229 
230 } // Intrepid Namespace
231 
232 
233 // include templated definitions
235 
236 #endif
int getNumPoints() const
Returns the number of cubature points.
void normalize()
Normalize CubatureLineSorted weights.
std::map< std::vector< Scalar >, int >::iterator begin()
Initiate iterator at the beginning of data.
Utilizes cubature (integration) rules contained in the library sandia_rules (John Burkardt...
Scalar getWeight(int node)
Get a specific weight described by the integer location.
std::vector< Scalar > getNode(typename std::map< std::vector< Scalar >, int >::iterator it)
Get a specific node described by the iterator location.
void update(Scalar alpha2, CubatureTensorSorted< Scalar > &cubRule2, Scalar alpha1)
Replace CubatureLineSorted values with &quot;this = alpha1*this+alpha2*cubRule2&quot;.
Header file for the Intrepid::CubatureLineSorted class.
std::vector< int > degree_
The degree of polynomials that are integrated exactly by this cubature rule.
int dimension_
Dimension of integration domain.
Definition file for the Intrepid::CubatureTensorSorted class.
void insert(typename std::map< std::vector< Scalar >, int >::iterator it, std::vector< Scalar > point, Scalar weight)
Insert a node and weight into data near the iterator position.
Defines the base class for cubature (integration) rules in Intrepid.
int getDimension() const
Returns dimension of domain of integration.
std::map< std::vector< Scalar >, int > points_
Contains nodes of this cubature rule.
void getAccuracy(std::vector< int > &accuracy) const
Returns max. degree of polynomials that are integrated exactly. The return vector has size 1...
Utilizes 1D cubature (integration) rules contained in the library sandia_rules (John Burkardt...
void getCubature(ArrayPoint &cubPoints, ArrayWeight &cubWeights) const
Returns cubature points and weights (return arrays must be pre-sized/pre-allocated).
std::map< std::vector< Scalar >, int >::iterator end()
Initiate iterator at the end of data.
int numPoints_
Contains the number of nodes for this cubature rule.
std::vector< Scalar > weights_
Contains weights of this cubature rule.