Panzer  Version of the Day
Panzer_Integrator_GradBasisCrossVector_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_GRADBASISCROSSVECTOR_IMPL_HPP
44 #define PANZER_EVALUATOR_GRADBASISCROSSVECTOR_IMPL_HPP
45 
46 #include "Intrepid2_FunctionSpaceTools.hpp"
48 #include "Panzer_BasisIRLayout.hpp"
50 #include "Kokkos_ViewFactory.hpp"
51 
52 namespace panzer {
53 
54 //**********************************************************************
55 template<typename EvalT, typename Traits>
58  const Teuchos::ParameterList& p):
59  _num_basis_nodes(0),
60  _num_quadrature_points(0),
61  _basis_index(-1)
62 {
63 
64  // TODO: Get all of the panzer modules to read in const integration rules and basisIRlayouts
67  Teuchos::RCP<PHX::DataLayout> vector_dl = p.get< Teuchos::RCP<PHX::DataLayout> >("Data Layout Vector");
68 
69  Teuchos::RCP<const PureBasis> basis = basis_layout->getBasis();
70  _basis_name = basis_layout->name();
71 
72  // Verify that this basis supports the gradient operation
73  TEUCHOS_TEST_FOR_EXCEPTION(!basis->supportsGrad(),std::logic_error,"Integrator_GradBasisCrossVector: Basis of type \"" << basis->name() << "\" does not support GRAD");
74 
75  // Setup residuals
76  _residuals.clear();
77  if (p.isType<Teuchos::RCP<const std::vector<std::string> > >("Residual Names")){
78  const std::vector<std::string> & names = *(p.get<Teuchos::RCP<const std::vector<std::string> > >("Residual Names"));
79  for(const auto & name : names){
80  PHX::MDField<ScalarT,Cell,BASIS> res(name, basis_layout->functional);
81  _residuals.push_back(res);
82  }
83  }
84  _num_dims = _residuals.size();
85 
86  // Currently we only allow for vectors of length 3
87  TEUCHOS_TEST_FOR_EXCEPTION(_num_dims != vector_dl->extent(2),std::logic_error,"Vector must be same length as number of residuals.");
88 
89  // Setup vector
90  _vector = PHX::MDField<const ScalarT,Cell,IP,Dim>(p.get<std::string>("Vector Name"), vector_dl);
91 
92  // TODO: I'm assuming that the gradient is held in the same number of dimensions as ir->dl_vector (dl_vector is initialized from mesh dimensions)
93  _num_grad_dims = ir->dl_vector->extent(2);
94 
95  // Vector must be at least the same dimension as the GRAD operation
96  TEUCHOS_TEST_FOR_EXCEPTION(_num_grad_dims > _num_dims,std::logic_error,"Vector size must have at least as many components as there are dimensions in the mesh.");
97 
98  // Read the scalar multiplier
99  _multiplier = ScalarT(p.get<double>("Multiplier"));
100 
101  // Setup field multipliers
102  if (p.isType<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers")){
103  const std::vector<std::string> & field_multiplier_names = *(p.get<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers"));
104  for (const std::string & name : field_multiplier_names){
105  _field_multipliers.push_back(PHX::MDField<const ScalarT,Cell,IP>(name, ir->dl_scalar));
106  }
107  }
108 
109  // Register evaluated fields
110  for (auto & residual : _residuals){
111  this->addEvaluatedField(residual);
112  }
113 
114  // Register dependent fields
115  this->addDependentField(_vector);
116  for (auto & field : _field_multipliers){
117  this->addDependentField(field);
118  }
119 
120  // Name the module
121  this->setName("Integrator_GradBasisCrossVector: " + _vector.fieldTag().name());
122 }
123 
124 //**********************************************************************
125 template<typename EvalT, typename Traits>
126 void
129  typename Traits::SetupData sd,
130  PHX::FieldManager<Traits>& /* fm */)
131 {
132  _num_basis_nodes = _residuals[0].extent(1);
133  _num_quadrature_points = _vector.extent(1);
134 
135  _basis_index = panzer::getBasisIndex(_basis_name, (*sd.worksets_)[0], this->wda);
136 
137  // TODO: figure out a clean way of cloning _vector
138  _tmp = Kokkos::createDynRankView(_residuals[0].get_static_view(),"tmp",_vector.extent(0), _num_quadrature_points, _num_dims);
139 }
140 
141 //**********************************************************************
142 template<typename EvalT, typename Traits>
143 void
146  typename Traits::EvalData workset)
147 {
148 
149  // Zero the residuals
150  for (int i(0); i < static_cast<int>(_num_dims); ++i)
151  Kokkos::deep_copy(_residuals[i].get_static_view(), ScalarT(0.0));
152 
153  // The curl operation will only do something if _num_dims == 3
154  if(_num_dims != 3){
155  return;
156  }
157 
158  // do a scaled copy to initialize _tmp
159  for (int i=0; i < _vector.extent_int(0); ++i)
160  for (int j=0; j < _vector.extent_int(1); ++j)
161  for (int k=0; k < _vector.extent_int(2); ++k)
162  _tmp(i,j,k) = _multiplier * _vector(i,j,k);
163 
164  // Apply the field multipliers
165  for(auto & field_data : _field_multipliers){
166  for (index_t cell = 0; cell < workset.num_cells; ++cell) {
167  for (std::size_t qp = 0; qp < _num_quadrature_points; ++qp) {
168  const ScalarT & tmpVar = field_data(cell,qp);
169  for (std::size_t dim = 0; dim < _num_dims; ++dim){
170  _tmp(cell,qp,dim) *= tmpVar;
171  }
172  }
173  }
174  }
175 
176  // const Kokkos::DynRankView<double,PHX::Device> & weighted_grad_basis = this->wda(workset).bases[basis_index]->weighted_grad_basis;
177  const BasisValues2<double> & bv = *this->wda(workset).bases[_basis_index];
178 
179  ScalarT r_0,r_1,r_2;
180 
181  // this part gets annoying since dimensionality of the weighted_grad_basis is defined by the mesh
182  // There might be a fix, but for now we will just make it ugly.
183  if(_num_grad_dims == 1){
184 
185  // Curl only exists if the vector is three dimensional
186  if(_num_dims == 3){
187  for (index_t cell = 0; cell < workset.num_cells; ++cell){
188  for (std::size_t basis = 0; basis < _num_basis_nodes; ++basis){
189  r_1=r_2=0;
190  for (std::size_t qp = 0; qp < _num_quadrature_points; ++qp){
191  r_1 += _tmp(cell,qp,2)*bv.weighted_grad_basis(cell,basis,qp,0);
192  r_2 += - _tmp(cell,qp,1)*bv.weighted_grad_basis(cell,basis,qp,0);
193  }
194  _residuals[1](cell,basis) += r_1;
195  _residuals[2](cell,basis) += r_2;
196  }
197  }
198  }
199 
200  } else if (_num_grad_dims == 2){
201 
202  // Curl only exists if the vector is three dimensional
203  if(_num_dims == 3){
204  for (index_t cell = 0; cell < workset.num_cells; ++cell){
205  for (std::size_t basis = 0; basis < _num_basis_nodes; ++basis){
206  r_0=r_1=r_2=0;
207  for (std::size_t qp = 0; qp < _num_quadrature_points; ++qp){
208  r_0 += - _tmp(cell,qp,2)*bv.weighted_grad_basis(cell,basis,qp,1);
209  r_1 += _tmp(cell,qp,2)*bv.weighted_grad_basis(cell,basis,qp,0);
210  r_2 += _tmp(cell,qp,0)*bv.weighted_grad_basis(cell,basis,qp,1) - _tmp(cell,qp,1)*bv.weighted_grad_basis(cell,basis,qp,0);
211  }
212  _residuals[0](cell,basis) += r_0;
213  _residuals[1](cell,basis) += r_1;
214  _residuals[2](cell,basis) += r_2;
215  }
216  }
217  }
218 
219  } else if (_num_grad_dims == 3){
220 
221  // Will have thrown an error if _num_grad_dims != _num_dims != 3
222  for (index_t cell = 0; cell < workset.num_cells; ++cell){
223  for (std::size_t basis = 0; basis < _num_basis_nodes; ++basis){
224  r_0=r_1=r_2=0;
225  for (std::size_t qp = 0; qp < _num_quadrature_points; ++qp){
226  r_0 += _tmp(cell,qp,1)*bv.weighted_grad_basis(cell,basis,qp,2) - _tmp(cell,qp,2)*bv.weighted_grad_basis(cell,basis,qp,1);
227  r_1 += _tmp(cell,qp,2)*bv.weighted_grad_basis(cell,basis,qp,0) - _tmp(cell,qp,0)*bv.weighted_grad_basis(cell,basis,qp,2);
228  r_2 += _tmp(cell,qp,0)*bv.weighted_grad_basis(cell,basis,qp,1) - _tmp(cell,qp,1)*bv.weighted_grad_basis(cell,basis,qp,0);
229  }
230  _residuals[0](cell,basis) += r_0;
231  _residuals[1](cell,basis) += r_1;
232  _residuals[2](cell,basis) += r_2;
233  }
234  }
235  }
236 
237 }
238 
239 //**********************************************************************
240 
241 }
242 
243 #endif
panzer::Integrator_GradBasisCrossVector::_num_grad_dims
std::size_t _num_grad_dims
Definition: Panzer_Integrator_GradBasisCrossVector_decl.hpp:95
panzer::Integrator_GradBasisCrossVector::Integrator_GradBasisCrossVector
Integrator_GradBasisCrossVector(const Teuchos::ParameterList &p)
Definition: Panzer_Integrator_GradBasisCrossVector_impl.hpp:57
panzer::Traits::SD::worksets_
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_
Definition: Panzer_Traits.hpp:124
panzer::Integrator_GradBasisCrossVector::evaluateFields
void evaluateFields(typename Traits::EvalData d)
Definition: Panzer_Integrator_GradBasisCrossVector_impl.hpp:145
panzer::Traits::SD
Definition: Panzer_Traits.hpp:123
panzer::PointRule::dl_scalar
Teuchos::RCP< PHX::DataLayout > dl_scalar
Data layout for scalar fields.
Definition: Panzer_PointRule.hpp:112
Panzer_BasisIRLayout.hpp
Teuchos::ParameterList::get
T & get(const std::string &name, T def_value)
panzer::Integrator_GradBasisCrossVector::_vector
PHX::MDField< const ScalarT, Cell, IP, Dim > _vector
Definition: Panzer_Integrator_GradBasisCrossVector_decl.hpp:82
panzer::BasisValues2::weighted_grad_basis
Array_CellBasisIPDim weighted_grad_basis
Definition: Panzer_BasisValues2.hpp:142
panzer::Workset
Definition: Panzer_Workset.hpp:186
Panzer_Workset_Utilities.hpp
_num_dims
int _num_dims
Definition: Panzer_IntegrationValues2.cpp:437
panzer::Integrator_GradBasisCrossVector::ScalarT
typename EvalT::ScalarT ScalarT
Definition: Panzer_Integrator_GradBasisCrossVector_decl.hpp:78
Teuchos::RCP
panzer::Integrator_GradBasisCrossVector::_basis_name
std::string _basis_name
Definition: Panzer_Integrator_GradBasisCrossVector_decl.hpp:101
panzer::Integrator_GradBasisCrossVector::_field_multipliers
std::vector< PHX::MDField< const ScalarT, Cell, IP > > _field_multipliers
Definition: Panzer_Integrator_GradBasisCrossVector_decl.hpp:83
panzer::getBasisIndex
std::vector< std::string >::size_type getBasisIndex(std::string basis_name, const panzer::Workset &workset, WorksetDetailsAccessor &wda)
Returns the index in the workset bases for a particular BasisIRLayout name.
Definition: Panzer_Workset_Utilities.cpp:80
Teuchos::ParameterList::isType
bool isType(const std::string &name) const
panzer::Integrator_GradBasisCrossVector::_multiplier
ScalarT _multiplier
Definition: Panzer_Integrator_GradBasisCrossVector_decl.hpp:98
panzer::BasisValues2
Definition: Panzer_BasisValues2.hpp:63
panzer::PointRule::dl_vector
Teuchos::RCP< PHX::DataLayout > dl_vector
Data layout for vector fields.
Definition: Panzer_PointRule.hpp:114
panzer::Workset::num_cells
index_t num_cells
Definition: Panzer_Workset.hpp:200
panzer::createDynRankView
Kokkos::DynRankView< typename InputArray::value_type, PHX::Device > createDynRankView(const InputArray &a, const std::string &name, const DimensionPack... dims)
Wrapper to simplify Panzer use of Sacado ViewFactory.
Definition: Panzer_ViewFactory.hpp:54
getConst
const T & getConst(T &t)
panzer
Definition: Panzer_BasisValues_Evaluator_decl.hpp:54
Teuchos::ParameterList
field
PHX::MDField< ScalarT, panzer::Cell, panzer::BASIS > field
A field to which we'll contribute, or in which we'll store, the result of computing this integral.
Definition: Panzer_Integrator_CurlBasisDotVector_impl.hpp:590
PHX::FieldManager
Definition: Panzer_BCStrategy_Base.hpp:53
Panzer_IntegrationRule.hpp
TEUCHOS_TEST_FOR_EXCEPTION
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
panzer::Integrator_GradBasisCrossVector::postRegistrationSetup
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
Definition: Panzer_Integrator_GradBasisCrossVector_impl.hpp:128
panzer::Integrator_GradBasisCrossVector::_num_dims
std::size_t _num_dims
Definition: Panzer_Integrator_GradBasisCrossVector_decl.hpp:92
panzer::Integrator_GradBasisCrossVector::_residuals
std::vector< PHX::MDField< ScalarT, Cell, BASIS > > _residuals
Definition: Panzer_Integrator_GradBasisCrossVector_decl.hpp:80