Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_GatherNormals_impl.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Panzer: A partial differential equation assembly
4// engine for strongly coupled complex multiphysics systems
5//
6// Copyright 2011 NTESS and the Panzer contributors.
7// SPDX-License-Identifier: BSD-3-Clause
8// *****************************************************************************
9// @HEADER
10
11#ifndef PANZER_GATHER_NORMALS_IMPL_HPP
12#define PANZER_GATHER_NORMALS_IMPL_HPP
13
14#include "Teuchos_Assert.hpp"
15#include "Phalanx_DataLayout.hpp"
16
17#include "Panzer_PureBasis.hpp"
18#include "Kokkos_ViewFactory.hpp"
19
20#include "Intrepid2_Kernels.hpp"
21#include "Intrepid2_OrientationTools.hpp"
22
23#include "Teuchos_FancyOStream.hpp"
24
25template<typename EvalT,typename Traits>
28 const Teuchos::ParameterList& p)
29{
30 dof_name_ = (p.get< std::string >("DOF Name"));
31
32 if(p.isType< Teuchos::RCP<PureBasis> >("Basis"))
33 basis_ = p.get< Teuchos::RCP<PureBasis> >("Basis");
34 else
35 basis_ = p.get< Teuchos::RCP<const PureBasis> >("Basis");
36
37 pointRule_ = p.get<Teuchos::RCP<const PointRule> >("Point Rule");
38
39 Teuchos::RCP<PHX::DataLayout> basis_layout = basis_->functional;
40 Teuchos::RCP<PHX::DataLayout> vector_layout_vector = basis_->functional_grad;
41
42 // some sanity checks
43 TEUCHOS_ASSERT(basis_->isVectorBasis());
44
45 // setup the orientation field
46 std::string orientationFieldName = basis_->name() + " Orientation";
47 // setup all fields to be evaluated and constructed
48 pointValues_ = panzer::PointValues2<double> (pointRule_->getName()+"_",false);
49 pointValues_.setupArrays(pointRule_);
50
51 // the field manager will allocate all of these field
52 constJac_ = pointValues_.jac;
53 this->addDependentField(constJac_);
54
55 gatherFieldNormals_ = PHX::MDField<ScalarT,Cell,NODE,Dim>(dof_name_+"_Normals",vector_layout_vector);
56 this->addEvaluatedField(gatherFieldNormals_);
57
58 this->setName("Gather Normals");
59}
60
61// **********************************************************************
62template<typename EvalT,typename Traits>
66{
67 auto orientations = d.orientations_;
68 orientations_ = Kokkos::View<Intrepid2::Orientation*>("orientations_",orientations->size());
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);
73
74 this->utils.setFieldData(pointValues_.jac,fm);
75
76 const shards::CellTopology & parentCell = *basis_->getCellTopology();
77 int sideDim = parentCell.getDimension()-1;
78 sideParam_ = Intrepid2::RefSubcellParametrization<PHX::Device>::get(sideDim, parentCell.getKey());
79
80 int numFaces = gatherFieldNormals_.extent(1);
81 keys_ = Kokkos::View<unsigned int*>("parentCell.keys",numFaces);
82 auto keys_host = Kokkos::create_mirror_view(keys_);
83 for (int i=0; i < numFaces; ++i)
84 keys_host(i) = parentCell.getKey(sideDim,i);
85 Kokkos::deep_copy(keys_,keys_host);
86
87 // allocate space that is sized correctly for AD
88 int cellDim = parentCell.getDimension();
89 refEdges_ = Kokkos::createDynRankViewWithType<Kokkos::DynRankView<ScalarT,PHX::Device>>(gatherFieldNormals_.get_static_view(),"ref_edges", (*d.worksets_)[0].num_cells, sideDim, cellDim);
90 phyEdges_ = Kokkos::createDynRankViewWithType<Kokkos::DynRankView<ScalarT,PHX::Device>>(gatherFieldNormals_.get_static_view(),"phy_edges", (*d.worksets_)[0].num_cells, sideDim, cellDim);
91}
92
93// **********************************************************************
94template<typename EvalT,typename Traits>
96evaluateFields(typename Traits::EvalData workset)
97{
98
99 if(workset.num_cells<=0)
100 return;
101
102 int numFaces = gatherFieldNormals_.extent(1);
103 const auto worksetJacobians = pointValues_.jac.get_view();
104 const auto cell_local_ids = workset.getLocalCellIDs();
105 auto gatherFieldNormals = gatherFieldNormals_;
106 auto sideParam = sideParam_;
107 auto keys = keys_;
108 auto orientations = orientations_;
109 auto refEdges = refEdges_;
110 auto phyEdges = phyEdges_;
111
112 // Loop over workset faces and edge points
113 Kokkos::parallel_for("panzer::GatherNormals",workset.num_cells,KOKKOS_LAMBDA(const int c){
114 int faceOrts[6] = {};
115 orientations(cell_local_ids(c)).getFaceOrientation(faceOrts, numFaces);
116
117 for(int pt = 0; pt < numFaces; pt++) {
118 auto ortEdgeTan_U = Kokkos::subview(refEdges, c, 0, Kokkos::ALL());
119 auto ortEdgeTan_V = Kokkos::subview(refEdges, c, 1, Kokkos::ALL());
120
121 auto tmpRefEdges = Kokkos::subview(refEdges, c, Kokkos::ALL(), Kokkos::ALL());
122 Intrepid2::Impl::OrientationTools::getRefSubcellTangents(tmpRefEdges,
123 sideParam,
124 keys(pt),
125 pt,
126 faceOrts[pt]);
127
128 auto phyEdgeTan_U = Kokkos::subview(phyEdges, c, 0, Kokkos::ALL());
129 auto phyEdgeTan_V = Kokkos::subview(phyEdges, c, 1, Kokkos::ALL());
130 auto J = Kokkos::subview(worksetJacobians, c, pt, Kokkos::ALL(), Kokkos::ALL());
131
132 Intrepid2::Kernels::Serial::matvec_product(phyEdgeTan_U, J, ortEdgeTan_U);
133 Intrepid2::Kernels::Serial::matvec_product(phyEdgeTan_V, J, ortEdgeTan_V);
134
135 // take the cross product of the two vectors
136 gatherFieldNormals(c,pt,0) = (phyEdgeTan_U(1)*phyEdgeTan_V(2) - phyEdgeTan_U(2)*phyEdgeTan_V(1));
137 gatherFieldNormals(c,pt,1) = (phyEdgeTan_U(2)*phyEdgeTan_V(0) - phyEdgeTan_U(0)*phyEdgeTan_V(2));
138 gatherFieldNormals(c,pt,2) = (phyEdgeTan_U(0)*phyEdgeTan_V(1) - phyEdgeTan_U(1)*phyEdgeTan_V(0));
139 }
140 });
141
142}
143
144#endif
void evaluateFields(typename Traits::EvalData d)
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &vm)
Kokkos::View< const int *, PHX::Device > getLocalCellIDs() const
Get the local cell IDs for the workset.
int num_cells
DEPRECATED - use: numCells()
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_
Teuchos::RCP< const std::vector< Intrepid2::Orientation > > orientations_