Panzer  Version of the Day
Panzer_PointValues2.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_POINT_VALUES2_HPP
44 #define PANZER_POINT_VALUES2_HPP
45 
46 #include "PanzerDiscFE_config.hpp"
47 
48 #include "Panzer_PointRule.hpp"
49 #include "Panzer_ArrayTraits.hpp"
50 #include "Panzer_Dimension.hpp"
51 
52 #include "Teuchos_RCP.hpp"
53 
54 namespace panzer {
55 
56  template <typename Scalar>
57  class PointValues2 {
58  public:
60 
61  template<typename SourceScalar>
64 
65  PointValues2(const std::string & pre="",
66  bool allocArrays=false)
67  : alloc_arrays_(allocArrays), prefix_(pre) {}
68 
69  PointValues2(const std::string & pre,
70  const std::vector<PHX::index_size_type> & ddims,
71  bool allocArrays=false)
72  : alloc_arrays_(allocArrays), prefix_(pre), ddims_(ddims) {}
73 
75  void setupArrays(const Teuchos::RCP<const panzer::PointRule>& pr);
76 
83  template <typename CoordinateArray,typename PointArray>
84  void evaluateValues(const CoordinateArray & node_coords,
85  const PointArray & in_point_coords)
86  { copyNodeCoords(node_coords);
87  copyPointCoords(in_point_coords);
88  evaluateValues(); }
89 
97  template <typename PointArray>
98  void evaluateValues(const PHX::MDField<Scalar,Cell,NODE,Dim> & node_coords,
99  const PointArray & in_point_coords,
100  bool shallow_copy_nodes)
101  { if(shallow_copy_nodes)
102  node_coordinates = node_coords;
103  else
104  copyNodeCoords(node_coords);
105  copyPointCoords(in_point_coords);
106  evaluateValues(); }
107 
109  PHX::MDField<Scalar,IP,Dim> & getRefCoordinates() const
110  { return coords_ref; }
111 
113  PHX::MDField<Scalar,Cell,NODE,Dim> & getVertexCoordinates() const
114  { return node_coordinates; }
115 
116  // input fields: both mutable because of getRefCoordinates/getVertexCoordinates
117  // Not sure this is the best design, but works for this iteration
118  mutable PHX::MDField<Scalar,IP,Dim> coords_ref; // <IP,Dim>
119  mutable PHX::MDField<Scalar,Cell,NODE,Dim> node_coordinates; // <Cell,NODE,Dim>
120 
121  // output fields
122  PHX::MDField<Scalar,Cell,IP,Dim,Dim> jac; // <Cell,IP,Dim,Dim>
123  PHX::MDField<Scalar,Cell,IP,Dim,Dim> jac_inv; // <Cell,IP,Dim,Dim>
124  PHX::MDField<Scalar,Cell,IP> jac_det; // <Cell,IP>
125  PHX::MDField<Scalar,Cell,IP,Dim> point_coords; // <Cell,IP,Dim> // cell points
126 
127  Teuchos::RCP<const panzer::PointRule> point_rule;
128 
129  private:
130  void evaluateValues();
131 
132  template <typename CoordinateArray>
133  void copyNodeCoords(const CoordinateArray& in_node_coords);
134 
135  template <typename CoordinateArray>
136  void copyPointCoords(const CoordinateArray& in_point_coords);
137 
139  std::string prefix_;
140  std::vector<PHX::index_size_type> ddims_;
141 
142  };
143 
144  template <typename Scalar>
145  template<typename SourceScalar>
148  {
149  // The separate template parameter for SourceScalar allows for
150  // assignment to a "const Scalar" from a non-const "Scalar", but
151  // we still need to enforce that the underlying scalar type is the
152  // same.
153  static_assert(std::is_same<typename std::decay<Scalar>::type,typename std::decay<SourceScalar>::type>::value,
154  "ERROR: PointValues assignment requires consistent scalar types!");
155 
156  coords_ref = source.coords_ref;
157  node_coordinates = source.node_coordinates;
158  jac = source.jac;
159  jac_inv = source.jac_inv;
160  jac_det = source.jac_det;
161  point_coords = source.point_coords;
162  point_rule = source.point_rule;
163  return *this;
164  }
165 
166 } // namespace panzer
167 
168 #endif
Teuchos::RCP< const panzer::PointRule > point_rule
std::vector< PHX::index_size_type > ddims_
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)
PHX::MDField< Scalar, Cell, NODE, Dim > & getVertexCoordinates() const
Return the vertex coordinates this class uses (Cell,NODE,Dim) sized.
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< ScalarT, PHX::MDField< ScalarT > >::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)
void copyPointCoords(const CoordinateArray &in_point_coords)
PHX::MDField< Scalar, Cell, NODE, Dim > node_coordinates
PHX::MDField< Scalar, Cell, IP, Dim, Dim > jac_inv
PHX::MDField< Scalar, IP, Dim > coords_ref
void evaluateValues(const CoordinateArray &node_coords, const PointArray &in_point_coords)
PointValues2(const std::string &pre="", bool allocArrays=false)