Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_DOF_PointField_decl.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_EVALUATOR_DOF_PointField_IMPL_HPP
12 #define PANZER_EVALUATOR_DOF_PointField_IMPL_HPP
13 
14 #include <string>
15 
16 #include "Phalanx_Evaluator_Macros.hpp"
17 #include "Phalanx_MDField.hpp"
18 #include "Phalanx_DataLayout.hpp"
19 #include "PanzerDiscFE_config.hpp"
20 #include "Panzer_PureBasis.hpp"
21 #include "Intrepid2_Basis.hpp"
23 
24 namespace panzer {
25 
27 template <typename EvalT, typename TRAITST>
29  : public panzer::EvaluatorWithBaseImpl<TRAITST>,
30  public PHX::EvaluatorDerived<EvalT, TRAITST> {
31 public:
32 
48  DOF_PointField(const std::string & postfixFieldName,
49  const std::string & fieldName,
50  const PureBasis & fieldBasis,
51  const std::string & coordinateName,
52  const Teuchos::RCP<PHX::DataLayout> & coordLayout,
53  const Teuchos::RCP<PHX::DataLayout> & quadLayout)
54  { initialize(fieldName,fieldBasis,coordinateName,coordLayout,quadLayout,postfixFieldName); }
55 
71  DOF_PointField(const std::string & fieldName,
72  const PureBasis & fieldBasis,
73  const std::string & coordinateName,
74  const Teuchos::RCP<PHX::DataLayout> & coordLayout,
75  const Teuchos::RCP<PHX::DataLayout> & quadLayout,
76  bool useCoordPostfix)
77  { std::string postfixFieldName = (useCoordPostfix ? coordinateName : "");
78  initialize(fieldName,fieldBasis,coordinateName,coordLayout,quadLayout,postfixFieldName); }
79 
80  void evaluateFields(typename TRAITST::EvalData workset);
81 
82 private:
83  typedef typename EvalT::ScalarT ScalarT;
84 
86  void initialize(const std::string & fieldName,
87  const PureBasis & fieldBasis,
88  const std::string & coordinateName,
89  const Teuchos::RCP<PHX::DataLayout> & coordLayout,
90  const Teuchos::RCP<PHX::DataLayout> & quadLayout,
91  const std::string & postfixFieldName);
92 
94  PHX::MDField<const ScalarT> dof_coeff; // coefficient values
95  PHX::MDField<ScalarT> dof_field; // evaluate field
96 
98  Kokkos::DynRankView<double,PHX::Device> intrpCoords, basisRef, basis;
99 };
100 
101 }
102 
103 #endif
Kokkos::DynRankView< double, PHX::Device > intrpCoords
void initialize(const std::string &fieldName, const PureBasis &fieldBasis, const std::string &coordinateName, const Teuchos::RCP< PHX::DataLayout > &coordLayout, const Teuchos::RCP< PHX::DataLayout > &quadLayout, const std::string &postfixFieldName)
Convenience initialization routine, see constructor above.
Interpolates basis DOF using reference coordinates defined by a field.
Kokkos::DynRankView< double, PHX::Device > basis
DOF_PointField(const std::string &postfixFieldName, const std::string &fieldName, const PureBasis &fieldBasis, const std::string &coordinateName, const Teuchos::RCP< PHX::DataLayout > &coordLayout, const Teuchos::RCP< PHX::DataLayout > &quadLayout)
void evaluateFields(typename TRAITST::EvalData workset)
Wrapper to PHX::EvaluatorWithBaseImpl that implements Panzer-specific helpers.
Kokkos::DynRankView< double, PHX::Device > basisRef
DOF_PointField(const std::string &fieldName, const PureBasis &fieldBasis, const std::string &coordinateName, const Teuchos::RCP< PHX::DataLayout > &coordLayout, const Teuchos::RCP< PHX::DataLayout > &quadLayout, bool useCoordPostfix)
PHX::MDField< ScalarT > dof_field
Teuchos::RCP< Intrepid2::Basis< PHX::exec_space, double, double > > intrepidBasis
PHX::MDField< const ScalarT, Point, Dim > coordinates
Description and data layouts associated with a particular basis.
PHX::MDField< const ScalarT > dof_coeff