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" 
   25 template<
typename EvalT, 
typename Traits> 
 
   34   const std::vector<std::string>& names = 
 
   41   gatherFields_.resize(names.size());
 
   42   stkFields_.resize(names.size());
 
   43   for (std::size_t fd = 0; fd < names.size(); ++fd) {
 
   46     this->addEvaluatedField(gatherFields_[fd]);
 
   49   this->setName(
"Gather STK Fields");
 
   53 template<
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);
 
   67                                  "panzer_stk::GatherFields: STK field " << 
"\"" << fieldName << 
"\" "  
   68                                  "not found.\n STK meta data follows: \n\n" << ss.str());
 
   74 template<
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);
 
T & get(const std::string &name, T def_value)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Teuchos::RCP< const PureBasis > getBasis() const 
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &vm)
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...
Teuchos::RCP< PHX::DataLayout > functional
<Cell,Basis> 
panzer_stk::STK_Interface::SolutionFieldType VariableField
void evaluateFields(typename Traits::EvalData d)