Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_BasisValues2.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Panzer: A partial differential equation assembly
4 // engine for strongly coupled complex multiphysics systems
5 //
6 // Copyright 2011 NTESS and the Panzer contributors.
7 // SPDX-License-Identifier: BSD-3-Clause
8 // *****************************************************************************
9 // @HEADER
10 
11 #ifndef PANZER_BASIS_VALUES2_HPP
12 #define PANZER_BASIS_VALUES2_HPP
13 
14 #include "Teuchos_RCP.hpp"
15 
16 #include "Intrepid2_Basis.hpp"
17 #include "Intrepid2_Orientation.hpp"
18 
19 #include "Panzer_BasisIRLayout.hpp"
20 #include "Panzer_ArrayTraits.hpp"
21 
22 namespace panzer {
23 
24  // For applying orientations with lazy evaluation
25  class OrientationsInterface;
26 
67  template <typename Scalar>
68  class BasisValues2 {
69  public:
71  using IntrepidBasis = Intrepid2::Basis<PHX::Device::execution_space,Scalar,Scalar>;
72 
80 
88 
89  //=======================================================================================================
90  // DEPRECATED Interface
91 
101  BasisValues2(const std::string & prefix="",
102  const bool allocArrays=false,
103  const bool buildWeighted=false);
104 
107  bool computeDerivatives=true);
108 
109  void evaluateValues(const PHX::MDField<Scalar,IP,Dim> & cub_points,
111  const PHX::MDField<Scalar,Cell,IP> & jac_det,
112  const PHX::MDField<Scalar,Cell,IP,Dim,Dim> & jac_inv,
113  const int in_num_cells = -1);
114 
115  void evaluateValues(const PHX::MDField<Scalar,IP,Dim> & cub_points,
117  const PHX::MDField<Scalar,Cell,IP> & jac_det,
118  const PHX::MDField<Scalar,Cell,IP,Dim,Dim> & jac_inv,
119  const PHX::MDField<Scalar,Cell,IP> & weighted_measure,
120  const PHX::MDField<Scalar,Cell,NODE,Dim> & node_coordinates,
121  bool use_node_coordinates=true,
122  const int in_num_cells = -1);
123 
124  void evaluateValuesCV(const PHX::MDField<Scalar,Cell,IP,Dim> & cell_cub_points,
126  const PHX::MDField<Scalar,Cell,IP> & jac_det,
127  const PHX::MDField<Scalar,Cell,IP,Dim,Dim> & jac_inv);
128 
129 
130  void evaluateValuesCV(const PHX::MDField<Scalar,Cell,IP,Dim> & cell_cub_points,
132  const PHX::MDField<Scalar,Cell,IP> & jac_det,
133  const PHX::MDField<Scalar,Cell,IP,Dim,Dim> & jac_inv,
134  const PHX::MDField<Scalar,Cell,NODE,Dim> & node_coordinates,
135  bool use_node_coordinates=true,
136  const int in_num_cells = -1);
137 
138 
139  void evaluateValues(const PHX::MDField<Scalar,Cell,IP,Dim> & cub_points,
141  const PHX::MDField<Scalar,Cell,IP> & jac_det,
142  const PHX::MDField<Scalar,Cell,IP,Dim,Dim> & jac_inv,
143  const PHX::MDField<Scalar,Cell,IP> & weighted_measure,
144  const PHX::MDField<Scalar,Cell,NODE,Dim> & node_coordinates,
145  bool use_node_coordinates=true,
146  const int in_num_cells = -1);
147 
148 
150  // some evaluators use this apply orientation (this will be deprecated)
151  void applyOrientations(const PHX::MDField<const Scalar,Cell,BASIS> & orientations);
152 
153  // this is used in workset factory
154  void applyOrientations(const std::vector<Intrepid2::Orientation> & orientations,
155  const int in_num_cells = -1);
156 
157  void setExtendedDimensions(const std::vector<PHX::index_size_type> & ddims)
158  { ddims_ = ddims; }
159 
161 
162  mutable Array_BasisIP basis_ref_scalar; // <BASIS,IP>
163  mutable Array_CellBasisIP basis_scalar; // <Cell,BASIS,IP>
164 
165  mutable Array_BasisIPDim basis_ref_vector; // <BASIS,IP,Dim>
166  mutable Array_CellBasisIPDim basis_vector; // <Cell,BASIS,IP,Dim>
167 
168  mutable Array_BasisIPDim grad_basis_ref; // <BASIS,IP,Dim>
169  mutable Array_CellBasisIPDim grad_basis; // <Cell,BASIS,IP,Dim>
170 
171  mutable Array_BasisIP curl_basis_ref_scalar; // <BASIS,IP> - 2D!
172  mutable Array_CellBasisIP curl_basis_scalar; // <Cell,BASIS,IP> - 2D!
173 
174  mutable Array_BasisIPDim curl_basis_ref_vector; // <BASIS,IP,Dim> - 3D!
175  mutable Array_CellBasisIPDim curl_basis_vector; // <Cell,BASIS,IP,Dim> - 3D!
176 
177  mutable Array_BasisIP div_basis_ref; // <BASIS,IP>
178  mutable Array_CellBasisIP div_basis; // <Cell,BASIS,IP>
179 
180  mutable Array_CellBasisIP weighted_basis_scalar; // <Cell,BASIS,IP>
181  mutable Array_CellBasisIPDim weighted_basis_vector; // <Cell,BASIS,IP,Dim>
182  mutable Array_CellBasisIPDim weighted_grad_basis; // <Cell,BASIS,IP,Dim>
183  mutable Array_CellBasisIP weighted_curl_basis_scalar; // <Cell,BASIS,IP>
184  mutable Array_CellBasisIPDim weighted_curl_basis_vector; // <Cell,BASIS,IP,Dim>
185  mutable Array_CellBasisIP weighted_div_basis; // <Cell,BASIS,IP>
186 
192  mutable Array_BasisDim basis_coordinates_ref; // <BASIS,Dim>
193 
199  mutable Array_CellBasisDim basis_coordinates; // <Cell,BASIS,Dim>
200 
202 
204 
208  std::string prefix;
209  std::vector<PHX::index_size_type> ddims_;
210 
211  bool orientationsApplied() const
212  {return orientations_applied_;}
213 
214  void evaluateBasisCoordinates(const PHX::MDField<Scalar,Cell,NODE,Dim> & node_coordinates,
215  const int in_num_cells = -1);
216 
217  private:
218 
220 
221  // Have orientations been applied
223 
224  //=======================================================================================================
225  // New Interface
226 
227  protected:
228 
229  // Reset internal data structures for lazy evaluation
230  void resetArrays();
231 
232  // Number of total cells (owned + ghost + virtual)
234 
235  // Number of cells to evaluate (HACK for backward compatability)
237 
238  // True if a uniform point space is used
240 
241  // Only valid for uniform reference space
243 
244  // For non-uniform reference space
246 
247  // Geometry objects that represent cubature/point values
252 
253  // Coordinates of the mesh nodes
255 
256  // Cell topology from the mesh
258 
259  // Number of cells to apply orientations to (required in situations where virtual cells exist)
261 
262  // Orientations object
263  std::vector<Intrepid2::Orientation> orientations_;
264 
271  mutable bool grad_basis_evaluated_;
277  mutable bool div_basis_evaluated_;
286 
287  public:
288 
302  void
306  PHX::MDField<const Scalar, Cell, IP> point_jacobian_determinant,
307  PHX::MDField<const Scalar, Cell, IP, Dim, Dim> point_jacobian_inverse,
308  const int num_evaluated_cells = -1);
309 
323  void
325  PHX::MDField<const Scalar, IP, Dim> reference_points,
327  PHX::MDField<const Scalar, Cell, IP> point_jacobian_determinant,
328  PHX::MDField<const Scalar, Cell, IP, Dim, Dim> point_jacobian_inverse,
329  const int num_evaluated_cells = -1);
330 
332  void
333  setOrientations(const std::vector<Intrepid2::Orientation> & orientations,
334  const int num_orientations_cells = -1);
335 
337  void
339 
342  void
345 
347  void
349 
351  bool
353  {return is_uniform_;}
354 
357 
359  const std::vector<PHX::index_size_type> &
361  {return ddims_;}
362 
363  //=======================================================================================================
364  // Reference space evaluations
365 
378  getBasisCoordinatesRef(const bool cache = true,
379  const bool force = false) const;
380 
393  getBasisValuesRef(const bool cache = true,
394  const bool force = false) const;
395 
408  getVectorBasisValuesRef(const bool cache = true,
409  const bool force = false) const;
410 
423  getGradBasisValuesRef(const bool cache = true,
424  const bool force = false) const;
425 
439  getCurl2DVectorBasisRef(const bool cache = true,
440  const bool force = false) const;
441 
455  getCurlVectorBasisRef(const bool cache = true,
456  const bool force = false) const;
457 
471  getDivVectorBasisRef(const bool cache = true,
472  const bool force = false) const;
473 
474  //=======================================================================================================
475  // Mesh evaluations
476 
485  getBasisCoordinates(const bool cache = true,
486  const bool force = false) const;
487 
500  getBasisValues(const bool weighted,
501  const bool cache = true,
502  const bool force = false) const;
503 
516  getVectorBasisValues(const bool weighted,
517  const bool cache = true,
518  const bool force = false) const;
519 
532  getGradBasisValues(const bool weighted,
533  const bool cache = true,
534  const bool force = false) const;
535 
549  getCurl2DVectorBasis(const bool weighted,
550  const bool cache = true,
551  const bool force = false) const;
552 
566  getCurlVectorBasis(const bool weighted,
567  const bool cache = true,
568  const bool force = false) const;
569 
583  getDivVectorBasis(const bool weighted,
584  const bool cache = true,
585  const bool force = false) const;
586 
587  //=======================================================================================================
588 
589  };
590 
591 } // namespace panzer
592 
593 #endif
Array_CellBasisIPDim weighted_curl_basis_vector
void setOrientations(const std::vector< Intrepid2::Orientation > &orientations, const int num_orientations_cells=-1)
Set the orientations object for applying orientations using the lazy evaluation path - required for c...
Teuchos::RCP< IntrepidBasis > intrepid_basis
ConstArray_BasisIP getCurl2DVectorBasisRef(const bool cache=true, const bool force=false) const
Get the curl of a vector basis evaluated at reference points.
void setCellVertexCoordinates(PHX::MDField< Scalar, Cell, NODE, Dim > vertex_coordinates)
Array_BasisIPDim curl_basis_ref_vector
void evaluateValuesCV(const PHX::MDField< Scalar, Cell, IP, Dim > &cell_cub_points, const PHX::MDField< Scalar, Cell, IP, Dim, Dim > &jac, const PHX::MDField< Scalar, Cell, IP > &jac_det, const PHX::MDField< Scalar, Cell, IP, Dim, Dim > &jac_inv)
PHX::MDField< const Scalar, Cell, BASIS, IP > ConstArray_CellBasisIP
bool basis_ref_scalar_evaluated_
Used to check if arrays have been cached.
ConstArray_BasisIPDim getVectorBasisValuesRef(const bool cache=true, const bool force=false) const
Get the vector basis values evaluated at reference points.
bool hasUniformReferenceSpace() const
Check if reference point space is uniform across all cells (faster evaluation)
PHX::MDField< const Scalar, Cell, IP > cubature_weights_
Array_CellBasisIPDim curl_basis_vector
PHX::MDField< const Scalar, IP, Dim > cubature_points_uniform_ref_
ConstArray_BasisIP getBasisValuesRef(const bool cache=true, const bool force=false) const
Get the basis values evaluated at reference points.
panzer::BasisDescriptor getBasisDescriptor() const
Return the basis descriptor.
ConstArray_CellBasisIPDim getCurlVectorBasis(const bool weighted, const bool cache=true, const bool force=false) const
Get the curl of a 3D vector basis evaluated at mesh points.
PHX::MDField< const Scalar, Cell, IP, Dim, Dim > cubature_jacobian_inverse_
Array_CellBasisIP weighted_basis_scalar
PHX::MDField< const Scalar, Cell, IP, Dim > cubature_points_ref_
PHX::MDField< const Scalar > ConstArrayDynamic
PHX::MDField< const Scalar, Cell, BASIS, Dim > ConstArray_CellBasisDim
Teuchos::RCP< const shards::CellTopology > cell_topology_
Array_CellBasisIP basis_scalar
std::vector< Intrepid2::Orientation > orientations_
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)
PHX::MDField< Scalar, Cell, BASIS, IP > Array_CellBasisIP
PHX::MDField< Scalar, Cell, BASIS, IP, Dim > Array_CellBasisIPDim
Array_CellBasisIPDim basis_vector
PHX::MDField< const Scalar, BASIS, IP > ConstArray_BasisIP
PHX::MDField< Scalar, Cell, BASIS, Dim > Array_CellBasisDim
PHX::MDField< const Scalar, Cell, IP, Dim, Dim > cubature_jacobian_
const std::vector< PHX::index_size_type > & getExtendedDimensions() const
Get the extended dimensions used by sacado AD allocations.
PHX::MDField< const Scalar, Cell, IP > cubature_jacobian_determinant_
Array_BasisIPDim basis_ref_vector
ConstArray_CellBasisIP getBasisValues(const bool weighted, const bool cache=true, const bool force=false) const
Get the basis values evaluated at mesh points.
std::vector< PHX::index_size_type > ddims_
ArrayTraits< Scalar, PHX::MDField< Scalar > >::size_type size_type
void setup(const Teuchos::RCP< const panzer::BasisIRLayout > &basis, PHX::MDField< const Scalar, Cell, IP, Dim > reference_points, PHX::MDField< const Scalar, Cell, IP, Dim, Dim > point_jacobian, PHX::MDField< const Scalar, Cell, IP > point_jacobian_determinant, PHX::MDField< const Scalar, Cell, IP, Dim, Dim > point_jacobian_inverse, const int num_evaluated_cells=-1)
Setup for lazy evaluation for non-uniform point layout.
ConstArray_CellBasisDim getBasisCoordinates(const bool cache=true, const bool force=false) const
Carterisan coordinates for basis coefficients in mesh space.
Teuchos::RCP< const panzer::BasisIRLayout > basis_layout
PureBasis::EElementSpace getElementSpace() const
ConstArray_BasisIP getDivVectorBasisRef(const bool cache=true, const bool force=false) const
Get the divergence of a vector basis evaluated at reference points.
Array_CellBasisIP weighted_curl_basis_scalar
Array_CellBasisDim basis_coordinates
Array_CellBasisIPDim weighted_grad_basis
Array_CellBasisIPDim weighted_basis_vector
ConstArray_CellBasisIP getDivVectorBasis(const bool weighted, const bool cache=true, const bool force=false) const
Get the divergence of a vector basis evaluated at mesh points.
void setupArrays(const Teuchos::RCP< const panzer::BasisIRLayout > &basis, bool computeDerivatives=true)
Sizes/allocates memory for arrays.
PHX::MDField< const Scalar, Cell, BASIS, IP, Dim > ConstArray_CellBasisIPDim
Intrepid2::Basis< PHX::Device::execution_space, Scalar, Scalar > IntrepidBasis
PHX::MDField< const Scalar, BASIS, IP, Dim > ConstArray_BasisIPDim
ConstArray_CellBasisIPDim getVectorBasisValues(const bool weighted, const bool cache=true, const bool force=false) const
Get the vector basis values evaluated at mesh points.
ConstArray_BasisIPDim getCurlVectorBasisRef(const bool cache=true, const bool force=false) const
Get the curl of a vector basis evaluated at reference points.
BasisValues2(const std::string &prefix="", const bool allocArrays=false, const bool buildWeighted=false)
Main constructor.
ConstArray_BasisDim getBasisCoordinatesRef(const bool cache=true, const bool force=false) const
Get the reference coordinates for basis.
void evaluateBasisCoordinates(const PHX::MDField< Scalar, Cell, NODE, Dim > &node_coordinates, const int in_num_cells=-1)
Array_CellBasisIPDim grad_basis
PHX::MDField< Scalar, BASIS, IP, Dim > Array_BasisIPDim
Array_CellBasisIP weighted_div_basis
PHX::MDField< Scalar, BASIS, Dim > Array_BasisDim
Array_BasisDim basis_coordinates_ref
Array_BasisIP curl_basis_ref_scalar
void setupUniform(const Teuchos::RCP< const panzer::BasisIRLayout > &basis, PHX::MDField< const Scalar, IP, Dim > reference_points, PHX::MDField< const Scalar, Cell, IP, Dim, Dim > point_jacobian, PHX::MDField< const Scalar, Cell, IP > point_jacobian_determinant, PHX::MDField< const Scalar, Cell, IP, Dim, Dim > point_jacobian_inverse, const int num_evaluated_cells=-1)
Setup for lazy evaluation for uniform point layout.
PHX::MDField< const Scalar, BASIS, Dim > ConstArray_BasisDim
PHX::MDField< const Scalar, Cell, NODE, Dim > cell_node_coordinates_
void setCellNodeCoordinates(PHX::MDField< Scalar, Cell, NODE, Dim > node_coordinates)
Set the cell node coordinates (required for getBasisCoordinates())
void evaluateValues(const PHX::MDField< Scalar, IP, Dim > &cub_points, const PHX::MDField< Scalar, Cell, IP, Dim, Dim > &jac, const PHX::MDField< Scalar, Cell, IP > &jac_det, const PHX::MDField< Scalar, Cell, IP, Dim, Dim > &jac_inv, const int in_num_cells=-1)
ConstArray_CellBasisIPDim getGradBasisValues(const bool weighted, const bool cache=true, const bool force=false) const
Get the gradient of the basis evaluated at mesh points.
void setWeightedMeasure(PHX::MDField< const Scalar, Cell, IP > weighted_measure)
Set the cubature weights (weighted measure) for the basis values object - required to get weighted ba...
PHX::MDField< Scalar, BASIS, IP > Array_BasisIP
ConstArray_CellBasisIP getCurl2DVectorBasis(const bool weighted, const bool cache=true, const bool force=false) const
Get the curl of a 2D vector basis evaluated at mesh points.
Array_BasisIPDim grad_basis_ref
ConstArray_BasisIPDim getGradBasisValuesRef(const bool cache=true, const bool force=false) const
Get the gradient of the basis evaluated at reference points.
Array_CellBasisIP curl_basis_scalar
PHX::MDField< Scalar > ArrayDynamic