Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_IntegrationRule.cpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Panzer: A partial differential equation assembly
4// engine for strongly coupled complex multiphysics systems
5//
6// Copyright 2011 NTESS and the Panzer contributors.
7// SPDX-License-Identifier: BSD-3-Clause
8// *****************************************************************************
9// @HEADER
10
12
13#include "Teuchos_ArrayRCP.hpp"
14#include "Teuchos_Assert.hpp"
15#include "Phalanx_DataLayout_MDALayout.hpp"
16#include "Intrepid2_DefaultCubatureFactory.hpp"
17#include "Intrepid2_CubatureControlVolume.hpp"
18#include "Intrepid2_CubatureControlVolumeSide.hpp"
19#include "Intrepid2_CubatureControlVolumeBoundary.hpp"
20#include "Panzer_Dimension.hpp"
21#include "Panzer_CellData.hpp"
22
23
25IntegrationRule(int in_cubature_degree, const panzer::CellData& cell_data)
27{
28 if(cell_data.isSide()){
29 IntegrationDescriptor::setup(in_cubature_degree, IntegrationDescriptor::SIDE, cell_data.side());
30 } else {
32 }
33 setup(in_cubature_degree,cell_data);
34}
35
37IntegrationRule(const panzer::CellData& cell_data, const std::string & in_cv_type)
39{
40 // Cubature orders are only used for indexing so we make them large enough not to interfere with other rules.
41 // TODO: This requirement (on arbitrary cubature order) will be dropped with the new Workset design (using descriptors to find integration rules)
42 if(in_cv_type == "volume"){
43 TEUCHOS_TEST_FOR_EXCEPT_MSG(cell_data.isSide(),
44 "IntegrationRule::IntegrationRule : Control Volume 'volume' type requested, but CellData is setup for sides.");
46 } else if(in_cv_type == "side"){
48 } else if(in_cv_type == "boundary"){
49 TEUCHOS_TEST_FOR_EXCEPT_MSG(not cell_data.isSide(),
50 "IntegrationRule::IntegrationRule : Control Volume 'boundary' type requested, but CellData is not setup for sides.");
52 } else {
53 TEUCHOS_ASSERT(false);
54 }
55 setup_cv(cell_data,in_cv_type);
56}
57
60 const Teuchos::RCP<const shards::CellTopology> & cell_topology,
61 const int num_cells,
62 const int num_faces)
63 : PointRule(), IntegrationDescriptor(description)
64{
65
66 TEUCHOS_ASSERT(description.getType() != panzer::IntegrationDescriptor::NONE);
67
68 cubature_degree = description.getOrder();
69 cv_type = "none";
70
71 panzer::CellData cell_data;
72 if(isSide()){
73 cell_data = panzer::CellData(num_cells, getSide(), cell_topology);
74 } else {
75 cell_data = panzer::CellData(num_cells, cell_topology);
76 }
77
79 setup(getOrder(), cell_data);
80 } else if(description.getType() == panzer::IntegrationDescriptor::SIDE){
81 setup(getOrder(), cell_data);
82 } else if(description.getType() == panzer::IntegrationDescriptor::SURFACE){
83 TEUCHOS_ASSERT(num_faces!=-1);
84 setup_surface(cell_topology, num_cells, num_faces);
85 } else if(description.getType() == panzer::IntegrationDescriptor::CV_VOLUME){
86 setup_cv(cell_data, "volume");
87 } else if(description.getType() == panzer::IntegrationDescriptor::CV_SIDE){
88 setup_cv(cell_data, "side");
89 } else if(description.getType() == panzer::IntegrationDescriptor::CV_BOUNDARY){
90 setup_cv(cell_data, "boundary");
91 } else {
92 TEUCHOS_ASSERT(false);
93 }
94
95}
96
97void panzer::IntegrationRule::setup(int in_cubature_degree, const panzer::CellData& cell_data)
98{
99 cubature_degree = in_cubature_degree;
100 cv_type = "none";
101 int spatialDimension = cell_data.baseCellDimension();
102
103 std::stringstream ss;
104 ss << "CubaturePoints (Degree=" << cubature_degree;
105
106 // Intrepid2 does not support a quadrature on a 0-dimensional object
107 // (which doesn't make much sense anyway) to work around this we
108 // will adjust the integration rule manually
109 if(cell_data.isSide() && spatialDimension==1) {
110 ss << ",side)";
111 PointRule::setup(ss.str(),1,cell_data);
112
113 return;
114 }
115
116 const shards::CellTopology & topo = *cell_data.getCellTopology();
117 Teuchos::RCP<shards::CellTopology> sideTopo = getSideTopology(cell_data);
118
119 Intrepid2::DefaultCubatureFactory cubature_factory;
120 Teuchos::RCP<Intrepid2::Cubature<PHX::Device::execution_space,double,double>> intrepid_cubature;
121
122 // get side topology
123 if (Teuchos::is_null(sideTopo)) {
124 ss << ",volume)";
125 intrepid_cubature = cubature_factory.create<PHX::Device::execution_space,double,double>(topo, cubature_degree);
126 }
127 else {
128 ss << ",side)";
129 intrepid_cubature = cubature_factory.create<PHX::Device::execution_space,double,double>(*sideTopo, cubature_degree);
130 }
131
132 PointRule::setup(ss.str(),intrepid_cubature->getNumPoints(),cell_data);
133}
134
135void panzer::IntegrationRule::setup_surface(const Teuchos::RCP<const shards::CellTopology> & cell_topology, const int num_cells, const int num_faces)
136{
137
138 const int cell_dim = cell_topology->getDimension();
139 const int subcell_dim = cell_dim-1;
140 const int num_faces_per_cell = cell_topology->getSubcellCount(subcell_dim);
141
142 panzer::CellData cell_data(num_cells, cell_topology);
143
144 std::string point_rule_name;
145 {
146 std::stringstream ss;
147 ss << "CubaturePoints (Degree=" << getOrder() << ",surface)";
148 point_rule_name = ss.str();
149 }
150
151 // We can skip some steps for 1D
152 if(cell_dim == 1){
153 const int num_points_per_cell = num_faces_per_cell;
154 const int num_points_per_face = 1;
155 PointRule::setup(point_rule_name, num_cells, num_points_per_cell, num_faces, num_points_per_face, cell_topology);
156 _point_offsets.resize(3,0);
157 _point_offsets[0] = 0;
158 _point_offsets[1] = num_points_per_face;
159 _point_offsets[2] = _point_offsets[1]+num_points_per_face;
160 return;
161 }
162
163 Intrepid2::DefaultCubatureFactory cubature_factory;
164
165 _point_offsets.resize(num_faces_per_cell+1,0);
166 int test_face_size = -1;
167 for(int subcell_index=0; subcell_index<num_faces_per_cell; ++subcell_index){
168 Teuchos::RCP<shards::CellTopology> face_topology = Teuchos::rcp(new shards::CellTopology(cell_topology->getCellTopologyData(subcell_dim,subcell_index)));
169 const auto & intrepid_cubature = cubature_factory.create<PHX::Device::execution_space,double,double>(*face_topology, getOrder());
170 const int num_face_points = intrepid_cubature->getNumPoints();
171 _point_offsets[subcell_index+1] = _point_offsets[subcell_index] + num_face_points;
172
173 // Right now we only support each face having the same number of points
174 if(test_face_size==-1){
175 test_face_size = num_face_points;
176 } else {
177 TEUCHOS_ASSERT(num_face_points == test_face_size);
178 }
179 }
180
181 const int num_points_per_cell = _point_offsets.back();
182 const int num_points_per_face = _point_offsets[1];
183
184 PointRule::setup(point_rule_name, num_cells, num_points_per_cell, num_faces, num_points_per_face, cell_topology);
185
186}
187
188void panzer::IntegrationRule::setup_cv(const panzer::CellData& cell_data, std::string in_cv_type)
189{
190 // set cubature degree to arbitrary constant for indexing
191 cubature_degree = 1;
192 cv_type = in_cv_type;
193 if (cv_type == "volume") {
194 cubature_degree = 75;
195 }
196 if (cv_type == "side") {
197 cubature_degree = 85;
198 }
199 if (cv_type == "boundary") {
200 cubature_degree = 95;
201 }
202
203 //int spatialDimension = cell_data.baseCellDimension();
204
205 std::stringstream ss;
206 ss << "CubaturePoints ControlVol (Index=" << cubature_degree;
207
208 const shards::CellTopology & topo = *cell_data.getCellTopology();
209
210 Teuchos::RCP<Intrepid2::Cubature<PHX::Device::execution_space,double,double> > intrepid_cubature;
211
212 int tmp_num_points = 0;
213 if (cv_type == "volume") {
214 ss << ",volume)";
215 intrepid_cubature = Teuchos::rcp(new Intrepid2::CubatureControlVolume<PHX::Device::execution_space,double,double>(topo));
216 tmp_num_points = intrepid_cubature->getNumPoints();
217 }
218 else if (cv_type == "side") {
219 ss << ",side)";
220 intrepid_cubature = Teuchos::rcp(new Intrepid2::CubatureControlVolumeSide<PHX::Device::execution_space,double,double>(topo));
221 tmp_num_points = intrepid_cubature->getNumPoints();
222 }
223 else if (cv_type == "boundary") {
224 ss << ",boundary)";
225 intrepid_cubature = Teuchos::rcp(new Intrepid2::CubatureControlVolumeBoundary<PHX::Device::execution_space,double,double>(topo,cell_data.side()));
226 tmp_num_points = intrepid_cubature->getNumPoints();
227 }
228
229 PointRule::setup(ss.str(),tmp_num_points,cell_data);
230}
231
233{ return cubature_degree; }
234
235
236int panzer::IntegrationRule::getPointOffset(const int subcell_index) const
237{
238 // Need to make sure this is a surface integrator
239 TEUCHOS_ASSERT(getType() == SURFACE);
240 return _point_offsets[subcell_index];
241}
242
243
244void panzer::IntegrationRule::print(std::ostream & os)
245{
246 os << "IntegrationRule ( "
247 << "Name = " << getName()
248 << ", Degree = " << cubature_degree
249 << ", Dimension = " << spatial_dimension
250 << ", Workset Size = " << workset_size
251 << ", Num Points = " << num_points
252 << ", Side = " << side
253 << " )";
254}
255
256void panzer::IntegrationRule::referenceCoordinates(Kokkos::DynRankView<double,PHX::Device> & cub_points)
257{
258 // build an interpid cubature rule
259 Teuchos::RCP< Intrepid2::Cubature<PHX::Device::execution_space,double,double> > intrepid_cubature;
260 Intrepid2::DefaultCubatureFactory cubature_factory;
261
262 if (!isSide())
263 intrepid_cubature = cubature_factory.create<PHX::Device::execution_space,double,double>(*(topology),cubature_degree);
264 else
265 intrepid_cubature = cubature_factory.create<PHX::Device::execution_space,double,double>(*(side_topology),cubature_degree);
266
267 int num_ip = intrepid_cubature->getNumPoints();
268 Kokkos::DynRankView<double,PHX::Device> cub_weights("cub_weights",num_ip);
269
270 // now compute weights (and throw them out) as well as reference points
271 if (!isSide()) {
272 cub_points = Kokkos::DynRankView<double,PHX::Device>("cub_points", num_ip, topology->getDimension());
273 intrepid_cubature->getCubature(cub_points, cub_weights);
274 }
275 else {
276 cub_points = Kokkos::DynRankView<double,PHX::Device>("cub_points", num_ip, side_topology->getDimension());
277 intrepid_cubature->getCubature(cub_points, cub_weights);
278 }
279}
Data for determining cell topology and dimensionality.
int baseCellDimension() const
Dimension of the base cell. NOT the dimension of the local side, even if the side() method returns tr...
Teuchos::RCP< const shards::CellTopology > getCellTopology() const
Get CellTopology for the base cell.
@ VOLUME
No integral specified - default state.
@ CV_VOLUME
Integral over a specific side of cells (side must be set)
@ CV_BOUNDARY
Control volume side integral.
@ SIDE
Integral over all sides of cells (closed surface integral)
const int & getOrder() const
Get order of integrator.
const int & getSide() const
Get side associated with integration - this is for backward compatibility.
const int & getType() const
Get type of integrator.
void setup(const int cubature_order, const int integration_type, const int side=-1)
Setup function.
IntegrationRule(int cubature_degree, const panzer::CellData &cell_data)
if side = -1 then we use the cell volume integration rule.
void setup_cv(const panzer::CellData &cell_data, std::string cv_type)
int getPointOffset(const int subcell_index) const
Returns the integration point offset for a given subcell_index (i.e. local face index)
void setup_surface(const Teuchos::RCP< const shards::CellTopology > &cell_topology, const int num_cells, const int num_faces)
Setup a surface integration.
virtual void print(std::ostream &os)
print information about the integration rule
void setup(int cubature_degree, const panzer::CellData &cell_data)
void referenceCoordinates(Kokkos::DynRankView< double, PHX::Device > &container)
Construct an array containing the reference coordinates.
int order() const
Returns the order of integration (cubature degree in intrepid lingo)
void setup(const std::string &ptName, int np, const panzer::CellData &cell_data)