Panzer Version of the Day
Loading...
Searching...
No Matches
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"
50#include "Panzer_Dimension.hpp"
51
52#include "Teuchos_RCP.hpp"
53
54namespace panzer {
55
56 template <typename Scalar>
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 const int in_num_cells = -1)
87 { copyNodeCoords(node_coords);
88 copyPointCoords(in_point_coords);
89 evaluateValues(in_num_cells); }
90
98 template <typename PointArray>
99 void evaluateValues(const PHX::MDField<Scalar,Cell,NODE,Dim> & node_coords,
100 const PointArray & in_point_coords,
101 bool shallow_copy_nodes,
102 const int in_num_cells = -1)
103 { if(shallow_copy_nodes)
104 node_coordinates = node_coords;
105 else
106 copyNodeCoords(node_coords);
107 copyPointCoords(in_point_coords);
108 evaluateValues(in_num_cells); }
109
111 const PHX::MDField<Scalar,IP,Dim> & getRefCoordinates() const
112 { return coords_ref; }
113
115 const PHX::MDField<Scalar,Cell,NODE,Dim> & getVertexCoordinates() const
116 { return node_coordinates; }
117
118 // input fields: both mutable because of getRefCoordinates/getVertexCoordinates
119 // Not sure this is the best design, but works for this iteration
120 mutable PHX::MDField<Scalar,IP,Dim> coords_ref; // <IP,Dim>
121 mutable PHX::MDField<Scalar,Cell,NODE,Dim> node_coordinates; // <Cell,NODE,Dim>
122
123 // output fields
124 PHX::MDField<Scalar,Cell,IP,Dim,Dim> jac; // <Cell,IP,Dim,Dim>
125 PHX::MDField<Scalar,Cell,IP,Dim,Dim> jac_inv; // <Cell,IP,Dim,Dim>
126 PHX::MDField<Scalar,Cell,IP> jac_det; // <Cell,IP>
127 PHX::MDField<Scalar,Cell,IP,Dim> point_coords; // <Cell,IP,Dim> // cell points
128
129 Teuchos::RCP<const panzer::PointRule> point_rule;
130
131 private:
132 void evaluateValues(const int in_num_cells);
133
134 template <typename CoordinateArray>
135 void copyNodeCoords(const CoordinateArray& in_node_coords);
136
137 template <typename CoordinateArray>
138 void copyPointCoords(const CoordinateArray& in_point_coords);
139
141 std::string prefix_;
142 std::vector<PHX::index_size_type> ddims_;
143
144 };
145
146 template <typename Scalar>
147 template<typename SourceScalar>
150 {
151 // The separate template parameter for SourceScalar allows for
152 // assignment to a "const Scalar" from a non-const "Scalar", but
153 // we still need to enforce that the underlying scalar type is the
154 // same.
155 static_assert(std::is_same<typename std::decay<Scalar>::type,typename std::decay<SourceScalar>::type>::value,
156 "ERROR: PointValues assignment requires consistent scalar types!");
157
158 coords_ref = source.coords_ref;
159 node_coordinates = source.node_coordinates;
160 jac = source.jac;
161 jac_inv = source.jac_inv;
162 jac_det = source.jac_det;
163 point_coords = source.point_coords;
164 point_rule = source.point_rule;
165 return *this;
166 }
167
168} // namespace panzer
169
170#endif
const PHX::MDField< Scalar, IP, Dim > & getRefCoordinates() const
Return reference cell coordinates this class uses (IP,Dim) sized.
PointValues2(const std::string &pre, const std::vector< PHX::index_size_type > &ddims, bool allocArrays=false)
ArrayTraits< Scalar, PHX::MDField< Scalar > >::size_type size_type
PHX::MDField< Scalar, IP, Dim > coords_ref
void copyPointCoords(const CoordinateArray &in_point_coords)
void evaluateValues(const CoordinateArray &node_coords, const PointArray &in_point_coords, const int in_num_cells=-1)
std::vector< PHX::index_size_type > ddims_
const PHX::MDField< Scalar, Cell, NODE, Dim > & getVertexCoordinates() const
Return the vertex coordinates this class uses (Cell,NODE,Dim) sized.
PHX::MDField< Scalar, Cell, NODE, Dim > node_coordinates
void copyNodeCoords(const CoordinateArray &in_node_coords)
PointValues2(const std::string &pre="", bool allocArrays=false)
PointValues2< Scalar > & operator=(const PointValues2< SourceScalar > &source)
PHX::MDField< Scalar, Cell, IP, Dim, Dim > jac
PHX::MDField< Scalar, Cell, IP, Dim > point_coords
Teuchos::RCP< const panzer::PointRule > point_rule
PHX::MDField< Scalar, Cell, IP > jac_det
void setupArrays(const Teuchos::RCP< const panzer::PointRule > &pr)
Sizes/allocates memory for arrays.
PHX::MDField< Scalar, Cell, IP, Dim, Dim > jac_inv
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)