11 #ifndef PANZER_EVALUATOR_DotProduct_IMPL_HPP
12 #define PANZER_EVALUATOR_DotProduct_IMPL_HPP
23 template <
typename EvalT,
typename TraitsT>
27 const std::string & vecA,
28 const std::string & vecB,
30 const std::string & fieldMultiplier)
33 pl.
set(
"Result Name",resultName);
34 pl.
set(
"Point Rule",Teuchos::rcpFromRef(pr));
35 pl.
set(
"Vector A Name",vecA);
36 pl.
set(
"Vector B Name",vecB);
37 pl.
set(
"Multiplier",multiplier);
38 pl.
set(
"Field Multiplier",fieldMultiplier);
44 template<
typename EvalT,
typename Traits>
48 : multiplier_field_on(false)
50 std::string result_name = p.
get<std::string>(
"Result Name");
51 std::string vec_a_name = p.
get<std::string>(
"Vector A Name");
52 std::string vec_b_name = p.
get<std::string>(
"Vector B Name");
54 std::string multiplier_name =
"";
55 if(p.
isType<std::string>(
"Field Multiplier"))
56 multiplier_name = p.
get<std::string>(
"Field Multiplier");
59 if(p.
isType<
double>(
"Multiplier"))
69 if(multiplier_name!=
"") {
75 this->addEvaluatedField(vec_a_dot_vec_b);
76 this->addDependentField(vec_a);
77 this->addDependentField(vec_b);
79 std::string n =
"DotProduct: " + result_name +
" = " + vec_a_name +
" . " + vec_b_name;
84 template<
typename EvalT,
typename Traits>
91 num_pts = vec_a.extent(1);
92 num_dim = vec_a.extent(2);
99 template<
typename EvalT,
typename Traits>
106 auto vec_a_v = vec_a.get_static_view();
107 auto vec_b_v = vec_b.get_static_view();
108 auto vec_a_dot_vec_b_v = vec_a_dot_vec_b.get_static_view();
109 auto multiplier_field_v = multiplier_field.get_static_view();
111 int l_num_pts = num_pts, l_num_dim = num_dim;
112 auto l_multiplier_field_on = multiplier_field_on;
113 auto l_multiplier_value = multiplier_value;
115 Kokkos::parallel_for (workset.
num_cells, KOKKOS_LAMBDA (
const int cell) {
116 for (
int p = 0; p < l_num_pts; ++p) {
117 vec_a_dot_vec_b_v(cell,p) =
ScalarT(0.0);
118 for (
int dim = 0; dim < l_num_dim; ++dim)
119 vec_a_dot_vec_b_v(cell,p) += vec_a_v(cell,p,dim) * vec_b_v(cell,p,dim);
121 if(l_multiplier_field_on)
122 vec_a_dot_vec_b_v(cell,p) *= l_multiplier_value*multiplier_field_v(cell,p);
124 vec_a_dot_vec_b_v(cell,p) *= l_multiplier_value;
PHX::MDField< const ScalarT > vec_a
int num_cells
DEPRECATED - use: numCells()
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
T & get(const std::string &name, T def_value)
typename EvalT::ScalarT ScalarT
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
Evaluates dot product at a set of points.
PHX::MDField< const ScalarT > vec_b
double multiplier
The scalar multiplier out in front of the integral ( ).
Teuchos::RCP< DotProduct< EvalT, TraitsT > > buildEvaluator_DotProduct(const std::string &resultName, const panzer::PointRule &pr, const std::string &vecA, const std::string &vecB, double multiplier=1, const std::string &fieldMultiplier="")
Build a dot product evaluator. Evaluates dot product at a set of points.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< PHX::DataLayout > dl_scalar
Data layout for scalar fields.
DotProduct(const Teuchos::ParameterList &p)
void evaluateFields(typename Traits::EvalData d)
Teuchos::RCP< PHX::DataLayout > dl_vector
Data layout for vector fields.
bool isType(const std::string &name) const
#define TEUCHOS_ASSERT(assertion_test)
PHX::MDField< const ScalarT > multiplier_field
PHX::MDField< ScalarT > vec_a_dot_vec_b