Panzer  Version of the Day
Panzer_GatherTangents_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_GATHER_TANGENTS_IMPL_HPP
44 #define PANZER_GATHER_TANGENTS_IMPL_HPP
45 
46 #include "Teuchos_Assert.hpp"
47 #include "Phalanx_DataLayout.hpp"
48 
49 #include "Intrepid2_Kernels.hpp"
50 #include "Intrepid2_OrientationTools.hpp"
51 
52 #include "Panzer_PureBasis.hpp"
53 #include "Kokkos_ViewFactory.hpp"
54 
55 #include "Teuchos_FancyOStream.hpp"
56 
57 template<typename EvalT,typename Traits>
60  const Teuchos::ParameterList& p)
61 {
62  dof_name = (p.get< std::string >("DOF Name"));
63 
64  if(p.isType< Teuchos::RCP<PureBasis> >("Basis"))
65  basis = p.get< Teuchos::RCP<PureBasis> >("Basis");
66  else
67  basis = p.get< Teuchos::RCP<const PureBasis> >("Basis");
68 
69  pointRule = p.get<Teuchos::RCP<const PointRule> >("Point Rule");
70 
71  Teuchos::RCP<PHX::DataLayout> basis_layout = basis->functional;
72  Teuchos::RCP<PHX::DataLayout> vector_layout_vector = basis->functional_grad;
73 
74  // some sanity checks
75  TEUCHOS_ASSERT(basis->isVectorBasis());
76 
77  // setup the orientation field
78  std::string orientationFieldName = basis->name() + " Orientation";
79  // setup all fields to be evaluated and constructed
80  pointValues = panzer::PointValues2<double> (pointRule->getName()+"_",false);
81  pointValues.setupArrays(pointRule);
82 
83  // the field manager will allocate all of these field
84  constJac_ = pointValues.jac;
85  this->addDependentField(constJac_);
86 
87  gatherFieldTangents = PHX::MDField<ScalarT,Cell,NODE,Dim>(dof_name+"_Tangents",vector_layout_vector);
88  this->addEvaluatedField(gatherFieldTangents);
89 
90  this->setName("Gather Tangents");
91 }
92 
93 // **********************************************************************
94 template<typename EvalT,typename Traits>
98 {
99  orientations = d.orientations_;
100  this->utils.setFieldData(pointValues.jac,fm);
101  const shards::CellTopology & parentCell = *basis->getCellTopology();
102  int edgeDim = 1;
103  edgeParam = Intrepid2::RefSubcellParametrization<PHX::Device>::get(edgeDim, parentCell.getKey());
104 }
105 
106 // **********************************************************************
107 template<typename EvalT,typename Traits>
110 {
111 
112  if(workset.num_cells<=0)
113  return;
114  else {
115  const shards::CellTopology & parentCell = *basis->getCellTopology();
116  int cellDim = parentCell.getDimension();
117  int edgeDim = 1;
118  int numEdges = gatherFieldTangents.extent(1);
119 
120  auto workspace = Kokkos::createDynRankView(gatherFieldTangents.get_static_view(),"workspace", 4, cellDim);
121 
122  const WorksetDetails & details = workset;
123 
124  const auto worksetJacobians = pointValues.jac.get_view();
125 
126  // Loop over workset faces and edge points
127  for(index_t c=0;c<workset.num_cells;c++) {
128 
129  int edgeOrts[12] = {};
130  orientations->at(details.cell_local_ids[c]).getEdgeOrientation(edgeOrts, numEdges);
131 
132  for(int pt = 0; pt < numEdges; pt++) {
133  auto phyEdgeTan = Kokkos::subview(gatherFieldTangents.get_static_view(), c, pt, Kokkos::ALL());
134  auto ortEdgeTan = Kokkos::subview(workspace, 0, Kokkos::ALL());
135 
136 
137  Intrepid2::Impl::OrientationTools::getRefSubcellTangents(
138  workspace,
139  edgeParam,
140  parentCell.getKey(edgeDim,pt),
141  pt,
142  edgeOrts[pt]);
143 
144  auto J = Kokkos::subview(worksetJacobians, c, pt, Kokkos::ALL(), Kokkos::ALL());
145  Intrepid2::Kernels::Serial::matvec_product(phyEdgeTan, J, ortEdgeTan);
146 
147  }// for pt
148  }// for pCell
149 
150  }
151 
152 }
153 
154 #endif
Teuchos::RCP< const std::vector< Intrepid2::Orientation > > orientations_
int num_cells
DEPRECATED - use: numCells()
void setupArrays(const Teuchos::RCP< const panzer::PointRule > &pr)
Sizes/allocates memory for arrays.
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.
void evaluateFields(typename Traits::EvalData d)
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &vm)