Panzer  Version of the Day
Panzer_IntegrationValues2.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_VALUES2_HPP
45 #define PANZER_INTEGRATION_VALUES2_HPP
46 
47 #include "Teuchos_RCP.hpp"
48 
49 #include "PanzerDiscFE_config.hpp"
51 #include "Panzer_ArrayTraits.hpp"
52 #include "Panzer_Dimension.hpp"
53 #include "Phalanx_MDField.hpp"
54 #include "Intrepid2_Cubature.hpp"
55 
56 namespace panzer {
57 
58  class SubcellConnectivity;
59 
71  template<typename Scalar,
72  typename T0,typename T1,typename T2,typename T3,
73  typename T4,typename T5,typename T6,typename T7>
74  KOKKOS_INLINE_FUNCTION
75  void swapQuadraturePoints(int cell,int a,int b,
76  T0& ref_ip_coordinates,
77  T1& ip_coordinates,
78  T2& weighted_measure,
79  T3& jac,
80  T4& jac_det,
81  T5& jac_inv,
82  T6& surface_normals,
83  T7& surface_rotation_matrices);
84 
85  template <typename Scalar>
87  public:
89 
90  typedef PHX::MDField<Scalar> ArrayDynamic;
91  typedef PHX::MDField<double> DblArrayDynamic;
92 
93  typedef PHX::MDField<Scalar,IP> Array_IP;
94  typedef PHX::MDField<Scalar,IP,Dim> Array_IPDim;
95 
96  typedef PHX::MDField<Scalar,Point> Array_Point;
97  typedef PHX::MDField<Scalar,Cell,IP> Array_CellIP;
98  typedef PHX::MDField<Scalar,Cell,IP,Dim> Array_CellIPDim;
99  typedef PHX::MDField<Scalar,Cell,IP,Dim,Dim> Array_CellIPDimDim;
100 
101  typedef PHX::MDField<Scalar,Cell,BASIS,Dim> Array_CellBASISDim;
102 
103  IntegrationValues2(const std::string & pre="",bool allocArrays=false)
104  : alloc_arrays(allocArrays), prefix(pre), ddims_(1,0) {}
105 
107  void setupArrays(const Teuchos::RCP<const panzer::IntegrationRule>& ir);
108 
109  void setupArraysForNodeRule(const Teuchos::RCP<const panzer::IntegrationRule>& ir);
110 
122  void evaluateValues(const PHX::MDField<Scalar,Cell,NODE,Dim> & vertex_coordinates,
123  const int num_cells = -1,
124  const Teuchos::RCP<const SubcellConnectivity> & face_connectivity = Teuchos::null);
125 
140  void evaluateValues(const PHX::MDField<Scalar,Cell,NODE,Dim> & vertex_coordinates,
141  const PHX::MDField<Scalar,Cell,IP,Dim> & other_ip_coordinates,
142  const int num_cells = -1);
143 
144  Array_IPDim cub_points; // <IP,Dim>
145  Array_IPDim side_cub_points; // <IP,Dim> points on face topology (dim-1)
148  Array_CellIPDimDim jac; // <Cell,IP,Dim,Dim>
149  Array_CellIPDimDim jac_inv; // <Cell,IP,Dim,Dim>
150  Array_CellIP jac_det; // <Cell,IP>
153 
156  // this (appears) is a matrix where the first row is the "normal" direction
157  // and the remaining two rows lie in the hyperplane
158 
159  Teuchos::RCP<const panzer::IntegrationRule> int_rule;
160 
161  Teuchos::RCP<Intrepid2::Cubature<PHX::Device::execution_space,double,double>> intrepid_cubature;
162 
163  // for Shakib stabilization <Cell,IP,Dim,Dim>
167 
168  // integration points
169  Array_CellIPDim ip_coordinates; // <Cell,IP,Dim>
170  Array_CellIPDim ref_ip_coordinates; // <Cell,IP,Dim> for Control Volumes or Surface integrals
171 
174 
175  Array_Point scratch_for_compute_side_measure; // <Point> size: span() == jac.span()
176 
178  void evaluateRemainingValues(const PHX::MDField<Scalar,Cell,NODE,Dim> & in_node_coordinates, const int in_num_cells);
179 
181  void getCubatureCV(const PHX::MDField<Scalar,Cell,NODE,Dim> & in_node_coordinates, const int in_num_cells);
182 
184  void generateSurfaceCubatureValues(const PHX::MDField<Scalar,Cell,NODE,Dim> & in_node_coordinates, const int in_num_cells,const SubcellConnectivity & face_connectivity);
185 
186  protected:
187 
188  // TODO: Make this a utility function that only exists in source file
189  Teuchos::RCP<Intrepid2::Cubature<PHX::Device::execution_space,double,double>> getIntrepidCubature(const panzer::IntegrationRule & ir) const;
190 
191  private:
193  std::string prefix;
194  std::vector<PHX::index_size_type> ddims_;
195 
196  void getCubature(const PHX::MDField<Scalar,Cell,NODE,Dim> & in_node_coordinates, const int in_num_cells);
197  void evaluateValuesCV(const PHX::MDField<Scalar,Cell,NODE,Dim> & vertex_coordinates,const int in_num_cells);
198  };
199 
200 } // namespace panzer
201 
202 #endif
void evaluateValues(const PHX::MDField< Scalar, Cell, NODE, Dim > &vertex_coordinates, const int num_cells=-1, const Teuchos::RCP< const SubcellConnectivity > &face_connectivity=Teuchos::null)
Evaluate basis values.
IntegrationValues2(const std::string &pre="", bool allocArrays=false)
void generateSurfaceCubatureValues(const PHX::MDField< Scalar, Cell, NODE, Dim > &in_node_coordinates, const int in_num_cells, const SubcellConnectivity &face_connectivity)
This should be a private method, but using lambdas on cuda forces this to be public.
PHX::MDField< Scalar, Cell, IP > Array_CellIP
PHX::MDField< Scalar, Point > Array_Point
void evaluateValuesCV(const PHX::MDField< Scalar, Cell, NODE, Dim > &vertex_coordinates, const int in_num_cells)
PHX::MDField< Scalar, IP > Array_IP
void getCubatureCV(const PHX::MDField< Scalar, Cell, NODE, Dim > &in_node_coordinates, const int in_num_cells)
This should be a private method, but using lambdas on cuda forces this to be public.
PHX::MDField< Scalar, Cell, IP, Dim, Dim > Array_CellIPDimDim
KOKKOS_INLINE_FUNCTION void swapQuadraturePoints(int cell, int a, int b, T0 &ref_ip_coordinates, T1 &ip_coordinates, T2 &weighted_measure, T3 &jac, T4 &jac_det, T5 &jac_inv, T6 &surface_normals, T7 &surface_rotation_matrices)
Swap the ordering of quadrature points in a specified cell.
Teuchos::RCP< Intrepid2::Cubature< PHX::Device::execution_space, double, double > > intrepid_cubature
Teuchos::RCP< const panzer::IntegrationRule > int_rule
Teuchos::RCP< Intrepid2::Cubature< PHX::Device::execution_space, double, double > > getIntrepidCubature(const panzer::IntegrationRule &ir) const
PHX::MDField< Scalar, Cell, BASIS, Dim > Array_CellBASISDim
void evaluateRemainingValues(const PHX::MDField< Scalar, Cell, NODE, Dim > &in_node_coordinates, const int in_num_cells)
This should be a private method, but using lambdas on cuda forces this to be public.
ArrayTraits< Scalar, PHX::MDField< Scalar > >::size_type size_type
void setupArrays(const Teuchos::RCP< const panzer::IntegrationRule > &ir)
Sizes/allocates memory for arrays.
void getCubature(const PHX::MDField< Scalar, Cell, NODE, Dim > &in_node_coordinates, const int in_num_cells)
PHX::MDField< Scalar, Cell, IP, Dim > Array_CellIPDim
PHX::MDField< double > DblArrayDynamic
PHX::MDField< Scalar, IP, Dim > Array_IPDim
std::vector< PHX::index_size_type > ddims_
void setupArraysForNodeRule(const Teuchos::RCP< const panzer::IntegrationRule > &ir)