Panzer  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Panzer_STK_GatherExodusCellDataToIP_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_GATHER_EXODUS_CELL_DATA_TO_IP_IMPL_HPP
12 #define PANZER_STK_GATHER_EXODUS_CELL_DATA_TO_IP_IMPL_HPP
13 
14 #include "Teuchos_Assert.hpp"
15 #include "Phalanx_DataLayout.hpp"
16 #include "Panzer_BasisIRLayout.hpp"
17 #include "Teuchos_FancyOStream.hpp"
18 #include <sstream>
19 
20 // **********************************************************************
21 // Specialization: Residual
22 // **********************************************************************
23 
24 template<typename EvalT, typename Traits>
27  const std::vector<std::string>& fieldNames,
28  const std::vector<std::string>& exodusNames,
29  const Teuchos::RCP<panzer::IntegrationRule>& integrationRule)
30  : mesh_(mesh),
31  exodusNames_(exodusNames)
32 {
33  using panzer::Cell;
34  using panzer::IP;
35 
36  TEUCHOS_ASSERT(fieldNames.size() == exodusNames.size());
37  TEUCHOS_ASSERT(nonnull(mesh));
38  TEUCHOS_ASSERT(nonnull(integrationRule));
39 
40  gatherFields_.resize(fieldNames.size());
41  stkFields_.resize(fieldNames.size());
42  for (std::size_t fd = 0; fd < fieldNames.size(); ++fd) {
43  gatherFields_[fd] = PHX::MDField<ScalarT,Cell,IP>(fieldNames[fd],integrationRule->dl_scalar);
44  this->addEvaluatedField(gatherFields_[fd]);
45  }
46 
47  std::ostringstream os;
48  for (size_t i=0; i < fieldNames.size(); ++i) {
49  os << fieldNames[i];
50  if (i < fieldNames.size()-1)
51  os << ",";
52  }
53 
54  this->setName("GatherExodusCellDataToIP: "+os.str());
55 }
56 
57 // **********************************************************************
58 template<typename EvalT, typename Traits>
60 postRegistrationSetup(typename Traits::SetupData /* d */,
61  PHX::FieldManager<Traits>& /* fm */)
62 {
63  for (std::size_t fd = 0; fd < gatherFields_.size(); ++fd) {
64  std::string fieldName = gatherFields_[fd].fieldTag().name();
65 
66  stkFields_[fd] = mesh_->getMetaData()->get_field<double>(stk::topology::ELEMENT_RANK, exodusNames_[fd]);
67 
68  if(stkFields_[fd]==0) {
69  std::stringstream ss;
70  mesh_->printMetaData(ss);
71  TEUCHOS_TEST_FOR_EXCEPTION(true,std::invalid_argument,
72  "panzer_stk::GatherExodusCellDataToIP: STK field " << "\"" << fieldName << "\" "
73  "not found.\n STK meta data follows: \n\n" << ss.str());
74  }
75  }
76 }
77 
78 // **********************************************************************
79 template<typename EvalT, typename Traits>
81 evaluateFields(typename Traits::EvalData workset)
82 {
83  const std::vector<stk::mesh::Entity>& localElementEntities = *mesh_->getElementsOrderedByLID();
84  const std::vector<std::size_t>& localWorksetCellIds = this->wda(workset).cell_local_ids;
85 
86  auto host_mirror = Kokkos::create_mirror_view(Kokkos::WithoutInitializing,gatherFields_[0].get_static_view());
87 
88  for (std::size_t fieldIndex=0; fieldIndex < gatherFields_.size(); ++fieldIndex) {
89  VariableField* field = stkFields_[fieldIndex];
90  const std::size_t numQuadPoints = gatherFields_[fieldIndex].extent(1);
91 
92  for(std::size_t worksetCellIndex=0;worksetCellIndex<localWorksetCellIds.size();++worksetCellIndex) {
93  std::size_t cellLocalId = localWorksetCellIds[worksetCellIndex];
94  for(std::size_t qp=0; qp < numQuadPoints; ++qp) {
95  double value = *stk::mesh::field_data(*field, localElementEntities[cellLocalId]);
96  host_mirror(worksetCellIndex,qp) = value;
97  }
98  }
99 
100  Kokkos::deep_copy(gatherFields_[fieldIndex].get_static_view(),host_mirror);
101  }
102 }
103 #endif
std::vector< PHX::MDField< ScalarT, panzer::Cell, panzer::IP > > gatherFields_
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
panzer_stk::STK_Interface::SolutionFieldType VariableField
Teuchos::RCP< PHX::DataLayout > dl_scalar
Data layout for scalar fields.
void postRegistrationSetup(typename Traits::SetupData d, PHX::FieldManager< Traits > &vm)
GatherExodusCellDataToIP(const Teuchos::RCP< const panzer_stk::STK_Interface > &mesh, const std::vector< std::string > &fieldNames, const std::vector< std::string > &exodusNames, const Teuchos::RCP< panzer::IntegrationRule > &integrationRule)
bool nonnull(const boost::shared_ptr< T > &p)
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...
#define TEUCHOS_ASSERT(assertion_test)