Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_Normals_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_NORMALS_IMPL_HPP
12#define PANZER_NORMALS_IMPL_HPP
13
14#include <algorithm>
17#include "Intrepid2_FunctionSpaceTools.hpp"
18#include "Intrepid2_CellTools.hpp"
19#include "Phalanx_Scratch_Utilities.hpp"
20
21namespace panzer {
22
23//**********************************************************************
24template<typename EvalT, typename Traits>
27 const Teuchos::ParameterList& p)
28 : normalize(true)
29{
30 // Read from parameters
31 const std::string name = p.get<std::string>("Name");
32 side_id = p.get<int>("Side ID");
33 Teuchos::RCP<panzer::IntegrationRule> quadRule
34 = p.get< Teuchos::RCP<panzer::IntegrationRule> >("IR");
35 if(p.isParameter("Normalize")) // set default
36 normalize = p.get<bool>("Normalize");
37
38 // grab information from quadrature rule
39 Teuchos::RCP<PHX::DataLayout> vector_dl = quadRule->dl_vector;
40 quad_order = quadRule->cubature_degree;
41
42 // build field, set as evaluated type
43 normals = PHX::MDField<ScalarT,Cell,Point,Dim>(name, vector_dl);
44 this->addEvaluatedField(normals);
45
46 std::string n = "Normals: " + name;
47 this->setName(n);
48}
49
50//**********************************************************************
51template<typename EvalT, typename Traits>
52void
55 typename Traits::SetupData sd,
57{
58 num_qp = normals.extent(1);
59 num_dim = normals.extent(2);
60
61 quad_index = panzer::getIntegrationRuleIndex(quad_order,(*sd.worksets_)[0], this->wda);
62}
63
64//**********************************************************************
65template<typename EvalT, typename Traits>
66void
69 typename Traits::EvalData workset)
70{
71 // ECC Fix: Get Physical Side Normals
72
73 if(workset.num_cells>0) {
74 Intrepid2::CellTools<PHX::exec_space>::getPhysicalSideNormals(normals.get_view(),
75 this->wda(workset).int_rules[quad_index]->jac.get_view(),
76 side_id, *this->wda(workset).int_rules[quad_index]->int_rule->topology);
77
78 if(normalize) {
79 // normalize vector: getPhysicalSideNormals does not
80 // return normalized vectors
81 auto local_normals = normals;
82 auto local_num_qp = num_qp;
83 auto local_num_dim = num_dim;
84
85 if (Sacado::IsADType<ScalarT>::value) {
86 scratch = PHX::View<ScalarT*>("normals:scratch",workset.num_cells,PHX::getFadSize(local_normals.get_static_view()));
87 }
88 else
89 scratch = PHX::View<ScalarT*>("normals:scratch",workset.num_cells);
90
91 auto norm = scratch;
92
93 Kokkos::parallel_for("normalize", workset.num_cells, KOKKOS_LAMBDA (index_t c) {
94 for(std::size_t q=0;q<local_num_qp;q++) {
95 norm(c) = 0.0;
96
97 // compute squared norm
98 for(std::size_t d=0;d<local_num_dim;d++)
99 norm(c) += local_normals(c,q,d)*local_normals(c,q,d);
100
101 // adjust for length of vector, now unit vectors
102 norm(c) = sqrt(norm(c));
103 for(std::size_t d=0;d<local_num_dim;d++)
104 local_normals(c,q,d) /= norm(c);
105 }
106 });
107 }
108 // else normals correspond to differential
109 }
110}
111
112//**********************************************************************
113
114}
115
116#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_