Panzer Version of the Day
Loading...
Searching...
No Matches
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
47//
48// Include Files
49//
51
52// Intrepid2
53#include "Intrepid2_FunctionSpaceTools.hpp"
54
55// Kokkos
56#include "Kokkos_ViewFactory.hpp"
57
58// Panzer
62
63namespace panzer
64{
66 //
67 // Constructor
68 //
70 template<typename EvalT, typename Traits>
74 const std::vector<std::string>& resNames,
75 const std::string& scalar,
76 const panzer::BasisIRLayout& basis,
78 const double& multiplier, /* = 1 */
79 const std::vector<std::string>& fmNames) /* =
80 std::vector<std::string>() */
81 :
82 evalStyle_(evalStyle),
83 multiplier_(multiplier),
84 numDims_(ir.dl_vector->extent(2)),
85 basisName_(basis.name())
86 {
87 using PHX::View;
88 using panzer::BASIS;
89 using panzer::Cell;
91 using panzer::IP;
93 using PHX::Device;
94 using PHX::DevLayout;
95 using PHX::MDField;
96 using std::invalid_argument;
97 using std::logic_error;
98 using std::string;
99 using Teuchos::RCP;
100
101 // Ensure the input makes sense.
102 for (const auto& name : resNames)
103 TEUCHOS_TEST_FOR_EXCEPTION(name == "", invalid_argument, "Error: " \
104 "Integrator_GradBasisTimesScalar called with an empty residual name.")
105 TEUCHOS_TEST_FOR_EXCEPTION(scalar == "", invalid_argument, "Error: " \
106 "Integrator_GradBasisTimesScalar called with an empty scalar name.")
107 TEUCHOS_TEST_FOR_EXCEPTION(numDims_ != static_cast<int>(resNames.size()),
108 logic_error, "Error: Integrator_GradBasisTimesScalar: The number of " \
109 "residuals must equal the number of gradient components (mesh " \
110 "dimensions).")
111 RCP<const PureBasis> tmpBasis = basis.getBasis();
112 TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->supportsGrad(), logic_error,
113 "Error: Integrator_GradBasisTimesScalar: Basis of type \""
114 << tmpBasis->name() << "\" does not support the gradient operation.")
115
116 // Create the field or the scalar function we're integrating.
117 scalar_ = MDField<const ScalarT, Cell, IP>(scalar, ir.dl_scalar);
118 this->addDependentField(scalar_);
119
120 // Create the fields that we're either contributing to or evaluating
121 // (storing).
122 fields_host_.resize(resNames.size());
123 fields_ = OuterView("Integrator_GradBasisCrossVector::fields_", resNames.size());
124 {
125 int i=0;
126 for (const auto& name : resNames)
127 fields_host_[i++] = PHX::MDField<ScalarT,Cell,BASIS>(name,basis.functional);
128 }
129 for (std::size_t i=0; i<fields_.extent(0); ++i) {
130 const auto& field = fields_host_[i];
132 this->addContributedField(field);
133 else // if (evalStyle_ == EvaluatorStyle::EVALUATES)
134 this->addEvaluatedField(field);
135 }
136
137 // Add the dependent field multipliers, if there are any.
138 int i = 0;
139 fieldMults_.resize(fmNames.size());
140 kokkosFieldMults_ = PHX::View<PHX::UnmanagedView<const ScalarT**>*>(
141 "GradBasisTimesScalar::KokkosFieldMultipliers", fmNames.size());
142 for (const auto& name : fmNames)
143 {
144 fieldMults_[i++] = MDField<const ScalarT, Cell, IP>(name, ir.dl_scalar);
145 this->addDependentField(fieldMults_[i - 1]);
146 } // end loop over the field multipliers
147
148 // Set the name of this object.
149 string n("Integrator_GradBasisTimesScalar (");
151 n += "CONTRIBUTES";
152 else // if (evalStyle_ == EvaluatorStyle::EVALUATES)
153 n += "EVALUATES";
154 n += "): {";
155 for (size_t j=0; j < fields_host_.size() - 1; ++j)
156 n += resNames[j] + ", ";
157 n += resNames[resNames.size()-1] + "}";
158 this->setName(n);
159 } // end of Constructor
160
162 //
163 // ParameterList Constructor
164 //
166 template<typename EvalT, typename Traits>
169 const Teuchos::ParameterList& p)
170 :
173 p.get<const std::vector<std::string>>("Residual Names"),
174 p.get<std::string>("Scalar Name"),
175 (*p.get<Teuchos::RCP<panzer::BasisIRLayout>>("Basis")),
176 (*p.get<Teuchos::RCP<panzer::IntegrationRule>>("IR")),
177 p.get<double>("Multiplier"),
178 p.isType<Teuchos::RCP<const std::vector<std::string>>>
179 ("Field Multipliers") ?
180 (*p.get<Teuchos::RCP<const std::vector<std::string>>>
181 ("Field Multipliers")) : std::vector<std::string>())
182 {
183 using Teuchos::ParameterList;
184 using Teuchos::RCP;
185
186 // Ensure that the input ParameterList didn't contain any bogus entries.
187 RCP<ParameterList> validParams = this->getValidParameters();
188 p.validateParameters(*validParams);
189 } // end of ParameterList Constructor
190
192 //
193 // postRegistrationSetup()
194 //
196 template<typename EvalT, typename Traits>
197 void
200 typename Traits::SetupData sd,
202 {
203 using Kokkos::createDynRankView;
205
206 // Get the PHX::Views of the fields.
207 auto fields_host_mirror = Kokkos::create_mirror_view(fields_);
208 for (size_t i(0); i < fields_host_.size(); ++i)
209 fields_host_mirror(i) = fields_host_[i].get_static_view();
210 Kokkos::deep_copy(fields_,fields_host_mirror);
211
212 // Get the PHX::Views of the field multipliers.
213 auto field_mults_host_mirror = Kokkos::create_mirror_view(kokkosFieldMults_);
214 for (size_t i(0); i < fieldMults_.size(); ++i)
215 field_mults_host_mirror(i) = fieldMults_[i].get_static_view();
216 Kokkos::deep_copy(kokkosFieldMults_,field_mults_host_mirror);
217
218 // Determine the index in the Workset bases for our particular basis name.
219 basisIndex_ = getBasisIndex(basisName_, (*sd.worksets_)[0], this->wda);
220 } // end of postRegistrationSetup()
221
222 template<typename EvalT, typename Traits>
223 template<int NUM_FIELD_MULT>
224 KOKKOS_INLINE_FUNCTION
225 void
228 const FieldMultTag<NUM_FIELD_MULT>& /* tag */,
229 const size_t& cell) const
230 {
232
233 // Initialize the evaluated fields.
234 const int numQP(scalar_.extent(1)), numBases(fields_[0].extent(1));
235 if (evalStyle_ == EvaluatorStyle::EVALUATES)
236 for (int dim(0); dim < numDims_; ++dim)
237 for (int basis(0); basis < numBases; ++basis)
238 fields_[dim](cell, basis) = 0.0;
239
240 // The following if-block is for the sake of optimization depending on the
241 // number of field multipliers.
242 ScalarT tmp;
243 if (NUM_FIELD_MULT == 0)
244 {
245 // Loop over the quadrature points, scale the integrand by the
246 // multiplier, and then perform the actual integration, looping over the
247 // bases and dimensions.
248 for (int qp(0); qp < numQP; ++qp)
249 {
250 tmp = multiplier_ * scalar_(cell, qp);
251 for (int basis(0); basis < numBases; ++basis)
252 for (int dim(0); dim < numDims_; ++dim)
253 fields_[dim](cell, basis) += basis_(cell, basis, qp, dim) * tmp;
254 } // end loop over the quadrature points
255 }
256 else if (NUM_FIELD_MULT == 1)
257 {
258 // Loop over the quadrature points, scale the integrand by the multiplier
259 // and the single field multiplier, and then perform the actual
260 // integration, looping over the bases and dimensions.
261 for (int qp(0); qp < numQP; ++qp)
262 {
263 tmp = multiplier_ * scalar_(cell, qp) * kokkosFieldMults_(0)(cell, qp);
264 for (int basis(0); basis < numBases; ++basis)
265 for (int dim(0); dim < numDims_; ++dim)
266 fields_[dim](cell, basis) += basis_(cell, basis, qp, dim) * tmp;
267 } // end loop over the quadrature points
268 }
269 else
270 {
271 // Loop over the quadrature points and pre-multiply all the field
272 // multipliers together, scale the integrand by the multiplier and
273 // the combination of the field multipliers, and then perform the actual
274 // integration, looping over the bases and dimensions.
275 const int numFieldMults(kokkosFieldMults_.extent(0));
276 for (int qp(0); qp < numQP; ++qp)
277 {
278 ScalarT fieldMultsTotal(1);
279 for (int fm(0); fm < numFieldMults; ++fm)
280 fieldMultsTotal *= kokkosFieldMults_(fm)(cell, qp);
281 tmp = multiplier_ * scalar_(cell, qp) * fieldMultsTotal;
282 for (int basis(0); basis < numBases; ++basis)
283 for (int dim(0); dim < numDims_; ++dim)
284 fields_[dim](cell, basis) += basis_(cell, basis, qp, dim) * tmp;
285 } // end loop over the quadrature points
286 } // end if (NUM_FIELD_MULT == something)
287 } // end of operator()()
288
290 //
291 // evaluateFields()
292 //
294 template<typename EvalT, typename Traits>
295 void
298 typename Traits::EvalData workset)
299 {
300 using Kokkos::parallel_for;
301 using Kokkos::RangePolicy;
302
303 // Grab the basis information.
304 basis_ = this->wda(workset).bases[basisIndex_]->weighted_grad_basis;
305
306 // The following if-block is for the sake of optimization depending on the
307 // number of field multipliers. The parallel_fors will loop over the cells
308 // in the Workset and execute operator()() above.
309 if (fieldMults_.size() == 0)
310 parallel_for(RangePolicy<FieldMultTag<0>>(0, workset.num_cells), *this);
311 else if (fieldMults_.size() == 1)
312 parallel_for(RangePolicy<FieldMultTag<1>>(0, workset.num_cells), *this);
313 else
314 parallel_for(RangePolicy<FieldMultTag<-1>>(0, workset.num_cells), *this);
315 } // end of evaluateFields()
316
318 //
319 // getValidParameters()
320 //
322 template<typename EvalT, typename TRAITS>
323 Teuchos::RCP<Teuchos::ParameterList>
326 {
329 using std::string;
330 using std::vector;
331 using Teuchos::ParameterList;
332 using Teuchos::RCP;
333 using Teuchos::rcp;
334
335 // Create a ParameterList with all the valid keys we support.
336 RCP<ParameterList> p = rcp(new ParameterList);
337
338 RCP<const vector<string>> resNames;
339 p->set("Residual Names", resNames);
340 p->set<string>("Scalar Name", "?");
341 RCP<BasisIRLayout> basis;
342 p->set("Basis", basis);
343 RCP<IntegrationRule> ir;
344 p->set("IR", ir);
345 p->set<double>("Multiplier", 1.0);
346 RCP<const vector<string>> fms;
347 p->set("Field Multipliers", fms);
348
349 return p;
350 } // end of getValidParameters()
351
352} // end of namespace panzer
353
354#endif // PANZER_EVALUATOR_GRADBASISTIMESSCALAR_IMPL_HPP
panzer::EvaluatorStyle evalStyle
The EvaluatorStyle of the parent Integrator_CurlBasisDotVector object.
double multiplier
The scalar multiplier out in front of the integral ( ).
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.
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis>
Teuchos::RCP< const PureBasis > getBasis() const
std::vector< PHX::MDField< ScalarT, Cell, BASIS > > fields_host_
The fields to which we'll contribute, or in which we'll store, the result of computing this integral.
int numDims_
The number of dimensions associated with the gradient.
PHX::MDField< const ScalarT, Cell, IP > scalar_
A field representing the scalar function we're integrating ( ).
Integrator_GradBasisTimesScalar(const panzer::EvaluatorStyle &evalStyle, const std::vector< std::string > &resNames, const std::string &scalar, const panzer::BasisIRLayout &basis, const panzer::IntegrationRule &ir, const double &multiplier=1, const std::vector< std::string > &fmNames=std::vector< std::string >())
Main Constructor.
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const std::size_t &cell) const
Perform the integration.
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
std::vector< PHX::MDField< const ScalarT, Cell, IP > > fieldMults_
The scalar multiplier out in front of the integral ( ).
PHX::View< PHX::UnmanagedView< const ScalarT ** > * > kokkosFieldMults_
The PHX::View representation of the (possibly empty) list of fields that are multipliers out in front...
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
void evaluateFields(typename Traits::EvalData d)
Evaluate Fields.
const panzer::EvaluatorStyle evalStyle_
An enum determining the behavior of this Evaluator.
Teuchos::RCP< PHX::DataLayout > dl_scalar
Data layout for scalar fields.
Description and data layouts associated with a particular basis.
int num_cells
DEPRECATED - use: numCells()
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.
EvaluatorStyle
An indication of how an Evaluator will behave.
This empty struct allows us to optimize operator()() depending on the number of field multipliers.
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_