Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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  template <typename Scalar>
60  public:
62 
63  typedef PHX::MDField<Scalar> ArrayDynamic;
64  typedef PHX::MDField<double> DblArrayDynamic;
65 
66  typedef PHX::MDField<Scalar,IP> Array_IP;
67  typedef PHX::MDField<Scalar,IP,Dim> Array_IPDim;
68 
69  typedef PHX::MDField<Scalar,Point> Array_Point;
70  typedef PHX::MDField<Scalar,Cell,IP> Array_CellIP;
71  typedef PHX::MDField<Scalar,Cell,IP,Dim> Array_CellIPDim;
72  typedef PHX::MDField<Scalar,Cell,IP,Dim,Dim> Array_CellIPDimDim;
73 
74  typedef PHX::MDField<Scalar,Cell,BASIS,Dim> Array_CellBASISDim;
75 
76  IntegrationValues2(const std::string & pre="",bool allocArrays=false)
77  : alloc_arrays(allocArrays), prefix(pre), ddims_(1,0) {}
78 
81 
83 
93  void evaluateValues(const PHX::MDField<Scalar,Cell,NODE,Dim> & vertex_coordinates,
94  const int num_cells = -1);
95 
110  void evaluateValues(const PHX::MDField<Scalar,Cell,NODE,Dim> & vertex_coordinates,
111  const PHX::MDField<Scalar,Cell,IP,Dim> & other_ip_coordinates,
112  const int num_cells = -1);
113 
114  Array_IPDim cub_points; // <IP,Dim>
115  Array_IPDim side_cub_points; // <IP,Dim> points on face topology (dim-1)
118  Array_CellIPDimDim jac; // <Cell,IP,Dim,Dim>
119  Array_CellIPDimDim jac_inv; // <Cell,IP,Dim,Dim>
120  Array_CellIP jac_det; // <Cell,IP>
123 
126  // this (appears) is a matrix where the first row is the "normal" direction
127  // and the remaining two rows lie in the hyperplane
128 
130 
132 
133  // for Shakib stabilization <Cell,IP,Dim,Dim>
137 
138  // integration points
139  Array_CellIPDim ip_coordinates; // <Cell,IP,Dim>
140  Array_CellIPDim ref_ip_coordinates; // <Cell,IP,Dim> for Control Volumes or Surface integrals
141 
144 
145  Array_Point scratch_for_compute_side_measure; // <Point> size: span() == jac.span()
146 
147 
160  static void uniqueCoordOrdering(Array_CellIPDim & coords,
161  int cell,
162  int offset,
163  std::vector<int> & order);
164 
172  void swapQuadraturePoints(int cell,int a,int b);
173 
174  static void convertNormalToRotationMatrix(const Scalar normal[3], Scalar transverse[3], Scalar binormal[3]);
175 
176  protected:
177 
178 
179 
180  // TODO: Make this a utility function that only exists in source file
182 
183 
184  private:
186  std::string prefix;
187  std::vector<PHX::index_size_type> ddims_;
188 
189  void generateSurfaceCubatureValues(const PHX::MDField<Scalar,Cell,NODE,Dim> & in_node_coordinates, const int in_num_cells);
190  void getCubature(const PHX::MDField<Scalar,Cell,NODE,Dim> & in_node_coordinates, const int in_num_cells);
191  void getCubatureCV(const PHX::MDField<Scalar,Cell,NODE,Dim> & in_node_coordinates, const int in_num_cells);
192  void evaluateRemainingValues(const PHX::MDField<Scalar,Cell,NODE,Dim> & in_node_coordinates, const int in_num_cells);
193  void evaluateValuesCV(const PHX::MDField<Scalar,Cell,NODE,Dim> & vertex_coordinates,const int in_num_cells);
194  };
195 
196 } // namespace panzer
197 
198 #endif
static void convertNormalToRotationMatrix(const Scalar normal[3], Scalar transverse[3], Scalar binormal[3])
static void uniqueCoordOrdering(Array_CellIPDim &coords, int cell, int offset, std::vector< int > &order)
Using coordinate build an arrray that specifies a unique ordering.
IntegrationValues2(const std::string &pre="", bool allocArrays=false)
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 generateSurfaceCubatureValues(const PHX::MDField< Scalar, Cell, NODE, Dim > &in_node_coordinates, const int in_num_cells)
void swapQuadraturePoints(int cell, int a, int b)
Swap the ordering of quadrature points in a specified cell.
void getCubatureCV(const PHX::MDField< Scalar, Cell, NODE, Dim > &in_node_coordinates, const int in_num_cells)
PHX::MDField< Scalar, Cell, IP, Dim, Dim > Array_CellIPDimDim
Teuchos::RCP< Intrepid2::Cubature< PHX::Device::execution_space, double, double > > intrepid_cubature
Teuchos::RCP< const panzer::IntegrationRule > int_rule
void evaluateValues(const PHX::MDField< Scalar, Cell, NODE, Dim > &vertex_coordinates, const int num_cells=-1)
Evaluate basis values.
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)
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
Teuchos::RCP< Intrepid2::Cubature< PHX::Device::execution_space, double, double > > getIntrepidCubature(const panzer::IntegrationRule &ir) const
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)