VertexCFD  0.0-dev
VertexCFD_Closure_FullInductionLocalTimeStepSize.hpp
1 #ifndef VERTEXCFD_CLOSURE_FULLINDUCTIONLOCALTIMESTEPSIZE_HPP
2 #define VERTEXCFD_CLOSURE_FULLINDUCTIONLOCALTIMESTEPSIZE_HPP
3 
4 #include "full_induction_mhd_solver/mhd_properties/VertexCFD_FullInductionMHDProperties.hpp"
5 
6 #include "utils/VertexCFD_Utils_MagneticDim.hpp"
7 
8 #include <Panzer_Dimension.hpp>
9 #include <Panzer_Evaluator_WithBaseImpl.hpp>
10 
11 #include <Phalanx_Evaluator_Derived.hpp>
12 #include <Phalanx_Evaluator_WithBaseImpl.hpp>
13 #include <Phalanx_FieldManager.hpp>
14 #include <Phalanx_config.hpp>
15 
16 #include <Kokkos_Core.hpp>
17 
18 namespace VertexCFD
19 {
20 namespace ClosureModel
21 {
22 //---------------------------------------------------------------------------//
23 // Compute local time step size based on the mesh size and maximum eigenvalue
24 // for the full induction MHD equation set (under incompressible assumption).
25 // The time step size used by the solver and based on the CFL condition is
26 // computed in the observer
27 // 'VertexCFD_TempusTimeStepControl_GlobalCFL_impl.hpp'
28 //---------------------------------------------------------------------------//
29 template<class EvalType, class Traits, int NumSpaceDim>
31  : public panzer::EvaluatorWithBaseImpl<Traits>,
32  public PHX::EvaluatorDerived<EvalType, Traits>
33 {
34  public:
35  using scalar_type = typename EvalType::ScalarT;
36  static constexpr int num_space_dim = NumSpaceDim;
37 
39  const panzer::IntegrationRule& ir,
41 
42  void evaluateFields(typename Traits::EvalData d) override;
43 
44  KOKKOS_INLINE_FUNCTION
45  void operator()(
46  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const;
47 
48  public:
49  PHX::MDField<scalar_type, panzer::Cell, panzer::Point> _local_dt;
50 
51  private:
52  double _magnetic_permeability;
53  double _c_h;
54 
55  PHX::MDField<const scalar_type, panzer::Cell, panzer::Point> _rho;
56  PHX::MDField<const double, panzer::Cell, panzer::Point, panzer::Dim>
57  _element_length;
58  Kokkos::Array<PHX::MDField<const scalar_type, panzer::Cell, panzer::Point>,
59  num_space_dim>
60  _velocity;
61  PHX::MDField<const scalar_type, panzer::Cell, panzer::Point, MagneticDim>
62  _total_magnetic_field;
63 };
64 
65 //---------------------------------------------------------------------------//
66 
67 } // end namespace ClosureModel
68 } // end namespace VertexCFD
69 
70 #endif // end VERTEXCFD_CLOSURE_FULLINDUCTIONLOCALTIMESTEPSIZE_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23
VertexCFD::MHDProperties::FullInductionMHDProperties
Definition: VertexCFD_FullInductionMHDProperties.hpp:17
VertexCFD::ClosureModel::FullInductionLocalTimeStepSize
Definition: VertexCFD_Closure_FullInductionLocalTimeStepSize.hpp:33