11 #ifndef PANZER_STK_SCATTER_FIELDS_IMPL_HPP
12 #define PANZER_STK_SCATTER_FIELDS_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"
27 namespace panzer_stk {
29 template <
typename EvalT,
typename TraitsT>
34 const std::vector<std::string> & names)
36 std::vector<double> scaling;
38 initialize(scatterName,mesh,basis,names,scaling);
41 template <
typename EvalT,
typename TraitsT>
46 const std::vector<std::string> & names,
47 const std::vector<double> & scaling)
49 initialize(scatterName,mesh,basis,names,scaling);
52 template <
typename EvalT,
typename TraitsT>
57 const std::vector<std::string> & names,
58 const std::vector<double> & scaling)
66 bool correctScaling = (names.size()==scaling.size()) || (scaling.size()==0);
68 "panzer_stk::ScatterFields evaluator requites a consistent number of scaling parameters (equal to the number of field names) "
69 "or an empty \"Field Scaling\" vector");
72 scatterFields_.resize(names.size());
73 for (std::size_t fd = 0; fd < names.size(); ++fd) {
76 this->addDependentField(scatterFields_[fd]);
84 this->addEvaluatedField(scatterHolder);
86 this->setName(scatterName+
": STK-Scatter Fields");
89 template <
typename EvalT,
typename TraitsT>
94 for (std::size_t fd = 0; fd < scatterFields_.size(); ++fd)
95 std::string fieldName = scatterFields_[fd].fieldTag().name();
98 template <
typename EvalT,
typename TraitsT>
110 const std::vector<std::size_t> & localCellIds = this->wda(workset).cell_local_ids;
111 std::string blockId = this->wda(workset).block_id;
113 for(std::size_t fieldIndex=0; fieldIndex<scatterFields_.size();fieldIndex++) {
115 double scaling = (scaling_.size()>0) ? scaling_[fieldIndex] : 1.0;
119 mesh_->setSolutionFieldData(scatterFields_[fieldIndex].fieldTag().name(),blockId,
120 localCellIds,scatterFields_[fieldIndex].get_view(),scaling);
122 mesh_->setCellFieldData(scatterFields_[fieldIndex].fieldTag().name(),blockId,
123 localCellIds,scatterFields_[fieldIndex].get_view(),scaling);
void initialize(const std::string &scatterName, const Teuchos::RCP< STK_Interface > mesh, const Teuchos::RCP< const panzer::PureBasis > &basis, const std::vector< std::string > &names, const std::vector< double > &scaling)
EElementSpace getElementSpace() const
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
void postRegistrationSetup(typename TraitsT::SetupData d, PHX::FieldManager< TraitsT > &fm)
#define TEUCHOS_ASSERT(assertion_test)
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis> or <Cell,Basis>
ScatterFields(const std::string &scatterName, const Teuchos::RCP< STK_Interface > mesh, const Teuchos::RCP< const panzer::PureBasis > &basis, const std::vector< std::string > &names)
void evaluateFields(typename TraitsT::EvalData d)