Panzer  Version of the Day
Panzer_CrossProduct_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_EVALUATOR_CrossProduct_IMPL_HPP
44 #define PANZER_EVALUATOR_CrossProduct_IMPL_HPP
45 
46 #include <string>
47 
48 #include "Panzer_PointRule.hpp"
50 
51 #include "Teuchos_RCP.hpp"
52 
53 namespace panzer {
54 
55 //**********************************************************************
56 PHX_EVALUATOR_CTOR(CrossProduct,p)
57 {
58  std::string result_name = p.get<std::string>("Result Name");
59  std::string vec_a_name = p.get<std::string>("Vector A Name");
60  std::string vec_b_name = p.get<std::string>("Vector B Name");
61 
62  const Teuchos::RCP<const panzer::PointRule> pr =
63  p.get< Teuchos::RCP<const panzer::PointRule> >("Point Rule");
64 
65  // use a scalar field only if dimension is 2D
66  useScalarField = (pr->spatial_dimension==2);
67 
68  if(!useScalarField)
69  vec_a_cross_vec_b = PHX::MDField<ScalarT>(result_name, pr->dl_vector);
70  else
71  vec_a_cross_vec_b = PHX::MDField<ScalarT>(result_name, pr->dl_scalar);
72 
73  vec_a = PHX::MDField<ScalarT>(vec_a_name, pr->dl_vector);
74  vec_b = PHX::MDField<ScalarT>(vec_b_name, pr->dl_vector);
75 
76  this->addEvaluatedField(vec_a_cross_vec_b);
77  this->addDependentField(vec_a);
78  this->addDependentField(vec_b);
79 
80  std::string n = "CrossProduct: " + result_name + " = " + vec_a_name + " . " + vec_b_name;
81  this->setName(n);
82 }
83 
84 //**********************************************************************
85 PHX_POST_REGISTRATION_SETUP(CrossProduct,sd,fm)
86 {
87  this->utils.setFieldData(vec_a_cross_vec_b,fm);
88  this->utils.setFieldData(vec_a,fm);
89  this->utils.setFieldData(vec_b,fm);
90 
91  num_pts = vec_a.dimension(1);
92  num_dim = vec_a.dimension(2);
93 
94  TEUCHOS_ASSERT(vec_a.dimension(1) == vec_b.dimension(1));
95  TEUCHOS_ASSERT(vec_a.dimension(2) == vec_b.dimension(2));
96 }
97 
98 //**********************************************************************
99 PHX_EVALUATE_FIELDS(CrossProduct,workset)
100 {
101  if(useScalarField) {
102  for (index_t cell = 0; cell < workset.num_cells; ++cell) {
103  for (int p = 0; p < num_pts; ++p) {
104  vec_a_cross_vec_b(cell,p) = vec_a(cell,p,0)*vec_b(cell,p,1)-vec_a(cell,p,1)*vec_b(cell,p,0);
105  }
106  }
107  }
108  else {
109  for (index_t cell = 0; cell < workset.num_cells; ++cell) {
110  for (int p = 0; p < num_pts; ++p) {
111  vec_a_cross_vec_b(cell,p,0) = vec_a(cell,p,1)*vec_b(cell,p,2)-vec_a(cell,p,2)*vec_b(cell,p,1);
112  vec_a_cross_vec_b(cell,p,1) = -(vec_a(cell,p,0)*vec_b(cell,p,2)-vec_a(cell,p,2)*vec_b(cell,p,0));
113  vec_a_cross_vec_b(cell,p,2) = vec_a(cell,p,0)*vec_b(cell,p,1)-vec_a(cell,p,1)*vec_b(cell,p,0);
114  }
115  }
116  }
117 }
118 
119 //**********************************************************************
120 
121 }
122 
123 #endif
PHX::MDField< ScalarT > vec_a_cross_vec_b
Evaluates cross product at a set of points.
PHX::MDField< ScalarT > vec_b
PHX_EVALUATOR_CTOR(BasisValues_Evaluator, p)
PHX::MDField< ScalarT > vec_a
PHX_EVALUATE_FIELDS(BasisValues_Evaluator, workset)
PHX_POST_REGISTRATION_SETUP(BasisValues_Evaluator, sd, fm)