Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_Normals_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_NORMALS_IMPL_HPP
44#define PANZER_NORMALS_IMPL_HPP
45
46#include <algorithm>
49#include "Intrepid2_FunctionSpaceTools.hpp"
50#include "Intrepid2_CellTools.hpp"
51#include "Phalanx_Scratch_Utilities.hpp"
52
53namespace panzer {
54
55//**********************************************************************
56template<typename EvalT, typename Traits>
59 const Teuchos::ParameterList& p)
60 : normalize(true)
61{
62 // Read from parameters
63 const std::string name = p.get<std::string>("Name");
64 side_id = p.get<int>("Side ID");
65 Teuchos::RCP<panzer::IntegrationRule> quadRule
66 = p.get< Teuchos::RCP<panzer::IntegrationRule> >("IR");
67 if(p.isParameter("Normalize")) // set default
68 normalize = p.get<bool>("Normalize");
69
70 // grab information from quadrature rule
71 Teuchos::RCP<PHX::DataLayout> vector_dl = quadRule->dl_vector;
72 quad_order = quadRule->cubature_degree;
73
74 // build field, set as evaluated type
75 normals = PHX::MDField<ScalarT,Cell,Point,Dim>(name, vector_dl);
76 this->addEvaluatedField(normals);
77
78 std::string n = "Normals: " + name;
79 this->setName(n);
80}
81
82//**********************************************************************
83template<typename EvalT, typename Traits>
84void
87 typename Traits::SetupData sd,
89{
90 num_qp = normals.extent(1);
91 num_dim = normals.extent(2);
92
93 quad_index = panzer::getIntegrationRuleIndex(quad_order,(*sd.worksets_)[0], this->wda);
94}
95
96//**********************************************************************
97template<typename EvalT, typename Traits>
98void
101 typename Traits::EvalData workset)
102{
103 // ECC Fix: Get Physical Side Normals
104
105 if(workset.num_cells>0) {
106 Intrepid2::CellTools<PHX::exec_space>::getPhysicalSideNormals(normals.get_view(),
107 this->wda(workset).int_rules[quad_index]->jac.get_view(),
108 side_id, *this->wda(workset).int_rules[quad_index]->int_rule->topology);
109
110 if(normalize) {
111 // normalize vector: getPhysicalSideNormals does not
112 // return normalized vectors
113 auto local_normals = normals;
114 auto local_num_qp = num_qp;
115 auto local_num_dim = num_dim;
116
117 if (Sacado::IsADType<ScalarT>::value) {
118 scratch = PHX::View<ScalarT*>("normals:scratch",workset.num_cells,PHX::getFadSize(local_normals.get_static_view()));
119 }
120 else
121 scratch = PHX::View<ScalarT*>("normals:scratch",workset.num_cells);
122
123 auto norm = scratch;
124
125 Kokkos::parallel_for("normalize", workset.num_cells, KOKKOS_LAMBDA (index_t c) {
126 for(std::size_t q=0;q<local_num_qp;q++) {
127 norm(c) = 0.0;
128
129 // compute squared norm
130 for(std::size_t d=0;d<local_num_dim;d++)
131 norm(c) += local_normals(c,q,d)*local_normals(c,q,d);
132
133 // adjust for length of vector, now unit vectors
134 norm(c) = sqrt(norm(c));
135 for(std::size_t d=0;d<local_num_dim;d++)
136 local_normals(c,q,d) /= norm(c);
137 }
138 });
139 }
140 // else normals correspond to differential
141 }
142}
143
144//**********************************************************************
145
146}
147
148#endif
Normals(const Teuchos::ParameterList &p)
PHX::MDField< ScalarT, Cell, Point, Dim > normals
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
void evaluateFields(typename Traits::EvalData d)
int num_cells
DEPRECATED - use: numCells()
std::vector< int >::size_type getIntegrationRuleIndex(int ir_degree, const panzer::Workset &workset, WorksetDetailsAccessor &wda)
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_