49 #include "Panzer_Workset_Builder.hpp"
55 #include <stk_mesh/base/Selector.hpp>
56 #include <stk_mesh/base/GetEntities.hpp>
57 #include <stk_mesh/base/GetBuckets.hpp>
58 #include <stk_mesh/base/CreateAdjacentEntities.hpp>
60 #include "Shards_CellTopology.hpp"
62 #include "Intrepid2_CellTools_Serial.hpp"
63 #include "Teuchos_Assert.hpp"
65 namespace panzer_stk {
69 const std::string& sidesetName,
70 const std::string& elementBlockName,
86 stk::mesh::Part * sidePart = mesh->
getSideset(sidesetName);
88 stk::mesh::Selector sideSelector = *sidePart;
89 stk::mesh::Selector blockSelector = *elmtPart;
90 stk::mesh::Selector mySelector = metaData->universal_part() & blockSelector & sideSelector;
91 std::vector<stk::mesh::Entity> sides;
92 stk::mesh::get_selected_entities(mySelector,bulkData->buckets(metaData->side_rank()),sides);
94 std::vector<std::size_t> localSideTopoIDs;
95 std::vector<stk::mesh::Entity> parentElements;
99 for (std::size_t i=0; i < localSideTopoIDs.size(); ++i) {
100 *pout <<
"parent element: "
101 <<
" gid(" << bulkData->identifier(parentElements[i]) <<
")"
102 <<
", local_face(" << localSideTopoIDs[i] <<
")"
111 std::unordered_map<unsigned,std::vector<double> > nodeNormals;
120 std::vector<stk::mesh::Entity>::const_iterator side = sides.begin();
121 std::vector<std::size_t>::const_iterator sideID = localSideTopoIDs.begin();
122 std::vector<stk::mesh::Entity>::const_iterator parentElement = parentElements.begin();
127 Intrepid2::Impl::CellTools::setSubcellParametrization();
128 Kokkos::DynRankView<double,Kokkos::HostSpace> normal_at_point(
"normal",3);
129 for ( ; sideID != localSideTopoIDs.end(); ++side,++sideID,++parentElement) {
131 std::vector<stk::mesh::Entity> elementEntities;
132 elementEntities.push_back(*parentElement);
133 PHX::MDField<double,panzer::Cell,panzer::NODE,panzer::Dim> vertices
146 auto jac_at_point = Kokkos::subview(iv.
jac.get_view(), 0, 0, Kokkos::ALL(), Kokkos::ALL());
148 CellTools::Serial::getPhysicalSideNormal(normal_at_point, jac_at_point, *sideID, *(ir->
topology));
154 *pout <<
"element normals: "
155 <<
"gid(" << bulkData->identifier(*parentElement) <<
")"
156 <<
", normal(" << normal_at_point(0) <<
"," << normal_at_point(1) <<
"," << normal_at_point(2) <<
")"
161 const size_t numNodes = bulkData->num_nodes(*side);
162 stk::mesh::Entity
const* nodeRelations = bulkData->begin_nodes(*side);
163 for (
size_t n=0; n<numNodes; ++n) {
164 stk::mesh::Entity node = nodeRelations[n];
165 for (
unsigned dim = 0; dim < parentTopology->getDimension(); ++dim) {
166 nodeNormals[bulkData->identifier(node)].push_back(normal_at_point(dim));
174 for (std::unordered_map<
unsigned,std::vector<double> >::const_iterator node = nodeNormals.begin(); node != nodeNormals.end(); ++node) {
176 TEUCHOS_ASSERT( (node->second.size() % parentTopology->getDimension()) == 0);
178 int numSurfaceContributions = node->second.size() / parentTopology->getDimension();
179 std::vector<double> contribArea(numSurfaceContributions);
180 double totalArea = 0.0;
181 for (
int surface = 0; surface < numSurfaceContributions; ++surface) {
184 for (
unsigned dim = 0; dim < parentTopology->getDimension(); ++dim)
185 sum += (node->second[surface*parentTopology->getDimension() + dim]) *
186 (node->second[surface*parentTopology->getDimension() + dim]);
188 contribArea[surface] = std::sqrt(sum);
190 totalArea += contribArea[surface];
194 for (std::size_t i = 0; i < contribArea.size(); ++i)
195 contribArea[i] /= totalArea;
198 normals[node->first].resize(parentTopology->getDimension());
199 for (
unsigned dim = 0; dim < parentTopology->getDimension(); ++dim) {
200 normals[node->first][dim] = 0.0;
201 for (
int surface = 0; surface < numSurfaceContributions; ++surface) {
202 normals[node->first][dim] += node->second[surface*parentTopology->getDimension() + dim] * contribArea[surface] / totalArea;
207 *pout <<
"surface normal before normalization: "
208 <<
"gid(" << node->first <<
")"
209 <<
", normal(" << normals[node->first][0] <<
"," << normals[node->first][1] <<
"," << normals[node->first][2] <<
")"
214 for (
unsigned dim = 0; dim < parentTopology->getDimension(); ++dim)
215 sum += normals[node->first][dim] * normals[node->first][dim];
217 for (
unsigned dim = 0; dim < parentTopology->getDimension(); ++dim)
218 normals[node->first][dim] /= std::sqrt(sum);
221 *pout <<
"surface normal after normalization: "
222 <<
"gid(" << node->first <<
")"
224 << normals[node->first][0] <<
","
225 << normals[node->first][1] <<
","
226 << normals[node->first][2] <<
")"
236 const std::string& sidesetName,
237 const std::string& elementBlockName,
243 std::unordered_map<unsigned,std::vector<double> > nodeEntityIdToNormals;
251 stk::mesh::Part * sidePart = mesh->
getSideset(sidesetName);
253 stk::mesh::Selector sideSelector = *sidePart;
254 stk::mesh::Selector blockSelector = *elmtPart;
255 stk::mesh::Selector mySelector = metaData->universal_part() & blockSelector & sideSelector;
256 std::vector<stk::mesh::Entity> sides;
257 stk::mesh::get_selected_entities(mySelector,bulkData->buckets(metaData->side_rank()),sides);
261 std::vector<std::size_t> localSideTopoIDs;
262 std::vector<stk::mesh::Entity> parentElements;
265 std::vector<stk::mesh::Entity>::const_iterator side = sides.begin();
266 std::vector<std::size_t>::const_iterator sideID = localSideTopoIDs.begin();
267 std::vector<stk::mesh::Entity>::const_iterator parentElement = parentElements.begin();
268 for ( ; sideID != localSideTopoIDs.end(); ++side,++sideID,++parentElement) {
271 const size_t numNodes = bulkData->num_nodes(*parentElement);
272 stk::mesh::Entity
const* nodeRelations = bulkData->begin_nodes(*parentElement);
274 normals[mesh->
elementLocalId(*parentElement)] = Kokkos::DynRankView<double,PHX::Device>(
"normals",numNodes,parentTopology->getDimension());
276 for (
size_t nodeIndex=0; nodeIndex<numNodes; ++nodeIndex) {
277 stk::mesh::Entity node = nodeRelations[nodeIndex];
280 if (nodeEntityIdToNormals.find(bulkData->identifier(node)) != nodeEntityIdToNormals.end()) {
281 for (
unsigned dim = 0; dim < parentTopology->getDimension(); ++dim) {
282 (normals[mesh->
elementLocalId(*parentElement)])(nodeIndex,dim) = (nodeEntityIdToNormals[bulkData->identifier(node)])[dim];
286 for (
unsigned dim = 0; dim < parentTopology->getDimension(); ++dim) {
287 (normals[mesh->
elementLocalId(*parentElement)])(nodeIndex,dim) = 0.0;
void computeSidesetNodeNormals(std::unordered_map< unsigned, std::vector< double > > &normals, const Teuchos::RCP< const panzer_stk::STK_Interface > &mesh, const std::string &sidesetName, const std::string &elementBlockName, std::ostream *, std::ostream *pout)
Computes the normals for all nodes associated with a sideset surface.
void getElementVerticesNoResize(const std::vector< std::size_t > &localIds, ArrayT &vertices) const
stk::mesh::Part * getElementBlockPart(const std::string &name) const
get the block count
PHX::MDField< Scalar, T0 > buildStaticArray(const std::string &str, int d0) const
stk::mesh::Part * getSideset(const std::string &name) const
std::size_t elementLocalId(stk::mesh::Entity elmt) const
unsigned getDimension() const
get the dimension
void getUniversalSubcellElements(const panzer_stk::STK_Interface &mesh, const std::string &blockId, const std::vector< stk::mesh::Entity > &entities, std::vector< std::size_t > &localEntityIds, std::vector< stk::mesh::Entity > &elements)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Data for determining cell topology and dimensionality.
Teuchos::RCP< stk::mesh::BulkData > getBulkData() const
void evaluateValues(const PHX::MDField< Scalar, Cell, NODE, Dim > &vertex_coordinates, const int num_cells=-1)
Evaluate basis values.
Teuchos::RCP< const shards::CellTopology > topology
void setupArrays(const Teuchos::RCP< const panzer::IntegrationRule > &ir)
Sizes/allocates memory for arrays.
Teuchos::RCP< stk::mesh::MetaData > getMetaData() const
#define TEUCHOS_ASSERT(assertion_test)
Teuchos::RCP< const shards::CellTopology > getCellTopology(const std::string &eBlock) const