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// Panzer: A partial differential equation assembly
4// engine for strongly coupled complex multiphysics systems
5//
6// Copyright 2011 NTESS and the Panzer contributors.
7// SPDX-License-Identifier: BSD-3-Clause
8// *****************************************************************************
9// @HEADER
10
11#ifndef __Panzer_Integerator_BasisTimesTensorTimesVector_impl_hpp__
12#define __Panzer_Integerator_BasisTimesTensorTimesVector_impl_hpp__
13
15//
16// Include Files
17//
19
20// Panzer
24
25namespace panzer
26{
28 //
29 // Main Constructor
30 //
32 template<typename EvalT, typename Traits>
36 const std::string& resName,
37 const std::string& valName,
38 const panzer::BasisIRLayout& basis,
40 const std::string& tensorName)
41 :
42 evalStyle_(evalStyle),
43 useDescriptors_(false),
44 basisName_(basis.name())
45 {
46 using PHX::View;
47 using panzer::BASIS;
48 using panzer::Cell;
50 using panzer::IP;
52 using PHX::MDField;
53 using PHX::print;
54 using std::invalid_argument;
55 using std::logic_error;
56 using std::string;
57 using Teuchos::RCP;
58
59 // Ensure the input makes sense.
60 TEUCHOS_TEST_FOR_EXCEPTION(resName == "", invalid_argument, "Error: " \
61 "Integrator_BasisTimesTensorTimesVector called with an empty residual name.")
62 TEUCHOS_TEST_FOR_EXCEPTION(valName == "", invalid_argument, "Error: " \
63 "Integrator_BasisTimesTensorTimesVector called with an empty value name.")
64 RCP<const PureBasis> tmpBasis = basis.getBasis();
65 TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->isVectorBasis(), logic_error,
66 "Error: Integrator_BasisTimesTensorTimesVector: Basis of type \""
67 << tmpBasis->name() << "\" is not a vector basis.");
68 TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->requiresOrientations(),
69 logic_error, "Integrator_BasisTimesTensorTimesVector: Basis of type \""
70 << tmpBasis->name() << "\" does not require orientations. This seems " \
71 "very strange, so I'm failing.");
72
73 // Create the field for the vector-valued quantity we're integrating.
74 vector_ = MDField<const ScalarT, Cell, IP, Dim>(valName, ir.dl_vector);
75 this->addDependentField(vector_);
76
77 // Create the field that we're either contributing to or evaluating
78 // (storing).
79 field_ = MDField<ScalarT, Cell, BASIS>(resName, basis.functional);
81 this->addContributedField(field_);
82 else // if (evalStyle == EvaluatorStyle::EVALUATES)
83 this->addEvaluatedField(field_);
84
85 // Add the tensor.
86 tensor_ = MDField<const ScalarT, Cell, IP, Dim, Dim>(tensorName, ir.dl_tensor);
87 this->addDependentField(tensor_);
88
89 // Set the name of this object.
90 string n("Integrator_BasisTimesTensorTimesVector (");
92 n += "Cont";
93 else // if (evalStyle == EvaluatorStyle::EVALUATES)
94 n += "Eval";
95 n += ", " + print<EvalT>() + "): " + field_.fieldTag().name();
96 this->setName(n);
97 } // end of Main Constructor
98
100 //
101 // Descriptor Constructor
102 //
104 template<typename EvalT, typename Traits>
108 const PHX::FieldTag& resTag,
109 const PHX::FieldTag& valTag,
110 const BasisDescriptor& bd,
111 const IntegrationDescriptor& id,
112 const PHX::FieldTag& tensorTag)
113 :
114 evalStyle_(evalStyle),
115 useDescriptors_(true),
116 bd_(bd),
117 id_(id)
118 {
119 using PHX::View;
120 using panzer::BASIS;
121 using panzer::Cell;
123 using panzer::IP;
124 using panzer::PureBasis;
125 using PHX::MDField;
126 using PHX::print;
127 using std::invalid_argument;
128 using std::logic_error;
129 using std::string;
130 using Teuchos::RCP;
131
132 // Ensure the input makes sense.
133 TEUCHOS_TEST_FOR_EXCEPTION(not ((bd_.getType() == "HCurl") or
134 (bd_.getType() == "HDiv")), logic_error, "Error: " \
135 "Integrator_BasisTimesTensorTimesVector: Basis of type \"" << bd_.getType()
136 << "\" is not a vector basis.")
137
138 // Create the field for the vector-valued quantity we're integrating.
139 vector_ = valTag;
140 this->addDependentField(vector_);
141
142 // Create the field that we're either contributing to or evaluating
143 // (storing).
144 field_ = resTag;
146 this->addContributedField(field_);
147 else // if (evalStyle == EvaluatorStyle::EVALUATES)
148 this->addEvaluatedField(field_);
149
150 // Add the tensor.
151 tensor_ = tensorTag;
152 this->addDependentField(tensor_);
153
154 // Set the name of this object.
155 string n("Integrator_BasisTimesTensorTimesVector (");
157 n += "Cont";
158 else // if (evalStyle == EvaluatorStyle::EVALUATES)
159 n += "Eval";
160 n += ", " + print<EvalT>() + "): " + field_.fieldTag().name();
161 this->setName(n);
162 } // end of Main Constructor
163
165 //
166 // ParameterList Constructor
167 //
169 template<typename EvalT, typename Traits>
172 const Teuchos::ParameterList& p)
173 :
176 p.get<std::string>("Residual Name"),
177 p.get<std::string>("Value Name"),
178 (*p.get<Teuchos::RCP<panzer::BasisIRLayout>>("Basis")),
179 (*p.get<Teuchos::RCP<panzer::IntegrationRule>>("IR")),
180 p.get<std::string>("Tensor Name"))
181 {
182 using Teuchos::ParameterList;
183 using Teuchos::RCP;
184
185 // Ensure that the input ParameterList didn't contain any bogus entries.
186 RCP<ParameterList> validParams = this->getValidParameters();
187 p.validateParameters(*validParams);
188 } // end of ParameterList Constructor
189
191 //
192 // postRegistrationSetup()
193 //
195 template<typename EvalT, typename Traits>
196 void
199 typename Traits::SetupData sd,
201 {
203 using std::size_t;
204
205 // Get the PHX::View of the tensor.
206 kokkosTensor_ = tensor_.get_static_view();
207
208 // Determine the number of quadrature points and the dimensionality of the
209 // vector that we're integrating.
210 numQP_ = vector_.extent(1);
211 numDim_ = vector_.extent(2);
212
213 // Determine the index in the Workset bases for our particular basis name.
214 if (not useDescriptors_)
215 basisIndex_ = getBasisIndex(basisName_, (*sd.worksets_)[0], this->wda);
216
217 if (Sacado::IsADType<ScalarT>::value) {
218 const auto fadSize = Kokkos::dimension_scalar(field_.get_view());
219 tmp_ = PHX::View<ScalarT*>("GradBasisDotVector::tmp_",field_.extent(0),fadSize);
220 } else {
221 tmp_ = PHX::View<ScalarT*>("GradBasisDotVector::tmp_",field_.extent(0));
222 }
223
224 } // end of postRegistrationSetup()
225
227 //
228 // operator()()
229 //
231 template<typename EvalT, typename Traits>
232 KOKKOS_INLINE_FUNCTION
233 void
236 const size_t& cell) const
237 {
239
240 // Initialize the evaluated field.
241 const int numBases(basis_.extent(1));
242 if (evalStyle_ == EvaluatorStyle::EVALUATES)
243 for (int basis(0); basis < numBases; ++basis)
244 field_(cell, basis) = 0.0;
245
246 // Loop over the quadrature points and dimensions of our vector fields,
247 // scale the integrand by the tensor,
248 // and then perform the actual integration, looping over the bases.
249 for (int qp(0); qp < numQP_; ++qp)
250 {
251 for (int dim(0); dim < numDim_; ++dim)
252 {
253 tmp_(cell) = 0.0;
254 for (int dim2(0); dim2 < numDim_; ++dim2)
255 tmp_(cell) += kokkosTensor_(cell, qp, dim, dim2) * vector_(cell, qp, dim2);
256 for (int basis(0); basis < numBases; ++basis)
257 field_(cell, basis) += basis_(cell, basis, qp, dim) * tmp_(cell);
258 } // end loop over the dimensions of the vector field
259 } // end loop over the quadrature points
260 } // end of operator()()
261
263 //
264 // evaluateFields()
265 //
267 template<typename EvalT, typename Traits>
268 void
271 typename Traits::EvalData workset)
272 {
273 using Kokkos::parallel_for;
274 using Kokkos::RangePolicy;
275
276 // Grab the basis information.
277 const panzer::BasisValues2<double>& bv = useDescriptors_ ?
278 this->wda(workset).getBasisValues(bd_,id_) :
279 *this->wda(workset).bases[basisIndex_];
281 basis_ = useDescriptors_ ? bv.getVectorBasisValues(true) : Array(bv.weighted_basis_vector);
282
283 parallel_for(RangePolicy<>(0, workset.num_cells), *this);
284 } // end of evaluateFields()
285
287 //
288 // getValidParameters()
289 //
291 template<typename EvalT, typename Traits>
292 Teuchos::RCP<Teuchos::ParameterList>
295 {
298 using std::string;
299 using std::vector;
300 using Teuchos::ParameterList;
301 using Teuchos::RCP;
302 using Teuchos::rcp;
303
304 // Create a ParameterList with all the valid keys we support.
305 RCP<ParameterList> p = rcp(new ParameterList);
306 p->set<string>("Residual Name", "?");
307 p->set<string>("Value Name", "?");
308 RCP<BasisIRLayout> basis;
309 p->set("Basis", basis);
310 RCP<IntegrationRule> ir;
311 p->set("IR", ir);
312 p->set<std::string>("Tensor Name", "?");
313 return p;
314 } // end of getValidParameters()
315
316} // end of namespace panzer
317
318#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_