Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 
52 namespace panzer {
53 
54 //**********************************************************************
55 template<typename EvalT, typename Traits>
58  const Teuchos::ParameterList& p)
59  : normalize(true)
60 {
61  // Read from parameters
62  const std::string name = p.get<std::string>("Name");
63  side_id = p.get<int>("Side ID");
66  if(p.isParameter("Normalize")) // set default
67  normalize = p.get<bool>("Normalize");
68 
69  // grab information from quadrature rule
70  Teuchos::RCP<PHX::DataLayout> vector_dl = quadRule->dl_vector;
71  quad_order = quadRule->cubature_degree;
72 
73  // build field, set as evaluated type
74  normals = PHX::MDField<ScalarT,Cell,Point,Dim>(name, vector_dl);
75  this->addEvaluatedField(normals);
76 
77  std::string n = "Normals: " + name;
78  this->setName(n);
79 }
80 
81 //**********************************************************************
82 template<typename EvalT, typename Traits>
83 void
86  typename Traits::SetupData sd,
87  PHX::FieldManager<Traits>& /* fm */)
88 {
89  num_qp = normals.extent(1);
90  num_dim = normals.extent(2);
91 
92  quad_index = panzer::getIntegrationRuleIndex(quad_order,(*sd.worksets_)[0], this->wda);
93 }
94 
95 //**********************************************************************
96 template<typename EvalT, typename Traits>
97 void
100  typename Traits::EvalData workset)
101 {
102  // ECC Fix: Get Physical Side Normals
103 
104  if(workset.num_cells>0) {
105  Intrepid2::CellTools<PHX::exec_space>::getPhysicalSideNormals(normals.get_view(),
106  this->wda(workset).int_rules[quad_index]->jac.get_view(),
107  side_id, *this->wda(workset).int_rules[quad_index]->int_rule->topology);
108 
109  if(normalize) {
110  // normalize vector: getPhysicalSideNormals does not
111  // return normalized vectors
112  for(index_t c=0;c<workset.num_cells;c++) {
113  for(std::size_t q=0;q<num_qp;q++) {
114  ScalarT norm = 0.0;
115 
116  // compute squared norm
117  for(std::size_t d=0;d<num_dim;d++)
118  norm += normals(c,q,d)*normals(c,q,d);
119 
120  // adjust for length of vector, now unit vectors
121  norm = sqrt(norm);
122  for(std::size_t d=0;d<num_dim;d++)
123  normals(c,q,d) /= norm;
124  }
125  }
126  }
127  // else normals correspond to differential
128  }
129 }
130 
131 //**********************************************************************
132 
133 }
134 
135 #endif
T & get(const std::string &name, T def_value)
Normals(const Teuchos::ParameterList &p)
std::vector< std::string >::size_type getIntegrationRuleIndex(int ir_degree, const panzer::Workset &workset, WorksetDetailsAccessor &wda)
bool isParameter(const std::string &name) const
PHX::MDField< ScalarT, Cell, Point, Dim > normals
typename EvalT::ScalarT ScalarT
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
void evaluateFields(typename Traits::EvalData d)
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_