VertexCFD  0.0-dev
VertexCFD_Closure_MagneticCorrectionDampingSource_impl.hpp
1 #ifndef VERTEXCFD_CLOSURE_MAGNETICCORRECTIONDAMPINGSOURCE_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_MAGNETICCORRECTIONDAMPINGSOURCE_IMPL_HPP
3 
4 #include <Panzer_HierarchicParallelism.hpp>
5 
6 namespace VertexCFD
7 {
8 namespace ClosureModel
9 {
10 //---------------------------------------------------------------------------//
11 template<class EvalType, class Traits>
12 MagneticCorrectionDampingSource<EvalType, Traits>::MagneticCorrectionDampingSource(
13  const panzer::IntegrationRule& ir,
14  const MHDProperties::FullInductionMHDProperties& mhd_props)
15  : _damping_potential_source("DAMPING_SOURCE_magnetic_correction_potential",
16  ir.dl_scalar)
17  , _alpha(mhd_props.magneticCorrectionDampingFactor())
18  , _scalar_magnetic_potential("scalar_magnetic_potential", ir.dl_scalar)
19 {
20  // Evaluated fields
21  this->addEvaluatedField(_damping_potential_source);
22 
23  // Dependent fields
24  this->addDependentField(_scalar_magnetic_potential);
25 
26  this->setName("Magnetic Correction Damping Source");
27 }
28 
29 //---------------------------------------------------------------------------//
30 template<class EvalType, class Traits>
31 void MagneticCorrectionDampingSource<EvalType, Traits>::evaluateFields(
32  typename Traits::EvalData workset)
33 {
34  auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
35  workset.num_cells);
36  Kokkos::parallel_for(this->getName(), policy, *this);
37 }
38 
39 //---------------------------------------------------------------------------//
40 template<class EvalType, class Traits>
41 KOKKOS_INLINE_FUNCTION void
42 MagneticCorrectionDampingSource<EvalType, Traits>::operator()(
43  const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) const
44 {
45  const int cell = team.league_rank();
46  const int num_point = _scalar_magnetic_potential.extent(1);
47 
48  Kokkos::parallel_for(
49  Kokkos::TeamThreadRange(team, 0, num_point), [&](const int point) {
50  _damping_potential_source(cell, point)
51  = -_alpha * _scalar_magnetic_potential(cell, point);
52  });
53 }
54 
55 //---------------------------------------------------------------------------//
56 
57 } // end namespace ClosureModel
58 } // end namespace VertexCFD
59 
60 #endif // end VERTEXCFD_CLOSURE_MAGNETICCORRECTIONDAMPINGSOURCE_IMPL_HPP
VertexCFD
Definition: tstMethodManufacturedSolutionBC.cpp:23