Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_PointValues_Evaluator_impl.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_PointValues_Evaluator_IMPL_HPP
12 #define PANZER_PointValues_Evaluator_IMPL_HPP
13 
14 #include <algorithm>
15 #include "Panzer_PointRule.hpp"
18 
19 namespace panzer {
20 
21 //**********************************************************************
22 template<typename EvalT, typename Traits>
25  const Teuchos::ParameterList& p)
26 {
27  basis_index = 0;
28 
30  = p.get< Teuchos::RCP<const panzer::PointRule> >("Point Rule");
33 
34  initialize(pointRule,userArray.ptr(),Teuchos::null);
35 }
36 
37 //**********************************************************************
38 template <typename EvalT, typename TRAITST>
40  const Kokkos::DynRankView<double,PHX::Device> & userArray)
41 {
42  basis_index = 0;
43 
44  initialize(pointRule,Teuchos::ptrFromRef(userArray),Teuchos::null);
45 }
46 
47 //**********************************************************************
48 template <typename EvalT, typename TRAITST>
51 {
52  basis_index = 0;
53 
54  initialize(pointRule,Teuchos::ptrFromRef(userArray),Teuchos::null);
55 }
56 
57 //**********************************************************************
58 template <typename EvalT, typename TRAITST>
60  const Teuchos::RCP<const panzer::PureBasis> & pureBasis)
61 {
62  basis_index = 0;
63 
65  initialize(pointRule,userArray,pureBasis);
66 }
67 
68 //**********************************************************************
69 template <typename EvalT, typename TRAITST>
70 template <typename ArrayT>
72  const Teuchos::Ptr<const ArrayT> & userArray,
73  const Teuchos::RCP<const panzer::PureBasis> & pureBasis)
74 {
75  basis = pureBasis;
76 
77  if(userArray!=Teuchos::null && basis==Teuchos::null)
78  useBasisValuesRefArray = false;
79  else if(userArray==Teuchos::null && basis!=Teuchos::null)
80  useBasisValuesRefArray = true;
81  else {
82  // this is a conflicting request, throw an exception
83  TEUCHOS_ASSERT(false);
84  }
85 
86  // copy user array data
87  if(userArray!=Teuchos::null) {
88  TEUCHOS_ASSERT(userArray->rank()==2);
89  MDFieldArrayFactory md_af("refPointArray",true);
90 
91  refPointArray = md_af.buildStaticArray<double,NODE,Dim>("refPointArray",userArray->extent(0),userArray->extent(1));
92  Kokkos::deep_copy(PHX::as_view(refPointArray), PHX::as_view(*userArray));
93 
94  }
95 
96  // setup all fields to be evaluated and constructed
97  pointValues = PointValues2<double>(pointRule->getName()+"_",false);
98  pointValues.setupArrays(pointRule);
99 
100  // the field manager will allocate all of these field
101  this->addEvaluatedField(pointValues.coords_ref);
102  this->addEvaluatedField(pointValues.node_coordinates);
103  this->addEvaluatedField(pointValues.jac);
104  this->addEvaluatedField(pointValues.jac_inv);
105  this->addEvaluatedField(pointValues.jac_det);
106  this->addEvaluatedField(pointValues.point_coords);
107 
108  std::string n = "PointValues_Evaluator: " + pointRule->getName();
109  this->setName(n);
110 }
111 
112 //**********************************************************************
113 template<typename EvalT, typename Traits>
114 void
117  typename Traits::SetupData sd,
119 {
120  // setup the pointers for evaluation
121  this->utils.setFieldData(pointValues.coords_ref,fm);
122  this->utils.setFieldData(pointValues.node_coordinates,fm);
123  this->utils.setFieldData(pointValues.jac,fm);
124  this->utils.setFieldData(pointValues.jac_inv,fm);
125  this->utils.setFieldData(pointValues.jac_det,fm);
126  this->utils.setFieldData(pointValues.point_coords,fm);
127 
128  if(useBasisValuesRefArray) {
129  basis_index = panzer::getPureBasisIndex(basis->name(), (*sd.worksets_)[0], this->wda);
130 
131  // basis better have coordinates if you want to use them! Assertion to protect
132  // a silent failure.
133  TEUCHOS_ASSERT(basis->supportsBasisCoordinates());
134  }
135 }
136 
137 //**********************************************************************
138 template<typename EvalT, typename Traits>
139 void
142  typename Traits::EvalData workset)
143 {
144  if(useBasisValuesRefArray) {
145  panzer::BasisValues2<double> & basisValues = *this->wda(workset).bases[basis_index];
146 
147  // evaluate the point values (construct jacobians etc...)
148  pointValues.evaluateValues(this->wda(workset).cell_node_coordinates,
149  basisValues.basis_coordinates_ref,
150  workset.num_cells);
151  }
152  else {
153  // evaluate the point values (construct jacobians etc...)
154  pointValues.evaluateValues(this->wda(workset).cell_node_coordinates,refPointArray,
155  workset.num_cells);
156  }
157 }
158 
159 //**********************************************************************
160 
161 }
162 
163 #endif
PointValues_Evaluator(const Teuchos::ParameterList &p)
int num_cells
DEPRECATED - use: numCells()
T & get(const std::string &name, T def_value)
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
const std::string & getName() const
void setupArrays(const Teuchos::RCP< const panzer::PointRule > &pr)
Sizes/allocates memory for arrays.
void initialize(const Teuchos::RCP< const panzer::PointRule > &pointRule, const Teuchos::Ptr< const ArrayT > &userArray, const Teuchos::RCP< const panzer::PureBasis > &pureBasis)
Initialization method to unify the constructors.
Ptr< T > ptr() const
void evaluateFields(typename Traits::EvalData d)
std::vector< std::string >::size_type getPureBasisIndex(std::string basis_name, const panzer::Workset &workset, WorksetDetailsAccessor &wda)
Returns the index in the workset bases for a particular PureBasis name.
Array_BasisDim basis_coordinates_ref
#define TEUCHOS_ASSERT(assertion_test)
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)
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_