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 
50 #include "Panzer_BasisIRLayout.hpp"
51 #include "Panzer_ArrayTraits.hpp"
52 
53 namespace panzer {
54 
61  template <typename Scalar>
62  class BasisValues2 {
63  public:
65 
66  typedef PHX::MDField<Scalar> ArrayDynamic;
67  typedef PHX::MDField<Scalar,BASIS,IP,void,void,void,void,void,void> Array_BasisIP;
68  typedef PHX::MDField<Scalar,Cell,BASIS,IP,void,void,void,void,void> Array_CellBasisIP;
69  typedef PHX::MDField<Scalar,BASIS,IP,Dim,void,void,void,void,void> Array_BasisIPDim;
70  typedef PHX::MDField<Scalar,Cell,BASIS,IP,Dim,void,void,void,void> Array_CellBasisIPDim;
71  typedef PHX::MDField<Scalar,BASIS,Dim,void,void,void,void,void,void> Array_BasisDim;
72  typedef PHX::MDField<Scalar,Cell,BASIS,Dim,void,void,void,void,void> Array_CellBasisDim;
73 
74  BasisValues2(const std::string & pre="",bool allocArrays=false,bool buildWeighted=false)
75  : build_weighted(buildWeighted), alloc_arrays(allocArrays), prefix(pre)
76  , ddims_(1,0), references_evaluated(false) {}
77 
80  bool computeDerivatives=true);
81 
82  void evaluateValues(const PHX::MDField<Scalar,IP,Dim,void,void,void,void,void,void> & cub_points,
83  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac,
84  const PHX::MDField<Scalar,Cell,IP,void,void,void,void,void,void> & jac_det,
85  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac_inv);
86 
87  void evaluateValues(const PHX::MDField<Scalar,IP,Dim,void,void,void,void,void,void> & cub_points,
88  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac,
89  const PHX::MDField<Scalar,Cell,IP,void,void,void,void,void,void> & jac_det,
90  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac_inv,
91  const PHX::MDField<Scalar,Cell,IP> & weighted_measure,
92  const PHX::MDField<Scalar,Cell,NODE,Dim> & vertex_coordinates,
93  bool use_vertex_coordinates=true);
94 
95  void evaluateValuesCV(const PHX::MDField<Scalar,Cell,IP,Dim,void,void,void,void,void> & cell_cub_points,
96  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac,
97  const PHX::MDField<Scalar,Cell,IP,void,void,void,void,void,void> & jac_det,
98  const PHX::MDField<Scalar,Cell,IP,Dim,Dim,void,void,void,void> & jac_inv);
99 
101  void applyOrientations(const PHX::MDField<Scalar,Cell,BASIS> & orientations);
102 
103  void setExtendedDimensions(const std::vector<PHX::index_size_type> & ddims)
104  { ddims_ = ddims; }
105 
107 
109  Array_CellBasisIP basis_scalar; // <Cell,BASIS,IP>
110 
111  Array_BasisIPDim basis_ref_vector; // <BASIS,IP,Dim>
112  Array_CellBasisIPDim basis_vector; // <Cell,BASIS,IP,Dim>
113 
114  Array_BasisIPDim grad_basis_ref; // <BASIS,IP,Dim>
115  Array_CellBasisIPDim grad_basis; // <Cell,BASIS,IP,Dim>
116 
118  Array_CellBasisIP curl_basis_scalar; // <Cell,BASIS,IP> - 2D!
119 
120  Array_BasisIPDim curl_basis_ref_vector; // <BASIS,IP,Dim> - 3D!
121  Array_CellBasisIPDim curl_basis_vector; // <Cell,BASIS,IP,Dim> - 3D!
122 
123  Array_BasisIP div_basis_ref; // <BASIS,IP>
124  Array_CellBasisIP div_basis; // <Cell,BASIS,IP>
125 
132 
139 
146 
148 
150 
154  std::string prefix;
155  std::vector<PHX::index_size_type> ddims_;
156 
157  private:
162  void evaluateReferenceValues(const PHX::MDField<Scalar,IP,Dim> & cub_points,bool compute_derivatives,bool use_vertex_coordinates);
163 
165  };
166 
167 } // namespace panzer
168 
169 #endif
Array_CellBasisIPDim weighted_curl_basis_vector
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
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
Array_CellBasisIP div_basis
void setExtendedDimensions(const std::vector< PHX::index_size_type > &ddims)
Array_CellBasisIPDim basis_vector
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
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.
Teuchos::RCP< Intrepid2::Basis< Scalar, ArrayDynamic > > intrepid_basis
PHX::MDField< Scalar, BASIS, IP, Dim, void, void, void, void, void > Array_BasisIPDim
Teuchos::RCP< const panzer::PureBasis > basis
Interpolates basis DOF values to IP DOF values.
Array_CellBasisIPDim grad_basis
Array_CellBasisIP weighted_div_basis
PHX::MDField< Scalar, Cell, BASIS, IP, void, void, void, void, void > Array_CellBasisIP
void applyOrientations(const PHX::MDField< Scalar, Cell, BASIS > &orientations)
Method to apply orientations to a basis values container.
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
Array_BasisIPDim grad_basis_ref
Array_CellBasisIP curl_basis_scalar
PHX::MDField< Scalar > ArrayDynamic