43 #ifndef PANZER_PRODUCT_IMPL_HPP 
   44 #define PANZER_PRODUCT_IMPL_HPP 
   53 template<
typename EvalT, 
typename Traits>
 
   59   std::string product_name = p.
get<std::string>(
"Product Name");
 
   65   if(p.
isType<
double>(
"Scaling"))
 
   68   product = PHX::MDField<ScalarT>(product_name, data_layout);
 
   70   this->addEvaluatedField(
product);
 
   72   values.resize(value_names->size());
 
   73   for (std::size_t i=0; i < value_names->size(); ++i) {
 
   74     values[i] = PHX::MDField<const ScalarT>( (*value_names)[i], data_layout);
 
   75     this->addDependentField(
values[i]);
 
   78   std::string n = 
"Product Evaluator";
 
   83 template<
int RANK, 
typename Scalar>
 
   86   const PHX::MDField<Scalar> 
base_;
 
   92   KOKKOS_INLINE_FUNCTION
 
   95     using idx_t = PHX::index_t;
 
  101       for (idx_t ind2=0; ind2 < static_cast<idx_t>(
base_.extent(1)); ind2++)
 
  105       for (idx_t ind2=0; ind2 < static_cast<idx_t>(
base_.extent(1)); ind2++)
 
  106         for (idx_t ind3=0; ind3 < static_cast<idx_t>(
base_.extent(2)); ind3++)
 
  110       for (idx_t ind2=0; ind2 < static_cast<idx_t>(
base_.extent(1)); ind2++)
 
  111         for (idx_t ind3=0; ind3 < static_cast<idx_t>(
base_.extent(2)); ind3++)
 
  112           for (idx_t ind4=0; ind4 < static_cast<idx_t>(
base_.extent(3)); ind4++)
 
  113             base_(ind1,ind2,ind3,ind4) = 
base_(ind1,ind2,ind3,ind4)*
source_(ind1,ind2,ind3,ind4);
 
  116       for (idx_t ind2=0; ind2 < static_cast<idx_t>(
base_.extent(1)); ind2++)
 
  117         for (idx_t ind3=0; ind3 < static_cast<idx_t>(
base_.extent(2)); ind3++)
 
  118           for (idx_t ind4=0; ind4 < static_cast<idx_t>(
base_.extent(3)); ind4++)
 
  119             for (idx_t ind5=0; ind5 < static_cast<idx_t>(
base_.extent(4)); ind5++)
 
  120               base_(ind1,ind2,ind3,ind4,ind5) = 
base_(ind1,ind2,ind3,ind4,ind5)*
source_(ind1,ind2,ind3,ind4,ind5);
 
  123       for (idx_t ind2=0; ind2 < static_cast<idx_t>(
base_.extent(1)); ind2++)
 
  124         for (idx_t ind3=0; ind3 < static_cast<idx_t>(
base_.extent(2)); ind3++)
 
  125           for (idx_t ind4=0; ind4 < static_cast<idx_t>(
base_.extent(3)); ind4++)
 
  126             for (idx_t ind5=0; ind5 < static_cast<idx_t>(
base_.extent(4)); ind5++)
 
  127               for (idx_t ind6=0; ind6 < static_cast<idx_t>(
base_.extent(5)); ind6++)
 
  128                 base_(ind1,ind2,ind3,ind4,ind5,ind6) = 
base_(ind1,ind2,ind3,ind4,ind5,ind6)*
source_(ind1,ind2,ind3,ind4,ind5,ind6);
 
  131       for (idx_t ind2=0; ind2 < static_cast<idx_t>(
base_.extent(1)); ind2++)
 
  132         for (idx_t ind3=0; ind3 < static_cast<idx_t>(
base_.extent(2)); ind3++)
 
  133           for (idx_t ind4=0; ind4 < static_cast<idx_t>(
base_.extent(3)); ind4++)
 
  134             for (idx_t ind5=0; ind5 < static_cast<idx_t>(
base_.extent(4)); ind5++)
 
  135               for (idx_t ind6=0; ind6 < static_cast<idx_t>(
base_.extent(5)); ind6++)
 
  136                 for (idx_t ind7=0; ind7 < static_cast<idx_t>(
base_.extent(6)); ind7++)
 
  137                   base_(ind1,ind2,ind3,ind4,ind5,ind6,ind7) = 
base_(ind1,ind2,ind3,ind4,ind5,ind6,ind7)*
source_(ind1,ind2,ind3,ind4,ind5,ind6,ind7);
 
  143 template<
typename EvalT, 
typename Traits>
 
  149   product.deep_copy(
ScalarT(scaling));
 
  151   for (std::size_t i = 0; i < values.size(); ++i) {
 
  152     const auto length = product.extent(0);
 
  153     if (product.rank() == 1){
 
  156     else if (product.rank() == 2){
 
  159     else if (product.rank() == 3){
 
  162     else if (product.rank() == 4){
 
  165     else if (product.rank() == 5){
 
  168     else if (product.rank() == 6){
 
  171     else if (product.rank() == 7){
 
KOKKOS_INLINE_FUNCTION void operator()(const PHX::index_t &ind1) const 
 
T & get(const std::string &name, T def_value)
 
const PHX::MDField< Scalar > base_
 
Product(const Teuchos::ParameterList &p)
 
V_MultiplyFunctor(PHX::MDField< Scalar > &base, const PHX::MDField< const Scalar > &source)
 
const PHX::MDField< const Scalar > source_
 
PHX::MDField< ScalarT > product
 
bool isType(const std::string &name) const 
 
std::vector< PHX::MDField< const ScalarT > > values
 
typename EvalT::ScalarT ScalarT
 
void evaluateFields(typename Traits::EvalData d)