Intrepid2
Intrepid2_CubatureControlVolume.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ************************************************************************
3 //
4 // Intrepid2 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 Kyungjoo Kim (kyukim@sandia.gov), or
38 // Mauro Perego (mperego@sandia.gov)
39 //
40 // ************************************************************************
41 // @HEADER
42 
49 #ifndef __INTREPID2_CUBATURE_CONTROLVOLUME_HPP__
50 #define __INTREPID2_CUBATURE_CONTROLVOLUME_HPP__
51 
52 #include "Intrepid2_ConfigDefs.hpp"
53 #include "Intrepid2_Cubature.hpp"
54 
55 #include "Shards_CellTopology.hpp"
56 #include "Intrepid2_CellTools.hpp"
58 
59 namespace Intrepid2{
60 
67  template<typename ExecSpaceType = void,
68  typename pointValueType = double,
69  typename weightValueType = double>
71  : public Cubature<ExecSpaceType,pointValueType,weightValueType> {
72  public:
73 
74  template<typename cubPointViewType,
75  typename cubWeightViewType,
76  typename subcvCoordViewType,
77  typename subcvWeightViewType,
78  typename jacDetViewType>
79  struct Functor {
80  cubPointViewType _cubPoints;
81  cubWeightViewType _cubWeights;
82  const subcvCoordViewType _subcvCoords;
83  const subcvWeightViewType _subcvWeights;
84  const jacDetViewType _jacDets;
85 
86  KOKKOS_INLINE_FUNCTION
87  Functor( cubPointViewType cubPoints_,
88  cubWeightViewType cubWeights_,
89  subcvCoordViewType subcvCoords_,
90  subcvWeightViewType subcvWeights_,
91  jacDetViewType jacDets_ )
92  : _cubPoints(cubPoints_), _cubWeights(cubWeights_),
93  _subcvCoords(subcvCoords_), _subcvWeights(subcvWeights_), _jacDets(jacDets_) {}
94 
95  KOKKOS_INLINE_FUNCTION
96  void operator()(const ordinal_type cell) const {
97  const ordinal_type numNodesPerCell = _subcvCoords.extent(1);
98  const ordinal_type numNodesPerSubcv = _subcvCoords.extent(2);
99  const ordinal_type spaceDim = _subcvCoords.extent(3);
100  const ordinal_type numSubcvPoints = _subcvWeights.extent(0);
101 
102  // compute subcv centers
103  for (ordinal_type node=0;node<numNodesPerCell;++node) {
104  typename cubPointViewType::value_type val[3] = {};
105  for (ordinal_type subcv=0;subcv<numNodesPerSubcv;++subcv) {
106  for (ordinal_type i=0;i<spaceDim;++i)
107  val[i] += _subcvCoords(cell, node, subcv, i);
108  }
109  for (ordinal_type i=0;i<spaceDim;++i)
110  _cubPoints(cell, node, i) = (val[i]/numNodesPerSubcv);
111  }
112 
113  // compute weights (area or volume)
114  for (ordinal_type node=0;node<numNodesPerCell;++node) {
115  typename cubWeightViewType::value_type val = 0;
116  for (ordinal_type pt=0;pt<numSubcvPoints;++pt)
117  val += _subcvWeights(pt)*_jacDets(cell, node, pt);
118  _cubWeights(cell, node) = val;
119  }
120  }
121  };
122 
123  protected:
124 
127  shards::CellTopology primaryCellTopo_;
128 
131  shards::CellTopology subcvCellTopo_;
132 
135  ordinal_type degree_;
136 
137  // cubature points and weights associated with sub-control volume.
138  Kokkos::DynRankView<pointValueType, ExecSpaceType> subcvCubaturePoints_;
139  Kokkos::DynRankView<weightValueType,ExecSpaceType> subcvCubatureWeights_;
140 
141  public:
142  typedef typename Cubature<ExecSpaceType,pointValueType,weightValueType>::PointViewType PointViewType;
143  typedef typename Cubature<ExecSpaceType,pointValueType,weightValueType>::weightViewType weightViewType;
144 
146 
154  virtual
155  void
156  getCubature( PointViewType cubPoints,
157  weightViewType cubWeights,
158  PointViewType cellCoords) const;
159 
162  virtual
163  ordinal_type
164  getNumPoints() const {
165  return primaryCellTopo_.getNodeCount();
166  }
167 
170  virtual
171  ordinal_type
172  getDimension() const {
173  return primaryCellTopo_.getDimension();
174  }
175 
178  virtual
179  const char*
180  getName() const {
181  return "CubatureControlVolume";
182  }
183 
187  CubatureControlVolume(const shards::CellTopology cellTopology);
188  virtual ~CubatureControlVolume() {}
189 
190  }; // end class CubatureControlVolume
191 
192 } // end namespace Intrepid2
193 
195 
196 #endif
197 
Header file for the abstract base class Intrepid2::DefaultCubatureFactory.
Defines the base class for cubature (integration) rules in Intrepid.
virtual void getCubature(PointViewType cubPoints, weightViewType cubWeights, PointViewType cellCoords) const
Returns cubature points and weights (return arrays must be pre-sized/pre-allocated).
Header file for the Intrepid2::CubatureControlVolume class.
ordinal_type degree_
The degree of the polynomials that are integrated exactly.
Header file for the Intrepid2::Cubature class.
Defines cubature (integration) rules over control volumes.
virtual ordinal_type getNumPoints() const
Returns the number of cubature points.
shards::CellTopology primaryCellTopo_
The topology of the primary cell.
virtual const char * getName() const
Returns cubature name.
shards::CellTopology subcvCellTopo_
The topology of the sub-control volume.
Header file for the Intrepid2::CellTools class.
CubatureControlVolume(const shards::CellTopology cellTopology)
virtual ordinal_type getDimension() const
Returns dimension of integration domain.