11#ifndef Panzer_Integrator_CurlBasisDotVector_impl_hpp
12#define Panzer_Integrator_CurlBasisDotVector_impl_hpp
21#include "Intrepid2_FunctionSpaceTools.hpp"
30#include "Phalanx_KokkosDeviceTypes.hpp"
39 template<
typename EvalT,
typename Traits>
43 const std::string& resName,
44 const std::string& valName,
48 const std::vector<std::string>& fmNames
52 useDescriptors_(false),
54 basisName_(basis.name())
64 using std::invalid_argument;
65 using std::logic_error;
70 TEUCHOS_TEST_FOR_EXCEPTION(resName ==
"", invalid_argument,
"Error: " \
71 "Integrator_CurlBasisDotVector called with an empty residual name.")
72 TEUCHOS_TEST_FOR_EXCEPTION(valName ==
"", invalid_argument,
"Error: " \
73 "Integrator_CurlBasisDotVector called with an empty value name.")
74 RCP<const PureBasis> tmpBasis = basis.
getBasis();
75 TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->isVectorBasis(), logic_error,
76 "Error: Integrator_CurlBasisDotVector: Basis of type \""
77 << tmpBasis->name() <<
"\" is not a vector basis.")
78 TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->requiresOrientations(),
79 logic_error,
"Error: Integrator_CurlBasisDotVector: Basis of type \""
80 << tmpBasis->name() <<
"\" does not require orientations, though it " \
81 "should for its use in this Evaluator.")
82 TEUCHOS_TEST_FOR_EXCEPTION(not tmpBasis->supportsCurl(), logic_error,
83 "Error: Integrator_CurlBasisDotVector: Basis of type \""
84 << tmpBasis->name() <<
"\" does not support curl.")
85 TEUCHOS_TEST_FOR_EXCEPTION(not ((tmpBasis->dimension() == 2) or
86 (tmpBasis->dimension() == 3)), logic_error,
87 "Error: Integrator_CurlBasisDotVector requires either a two- or " \
88 "three-dimensional basis. The basis \"" << tmpBasis->name()
110 this->addContributedField(
field_);
112 this->addEvaluatedField(
field_);
117 kokkosFieldMults_ = View<Kokkos::View<const ScalarT**, typename PHX::DevLayout<ScalarT>::type, Kokkos::MemoryUnmanaged>*>(
118 "CurlBasisDotVector::KokkosFieldMultipliers", fmNames.size());
119 for (
const auto& name : fmNames)
126 string n(
"Integrator_CurlBasisDotVector (");
131 n +=
"): " +
field_.fieldTag().name();
140 template<
typename EvalT,
typename Traits>
143 const Teuchos::ParameterList& p)
147 p.get<
std::string>(
"Residual Name"),
148 p.get<
std::string>(
"Value Name"),
151 p.get<double>(
"Multiplier"),
153 (
"Field Multipliers") ?
155 (
"Field Multipliers")) :
std::vector<
std::string>())
157 using Teuchos::ParameterList;
162 p.validateParameters(*validParams);
170 template<
typename EvalT,
typename TRAITS>
174 const PHX::FieldTag& resTag,
175 const PHX::FieldTag& valTag,
180 const std::vector<PHX::FieldTag>& multipliers
184 useDescriptors_(true),
192 using std::logic_error;
196 TEUCHOS_TEST_FOR_EXCEPTION(
bd_.
getType() !=
"HCurl", logic_error,
197 "Error: Integrator_CurlBasisDotVector: Basis of type \""
198 <<
bd_.
getType() <<
"\" does not support curl.")
199 TEUCHOS_TEST_FOR_EXCEPTION(not ((spaceDim == 2) or (spaceDim == 3)),
200 logic_error,
"Error: Integrator_CurlBasisDotVector works on either " \
201 "two- or three-dimensional problems. You've provided spaceDim = "
220 this->addContributedField(
field_);
222 this->addEvaluatedField(
field_);
227 kokkosFieldMults_ = View<Kokkos::View<const ScalarT**, typename PHX::DevLayout<ScalarT>::type, Kokkos::MemoryUnmanaged>*>(
228 "CurlBasisDotVector::KokkosFieldMultipliers", multipliers.size());
229 for (
const auto& fm : multipliers)
236 string n(
"Integrator_CurlBasisDotVector (");
241 n +=
"): " +
field_.fieldTag().name();
250 template<
typename EvalT,
typename Traits>
262 auto kokkosFieldMults_h = Kokkos::create_mirror_view(kokkosFieldMults_);
263 for (
size_t i(0); i < fieldMults_.size(); ++i)
264 kokkosFieldMults_h(i) = fieldMults_[i].get_static_view();
265 Kokkos::deep_copy(kokkosFieldMults_, kokkosFieldMults_h);
268 if (not useDescriptors_)
274 fm.template getKokkosExtendedDataTypeDimensions<EvalT>(),
true);
277 "Integrator_CurlBasisDotVector: 2-D Result", vector2D_.extent(0),
278 vector2D_.extent(1));
281 "Integrator_CurlBasisDotVector: 3-D Result", vector3D_.extent(0),
282 vector3D_.extent(1), 3);
303 template<
typename ScalarT>
312 struct ScalarMultiplierTag
320 struct FieldMultiplierTag
335 KOKKOS_INLINE_FUNCTION
336 void operator()(
const ScalarMultiplierTag,
const unsigned cell)
const
338 int numQP(
result.extent(1));
339 for (
int qp(0); qp < numQP; ++qp)
355 KOKKOS_INLINE_FUNCTION
356 void operator()(
const FieldMultiplierTag,
const unsigned cell)
const
358 int numQP(
result.extent(1));
359 for (
int qp(0); qp < numQP; ++qp)
367 using execution_space = PHX::Device;
373 PHX::MDField<ScalarT, panzer::Cell, panzer::IP>
result;
385 PHX::MDField<const ScalarT, panzer::Cell, panzer::IP>
vectorField;
391 PHX::MDField<const ScalarT, panzer::Cell, panzer::IP>
fieldMult;
405 template<
typename ScalarT>
414 struct ScalarMultiplierTag
422 struct FieldMultiplierTag
437 KOKKOS_INLINE_FUNCTION
438 void operator()(
const ScalarMultiplierTag,
const unsigned cell)
const
441 for (
int qp(0); qp < numQP; ++qp)
442 for (
int dim(0); dim < numDim; ++dim)
458 KOKKOS_INLINE_FUNCTION
459 void operator()(
const FieldMultiplierTag,
const unsigned cell)
const
462 for (
int qp(0); qp < numQP; ++qp)
463 for (
int dim(0); dim < numDim; ++dim)
471 using execution_space = PHX::Device;
477 PHX::MDField<ScalarT, panzer::Cell, panzer::IP, panzer::Dim>
result;
489 PHX::MDField<const ScalarT, panzer::Cell, panzer::IP, panzer::Dim>
496 PHX::MDField<const ScalarT, panzer::Cell, panzer::IP>
fieldMult;
510 template<
typename ScalarT,
int spaceDim>
527 KOKKOS_INLINE_FUNCTION
528 void operator()(
const unsigned cell)
const
533 for (
int basis(0); basis < numBases; ++basis)
536 field(cell, basis) = 0.0;
537 for (
int qp(0); qp < numQP; ++qp)
538 for (
int dim(0); dim < spaceDim; ++dim)
548 using execution_space = PHX::Device;
553 PHX::MDField<ScalarT, panzer::Cell, panzer::IP, panzer::Dim>
result;
559 PHX::MDField<ScalarT, panzer::Cell, panzer::BASIS>
field;
585 template<
typename ScalarT>
601 KOKKOS_INLINE_FUNCTION
602 void operator()(
const unsigned cell)
const
607 for (
int basis(0); basis < numBases; ++basis)
610 field(cell, basis) = 0.0;
611 for (
int qp(0); qp < numQP; ++qp)
621 using execution_space = PHX::Device;
626 PHX::MDField<ScalarT, panzer::Cell, panzer::IP>
result;
632 PHX::MDField<ScalarT, panzer::Cell, panzer::BASIS>
field;
637 PHX::MDField<const double, panzer::Cell, panzer::BASIS, panzer::IP>
654 template<
typename EvalT,
typename Traits>
660 using Kokkos::parallel_for;
661 using Kokkos::RangePolicy;
670 *this->wda(workset).bases[basisIndex_];
677 using PreMultiply = PreMultiply2D<ScalarT>;
678 using ScalarMultiplierTag =
typename PreMultiply::ScalarMultiplierTag;
679 using FieldMultiplierTag =
typename PreMultiply::FieldMultiplierTag;
680 PreMultiply preMultiply;
681 preMultiply.result = result2D_;
682 preMultiply.multiplier = multiplier_;
683 preMultiply.vectorField = vector2D_;
687 parallel_for(RangePolicy<Device, ScalarMultiplierTag>(0,
692 for (
const auto&
field : fieldMults_)
694 preMultiply.fieldMult =
field;
695 parallel_for(RangePolicy<Device, FieldMultiplierTag>(0,
700 Integrate2D<ScalarT> integrate;
701 integrate.result = result2D_;
702 integrate.field = field_;
705 integrate.evalStyle = evalStyle_;
706 parallel_for(workset.
num_cells, integrate);
712 using PreMultiply = PreMultiply3D<ScalarT>;
713 using ScalarMultiplierTag =
typename PreMultiply::ScalarMultiplierTag;
714 using FieldMultiplierTag =
typename PreMultiply::FieldMultiplierTag;
715 PreMultiply preMultiply;
716 preMultiply.result = result3D_;
717 preMultiply.multiplier = multiplier_;
718 preMultiply.vectorField = vector3D_;
722 parallel_for(RangePolicy<Device, ScalarMultiplierTag>(0,
727 for (
const auto&
field : fieldMults_)
729 preMultiply.fieldMult =
field;
730 parallel_for(RangePolicy<Device, FieldMultiplierTag>(0,
735 Integrate3D<ScalarT, 3> integrate;
736 integrate.result = result3D_;
737 integrate.field = field_;
740 integrate.evalStyle = evalStyle_;
741 parallel_for(workset.
num_cells, integrate);
750 template<
typename EvalT,
typename TRAITS>
751 Teuchos::RCP<Teuchos::ParameterList>
759 using Teuchos::ParameterList;
763 RCP<ParameterList> p = rcp(
new ParameterList);
764 p->set<
string>(
"Residual Name",
"?");
765 p->set<
string>(
"Value Name",
"?");
766 RCP<BasisIRLayout> basis;
767 p->set(
"Basis", basis);
768 RCP<IntegrationRule> ir;
770 p->set<
double>(
"Multiplier", 1.0);
771 RCP<const vector<string>> fms;
772 p->set(
"Field Multipliers", fms);
PHX::MDField< const double, panzer::Cell, panzer::BASIS, panzer::IP, panzer::Dim > weightedCurlBasis
The vector basis information necessary for integration.
PHX::MDField< const ScalarT, panzer::Cell, panzer::IP > vectorField
A field representing the vector-valued function we're integrating ( ).
PHX::MDField< ScalarT, panzer::Cell, panzer::IP > result
A field that will be used to build up the result of the integral we're performing.
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.
PHX::MDField< const ScalarT, panzer::Cell, panzer::IP > fieldMult
One of the field multipliers ( , , etc.) 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
PHX::MDField< const Scalar, Cell, BASIS, IP > ConstArray_CellBasisIP
Array_CellBasisIPDim weighted_curl_basis_vector
Array_CellBasisIP weighted_curl_basis_scalar
ConstArray_CellBasisIPDim getCurlVectorBasis(const bool weighted, const bool cache=true, const bool force=false) const
Get the curl of a 3D vector basis evaluated at mesh points.
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
ConstArray_CellBasisIP getCurl2DVectorBasis(const bool weighted, const bool cache=true, const bool force=false) const
Get the curl of a 2D vector basis evaluated at mesh points.
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, panzer::Dim > vector3D_
A field representing the vector-valued function we're integrating ( ) in a three-dimensional problem.
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::MDField< const ScalarT, panzer::Cell, panzer::IP > vector2D_
A field representing the vector-valued function we're integrating ( ) in a two-dimensional problem.
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
void postRegistrationSetup(typename Traits::SetupData sd, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
panzer::BasisDescriptor bd_
The BasisDescriptor for the basis to use.
int spaceDim_
The spatial dimension of the vector-valued function we're integrating, either 2 or 3.
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.
typename EvalT::ScalarT ScalarT
The scalar type.
void evaluateFields(typename Traits::EvalData d)
Evaluate Fields.
Integrator_CurlBasisDotVector(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< Scalar, T0 > buildStaticArray(const std::string &str, int d0) const
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.
Teuchos::RCP< const std::vector< panzer::Workset > > worksets_