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//
4// Panzer: A partial differential equation assembly
5// engine for strongly coupled complex multiphysics systems
6// Copyright (2011) Sandia Corporation
7//
8// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9// the U.S. Government retains certain rights in this software.
10//
11// Redistribution and use in source and binary forms, with or without
12// modification, are permitted provided that the following conditions are
13// met:
14//
15// 1. Redistributions of source code must retain the above copyright
16// notice, this list of conditions and the following disclaimer.
17//
18// 2. Redistributions in binary form must reproduce the above copyright
19// notice, this list of conditions and the following disclaimer in the
20// documentation and/or other materials provided with the distribution.
21//
22// 3. Neither the name of the Corporation nor the names of the
23// contributors may be used to endorse or promote products derived from
24// this software without specific prior written permission.
25//
26// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37//
38// Questions? Contact Roger P. Pawlowski (rppawlo@sandia.gov) and
39// Eric C. Cyr (eccyr@sandia.gov)
40// ***********************************************************************
41// @HEADER
42
43#ifndef PANZER_DIRICHLET_RESIDUAL_FACEBASIS_IMPL_HPP
44#define PANZER_DIRICHLET_RESIDUAL_FACEBASIS_IMPL_HPP
45
46#include <cstddef>
47#include <string>
48#include <vector>
49
50#include "Intrepid2_CellTools.hpp"
51
52#include "Kokkos_ViewFactory.hpp"
53
54namespace panzer {
55
56//**********************************************************************
57template<typename EvalT, typename Traits>
60 const Teuchos::ParameterList& p)
61{
62 std::string residual_name = p.get<std::string>("Residual Name");
63
64 basis = p.get<Teuchos::RCP<const panzer::PureBasis> >("Basis");
65 pointRule = p.get<Teuchos::RCP<const panzer::PointRule> >("Point Rule");
66
67 std::string field_name = p.get<std::string>("DOF Name");
68 std::string dof_name = field_name+"_"+pointRule->getName();
69 std::string value_name = p.get<std::string>("Value Name");
70
71 Teuchos::RCP<PHX::DataLayout> basis_layout = basis->functional;
72 Teuchos::RCP<PHX::DataLayout> vector_layout_dof = pointRule->dl_vector;
73 Teuchos::RCP<PHX::DataLayout> vector_layout_vector = basis->functional_grad;
74
75 // some sanity checks
76 TEUCHOS_ASSERT(basis->isVectorBasis());
77 TEUCHOS_ASSERT(basis_layout->extent(0)==vector_layout_dof->extent(0));
78 TEUCHOS_ASSERT(basis_layout->extent(1)==vector_layout_dof->extent(1));
79 TEUCHOS_ASSERT(basis->dimension()==vector_layout_dof->extent_int(2));
80 TEUCHOS_ASSERT(vector_layout_vector->extent(0)==vector_layout_dof->extent(0));
81 TEUCHOS_ASSERT(vector_layout_vector->extent(1)==vector_layout_dof->extent(1));
82 TEUCHOS_ASSERT(vector_layout_vector->extent(2)==vector_layout_dof->extent(2));
83
84 residual = PHX::MDField<ScalarT,Cell,BASIS>(residual_name, basis_layout);
85 dof = PHX::MDField<const ScalarT,Cell,Point,Dim>(dof_name, vector_layout_dof);
86 value = PHX::MDField<const ScalarT,Cell,Point,Dim>(value_name, vector_layout_vector);
87
88 // setup all fields to be evaluated and constructed
89 pointValues = PointValues2<double> (pointRule->getName()+"_",false);
90 pointValues.setupArrays(pointRule);
91
92 // the field manager will allocate all of these field
93 constJac_ = pointValues.jac;
94 this->addDependentField(constJac_);
95
96
97 this->addEvaluatedField(residual);
98 this->addDependentField(dof);
99 this->addDependentField(value);
100
101 std::string n = "Dirichlet Residual Face Basis Evaluator";
102 this->setName(n);
103}
104
105//**********************************************************************
106template<typename EvalT, typename Traits>
107void
110 typename Traits::SetupData sd,
112{
113 orientations = sd.orientations_;
114 this->utils.setFieldData(pointValues.jac,fm);
115 faceNormal = Kokkos::createDynRankView(residual.get_static_view(),"faceNormal",dof.extent(0),dof.extent(1),dof.extent(2));
116}
117
118//**********************************************************************
119template<typename EvalT, typename Traits>
120void
123 typename Traits::EvalData workset)
124{
125 // basic cell topology data
126 const shards::CellTopology & parentCell = *basis->getCellTopology();
127 const int cellDim = parentCell.getDimension();
128
129 // face only, subcellOrd does not include edge count
130 const int subcellDim = 2;
131 const int subcellOrd = this->wda(workset).subcell_index;
132
133 const int numFaces = parentCell.getSubcellCount(subcellDim);
134 const int numFaceDofs = dof.extent_int(1);
135
136 TEUCHOS_ASSERT(cellDim == dof.extent_int(2));
137
138 if(workset.num_cells<=0)
139 return;
140 else {
141 Intrepid2::CellTools<PHX::exec_space>::getPhysicalFaceNormals(faceNormal,
142 pointValues.jac.get_view(),
143 subcellOrd,
144 parentCell);
145 PHX::Device().fence();
146
147 const auto subcellTopo = shards::CellTopology(parentCell.getBaseCellTopologyData(subcellDim, subcellOrd));
148 TEUCHOS_ASSERT(subcellTopo.getBaseKey() == shards::Triangle<>::key ||
149 subcellTopo.getBaseKey() == shards::Quadrilateral<>::key);
150
151 const WorksetDetails & details = workset;
152 const auto subcellVertexCount = static_cast<Intrepid2::ordinal_type>(subcellTopo.getVertexCount());
153
154 // Copy orientations to device.
155 Kokkos::View<Intrepid2::Orientation*,PHX::Device> orientations_device(Kokkos::view_alloc("orientations_device",Kokkos::WithoutInitializing),orientations->size());
156 auto orientations_host = Kokkos::create_mirror_view(Kokkos::WithoutInitializing,orientations_device);
157 for (size_t i=0; i < orientations_host.extent(0); ++i)
158 orientations_host(i) = orientations->at(i);
159 Kokkos::deep_copy(orientations_device,orientations_host);
160
161 // Local temporaries for device lambda capture
162 const auto cell_local_ids_k = details.cell_local_ids_k;
163 auto residual_local = residual;
164 auto dof_local = dof;
165 auto value_local = value;
166 auto faceNormal_local = faceNormal;
167
168 Kokkos::parallel_for("panzer::DirichletRsidual_FaceBasis::evalauteFields",
169 workset.num_cells,
170 KOKKOS_LAMBDA(const index_t c) {
171 const auto ort = orientations_device(cell_local_ids_k[c]);
172 Intrepid2::ordinal_type faceOrts[6] = {};
173 ort.getFaceOrientation(faceOrts, numFaces);
174
175 // vertex count represent rotation count before it flips
176 const double ortVal = faceOrts[subcellOrd] < subcellVertexCount ? 1.0 : -1.0;
177 for(int b=0;b<numFaceDofs;b++) {
178 residual_local(c,b) = 0.0;
179 for(int d=0;d<cellDim;d++)
180 residual_local(c,b) += (dof_local(c,b,d)-value_local(c,b,d))*faceNormal_local(c,b,d);
181 residual_local(c,b) *= ortVal;
182 }
183 });
184 }
185}
186
187//**********************************************************************
188
189}
190
191#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_