| Intrepid
    | 
Utilizes 1D cubature (integration) rules contained in the library sandia_rules (John Burkardt, Scientific Computing, Florida State University) within Intrepid. More...
#include <Intrepid_CubatureTensorSorted.hpp>
 
  
 | Public Member Functions | |
| CubatureTensorSorted (int numPoints=0, int dimension=1) | |
| CubatureTensorSorted (CubatureLineSorted< Scalar > &cubLine) | |
| Constructor.  More... | |
| CubatureTensorSorted (int dimension, std::vector< int > numPoints1D, std::vector< EIntrepidBurkardt > rule1D, bool isNormalized) | |
| Constructor.  More... | |
| CubatureTensorSorted (int dimension, std::vector< int > numPoints1D, std::vector< EIntrepidBurkardt > rule1D, std::vector< EIntrepidGrowth > growth1D, bool isNormalized) | |
| Constructor.  More... | |
| CubatureTensorSorted (int dimension, int maxNumPoints, std::vector< EIntrepidBurkardt > rule1D, std::vector< EIntrepidGrowth > growth1D, bool isNormalized) | |
| Constructor.  More... | |
| void | getCubature (ArrayPoint &cubPoints, ArrayWeight &cubWeights) const | 
| Returns cubature points and weights (return arrays must be pre-sized/pre-allocated).  More... | |
| void | getCubature (ArrayPoint &cubPoints, ArrayWeight &cubWeights, ArrayPoint &cellCoords) const | 
| Returns cubature points and weights. Method for physical space cubature, throws an exception.  More... | |
| int | getNumPoints () const | 
| Returns the number of cubature points. | |
| void | getAccuracy (std::vector< int > &accuracy) const | 
| Returns max. degree of polynomials that are integrated exactly. The return vector has size 1. | |
| int | getDimension () const | 
| Returns dimension of domain of integration. | |
| std::map< std::vector< Scalar > , int >::iterator | begin () | 
| Initiate iterator at the beginning of data. | |
| std::map< std::vector< Scalar > , int >::iterator | end () | 
| Initiate iterator at the end of data. | |
| 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. | |
| std::vector< Scalar > | getNode (typename std::map< std::vector< Scalar >, int >::iterator it) | 
| Get a specific node described by the iterator location. | |
| Scalar | getWeight (int node) | 
| Get a specific weight described by the integer location. | |
| Scalar | getWeight (std::vector< Scalar > point) | 
| Get a specific weight described by the iterator location. | |
| void | update (Scalar alpha2, CubatureTensorSorted< Scalar > &cubRule2, Scalar alpha1) | 
| Replace CubatureLineSorted values with "this = alpha1*this+alpha2*cubRule2". | |
| void | normalize () | 
| Normalize CubatureLineSorted weights. | |
| Private Attributes | |
| std::map< std::vector< Scalar > , int > | points_ | 
| Contains nodes of this cubature rule. | |
| std::vector< Scalar > | weights_ | 
| Contains weights of this cubature rule. | |
| int | numPoints_ | 
| Contains the number of nodes for this cubature rule. | |
| std::vector< int > | degree_ | 
| The degree of polynomials that are integrated exactly by this cubature rule. | |
| int | dimension_ | 
| Dimension of integration domain. | |
Utilizes 1D cubature (integration) rules contained in the library sandia_rules (John Burkardt, Scientific Computing, Florida State University) within Intrepid.
This class contains ten rules:
Definition at line 90 of file Intrepid_CubatureTensorSorted.hpp.
| Intrepid::CubatureTensorSorted< Scalar, ArrayPoint, ArrayWeight >::CubatureTensorSorted | ( | CubatureLineSorted< Scalar > & | cubLine | ) | 
Constructor.
| cubLine | [in] - 1D cubature rule. | 
Definition at line 64 of file Intrepid_CubatureTensorSortedDef.hpp.
References Intrepid::CubatureLineSorted< Scalar, ArrayPoint, ArrayWeight >::begin(), Intrepid::CubatureLineSorted< Scalar, ArrayPoint, ArrayWeight >::end(), Intrepid::CubatureLineSorted< Scalar, ArrayPoint, ArrayWeight >::getAccuracy(), Intrepid::CubatureLineSorted< Scalar, ArrayPoint, ArrayWeight >::getNode(), Intrepid::CubatureLineSorted< Scalar, ArrayPoint, ArrayWeight >::getNumPoints(), and Intrepid::CubatureLineSorted< Scalar, ArrayPoint, ArrayWeight >::getWeight().
| Intrepid::CubatureTensorSorted< Scalar, ArrayPoint, ArrayWeight >::CubatureTensorSorted | ( | int | dimension, | 
| std::vector< int > | numPoints1D, | ||
| std::vector< EIntrepidBurkardt > | rule1D, | ||
| bool | isNormalized | ||
| ) | 
Constructor.
| dimension | [in] - The number of spatial dimensions. | 
| numPoints1D | [in] - The number of cubature points in each direction. | 
| rule1D | [in] - The cubature rule for each direction. | 
| isNormalized | [in] - Flag whether or not to normalize the cubature weights. | 
Definition at line 87 of file Intrepid_CubatureTensorSortedDef.hpp.
References Intrepid::CubatureTensorSorted< Scalar, ArrayPoint, ArrayWeight >::begin(), Intrepid::CubatureTensorSorted< Scalar, ArrayPoint, ArrayWeight >::end(), Intrepid::CubatureLineSorted< Scalar, ArrayPoint, ArrayWeight >::getAccuracy(), Intrepid::CubatureTensorSorted< Scalar, ArrayPoint, ArrayWeight >::getNumPoints(), and Intrepid::CubatureTensorSorted< Scalar, ArrayPoint, ArrayWeight >::getWeight().
| Intrepid::CubatureTensorSorted< Scalar, ArrayPoint, ArrayWeight >::CubatureTensorSorted | ( | int | dimension, | 
| std::vector< int > | numPoints1D, | ||
| std::vector< EIntrepidBurkardt > | rule1D, | ||
| std::vector< EIntrepidGrowth > | growth1D, | ||
| bool | isNormalized | ||
| ) | 
Constructor.
| dimension | [in] - The number of spatial dimensions. | 
| numPoints1D | [in] - The number of cubature points in each direction. | 
| rule1D | [in] - The cubature rule for each direction. | 
| isNormalized | [in] - Flag whether or not to normalize the cubature weights. | 
Definition at line 123 of file Intrepid_CubatureTensorSortedDef.hpp.
References Intrepid::CubatureTensorSorted< Scalar, ArrayPoint, ArrayWeight >::begin(), Intrepid::CubatureTensorSorted< Scalar, ArrayPoint, ArrayWeight >::end(), Intrepid::CubatureLineSorted< Scalar, ArrayPoint, ArrayWeight >::getAccuracy(), Intrepid::CubatureTensorSorted< Scalar, ArrayPoint, ArrayWeight >::getNumPoints(), and Intrepid::CubatureTensorSorted< Scalar, ArrayPoint, ArrayWeight >::getWeight().
| Intrepid::CubatureTensorSorted< Scalar, ArrayPoint, ArrayWeight >::CubatureTensorSorted | ( | int | dimension, | 
| int | maxNumPoints, | ||
| std::vector< EIntrepidBurkardt > | rule1D, | ||
| std::vector< EIntrepidGrowth > | growth1D, | ||
| bool | isNormalized | ||
| ) | 
Constructor.
| dimension | [in] - The number of spatial dimensions. | 
| maxNumPoints | [in] - The maximum number of cubature points in each direction. | 
| rule1D | [in] - The cubature rule for each direction. | 
| growth1D | [in] - Growth rule for each direction. | 
| isNormalized | [in] - Flag whether or not to normalize the cubature weights. | 
Definition at line 163 of file Intrepid_CubatureTensorSortedDef.hpp.
References Intrepid::CubatureTensorSorted< Scalar, ArrayPoint, ArrayWeight >::begin(), Intrepid::CubatureTensorSorted< Scalar, ArrayPoint, ArrayWeight >::end(), Intrepid::CubatureLineSorted< Scalar, ArrayPoint, ArrayWeight >::getAccuracy(), Intrepid::CubatureTensorSorted< Scalar, ArrayPoint, ArrayWeight >::getNumPoints(), and Intrepid::CubatureTensorSorted< Scalar, ArrayPoint, ArrayWeight >::getWeight().
| 
 | virtual | 
Returns cubature points and weights (return arrays must be pre-sized/pre-allocated).
| cubPoints | [out] - Array containing the cubature points. | 
| cubWeights | [out] - Array of corresponding cubature weights. | 
Implements Intrepid::Cubature< Scalar, ArrayPoint, ArrayWeight >.
Definition at line 216 of file Intrepid_CubatureTensorSortedDef.hpp.
| 
 | virtual | 
Returns cubature points and weights. Method for physical space cubature, throws an exception.
| cubPoints | [out] - Array containing the cubature points. | 
| cubWeights | [out] - Array of corresponding cubature weights. | 
| cellCoords | [in] - Array of cell coordinates | 
Implements Intrepid::Cubature< Scalar, ArrayPoint, ArrayWeight >.
Definition at line 241 of file Intrepid_CubatureTensorSortedDef.hpp.
 1.8.5
 1.8.5