Panzer  Version of the Day
Panzer_BasisValues2.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 #ifndef PANZER_BASIS_VALUES2_HPP
44 #define PANZER_BASIS_VALUES2_HPP
45 
46 #include "Teuchos_RCP.hpp"
47 
48 #include "Intrepid2_Basis.hpp"
49 #include "Intrepid2_Orientation.hpp"
50 
51 #include "Panzer_BasisIRLayout.hpp"
52 #include "Panzer_ArrayTraits.hpp"
53 
54 namespace panzer {
55 
62  template <typename Scalar>
63  class BasisValues2 {
64  public:
66 
67  typedef PHX::MDField<Scalar> ArrayDynamic;
68  typedef PHX::MDField<Scalar,BASIS,IP,void,void,void,void,void,void> Array_BasisIP;
69  typedef PHX::MDField<Scalar,Cell,BASIS,IP,void,void,void,void,void> Array_CellBasisIP;
70  typedef PHX::MDField<Scalar,BASIS,IP,Dim,void,void,void,void,void> Array_BasisIPDim;
71  typedef PHX::MDField<Scalar,Cell,BASIS,IP,Dim,void,void,void,void> Array_CellBasisIPDim;
72  typedef PHX::MDField<Scalar,BASIS,Dim,void,void,void,void,void,void> Array_BasisDim;
73  typedef PHX::MDField<Scalar,Cell,BASIS,Dim,void,void,void,void,void> Array_CellBasisDim;
74 
75  BasisValues2(const std::string & pre="",bool allocArrays=false,bool buildWeighted=false)
76  : build_weighted(buildWeighted), alloc_arrays(allocArrays), prefix(pre)
77  , ddims_(1,0), references_evaluated(false) {}
78 
80  void setupArrays(const Teuchos::RCP<const panzer::BasisIRLayout>& basis,
81  bool computeDerivatives=true);
82 
83  void evaluateValues(const PHX::MDField<Scalar,IP,Dim,void,void,void,void,void,void> & cub_points,
84  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac,
85  const PHX::MDField<Scalar,Cell,IP,void,void,void,void,void,void> & jac_det,
86  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac_inv);
87 
88  void evaluateValues(const PHX::MDField<Scalar,IP,Dim,void,void,void,void,void,void> & cub_points,
89  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac,
90  const PHX::MDField<Scalar,Cell,IP,void,void,void,void,void,void> & jac_det,
91  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac_inv,
92  const PHX::MDField<Scalar,Cell,IP> & weighted_measure,
93  const PHX::MDField<Scalar,Cell,NODE,Dim> & vertex_coordinates,
94  bool use_vertex_coordinates=true);
95 
96  void evaluateValuesCV(const PHX::MDField<Scalar,Cell,IP,Dim,void,void,void,void,void> & cell_cub_points,
97  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac,
98  const PHX::MDField<Scalar,Cell,IP,void,void,void,void,void,void> & jac_det,
99  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac_inv);
100 
101  void evaluateValues(const PHX::MDField<Scalar,Cell,IP,Dim,void,void,void,void,void> & cub_points,
102  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac,
103  const PHX::MDField<Scalar,Cell,IP,void,void,void,void,void,void> & jac_det,
104  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac_inv,
105  const PHX::MDField<Scalar,Cell,IP> & weighted_measure,
106  const PHX::MDField<Scalar,Cell,NODE,Dim> & vertex_coordinates,
107  bool use_vertex_coordinates=true);
108 
109 
111  // some evaluators use this apply orientation (this will be deprecated)
112  void applyOrientations(const PHX::MDField<const Scalar,Cell,BASIS> & orientations);
113 
114  // this is used in workset factory
115  void applyOrientations(const std::vector<Intrepid2::Orientation> & orientations);
116 
117  void setExtendedDimensions(const std::vector<PHX::index_size_type> & ddims)
118  { ddims_ = ddims; }
119 
121 
123  Array_CellBasisIP basis_scalar; // <Cell,BASIS,IP>
124 
125  Array_BasisIPDim basis_ref_vector; // <BASIS,IP,Dim>
126  Array_CellBasisIPDim basis_vector; // <Cell,BASIS,IP,Dim>
127 
128  Array_BasisIPDim grad_basis_ref; // <BASIS,IP,Dim>
129  Array_CellBasisIPDim grad_basis; // <Cell,BASIS,IP,Dim>
130 
132  Array_CellBasisIP curl_basis_scalar; // <Cell,BASIS,IP> - 2D!
133 
134  Array_BasisIPDim curl_basis_ref_vector; // <BASIS,IP,Dim> - 3D!
135  Array_CellBasisIPDim curl_basis_vector; // <Cell,BASIS,IP,Dim> - 3D!
136 
137  Array_BasisIP div_basis_ref; // <BASIS,IP>
138  Array_CellBasisIP div_basis; // <Cell,BASIS,IP>
139 
146 
153 
160 
161  Teuchos::RCP<const panzer::BasisIRLayout> basis_layout;
162 
163  Teuchos::RCP<Intrepid2::Basis<PHX::Device::execution_space,Scalar,Scalar> > intrepid_basis;
164 
168  std::string prefix;
169  std::vector<PHX::index_size_type> ddims_;
170 
171  protected:
172 
173 
174  void evaluateValues_Const(const PHX::MDField<Scalar,Cell,IP,Dim,void,void,void,void,void> & cub_points,
175  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac_inv,
176  const PHX::MDField<Scalar,Cell,IP> & weighted_measure);
177 
178  void evaluateValues_HGrad(const PHX::MDField<Scalar,Cell,IP,Dim,void,void,void,void,void> & cub_points,
179  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac_inv,
180  const PHX::MDField<Scalar,Cell,IP> & weighted_measure);
181 
182  void evaluateValues_HCurl(const PHX::MDField<Scalar,Cell,IP,Dim,void,void,void,void,void> & cub_points,
183  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac,
184  const PHX::MDField<Scalar,Cell,IP,void,void,void,void,void,void> & jac_det,
185  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac_inv,
186  const PHX::MDField<Scalar,Cell,IP> & weighted_measure);
187 
188  void evaluateValues_HDiv(const PHX::MDField<Scalar,Cell,IP,Dim,void,void,void,void,void> & cub_points,
189  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac,
190  const PHX::MDField<Scalar,Cell,IP,void,void,void,void,void,void> & jac_det,
191  const PHX::MDField<Scalar,Cell,IP> & weighted_measure);
192 
193  private:
194 
195  void evaluateBasisCoordinates(const PHX::MDField<Scalar,Cell,NODE,Dim> & vertex_coordinates);
196 
201  void evaluateReferenceValues(const PHX::MDField<Scalar,IP,Dim> & cub_points,bool compute_derivatives,bool use_vertex_coordinates);
202 
204  };
205 
206 } // namespace panzer
207 
208 #endif
Array_CellBasisIPDim weighted_curl_basis_vector
void evaluateValues_HGrad(const PHX::MDField< Scalar, Cell, IP, Dim, void, void, void, void, void > &cub_points, const PHX::MDField< Scalar, Cell, IP, Dim, Dim, void, void, void, void > &jac_inv, const PHX::MDField< Scalar, Cell, IP > &weighted_measure)
PHX::MDField< Scalar, Cell, BASIS, IP, Dim, void, void, void, void > Array_CellBasisIPDim
Array_BasisIPDim curl_basis_ref_vector
PHX::MDField< Scalar, Cell, BASIS, Dim, void, void, void, void, void > Array_CellBasisDim
void evaluateValues(const PHX::MDField< Scalar, IP, Dim, void, void, void, void, void, void > &cub_points, const PHX::MDField< Scalar, Cell, IP, Dim, Dim, void, void, void, void > &jac, const PHX::MDField< Scalar, Cell, IP, void, void, void, void, void, void > &jac_det, const PHX::MDField< Scalar, Cell, IP, Dim, Dim, void, void, void, void > &jac_inv)
void evaluateReferenceValues(const PHX::MDField< Scalar, IP, Dim > &cub_points, bool compute_derivatives, bool use_vertex_coordinates)
ArrayTraits< Scalar, PHX::MDField< Scalar, void, void, void, void, void, void, void, void > >::size_type size_type
void evaluateBasisCoordinates(const PHX::MDField< Scalar, Cell, NODE, Dim > &vertex_coordinates)
Array_CellBasisIPDim curl_basis_vector
BasisValues2(const std::string &pre="", bool allocArrays=false, bool buildWeighted=false)
Array_CellBasisIP weighted_basis_scalar
Array_CellBasisIP basis_scalar
void evaluateValues_Const(const PHX::MDField< Scalar, Cell, IP, Dim, void, void, void, void, void > &cub_points, const PHX::MDField< Scalar, Cell, IP, Dim, Dim, void, void, void, void > &jac_inv, const PHX::MDField< Scalar, Cell, IP > &weighted_measure)
Array_CellBasisIP div_basis
void applyOrientations(const PHX::MDField< const Scalar, Cell, BASIS > &orientations)
Method to apply orientations to a basis values container.
void setExtendedDimensions(const std::vector< PHX::index_size_type > &ddims)
Array_CellBasisIPDim basis_vector
void evaluateValues_HCurl(const PHX::MDField< Scalar, Cell, IP, Dim, void, void, void, void, void > &cub_points, const PHX::MDField< Scalar, Cell, IP, Dim, Dim, void, void, void, void > &jac, const PHX::MDField< Scalar, Cell, IP, void, void, void, void, void, void > &jac_det, const PHX::MDField< Scalar, Cell, IP, Dim, Dim, void, void, void, void > &jac_inv, const PHX::MDField< Scalar, Cell, IP > &weighted_measure)
PHX::MDField< Scalar, BASIS, IP, void, void, void, void, void, void > Array_BasisIP
Array_BasisIPDim basis_ref_vector
std::vector< PHX::index_size_type > ddims_
PHX::MDField< Scalar, BASIS, Dim, void, void, void, void, void, void > Array_BasisDim
Teuchos::RCP< const panzer::BasisIRLayout > basis_layout
Teuchos::RCP< Intrepid2::Basis< PHX::Device::execution_space, Scalar, Scalar > > intrepid_basis
Array_CellBasisIP weighted_curl_basis_scalar
Array_CellBasisDim basis_coordinates
Array_CellBasisIPDim weighted_grad_basis
Array_CellBasisIPDim weighted_basis_vector
void setupArrays(const Teuchos::RCP< const panzer::BasisIRLayout > &basis, bool computeDerivatives=true)
Sizes/allocates memory for arrays.
PHX::MDField< Scalar, BASIS, IP, Dim, void, void, void, void, void > Array_BasisIPDim
Array_CellBasisIPDim grad_basis
Array_CellBasisIP weighted_div_basis
PHX::MDField< Scalar, Cell, BASIS, IP, void, void, void, void, void > Array_CellBasisIP
Array_BasisDim basis_coordinates_ref
void evaluateValuesCV(const PHX::MDField< Scalar, Cell, IP, Dim, void, void, void, void, void > &cell_cub_points, const PHX::MDField< Scalar, Cell, IP, Dim, Dim, void, void, void, void > &jac, const PHX::MDField< Scalar, Cell, IP, void, void, void, void, void, void > &jac_det, const PHX::MDField< Scalar, Cell, IP, Dim, Dim, void, void, void, void > &jac_inv)
Array_BasisIP curl_basis_ref_scalar
PureBasis::EElementSpace getElementSpace() const
void evaluateValues_HDiv(const PHX::MDField< Scalar, Cell, IP, Dim, void, void, void, void, void > &cub_points, const PHX::MDField< Scalar, Cell, IP, Dim, Dim, void, void, void, void > &jac, const PHX::MDField< Scalar, Cell, IP, void, void, void, void, void, void > &jac_det, const PHX::MDField< Scalar, Cell, IP > &weighted_measure)
Array_BasisIPDim grad_basis_ref
Array_CellBasisIP curl_basis_scalar
PHX::MDField< Scalar > ArrayDynamic