Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_Integrator_BasisTimesVector_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_Integerator_BasisTimesVector_impl_hpp__
44#define __Panzer_Integerator_BasisTimesVector_impl_hpp__
45
47//
48// Include Files
49//
51
52// Panzer
56
57namespace panzer
58{
60 //
61 // Main Constructor
62 //
64 template<typename EvalT, typename Traits>
68 const std::string& resName,
69 const std::string& valName,
70 const panzer::BasisIRLayout& basis,
72 const double& multiplier, /* = 1 */
73 const std::vector<std::string>& fmNames /* =
74 std::vector<std::string>() */)
75 :
76 evalStyle_(evalStyle),
77 useDescriptors_(false),
78 multiplier_(multiplier),
79 basisName_(basis.name())
80 {
81 using PHX::View;
82 using panzer::BASIS;
83 using panzer::Cell;
85 using panzer::IP;
87 using PHX::MDField;
88 using PHX::print;
89 using std::invalid_argument;
90 using std::logic_error;
91 using std::string;
92 using Teuchos::RCP;
93
94 // Ensure the input makes sense.
95 TEUCHOS_TEST_FOR_EXCEPTION(resName == "", invalid_argument, "Error: " \
96 "Integrator_BasisTimesVector called with an empty residual name.")
97 TEUCHOS_TEST_FOR_EXCEPTION(valName == "", invalid_argument, "Error: " \
98 "Integrator_BasisTimesVector called with an empty value name.")
99 RCP<const PureBasis> tmpBasis = basis.getBasis();
100 TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->isVectorBasis(), logic_error,
101 "Error: Integrator_BasisTimesVector: Basis of type \""
102 << tmpBasis->name() << "\" is not a vector basis.");
103 TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->requiresOrientations(),
104 logic_error, "Integrator_BasisTimesVector: Basis of type \""
105 << tmpBasis->name() << "\" does not require orientations. This seems " \
106 "very strange, so I'm failing.");
107
108 // Create the field for the vector-valued quantity we're integrating.
109 vector_ = MDField<const ScalarT, Cell, IP, Dim>(valName, ir.dl_vector);
110 this->addDependentField(vector_);
111
112 // Create the field that we're either contributing to or evaluating
113 // (storing).
114 field_ = MDField<ScalarT, Cell, BASIS>(resName, basis.functional);
116 this->addContributedField(field_);
117 else // if (evalStyle == EvaluatorStyle::EVALUATES)
118 this->addEvaluatedField(field_);
119
120 // Add the dependent field multipliers, if there are any.
121 int i(0);
122 fieldMults_.resize(fmNames.size());
124 View<Kokkos::View<const ScalarT**, typename PHX::DevLayout<ScalarT>::type, Kokkos::MemoryUnmanaged>*>("BasisTimesVector::KokkosFieldMultipliers",
125 fmNames.size());
126 for (const auto& name : fmNames)
127 {
128 fieldMults_[i++] = MDField<const ScalarT, Cell, IP>(name, ir.dl_scalar);
129 this->addDependentField(fieldMults_[i - 1]);
130 } // end loop over the field multipliers
131
132 // Set the name of this object.
133 string n("Integrator_BasisTimesVector (");
135 n += "Cont";
136 else // if (evalStyle == EvaluatorStyle::EVALUATES)
137 n += "Eval";
138 n += ", " + print<EvalT>() + "): " + field_.fieldTag().name();
139 this->setName(n);
140 } // end of Main Constructor
141
143 //
144 // Descriptor Constructor
145 //
147 template<typename EvalT, typename Traits>
151 const PHX::FieldTag& resTag,
152 const PHX::FieldTag& valTag,
153 const BasisDescriptor& bd,
154 const IntegrationDescriptor& id,
155 const double& multiplier, /* = 1 */
156 const std::vector<PHX::FieldTag>& multipliers /* =
157 std::vector<PHX::FieldTag>() */)
158 :
159 evalStyle_(evalStyle),
160 useDescriptors_(true),
161 bd_(bd),
162 id_(id),
163 multiplier_(multiplier)
164 {
165 using PHX::View;
166 using panzer::BASIS;
167 using panzer::Cell;
169 using panzer::IP;
170 using panzer::PureBasis;
171 using PHX::MDField;
172 using PHX::print;
173 using std::invalid_argument;
174 using std::logic_error;
175 using std::string;
176 using Teuchos::RCP;
177
178 // Ensure the input makes sense.
179 TEUCHOS_TEST_FOR_EXCEPTION(not ((bd_.getType() == "HCurl") or
180 (bd_.getType() == "HDiv")), logic_error, "Error: " \
181 "Integrator_BasisTimesVector: Basis of type \"" << bd_.getType()
182 << "\" is not a vector basis.")
183
184 // Create the field for the vector-valued quantity we're integrating.
185 vector_ = valTag;
186 this->addDependentField(vector_);
187
188 // Create the field that we're either contributing to or evaluating
189 // (storing).
190 field_ = resTag;
192 this->addContributedField(field_);
193 else // if (evalStyle == EvaluatorStyle::EVALUATES)
194 this->addEvaluatedField(field_);
195
196 // Add the dependent field multipliers, if there are any.
197 int i(0);
198 fieldMults_.resize(multipliers.size());
200 View<Kokkos::View<const ScalarT**, typename PHX::DevLayout<ScalarT>::type, Kokkos::MemoryUnmanaged>*>("BasisTimesVector::KokkosFieldMultipliers",
201 multipliers.size());
202 for (const auto& fm : multipliers)
203 {
204 fieldMults_[i++] = fm;
205 this->addDependentField(fieldMults_[i - 1]);
206 } // end loop over the field multipliers
207
208 // Set the name of this object.
209 string n("Integrator_BasisTimesVector (");
211 n += "Cont";
212 else // if (evalStyle == EvaluatorStyle::EVALUATES)
213 n += "Eval";
214 n += ", " + print<EvalT>() + "): " + field_.fieldTag().name();
215 this->setName(n);
216 } // end of Main Constructor
217
219 //
220 // ParameterList Constructor
221 //
223 template<typename EvalT, typename Traits>
226 const Teuchos::ParameterList& p)
227 :
230 p.get<std::string>("Residual Name"),
231 p.get<std::string>("Value Name"),
232 (*p.get<Teuchos::RCP<panzer::BasisIRLayout>>("Basis")),
233 (*p.get<Teuchos::RCP<panzer::IntegrationRule>>("IR")),
234 p.get<double>("Multiplier"),
235 p.isType<Teuchos::RCP<const std::vector<std::string>>>
236 ("Field Multipliers") ?
237 (*p.get<Teuchos::RCP<const std::vector<std::string>>>
238 ("Field Multipliers")) : std::vector<std::string>())
239 {
240 using Teuchos::ParameterList;
241 using Teuchos::RCP;
242
243 // Ensure that the input ParameterList didn't contain any bogus entries.
244 RCP<ParameterList> validParams = this->getValidParameters();
245 p.validateParameters(*validParams);
246 } // end of ParameterList Constructor
247
249 //
250 // postRegistrationSetup()
251 //
253 template<typename EvalT, typename Traits>
254 void
257 typename Traits::SetupData sd,
259 {
261 using std::size_t;
262
263 // Get the PHX::Views of the field multipliers.
264 auto kokkosFieldMults_h = Kokkos::create_mirror_view(kokkosFieldMults_);
265 for (size_t i(0); i < fieldMults_.size(); ++i)
266 kokkosFieldMults_h(i) = fieldMults_[i].get_static_view();
267 Kokkos::deep_copy(kokkosFieldMults_, kokkosFieldMults_h);
268
269 // Determine the number of quadrature points and the dimensionality of the
270 // vector that we're integrating.
271 numQP_ = vector_.extent(1);
272 numDim_ = vector_.extent(2);
273
274 // Determine the index in the Workset bases for our particular basis name.
275 if (not useDescriptors_)
276 basisIndex_ = getBasisIndex(basisName_, (*sd.worksets_)[0], this->wda);
277
278 // Allocate temporary memory
279 if (Sacado::IsADType<ScalarT>::value) {
280 const auto fadSize = Kokkos::dimension_scalar(field_.get_static_view());
281 tmp_ = PHX::View<ScalarT*>("panzer::Integrator::BasisTimesVector::tmp_",field_.extent(0),fadSize);
282 } else {
283 tmp_ = PHX::View<ScalarT*>("panzer::Integrator::BasisTimesVector::tmp_",field_.extent(0));
284 }
285
286 } // end of postRegistrationSetup()
287
289 //
290 // operator()()
291 //
293 template<typename EvalT, typename Traits>
294 template<int NUM_FIELD_MULT>
295 KOKKOS_INLINE_FUNCTION
296 void
299 const FieldMultTag<NUM_FIELD_MULT>& /* tag */,
300 const size_t& cell) const
301 {
303
304 // Initialize the evaluated field.
305 const int numBases(basis_.extent(1));
306 if (evalStyle_ == EvaluatorStyle::EVALUATES)
307 for (int basis(0); basis < numBases; ++basis)
308 field_(cell, basis) = 0.0;
309
310 // The following if-block is for the sake of optimization depending on the
311 // number of field multipliers.
312 if (NUM_FIELD_MULT == 0)
313 {
314 // Loop over the quadrature points and dimensions of our vector fields,
315 // scale the integrand by the multiplier, and then perform the
316 // integration, looping over the bases.
317 for (int qp(0); qp < numQP_; ++qp)
318 {
319 for (int dim(0); dim < numDim_; ++dim)
320 {
321 for (int basis(0); basis < numBases; ++basis)
322 field_(cell, basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim);
323 } // end loop over the dimensions of the vector field
324 } // end loop over the quadrature points
325 }
326 else if (NUM_FIELD_MULT == 1)
327 {
328 // Loop over the quadrature points and dimensions of our vector fields,
329 // scale the integrand by the multiplier and the single field multiplier,
330 // and then perform the actual integration, looping over the bases.
331 for (int qp(0); qp < numQP_; ++qp)
332 {
333 for (int dim(0); dim < numDim_; ++dim)
334 {
335 for (int basis(0); basis < numBases; ++basis)
336 field_(cell, basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim) * kokkosFieldMults_(0)(cell, qp);
337 } // end loop over the dimensions of the vector field
338 } // end loop over the quadrature points
339 }
340 else
341 {
342 // Loop over the quadrature points and pre-multiply all the field
343 // multipliers together. Then loop over the dimensions of our vector
344 // fields, scale the integrand by the multiplier and the combination of
345 // the field multipliers, and then perform the actual integration,
346 // looping over the bases.
347 const int numFieldMults(kokkosFieldMults_.extent(0));
348 for (int qp(0); qp < numQP_; ++qp)
349 {
350 tmp_(cell) = 1.0;
351 for (int fm(0); fm < numFieldMults; ++fm)
352 tmp_(cell) *= kokkosFieldMults_(fm)(cell, qp);
353 for (int dim(0); dim < numDim_; ++dim)
354 {
355 for (int basis(0); basis < numBases; ++basis)
356 field_(cell, basis) += basis_(cell, basis, qp, dim) * multiplier_ * vector_(cell, qp, dim) * tmp_(cell);
357 } // end loop over the dimensions of the vector field
358 } // end loop over the quadrature points
359 } // end if (NUM_FIELD_MULT == something)
360 } // end of operator()()
361
363 //
364 // evaluateFields()
365 //
367 template<typename EvalT, typename Traits>
368 void
371 typename Traits::EvalData workset)
372 {
373 using Kokkos::parallel_for;
374 using Kokkos::RangePolicy;
375
376 // Grab the basis information.
377 const panzer::BasisValues2<double>& bv = useDescriptors_ ?
378 this->wda(workset).getBasisValues(bd_,id_) :
379 *this->wda(workset).bases[basisIndex_];
381 basis_ = useDescriptors_ ? bv.getVectorBasisValues(true) : Array(bv.weighted_basis_vector);
382
383 // The following if-block is for the sake of optimization depending on the
384 // number of field multipliers. The parallel_fors will loop over the cells
385 // in the Workset and execute operator()() above.
386 if (fieldMults_.size() == 0)
387 parallel_for(RangePolicy<FieldMultTag<0>>(0, workset.num_cells), *this);
388 else if (fieldMults_.size() == 1)
389 parallel_for(RangePolicy<FieldMultTag<1>>(0, workset.num_cells), *this);
390 else
391 parallel_for(RangePolicy<FieldMultTag<-1>>(0, workset.num_cells), *this);
392 } // end of evaluateFields()
393
395 //
396 // getValidParameters()
397 //
399 template<typename EvalT, typename Traits>
400 Teuchos::RCP<Teuchos::ParameterList>
403 {
406 using std::string;
407 using std::vector;
408 using Teuchos::ParameterList;
409 using Teuchos::RCP;
410 using Teuchos::rcp;
411
412 // Create a ParameterList with all the valid keys we support.
413 RCP<ParameterList> p = rcp(new ParameterList);
414 p->set<string>("Residual Name", "?");
415 p->set<string>("Value Name", "?");
416 RCP<BasisIRLayout> basis;
417 p->set("Basis", basis);
418 RCP<IntegrationRule> ir;
419 p->set("IR", ir);
420 p->set<double>("Multiplier", 1.0);
421 RCP<const vector<string> > fms;
422 p->set("Field Multipliers", fms);
423 return p;
424 } // end of getValidParameters()
425
426} // end of namespace panzer
427
428#endif // __Panzer_Integerator_BasisTimesVector_impl_hpp__
panzer::EvaluatorStyle evalStyle
The EvaluatorStyle of the parent Integrator_CurlBasisDotVector object.
double multiplier
The scalar multiplier out in front of the integral ( ).
const std::string & getType() const
Get type of basis.
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis>
Teuchos::RCP< const PureBasis > getBasis() const
ConstArray_CellBasisIP getBasisValues(const bool weighted, const bool cache=true, const bool force=false) const
Get the basis values evaluated at mesh points.
PHX::MDField< const Scalar, Cell, BASIS, IP, Dim > ConstArray_CellBasisIPDim
Array_CellBasisIPDim weighted_basis_vector
ConstArray_CellBasisIPDim getVectorBasisValues(const bool weighted, const bool cache=true, const bool force=false) const
Get the vector basis values evaluated at mesh points.
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.
Integrator_BasisTimesVector(const panzer::EvaluatorStyle &evalStyle, const std::string &resName, const std::string &valName, const panzer::BasisIRLayout &basis, const panzer::IntegrationRule &ir, const double &multiplier=1, const std::vector< std::string > &fmNames=std::vector< std::string >())
Main Constructor.
PHX::MDField< const ScalarT, panzer::Cell, panzer::IP, panzer::Dim > vector_
A field representing the vector-valued function we're integrating ( ).
BasisDescriptor bd_
The BasisDescriptor for the basis to use.
void postRegistrationSetup(typename Traits::SetupData sd, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const std::size_t &cell) const
Perform the integration.
std::vector< PHX::MDField< const ScalarT, panzer::Cell, panzer::IP > > fieldMults_
The (possibly empty) list of fields that are multipliers out in front of the integral ( ,...
PHX::View< Kokkos::View< const ScalarT **, typename PHX::DevLayout< ScalarT >::type, Kokkos::MemoryUnmanaged > * > kokkosFieldMults_
The PHX::View representation of the (possibly empty) list of fields that are multipliers out in front...
void evaluateFields(typename Traits::EvalData workset)
Evaluate Fields.
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
Teuchos::RCP< PHX::DataLayout > dl_vector
Data layout for vector fields.
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_