43 #ifndef PANZER_GATHER_NORMALS_IMPL_HPP
44 #define PANZER_GATHER_NORMALS_IMPL_HPP
46 #include "Teuchos_Assert.hpp"
47 #include "Phalanx_DataLayout.hpp"
50 #include "Kokkos_ViewFactory.hpp"
52 #include "Intrepid2_Kernels.hpp"
53 #include "Intrepid2_CellTools.hpp"
54 #include "Intrepid2_OrientationTools.hpp"
56 #include "Teuchos_FancyOStream.hpp"
58 template<
typename EvalT,
typename Traits>
63 dof_name = (p.
get< std::string >(
"DOF Name"));
79 std::string orientationFieldName = basis->name() +
" Orientation";
85 constJac_ = pointValues.jac;
86 this->addDependentField(constJac_);
88 gatherFieldNormals = PHX::MDField<ScalarT,Cell,NODE,Dim>(dof_name+
"_Normals",vector_layout_vector);
89 this->addEvaluatedField(gatherFieldNormals);
91 this->setName(
"Gather Normals");
95 template<
typename EvalT,
typename Traits>
101 this->utils.setFieldData(pointValues.jac,fm);
104 gatherFieldNormals.extent(0),
105 gatherFieldNormals.extent(1),
106 gatherFieldNormals.extent(2));
110 template<
typename EvalT,
typename Traits>
118 const shards::CellTopology & parentCell = *basis->getCellTopology();
119 int cellDim = parentCell.getDimension();
120 int numFaces = gatherFieldNormals.extent(1);
127 const auto worksetJacobians = pointValues.jac.get_view();
130 for(index_t c=0;c<workset.
num_cells;c++) {
131 int faceOrts[6] = {};
132 orientations->at(details.cell_local_ids[c]).getFaceOrientation(faceOrts, numFaces);
134 for(
int pt = 0; pt < numFaces; pt++) {
135 auto ortEdgeTan_U = Kokkos::subview(refEdges, 0, Kokkos::ALL());
136 auto ortEdgeTan_V = Kokkos::subview(refEdges, 1, Kokkos::ALL());
139 Intrepid2::Orientation::getReferenceFaceTangents(ortEdgeTan_U,
145 auto phyEdgeTan_U = Kokkos::subview(phyEdges, 0, Kokkos::ALL());
146 auto phyEdgeTan_V = Kokkos::subview(phyEdges, 1, Kokkos::ALL());
147 auto J = Kokkos::subview(worksetJacobians, c, pt, Kokkos::ALL(), Kokkos::ALL());
149 Intrepid2::Kernels::Serial::matvec_product(phyEdgeTan_U, J, ortEdgeTan_U);
150 Intrepid2::Kernels::Serial::matvec_product(phyEdgeTan_V, J, ortEdgeTan_V);
153 gatherFieldNormals(c,pt,0) = (phyEdgeTan_U(1)*phyEdgeTan_V(2) - phyEdgeTan_U(2)*phyEdgeTan_V(1));
154 gatherFieldNormals(c,pt,1) = (phyEdgeTan_U(2)*phyEdgeTan_V(0) - phyEdgeTan_U(0)*phyEdgeTan_V(2));
155 gatherFieldNormals(c,pt,2) = (phyEdgeTan_U(0)*phyEdgeTan_V(1) - phyEdgeTan_U(1)*phyEdgeTan_V(0));
Kokkos::DynRankView< typename InputArray::value_type, PHX::Device > createDynRankView(const InputArray &a, const std::string &name, const DimensionPack...dims)
Wrapper to simplify Panzer use of Sacado ViewFactory.
Teuchos::RCP< const std::vector< Intrepid2::Orientation > > orientations_
T & get(const std::string &name, T def_value)
void setupArrays(const Teuchos::RCP< const panzer::PointRule > &pr)
Sizes/allocates memory for arrays.
void evaluateFields(typename Traits::EvalData d)
bool isType(const std::string &name) const
#define TEUCHOS_ASSERT(assertion_test)
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &vm)