Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_ProjectToEdges_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_EDGES_IMPL_HPP
12#define PANZER_PROJECT_TO_EDGES_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_OrientationTools.hpp"
21
22#include "Panzer_PureBasis.hpp"
25#include "Kokkos_ViewFactory.hpp"
26
27#include "Teuchos_FancyOStream.hpp"
28
29template<typename EvalT,typename Traits>
32 const Teuchos::ParameterList& p)
33{
34 dof_name_ = (p.get< std::string >("DOF Name"));
35
36 if(p.isType< Teuchos::RCP<PureBasis> >("Basis"))
37 basis_ = p.get< Teuchos::RCP<PureBasis> >("Basis");
38 else
39 basis_ = p.get< Teuchos::RCP<const PureBasis> >("Basis");
40
41 Teuchos::RCP<PHX::DataLayout> basis_layout = basis_->functional;
42 Teuchos::RCP<PHX::DataLayout> vector_layout = basis_->functional_grad;
43
44 // some sanity checks
45 TEUCHOS_ASSERT(basis_->isVectorBasis());
46
47 result_ = PHX::MDField<ScalarT,Cell,BASIS>(dof_name_,basis_layout);
48 this->addEvaluatedField(result_);
49
50 tangents_ = PHX::MDField<const ScalarT,Cell,BASIS,Dim>(dof_name_+"_Tangents",vector_layout);
51 this->addDependentField(tangents_);
52
53 vector_values_ = PHX::MDField<const ScalarT,Cell,BASIS,Dim>(dof_name_+"_Vector",vector_layout);
54 this->addDependentField(vector_values_);
55
56 this->setName("Project To Edges");
57}
58
59// **********************************************************************
60template<typename EvalT,typename Traits>
64{
65 num_edges_ = vector_values_.extent(1);
66 num_dim_ = vector_values_.extent(2);
67
68 TEUCHOS_ASSERT(vector_values_.extent(1) == tangents_.extent(1));
69 TEUCHOS_ASSERT(vector_values_.extent(2) == tangents_.extent(2));
70}
71
72// **********************************************************************
73template<typename EvalT,typename Traits>
75evaluateFields(typename Traits::EvalData workset)
76{
77 // Restricting HCurl field, multiplied by the tangent to the edge, into HVol on the edges.
78 // This code assumes affine mapping and the projection into 1 quadrature point for each edge,
79 // which is identified with the edge. This makes sense only for low order bases, for which
80 // HVol is constant
81
82 //TODO: make this work w/ high order basis
83 const int intDegree = basis_->order();
84 TEUCHOS_ASSERT(intDegree == 1);
85
86 auto result = result_.get_static_view();
87 auto vector_values = vector_values_.get_static_view();
88 auto tangents = tangents_.get_static_view();
89 auto num_edges = num_edges_;
90 auto num_dim = num_dim_;
91
92 auto policy = panzer::HP::inst().teamPolicy<ScalarT,PHX::exec_space>(workset.num_cells);
93 Kokkos::parallel_for("panzer::ProjectToEdges",policy,KOKKOS_LAMBDA(const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) {
94 const auto cell = team.league_rank();
95 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,num_edges),[&] (const int p) {
96 result(cell,p) = ScalarT(0.0);
97 for (int dim = 0; dim < num_dim; ++dim)
98 result(cell,p) += vector_values(cell,p,dim) * tangents(cell,p,dim);
99 });
100 });
101}
102
103#endif
PHX::MDField< ScalarT, panzer::Cell, panzer::IP > result
A field that will be used to build up the result of the integral we're performing.
static HP & inst()
Private ctor.
Kokkos::TeamPolicy< TeamPolicyProperties... > teamPolicy(const int &league_size)
Returns a TeamPolicy for hierarchic parallelism.
ProjectToEdges(const Teuchos::ParameterList &p)
void evaluateFields(typename Traits::EvalData d)
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &vm)
int num_cells
DEPRECATED - use: numCells()