43 #ifndef Panzer_Integrator_CurlBasisDotVector_impl_hpp
44 #define Panzer_Integrator_CurlBasisDotVector_impl_hpp
53 #include "Intrepid2_FunctionSpaceTools.hpp"
62 #include "Phalanx_KokkosDeviceTypes.hpp"
71 template<
typename EvalT,
typename Traits>
75 const std::string& resName,
76 const std::string& valName,
80 const std::vector<std::string>& fmNames
83 evalStyle_(evalStyle),
84 useDescriptors_(false),
85 multiplier_(multiplier),
86 basisName_(basis.name())
96 using std::invalid_argument;
97 using std::logic_error;
103 "Integrator_CurlBasisDotVector called with an empty residual name.")
105 "Integrator_CurlBasisDotVector called with an empty value name.")
108 "Error: Integrator_CurlBasisDotVector: Basis of type \""
109 << tmpBasis->name() <<
"\" is not a vector basis.")
111 logic_error,
"Error: Integrator_CurlBasisDotVector: Basis of type \""
112 << tmpBasis->name() <<
"\" does not require orientations, though it " \
113 "should for its use in this Evaluator.")
115 "Error: Integrator_CurlBasisDotVector: Basis of type \""
116 << tmpBasis->name() <<
"\" does not support curl.")
118 (tmpBasis->dimension() == 3)), logic_error,
119 "Error: Integrator_CurlBasisDotVector requires either a two- or " \
120 "three-dimensional basis. The basis \"" << tmpBasis->name()
142 this->addContributedField(
field_);
144 this->addEvaluatedField(
field_);
149 kokkosFieldMults_ = View<View<const ScalarT**,typename PHX::DevLayout<ScalarT>::type,PHX::Device>*>(
150 "CurlBasisDotVector::KokkosFieldMultipliers", fmNames.size());
151 for (
const auto& name : fmNames)
158 string n(
"Integrator_CurlBasisDotVector (");
163 n +=
"): " +
field_.fieldTag().name();
172 template<
typename EvalT,
typename Traits>
179 p.get<std::string>(
"Residual Name"),
180 p.get<std::string>(
"Value Name"),
183 p.get<double>(
"Multiplier"),
184 p.isType<Teuchos::
RCP<
const std::vector<std::string>>>
185 (
"Field Multipliers") ?
186 (*p.get<Teuchos::
RCP<
const std::vector<std::string>>>
187 (
"Field Multipliers")) : std::vector<std::string>())
202 template<
typename EvalT,
typename TRAITS>
206 const PHX::FieldTag& resTag,
207 const PHX::FieldTag& valTag,
212 const std::vector<PHX::FieldTag>& multipliers
215 evalStyle_(evalStyle),
216 useDescriptors_(true),
219 multiplier_(multiplier),
224 using std::logic_error;
229 "Error: Integrator_CurlBasisDotVector: Basis of type \""
230 <<
bd_.
getType() <<
"\" does not support curl.")
232 logic_error,
"Error: Integrator_CurlBasisDotVector works on either " \
233 "two- or three-dimensional problems. You've provided spaceDim = "
252 this->addContributedField(
field_);
254 this->addEvaluatedField(
field_);
259 kokkosFieldMults_ = View<View<const ScalarT**,typename PHX::DevLayout<ScalarT>::type,PHX::Device>*>(
260 "CurlBasisDotVector::KokkosFieldMultipliers", multipliers.size());
261 for (
const auto& fm : multipliers)
268 string n(
"Integrator_CurlBasisDotVector (");
273 n +=
"): " +
field_.fieldTag().name();
282 template<
typename EvalT,
typename Traits>
294 for (
size_t i(0); i < fieldMults_.size(); ++i)
295 kokkosFieldMults_(i) = fieldMults_[i].get_static_view();
299 if (not useDescriptors_)
305 fm.template getKokkosExtendedDataTypeDimensions<EvalT>(),
true);
308 "Integrator_CurlBasisDotVector: 2-D Result", vector2D_.extent(0),
309 vector2D_.extent(1));
312 "Integrator_CurlBasisDotVector: 3-D Result", vector3D_.extent(0),
313 vector3D_.extent(1), 3);
334 template<
typename ScalarT>
343 struct ScalarMultiplierTag
351 struct FieldMultiplierTag
366 KOKKOS_INLINE_FUNCTION
367 void operator()(
const ScalarMultiplierTag,
const unsigned cell)
const
369 int numQP(
result.extent(1));
370 for (
int qp(0); qp < numQP; ++qp)
386 KOKKOS_INLINE_FUNCTION
387 void operator()(
const FieldMultiplierTag,
const unsigned cell)
const
389 int numQP(
result.extent(1));
390 for (
int qp(0); qp < numQP; ++qp)
398 using execution_space = PHX::Device;
404 PHX::MDField<ScalarT, panzer::Cell, panzer::IP>
result;
416 PHX::MDField<const ScalarT, panzer::Cell, panzer::IP>
vectorField;
422 PHX::MDField<const ScalarT, panzer::Cell, panzer::IP>
fieldMult;
436 template<
typename ScalarT>
445 struct ScalarMultiplierTag
453 struct FieldMultiplierTag
468 KOKKOS_INLINE_FUNCTION
469 void operator()(
const ScalarMultiplierTag,
const unsigned cell)
const
472 for (
int qp(0); qp < numQP; ++qp)
473 for (
int dim(0); dim < numDim; ++dim)
489 KOKKOS_INLINE_FUNCTION
490 void operator()(
const FieldMultiplierTag,
const unsigned cell)
const
493 for (
int qp(0); qp < numQP; ++qp)
494 for (
int dim(0); dim < numDim; ++dim)
502 using execution_space = PHX::Device;
508 PHX::MDField<ScalarT, panzer::Cell, panzer::IP, panzer::Dim>
result;
520 PHX::MDField<const ScalarT, panzer::Cell, panzer::IP, panzer::Dim>
527 PHX::MDField<const ScalarT, panzer::Cell, panzer::IP>
fieldMult;
541 template<
typename ScalarT,
int spaceDim>
558 KOKKOS_INLINE_FUNCTION
559 void operator()(
const unsigned cell)
const
564 for (
int basis(0); basis < numBases; ++basis)
567 field(cell, basis) = 0.0;
568 for (
int qp(0); qp < numQP; ++qp)
569 for (
int dim(0); dim < spaceDim; ++dim)
579 using execution_space = PHX::Device;
584 PHX::MDField<ScalarT, panzer::Cell, panzer::IP, panzer::Dim>
result;
590 PHX::MDField<ScalarT, panzer::Cell, panzer::BASIS>
field;
616 template<
typename ScalarT>
632 KOKKOS_INLINE_FUNCTION
633 void operator()(
const unsigned cell)
const
638 for (
int basis(0); basis < numBases; ++basis)
641 field(cell, basis) = 0.0;
642 for (
int qp(0); qp < numQP; ++qp)
652 using execution_space = PHX::Device;
657 PHX::MDField<ScalarT, panzer::Cell, panzer::IP>
result;
663 PHX::MDField<ScalarT, panzer::Cell, panzer::BASIS>
field;
668 PHX::MDField<double, panzer::Cell, panzer::BASIS, panzer::IP>
685 template<
typename EvalT,
typename Traits>
691 using Kokkos::parallel_for;
692 using Kokkos::RangePolicy;
700 this->wda(workset).getBasisValues(bd_, id_) :
701 *this->wda(workset).bases[basisIndex_];
708 using PreMultiply = PreMultiply2D<ScalarT>;
709 using ScalarMultiplierTag =
typename PreMultiply::ScalarMultiplierTag;
710 using FieldMultiplierTag =
typename PreMultiply::FieldMultiplierTag;
711 PreMultiply preMultiply;
712 preMultiply.result = result2D_;
713 preMultiply.multiplier = multiplier_;
714 preMultiply.vectorField = vector2D_;
718 parallel_for(RangePolicy<Device, ScalarMultiplierTag>(0,
723 for (
const auto&
field : fieldMults_)
725 preMultiply.fieldMult =
field;
726 parallel_for(RangePolicy<Device, FieldMultiplierTag>(0,
731 Integrate2D<ScalarT> integrate;
732 integrate.result = result2D_;
733 integrate.field = field_;
735 integrate.evalStyle = evalStyle_;
736 parallel_for(workset.
num_cells, integrate);
742 using PreMultiply = PreMultiply3D<ScalarT>;
743 using ScalarMultiplierTag =
typename PreMultiply::ScalarMultiplierTag;
744 using FieldMultiplierTag =
typename PreMultiply::FieldMultiplierTag;
745 PreMultiply preMultiply;
746 preMultiply.result = result3D_;
747 preMultiply.multiplier = multiplier_;
748 preMultiply.vectorField = vector3D_;
752 parallel_for(RangePolicy<Device, ScalarMultiplierTag>(0,
757 for (
const auto&
field : fieldMults_)
759 preMultiply.fieldMult =
field;
760 parallel_for(RangePolicy<Device, FieldMultiplierTag>(0,
765 Integrate3D<ScalarT, 3> integrate;
766 integrate.result = result3D_;
767 integrate.field = field_;
769 integrate.evalStyle = evalStyle_;
770 parallel_for(workset.
num_cells, integrate);
779 template<
typename EvalT,
typename TRAITS>
793 p->set<
string>(
"Residual Name",
"?");
794 p->set<
string>(
"Value Name",
"?");
796 p->set(
"Basis", basis);
799 p->set<
double>(
"Multiplier", 1.0);
801 p->set(
"Field Multipliers", fms);
807 #endif // Panzer_Integrator_CurlBasisDotVector_impl_hpp
Array_CellBasisIPDim weighted_curl_basis_vector
Kokkos::View< Kokkos::View< const ScalarT **, typename PHX::DevLayout< ScalarT >::type, PHX::Device > * > kokkosFieldMults_
The Kokkos::View representation of the (possibly empty) list of fields that are multipliers out in fr...
void evaluateFields(typename Traits::EvalData d)
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 ( ...
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
PHX::MDField< const ScalarT, panzer::Cell, panzer::IP > vectorField
A field representing the vector-valued function we're integrating ( ).
Teuchos::RCP< Teuchos::ParameterList > getValidParameters() const
Get Valid Parameters.
EvaluatorStyle
An indication of how an Evaluator will behave.
PHX::MDField< Scalar, T0 > buildStaticArray(const std::string &str, int d0) const
Teuchos::RCP< const PureBasis > getBasis() const
panzer::BasisDescriptor bd_
The BasisDescriptor for the basis to use.
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.
double multiplier
The scalar multiplier out in front of the integral ( ).
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
const std::string & getType() const
Get type of basis.
Teuchos::RCP< PHX::DataLayout > dl_scalar
Data layout for scalar fields.
panzer::EvaluatorStyle evalStyle
The EvaluatorStyle of the parent Integrator_CurlBasisDotVector object.
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...
int spaceDim_
The spatial dimension of the vector-valued function we're integrating, either 2 or 3...
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.
PHX::MDField< double, panzer::Cell, panzer::BASIS, panzer::IP, panzer::Dim > weightedCurlBasis
The vector basis information necessary for integration.
Array_CellBasisIP weighted_curl_basis_scalar
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
void postRegistrationSetup(typename Traits::SetupData sd, PHX::FieldManager< Traits > &fm)
Post-Registration Setup.
Teuchos::RCP< PHX::DataLayout > dl_vector
Data layout for vector fields.
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< 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_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.
const std::vector< std::pair< int, LocalOrdinal > > &pid_and_lid const
Description and data layouts associated with a particular basis.
PHX::MDField< const ScalarT, panzer::Cell, panzer::IP > fieldMult
One of the field multipliers ( , , etc.) out in front of the integral.
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis>
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< const std::vector< panzer::Workset > > worksets_
typename EvalT::ScalarT ScalarT
The scalar type.