1 #ifndef VERTEXCFD_CLOSURE_METRICTENSORELEMENTLENGTH_IMPL_HPP
2 #define VERTEXCFD_CLOSURE_METRICTENSORELEMENTLENGTH_IMPL_HPP
4 #include <Panzer_HierarchicParallelism.hpp>
5 #include <Panzer_Workset_Utilities.hpp>
11 namespace ClosureModel
14 template<
class EvalType,
class Traits>
15 MetricTensorElementLength<EvalType, Traits>::MetricTensorElementLength(
16 const panzer::IntegrationRule& ir,
const std::string& prefix)
17 : _element_length(prefix +
"element_length", ir.dl_vector)
18 , _metric_tensor(
"metric_tensor", ir.dl_tensor)
20 this->addEvaluatedField(_element_length);
21 this->addDependentField(_metric_tensor);
22 this->setName(
"Metric Tensor Element Length");
26 template<
class EvalType,
class Traits>
27 void MetricTensorElementLength<EvalType, Traits>::evaluateFields(
28 typename Traits::EvalData workset)
30 auto policy = panzer::HP::inst().teamPolicy<scalar_type, PHX::Device>(
32 Kokkos::parallel_for(this->getName(), policy, *
this);
36 template<
class EvalType,
class Traits>
37 KOKKOS_INLINE_FUNCTION
void
38 MetricTensorElementLength<EvalType, Traits>::operator()(
39 const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team)
const
41 const int cell = team.league_rank();
42 const int num_point = _metric_tensor.extent(1);
43 const int num_space_dim = _metric_tensor.extent(2);
45 Kokkos::parallel_for(Kokkos::TeamThreadRange(team, 0, num_point),
46 [&](
const int point) {
48 for (
int i = 0; i < num_space_dim; ++i)
50 _element_length(cell, point, i)
51 = sqrt(_metric_tensor(cell, point, i, i));
61 #endif // end VERTEXCFD_CLOSURE_METRICTENSORELEMENTLENGTH_IMPL_HPP