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"
66template<
typename EvalT,
typename Traits>
70 dof_name_ = (p.get< std::string >(
"DOF Name"));
72 if(p.isType< Teuchos::RCP<PureBasis> >(
"Basis"))
73 basis_ = p.get< Teuchos::RCP<PureBasis> >(
"Basis");
75 basis_ = p.get< Teuchos::RCP<const PureBasis> >(
"Basis");
77 Teuchos::RCP<PHX::DataLayout> basis_layout = basis_->functional;
78 Teuchos::RCP<PHX::DataLayout> vector_layout = basis_->functional_grad;
81 TEUCHOS_ASSERT(basis_->isVectorBasis());
83 result_ = PHX::MDField<ScalarT,Cell,BASIS>(dof_name_,basis_layout);
84 this->addEvaluatedField(result_);
86 normals_ = PHX::MDField<const ScalarT,Cell,BASIS,Dim>(dof_name_+
"_Normals",vector_layout);
87 this->addDependentField(normals_);
89 vector_values_ = PHX::MDField<const ScalarT,Cell,BASIS,Dim>(dof_name_+
"_Vector",vector_layout);
90 this->addDependentField(vector_values_);
92 this->setName(
"Project To Faces");
96template<
typename EvalT,
typename Traits>
101 num_faces_ = result_.extent(1);
102 num_dim_ = vector_values_.extent(2);
104 TEUCHOS_ASSERT(result_.extent(1) == normals_.extent(1));
105 TEUCHOS_ASSERT(vector_values_.extent(2) == normals_.extent(2));
109template<
typename EvalT,
typename Traits>
120 const int intDegree = basis_->order();
121 TEUCHOS_ASSERT(intDegree == 1);
123 auto result = result_.get_static_view();
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);
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.
static HP & inst()
Private ctor.
Kokkos::TeamPolicy< TeamPolicyProperties... > teamPolicy(const int &league_size)
Returns a TeamPolicy for hierarchic parallelism.
ProjectToFaces(const Teuchos::ParameterList &p)
void evaluateFields(typename Traits::EvalData d)
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &vm)
int num_cells
DEPRECATED - use: numCells()