Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_ProjectToFaces_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_PROJECT_TO_FACES_IMPL_HPP
12 #define PANZER_PROJECT_TO_FACES_IMPL_HPP
13 
14 #include "Teuchos_Assert.hpp"
15 #include "Phalanx_DataLayout.hpp"
16 
17 #include "Intrepid2_Cubature.hpp"
18 #include "Intrepid2_DefaultCubatureFactory.hpp"
19 #include "Intrepid2_FunctionSpaceTools.hpp"
20 #include "Intrepid2_Kernels.hpp"
21 #include "Intrepid2_CellTools.hpp"
22 #include "Intrepid2_OrientationTools.hpp"
23 
25 #include "Panzer_PureBasis.hpp"
28 #include "Kokkos_ViewFactory.hpp"
29 
30 #include "Teuchos_FancyOStream.hpp"
31 
32 #include <cstring>
33 
34 template<typename EvalT,typename Traits>
37 {
38  dof_name_ = (p.get< std::string >("DOF Name"));
39 
40  if(p.isType< Teuchos::RCP<PureBasis> >("Basis"))
41  basis_ = p.get< Teuchos::RCP<PureBasis> >("Basis");
42  else
43  basis_ = p.get< Teuchos::RCP<const PureBasis> >("Basis");
44 
45  Teuchos::RCP<PHX::DataLayout> basis_layout = basis_->functional;
46  Teuchos::RCP<PHX::DataLayout> vector_layout = basis_->functional_grad;
47 
48  // some sanity checks
49  TEUCHOS_ASSERT(basis_->isVectorBasis());
50 
51  result_ = PHX::MDField<ScalarT,Cell,BASIS>(dof_name_,basis_layout);
52  this->addEvaluatedField(result_);
53 
54  normals_ = PHX::MDField<const ScalarT,Cell,BASIS,Dim>(dof_name_+"_Normals",vector_layout);
55  this->addDependentField(normals_);
56 
57  vector_values_ = PHX::MDField<const ScalarT,Cell,BASIS,Dim>(dof_name_+"_Vector",vector_layout);
58  this->addDependentField(vector_values_);
59 
60  this->setName("Project To Faces");
61 }
62 
63 // **********************************************************************
64 template<typename EvalT,typename Traits>
67  PHX::FieldManager<Traits>& /* fm */)
68 {
69  num_faces_ = result_.extent(1);
70  num_dim_ = vector_values_.extent(2);
71 
72  TEUCHOS_ASSERT(result_.extent(1) == normals_.extent(1));
73  TEUCHOS_ASSERT(vector_values_.extent(2) == normals_.extent(2));
74 }
75 
76 // **********************************************************************
77 template<typename EvalT,typename Traits>
80 {
81 
82  // Restricting HDiv field, multiplied by the normal to the faces, into HVol on the faces.
83  // This code assumes affine mapping and the projection into 1 quadrature point for each face,
84  // which is identified with the face. This makes sense only for low order bases, for which
85  // HVol is constant
86 
87  //TODO: make this work w/ high order basis
88  const int intDegree = basis_->order();
89  TEUCHOS_ASSERT(intDegree == 1);
90 
91  auto result = result_.get_static_view();
92  auto vector_values = vector_values_.get_static_view();
93  auto normals = normals_.get_static_view();
94  auto num_faces = num_faces_;
95  auto num_dim = num_dim_;
96 
97  auto policy = panzer::HP::inst().teamPolicy<ScalarT,PHX::exec_space>(workset.num_cells);
98  Kokkos::parallel_for("panzer::ProjectToFaces",policy,KOKKOS_LAMBDA(const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) {
99  const auto cell = team.league_rank();
100  Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,num_faces),[&] (const int p) {
101  result(cell,p) = ScalarT(0.0);
102  for (int dim = 0; dim < num_dim; ++dim)
103  result(cell,p) += vector_values(cell,p,dim) * normals(cell,p,dim);
104  });
105  });
106 }
107 
108 
109 #endif
int num_cells
DEPRECATED - use: numCells()
void evaluateFields(typename Traits::EvalData d)
T & get(const std::string &name, T def_value)
Kokkos::TeamPolicy< TeamPolicyProperties...> teamPolicy(const int &league_size)
Returns a TeamPolicy for hierarchic parallelism.
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &vm)
PHX::MDField< ScalarT, panzer::Cell, panzer::IP > result
A field that will be used to build up the result of the integral we&#39;re performing.
KOKKOS_FORCEINLINE_FUNCTION array_type get_static_view()
ProjectToFaces(const Teuchos::ParameterList &p)
bool isType(const std::string &name) const
static HP & inst()
Private ctor.
#define TEUCHOS_ASSERT(assertion_test)