Panzer
Version of the Day
Toggle main menu visibility
Loading...
Searching...
No Matches
disc-fe
src
Panzer_IntegrationRule.hpp
Go to the documentation of this file.
1
// @HEADER
2
// ***********************************************************************
3
//
4
// Panzer: A partial differential equation assembly
5
// engine for strongly coupled complex multiphysics systems
6
// Copyright (2011) Sandia Corporation
7
//
8
// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9
// the U.S. Government retains certain rights in this software.
10
//
11
// Redistribution and use in source and binary forms, with or without
12
// modification, are permitted provided that the following conditions are
13
// met:
14
//
15
// 1. Redistributions of source code must retain the above copyright
16
// notice, this list of conditions and the following disclaimer.
17
//
18
// 2. Redistributions in binary form must reproduce the above copyright
19
// notice, this list of conditions and the following disclaimer in the
20
// documentation and/or other materials provided with the distribution.
21
//
22
// 3. Neither the name of the Corporation nor the names of the
23
// contributors may be used to endorse or promote products derived from
24
// this software without specific prior written permission.
25
//
26
// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27
// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28
// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29
// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30
// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31
// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32
// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33
// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34
// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35
// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36
// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37
//
38
// Questions? Contact Roger P. Pawlowski (rppawlo@sandia.gov) and
39
// Eric C. Cyr (eccyr@sandia.gov)
40
// ***********************************************************************
41
// @HEADER
42
43
44
#ifndef PANZER_INTEGRATION_RULE_HPP
45
#define PANZER_INTEGRATION_RULE_HPP
46
47
#include "Teuchos_RCP.hpp"
48
//#include "Teuchos_ArrayRCP.hpp"
49
#include "Shards_CellTopology.hpp"
50
//#include "Phalanx_DataLayout.hpp"
51
52
#include "
Panzer_PointRule.hpp
"
53
//
54
//#include "Intrepid2_DefaultCubatureFactory.hpp"
55
//#include "Intrepid2_CubatureControlVolume.hpp"
56
//#include "Intrepid2_CubatureControlVolumeSide.hpp"
57
//#include "Intrepid2_CubatureControlVolumeBoundary.hpp"
58
#include "Kokkos_DynRankView.hpp"
59
#include "
Panzer_IntegrationDescriptor.hpp
"
60
61
#include <ostream>
62
#include <string>
63
64
namespace
panzer
{
65
66
class
CellData
;
67
73
class
IntegrationRule
:
public
PointRule
,
public
IntegrationDescriptor
{
74
public
:
75
77
//TEUCHOS_DEPRECATED
78
IntegrationRule
(
int
cubature_degree
,
const
panzer::CellData
& cell_data);
79
80
//TEUCHOS_DEPRECATED
81
IntegrationRule
(
const
panzer::CellData
& cell_data,
const
std::string &
cv_type
);
82
83
IntegrationRule
(
const
panzer::IntegrationDescriptor
& description,
84
const
Teuchos::RCP<const shards::CellTopology> & cell_topology,
85
const
int
num_cells,
86
const
int
num_faces=-1);
87
88
// TODO: Move to protected
89
void
setup
(
int
cubature_degree
,
const
panzer::CellData
& cell_data);
90
91
// TODO: Move to protected
92
void
setup_cv
(
const
panzer::CellData
& cell_data, std::string
cv_type
);
93
95
// Use getOrder() from base class
96
//TEUCHOS_DEPRECATED
97
int
order
()
const
;
98
99
// Use _cubature_order if inside class, use getOrder() if outside
100
//TEUCHOS_DEPRECATED
101
int
cubature_degree
;
102
103
// Use _type if inside class, use getType() if outside class
104
//TEUCHOS_DEPRECATED
105
std::string
cv_type
;
106
108
virtual
void
print
(std::ostream & os);
109
111
void
referenceCoordinates
(Kokkos::DynRankView<double,PHX::Device> & container);
112
114
int
getPointOffset
(
const
int
subcell_index)
const
;
115
116
protected
:
117
119
void
setup_surface
(
const
Teuchos::RCP<const shards::CellTopology> & cell_topology,
120
const
int
num_cells,
121
const
int
num_faces);
122
123
std::vector<int>
_point_offsets
;
124
125
private
:
126
127
};
128
129
}
130
131
#endif
Panzer_IntegrationDescriptor.hpp
Panzer_PointRule.hpp
panzer::CellData
Data for determining cell topology and dimensionality.
Definition
Panzer_CellData.hpp:65
panzer::IntegrationDescriptor
Definition
Panzer_IntegrationDescriptor.hpp:52
panzer::IntegrationDescriptor::IntegrationDescriptor
IntegrationDescriptor()
Constructor for empty integrator.
Definition
Panzer_IntegrationDescriptor.cpp:52
panzer::IntegrationRule::IntegrationRule
IntegrationRule(int cubature_degree, const panzer::CellData &cell_data)
if side = -1 then we use the cell volume integration rule.
Definition
Panzer_IntegrationRule.cpp:58
panzer::IntegrationRule::setup_cv
void setup_cv(const panzer::CellData &cell_data, std::string cv_type)
Definition
Panzer_IntegrationRule.cpp:221
panzer::IntegrationRule::getPointOffset
int getPointOffset(const int subcell_index) const
Returns the integration point offset for a given subcell_index (i.e. local face index).
Definition
Panzer_IntegrationRule.cpp:269
panzer::IntegrationRule::setup_surface
void setup_surface(const Teuchos::RCP< const shards::CellTopology > &cell_topology, const int num_cells, const int num_faces)
Setup a surface integration.
Definition
Panzer_IntegrationRule.cpp:168
panzer::IntegrationRule::print
virtual void print(std::ostream &os)
print information about the integration rule
Definition
Panzer_IntegrationRule.cpp:277
panzer::IntegrationRule::cubature_degree
int cubature_degree
Definition
Panzer_IntegrationRule.hpp:101
panzer::IntegrationRule::setup
void setup(int cubature_degree, const panzer::CellData &cell_data)
Definition
Panzer_IntegrationRule.cpp:130
panzer::IntegrationRule::referenceCoordinates
void referenceCoordinates(Kokkos::DynRankView< double, PHX::Device > &container)
Construct an array containing the reference coordinates.
Definition
Panzer_IntegrationRule.cpp:289
panzer::IntegrationRule::order
int order() const
Returns the order of integration (cubature degree in intrepid lingo).
Definition
Panzer_IntegrationRule.cpp:265
panzer::IntegrationRule::_point_offsets
std::vector< int > _point_offsets
Definition
Panzer_IntegrationRule.hpp:123
panzer::IntegrationRule::cv_type
std::string cv_type
Definition
Panzer_IntegrationRule.hpp:105
panzer::PointRule::PointRule
PointRule(const std::string &ptName, int np, const panzer::CellData &cell_data)
Definition
Panzer_PointRule.cpp:54
panzer
Computes .
Definition
Panzer_BasisValues_Evaluator_decl.hpp:54
Generated by
1.17.0