11 #ifndef __Panzer_STK_Interface_hpp__
12 #define __Panzer_STK_Interface_hpp__
17 #include <stk_mesh/base/Types.hpp>
18 #include <stk_mesh/base/MetaData.hpp>
19 #include <stk_mesh/base/BulkData.hpp>
20 #include <stk_mesh/base/Field.hpp>
21 #include <stk_mesh/base/FieldBase.hpp>
23 #include "Kokkos_Core.hpp"
25 #include <Shards_CellTopology.hpp>
26 #include <Shards_CellTopologyData.h>
28 #include <PanzerAdaptersSTK_config.hpp>
29 #include <Kokkos_ViewFactory.hpp>
31 #include <unordered_map>
33 #ifdef PANZER_HAVE_IOSS
34 #include <stk_io/StkMeshIoBroker.hpp>
37 #ifdef PANZER_HAVE_PERCEPT
40 class URP_Heterogeneous_3D;
44 namespace panzer_stk {
46 class PeriodicBC_MatcherBase;
53 ElementDescriptor(stk::mesh::EntityId gid,
const std::vector<stk::mesh::EntityId> & nodes);
60 std::vector<stk::mesh::EntityId>
nodes_;
103 void addElementBlock(
const std::string & name,
const CellTopologyData * ctData);
109 const std::string & edgeBlockName,
110 const stk::topology & topology);
112 const std::string & edgeBlockName,
113 const CellTopologyData * ctData);
118 const std::string & faceBlockName,
119 const stk::topology & topology);
121 const std::string & faceBlockName,
122 const CellTopologyData * ctData);
126 void addSideset(
const std::string & name,
const CellTopologyData * ctData);
134 void addSolutionField(
const std::string & fieldName,
const std::string & blockId);
138 void addCellField(
const std::string & fieldName,
const std::string & blockId);
142 void addEdgeField(
const std::string & fieldName,
const std::string & blockId);
146 void addFaceField(
const std::string & fieldName,
const std::string & blockId);
157 const std::vector<std::string> & coordField,
158 const std::string & dispPrefix);
179 void initialize(stk::ParallelMachine parallelMach,
bool setupIO=
true,
180 const bool buildRefinementSupport =
false);
205 void addNode(stk::mesh::EntityId gid,
const std::vector<double> & coord);
262 std::vector<stk::mesh::EntityId> & subcellIds)
const;
266 void getMyElements(std::vector<stk::mesh::Entity> & elements)
const;
270 void getMyElements(
const std::string & blockID,std::vector<stk::mesh::Entity> & elements)
const;
279 void getNeighborElements(
const std::string & blockID,std::vector<stk::mesh::Entity> & elements)
const;
283 void getMyEdges(std::vector<stk::mesh::Entity> & edges)
const;
292 void getMyEdges(
const std::string & edgeBlockName,std::vector<stk::mesh::Entity> & edges)
const;
302 void getMyEdges(
const std::string & edgeBlockName,
const std::string & blockName,std::vector<stk::mesh::Entity> & edges)
const;
311 void getAllEdges(
const std::string & edgeBlockName,std::vector<stk::mesh::Entity> & edges)
const;
321 void getAllEdges(
const std::string & edgeBlockName,
const std::string & blockName,std::vector<stk::mesh::Entity> & edges)
const;
325 void getMyFaces(std::vector<stk::mesh::Entity> & faces)
const;
334 void getMyFaces(
const std::string & faceBlockName,std::vector<stk::mesh::Entity> & faces)
const;
344 void getMyFaces(
const std::string & faceBlockName,
const std::string & blockName,std::vector<stk::mesh::Entity> & faces)
const;
353 void getAllFaces(
const std::string & faceBlockName,std::vector<stk::mesh::Entity> & faces)
const;
363 void getAllFaces(
const std::string & faceBlockName,
const std::string & blockName,std::vector<stk::mesh::Entity> & faces)
const;
372 void getMySides(
const std::string & sideName,std::vector<stk::mesh::Entity> & sides)
const;
382 void getMySides(
const std::string & sideName,
const std::string & blockName,std::vector<stk::mesh::Entity> & sides)
const;
391 void getAllSides(
const std::string & sideName,std::vector<stk::mesh::Entity> & sides)
const;
401 void getAllSides(
const std::string & sideName,
const std::string & blockName,std::vector<stk::mesh::Entity> & sides)
const;
412 void getMyNodes(
const std::string & sideName,
const std::string & blockName,std::vector<stk::mesh::Entity> & nodes)
const;
422 stk::mesh::Entity
findConnectivityById(stk::mesh::Entity src, stk::mesh::EntityRank tgt_rank,
unsigned rel_id)
const;
437 const bool append =
false);
458 const bool append =
false,
459 const bool append_after_restart_time =
false,
460 const double restart_time = 0.0);
464 const std::vector<Ioss::Property>& ioss_properties,
465 const bool append =
false,
466 const bool append_after_restart_time =
false,
467 const double restart_time = 0.0);
506 const std::string& key,
530 const std::string& key,
531 const double& value);
554 const std::string& key,
555 const std::vector<int>& value);
578 const std::string& key,
579 const std::vector<double>& value);
590 #ifdef PANZER_HAVE_PERCEPT
599 {
if(
bulkData_==Teuchos::null)
return false;
600 return bulkData_->in_modifiable_state(); }
662 std::map<std::string, stk::mesh::Part*>::const_iterator itr =
elementBlocks_.find(name);
670 std::map<std::string, stk::mesh::Part*>::const_iterator itr =
edgeBlocks_.find(name);
678 std::map<std::string, stk::mesh::Part*>::const_iterator itr =
faceBlocks_.find(name);
690 return (itr !=
sidesets_.end()) ? itr->second :
nullptr;
700 return (itr !=
nodesets_.end()) ? itr->second :
nullptr;
716 void getNodeIdsForElement(stk::mesh::Entity element,std::vector<stk::mesh::EntityId> & nodeIds)
const;
721 std::vector<int> & relIds)
const;
726 std::vector<int> & relIds,
unsigned int matchType)
const;
730 void getElementsSharingNodes(
const std::vector<stk::mesh::EntityId> nodeId,std::vector<stk::mesh::Entity> & elements)
const;
763 std::size_t
edgeLocalId(stk::mesh::Entity elmt)
const;
767 std::size_t
edgeLocalId(stk::mesh::EntityId gid)
const;
789 std::size_t
faceLocalId(stk::mesh::Entity elmt)
const;
793 std::size_t
faceLocalId(stk::mesh::EntityId gid)
const;
808 {
return bulkData_->parallel_owner_rank(entity); }
812 inline bool isValid(stk::mesh::Entity entity)
const
824 const std::string & blockId)
const;
830 stk::mesh::Field<double> *
getCellField(
const std::string & fieldName,
831 const std::string & blockId)
const;
837 stk::mesh::Field<double> *
getEdgeField(
const std::string & fieldName,
838 const std::string & blockId)
const;
844 stk::mesh::Field<double> *
getFaceField(
const std::string & fieldName,
845 const std::string & blockId)
const;
871 template <
typename ArrayT>
873 const std::vector<std::size_t> & localElementIds,
const ArrayT & solutionValues,
double scaleValue=1.0);
889 template <
typename ArrayT>
891 const std::vector<std::size_t> & localElementIds,ArrayT & solutionValues)
const;
907 template <
typename ArrayT>
908 void setCellFieldData(
const std::string & fieldName,
const std::string & blockId,
909 const std::vector<std::size_t> & localElementIds,
const ArrayT & solutionValues,
double scaleValue=1.0);
935 template <
typename ArrayT>
936 void setEdgeFieldData(
const std::string & fieldName,
const std::string & blockId,
937 const std::vector<std::size_t> & localEdgeIds,
const ArrayT & edgeValues,
double scaleValue=1.0);
953 template <
typename ArrayT>
954 void setFaceFieldData(
const std::string & fieldName,
const std::string & blockId,
955 const std::vector<std::size_t> & localFaceIds,
const ArrayT & faceValues,
double scaleValue=1.0);
967 template <
typename ArrayT>
968 void getElementVertices(
const std::vector<std::size_t> & localIds, ArrayT & vertices)
const;
978 template <
typename ArrayT>
979 void getElementVertices(
const std::vector<stk::mesh::Entity> & elements, ArrayT & vertices)
const;
990 template <
typename ArrayT>
991 void getElementVertices(
const std::vector<std::size_t> & localIds,
const std::string & eBlock, ArrayT & vertices)
const;
1002 template <
typename ArrayT>
1003 void getElementVertices(
const std::vector<stk::mesh::Entity> & elements,
const std::string & eBlock, ArrayT & vertices)
const;
1013 template <
typename ArrayT>
1024 template <
typename ArrayT>
1036 template <
typename ArrayT>
1037 void getElementVerticesNoResize(
const std::vector<std::size_t> & localIds,
const std::string & eBlock, ArrayT & vertices)
const;
1048 template <
typename ArrayT>
1049 void getElementVerticesNoResize(
const std::vector<stk::mesh::Entity> & elements,
const std::string & eBlock, ArrayT & vertices)
const;
1062 template <
typename ArrayT>
1063 void getElementNodes(
const std::vector<std::size_t> & localIds, ArrayT & nodes)
const;
1073 template <
typename ArrayT>
1074 void getElementNodes(
const std::vector<stk::mesh::Entity> & elements, ArrayT & nodes)
const;
1085 template <
typename ArrayT>
1086 void getElementNodes(
const std::vector<std::size_t> & localIds,
const std::string & eBlock, ArrayT & nodes)
const;
1097 template <
typename ArrayT>
1098 void getElementNodes(
const std::vector<stk::mesh::Entity> & elements,
const std::string & eBlock, ArrayT & nodes)
const;
1108 template <
typename ArrayT>
1119 template <
typename ArrayT>
1131 template <
typename ArrayT>
1132 void getElementNodesNoResize(
const std::vector<std::size_t> & localIds,
const std::string & eBlock, ArrayT & nodes)
const;
1143 template <
typename ArrayT>
1144 void getElementNodesNoResize(
const std::vector<stk::mesh::Entity> & elements,
const std::string & eBlock, ArrayT & nodes)
const;
1151 stk::mesh::EntityRank
getFaceRank()
const {
return stk::topology::FACE_RANK; }
1152 stk::mesh::EntityRank
getEdgeRank()
const {
return stk::topology::EDGE_RANK; }
1153 stk::mesh::EntityRank
getNodeRank()
const {
return stk::topology::NODE_RANK; }
1173 const std::vector<Teuchos::RCP<const PeriodicBC_MatcherBase> > &
1179 std::vector<Teuchos::RCP<const PeriodicBC_MatcherBase> > &
1219 void print(std::ostream & os)
const;
1291 template <
typename ArrayT>
1292 void getElementVertices_FromField(
const std::vector<stk::mesh::Entity> & elements,
const std::string & eBlock, ArrayT & vertices)
const;
1294 template <
typename ArrayT>
1296 const std::string & eBlock, ArrayT & vertices)
const;
1307 template <
typename ArrayT>
1310 template <
typename ArrayT>
1325 template <
typename ArrayT>
1326 void getElementNodes_FromField(
const std::vector<stk::mesh::Entity> & elements,
const std::string & eBlock, ArrayT & nodes)
const;
1328 template <
typename ArrayT>
1330 const std::string & eBlock, ArrayT & nodes)
const;
1341 template <
typename ArrayT>
1344 template <
typename ArrayT>
1352 void refineMesh(
const int numberOfLevels,
const bool deleteParentElements);
1394 bool isMeshCoordField(
const std::string & eBlock,
const std::string & fieldName,
int & axis)
const;
1411 template <
typename ArrayT>
1412 void setDispFieldData(
const std::string & fieldName,
const std::string & blockId,
int axis,
1413 const std::vector<std::size_t> & localElementIds,
const ArrayT & solutionValues);
1420 #ifdef PANZER_HAVE_PERCEPT
1474 #ifdef PANZER_HAVE_IOSS
1483 enum class GlobalVariable
1507 const GlobalVariable& flag);
1556 template <
typename ArrayT>
1558 const std::vector<std::size_t> & localElementIds,
const ArrayT & solutionValues,
double scaleValue)
1561 auto solutionValues_h = Kokkos::create_mirror_view(solutionValues);
1562 Kokkos::deep_copy(solutionValues_h, solutionValues);
1564 int field_axis = -1;
1566 setDispFieldData(fieldName,blockId,field_axis,localElementIds,solutionValues_h);
1572 for(std::size_t cell=0;cell<localElementIds.size();cell++) {
1573 std::size_t localId = localElementIds[cell];
1574 stk::mesh::Entity element = elements[localId];
1577 const size_t num_nodes =
bulkData_->num_nodes(element);
1578 stk::mesh::Entity
const* nodes =
bulkData_->begin_nodes(element);
1579 for(std::size_t i=0; i<num_nodes; ++i) {
1580 stk::mesh::Entity node = nodes[i];
1582 double * solnData = stk::mesh::field_data(*field,node);
1584 solnData[0] = scaleValue*solutionValues_h(cell,i);
1589 template <
typename ArrayT>
1591 const std::vector<std::size_t> & localElementIds,
const ArrayT & dispValues)
1600 for(std::size_t cell=0;cell<localElementIds.size();cell++) {
1601 std::size_t localId = localElementIds[cell];
1602 stk::mesh::Entity element = elements[localId];
1605 const size_t num_nodes =
bulkData_->num_nodes(element);
1606 stk::mesh::Entity
const* nodes =
bulkData_->begin_nodes(element);
1607 for(std::size_t i=0; i<num_nodes; ++i) {
1608 stk::mesh::Entity node = nodes[i];
1610 double * solnData = stk::mesh::field_data(*field,node);
1611 double * coordData = stk::mesh::field_data(coord_field,node);
1613 solnData[0] = dispValues(cell,i)-coordData[axis];
1618 template <
typename ArrayT>
1620 const std::vector<std::size_t> & localElementIds,ArrayT & solutionValues)
const
1626 localElementIds.size(),
1627 bulkData_->num_nodes(elements[localElementIds[0]]));
1631 for(std::size_t cell=0;cell<localElementIds.size();cell++) {
1632 std::size_t localId = localElementIds[cell];
1633 stk::mesh::Entity element = elements[localId];
1636 const size_t num_nodes =
bulkData_->num_nodes(element);
1637 stk::mesh::Entity
const* nodes =
bulkData_->begin_nodes(element);
1638 for(std::size_t i=0; i<num_nodes; ++i) {
1639 stk::mesh::Entity node = nodes[i];
1641 double * solnData = stk::mesh::field_data(*field,node);
1643 solutionValues(cell,i) = solnData[0];
1648 template <
typename ArrayT>
1650 const std::vector<std::size_t> & localElementIds,
const ArrayT & solutionValues,
double scaleValue)
1656 auto solutionValues_h = Kokkos::create_mirror_view(solutionValues);
1657 Kokkos::deep_copy(solutionValues_h, solutionValues);
1659 for(std::size_t cell=0;cell<localElementIds.size();cell++) {
1660 std::size_t localId = localElementIds[cell];
1661 stk::mesh::Entity element = elements[localId];
1663 double * solnData = stk::mesh::field_data(*field,element);
1665 solnData[0] = scaleValue*solutionValues_h.access(cell,0);
1669 template <
typename ArrayT>
1671 const std::vector<std::size_t> & localEdgeIds,
const ArrayT & edgeValues,
double scaleValue)
1677 auto edgeValues_h = Kokkos::create_mirror_view(edgeValues);
1678 Kokkos::deep_copy(edgeValues_h, edgeValues);
1680 for(std::size_t idx=0;idx<localEdgeIds.size();idx++) {
1681 std::size_t localId = localEdgeIds[idx];
1682 stk::mesh::Entity edge = edges[localId];
1684 double * solnData = stk::mesh::field_data(*field,edge);
1686 solnData[0] = scaleValue*edgeValues_h.access(idx,0);
1690 template <
typename ArrayT>
1692 const std::vector<std::size_t> & localFaceIds,
const ArrayT & faceValues,
double scaleValue)
1698 auto faceValues_h = Kokkos::create_mirror_view(faceValues);
1699 Kokkos::deep_copy(faceValues_h, faceValues);
1701 for(std::size_t idx=0;idx<localFaceIds.size();idx++) {
1702 std::size_t localId = localFaceIds[idx];
1703 stk::mesh::Entity face = faces[localId];
1705 double * solnData = stk::mesh::field_data(*field,face);
1707 solnData[0] = scaleValue*faceValues_h.access(idx,0);
1713 template <
typename ArrayT>
1724 std::vector<stk::mesh::Entity> selected_elements;
1725 for(std::size_t cell=0;cell<localElementIds.size();cell++)
1726 selected_elements.push_back(elements[localElementIds[cell]]);
1732 "STK_Interface::getElementVertices: Cannot call this method when field coordinates are used "
1733 "without specifying an element block.");
1737 template <
typename ArrayT>
1745 "STK_Interface::getElementVertices: Cannot call this method when field coordinates are used "
1746 "without specifying an element block.");
1750 template <
typename ArrayT>
1761 template <
typename ArrayT>
1767 std::vector<stk::mesh::Entity> selected_elements;
1768 for(std::size_t cell=0;cell<localElementIds.size();cell++)
1769 selected_elements.push_back(elements[localElementIds[cell]]);
1779 template <
typename ArrayT>
1790 std::vector<stk::mesh::Entity> selected_elements;
1791 for(std::size_t cell=0;cell<localElementIds.size();cell++)
1792 selected_elements.push_back(elements[localElementIds[cell]]);
1798 "STK_Interface::getElementVerticesNoResize: Cannot call this method when field coordinates are used "
1799 "without specifying an element block.");
1803 template <
typename ArrayT>
1811 "STK_Interface::getElementVerticesNoResize: Cannot call this method when field coordinates are used "
1812 "without specifying an element block.");
1816 template <
typename ArrayT>
1827 template <
typename ArrayT>
1833 std::vector<stk::mesh::Entity> selected_elements;
1834 for(std::size_t cell=0;cell<localElementIds.size();cell++)
1835 selected_elements.push_back(elements[localElementIds[cell]]);
1845 template <
typename ArrayT>
1849 if(elements.size() == 0) {
1859 const auto masterVertexCount
1860 = stk::mesh::get_cell_topology(
bulkData_->bucket(elements[0]).topology()).getCellTopologyData()->vertex_count;
1864 auto vertices_h = Kokkos::create_mirror_view(vertices);
1865 Kokkos::deep_copy(vertices_h, vertices);
1869 for(std::size_t cell = 0; cell < elements.size(); cell++) {
1870 const auto element = elements[cell];
1873 const auto vertexCount
1874 = stk::mesh::get_cell_topology(
bulkData_->bucket(element).topology()).getCellTopologyData()->vertex_count;
1876 "In call to STK_Interface::getElementVertices all elements "
1877 "must have the same vertex count!");
1880 const size_t num_nodes =
bulkData_->num_nodes(element);
1881 auto const* nodes =
bulkData_->begin_nodes(element);
1883 "In call to STK_Interface::getElementVertices cardinality of "
1884 "element node relations must be the vertex count!");
1885 for(std::size_t node = 0; node < num_nodes; ++node) {
1889 for(
unsigned d=0;d<dim;d++)
1890 vertices_h(cell,node,d) = coord[d];
1893 Kokkos::deep_copy(vertices, vertices_h);
1896 template <
typename ArrayT>
1900 if(elements.size()==0) {
1909 unsigned masterVertexCount
1910 = stk::mesh::get_cell_topology(
bulkData_->bucket(elements[0]).topology()).getCellTopologyData()->vertex_count;
1914 auto vertices_h = Kokkos::create_mirror_view(vertices);
1915 for(std::size_t cell=0;cell<elements.size();cell++) {
1916 stk::mesh::Entity element = elements[cell];
1919 unsigned vertexCount
1920 = stk::mesh::get_cell_topology(
bulkData_->bucket(element).topology()).getCellTopologyData()->vertex_count;
1922 "In call to STK_Interface::getElementVertices all elements "
1923 "must have the same vertex count!");
1926 const size_t num_nodes =
bulkData_->num_nodes(element);
1927 stk::mesh::Entity
const* nodes =
bulkData_->begin_nodes(element);
1929 "In call to STK_Interface::getElementVertices cardinality of "
1930 "element node relations must be the vertex count!");
1931 for(std::size_t node=0; node<num_nodes; ++node) {
1935 for(
unsigned d=0;d<dim;d++)
1936 vertices_h(cell,node,d) = coord[d];
1939 Kokkos::deep_copy(vertices, vertices_h);
1942 template <
typename ArrayT>
1948 if(elements.size()==0) {
1954 unsigned masterVertexCount
1955 = stk::mesh::get_cell_topology(
bulkData_->bucket(elements[0]).topology()).getCellTopologyData()->vertex_count;
1959 auto vertices_h = Kokkos::create_mirror_view(vertices);
1960 std::map<std::string,std::vector<std::string> >::const_iterator itr =
meshCoordFields_.find(eBlock);
1966 const std::vector<std::string> & coordField = itr->second;
1967 std::vector<SolutionFieldType*> fields(
getDimension());
1968 for(std::size_t d=0;d<fields.size();d++) {
1972 for(std::size_t cell=0;cell<elements.size();cell++) {
1973 stk::mesh::Entity element = elements[cell];
1976 const size_t num_nodes =
bulkData_->num_nodes(element);
1977 stk::mesh::Entity
const* nodes =
bulkData_->begin_nodes(element);
1978 for(std::size_t i=0; i<num_nodes; ++i) {
1979 stk::mesh::Entity node = nodes[i];
1984 double * solnData = stk::mesh::field_data(*fields[d],node);
1988 vertices_h(cell,i,d) = solnData[0]+coord[d];
1992 Kokkos::deep_copy(vertices, vertices_h);
1995 template <
typename ArrayT>
1997 const std::string & eBlock, ArrayT & vertices)
const
2002 if(elements.size()==0) {
2006 std::map<std::string,std::vector<std::string> >::const_iterator itr =
meshCoordFields_.find(eBlock);
2012 const std::vector<std::string> & coordField = itr->second;
2013 std::vector<SolutionFieldType*> fields(
getDimension());
2014 for(std::size_t d=0;d<fields.size();d++) {
2018 for(std::size_t cell=0;cell<elements.size();cell++) {
2019 stk::mesh::Entity element = elements[cell];
2022 const size_t num_nodes =
bulkData_->num_nodes(element);
2023 stk::mesh::Entity
const* nodes =
bulkData_->begin_nodes(element);
2024 for(std::size_t i=0; i<num_nodes; ++i) {
2025 stk::mesh::Entity node = nodes[i];
2030 double * solnData = stk::mesh::field_data(*fields[d],node);
2034 vertices(cell,i,d) = solnData[0]+coord[d];
2042 template <
typename ArrayT>
2053 std::vector<stk::mesh::Entity> selected_elements;
2054 for(std::size_t cell=0;cell<localElementIds.size();cell++)
2055 selected_elements.push_back(elements[localElementIds[cell]]);
2061 "STK_Interface::getElementNodes: Cannot call this method when field coordinates are used "
2062 "without specifying an element block.");
2066 template <
typename ArrayT>
2074 "STK_Interface::getElementNodes: Cannot call this method when field coordinates are used "
2075 "without specifying an element block.");
2079 template <
typename ArrayT>
2090 template <
typename ArrayT>
2096 std::vector<stk::mesh::Entity> selected_elements;
2097 for(std::size_t cell=0;cell<localElementIds.size();cell++)
2098 selected_elements.push_back(elements[localElementIds[cell]]);
2108 template <
typename ArrayT>
2119 std::vector<stk::mesh::Entity> selected_elements;
2120 for(std::size_t cell=0;cell<localElementIds.size();cell++)
2121 selected_elements.push_back(elements[localElementIds[cell]]);
2127 "STK_Interface::getElementNodesNoResize: Cannot call this method when field coordinates are used "
2128 "without specifying an element block.");
2132 template <
typename ArrayT>
2140 "STK_Interface::getElementNodesNoResize: Cannot call this method when field coordinates are used "
2141 "without specifying an element block.");
2145 template <
typename ArrayT>
2156 template <
typename ArrayT>
2162 std::vector<stk::mesh::Entity> selected_elements;
2163 for(std::size_t cell=0;cell<localElementIds.size();cell++)
2164 selected_elements.push_back(elements[localElementIds[cell]]);
2174 template <
typename ArrayT>
2178 if(elements.size() == 0) {
2188 const auto masterNodeCount
2189 = stk::mesh::get_cell_topology(
bulkData_->bucket(elements[0]).topology()).getCellTopologyData()->node_count;
2193 auto nodes_h = Kokkos::create_mirror_view(nodes);
2194 Kokkos::deep_copy(nodes_h, nodes);
2198 for(std::size_t cell = 0; cell < elements.size(); cell++) {
2199 const auto element = elements[cell];
2202 const auto nodeCount
2203 = stk::mesh::get_cell_topology(
bulkData_->bucket(element).topology()).getCellTopologyData()->node_count;
2205 "In call to STK_Interface::getElementNodes all elements "
2206 "must have the same node count!");
2209 const size_t num_nodes =
bulkData_->num_nodes(element);
2210 auto const* elem_nodes =
bulkData_->begin_nodes(element);
2212 "In call to STK_Interface::getElementNodes cardinality of "
2213 "element node relations must be the node count!");
2214 for(std::size_t node = 0; node < num_nodes; ++node) {
2218 for(
unsigned d=0;d<dim;d++)
2219 nodes_h(cell,node,d) = coord[d];
2222 Kokkos::deep_copy(nodes, nodes_h);
2225 template <
typename ArrayT>
2229 if(elements.size()==0) {
2238 unsigned masterNodeCount
2239 = stk::mesh::get_cell_topology(
bulkData_->bucket(elements[0]).topology()).getCellTopologyData()->node_count;
2243 auto nodes_h = Kokkos::create_mirror_view(nodes);
2244 for(std::size_t cell=0;cell<elements.size();cell++) {
2245 stk::mesh::Entity element = elements[cell];
2249 = stk::mesh::get_cell_topology(
bulkData_->bucket(element).topology()).getCellTopologyData()->node_count;
2251 "In call to STK_Interface::getElementNodes all elements "
2252 "must have the same node count!");
2255 const size_t num_nodes =
bulkData_->num_nodes(element);
2256 stk::mesh::Entity
const* elem_nodes =
bulkData_->begin_nodes(element);
2258 "In call to STK_Interface::getElementNodes cardinality of "
2259 "element node relations must be the node count!");
2260 for(std::size_t node=0; node<num_nodes; ++node) {
2264 for(
unsigned d=0;d<dim;d++)
2265 nodes_h(cell,node,d) = coord[d];
2268 Kokkos::deep_copy(nodes, nodes_h);
2271 template <
typename ArrayT>
2277 if(elements.size()==0) {
2283 unsigned masterNodeCount
2284 = stk::mesh::get_cell_topology(
bulkData_->bucket(elements[0]).topology()).getCellTopologyData()->node_count;
2288 auto nodes_h = Kokkos::create_mirror_view(nodes);
2289 std::map<std::string,std::vector<std::string> >::const_iterator itr =
meshCoordFields_.find(eBlock);
2295 const std::vector<std::string> & coordField = itr->second;
2296 std::vector<SolutionFieldType*> fields(
getDimension());
2297 for(std::size_t d=0;d<fields.size();d++) {
2302 for(std::size_t cell=0;cell<elements.size();cell++) {
2303 stk::mesh::Entity element = elements[cell];
2306 const size_t num_nodes =
bulkData_->num_nodes(element);
2307 stk::mesh::Entity
const* elem_nodes =
bulkData_->begin_nodes(element);
2308 for(std::size_t i=0; i<num_nodes; ++i) {
2309 stk::mesh::Entity node = elem_nodes[i];
2314 double * solnData = stk::mesh::field_data(*fields[d],node);
2318 nodes_h(cell,i,d) = solnData[0]+coord[d];
2322 Kokkos::deep_copy(nodes, nodes_h);
2325 template <
typename ArrayT>
2327 const std::string & eBlock, ArrayT & nodes)
const
2332 if(elements.size()==0) {
2336 std::map<std::string,std::vector<std::string> >::const_iterator itr =
meshCoordFields_.find(eBlock);
2342 const std::vector<std::string> & coordField = itr->second;
2343 std::vector<SolutionFieldType*> fields(
getDimension());
2344 for(std::size_t d=0;d<fields.size();d++) {
2349 for(std::size_t cell=0;cell<elements.size();cell++) {
2350 stk::mesh::Entity element = elements[cell];
2353 const size_t num_nodes =
bulkData_->num_nodes(element);
2354 stk::mesh::Entity
const* elem_nodes =
bulkData_->begin_nodes(element);
2355 for(std::size_t i=0; i<num_nodes; ++i) {
2356 stk::mesh::Entity node = elem_nodes[i];
2361 double * solnData = stk::mesh::field_data(*fields[d],node);
2365 nodes(cell,i,d) = solnData[0]+coord[d];
void initializeFromMetaData()
Teuchos::RCP< ElementDescriptor > buildElementDescriptor(stk::mesh::EntityId elmtId, std::vector< stk::mesh::EntityId > &nodes)
Teuchos::RCP< Teuchos::MpiComm< int > > getSafeCommunicator(stk::ParallelMachine parallelMach) const
void getSidesetNames(std::vector< std::string > &name) const
std::pair< Teuchos::RCP< std::vector< std::pair< std::size_t, std::size_t > > >, Teuchos::RCP< std::vector< unsigned int > > > getPeriodicNodePairing() const
const bool & useBoundingBoxSearch() const
Teuchos::RCP< std::vector< stk::mesh::Entity > > orderedEdgeVector_
Kokkos::DynRankView< typename InputArray::value_type, PHX::Device > createDynRankView(const InputArray &a, const std::string &name, const DimensionPack...dims)
Wrapper to simplify Panzer use of Sacado ViewFactory.
ElementBlockException(const std::string &what)
void addNodeset(const std::string &name)
EdgeBlockException(const std::string &what)
Teuchos::RCP< stk::mesh::BulkData > bulkData_
void getFaceBlockNames(std::vector< std::string > &names) const
void print(std::ostream &os) const
void getElementNodesNoResize(const std::vector< std::size_t > &localIds, ArrayT &nodes) const
unsigned entityOwnerRank(stk::mesh::Entity entity) const
void getMyEdges(std::vector< stk::mesh::Entity > &edges) const
void setEdgeFieldData(const std::string &fieldName, const std::string &blockId, const std::vector< std::size_t > &localEdgeIds, const ArrayT &edgeValues, double scaleValue=1.0)
stk::mesh::Field< double > VectorFieldType
Teuchos::RCP< const std::vector< stk::mesh::Entity > > getElementsOrderedByLID() const
void setSolutionFieldData(const std::string &fieldName, const std::string &blockId, const std::vector< std::size_t > &localElementIds, const ArrayT &solutionValues, double scaleValue=1.0)
std::vector< stk::mesh::Part * > facesPartVec_
std::vector< stk::mesh::EntityId > nodes_
std::map< std::pair< std::string, std::string >, SolutionFieldType * > fieldNameToSolution_
void getElementVerticesNoResize(const std::vector< std::size_t > &localIds, ArrayT &vertices) const
void addEntityToNodeset(stk::mesh::Entity entity, stk::mesh::Part *nodeset)
void addGlobalToExodus(const std::string &key, const int &value)
Add an int global variable to the information to be written to the Exodus output file.
void setDispFieldData(const std::string &fieldName, const std::string &blockId, int axis, const std::vector< std::size_t > &localElementIds, const ArrayT &solutionValues)
void getElementNodes_FromCoordsNoResize(const std::vector< stk::mesh::Entity > &elements, ArrayT &nodes) const
void getElementBlockNames(std::vector< std::string > &names) const
void addPeriodicBC(const Teuchos::RCP< const PeriodicBC_MatcherBase > &bc)
void getSubcellIndices(unsigned entityRank, stk::mesh::EntityId elementId, std::vector< stk::mesh::EntityId > &subcellIds) const
std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > & getPeriodicBCVector()
void getElementsSharingNodes(const std::vector< stk::mesh::EntityId > nodeId, std::vector< stk::mesh::Entity > &elements) const
get a set of elements sharing multiple nodes
stk::mesh::Field< double > SolutionFieldType
stk::mesh::Part * getElementBlockPart(const std::string &name) const
get the block part
std::vector< stk::mesh::Part * > edgesPartVec_
bool nonnull(const std::shared_ptr< T > &p)
bool isFaceLocal(stk::mesh::Entity face) const
std::map< std::string, double > blockWeights_
Teuchos::RCP< const std::vector< stk::mesh::Entity > > getEdgesOrderedByLID() const
VectorFieldType * coordinatesField_
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
std::size_t getEntityCounts(unsigned entityRank) const
get the global counts for the entity of specified rank
void printMetaData(std::ostream &os) const
void getElementNodes_FromFieldNoResize(const std::vector< stk::mesh::Entity > &elements, const std::string &eBlock, ArrayT &nodes) const
void getElementVertices(const std::vector< std::size_t > &localIds, ArrayT &vertices) const
void addEntityToSideset(stk::mesh::Entity entity, stk::mesh::Part *sideset)
std::map< std::string, stk::mesh::Part * > faceBlocks_
stk::mesh::Part * nodesPart_
stk::mesh::EntityId edgeGlobalId(std::size_t lid) const
void getElementVertices_FromCoords(const std::vector< stk::mesh::Entity > &elements, ArrayT &vertices) const
void getElementNodes(const std::vector< std::size_t > &localIds, ArrayT &nodes) const
void addSolutionField(const std::string &fieldName, const std::string &blockId)
std::map< std::string, Teuchos::RCP< shards::CellTopology > > elementBlockCT_
void addEntityToFaceBlock(stk::mesh::Entity entity, stk::mesh::Part *faceblock)
stk::mesh::EntityId getMaxEntityId(unsigned entityRank) const
get max entity ID of type entityRank
std::size_t currentLocalId_
std::map< std::pair< std::string, std::string >, SolutionFieldType * > fieldNameToFaceField_
std::map< std::string, stk::mesh::Part * > sidesets_
std::vector< std::size_t > entityCounts_
bool validBlockId(const std::string &blockId) const
stk::mesh::Part * getNodeset(const std::string &name) const
void getEdgeBlockNames(std::vector< std::string > &names) const
stk::mesh::Part * getFaceBlock(const std::string &name) const
get the block part
void addInformationRecords(const std::vector< std::string > &info_records)
stk::mesh::EntityId faceGlobalId(std::size_t lid) const
void getAllFaces(const std::string &faceBlockName, std::vector< stk::mesh::Entity > &faces) const
stk::mesh::Field< double > * getSolutionField(const std::string &fieldName, const std::string &blockId) const
SolutionFieldType * loadBalField_
void addEdgeBlock(const std::string &elemBlockName, const std::string &edgeBlockName, const stk::topology &topology)
const double * getNodeCoordinates(stk::mesh::EntityId nodeId) const
stk::mesh::Part * getSideset(const std::string &name) const
static const std::string nodesString
void instantiateBulkData(stk::ParallelMachine parallelMach)
void getElementVertices_FromFieldNoResize(const std::vector< stk::mesh::Entity > &elements, const std::string &eBlock, ArrayT &vertices) const
std::set< std::string > informationRecords_
void getMySides(const std::string &sideName, std::vector< stk::mesh::Entity > &sides) const
std::size_t getNumSidesets() const
get the side set count
stk::mesh::EntityId elementGlobalId(stk::mesh::Entity elmt) const
stk::mesh::Field< double > * getEdgeField(const std::string &fieldName, const std::string &blockId) const
std::size_t elementLocalId(stk::mesh::Entity elmt) const
bool useFieldCoordinates_
stk::mesh::EntityId faceGlobalId(stk::mesh::Entity face) const
void getElementVertices_FromCoordsNoResize(const std::vector< stk::mesh::Entity > &elements, ArrayT &vertices) const
unsigned getDimension() const
get the dimension
void addElement(const Teuchos::RCP< ElementDescriptor > &ed, stk::mesh::Part *block)
std::map< std::string, std::vector< std::string > > meshDispFields_
void setupExodusFile(const std::string &filename, const bool append=false, const bool append_after_restart_time=false, const double restart_time=0.0)
Set up an output Exodus file for writing results.
VectorFieldType * edgesField_
void setUseFieldCoordinates(bool useFieldCoordinates)
stk::mesh::EntityRank getSideRank() const
void getAllEdges(const std::string &edgeBlockName, std::vector< stk::mesh::Entity > &edges) const
void initialize(stk::ParallelMachine parallelMach, bool setupIO=true, const bool buildRefinementSupport=false)
FaceBlockException(const std::string &what)
static const std::string coordsString
void setFaceFieldData(const std::string &fieldName, const std::string &blockId, const std::vector< std::size_t > &localFaceIds, const ArrayT &faceValues, double scaleValue=1.0)
std::size_t getNumElementBlocks() const
get the block count
void getElementNodes_FromField(const std::vector< stk::mesh::Entity > &elements, const std::string &eBlock, ArrayT &nodes) const
Teuchos::RCP< stk::mesh::BulkData > getBulkData() const
std::map< std::string, stk::mesh::Part * > edgeBlocks_
void getElementVertices_FromField(const std::vector< stk::mesh::Entity > &elements, const std::string &eBlock, ArrayT &vertices) const
bool getUseFieldCoordinates() const
stk::mesh::EntityId elementGlobalId(std::size_t lid) const
std::map< std::string, stk::mesh::Part * > elementBlocks_
std::vector< stk::mesh::Part * > nodesPartVec_
static const std::string edgeBlockString
static const std::string edgesString
SidesetException(const std::string &what)
std::string containingBlockId(stk::mesh::Entity elmt) const
stk::mesh::EntityRank getFaceRank() const
std::size_t edgeLocalId(stk::mesh::Entity elmt) const
const STK_Interface * mesh_
std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > periodicBCs_
void getOwnedElementsSharingNode(stk::mesh::Entity node, std::vector< stk::mesh::Entity > &elements, std::vector< int > &relIds) const
void addNode(stk::mesh::EntityId gid, const std::vector< double > &coord)
void addEntitiesToEdgeBlock(std::vector< stk::mesh::Entity > entities, stk::mesh::Part *edgeblock)
void setBlockWeight(const std::string &blockId, double weight)
void buildSubcells()
force the mesh to build subcells: edges and faces
Teuchos::RCP< Teuchos::MpiComm< int > > mpiComm_
LocalIdCompare(const STK_Interface *mesh)
void setInitialStateTime(double value)
Teuchos::RCP< stk::mesh::MetaData > metaData_
Teuchos::RCP< stk::mesh::MetaData > getMetaData() const
virtual ~ElementDescriptor()
Teuchos::RCP< std::vector< stk::mesh::Entity > > orderedElementVector_
std::size_t getNumNodesets() const
get the side set count
double getInitialStateTime() const
bool isInitialized() const
Has initialize been called on this mesh object?
void writeToExodus(const std::string &filename, const bool append=false)
Write this mesh and associated fields to the given output file.
void getMyNodes(const std::string &sideName, const std::string &blockName, std::vector< stk::mesh::Entity > &nodes) const
void getSolutionFieldData(const std::string &fieldName, const std::string &blockId, const std::vector< std::size_t > &localElementIds, ArrayT &solutionValues) const
void addEdgeField(const std::string &fieldName, const std::string &blockId)
void addMeshCoordFields(const std::string &blockId, const std::vector< std::string > &coordField, const std::string &dispPrefix)
void applyElementLoadBalanceWeights()
stk::mesh::Entity findConnectivityById(stk::mesh::Entity src, stk::mesh::EntityRank tgt_rank, unsigned rel_id) const
void setUseLowerCaseForIO(bool useLowerCase)
bool isMeshCoordField(const std::string &eBlock, const std::string &fieldName, int &axis) const
void getAllSides(const std::string &sideName, std::vector< stk::mesh::Entity > &sides) const
void addSideset(const std::string &name, const CellTopologyData *ctData)
static const std::string faceBlockString
std::unordered_map< stk::mesh::EntityId, std::size_t > localFaceIDHash_
void addEntityToEdgeBlock(stk::mesh::Entity entity, stk::mesh::Part *edgeblock)
std::map< std::pair< std::string, std::string >, SolutionFieldType * > fieldNameToEdgeField_
stk::mesh::Part * getEdgeBlock(const std::string &name) const
get the block part
void setCellFieldData(const std::string &fieldName, const std::string &blockId, const std::vector< std::size_t > &localElementIds, const ArrayT &solutionValues, double scaleValue=1.0)
const VectorFieldType & getEdgesField() const
std::unordered_map< stk::mesh::EntityId, std::size_t > localEdgeIDHash_
ProcIdFieldType * getProcessorIdField()
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 rebalance(const Teuchos::ParameterList ¶ms)
std::map< std::string, std::vector< std::string > > meshCoordFields_
stk::mesh::EntityRank getEdgeRank() const
static const std::string facesString
void addEntitiesToFaceBlock(std::vector< stk::mesh::Entity > entities, stk::mesh::Part *faceblock)
void initializeFieldsInSTK(const std::map< std::pair< std::string, std::string >, SolutionFieldType * > &nameToField, bool setupIO)
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
get the comm associated with this mesh
const std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > & getPeriodicBCVector() const
const VectorFieldType & getFacesField() const
void refineMesh(const int numberOfLevels, const bool deleteParentElements)
stk::mesh::Field< double > * getFaceField(const std::string &fieldName, const std::string &blockId) const
void buildLocalElementIDs()
std::size_t faceLocalId(stk::mesh::Entity elmt) const
bool isModifiable() const
void addFaceBlock(const std::string &elemBlockName, const std::string &faceBlockName, const stk::topology &topology)
Teuchos::RCP< std::vector< stk::mesh::Entity > > orderedFaceVector_
std::unordered_map< stk::mesh::EntityId, std::size_t > localIDHash_
bool isValid(stk::mesh::Entity entity) const
std::vector< stk::mesh::EntityId > maxEntityId_
stk::mesh::Field< ProcIdData > ProcIdFieldType
void addPeriodicBCs(const std::vector< Teuchos::RCP< const PeriodicBC_MatcherBase > > &bc_vec)
void getElementsSharingNode(stk::mesh::EntityId nodeId, std::vector< stk::mesh::Entity > &elements) const
get a set of elements sharing a single node
void addFaceField(const std::string &fieldName, const std::string &blockId)
#define TEUCHOS_ASSERT(assertion_test)
void getMyElements(std::vector< stk::mesh::Entity > &elements) const
stk::mesh::EntityRank getNodeRank() const
void getMyFaces(std::vector< stk::mesh::Entity > &faces) const
void addCellField(const std::string &fieldName, const std::string &blockId)
void getNodeIdsForElement(stk::mesh::Entity element, std::vector< stk::mesh::EntityId > &nodeIds) const
get a list of node ids for nodes connected to an element
Teuchos::RCP< const std::vector< stk::mesh::Entity > > getFacesOrderedByLID() const
stk::mesh::Part * edgesPart_
bool operator()(stk::mesh::Entity a, stk::mesh::Entity b)
stk::mesh::EntityRank getElementRank() const
void getElementNodes_FromCoords(const std::vector< stk::mesh::Entity > &elements, ArrayT &nodes) const
stk::mesh::EntityId edgeGlobalId(stk::mesh::Entity edge) const
stk::mesh::Part * facesPart_
ProcIdFieldType * processorIdField_
void addElementBlock(const std::string &name, const CellTopologyData *ctData)
const std::vector< stk::mesh::EntityId > & getNodes() const
stk::mesh::Field< double > * getCellField(const std::string &fieldName, const std::string &blockId) const
const VectorFieldType & getCoordinatesField() const
stk::mesh::Part * getOwnedPart() const
Get a pointer to the locally owned part.
stk::mesh::EntityId getGID() const
void setBoundingBoxSearchFlag(const bool &searchFlag)
bool getUseLowerCaseForIO() const
bool isEdgeLocal(stk::mesh::Entity edge) const
double getCurrentStateTime() const
void getNeighborElements(std::vector< stk::mesh::Entity > &elements) const
void getNodesetNames(std::vector< std::string > &name) const
Teuchos::RCP< const shards::CellTopology > getCellTopology(const std::string &eBlock) const
std::map< std::pair< std::string, std::string >, SolutionFieldType * > fieldNameToCellField_
std::map< std::string, stk::mesh::Part * > nodesets_
VectorFieldType * facesField_