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 PHX_EVALUATOR_CTOR(Integrator_GradBasisCrossVector,p):
58  _basis_index(-1)
59 {
60 
61  // TODO: Get all of the panzer modules to read in const integration rules and basisIRlayouts
64  Teuchos::RCP<PHX::DataLayout> vector_dl = p.get< Teuchos::RCP<PHX::DataLayout> >("Data Layout Vector");
65 
67  _basis_name = basis_layout->name();
68 
69  // Verify that this basis supports the gradient operation
70  TEUCHOS_TEST_FOR_EXCEPTION(!basis->supportsGrad(),std::logic_error,"Integrator_GradBasisCrossVector: Basis of type \"" << basis->name() << "\" does not support GRAD");
71 
72  // Setup residuals
73  _residuals.clear();
74  if (p.isType<Teuchos::RCP<const std::vector<std::string> > >("Residual Names")){
75  const std::vector<std::string> & names = *(p.get<Teuchos::RCP<const std::vector<std::string> > >("Residual Names"));
76  for(const auto & name : names){
77  PHX::MDField<ScalarT,Cell,BASIS> res(name, basis_layout->functional);
78  _residuals.push_back(res);
79  }
80  }
81  _num_dims = _residuals.size();
82 
83  // Currently we only allow for vectors of length 3
84  TEUCHOS_TEST_FOR_EXCEPTION(_num_dims != vector_dl->dimension(2),std::logic_error,"Vector must be same length as number of residuals.");
85 
86  // Setup vector
87  _vector = PHX::MDField<const ScalarT,Cell,IP,Dim>(p.get<std::string>("Vector Name"), vector_dl);
88 
89  // 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)
90  _num_grad_dims = ir->dl_vector->dimension(2);
91 
92  // Vector must be at least the same dimension as the GRAD operation
93  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.");
94 
95  // Read the scalar multiplier
96  _multiplier = ScalarT(p.get<double>("Multiplier"));
97 
98  // Setup field multipliers
99  if (p.isType<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers")){
100  const std::vector<std::string> & field_multiplier_names = *(p.get<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers"));
101  for (const std::string & name : field_multiplier_names){
102  _field_multipliers.push_back(PHX::MDField<const ScalarT,Cell,IP>(name, ir->dl_scalar));
103  }
104  }
105 
106  // Register evaluated fields
107  for (auto & residual : _residuals){
108  this->addEvaluatedField(residual);
109  }
110 
111  // Register dependent fields
112  this->addDependentField(_vector);
113  for (auto & field : _field_multipliers){
114  this->addDependentField(field);
115  }
116 
117  // Name the module
118  this->setName("Integrator_GradBasisCrossVector: " + _vector.fieldTag().name());
119 }
120 
121 //**********************************************************************
122 PHX_POST_REGISTRATION_SETUP(Integrator_GradBasisCrossVector,sd,fm)
123 {
124 
125  for (auto & residual : _residuals){
126  this->utils.setFieldData(residual,fm);
127  }
128 
129  this->utils.setFieldData(_vector,fm);
130 
131  for (auto & field : _field_multipliers){
132  this->utils.setFieldData(field,fm);
133  }
134 
135  _num_basis_nodes = _residuals[0].dimension(1);
136  _num_quadrature_points = _vector.dimension(1);
137 
138  _basis_index = panzer::getBasisIndex(_basis_name, (*sd.worksets_)[0], this->wda);
139 
140  // TODO: figure out a clean way of cloning _vector
141  _tmp = Kokkos::createDynRankView(_residuals[0].get_static_view(),"tmp",_vector.dimension(0), _num_quadrature_points, _num_dims);
142 }
143 
144 //**********************************************************************
145 PHX_EVALUATE_FIELDS(Integrator_GradBasisCrossVector,workset)
146 {
147 
148  // Zero the residuals
149  for(int i=0; i<_num_dims; ++i){
150  Kokkos::deep_copy(_residuals[i].get_static_view(), ScalarT(0.0));
151  }
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
PHX::MDField< ScalarT > residual
Evaluates a Dirichlet BC residual corresponding to a field value.
std::vector< std::string >::size_type getBasisIndex(std::string basis_name, panzer::Workset &workset, WorksetDetailsAccessor &wda)
Returns the index in the workset bases for a particular BasisIRLayout name.
std::string name() const
A unique key that is the combination of the basis type and basis order.
Kokkos::DynRankView< ScalarT, PHX::Device > _tmp
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
std::vector< PHX::MDField< ScalarT, Cell, BASIS > > _residuals
Teuchos::RCP< const PureBasis > getBasis() const
bool supportsGrad() const
T * get() const
PHX::MDField< const ScalarT, Cell, IP, Dim > _vector
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.
Teuchos::RCP< PHX::DataLayout > dl_scalar
Data layout for scalar fields.
PHX_EVALUATOR_CTOR(BasisValues_Evaluator, p)
PHX::MDField< const ScalarT, Cell, IP > field
Array_CellBasisIPDim weighted_grad_basis
PHX_EVALUATE_FIELDS(BasisValues_Evaluator, workset)
Teuchos::RCP< PHX::DataLayout > dl_vector
Data layout for vector fields.
std::string name() const
Unique key for workset indexing composed of basis name and point rule name.
Teuchos::RCP< const panzer::PureBasis > basis
Interpolates basis DOF values to IP DOF values.
PHX_POST_REGISTRATION_SETUP(BasisValues_Evaluator, sd, fm)
const T & getConst(T &t)
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis>
std::vector< PHX::MDField< const ScalarT, Cell, IP > > _field_multipliers