Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_Integrator_BasisTimesScalar_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_INTEGRATOR_BASISTIMESSCALAR_IMPL_HPP
44#define PANZER_INTEGRATOR_BASISTIMESSCALAR_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 multiplier_(multiplier),
78 basisName_(basis.name())
79 {
80 using PHX::View;
81 using panzer::BASIS;
82 using panzer::Cell;
84 using panzer::IP;
86 using PHX::MDField;
87 using std::invalid_argument;
88 using std::logic_error;
89 using std::string;
90 using Teuchos::RCP;
91
92 // Ensure the input makes sense.
93 TEUCHOS_TEST_FOR_EXCEPTION(resName == "", invalid_argument, "Error: " \
94 "Integrator_BasisTimesScalar called with an empty residual name.")
95 TEUCHOS_TEST_FOR_EXCEPTION(valName == "", invalid_argument, "Error: " \
96 "Integrator_BasisTimesScalar called with an empty value name.")
97 RCP<const PureBasis> tmpBasis = basis.getBasis();
98 TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->isScalarBasis(), logic_error,
99 "Error: Integrator_BasisTimesScalar: Basis of type \""
100 << tmpBasis->name() << "\" is not a scalar basis.")
101
102 // Create the field for the scalar quantity we're integrating.
103 scalar_ = MDField<const ScalarT, Cell, IP>(valName, ir.dl_scalar);
104 this->addDependentField(scalar_);
105
106 // Create the field that we're either contributing to or evaluating
107 // (storing).
108 field_ = MDField<ScalarT, Cell, BASIS>(resName, basis.functional);
110 this->addContributedField(field_);
111 else // if (evalStyle == EvaluatorStyle::EVALUATES)
112 this->addEvaluatedField(field_);
113
114 // Add the dependent field multipliers, if there are any.
115 int i(0);
116 fieldMults_.resize(fmNames.size());
118 View<Kokkos::View<const ScalarT**, typename PHX::DevLayout<ScalarT>::type, Kokkos::MemoryUnmanaged>*>("BasisTimesScalar::KokkosFieldMultipliers",
119 fmNames.size());
120 for (const auto& name : fmNames)
121 {
122 fieldMults_[i++] = MDField<const ScalarT, Cell, IP>(name, ir.dl_scalar);
123 this->addDependentField(fieldMults_[i - 1]);
124 } // end loop over the field multipliers
125
126 // Set the name of this object.
127 string n("Integrator_BasisTimesScalar (");
129 n += "CONTRIBUTES";
130 else // if (evalStyle_ == EvaluatorStyle::EVALUATES)
131 n += "EVALUATES";
132 n += "): " + field_.fieldTag().name();
133 this->setName(n);
134 } // end of Main Constructor
135
137 //
138 // ParameterList Constructor
139 //
141 template<typename EvalT, typename Traits>
144 const Teuchos::ParameterList& p)
145 :
148 p.get<std::string>("Residual Name"),
149 p.get<std::string>("Value Name"),
150 (*p.get<Teuchos::RCP<panzer::BasisIRLayout>>("Basis")),
151 (*p.get<Teuchos::RCP<panzer::IntegrationRule>>("IR")),
152 p.get<double>("Multiplier"),
153 p.isType<Teuchos::RCP<const std::vector<std::string>>>
154 ("Field Multipliers") ?
155 (*p.get<Teuchos::RCP<const std::vector<std::string>>>
156 ("Field Multipliers")) : std::vector<std::string>())
157 {
158 using Teuchos::ParameterList;
159 using Teuchos::RCP;
160
161 // Ensure that the input ParameterList didn't contain any bogus entries.
162 RCP<ParameterList> validParams = this->getValidParameters();
163 p.validateParameters(*validParams);
164 } // end of ParameterList Constructor
165
167 //
168 // postRegistrationSetup()
169 //
171 template<typename EvalT, typename Traits>
172 void
175 typename Traits::SetupData sd,
177 {
179 using std::size_t;
180
181 auto kokkosFieldMults_h = Kokkos::create_mirror_view(kokkosFieldMults_);
182
183 // Get the PHX::Views of the field multipliers.
184 for (size_t i(0); i < fieldMults_.size(); ++i)
185 kokkosFieldMults_h(i) = fieldMults_[i].get_static_view();
186
187 Kokkos::deep_copy(kokkosFieldMults_, kokkosFieldMults_h);
188
189 // Determine the index in the Workset bases for our particular basis name.
190 basisIndex_ = getBasisIndex(basisName_, (*sd.worksets_)[0], this->wda);
191
192 // Allocate temporary
193 if (Sacado::IsADType<ScalarT>::value) {
194 const auto fadSize = Kokkos::dimension_scalar(field_.get_view());
195 tmp_ = PHX::View<ScalarT*>("IntegratorBasisTimesScalar::tmp_",field_.extent(0),fadSize);
196 if (fieldMults_.size() > 1)
197 tmp2_ = PHX::View<ScalarT*>("IntegratorBasisTimesScalar::tmp_",field_.extent(0),fadSize);
198 } else {
199 tmp_ = PHX::View<ScalarT*>("IntegratorBasisTimesScalar::tmp_",field_.extent(0));
200 if (fieldMults_.size() > 1)
201 tmp2_ = PHX::View<ScalarT*>("IntegratorBasisTimesScalar::tmp_",field_.extent(0));
202 }
203 } // end of postRegistrationSetup()
204
206 //
207 // operator()()
208 //
210 template<typename EvalT, typename Traits>
211 template<int NUM_FIELD_MULT>
212 KOKKOS_INLINE_FUNCTION
213 void
216 const FieldMultTag<NUM_FIELD_MULT>& /* tag */,
217 const size_t& cell) const
218 {
220
221 // Initialize the evaluated field.
222 const int numQP(scalar_.extent(1)), numBases(basis_.extent(1));
223 if (evalStyle_ == EvaluatorStyle::EVALUATES)
224 for (int basis(0); basis < numBases; ++basis)
225 field_(cell, basis) = 0.0;
226
227 // The following if-block is for the sake of optimization depending on the
228 // number of field multipliers.
229 if (NUM_FIELD_MULT == 0)
230 {
231 // Loop over the quadrature points, scale the integrand by the
232 // multiplier, and then perform the actual integration, looping over the
233 // bases.
234 for (int qp(0); qp < numQP; ++qp)
235 {
236 tmp_(cell) = multiplier_ * scalar_(cell, qp);
237 for (int basis(0); basis < numBases; ++basis)
238 field_(cell, basis) += basis_(cell, basis, qp) * tmp_(cell);
239 } // end loop over the quadrature points
240 }
241 else if (NUM_FIELD_MULT == 1)
242 {
243 // Loop over the quadrature points, scale the integrand by the multiplier
244 // and the single field multiplier, and then perform the actual
245 // integration, looping over the bases.
246 for (int qp(0); qp < numQP; ++qp)
247 {
248 tmp_(cell) = multiplier_ * scalar_(cell, qp) * kokkosFieldMults_(0)(cell, qp);
249 for (int basis(0); basis < numBases; ++basis)
250 field_(cell, basis) += basis_(cell, basis, qp) * tmp_(cell);
251 } // end loop over the quadrature points
252 }
253 else
254 {
255 // Loop over the quadrature points and pre-multiply all the field
256 // multipliers together, scale the integrand by the multiplier and
257 // the combination of the field multipliers, and then perform the actual
258 // integration, looping over the bases.
259 const int numFieldMults(kokkosFieldMults_.extent(0));
260 for (int qp(0); qp < numQP; ++qp)
261 {
262 tmp2_(cell) = 1.0;
263 for (int fm(0); fm < numFieldMults; ++fm)
264 tmp2_(cell) *= kokkosFieldMults_(fm)(cell, qp);
265 tmp_(cell) = multiplier_ * scalar_(cell, qp) * tmp2_(cell);
266 for (int basis(0); basis < numBases; ++basis)
267 field_(cell, basis) += basis_(cell, basis, qp) * tmp_(cell);
268 } // end loop over the quadrature points
269 } // end if (NUM_FIELD_MULT == something)
270 } // end of operator()()
271
273 //
274 // evaluateFields()
275 //
277 template<typename EvalT, typename Traits>
278 void
281 typename Traits::EvalData workset)
282 {
283 using Kokkos::parallel_for;
284 using Kokkos::RangePolicy;
285
286 // Grab the basis information.
287 basis_ = this->wda(workset).bases[basisIndex_]->weighted_basis_scalar;
288
289 // The following if-block is for the sake of optimization depending on the
290 // number of field multipliers. The parallel_fors will loop over the cells
291 // in the Workset and execute operator()() above.
292 if (fieldMults_.size() == 0)
293 parallel_for(RangePolicy<FieldMultTag<0>>(0, workset.num_cells), *this);
294 else if (fieldMults_.size() == 1)
295 parallel_for(RangePolicy<FieldMultTag<1>>(0, workset.num_cells), *this);
296 else
297 parallel_for(RangePolicy<FieldMultTag<-1>>(0, workset.num_cells), *this);
298 } // end of evaluateFields()
299
301 //
302 // getValidParameters()
303 //
305 template<typename EvalT, typename TRAITS>
306 Teuchos::RCP<Teuchos::ParameterList>
309 {
312 using std::string;
313 using std::vector;
314 using Teuchos::ParameterList;
315 using Teuchos::RCP;
316 using Teuchos::rcp;
317
318 // Create a ParameterList with all the valid keys we support.
319 RCP<ParameterList> p = rcp(new ParameterList);
320 p->set<string>("Residual Name", "?");
321 p->set<string>("Value Name", "?");
322 RCP<BasisIRLayout> basis;
323 p->set("Basis", basis);
324 RCP<IntegrationRule> ir;
325 p->set("IR", ir);
326 p->set<double>("Multiplier", 1.0);
327 RCP<const vector<string>> fms;
328 p->set("Field Multipliers", fms);
329 return p;
330 } // end of getValidParameters()
331
332} // end of namespace panzer
333
334#endif // PANZER_INTEGRATOR_BASISTIMESSCALAR_IMPL_HPP
panzer::EvaluatorStyle evalStyle
The EvaluatorStyle of the parent Integrator_CurlBasisDotVector object.
double multiplier
The scalar multiplier out in front of the integral ( ).
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis>
Teuchos::RCP< const PureBasis > getBasis() const
void evaluateFields(typename Traits::EvalData workset)
Evaluate Fields.
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 ( ,...
KOKKOS_INLINE_FUNCTION void operator()(const FieldMultTag< NUM_FIELD_MULT > &tag, const std::size_t &cell) const
Perform the integration.
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...
PHX::MDField< const ScalarT, panzer::Cell, panzer::IP > scalar_
A field representing the scalar function we're integrating ( ).
void postRegistrationSetup(typename Traits::SetupData sd, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
const panzer::EvaluatorStyle evalStyle_
An enum determining the behavior of this Evaluator.
Integrator_BasisTimesScalar(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< 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 > 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_