VertexCFD  0.0-dev
VertexCFD_Integrator_BoundaryGradBasisDotVector.hpp
1 #ifndef VERTEXCFD_INTEGRATOR_BOUNDARYGRADBASISDOTVECTOR_HPP
2 #define VERTEXCFD_INTEGRATOR_BOUNDARYGRADBASISDOTVECTOR_HPP
3 
4 #include <Panzer_EvaluatorStyle.hpp>
5 #include <Panzer_Evaluator_WithBaseImpl.hpp>
6 
7 #include <Phalanx_Evaluator_Derived.hpp>
8 #include <Phalanx_MDField.hpp>
9 
10 #include <string>
11 
12 //---------------------------------------------------------------------------//
13 // Special boundary integration operator for the penalty term of the boundary
14 // flux. Only the gradient tangential to the surface of the boundary side is
15 // used.
16 //
17 // NOTE: The normal vectors used by this evaluator to compute the tangential
18 // gradient are expected to be unit normals.
19 //---------------------------------------------------------------------------//
20 namespace VertexCFD
21 {
22 namespace Integrator
23 {
24 template<typename EvalType, typename Traits>
26  : public panzer::EvaluatorWithBaseImpl<Traits>,
27  public PHX::EvaluatorDerived<EvalType, Traits>
28 {
29  public:
30  BoundaryGradBasisDotVector(const panzer::EvaluatorStyle& eval_style,
31  const std::string& res_name,
32  const std::string& flux_name,
33  const panzer::BasisIRLayout& basis,
34  const panzer::IntegrationRule& ir,
35  const double& multiplier = 1,
36  const std::vector<std::string>& fm_names
37  = std::vector<std::string>{});
38 
39  void postRegistrationSetup(typename Traits::SetupData d,
40  PHX::FieldManager<Traits>& fm);
41 
42  void evaluateFields(typename Traits::EvalData d);
43 
44  // Regular memory version.
45  template<int NUM_FIELD_MULT>
46  struct FieldMultTag
47  {
48  };
49 
50  template<int NUM_FIELD_MULT>
51  KOKKOS_INLINE_FUNCTION void operator()(
53  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const;
54 
55  // Shared memory version.
56  template<int NUM_FIELD_MULT>
58  {
59  };
60 
61  template<int NUM_FIELD_MULT>
62  KOKKOS_INLINE_FUNCTION void operator()(
64  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const;
65 
66  private:
67  Teuchos::RCP<Teuchos::ParameterList> getValidParameters() const;
68 
69  using ScalarT = typename EvalType::ScalarT;
70  using scratch_view
71  = Kokkos::View<ScalarT*,
72  typename PHX::DevLayout<ScalarT>::type,
73  typename PHX::exec_space::scratch_memory_space,
74  Kokkos::MemoryUnmanaged>;
75 
76  const panzer::EvaluatorStyle _eval_style;
77  double _multiplier;
78  std::string _basis_name;
79  std::size_t _basis_index;
80 
81  PHX::MDField<ScalarT, panzer::Cell, panzer::BASIS> _field;
82  PHX::MDField<ScalarT, panzer::Cell, panzer::BASIS, panzer::Point, panzer::Dim>
83  _boundary_grad_basis;
84 
85  std::vector<PHX::MDField<const ScalarT, panzer::Cell, panzer::Point>>
86  _field_mults;
87  Kokkos::View<Kokkos::View<const ScalarT**,
88  typename PHX::DevLayout<ScalarT>::type,
89  PHX::Device>*>
90  _kokkos_field_mults;
91 
92  PHX::MDField<const ScalarT, panzer::Cell, panzer::Point, panzer::Dim> _vector;
93  PHX::MDField<const ScalarT, panzer::Cell, panzer::Point, panzer::Dim> _normals;
94  PHX::MDField<const double, panzer::Cell, panzer::BASIS, panzer::Point, panzer::Dim>
95  _grad_basis;
96 };
97 
98 //---------------------------------------------------------------------------//
99 
100 } // end namespace Integrator
101 } // end namespace VertexCFD
102 
103 #endif // VERTEXCFD_INTEGRATOR_BOUNDARYGRADBASISDOTVECTOR_HPP
VertexCFD::Integrator::BoundaryGradBasisDotVector
Definition: VertexCFD_Integrator_BoundaryGradBasisDotVector.hpp:28
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23
VertexCFD::Integrator::BoundaryGradBasisDotVector::SharedFieldMultTag
Definition: VertexCFD_Integrator_BoundaryGradBasisDotVector.hpp:58
VertexCFD::Integrator::BoundaryGradBasisDotVector::FieldMultTag
Definition: VertexCFD_Integrator_BoundaryGradBasisDotVector.hpp:47