Panzer Version of the Day
Loading...
Searching...
No Matches
Panzer_DotProduct_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_EVALUATOR_DotProduct_IMPL_HPP
12#define PANZER_EVALUATOR_DotProduct_IMPL_HPP
13
14#include <string>
15
16#include "Panzer_PointRule.hpp"
18
19#include "Teuchos_RCP.hpp"
20
21namespace panzer {
22
23template <typename EvalT,typename TraitsT>
24Teuchos::RCP<DotProduct<EvalT,TraitsT> >
25buildEvaluator_DotProduct(const std::string & resultName,
26 const panzer::PointRule & pr,
27 const std::string & vecA,
28 const std::string & vecB,
29 double multiplier,
30 const std::string & fieldMultiplier)
31{
32 Teuchos::ParameterList pl;
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);
39
40 return Teuchos::rcp(new DotProduct<EvalT,TraitsT>(pl));
41}
42
43//**********************************************************************
44template<typename EvalT, typename Traits>
47 const Teuchos::ParameterList& p)
48 : multiplier_field_on(false)
49{
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");
53
54 std::string multiplier_name = "";
55 if(p.isType<std::string>("Field Multiplier"))
56 multiplier_name = p.get<std::string>("Field Multiplier");
57
58 multiplier_value = 1.0;
59 if(p.isType<double>("Multiplier"))
60 multiplier_value = p.get<double>("Multiplier");
61
62 const Teuchos::RCP<const panzer::PointRule> pr =
63 p.get< Teuchos::RCP<const panzer::PointRule> >("Point Rule");
64
65 vec_a_dot_vec_b = PHX::MDField<ScalarT>(result_name, pr->dl_scalar);
66 vec_a = PHX::MDField<const ScalarT>(vec_a_name, pr->dl_vector);
67 vec_b = PHX::MDField<const ScalarT>(vec_b_name, pr->dl_vector);
68
69 if(multiplier_name!="") {
70 multiplier_field = PHX::MDField<const ScalarT>(multiplier_name,pr->dl_scalar);
72 this->addDependentField(multiplier_field);
73 }
74
75 this->addEvaluatedField(vec_a_dot_vec_b);
76 this->addDependentField(vec_a);
77 this->addDependentField(vec_b);
78
79 std::string n = "DotProduct: " + result_name + " = " + vec_a_name + " . " + vec_b_name;
80 this->setName(n);
81}
82
83//**********************************************************************
84template<typename EvalT, typename Traits>
85void
88 typename Traits::SetupData /* sd */,
90{
91 num_pts = vec_a.extent(1);
92 num_dim = vec_a.extent(2);
93
94 TEUCHOS_ASSERT(vec_a.extent(1) == vec_b.extent(1));
95 TEUCHOS_ASSERT(vec_a.extent(2) == vec_b.extent(2));
96}
97
98//**********************************************************************
99template<typename EvalT, typename Traits>
100void
103 typename Traits::EvalData workset)
104{
105
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();
110
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;
114
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);
120
121 if(l_multiplier_field_on)
122 vec_a_dot_vec_b_v(cell,p) *= l_multiplier_value*multiplier_field_v(cell,p);
123 else
124 vec_a_dot_vec_b_v(cell,p) *= l_multiplier_value;
125 }
126 });
127 Kokkos::fence();
128}
129
130//**********************************************************************
131
132}
133
134#endif
double multiplier
The scalar multiplier out in front of the integral ( ).
Evaluates dot product at a set of points.
PHX::MDField< ScalarT > vec_a_dot_vec_b
PHX::MDField< const ScalarT > vec_a
PHX::MDField< const ScalarT > multiplier_field
PHX::MDField< const ScalarT > vec_b
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
void evaluateFields(typename Traits::EvalData d)
typename EvalT::ScalarT ScalarT
DotProduct(const Teuchos::ParameterList &p)
int num_cells
DEPRECATED - use: numCells()
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.