Panzer  Version of the Day
Panzer_Integrator_GradBasisTimesScalar_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_GRADBASISTIMESSCALAR_IMPL_HPP
44 #define PANZER_EVALUATOR_GRADBASISTIMESSCALAR_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  _num_dims = ir->dl_vector->extent(2);
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_GradBasisTimesScalar: 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 
85  // Currently we only allow for vectors of length 3
86  TEUCHOS_TEST_FOR_EXCEPTION(_num_dims != _residuals.size(),std::logic_error,"Number of residuals must equal number of gradient components (mesh dimensions).");
87 
88  // Setup vector
89  _scalar = PHX::MDField<const ScalarT,Cell,IP>(p.get<std::string>("Scalar Name"), ir->dl_scalar);
90 
91  // Read the scalar multiplier
92  _multiplier = ScalarT(p.get<double>("Multiplier"));
93 
94  // Setup field multipliers
95  if (p.isType<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers")){
96  const std::vector<std::string> & field_multiplier_names = *(p.get<Teuchos::RCP<const std::vector<std::string> > >("Field Multipliers"));
97  for (const std::string & name : field_multiplier_names){
98  _field_multipliers.push_back(PHX::MDField<const ScalarT,Cell,IP>(name, ir->dl_scalar));
99  }
100  }
101 
102  // Register evaluated fields
103  for (auto & residual : _residuals){
104  this->addEvaluatedField(residual);
105  }
106 
107  // Register dependent fields
108  this->addDependentField(_scalar);
109  for (auto & field : _field_multipliers){
110  this->addDependentField(field);
111  }
112 
113  // Name the module
114  this->setName("Integrator_GradBasisTimesScalar: " + _scalar.fieldTag().name());
115 }
116 
117 //**********************************************************************
118 template<typename EvalT, typename Traits>
119 void
122  typename Traits::SetupData sd,
123  PHX::FieldManager<Traits>& /* fm */)
124 {
125  _num_basis_nodes = _residuals[0].extent(1);
126  _num_quadrature_points = _scalar.extent(1);
127 
128  _basis_index = panzer::getBasisIndex(_basis_name, (*sd.worksets_)[0], this->wda);
129 
130  // TODO: figure out a clean way of cloning _vector
131  _tmp = Kokkos::createDynRankView(_residuals[0].get_static_view(),"tmp",_scalar.extent(0), _num_quadrature_points);
132 }
133 
134 //**********************************************************************
135 template<typename EvalT, typename Traits>
136 void
139  typename Traits::EvalData workset)
140 {
141 
142  // Zero the residuals
143  for (int i(0); i < static_cast<int>(_num_dims); ++i)
144  Kokkos::deep_copy(_residuals[i].get_static_view(), ScalarT(0.0));
145 
146  // do a scaled copy to initialize _tmp
147  for (int i=0; i < _scalar.extent_int(0); ++i)
148  for (int j=0; j < _scalar.extent_int(1); ++j)
149  _tmp(i,j) = _multiplier * _scalar(i,j);
150 
151  // Apply the field multipliers
152  for(auto & field_data : _field_multipliers){
153  for (index_t cell = 0; cell < workset.num_cells; ++cell) {
154  for (std::size_t qp = 0; qp < _num_quadrature_points; ++qp) {
155  _tmp(cell,qp) *= field_data(cell,qp);
156  }
157  }
158  }
159 
160  const BasisValues2<double> & bv = *this->wda(workset).bases[_basis_index];
161 
162  for (index_t cell = 0; cell < workset.num_cells; ++cell){
163  for (std::size_t basis = 0; basis < _num_basis_nodes; ++basis){
164  for (std::size_t qp = 0; qp < _num_quadrature_points; ++qp){
165  for(std::size_t dim = 0; dim < _num_dims; ++dim){
166  _residuals[dim](cell,basis) += _tmp(cell,qp)*bv.weighted_grad_basis(cell,basis,qp,dim);
167  }
168  }
169  }
170  }
171 
172 }
173 
174 //**********************************************************************
175 
176 }
177 
178 #endif
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
T & get(const std::string &name, T def_value)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
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.
bool isType(const std::string &name) const
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.
std::vector< PHX::MDField< ScalarT, Cell, BASIS > > _residuals
Array_CellBasisIPDim weighted_grad_basis
Teuchos::RCP< PHX::DataLayout > dl_vector
Data layout for vector fields.
PHX::MDField< ScalarT, panzer::Cell, panzer::BASIS > field
A field to which we&#39;ll contribute, or in which we&#39;ll store, the result of computing this integral...
const T & getConst(T &t)
std::vector< PHX::MDField< const ScalarT, Cell, IP > > _field_multipliers
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_