11 #ifndef PANZER_STK_SCATTER_CELL_AVG_QUANTITY_IMPL_HPP 
   12 #define PANZER_STK_SCATTER_CELL_AVG_QUANTITY_IMPL_HPP 
   14 #include "Teuchos_Assert.hpp" 
   16 #include "Phalanx_config.hpp" 
   17 #include "Phalanx_Evaluator_Macros.hpp" 
   18 #include "Phalanx_MDField.hpp" 
   19 #include "Phalanx_DataLayout.hpp" 
   20 #include "Phalanx_DataLayout_MDALayout.hpp" 
   25 #include "Teuchos_FancyOStream.hpp" 
   26 #include "Teuchos_ArrayRCP.hpp" 
   28 namespace panzer_stk {
 
   30 template<
typename EvalT, 
typename Traits>
 
   39   std::string scatterName = p.
get<std::string>(
"Scatter Name");
 
   45   const std::vector<std::string> & names =
 
   54   for (std::size_t fd = 0; fd < names.size(); ++fd) {
 
   62   this->addEvaluatedField(scatterHolder);
 
   64   this->setName(scatterName+
": STK-Scatter Cell Fields");
 
   67 template<
typename EvalT, 
typename Traits>
 
   71   typename Traits::SetupData ,
 
   74   for (std::size_t fd = 0; fd < scatterFields_.size(); ++fd) {
 
   75     std::string fieldName = scatterFields_[fd].fieldTag().name();
 
   77     stkFields_[fd] = mesh_->getMetaData()->get_field<
double>(stk::topology::ELEMENT_RANK, fieldName);
 
   81 template<
typename EvalT, 
typename Traits>
 
   85   typename Traits::EvalData workset)
 
   90    const std::vector<std::size_t> & localCellIds = this->wda(workset).cell_local_ids;
 
   91    std::string blockId = this->wda(workset).block_id;
 
   93    for(std::size_t fieldIndex=0; fieldIndex<scatterFields_.size();fieldIndex++) {
 
   97      Kokkos::parallel_for(
"ScatterCellAvgQuantity",
field.extent(0), KOKKOS_LAMBDA(
int i) {
 
   98        for(
unsigned j=0; j<
field.extent(1);j++)
 
   99    average(i,0) += Sacado::scalarValue(
field(i,j));
 
  100        average(i,0) /= 
field.extent(1);
 
  104      std::string varname = scatterFields_[fieldIndex].fieldTag().name();
 
  107      if (!varScaleFactors_.is_null())
 
  109        std::map<std::string,double> *tmp_sfs = varScaleFactors_.get();
 
  110        if(tmp_sfs->find(varname) != tmp_sfs->end())
 
  111          scalef = (*tmp_sfs)[varname];
 
  114      mesh_->setCellFieldData(varname,blockId,localCellIds,average,scalef);
 
T & get(const std::string &name, T def_value)
PHX::MDField< Scalar, T0 > buildStaticArray(const std::string &str, int d0) const 
std::vector< PHX::MDField< const ScalarT, panzer::Cell, panzer::Point > > scatterFields_
bool isParameter(const std::string &name) const 
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< PHX::DataLayout > dl_scalar
Data layout for scalar fields. 
KOKKOS_FORCEINLINE_FUNCTION array_type get_static_view()
Teuchos::RCP< std::map< std::string, double > > varScaleFactors_
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...
void evaluateFields(typename Traits::EvalData d)
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
ScatterCellAvgQuantity(const Teuchos::ParameterList &p)
std::vector< VariableField * > stkFields_