Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_GatherTangents_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_TANGENTS_IMPL_HPP
12 #define PANZER_GATHER_TANGENTS_IMPL_HPP
13 
14 #include "Teuchos_Assert.hpp"
15 #include "Phalanx_DataLayout.hpp"
16 
17 #include "Intrepid2_Kernels.hpp"
18 #include "Intrepid2_OrientationTools.hpp"
19 
20 #include "Panzer_PureBasis.hpp"
21 #include "Kokkos_ViewFactory.hpp"
22 
23 #include "Teuchos_FancyOStream.hpp"
24 
25 template<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  gatherFieldTangents_ = PHX::MDField<ScalarT,Cell,NODE,Dim>(dof_name_+"_Tangents",vector_layout_vector);
56  this->addEvaluatedField(gatherFieldTangents_);
57 
58  this->setName("Gather Tangents");
59 }
60 
61 // **********************************************************************
62 template<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  const shards::CellTopology & parentCell = *basis_->getCellTopology();
76  const int edgeDim = 1;
77  edgeParam_ = Intrepid2::RefSubcellParametrization<PHX::Device>::get(edgeDim, parentCell.getKey());
78 
79  const int numEdges = gatherFieldTangents_.extent(1);
80  keys_ = Kokkos::View<unsigned int*>("parentCell.keys",numEdges);
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);
85 }
86 
87 // **********************************************************************
88 template<typename EvalT,typename Traits>
91 {
92 
93  if(workset.num_cells<=0)
94  return;
95  else {
96  const shards::CellTopology & parentCell = *basis_->getCellTopology();
97  const int cellDim = parentCell.getDimension();
98  const int numEdges = gatherFieldTangents_.extent(1);
99 
100  auto workspace_tmp = Kokkos::createDynRankView(gatherFieldTangents_.get_static_view(),"workspace", workset.num_cells,4, cellDim);
101  const auto worksetJacobians = pointValues_.jac.get_view();
102  const auto cell_local_ids = workset.getLocalCellIDs();
103  auto gatherFieldTangents = gatherFieldTangents_.get_static_view();
104  auto keys = keys_;
105  auto orientations = orientations_;
106  auto edgeParam = edgeParam_;
107 
108  // Loop over workset faces and edge points
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);
112 
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());
117 
118  Intrepid2::Impl::OrientationTools::getRefSubcellTangents(
119  workspace,
120  edgeParam,
121  keys(pt),
122  pt,
123  edgeOrts[pt]);
124 
125  auto J = Kokkos::subview(worksetJacobians, c, pt, Kokkos::ALL(), Kokkos::ALL());
126  Intrepid2::Kernels::Serial::matvec_product(phyEdgeTan, J, ortEdgeTan);
127 
128  }// for pt
129  });// for pCell
130 
131  }
132 
133 }
134 
135 #endif
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)