11#ifndef PANZER_STK_GATHER_FIELDS_IMPL_HPP
12#define PANZER_STK_GATHER_FIELDS_IMPL_HPP
14#include "Teuchos_Assert.hpp"
15#include "Phalanx_DataLayout.hpp"
19#include "Teuchos_FancyOStream.hpp"
25template<
typename EvalT,
typename Traits>
27 GatherFields(
const Teuchos::RCP<const STK_Interface> & mesh,
const Teuchos::ParameterList& p)
34 const std::vector<std::string>& names =
35 *(p.get< Teuchos::RCP< std::vector<std::string> > >(
"Field Names"));
37 Teuchos::RCP<panzer::BasisIRLayout> basis =
38 p.get< Teuchos::RCP<panzer::BasisIRLayout> >(
"Basis");
41 gatherFields_.resize(names.size());
42 stkFields_.resize(names.size());
43 for (std::size_t fd = 0; fd < names.size(); ++fd) {
45 PHX::MDField<ScalarT,Cell,NODE>(names[fd],basis->functional);
46 this->addEvaluatedField(gatherFields_[fd]);
49 this->setName(
"Gather STK Fields");
53template<
typename EvalT,
typename Traits>
58 for (std::size_t fd = 0; fd < gatherFields_.size(); ++fd) {
59 std::string fieldName = gatherFields_[fd].fieldTag().name();
61 stkFields_[fd] = mesh_->getMetaData()->get_field<
double>(stk::topology::NODE_RANK, fieldName);
63 if(stkFields_[fd]==0) {
65 mesh_->printMetaData(ss);
66 TEUCHOS_TEST_FOR_EXCEPTION(
true,std::invalid_argument,
67 "panzer_stk::GatherFields: STK field " <<
"\"" << fieldName <<
"\" "
68 "not found.\n STK meta data follows: \n\n" << ss.str());
74template<
typename EvalT,
typename Traits>
78 const std::vector<stk::mesh::Entity> & localElements = *mesh_->getElementsOrderedByLID();
81 const std::vector<std::size_t> & localCellIds = this->wda(workset).cell_local_ids;
84 for(std::size_t worksetCellIndex=0;worksetCellIndex<localCellIds.size();++worksetCellIndex) {
85 std::size_t cellLocalId = localCellIds[worksetCellIndex];
86 stk::mesh::Entity
const* relations = mesh_->getBulkData()->begin_nodes(localElements[cellLocalId]);
89 for (std::size_t fieldIndex=0; fieldIndex<gatherFields_.size();fieldIndex++) {
92 std::size_t basisCnt = gatherFields_[fieldIndex].extent(1);
96 (gatherFields_[fieldIndex])(worksetCellIndex,0) = *stk::mesh::field_data(*
field, localElements[cellLocalId]);
100 for(std::size_t basis=0;basis<basisCnt;basis++) {
101 stk::mesh::Entity node = relations[basis];
102 (gatherFields_[fieldIndex])(worksetCellIndex,basis) = *stk::mesh::field_data(*
field, node);
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)
panzer_stk::STK_Interface::SolutionFieldType VariableField
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &vm)