From: Konstantin Leontev Date: Thu, 1 Aug 2024 19:11:10 +0000 (+0100) Subject: [bos #42217][EDF 28921] Horseshoe with bodyfitting. X-Git-Url: http://git.salome-platform.org/gitweb/?a=commitdiff_plain;h=refs%2Fheads%2Fkleontev%2F42217_horseshoe_with_bodyfitting;p=modules%2Fsmesh.git [bos #42217][EDF 28921] Horseshoe with bodyfitting. Fix for prev commit: crash on parallel execution caused reallocation of vector after push_back(). --- diff --git a/src/StdMeshers/StdMeshers_Cartesian_3D.cxx b/src/StdMeshers/StdMeshers_Cartesian_3D.cxx index 34bf2fb61..26b748e07 100644 --- a/src/StdMeshers/StdMeshers_Cartesian_3D.cxx +++ b/src/StdMeshers/StdMeshers_Cartesian_3D.cxx @@ -2825,23 +2825,24 @@ namespace bool Hexahedron::isSplittedLink(const Solid* solid, SMESH_MesherHelper& helper, const Hexahedron::_Link& linkIn) const { _Link split; - std::vector<_Node> intNodes; _Link link = linkIn; link._fIntNodes.clear(); link._fIntNodes.reserve( link._fIntPoints.size() ); + std::vector<_Node> intNodes; + intNodes.reserve(link._fIntPoints.size()); + for ( size_t i = 0; i < link._fIntPoints.size(); ++i ) if ( solid->ContainsAny( link._fIntPoints[i]->_faceIDs )) { - intNodes.push_back( _Node( 0, link._fIntPoints[i] )); + intNodes.emplace_back(nullptr, link._fIntPoints[i]); link._fIntNodes.push_back( & intNodes.back() ); } - link._splits.clear(); split._nodes[ 0 ] = link._nodes[0]; bool isOut = ( ! link._nodes[0]->Node() ); - bool checkTransition; + bool checkTransition = false; for ( size_t i = 0; i < link._fIntNodes.size(); ++i ) { const bool isGridNode = ( ! link._fIntNodes[i]->Node() ); @@ -4077,12 +4078,14 @@ bool Hexahedron::createVolume(const Solid* solid) { // A simple case - just one task executed in one thread. // TODO: check if it's faster to do it sequentially - // threads.reserve(numTasksTotal); - // for (; it < last; ++it) - // { - // threads.emplace_back(f, std::ref(*it)); - // } - std::for_each(it, last, f); + threads.reserve(numTasksTotal); + for (; it < last; ++it) + { + threads.emplace_back(f, std::ref(*it)); + } + + // This line for debug in sequential mode + // std::for_each(it, last, f); } else {