Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_PointValues2.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_POINT_VALUES2_HPP
12 #define PANZER_POINT_VALUES2_HPP
13 
14 #include "PanzerDiscFE_config.hpp"
15 
16 #include "Panzer_PointRule.hpp"
17 #include "Panzer_ArrayTraits.hpp"
18 #include "Panzer_Dimension.hpp"
19 
20 #include "Teuchos_RCP.hpp"
21 
22 namespace panzer {
23 
24  template <typename Scalar>
25  class PointValues2 {
26  public:
28 
29  template<typename SourceScalar>
32 
33  PointValues2(const std::string & pre="",
34  bool allocArrays=false)
35  : alloc_arrays_(allocArrays), prefix_(pre) {}
36 
37  PointValues2(const std::string & pre,
38  const std::vector<PHX::index_size_type> & ddims,
39  bool allocArrays=false)
40  : alloc_arrays_(allocArrays), prefix_(pre), ddims_(ddims) {}
41 
44 
51  template <typename CoordinateArray,typename PointArray>
52  void evaluateValues(const CoordinateArray & node_coords,
53  const PointArray & in_point_coords,
54  const int in_num_cells = -1)
55  { copyNodeCoords(node_coords);
56  copyPointCoords(in_point_coords);
57  evaluateValues(in_num_cells); }
58 
66  template <typename PointArray>
68  const PointArray & in_point_coords,
69  bool shallow_copy_nodes,
70  const int in_num_cells = -1)
71  { if(shallow_copy_nodes)
72  node_coordinates = node_coords;
73  else
74  copyNodeCoords(node_coords);
75  copyPointCoords(in_point_coords);
76  evaluateValues(in_num_cells); }
77 
80  { return coords_ref; }
81 
85  { return node_coordinates; }
87 
90  { return node_coordinates; }
91 
92  // input fields: both mutable because of getRefCoordinates/getNodeCoordinates
93  // Not sure this is the best design, but works for this iteration
96 
97  // output fields
101  PHX::MDField<Scalar,Cell,IP,Dim> point_coords; // <Cell,IP,Dim> // cell points
102 
104 
105  private:
106  void evaluateValues(const int in_num_cells);
107 
108  template <typename CoordinateArray>
109  void copyNodeCoords(const CoordinateArray& in_node_coords);
110 
111  template <typename CoordinateArray>
112  void copyPointCoords(const CoordinateArray& in_point_coords);
113 
115  std::string prefix_;
116  std::vector<PHX::index_size_type> ddims_;
117 
118  };
119 
120  template <typename Scalar>
121  template<typename SourceScalar>
124  {
125  // The separate template parameter for SourceScalar allows for
126  // assignment to a "const Scalar" from a non-const "Scalar", but
127  // we still need to enforce that the underlying scalar type is the
128  // same.
129  static_assert(std::is_same<typename std::decay<Scalar>::type,typename std::decay<SourceScalar>::type>::value,
130  "ERROR: PointValues assignment requires consistent scalar types!");
131 
132  coords_ref = source.coords_ref;
133  node_coordinates = source.node_coordinates;
134  jac = source.jac;
135  jac_inv = source.jac_inv;
136  jac_det = source.jac_det;
137  point_coords = source.point_coords;
138  point_rule = source.point_rule;
139  return *this;
140  }
141 
142 } // namespace panzer
143 
144 #endif
Teuchos::RCP< const panzer::PointRule > point_rule
std::vector< PHX::index_size_type > ddims_
void evaluateValues(const CoordinateArray &node_coords, const PointArray &in_point_coords, const int in_num_cells=-1)
PHX::MDField< Scalar, Cell, IP, Dim, Dim > jac
PointValues2(const std::string &pre, const std::vector< PHX::index_size_type > &ddims, bool allocArrays=false)
PHX::MDField< Scalar, Cell, IP > jac_det
void setupArrays(const Teuchos::RCP< const panzer::PointRule > &pr)
Sizes/allocates memory for arrays.
ArrayTraits< const double, PHX::MDField< const double > >::size_type size_type
PointValues2< Scalar > & operator=(const PointValues2< SourceScalar > &source)
PHX::MDField< Scalar, Cell, IP, Dim > point_coords
void copyNodeCoords(const CoordinateArray &in_node_coords)
const PHX::MDField< Scalar, Cell, NODE, Dim > & getNodeCoordinates() const
Return the node coordinates this class uses (Cell,NODE,Dim) sized.
void copyPointCoords(const CoordinateArray &in_point_coords)
PHX::MDField< Scalar, Cell, NODE, Dim > node_coordinates
const PHX::MDField< Scalar, IP, Dim > & getRefCoordinates() const
Return reference cell coordinates this class uses (IP,Dim) sized.
void evaluateValues(const PHX::MDField< Scalar, Cell, NODE, Dim > &node_coords, const PointArray &in_point_coords, bool shallow_copy_nodes, const int in_num_cells=-1)
PHX::MDField< Scalar, Cell, IP, Dim, Dim > jac_inv
const PHX::MDField< Scalar, Cell, NODE, Dim > & getVertexCoordinates() const
PHX::MDField< Scalar, IP, Dim > coords_ref
PointValues2(const std::string &pre="", bool allocArrays=false)