Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_Integrator_BasisTimesTensorTimesVector_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_BasisTimesTensorTimesVector_impl_hpp__
44#define __Panzer_Integerator_BasisTimesTensorTimesVector_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 std::string& tensorName)
73 :
74 evalStyle_(evalStyle),
75 useDescriptors_(false),
76 basisName_(basis.name())
77 {
78 using PHX::View;
79 using panzer::BASIS;
80 using panzer::Cell;
82 using panzer::IP;
84 using PHX::MDField;
85 using PHX::print;
86 using std::invalid_argument;
87 using std::logic_error;
88 using std::string;
89 using Teuchos::RCP;
90
91 // Ensure the input makes sense.
92 TEUCHOS_TEST_FOR_EXCEPTION(resName == "", invalid_argument, "Error: " \
93 "Integrator_BasisTimesTensorTimesVector called with an empty residual name.")
94 TEUCHOS_TEST_FOR_EXCEPTION(valName == "", invalid_argument, "Error: " \
95 "Integrator_BasisTimesTensorTimesVector called with an empty value name.")
96 RCP<const PureBasis> tmpBasis = basis.getBasis();
97 TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->isVectorBasis(), logic_error,
98 "Error: Integrator_BasisTimesTensorTimesVector: Basis of type \""
99 << tmpBasis->name() << "\" is not a vector basis.");
100 TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->requiresOrientations(),
101 logic_error, "Integrator_BasisTimesTensorTimesVector: Basis of type \""
102 << tmpBasis->name() << "\" does not require orientations. This seems " \
103 "very strange, so I'm failing.");
104
105 // Create the field for the vector-valued quantity we're integrating.
106 vector_ = MDField<const ScalarT, Cell, IP, Dim>(valName, ir.dl_vector);
107 this->addDependentField(vector_);
108
109 // Create the field that we're either contributing to or evaluating
110 // (storing).
111 field_ = MDField<ScalarT, Cell, BASIS>(resName, basis.functional);
113 this->addContributedField(field_);
114 else // if (evalStyle == EvaluatorStyle::EVALUATES)
115 this->addEvaluatedField(field_);
116
117 // Add the tensor.
118 tensor_ = MDField<const ScalarT, Cell, IP, Dim, Dim>(tensorName, ir.dl_tensor);
119 this->addDependentField(tensor_);
120
121 // Set the name of this object.
122 string n("Integrator_BasisTimesTensorTimesVector (");
124 n += "Cont";
125 else // if (evalStyle == EvaluatorStyle::EVALUATES)
126 n += "Eval";
127 n += ", " + print<EvalT>() + "): " + field_.fieldTag().name();
128 this->setName(n);
129 } // end of Main Constructor
130
132 //
133 // Descriptor Constructor
134 //
136 template<typename EvalT, typename Traits>
140 const PHX::FieldTag& resTag,
141 const PHX::FieldTag& valTag,
142 const BasisDescriptor& bd,
143 const IntegrationDescriptor& id,
144 const PHX::FieldTag& tensorTag)
145 :
146 evalStyle_(evalStyle),
147 useDescriptors_(true),
148 bd_(bd),
149 id_(id)
150 {
151 using PHX::View;
152 using panzer::BASIS;
153 using panzer::Cell;
155 using panzer::IP;
156 using panzer::PureBasis;
157 using PHX::MDField;
158 using PHX::print;
159 using std::invalid_argument;
160 using std::logic_error;
161 using std::string;
162 using Teuchos::RCP;
163
164 // Ensure the input makes sense.
165 TEUCHOS_TEST_FOR_EXCEPTION(not ((bd_.getType() == "HCurl") or
166 (bd_.getType() == "HDiv")), logic_error, "Error: " \
167 "Integrator_BasisTimesTensorTimesVector: Basis of type \"" << bd_.getType()
168 << "\" is not a vector basis.")
169
170 // Create the field for the vector-valued quantity we're integrating.
171 vector_ = valTag;
172 this->addDependentField(vector_);
173
174 // Create the field that we're either contributing to or evaluating
175 // (storing).
176 field_ = resTag;
178 this->addContributedField(field_);
179 else // if (evalStyle == EvaluatorStyle::EVALUATES)
180 this->addEvaluatedField(field_);
181
182 // Add the tensor.
183 tensor_ = tensorTag;
184 this->addDependentField(tensor_);
185
186 // Set the name of this object.
187 string n("Integrator_BasisTimesTensorTimesVector (");
189 n += "Cont";
190 else // if (evalStyle == EvaluatorStyle::EVALUATES)
191 n += "Eval";
192 n += ", " + print<EvalT>() + "): " + field_.fieldTag().name();
193 this->setName(n);
194 } // end of Main Constructor
195
197 //
198 // ParameterList Constructor
199 //
201 template<typename EvalT, typename Traits>
204 const Teuchos::ParameterList& p)
205 :
208 p.get<std::string>("Residual Name"),
209 p.get<std::string>("Value Name"),
210 (*p.get<Teuchos::RCP<panzer::BasisIRLayout>>("Basis")),
211 (*p.get<Teuchos::RCP<panzer::IntegrationRule>>("IR")),
212 p.get<std::string>("Tensor Name"))
213 {
214 using Teuchos::ParameterList;
215 using Teuchos::RCP;
216
217 // Ensure that the input ParameterList didn't contain any bogus entries.
218 RCP<ParameterList> validParams = this->getValidParameters();
219 p.validateParameters(*validParams);
220 } // end of ParameterList Constructor
221
223 //
224 // postRegistrationSetup()
225 //
227 template<typename EvalT, typename Traits>
228 void
231 typename Traits::SetupData sd,
233 {
235 using std::size_t;
236
237 // Get the PHX::View of the tensor.
238 kokkosTensor_ = tensor_.get_static_view();
239
240 // Determine the number of quadrature points and the dimensionality of the
241 // vector that we're integrating.
242 numQP_ = vector_.extent(1);
243 numDim_ = vector_.extent(2);
244
245 // Determine the index in the Workset bases for our particular basis name.
246 if (not useDescriptors_)
247 basisIndex_ = getBasisIndex(basisName_, (*sd.worksets_)[0], this->wda);
248
249 if (Sacado::IsADType<ScalarT>::value) {
250 const auto fadSize = Kokkos::dimension_scalar(field_.get_view());
251 tmp_ = PHX::View<ScalarT*>("GradBasisDotVector::tmp_",field_.extent(0),fadSize);
252 } else {
253 tmp_ = PHX::View<ScalarT*>("GradBasisDotVector::tmp_",field_.extent(0));
254 }
255
256 } // end of postRegistrationSetup()
257
259 //
260 // operator()()
261 //
263 template<typename EvalT, typename Traits>
264 KOKKOS_INLINE_FUNCTION
265 void
268 const size_t& cell) const
269 {
271
272 // Initialize the evaluated field.
273 const int numBases(basis_.extent(1));
274 if (evalStyle_ == EvaluatorStyle::EVALUATES)
275 for (int basis(0); basis < numBases; ++basis)
276 field_(cell, basis) = 0.0;
277
278 // Loop over the quadrature points and dimensions of our vector fields,
279 // scale the integrand by the tensor,
280 // and then perform the actual integration, looping over the bases.
281 for (int qp(0); qp < numQP_; ++qp)
282 {
283 for (int dim(0); dim < numDim_; ++dim)
284 {
285 tmp_(cell) = 0.0;
286 for (int dim2(0); dim2 < numDim_; ++dim2)
287 tmp_(cell) += kokkosTensor_(cell, qp, dim, dim2) * vector_(cell, qp, dim2);
288 for (int basis(0); basis < numBases; ++basis)
289 field_(cell, basis) += basis_(cell, basis, qp, dim) * tmp_(cell);
290 } // end loop over the dimensions of the vector field
291 } // end loop over the quadrature points
292 } // end of operator()()
293
295 //
296 // evaluateFields()
297 //
299 template<typename EvalT, typename Traits>
300 void
303 typename Traits::EvalData workset)
304 {
305 using Kokkos::parallel_for;
306 using Kokkos::RangePolicy;
307
308 // Grab the basis information.
309 const panzer::BasisValues2<double>& bv = useDescriptors_ ?
310 this->wda(workset).getBasisValues(bd_,id_) :
311 *this->wda(workset).bases[basisIndex_];
313 basis_ = useDescriptors_ ? bv.getVectorBasisValues(true) : Array(bv.weighted_basis_vector);
314
315 parallel_for(RangePolicy<>(0, workset.num_cells), *this);
316 } // end of evaluateFields()
317
319 //
320 // getValidParameters()
321 //
323 template<typename EvalT, typename Traits>
324 Teuchos::RCP<Teuchos::ParameterList>
327 {
330 using std::string;
331 using std::vector;
332 using Teuchos::ParameterList;
333 using Teuchos::RCP;
334 using Teuchos::rcp;
335
336 // Create a ParameterList with all the valid keys we support.
337 RCP<ParameterList> p = rcp(new ParameterList);
338 p->set<string>("Residual Name", "?");
339 p->set<string>("Value Name", "?");
340 RCP<BasisIRLayout> basis;
341 p->set("Basis", basis);
342 RCP<IntegrationRule> ir;
343 p->set("IR", ir);
344 p->set<std::string>("Tensor Name", "?");
345 return p;
346 } // end of getValidParameters()
347
348} // end of namespace panzer
349
350#endif // __Panzer_Integerator_BasisTimesTensorTimesVector_impl_hpp__
panzer::EvaluatorStyle evalStyle
The EvaluatorStyle of the parent Integrator_CurlBasisDotVector object.
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.
void postRegistrationSetup(typename Traits::SetupData sd, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
KOKKOS_INLINE_FUNCTION void operator()(const std::size_t &cell) const
Perform the integration.
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_BasisTimesTensorTimesVector(const panzer::EvaluatorStyle &evalStyle, const std::string &resName, const std::string &valName, const panzer::BasisIRLayout &basis, const panzer::IntegrationRule &ir, const std::string &tensorName)
Main Constructor.
void evaluateFields(typename Traits::EvalData workset)
Evaluate Fields.
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.
PHX::MDField< const ScalarT, panzer::Cell, panzer::IP, panzer::Dim, panzer::Dim > tensor_
The tensor field( ).
Teuchos::RCP< PHX::DataLayout > dl_vector
Data layout for vector fields.
Teuchos::RCP< PHX::DataLayout > dl_tensor
Data layout for rank-2 tensor 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.
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_