43 #ifndef PANZER_PROJECT_TO_FACES_IMPL_HPP
44 #define PANZER_PROJECT_TO_FACES_IMPL_HPP
46 #include "Teuchos_Assert.hpp"
47 #include "Phalanx_DataLayout.hpp"
49 #include "Intrepid2_Cubature.hpp"
50 #include "Intrepid2_DefaultCubatureFactory.hpp"
51 #include "Intrepid2_FunctionSpaceTools.hpp"
52 #include "Intrepid2_Kernels.hpp"
53 #include "Intrepid2_CellTools.hpp"
54 #include "Intrepid2_OrientationTools.hpp"
60 #include "Kokkos_ViewFactory.hpp"
62 #include "Teuchos_FancyOStream.hpp"
66 template<
typename EvalT,
typename Traits>
70 dof_name_ = (p.
get< std::string >(
"DOF Name"));
84 this->addEvaluatedField(result_);
87 this->addDependentField(normals_);
90 this->addDependentField(vector_values_);
92 this->setName(
"Project To Faces");
96 template<
typename EvalT,
typename Traits>
101 num_faces_ = result_.extent(1);
102 num_dim_ = vector_values_.extent(2);
109 template<
typename EvalT,
typename Traits>
120 const int intDegree = basis_->order();
124 auto vector_values = vector_values_.get_static_view();
125 auto normals = normals_.get_static_view();
126 auto num_faces = num_faces_;
127 auto num_dim = num_dim_;
130 Kokkos::parallel_for(
"panzer::ProjectToFaces",policy,KOKKOS_LAMBDA(
const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) {
131 const auto cell = team.league_rank();
132 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,num_faces),[&] (
const int p) {
134 for (
int dim = 0; dim < num_dim; ++dim)
135 result(cell,p) += vector_values(cell,p,dim) * normals(cell,p,dim);
int num_cells
DEPRECATED - use: numCells()
void evaluateFields(typename Traits::EvalData d)
T & get(const std::string &name, T def_value)
Kokkos::TeamPolicy< TeamPolicyProperties...> teamPolicy(const int &league_size)
Returns a TeamPolicy for hierarchic parallelism.
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &vm)
PHX::MDField< ScalarT, panzer::Cell, panzer::IP > result
A field that will be used to build up the result of the integral we're performing.
ProjectToFaces(const Teuchos::ParameterList &p)
KOKKOS_FORCEINLINE_FUNCTION array_type get_static_view()
bool isType(const std::string &name) const
static HP & inst()
Private ctor.
#define TEUCHOS_ASSERT(assertion_test)