Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_STK_ScatterCellAvgVector_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_STK_SCATTER_CELL_AVG_VECTOR_IMPL_HPP
12 #define PANZER_STK_SCATTER_CELL_AVG_VECTOR_IMPL_HPP
13 
14 #include "Teuchos_Assert.hpp"
15 
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"
21 
24 
25 #include "Teuchos_FancyOStream.hpp"
26 #include "Teuchos_ArrayRCP.hpp"
27 
28 namespace panzer_stk {
29 
30 template<typename EvalT, typename Traits>
33  const Teuchos::ParameterList& p) :
34  mesh_(p.get<Teuchos::RCP<STK_Interface> >("Mesh"))
35 {
36  using panzer::Cell;
37  using panzer::Point;
38  using panzer::Dim;
39 
40  std::string scatterName = p.get<std::string>("Scatter Name");
41 
42  if (p.isParameter("Variable Scale Factors Map"))
43  varScaleFactors_ = p.get<Teuchos::RCP<std::map<std::string,double>>>("Variable Scale Factors Map");
44 
45  const std::vector<std::string> & names =
46  *(p.get< Teuchos::RCP< std::vector<std::string> > >("Field Names"));
47 
50 
51  // build dependent fields
52  scatterFields_.resize(names.size());
53  stkFields_.resize(names.size());
54  for (std::size_t fd = 0; fd < names.size(); ++fd)
55  {
57  this->addDependentField(scatterFields_[fd]);
58  }
59 
60  // setup a dummy field to evaluate
61  PHX::Tag<ScalarT> scatterHolder(scatterName,Teuchos::rcp(new PHX::MDALayout<panzer::Dummy>(0)));
62  this->addEvaluatedField(scatterHolder);
63 
64  this->setName(scatterName+": PanzerSTK::ScatterCellAvgVectors");
65 }
66 
67 
68 template<typename EvalT, typename Traits>
69 void
72  typename Traits::SetupData /* d */,
73  PHX::FieldManager<Traits>& /* fm */)
74 {
75  for (std::size_t fd = 0; fd < scatterFields_.size(); ++fd)
76  {
77  std::string fieldName = scatterFields_[fd].fieldTag().name();
78 
79  stkFields_[fd] = mesh_->getMetaData()->get_field<double>(stk::topology::ELEMENT_RANK, fieldName);
80  }
81 }
82 
83 
84 template<typename EvalT, typename Traits>
85 void
88  typename Traits::EvalData workset)
89 {
90  panzer::MDFieldArrayFactory af("",true);
91 
92  // for convenience pull out some objects from workset
93  const std::vector<std::size_t> & localCellIds = this->wda(workset).cell_local_ids;
94  std::string blockId = this->wda(workset).block_id;
95  std::string d_mod[3] = {"X","Y","Z"};
96 
97  // loop over the number of vector fields requested for exodus output
98  for(std::size_t fieldIndex = 0; fieldIndex < scatterFields_.size(); fieldIndex++)
99  {
101  std::string fieldName = field.fieldTag().name();
102  const int numCells = workset.num_cells;
103  const int numPoints = field.extent(1);
104  const int numDims = field.extent(2);
105 
106  for (int dim = 0; dim < numDims; dim++)
107  {
108  // std::vector<double> average(numCells,0.0);
110 
111  // write to double field
112  Kokkos::parallel_for("ScatterCellAvgVector",numCells,KOKKOS_LAMBDA(const int i){
113  average(i,0) = 0.0;
114  for(int j = 0; j < numPoints; j++) { // loop over IPs
115  average(i,0) += Sacado::scalarValue(field(i,j,dim));
116  }
117  average(i,0) /= numPoints;
118  });
119  double scalef = 1.0;
120 
121  if (!varScaleFactors_.is_null())
122  {
123  std::map<std::string,double> *tmp_sfs = varScaleFactors_.get();
124  if(tmp_sfs->find(fieldName) != tmp_sfs->end())
125  scalef = (*tmp_sfs)[fieldName];
126  }
127 
128  PHX::Device().fence();
129  mesh_->setCellFieldData(fieldName+d_mod[dim],blockId,localCellIds,average.get_view(),scalef);
130  }
131  }
132 }
133 
134 } // end panzer_stk
135 
136 #endif
T & get(const std::string &name, T def_value)
PHX::MDField< Scalar, T0 > buildStaticArray(const std::string &str, int d0) const
int numPoints
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)
ScatterCellAvgVector(const Teuchos::ParameterList &p)
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_vector
Data layout for vector fields.
PHX::MDField< ScalarT, panzer::Cell, panzer::BASIS > field
A field to which we&#39;ll contribute, or in which we&#39;ll store, the result of computing this integral...
std::vector< PHX::MDField< const ScalarT, panzer::Cell, panzer::Point, panzer::Dim > > scatterFields_
Teuchos::RCP< std::map< std::string, double > > varScaleFactors_