43 #ifndef PANZER_EVALUATOR_DotProduct_IMPL_HPP 
   44 #define PANZER_EVALUATOR_DotProduct_IMPL_HPP 
   55 template <
typename EvalT,
typename TraitsT>
 
   59                           const std::string & vecA,
 
   60                           const std::string & vecB,
 
   62                           const std::string & fieldMultiplier)
 
   65   pl.
set(
"Result Name",resultName);
 
   66   pl.
set(
"Point Rule",Teuchos::rcpFromRef(pr));
 
   67   pl.
set(
"Vector A Name",vecA);
 
   68   pl.
set(
"Vector B Name",vecB);
 
   69   pl.
set(
"Multiplier",multiplier);
 
   70   pl.
set(
"Field Multiplier",fieldMultiplier);
 
   76 template<
typename EvalT, 
typename Traits>
 
   80   : multiplier_field_on(false)
 
   82   std::string result_name = p.
get<std::string>(
"Result Name");
 
   83   std::string vec_a_name = p.
get<std::string>(
"Vector A Name");
 
   84   std::string vec_b_name = p.
get<std::string>(
"Vector B Name");
 
   86   std::string multiplier_name = 
"";
 
   87   if(p.
isType<std::string>(
"Field Multiplier"))
 
   88     multiplier_name = p.
get<std::string>(
"Field Multiplier");
 
   91   if(p.
isType<
double>(
"Multiplier"))
 
  101   if(multiplier_name!=
"") {
 
  107   this->addEvaluatedField(vec_a_dot_vec_b);
 
  108   this->addDependentField(vec_a);
 
  109   this->addDependentField(vec_b);
 
  111   std::string n = 
"DotProduct: " + result_name + 
" = " + vec_a_name + 
" . " + vec_b_name;
 
  116 template<
typename EvalT, 
typename Traits>
 
  123   num_pts = vec_a.extent(1);
 
  124   num_dim = vec_a.extent(2);
 
  131 template<
typename EvalT, 
typename Traits>
 
  137   for (index_t cell = 0; cell < workset.
num_cells; ++cell) {
 
  138     for (
int p = 0; p < num_pts; ++p) {
 
  139       vec_a_dot_vec_b(cell,p) = 
ScalarT(0.0);
 
  140       for (
int dim = 0; dim < num_dim; ++dim)
 
  141   vec_a_dot_vec_b(cell,p) += vec_a(cell,p,dim) * vec_b(cell,p,dim); 
 
  143       if(multiplier_field_on)
 
  144         vec_a_dot_vec_b(cell,p) *= multiplier_value*multiplier_field(cell,p);
 
  146         vec_a_dot_vec_b(cell,p) *= multiplier_value;
 
PHX::MDField< const ScalarT > vec_a
 
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
 
T & get(const std::string &name, T def_value)
 
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
 
typename EvalT::ScalarT ScalarT
 
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