1 #ifndef VERTEXCFD_BOUNDARYSTATE_INCOMPRESSIBLEROTATINGWALL_HPP
2 #define VERTEXCFD_BOUNDARYSTATE_INCOMPRESSIBLEROTATINGWALL_HPP
4 #include <Panzer_Dimension.hpp>
5 #include <Panzer_Evaluator_WithBaseImpl.hpp>
7 #include <Phalanx_Evaluator_Derived.hpp>
8 #include <Phalanx_Evaluator_WithBaseImpl.hpp>
9 #include <Phalanx_FieldManager.hpp>
10 #include <Phalanx_config.hpp>
12 #include <Teuchos_ParameterList.hpp>
14 #include <Kokkos_Core.hpp>
18 namespace BoundaryCondition
23 template<
class EvalType,
class Traits,
int NumSpaceDim>
25 :
public panzer::EvaluatorWithBaseImpl<Traits>,
26 public PHX::EvaluatorDerived<EvalType, Traits>
29 using scalar_type =
typename EvalType::ScalarT;
30 static constexpr
int num_space_dim = NumSpaceDim;
33 const Teuchos::ParameterList& fluid_params,
34 const Teuchos::ParameterList& bc_params,
37 void postRegistrationSetup(
typename Traits::SetupData sd,
38 PHX::FieldManager<Traits>& fm)
override;
40 void evaluateFields(
typename Traits::EvalData workset)
override;
42 KOKKOS_INLINE_FUNCTION
44 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const;
47 PHX::MDField<scalar_type, panzer::Cell, panzer::Point>
48 _boundary_lagrange_pressure;
49 PHX::MDField<scalar_type, panzer::Cell, panzer::Point, panzer::Dim>
50 _boundary_grad_lagrange_pressure;
51 PHX::MDField<scalar_type, panzer::Cell, panzer::Point> _boundary_temperature;
52 Kokkos::Array<PHX::MDField<scalar_type, panzer::Cell, panzer::Point>,
57 PHX::MDField<scalar_type, panzer::Cell, panzer::Point, panzer::Dim>,
59 _boundary_grad_velocity;
60 PHX::MDField<scalar_type, panzer::Cell, panzer::Point, panzer::Dim>
61 _boundary_grad_temperature;
65 bool _set_lagrange_pressure;
68 double _a_vel, _b_vel;
69 double _time_init, _time_final;
70 double _angular_velocity;
75 PHX::MDField<const scalar_type, panzer::Cell, panzer::Point> _lagrange_pressure;
76 PHX::MDField<const scalar_type, panzer::Cell, panzer::Point, panzer::Dim>
77 _grad_lagrange_pressure;
80 PHX::MDField<const scalar_type, panzer::Cell, panzer::Point, panzer::Dim>,
84 PHX::MDField<const scalar_type, panzer::Cell, panzer::Point, panzer::Dim>
87 PHX::MDField<const scalar_type, panzer::Cell, panzer::Point, panzer::Dim>
90 PHX::MDField<const double, panzer::Cell, panzer::Point, panzer::Dim> _ip_coords;
98 #endif // VERTEXCFD_BOUNDARYSTATE_INCOMPRESSIBLEROTATINGWALL_HPP