Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_ProjectToEdges_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_PROJECT_TO_EDGES_IMPL_HPP
44#define PANZER_PROJECT_TO_EDGES_IMPL_HPP
45
46#include "Teuchos_Assert.hpp"
47#include "Phalanx_DataLayout.hpp"
48
49#include "Intrepid2_Cubature.hpp"
50#include "Intrepid2_DefaultCubatureFactory.hpp"
51#include "Intrepid2_FunctionSpaceTools.hpp"
52#include "Intrepid2_OrientationTools.hpp"
53
54#include "Panzer_PureBasis.hpp"
57#include "Kokkos_ViewFactory.hpp"
58
59#include "Teuchos_FancyOStream.hpp"
60
61template<typename EvalT,typename Traits>
64 const Teuchos::ParameterList& p)
65{
66 dof_name_ = (p.get< std::string >("DOF Name"));
67
68 if(p.isType< Teuchos::RCP<PureBasis> >("Basis"))
69 basis_ = p.get< Teuchos::RCP<PureBasis> >("Basis");
70 else
71 basis_ = p.get< Teuchos::RCP<const PureBasis> >("Basis");
72
73 Teuchos::RCP<PHX::DataLayout> basis_layout = basis_->functional;
74 Teuchos::RCP<PHX::DataLayout> vector_layout = basis_->functional_grad;
75
76 // some sanity checks
77 TEUCHOS_ASSERT(basis_->isVectorBasis());
78
79 result_ = PHX::MDField<ScalarT,Cell,BASIS>(dof_name_,basis_layout);
80 this->addEvaluatedField(result_);
81
82 tangents_ = PHX::MDField<const ScalarT,Cell,BASIS,Dim>(dof_name_+"_Tangents",vector_layout);
83 this->addDependentField(tangents_);
84
85 vector_values_ = PHX::MDField<const ScalarT,Cell,BASIS,Dim>(dof_name_+"_Vector",vector_layout);
86 this->addDependentField(vector_values_);
87
88 this->setName("Project To Edges");
89}
90
91// **********************************************************************
92template<typename EvalT,typename Traits>
96{
97 num_edges_ = vector_values_.extent(1);
98 num_dim_ = vector_values_.extent(2);
99
100 TEUCHOS_ASSERT(vector_values_.extent(1) == tangents_.extent(1));
101 TEUCHOS_ASSERT(vector_values_.extent(2) == tangents_.extent(2));
102}
103
104// **********************************************************************
105template<typename EvalT,typename Traits>
107evaluateFields(typename Traits::EvalData workset)
108{
109 // Restricting HCurl field, multiplied by the tangent to the edge, into HVol on the edges.
110 // This code assumes affine mapping and the projection into 1 quadrature point for each edge,
111 // which is identified with the edge. This makes sense only for low order bases, for which
112 // HVol is constant
113
114 //TODO: make this work w/ high order basis
115 const int intDegree = basis_->order();
116 TEUCHOS_ASSERT(intDegree == 1);
117
118 auto result = result_.get_static_view();
119 auto vector_values = vector_values_.get_static_view();
120 auto tangents = tangents_.get_static_view();
121 auto num_edges = num_edges_;
122 auto num_dim = num_dim_;
123
124 auto policy = panzer::HP::inst().teamPolicy<ScalarT,PHX::exec_space>(workset.num_cells);
125 Kokkos::parallel_for("panzer::ProjectToEdges",policy,KOKKOS_LAMBDA(const Kokkos::TeamPolicy<PHX::exec_space>::member_type& team) {
126 const auto cell = team.league_rank();
127 Kokkos::parallel_for(Kokkos::TeamThreadRange(team,0,num_edges),[&] (const int p) {
128 result(cell,p) = ScalarT(0.0);
129 for (int dim = 0; dim < num_dim; ++dim)
130 result(cell,p) += vector_values(cell,p,dim) * tangents(cell,p,dim);
131 });
132 });
133}
134
135#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()