Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_Dirichlet_Residual_FaceBasis_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_DIRICHLET_RESIDUAL_FACEBASIS_IMPL_HPP
12#define PANZER_DIRICHLET_RESIDUAL_FACEBASIS_IMPL_HPP
13
14#include <cstddef>
15#include <string>
16#include <vector>
17
18#include "Intrepid2_CellTools.hpp"
19
20#include "Kokkos_ViewFactory.hpp"
21
22namespace panzer {
23
24//**********************************************************************
25template<typename EvalT, typename Traits>
28 const Teuchos::ParameterList& p)
29{
30 std::string residual_name = p.get<std::string>("Residual Name");
31
32 basis = p.get<Teuchos::RCP<const panzer::PureBasis> >("Basis");
33 pointRule = p.get<Teuchos::RCP<const panzer::PointRule> >("Point Rule");
34
35 std::string field_name = p.get<std::string>("DOF Name");
36 std::string dof_name = field_name+"_"+pointRule->getName();
37 std::string value_name = p.get<std::string>("Value Name");
38
39 Teuchos::RCP<PHX::DataLayout> basis_layout = basis->functional;
40 Teuchos::RCP<PHX::DataLayout> vector_layout_dof = pointRule->dl_vector;
41 Teuchos::RCP<PHX::DataLayout> vector_layout_vector = basis->functional_grad;
42
43 // some sanity checks
44 TEUCHOS_ASSERT(basis->isVectorBasis());
45 TEUCHOS_ASSERT(basis_layout->extent(0)==vector_layout_dof->extent(0));
46 TEUCHOS_ASSERT(basis_layout->extent(1)==vector_layout_dof->extent(1));
47 TEUCHOS_ASSERT(basis->dimension()==vector_layout_dof->extent_int(2));
48 TEUCHOS_ASSERT(vector_layout_vector->extent(0)==vector_layout_dof->extent(0));
49 TEUCHOS_ASSERT(vector_layout_vector->extent(1)==vector_layout_dof->extent(1));
50 TEUCHOS_ASSERT(vector_layout_vector->extent(2)==vector_layout_dof->extent(2));
51
52 residual = PHX::MDField<ScalarT,Cell,BASIS>(residual_name, basis_layout);
53 dof = PHX::MDField<const ScalarT,Cell,Point,Dim>(dof_name, vector_layout_dof);
54 value = PHX::MDField<const ScalarT,Cell,Point,Dim>(value_name, vector_layout_vector);
55
56 // setup all fields to be evaluated and constructed
57 pointValues = PointValues2<double> (pointRule->getName()+"_",false);
58 pointValues.setupArrays(pointRule);
59
60 // the field manager will allocate all of these field
61 constJac_ = pointValues.jac;
62 this->addDependentField(constJac_);
63
64
65 this->addEvaluatedField(residual);
66 this->addDependentField(dof);
67 this->addDependentField(value);
68
69 std::string n = "Dirichlet Residual Face Basis Evaluator";
70 this->setName(n);
71}
72
73//**********************************************************************
74template<typename EvalT, typename Traits>
75void
78 typename Traits::SetupData sd,
80{
81 orientations = sd.orientations_;
82 this->utils.setFieldData(pointValues.jac,fm);
83 faceNormal = Kokkos::createDynRankView(residual.get_static_view(),"faceNormal",dof.extent(0),dof.extent(1),dof.extent(2));
84}
85
86//**********************************************************************
87template<typename EvalT, typename Traits>
88void
91 typename Traits::EvalData workset)
92{
93 // basic cell topology data
94 const shards::CellTopology & parentCell = *basis->getCellTopology();
95 const int cellDim = parentCell.getDimension();
96
97 // face only, subcellOrd does not include edge count
98 const int subcellDim = 2;
99 const int subcellOrd = this->wda(workset).subcell_index;
100
101 const int numFaces = parentCell.getSubcellCount(subcellDim);
102 const int numFaceDofs = dof.extent_int(1);
103
104 TEUCHOS_ASSERT(cellDim == dof.extent_int(2));
105
106 if(workset.num_cells<=0)
107 return;
108 else {
109 Intrepid2::CellTools<PHX::exec_space>::getPhysicalFaceNormals(faceNormal,
110 pointValues.jac.get_view(),
111 subcellOrd,
112 parentCell);
113 PHX::Device().fence();
114
115 const auto subcellBaseTopo = shards::CellTopology(parentCell.getBaseCellTopologyData(subcellDim, subcellOrd));
116 TEUCHOS_ASSERT(subcellBaseTopo.getBaseKey() == shards::Triangle<>::key ||
117 subcellBaseTopo.getBaseKey() == shards::Quadrilateral<>::key);
118
119 const WorksetDetails & details = workset;
120 const auto subcellVertexCount = static_cast<Intrepid2::ordinal_type>(subcellBaseTopo.getVertexCount());
121
122 // Copy orientations to device.
123 Kokkos::View<Intrepid2::Orientation*,PHX::Device> orientations_device(Kokkos::view_alloc("orientations_device",Kokkos::WithoutInitializing),orientations->size());
124 auto orientations_host = Kokkos::create_mirror_view(Kokkos::WithoutInitializing,orientations_device);
125 for (size_t i=0; i < orientations_host.extent(0); ++i)
126 orientations_host(i) = orientations->at(i);
127 Kokkos::deep_copy(orientations_device,orientations_host);
128
129 // Local temporaries for device lambda capture
130 const auto cell_local_ids_k = details.cell_local_ids_k;
131 auto residual_local = residual;
132 auto dof_local = dof;
133 auto value_local = value;
134 auto faceNormal_local = faceNormal;
135
136 Kokkos::parallel_for("panzer::DirichletRsidual_FaceBasis::evalauteFields",
137 workset.num_cells,
138 KOKKOS_LAMBDA(const index_t c) {
139 const auto ort = orientations_device(cell_local_ids_k[c]);
140 Intrepid2::ordinal_type faceOrts[6] = {};
141 ort.getFaceOrientation(faceOrts, numFaces);
142
143 // vertex count represent rotation count before it flips
144 const double ortVal = faceOrts[subcellOrd] < subcellVertexCount ? 1.0 : -1.0;
145 for(int b=0;b<numFaceDofs;b++) {
146 residual_local(c,b) = 0.0;
147 for(int d=0;d<cellDim;d++)
148 residual_local(c,b) += (dof_local(c,b,d)-value_local(c,b,d))*faceNormal_local(c,b,d);
149 residual_local(c,b) *= ortVal;
150 }
151 });
152 }
153}
154
155//**********************************************************************
156
157}
158
159#endif
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
int num_cells
DEPRECATED - use: numCells()
Teuchos::RCP< const std::vector< Intrepid2::Orientation > > orientations_