43 #ifndef PANZER_STK_SCATTER_CELL_QUANTITY_IMPL_HPP 44 #define PANZER_STK_SCATTER_CELL_QUANTITY_IMPL_HPP 46 #include "Teuchos_Assert.hpp" 48 #include "Phalanx_config.hpp" 49 #include "Phalanx_Evaluator_Macros.hpp" 50 #include "Phalanx_MDField.hpp" 51 #include "Phalanx_DataLayout.hpp" 52 #include "Phalanx_DataLayout_MDALayout.hpp" 57 #include "Teuchos_FancyOStream.hpp" 58 #include "Teuchos_ArrayRCP.hpp" 62 template<
typename EvalT,
typename Traits>
65 const Teuchos::ParameterList& p) :
71 std::string scatterName = p.get<std::string>(
"Scatter Name");
72 int worksetSize = p.get<
int>(
"Workset Size");
75 if (p.isParameter(
"Variable Scale Factors Map"))
76 varScaleFactors_ = p.get<Teuchos::RCP<std::map<std::string,double>>>(
"Variable Scale Factors Map");
78 const std::vector<std::string> & names =
79 *(p.get< Teuchos::RCP< std::vector<std::string> > >(
"Field Names"));
81 Teuchos::RCP<PHX::DataLayout> dl_cell = Teuchos::rcp(
new PHX::MDALayout<Cell>(worksetSize));
85 for (std::size_t fd = 0; fd < names.size(); ++fd) {
86 scatterFields_[fd] = PHX::MDField<const ScalarT,Cell>(names[fd],dl_cell);
91 PHX::Tag<ScalarT> scatterHolder(scatterName,Teuchos::rcp(
new PHX::MDALayout<panzer::Dummy>(0)));
92 this->addEvaluatedField(scatterHolder);
94 this->setName(scatterName+
": STK-Scatter Cell Quantity Fields");
97 template<
typename EvalT,
typename Traits>
101 typename Traits::SetupData ,
104 for (std::size_t fd = 0; fd < scatterFields_.size(); ++fd)
105 std::string fieldName = scatterFields_[fd].fieldTag().name();
108 template<
typename EvalT,
typename Traits>
112 typename Traits::EvalData workset)
117 const std::vector<std::size_t> & localCellIds = this->wda(workset).cell_local_ids;
118 std::string blockId = this->wda(workset).block_id;
120 for(std::size_t fieldIndex=0; fieldIndex<scatterFields_.size();fieldIndex++) {
121 PHX::MDField<const ScalarT,panzer::Cell> &
field = scatterFields_[fieldIndex];
126 for(
unsigned i=0; i<
field.extent(0);i++)
127 value(i,0) = Sacado::scalarValue(
field(i));
129 std::string varname =
field.fieldTag().name();
132 if (!varScaleFactors_.is_null())
134 std::map<std::string,double> *tmp_sfs = varScaleFactors_.get();
135 if(tmp_sfs->find(varname) != tmp_sfs->end())
136 scalef = (*tmp_sfs)[varname];
139 mesh_->setCellFieldData(
field.fieldTag().name(),blockId,localCellIds,value.get_view(),scalef);
PHX::MDField< Scalar, T0 > buildStaticArray(const std::string &str, int d0) const
ScatterCellQuantity(const Teuchos::ParameterList &p)
std::vector< PHX::MDField< const ScalarT, panzer::Cell > > scatterFields_
Teuchos::RCP< std::map< std::string, double > > varScaleFactors_
void evaluateFields(typename Traits::EvalData d)
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 postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &fm)