11 #ifndef PANZER_GATHER_TANGENTS_IMPL_HPP
12 #define PANZER_GATHER_TANGENTS_IMPL_HPP
14 #include "Teuchos_Assert.hpp"
15 #include "Phalanx_DataLayout.hpp"
17 #include "Intrepid2_Kernels.hpp"
18 #include "Intrepid2_OrientationTools.hpp"
21 #include "Kokkos_ViewFactory.hpp"
23 #include "Teuchos_FancyOStream.hpp"
25 template<
typename EvalT,
typename Traits>
30 dof_name_ = (p.
get< std::string >(
"DOF Name"));
46 std::string orientationFieldName = basis_->name() +
" Orientation";
52 constJac_ = pointValues_.jac;
53 this->addDependentField(constJac_);
56 this->addEvaluatedField(gatherFieldTangents_);
58 this->setName(
"Gather Tangents");
62 template<
typename EvalT,
typename Traits>
69 auto orientations_host = Kokkos::create_mirror_view(orientations_);
70 for (
size_t i=0; i < orientations->size(); ++i)
71 orientations_host(i) = (*orientations)[i];
72 Kokkos::deep_copy(orientations_,orientations_host);
74 this->utils.setFieldData(pointValues_.jac,fm);
75 const shards::CellTopology & parentCell = *basis_->getCellTopology();
76 const int edgeDim = 1;
77 edgeParam_ = Intrepid2::RefSubcellParametrization<PHX::Device>::get(edgeDim, parentCell.getKey());
79 const int numEdges = gatherFieldTangents_.extent(1);
81 auto keys_host = Kokkos::create_mirror_view(keys_);
82 for (
int i=0; i < numEdges; ++i)
83 keys_host(i) = parentCell.getKey(edgeDim,i);
84 Kokkos::deep_copy(keys_,keys_host);
88 template<
typename EvalT,
typename Traits>
96 const shards::CellTopology & parentCell = *basis_->getCellTopology();
97 const int cellDim = parentCell.getDimension();
98 const int numEdges = gatherFieldTangents_.extent(1);
101 const auto worksetJacobians = pointValues_.jac.get_view();
103 auto gatherFieldTangents = gatherFieldTangents_.get_static_view();
105 auto orientations = orientations_;
106 auto edgeParam = edgeParam_;
109 Kokkos::parallel_for(
"panzer::GatherTangets",workset.
num_cells,KOKKOS_LAMBDA(
const int c){
110 int edgeOrts[12] = {};
111 orientations(cell_local_ids(c)).getEdgeOrientation(edgeOrts, numEdges);
113 auto workspace = Kokkos::subview(workspace_tmp,c,Kokkos::ALL(),Kokkos::ALL());
114 for(
int pt = 0; pt < numEdges; pt++) {
115 auto phyEdgeTan = Kokkos::subview(gatherFieldTangents, c, pt, Kokkos::ALL());
116 auto ortEdgeTan = Kokkos::subview(workspace_tmp,c,0,Kokkos::ALL());
118 Intrepid2::Impl::OrientationTools::getRefSubcellTangents(
125 auto J = Kokkos::subview(worksetJacobians, c, pt, Kokkos::ALL(), Kokkos::ALL());
126 Intrepid2::Kernels::Serial::matvec_product(phyEdgeTan, J, ortEdgeTan);
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_
int num_cells
DEPRECATED - use: numCells()
T & get(const std::string &name, T def_value)
GatherTangents(const Teuchos::ParameterList &p)
void setupArrays(const Teuchos::RCP< const panzer::PointRule > &pr)
Sizes/allocates memory for arrays.
Kokkos::View< const int *, PHX::Device > getLocalCellIDs() const
Get the local cell IDs for the workset.
bool isType(const std::string &name) const
#define TEUCHOS_ASSERT(assertion_test)
void evaluateFields(typename Traits::EvalData d)
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &vm)