1 // Copyright (C) 2007-2023 CEA, EDF
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License, or (at your option) any later version.
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
11 // Lesser General Public License for more details.
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 // File : HYBRIDPlugin_HYBRID.cxx
22 // Author : Christian VAN WAMBEKE (CEA) (from GHS3D plugin V730)
25 #include "HYBRIDPlugin_HYBRID.hxx"
26 #include "HYBRIDPlugin_Hypothesis.hxx"
27 #include "MG_HYBRID_API.hxx"
29 #include <SMDS_FaceOfNodes.hxx>
30 #include <SMDS_LinearEdge.hxx>
31 #include <SMDS_VolumeOfNodes.hxx>
32 #include <SMESHDS_Group.hxx>
33 #include <SMESHDS_Mesh.hxx>
34 #include <SMESH_Comment.hxx>
35 #include <SMESH_File.hxx>
36 #include <SMESH_Group.hxx>
37 #include <SMESH_HypoFilter.hxx>
38 #include <SMESH_Mesh.hxx>
39 #include <SMESH_MeshAlgos.hxx>
40 #include <SMESH_MeshEditor.hxx>
41 #include <SMESH_MesherHelper.hxx>
42 #include <SMESH_ProxyMesh.hxx>
43 #include <SMESH_subMeshEventListener.hxx>
45 #include <BRepAdaptor_Surface.hxx>
46 #include <BRepBndLib.hxx>
47 #include <BRepBuilderAPI_MakeVertex.hxx>
48 #include <BRepClass3d.hxx>
49 #include <BRepClass3d_SolidClassifier.hxx>
50 #include <BRepExtrema_DistShapeShape.hxx>
51 #include <BRepGProp.hxx>
52 #include <BRepTools.hxx>
53 #include <BRep_Tool.hxx>
54 #include <Bnd_Box.hxx>
55 #include <GProp_GProps.hxx>
56 #include <GeomAPI_ProjectPointOnSurf.hxx>
57 #include <Precision.hxx>
58 #include <Standard_ErrorHandler.hxx>
59 #include <Standard_Failure.hxx>
60 #include <Standard_ProgramError.hxx>
62 #include <TopExp_Explorer.hxx>
63 #include <TopTools_IndexedMapOfShape.hxx>
64 #include <TopTools_ListIteratorOfListOfShape.hxx>
65 #include <TopTools_MapOfShape.hxx>
67 #include <TopoDS_Shell.hxx>
68 #include <TopoDS_Solid.hxx>
70 #include <Basics_Utils.hxx>
71 #include <utilities.h>
73 #include <boost/filesystem.hpp>
74 namespace boofs = boost::filesystem;
79 #define castToNode(n) static_cast<const SMDS_MeshNode *>( n );
82 #define GMFVERSION GmfDouble
84 #define GMFDIMENSION 3
88 typedef const std::list<const SMDS_MeshFace*> TTriaList;
90 static const char theDomainGroupNamePrefix[] = "Domain_";
92 static void removeFile( const std::string& fileName )
95 SMESH_File( fileName ).remove();
98 MESSAGE("Can't remove file: " << fileName << " ; file does not exist or permission denied");
102 // change results files permissions to user only (using boost to be used without C++17)
103 static void chmodUserOnly(const char* filename)
105 if (boofs::exists(filename))
106 boofs::permissions(filename, boofs::remove_perms | boofs::group_all | boofs::others_all );
109 //=============================================================================
113 //=============================================================================
115 HYBRIDPlugin_HYBRID::HYBRIDPlugin_HYBRID(int hypId, SMESH_Gen* gen)
116 : SMESH_3D_Algo(hypId, gen)
119 _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID);// 1 bit /shape type
120 _onlyUnaryInput = true; // Compute() will be called on each solid
123 _compatibleHypothesis.push_back( HYBRIDPlugin_Hypothesis::GetHypType());
124 _requireShape = false; // can work without shape
125 _computeCanceled = false;
128 //=============================================================================
132 //=============================================================================
134 HYBRIDPlugin_HYBRID::~HYBRIDPlugin_HYBRID()
138 //=============================================================================
142 //=============================================================================
144 bool HYBRIDPlugin_HYBRID::CheckHypothesis ( SMESH_Mesh& aMesh,
145 const TopoDS_Shape& aShape,
146 Hypothesis_Status& aStatus )
148 aStatus = SMESH_Hypothesis::HYP_OK;
152 _removeLogOnSuccess = true;
153 _logInStandardOutput = false;
155 const std::list <const SMESHDS_Hypothesis * >& hyps =
156 GetUsedHypothesis(aMesh, aShape, /*ignoreAuxiliary=*/false);
157 std::list <const SMESHDS_Hypothesis* >::const_iterator h = hyps.begin();
158 for ( ; h != hyps.end(); ++h )
161 _hyp = dynamic_cast< const HYBRIDPlugin_Hypothesis*> ( *h );
165 _keepFiles = _hyp->GetKeepFiles();
166 _removeLogOnSuccess = _hyp->GetRemoveLogOnSuccess();
167 _logInStandardOutput = _hyp->GetStandardOutputLog();
174 //=======================================================================
175 //function : entryToShape
177 //=======================================================================
179 TopoDS_Shape HYBRIDPlugin_HYBRID::entryToShape(std::string entry)
181 if ( SMESH_Gen_i::GetSMESHGen()->getStudyServant()->_is_nil() )
182 throw SALOME_Exception("MG-HYBRID plugin can't work w/o publishing in the study");
183 GEOM::GEOM_Object_var aGeomObj;
184 TopoDS_Shape S = TopoDS_Shape();
185 SALOMEDS::SObject_var aSObj = SMESH_Gen_i::GetSMESHGen()->getStudyServant()->FindObjectID( entry.c_str() );
186 if (!aSObj->_is_nil() ) {
187 CORBA::Object_var obj = aSObj->GetObject();
188 aGeomObj = GEOM::GEOM_Object::_narrow(obj);
191 if ( !aGeomObj->_is_nil() )
192 S = SMESH_Gen_i::GetSMESHGen()->GeomObjectToShape( aGeomObj.in() );
196 //=======================================================================
197 //function : addElemInMeshGroup
198 //purpose : Update or create groups in mesh
199 //=======================================================================
201 static void addElemInMeshGroup(SMESH_Mesh* theMesh,
202 const SMDS_MeshElement* anElem,
203 std::string& groupName,
204 std::set<std::string>& /*groupsToRemove*/)
206 if ( !anElem ) return; // issue 0021776
208 bool groupDone = false;
209 SMESH_Mesh::GroupIteratorPtr grIt = theMesh->GetGroups();
210 while (grIt->more()) {
211 SMESH_Group * group = grIt->next();
212 if ( !group ) continue;
213 SMESHDS_GroupBase* groupDS = group->GetGroupDS();
214 if ( !groupDS ) continue;
215 if ( groupDS->GetType()==anElem->GetType() &&groupName.compare(group->GetName())==0) {
216 SMESHDS_Group* aGroupDS = static_cast<SMESHDS_Group*>( groupDS );
217 aGroupDS->SMDSGroup().Add(anElem);
225 SMESH_Group* aGroup = theMesh->AddGroup(anElem->GetType(), groupName.c_str());
226 aGroup->SetName( groupName.c_str() );
227 SMESHDS_Group* aGroupDS = static_cast<SMESHDS_Group*>( aGroup->GetGroupDS() );
228 aGroupDS->SMDSGroup().Add(anElem);
232 throw SALOME_Exception(LOCALIZED("A given element was not added to a group"));
236 //=======================================================================
237 //function : updateMeshGroups
238 //purpose : Update or create groups in mesh
239 //=======================================================================
241 static void updateMeshGroups(SMESH_Mesh* theMesh, std::set<std::string> groupsToRemove)
243 SMESH_Mesh::GroupIteratorPtr grIt = theMesh->GetGroups();
244 while (grIt->more()) {
245 SMESH_Group * group = grIt->next();
246 if ( !group ) continue;
247 SMESHDS_GroupBase* groupDS = group->GetGroupDS();
248 if ( !groupDS ) continue;
249 std::string currentGroupName = (std::string)group->GetName();
250 if (groupDS->IsEmpty() && groupsToRemove.find(currentGroupName) != groupsToRemove.end()) {
251 // Previous group created by enforced elements
252 theMesh->RemoveGroup(groupDS->GetID());
257 //=======================================================================
258 //function : removeEmptyGroupsOfDomains
259 //purpose : remove empty groups named "Domain_nb" created due to
260 // "To make groups of domains" option.
261 //=======================================================================
263 static void removeEmptyGroupsOfDomains(SMESH_Mesh* mesh,
264 bool notEmptyAsWell = false)
266 const char* refName = theDomainGroupNamePrefix;
267 const size_t refLen = strlen( theDomainGroupNamePrefix );
269 std::list<int> groupIDs = mesh->GetGroupIds();
270 std::list<int>::const_iterator id = groupIDs.begin();
271 for ( ; id != groupIDs.end(); ++id )
273 SMESH_Group* group = mesh->GetGroup( *id );
274 if ( !group || ( !group->GetGroupDS()->IsEmpty() && !notEmptyAsWell ))
276 const char* name = group->GetName();
279 if ( strncmp( name, refName, refLen ) == 0 && // starts from refName;
280 isdigit( *( name + refLen )) && // refName is followed by a digit;
281 strtol( name + refLen, &end, 10) >= 0 && // there are only digits ...
282 *end == '\0') // ... till a string end.
284 mesh->RemoveGroup( *id );
289 //================================================================================
291 * \brief Create the groups corresponding to domains
293 //================================================================================
295 static void makeDomainGroups( std::vector< std::vector< const SMDS_MeshElement* > >& elemsOfDomain,
296 SMESH_MesherHelper* theHelper)
298 for ( size_t iDomain = 0; iDomain < elemsOfDomain.size(); ++iDomain )
300 std::vector< const SMDS_MeshElement* > & elems = elemsOfDomain[ iDomain ];
301 if ( elems.empty() ) continue;
303 // find existing groups
304 std::vector< SMESH_Group* > groupOfType( SMDSAbs_NbElementTypes, (SMESH_Group*)NULL );
305 const std::string domainName = ( SMESH_Comment( theDomainGroupNamePrefix ) << iDomain );
306 SMESH_Mesh::GroupIteratorPtr groupIt = theHelper->GetMesh()->GetGroups();
307 while ( groupIt->more() )
309 SMESH_Group* group = groupIt->next();
310 if ( domainName == group->GetName() &&
311 dynamic_cast< SMESHDS_Group* >( group->GetGroupDS()) )
312 groupOfType[ group->GetGroupDS()->GetType() ] = group;
314 // create and fill the groups
318 SMESH_Group* group = groupOfType[ elems[ iElem ]->GetType() ];
320 group = theHelper->GetMesh()->AddGroup( elems[ iElem ]->GetType(),
321 domainName.c_str() );
322 SMDS_MeshGroup& groupDS =
323 static_cast< SMESHDS_Group* >( group->GetGroupDS() )->SMDSGroup();
325 while ( iElem < elems.size() && groupDS.Add( elems[iElem] ))
328 } while ( iElem < elems.size() );
332 //=======================================================================
333 //function : readGMFFile
334 //purpose : read GMF file w/o geometry associated to mesh
335 //=======================================================================
337 static bool readGMFFile(MG_HYBRID_API* MGOutput,
339 HYBRIDPlugin_HYBRID* theAlgo,
340 SMESH_MesherHelper* theHelper,
341 std::vector <const SMDS_MeshNode*> & theNodeByHybridId,
342 std::vector <const SMDS_MeshElement*> & theFaceByHybridId,
343 std::map<const SMDS_MeshNode*,int> & /*theNodeToHybridIdMap*/,
344 std::vector<std::string> & aNodeGroupByHybridId,
345 std::vector<std::string> & anEdgeGroupByHybridId,
346 std::vector<std::string> & aFaceGroupByHybridId,
347 std::set<std::string> & groupsToRemove,
348 bool toMakeGroupsOfDomains=false,
349 bool /*toMeshHoles*/=true)
352 SMESHDS_Mesh* theMeshDS = theHelper->GetMeshDS();
353 const bool hasGeom = ( theHelper->GetMesh()->HasShapeToMesh() );
355 // if imprinting, the original mesh faces are modified
356 // => we clear all the faces to retrieve them from Hybrid output mesh.
357 std::vector<int> facesWithImprinting;
358 if (theAlgo->getHyp())
359 facesWithImprinting = theAlgo->getHyp()->GetFacesWithImprinting();
361 if ( ! facesWithImprinting.empty() ) {
363 std::cout << "Imprinting => Clear original mesh" << std::endl;
365 SMESH_subMesh* smOfSolid =
366 theHelper->GetMesh()->GetSubMesh( theHelper->GetSubShape() );
367 SMESH_subMeshIteratorPtr smIt =
368 smOfSolid->getDependsOnIterator(/*includeSelf=*/false, /*complexShapeFirst=*/true);
369 while ( smIt->more() )
371 SMESH_subMesh* sm = smIt->next();
372 if ( SMESHDS_SubMesh * smDS = sm->GetSubMeshDS() )
374 SMDS_ElemIteratorPtr eIt = smDS->GetElements();
377 theMeshDS->RemoveFreeElement( eIt->next(), smDS );
379 SMDS_NodeIteratorPtr nIt = smDS->GetNodes();
380 while ( nIt->more() )
382 const SMDS_MeshNode* n = nIt->next();
383 if ( n->NbInverseElements() == 0 )
384 theMeshDS->RemoveFreeNode( n, smDS );
388 theNodeByHybridId.clear();
389 theFaceByHybridId.clear();
392 int nbInitialNodes = theNodeByHybridId.size();
395 int nbMeshNodes = theMeshDS->NbNodes();
396 const bool isQuadMesh =
397 theHelper->GetMesh()->NbEdges( ORDER_QUADRATIC ) ||
398 theHelper->GetMesh()->NbFaces( ORDER_QUADRATIC ) ||
399 theHelper->GetMesh()->NbVolumes( ORDER_QUADRATIC );
401 std::cout << "theNodeByHybridId.size(): " << nbInitialNodes << std::endl;
402 std::cout << "theHelper->GetMesh()->NbNodes(): " << nbMeshNodes << std::endl;
403 std::cout << "isQuadMesh: " << isQuadMesh << std::endl;
406 // ---------------------------------
407 // Read generated elements and nodes
408 // ---------------------------------
410 int nbElem = 0, nbRef = 0;
412 std::vector< const SMDS_MeshNode* > GMFNode;
414 std::map<int, std::set<int> > subdomainId2tetraId;
416 std::map <GmfKwdCod,int> tabRef;
417 const bool force3d = !hasGeom;
420 tabRef[GmfVertices] = 3; // for new nodes and enforced nodes
421 tabRef[GmfCorners] = 1;
422 tabRef[GmfEdges] = 2; // for enforced edges
423 tabRef[GmfRidges] = 1;
424 tabRef[GmfTriangles] = 3; // for enforced faces
425 tabRef[GmfQuadrilaterals] = 4;
426 tabRef[GmfTetrahedra] = 4; // for new tetras
427 tabRef[GmfPyramids] = 5; // for new pyramids
428 tabRef[GmfPrisms] = 6; // for new prisms
429 tabRef[GmfHexahedra] = 8;
432 int InpMsh = MGOutput->GmfOpenMesh(theFile, GmfRead, &ver, &dim);
436 // Hybrid is not multi-domain => We can't (and don't need to) read ids of domains in ouput file like in GHS3DPlugin
437 // We just need to get the id of the one and only solid
441 if ( theHelper->GetSubShape().ShapeType() == TopAbs_SOLID )
442 solidID = theHelper->GetSubShapeID();
444 solidID = theMeshDS->ShapeToIndex
445 ( TopExp_Explorer( theHelper->GetSubShape(), TopAbs_SOLID ).Current() );
448 // Issue 0020682. Avoid creating nodes and tetras at place where
449 // volumic elements already exist
450 SMESH_ElementSearcher* elemSearcher = 0;
451 std::vector< const SMDS_MeshElement* > foundVolumes;
452 if ( !hasGeom && theHelper->GetMesh()->NbVolumes() > 0 )
453 elemSearcher = SMESH_MeshAlgos::GetElementSearcher( *theMeshDS );
454 SMESHUtils::Deleter< SMESH_ElementSearcher > elemSearcherDeleter( elemSearcher );
456 // IMP 0022172: [CEA 790] create the groups corresponding to domains
457 std::vector< std::vector< const SMDS_MeshElement* > > elemsOfDomain;
459 int nbVertices = MGOutput->GmfStatKwd(InpMsh, GmfVertices) - nbInitialNodes;
460 if ( nbVertices < 0 )
462 GMFNode.resize( nbVertices + 1 );
464 std::map <GmfKwdCod,int>::const_iterator it = tabRef.begin();
465 for ( ; it != tabRef.end() ; ++it)
467 if(theAlgo->computeCanceled()) {
468 MGOutput->GmfCloseMesh(InpMsh);
472 GmfKwdCod token = it->first;
475 nbElem = MGOutput->GmfStatKwd(InpMsh, token);
477 MGOutput->GmfGotoKwd(InpMsh, token);
478 std::cout << "Read " << nbElem;
483 std::vector<int> id (nbElem*tabRef[token]); // node ids
484 std::vector<int> domainID( nbElem ); // domain
486 if (token == GmfVertices) {
487 (nbElem <= 1) ? tmpStr = " vertex" : tmpStr = " vertices";
492 const SMDS_MeshNode * aGMFNode;
494 for ( int iElem = 0; iElem < nbElem; iElem++ ) {
495 if(theAlgo->computeCanceled()) {
496 MGOutput->GmfCloseMesh(InpMsh);
499 if (ver == GmfFloat) {
500 MGOutput->GmfGetLin(InpMsh, token, &VerTab_f[0], &VerTab_f[1], &VerTab_f[2], &dummy);
506 MGOutput->GmfGetLin(InpMsh, token, &x, &y, &z, &dummy);
508 if (iElem >= nbInitialNodes) {
510 elemSearcher->FindElementsByPoint( gp_Pnt(x,y,z), SMDSAbs_Volume, foundVolumes))
513 aGMFNode = theHelper->AddNode(x, y, z);
515 aGMFID = iElem -nbInitialNodes +1;
516 GMFNode[ aGMFID ] = aGMFNode;
517 if (aGMFID-1 < (int)aNodeGroupByHybridId.size() && !aNodeGroupByHybridId.at(aGMFID-1).empty())
518 addElemInMeshGroup(theHelper->GetMesh(), aGMFNode, aNodeGroupByHybridId.at(aGMFID-1), groupsToRemove);
522 else if (token == GmfCorners && nbElem > 0) {
523 (nbElem <= 1) ? tmpStr = " corner" : tmpStr = " corners";
524 for ( int iElem = 0; iElem < nbElem; iElem++ )
525 MGOutput->GmfGetLin(InpMsh, token, &id[iElem*tabRef[token]]);
527 else if (token == GmfRidges && nbElem > 0) {
528 (nbElem <= 1) ? tmpStr = " ridge" : tmpStr = " ridges";
529 for ( int iElem = 0; iElem < nbElem; iElem++ )
530 MGOutput->GmfGetLin(InpMsh, token, &id[iElem*tabRef[token]]);
532 else if (token == GmfEdges && nbElem > 0) {
533 (nbElem <= 1) ? tmpStr = " edge" : tmpStr = " edges";
534 for ( int iElem = 0; iElem < nbElem; iElem++ )
535 MGOutput->GmfGetLin(InpMsh, token, &id[iElem*tabRef[token]], &id[iElem*tabRef[token]+1], &domainID[iElem]);
537 else if (token == GmfTriangles && nbElem > 0) {
538 (nbElem <= 1) ? tmpStr = " triangle" : tmpStr = " triangles";
539 for ( int iElem = 0; iElem < nbElem; iElem++ )
540 MGOutput->GmfGetLin(InpMsh, token, &id[iElem*tabRef[token]], &id[iElem*tabRef[token]+1], &id[iElem*tabRef[token]+2], &domainID[iElem]);
542 else if (token == GmfQuadrilaterals && nbElem > 0) {
543 (nbElem <= 1) ? tmpStr = " Quadrilateral" : tmpStr = " Quadrilaterals";
544 for ( int iElem = 0; iElem < nbElem; iElem++ )
545 MGOutput->GmfGetLin(InpMsh, token, &id[iElem*tabRef[token]], &id[iElem*tabRef[token]+1], &id[iElem*tabRef[token]+2], &id[iElem*tabRef[token]+3], &domainID[iElem]);
547 else if (token == GmfTetrahedra && nbElem > 0) {
548 (nbElem <= 1) ? tmpStr = " Tetrahedron" : tmpStr = " Tetrahedra";
549 for ( int iElem = 0; iElem < nbElem; iElem++ ) {
550 MGOutput->GmfGetLin(InpMsh, token, &id[iElem*tabRef[token]], &id[iElem*tabRef[token]+1], &id[iElem*tabRef[token]+2], &id[iElem*tabRef[token]+3], &domainID[iElem]);
552 subdomainId2tetraId[dummy].insert(iElem+1);
556 else if (token == GmfPyramids && nbElem > 0) {
557 (nbElem <= 1) ? tmpStr = " Pyramid" : tmpStr = " Pyramids";
558 for ( int iElem = 0; iElem < nbElem; iElem++ )
559 MGOutput->GmfGetLin(InpMsh, token, &id[iElem*tabRef[token]], &id[iElem*tabRef[token]+1], &id[iElem*tabRef[token]+2], &id[iElem*tabRef[token]+3],
560 &id[iElem*tabRef[token]+4], &domainID[iElem]);
562 else if (token == GmfPrisms && nbElem > 0) {
563 (nbElem <= 1) ? tmpStr = " Prism" : tmpStr = " Prisms";
564 for ( int iElem = 0; iElem < nbElem; iElem++ )
565 MGOutput->GmfGetLin(InpMsh, token, &id[iElem*tabRef[token]], &id[iElem*tabRef[token]+1], &id[iElem*tabRef[token]+2], &id[iElem*tabRef[token]+3],
566 &id[iElem*tabRef[token]+4], &id[iElem*tabRef[token]+5], &domainID[iElem]);
568 else if (token == GmfHexahedra && nbElem > 0) {
569 (nbElem <= 1) ? tmpStr = " Hexahedron" : tmpStr = " Hexahedra";
570 for ( int iElem = 0; iElem < nbElem; iElem++ )
571 MGOutput->GmfGetLin(InpMsh, token, &id[iElem*tabRef[token]], &id[iElem*tabRef[token]+1], &id[iElem*tabRef[token]+2], &id[iElem*tabRef[token]+3],
572 &id[iElem*tabRef[token]+4], &id[iElem*tabRef[token]+5], &id[iElem*tabRef[token]+6], &id[iElem*tabRef[token]+7], &domainID[iElem]);
574 std::cout << tmpStr << std::endl;
575 //std::cout << std::endl;
582 case GmfQuadrilaterals:
588 std::vector< const SMDS_MeshNode* > node( nbRef );
589 std::vector< int > nodeID( nbRef );
590 std::vector< SMDS_MeshNode* > enfNode( nbRef );
591 const SMDS_MeshElement* aCreatedElem;
593 for ( int iElem = 0; iElem < nbElem; iElem++ )
595 if(theAlgo->computeCanceled()) {
596 MGOutput->GmfCloseMesh(InpMsh);
599 // Check if elem is already in input mesh. If yes => skip
600 bool fullyCreatedElement = false; // if at least one of the nodes was created
601 for ( int iRef = 0; iRef < nbRef; iRef++ )
603 aGMFNodeID = id[iElem*tabRef[token]+iRef]; // read nbRef aGMFNodeID
604 if (aGMFNodeID <= nbInitialNodes) // input nodes
607 node[ iRef ] = theNodeByHybridId[aGMFNodeID];
611 fullyCreatedElement = true;
612 aGMFNodeID -= nbInitialNodes;
613 nodeID[ iRef ] = aGMFNodeID ;
614 node [ iRef ] = GMFNode[ aGMFNodeID ];
621 if (fullyCreatedElement) {
622 aCreatedElem = theHelper->AddEdge( node[0], node[1], noID, force3d );
623 if ( !anEdgeGroupByHybridId.empty() && !anEdgeGroupByHybridId[iElem].empty())
624 addElemInMeshGroup(theHelper->GetMesh(), aCreatedElem, anEdgeGroupByHybridId[iElem], groupsToRemove);
628 if (fullyCreatedElement) {
629 aCreatedElem = theHelper->AddFace( node[0], node[1], node[2], noID, force3d );
630 // add iElem < aFaceGroupByHybridId.size() to avoid crash if imprinting with hexa core with MeshGems <= 2.4-5
631 if ( iElem < (int)aFaceGroupByHybridId.size() && !aFaceGroupByHybridId[iElem].empty() ) {
632 addElemInMeshGroup(theHelper->GetMesh(), aCreatedElem, aFaceGroupByHybridId[iElem], groupsToRemove);
634 // add element in shape for groups on geom to work
635 if ( domainID[iElem] > 0 )
637 theMeshDS->SetMeshElementOnShape( aCreatedElem, domainID[iElem] );
638 for ( int iN = 0; iN < 3; ++iN )
639 if ( node[iN]->getshapeId() < 1 )
640 theMeshDS->SetNodeOnFace( node[iN], domainID[iElem] );
644 case GmfQuadrilaterals:
645 if (fullyCreatedElement) {
646 aCreatedElem = theHelper->AddFace( node[0], node[1], node[2], node[3], noID, force3d );
647 // add element in shape for groups on geom to work
648 if ( domainID[iElem] > 0 )
650 theMeshDS->SetMeshElementOnShape( aCreatedElem, domainID[iElem] );
651 for ( int iN = 0; iN < 3; ++iN )
652 if ( node[iN]->getshapeId() < 1 )
653 theMeshDS->SetNodeOnFace( node[iN], domainID[iElem] );
660 if ( solidID != HOLE_ID )
662 aCreatedElem = theHelper->AddVolume( node[1], node[0], node[2], node[3],
664 theMeshDS->SetMeshElementOnShape( aCreatedElem, solidID );
665 for ( int iN = 0; iN < 4; ++iN )
666 if ( node[iN]->getshapeId() < 1 )
667 theMeshDS->SetNodeInVolume( node[iN], solidID );
672 if ( elemSearcher ) {
673 // Issue 0020682. Avoid creating nodes and tetras at place where
674 // volumic elements already exist
675 if ( !node[1] || !node[0] || !node[2] || !node[3] )
677 if ( elemSearcher->FindElementsByPoint((SMESH_TNodeXYZ(node[0]) +
678 SMESH_TNodeXYZ(node[1]) +
679 SMESH_TNodeXYZ(node[2]) +
680 SMESH_TNodeXYZ(node[3]) ) / 4.,
681 SMDSAbs_Volume, foundVolumes ))
684 aCreatedElem = theHelper->AddVolume( node[1], node[0], node[2], node[3],
691 if ( solidID != HOLE_ID )
693 aCreatedElem = theHelper->AddVolume( node[3], node[2], node[1],
696 theMeshDS->SetMeshElementOnShape( aCreatedElem, solidID );
697 for ( int iN = 0; iN < 5; ++iN )
698 if ( node[iN]->getshapeId() < 1 )
699 theMeshDS->SetNodeInVolume( node[iN], solidID );
704 if ( elemSearcher ) {
705 // Issue 0020682. Avoid creating nodes and tetras at place where
706 // volumic elements already exist
707 if ( !node[1] || !node[0] || !node[2] || !node[3] || !node[4] || !node[5] )
709 if ( elemSearcher->FindElementsByPoint((SMESH_TNodeXYZ(node[0]) +
710 SMESH_TNodeXYZ(node[1]) +
711 SMESH_TNodeXYZ(node[2]) +
712 SMESH_TNodeXYZ(node[3]) +
713 SMESH_TNodeXYZ(node[4])) / 5.,
714 SMDSAbs_Volume, foundVolumes ))
717 aCreatedElem = theHelper->AddVolume( node[3], node[2], node[1],
725 if ( solidID != HOLE_ID )
727 aCreatedElem = theHelper->AddVolume( node[0], node[2], node[1],
728 node[3], node[5], node[4],
730 theMeshDS->SetMeshElementOnShape( aCreatedElem, solidID );
731 for ( int iN = 0; iN < 6; ++iN )
732 if ( node[iN]->getshapeId() < 1 )
733 theMeshDS->SetNodeInVolume( node[iN], solidID );
738 if ( elemSearcher ) {
739 // Issue 0020682. Avoid creating nodes and tetras at place where
740 // volumic elements already exist
741 if ( !node[1] || !node[0] || !node[2] || !node[3] || !node[4] || !node[5] )
743 if ( elemSearcher->FindElementsByPoint((SMESH_TNodeXYZ(node[0]) +
744 SMESH_TNodeXYZ(node[1]) +
745 SMESH_TNodeXYZ(node[2]) +
746 SMESH_TNodeXYZ(node[3]) +
747 SMESH_TNodeXYZ(node[4]) +
748 SMESH_TNodeXYZ(node[5])) / 6.,
749 SMDSAbs_Volume, foundVolumes ))
752 aCreatedElem = theHelper->AddVolume( node[0], node[2], node[1],
753 node[3], node[5], node[4],
760 if ( solidID != HOLE_ID )
762 aCreatedElem = theHelper->AddVolume( node[0], node[3], node[2], node[1],
763 node[4], node[7], node[6], node[5],
765 theMeshDS->SetMeshElementOnShape( aCreatedElem, solidID );
766 for ( int iN = 0; iN < 8; ++iN )
767 if ( node[iN]->getshapeId() < 1 )
768 theMeshDS->SetNodeInVolume( node[iN], solidID );
773 if ( elemSearcher ) {
774 // Issue 0020682. Avoid creating nodes and tetras at place where
775 // volumic elements already exist
776 if ( !node[1] || !node[0] || !node[2] || !node[3] || !node[4] || !node[5] || !node[6] || !node[7])
778 if ( elemSearcher->FindElementsByPoint((SMESH_TNodeXYZ(node[0]) +
779 SMESH_TNodeXYZ(node[1]) +
780 SMESH_TNodeXYZ(node[2]) +
781 SMESH_TNodeXYZ(node[3]) +
782 SMESH_TNodeXYZ(node[4]) +
783 SMESH_TNodeXYZ(node[5]) +
784 SMESH_TNodeXYZ(node[6]) +
785 SMESH_TNodeXYZ(node[7])) / 8.,
786 SMDSAbs_Volume, foundVolumes ))
789 aCreatedElem = theHelper->AddVolume( node[0], node[3], node[2], node[1],
790 node[4], node[7], node[6], node[5],
797 if ( aCreatedElem && toMakeGroupsOfDomains )
799 if ( domainID[iElem] >= (int) elemsOfDomain.size() )
800 elemsOfDomain.resize( domainID[iElem] + 1 );
801 elemsOfDomain[ domainID[iElem] ].push_back( aCreatedElem );
803 } // loop on elements of one type
811 // remove nodes in holes
814 for ( int i = 1; i <= nbVertices; ++i )
815 if ( GMFNode[i]->NbInverseElements() == 0 )
816 theMeshDS->RemoveFreeNode( GMFNode[i], /*sm=*/0, /*fromGroups=*/false );
819 MGOutput->GmfCloseMesh(InpMsh);
821 // 0022172: [CEA 790] create the groups corresponding to domains
822 if ( toMakeGroupsOfDomains )
823 makeDomainGroups( elemsOfDomain, theHelper );
826 std::map<int, std::set<int> >::const_iterator subdomainIt = subdomainId2tetraId.begin();
827 std::string aSubdomainFileName = theFile;
828 aSubdomainFileName = aSubdomainFileName + ".subdomain";
829 std::ofstream aSubdomainFile ( aSubdomainFileName , ios::out);
831 aSubdomainFile << "Nb subdomains " << subdomainId2tetraId.size() << std::endl;
832 for(;subdomainIt != subdomainId2tetraId.end() ; ++subdomainIt) {
833 int subdomainId = subdomainIt->first;
834 std::set<int> tetraIds = subdomainIt->second;
835 std::set<int>::const_iterator tetraIdsIt = tetraIds.begin();
836 aSubdomainFile << subdomainId << std::endl;
837 for(;tetraIdsIt != tetraIds.end() ; ++tetraIdsIt) {
838 aSubdomainFile << (*tetraIdsIt) << " ";
840 aSubdomainFile << std::endl;
842 aSubdomainFile.close();
849 static bool writeGMFFile(MG_HYBRID_API* MGInput,
850 const char* theMeshFileName,
851 const char* theRequiredFileName,
852 const char* theSolFileName,
853 const SMESH_ProxyMesh& theProxyMesh,
854 SMESH_MesherHelper& theHelper,
855 std::vector <const SMDS_MeshNode*> & theNodeByHybridId,
856 std::vector <const SMDS_MeshElement*> & theFaceByHybridId,
857 std::map<const SMDS_MeshNode*,int> & aNodeToHybridIdMap,
858 std::vector<std::string> & aNodeGroupByHybridId,
859 std::vector<std::string> & anEdgeGroupByHybridId,
860 std::vector<std::string> & aFaceGroupByHybridId,
861 HYBRIDPlugin_Hypothesis::TIDSortedNodeGroupMap & theEnforcedNodes,
862 HYBRIDPlugin_Hypothesis::TIDSortedElemGroupMap & theEnforcedEdges,
863 HYBRIDPlugin_Hypothesis::TIDSortedElemGroupMap & theEnforcedTriangles,
864 std::map<std::vector<double>, std::string> & enfVerticesWithGroup,
865 HYBRIDPlugin_Hypothesis::THYBRIDEnforcedVertexCoordsValues & theEnforcedVertices)
868 int idx, idxRequired = 0, idxSol = 0;
870 //const int dummyint = 0;
871 const int dummyint1 = 1;
872 const int dummyint2 = 2;
873 const int dummyint3 = 3;
874 const int dummyint4 = 4;
875 const int enforcedTag = HYBRIDPlugin_Hypothesis::EnforcedTag();
876 //const int dummyint6 = 6; //are interesting for layers
877 HYBRIDPlugin_Hypothesis::THYBRIDEnforcedVertexCoordsValues::const_iterator vertexIt;
878 std::vector<double> enfVertexSizes;
879 const SMDS_MeshElement* elem;
880 TIDSortedElemSet anElemSetTri, anElemSetQuad, theKeptEnforcedEdges, theKeptEnforcedTriangles;
881 SMDS_ElemIteratorPtr nodeIt;
882 std::vector <const SMDS_MeshNode*> theEnforcedNodeByHybridId;
883 std::map<const SMDS_MeshNode*,int> anEnforcedNodeToHybridIdMap, anExistingEnforcedNodeToHybridIdMap;
884 std::vector< const SMDS_MeshElement* > foundElems;
885 std::map<const SMDS_MeshNode*,TopAbs_State> aNodeToTopAbs_StateMap;
887 HYBRIDPlugin_Hypothesis::TIDSortedElemGroupMap::iterator elemIt;
888 TIDSortedElemSet::iterator elemSetIt;
890 SMESH_Mesh* theMesh = theHelper.GetMesh();
891 const bool hasGeom = theMesh->HasShapeToMesh();
892 SMESHUtils::Deleter< SMESH_ElementSearcher > pntCls
893 ( SMESH_MeshAlgos::GetElementSearcher(*theMesh->GetMeshDS()));
895 int nbEnforcedVertices = theEnforcedVertices.size();
898 int nbFaces = theProxyMesh.NbFaces();
900 theFaceByHybridId.reserve( nbFaces );
903 int usedEnforcedNodes = 0;
909 idx = MGInput->GmfOpenMesh(theMeshFileName, GmfWrite, GMFVERSION, GMFDIMENSION);
913 // ========================== FACES ==========================
914 // TRIANGLES ==========================
915 SMDS_ElemIteratorPtr eIt =
916 hasGeom ? theProxyMesh.GetFaces( theHelper.GetSubShape()) : theProxyMesh.GetFaces();
917 while ( eIt->more() )
920 nodeIt = elem->nodesIterator();
921 nbNodes = elem->NbCornerNodes();
923 anElemSetTri.insert(elem);
924 else if (nbNodes == 4)
925 anElemSetQuad.insert(elem);
928 std::cout << "Unexpected number of nodes: " << nbNodes << std::endl;
929 throw ("Unexpected number of nodes" );
931 while ( nodeIt->more() && nbNodes--)
934 const SMDS_MeshNode* node = castToNode( nodeIt->next() );
935 int newId = aNodeToHybridIdMap.size() + 1; // hybrid ids count from 1
936 aNodeToHybridIdMap.insert( std::make_pair( node, newId ));
940 //EDGES ==========================
942 // Iterate over the enforced edges
943 for(elemIt = theEnforcedEdges.begin() ; elemIt != theEnforcedEdges.end() ; ++elemIt) {
944 elem = elemIt->first;
946 nodeIt = elem->nodesIterator();
948 while ( nodeIt->more() && nbNodes-- ) {
950 const SMDS_MeshNode* node = castToNode( nodeIt->next() );
951 // Test if point is inside shape to mesh
952 gp_Pnt myPoint(node->X(),node->Y(),node->Z());
953 TopAbs_State result = pntCls->GetPointState( myPoint );
954 if ( result == TopAbs_OUT ) {
958 aNodeToTopAbs_StateMap.insert( std::make_pair( node, result ));
961 nodeIt = elem->nodesIterator();
964 while ( nodeIt->more() && nbNodes-- ) {
966 const SMDS_MeshNode* node = castToNode( nodeIt->next() );
967 gp_Pnt myPoint(node->X(),node->Y(),node->Z());
968 nbFoundElems = pntCls->FindElementsByPoint(myPoint, SMDSAbs_Node, foundElems);
970 std::cout << "Node at "<<node->X()<<", "<<node->Y()<<", "<<node->Z()<<std::endl;
971 std::cout << "Nb nodes found : "<<nbFoundElems<<std::endl;
973 if (nbFoundElems ==0) {
974 if ((*aNodeToTopAbs_StateMap.find(node)).second == TopAbs_IN) {
975 newId = aNodeToHybridIdMap.size() + anEnforcedNodeToHybridIdMap.size() + 1; // hybrid ids count from 1
976 anEnforcedNodeToHybridIdMap.insert( std::make_pair( node, newId ));
979 else if (nbFoundElems ==1) {
980 const SMDS_MeshNode* existingNode = (SMDS_MeshNode*) foundElems.at(0);
981 newId = (*aNodeToHybridIdMap.find(existingNode)).second;
982 anExistingEnforcedNodeToHybridIdMap.insert( std::make_pair( node, newId ));
987 std::cout << "HYBRID node ID: "<<newId<<std::endl;
991 theKeptEnforcedEdges.insert(elem);
995 //ENFORCED TRIANGLES ==========================
997 // Iterate over the enforced triangles
998 for(elemIt = theEnforcedTriangles.begin() ; elemIt != theEnforcedTriangles.end() ; ++elemIt) {
999 elem = elemIt->first;
1001 nodeIt = elem->nodesIterator();
1003 while ( nodeIt->more() && nbNodes--) {
1005 const SMDS_MeshNode* node = castToNode( nodeIt->next() );
1006 // Test if point is inside shape to mesh
1007 gp_Pnt myPoint(node->X(),node->Y(),node->Z());
1008 TopAbs_State result = pntCls->GetPointState( myPoint );
1009 if ( result == TopAbs_OUT ) {
1013 aNodeToTopAbs_StateMap.insert( std::make_pair( node, result ));
1016 nodeIt = elem->nodesIterator();
1019 while ( nodeIt->more() && nbNodes--) {
1021 const SMDS_MeshNode* node = castToNode( nodeIt->next() );
1022 gp_Pnt myPoint(node->X(),node->Y(),node->Z());
1023 nbFoundElems = pntCls->FindElementsByPoint(myPoint, SMDSAbs_Node, foundElems);
1025 std::cout << "Nb nodes found : "<<nbFoundElems<<std::endl;
1027 if (nbFoundElems ==0) {
1028 if ((*aNodeToTopAbs_StateMap.find(node)).second == TopAbs_IN) {
1029 newId = aNodeToHybridIdMap.size() + anEnforcedNodeToHybridIdMap.size() + 1; // hybrid ids count from 1
1030 anEnforcedNodeToHybridIdMap.insert( std::make_pair( node, newId ));
1033 else if (nbFoundElems ==1) {
1034 const SMDS_MeshNode* existingNode = (SMDS_MeshNode*) foundElems.at(0);
1035 newId = (*aNodeToHybridIdMap.find(existingNode)).second;
1036 anExistingEnforcedNodeToHybridIdMap.insert( std::make_pair( node, newId ));
1041 std::cout << "HYBRID node ID: "<<newId<<std::endl;
1045 theKeptEnforcedTriangles.insert(elem);
1049 // put nodes to theNodeByHybridId vector
1051 std::cout << "aNodeToHybridIdMap.size(): "<<aNodeToHybridIdMap.size()<<std::endl;
1053 theNodeByHybridId.resize( aNodeToHybridIdMap.size() );
1054 std::map<const SMDS_MeshNode*,int>::const_iterator n2id = aNodeToHybridIdMap.begin();
1055 for ( ; n2id != aNodeToHybridIdMap.end(); ++ n2id)
1057 // std::cout << "n2id->first: "<<n2id->first<<std::endl;
1058 theNodeByHybridId[ n2id->second - 1 ] = n2id->first; // hybrid ids count from 1
1061 // put nodes to anEnforcedNodeToHybridIdMap vector
1063 std::cout << "anEnforcedNodeToHybridIdMap.size(): "<<anEnforcedNodeToHybridIdMap.size()<<std::endl;
1065 theEnforcedNodeByHybridId.resize( anEnforcedNodeToHybridIdMap.size());
1066 n2id = anEnforcedNodeToHybridIdMap.begin();
1067 for ( ; n2id != anEnforcedNodeToHybridIdMap.end(); ++ n2id)
1069 if (n2id->second > (int)aNodeToHybridIdMap.size()) {
1070 theEnforcedNodeByHybridId[ n2id->second - aNodeToHybridIdMap.size() - 1 ] = n2id->first; // hybrid ids count from 1
1075 //========================== NODES ==========================
1076 std::vector<const SMDS_MeshNode*> theOrderedNodes, theRequiredNodes;
1077 std::set< std::vector<double> > nodesCoords;
1078 std::vector<const SMDS_MeshNode*>::const_iterator hybridNodeIt = theNodeByHybridId.begin();
1079 std::vector<const SMDS_MeshNode*>::const_iterator after = theNodeByHybridId.end();
1081 (theNodeByHybridId.size() <= 1) ? tmpStr = " node" : " nodes";
1082 std::cout << theNodeByHybridId.size() << tmpStr << " from mesh ..." << std::endl;
1083 for ( ; hybridNodeIt != after; ++hybridNodeIt )
1085 const SMDS_MeshNode* node = *hybridNodeIt;
1086 std::vector<double> coords;
1087 coords.push_back(node->X());
1088 coords.push_back(node->Y());
1089 coords.push_back(node->Z());
1090 nodesCoords.insert(coords);
1091 theOrderedNodes.push_back(node);
1094 // Iterate over the enforced nodes given by enforced elements
1095 hybridNodeIt = theEnforcedNodeByHybridId.begin();
1096 after = theEnforcedNodeByHybridId.end();
1097 (theEnforcedNodeByHybridId.size() <= 1) ? tmpStr = " node" : " nodes";
1098 std::cout << theEnforcedNodeByHybridId.size() << tmpStr << " from enforced elements ..." << std::endl;
1099 for ( ; hybridNodeIt != after; ++hybridNodeIt )
1101 const SMDS_MeshNode* node = *hybridNodeIt;
1102 std::vector<double> coords;
1103 coords.push_back(node->X());
1104 coords.push_back(node->Y());
1105 coords.push_back(node->Z());
1107 std::cout << "Node at " << node->X()<<", " <<node->Y()<<", " <<node->Z();
1110 if (nodesCoords.find(coords) != nodesCoords.end()) {
1111 // node already exists in original mesh
1113 std::cout << " found" << std::endl;
1118 if (theEnforcedVertices.find(coords) != theEnforcedVertices.end()) {
1119 // node already exists in enforced vertices
1121 std::cout << " found" << std::endl;
1127 std::cout << " not found" << std::endl;
1130 nodesCoords.insert(coords);
1131 theOrderedNodes.push_back(node);
1132 // theRequiredNodes.push_back(node);
1136 // Iterate over the enforced nodes
1137 HYBRIDPlugin_Hypothesis::TIDSortedNodeGroupMap::const_iterator enfNodeIt;
1138 (theEnforcedNodes.size() <= 1) ? tmpStr = " node" : " nodes";
1139 std::cout << theEnforcedNodes.size() << tmpStr << " from enforced nodes ..." << std::endl;
1140 for(enfNodeIt = theEnforcedNodes.begin() ; enfNodeIt != theEnforcedNodes.end() ; ++enfNodeIt)
1142 const SMDS_MeshNode* node = enfNodeIt->first;
1143 std::vector<double> coords;
1144 coords.push_back(node->X());
1145 coords.push_back(node->Y());
1146 coords.push_back(node->Z());
1148 std::cout << "Node at " << node->X()<<", " <<node->Y()<<", " <<node->Z();
1151 // Test if point is inside shape to mesh
1152 gp_Pnt myPoint(node->X(),node->Y(),node->Z());
1153 TopAbs_State result = pntCls->GetPointState( myPoint );
1154 if ( result == TopAbs_OUT ) {
1156 std::cout << " out of volume" << std::endl;
1161 if (nodesCoords.find(coords) != nodesCoords.end()) {
1163 std::cout << " found in nodesCoords" << std::endl;
1165 // theRequiredNodes.push_back(node);
1169 if (theEnforcedVertices.find(coords) != theEnforcedVertices.end()) {
1171 std::cout << " found in theEnforcedVertices" << std::endl;
1177 std::cout << " not found" << std::endl;
1179 nodesCoords.insert(coords);
1180 // theOrderedNodes.push_back(node);
1181 theRequiredNodes.push_back(node);
1183 int requiredNodes = theRequiredNodes.size();
1186 std::vector<std::vector<double> > ReqVerTab;
1187 if (nbEnforcedVertices) {
1188 (nbEnforcedVertices <= 1) ? tmpStr = " node" : " nodes";
1189 std::cout << nbEnforcedVertices << tmpStr << " from enforced vertices ..." << std::endl;
1190 // Iterate over the enforced vertices
1191 for(vertexIt = theEnforcedVertices.begin() ; vertexIt != theEnforcedVertices.end() ; ++vertexIt) {
1192 double x = vertexIt->first[0];
1193 double y = vertexIt->first[1];
1194 double z = vertexIt->first[2];
1195 // Test if point is inside shape to mesh
1196 gp_Pnt myPoint(x,y,z);
1197 TopAbs_State result = pntCls->GetPointState( myPoint );
1198 if ( result == TopAbs_OUT )
1200 std::vector<double> coords;
1201 coords.push_back(x);
1202 coords.push_back(y);
1203 coords.push_back(z);
1204 ReqVerTab.push_back(coords);
1205 enfVertexSizes.push_back(vertexIt->second);
1212 std::cout << "Begin writing required nodes in GmfVertices" << std::endl;
1213 std::cout << "Nb vertices: " << theOrderedNodes.size() << std::endl;
1214 MGInput->GmfSetKwd(idx, GmfVertices, theOrderedNodes.size());
1215 for (hybridNodeIt = theOrderedNodes.begin();hybridNodeIt != theOrderedNodes.end();++hybridNodeIt) {
1216 MGInput->GmfSetLin(idx, GmfVertices, (*hybridNodeIt)->X(), (*hybridNodeIt)->Y(), (*hybridNodeIt)->Z(), dummyint1);
1219 std::cout << "End writing required nodes in GmfVertices" << std::endl;
1221 if (requiredNodes + solSize) {
1222 std::cout << "Begin writing in req and sol file" << std::endl;
1223 aNodeGroupByHybridId.resize( requiredNodes + solSize );
1224 idxRequired = MGInput->GmfOpenMesh(theRequiredFileName, GmfWrite, GMFVERSION, GMFDIMENSION);
1226 MGInput->GmfCloseMesh(idx);
1229 idxSol = MGInput->GmfOpenMesh(theSolFileName, GmfWrite, GMFVERSION, GMFDIMENSION);
1231 MGInput->GmfCloseMesh(idx);
1233 MGInput->GmfCloseMesh(idxRequired);
1236 int TypTab[] = {GmfSca};
1237 double ValTab[] = {0.0};
1238 MGInput->GmfSetKwd(idxRequired, GmfVertices, requiredNodes + solSize);
1239 MGInput->GmfSetKwd(idxSol, GmfSolAtVertices, requiredNodes + solSize, 1, TypTab);
1240 for (hybridNodeIt = theRequiredNodes.begin();hybridNodeIt != theRequiredNodes.end();++hybridNodeIt) {
1241 MGInput->GmfSetLin(idxRequired, GmfVertices, (*hybridNodeIt)->X(), (*hybridNodeIt)->Y(), (*hybridNodeIt)->Z(), dummyint2);
1242 MGInput->GmfSetLin(idxSol, GmfSolAtVertices, ValTab);
1243 if (theEnforcedNodes.find((*hybridNodeIt)) != theEnforcedNodes.end())
1244 gn = theEnforcedNodes.find((*hybridNodeIt))->second;
1245 aNodeGroupByHybridId[usedEnforcedNodes] = gn;
1246 usedEnforcedNodes++;
1249 for (int i=0;i<solSize;i++) {
1250 std::cout << ReqVerTab[i][0] <<" "<< ReqVerTab[i][1] << " "<< ReqVerTab[i][2] << std::endl;
1252 std::cout << "enfVertexSizes.at("<<i<<"): " << enfVertexSizes.at(i) << std::endl;
1254 double solTab[] = {enfVertexSizes.at(i)};
1255 MGInput->GmfSetLin(idxRequired, GmfVertices, ReqVerTab[i][0], ReqVerTab[i][1], ReqVerTab[i][2], dummyint3);
1256 MGInput->GmfSetLin(idxSol, GmfSolAtVertices, solTab);
1257 aNodeGroupByHybridId[usedEnforcedNodes] = enfVerticesWithGroup.find(ReqVerTab[i])->second;
1259 std::cout << "aNodeGroupByHybridId["<<usedEnforcedNodes<<"] = \""<<aNodeGroupByHybridId[usedEnforcedNodes]<<"\""<<std::endl;
1261 usedEnforcedNodes++;
1263 std::cout << "End writing in req and sol file" << std::endl;
1266 int nedge[2], ntri[3], nquad[4];
1269 int usedEnforcedEdges = 0;
1270 if (theKeptEnforcedEdges.size()) {
1271 anEdgeGroupByHybridId.resize( theKeptEnforcedEdges.size() );
1272 MGInput->GmfSetKwd(idx, GmfEdges, theKeptEnforcedEdges.size());
1273 for(elemSetIt = theKeptEnforcedEdges.begin() ; elemSetIt != theKeptEnforcedEdges.end() ; ++elemSetIt) {
1274 elem = (*elemSetIt);
1275 nodeIt = elem->nodesIterator();
1277 while ( nodeIt->more() ) {
1279 const SMDS_MeshNode* node = castToNode( nodeIt->next() );
1280 std::map< const SMDS_MeshNode*,int >::iterator it = anEnforcedNodeToHybridIdMap.find(node);
1281 if (it == anEnforcedNodeToHybridIdMap.end()) {
1282 it = anExistingEnforcedNodeToHybridIdMap.find(node);
1283 if (it == anEnforcedNodeToHybridIdMap.end())
1284 throw "Node not found";
1286 nedge[index] = it->second;
1289 MGInput->GmfSetLin(idx, GmfEdges, nedge[0], nedge[1], dummyint4);
1290 anEdgeGroupByHybridId[usedEnforcedEdges] = theEnforcedEdges.find(elem)->second;
1291 usedEnforcedEdges++;
1296 if (usedEnforcedEdges) {
1297 MGInput->GmfSetKwd(idx, GmfRequiredEdges, usedEnforcedEdges);
1298 for (int enfID=1;enfID<=usedEnforcedEdges;enfID++) {
1299 MGInput->GmfSetLin(idx, GmfRequiredEdges, enfID);
1304 int usedEnforcedTriangles = 0;
1305 if (anElemSetTri.size()+theKeptEnforcedTriangles.size())
1307 aFaceGroupByHybridId.resize( anElemSetTri.size()+theKeptEnforcedTriangles.size() );
1308 MGInput->GmfSetKwd(idx, GmfTriangles, anElemSetTri.size()+theKeptEnforcedTriangles.size());
1310 for(elemSetIt = anElemSetTri.begin() ; elemSetIt != anElemSetTri.end() ; ++elemSetIt,++k)
1312 elem = (*elemSetIt);
1313 theFaceByHybridId.push_back( elem );
1314 nodeIt = elem->nodesIterator();
1316 for ( int j = 0; j < 3; ++j )
1319 const SMDS_MeshNode* node = castToNode( nodeIt->next() );
1320 std::map< const SMDS_MeshNode*,int >::iterator it = aNodeToHybridIdMap.find(node);
1321 if (it == aNodeToHybridIdMap.end())
1322 throw "Node not found";
1323 ntri[index] = it->second;
1326 MGInput->GmfSetLin(idx, GmfTriangles, ntri[0], ntri[1], ntri[2], /*tag=*/elem->getshapeId() );
1327 aFaceGroupByHybridId[k] = "";
1330 if ( !theHelper.GetMesh()->HasShapeToMesh() ) SMESHUtils::FreeVector( theFaceByHybridId );
1331 std::cout << "Enforced triangles size " << theKeptEnforcedTriangles.size() << std::endl;
1332 if (theKeptEnforcedTriangles.size())
1334 for(elemSetIt = theKeptEnforcedTriangles.begin() ; elemSetIt != theKeptEnforcedTriangles.end() ; ++elemSetIt,++k)
1336 elem = (*elemSetIt);
1337 nodeIt = elem->nodesIterator();
1339 for ( int j = 0; j < 3; ++j )
1342 const SMDS_MeshNode* node = castToNode( nodeIt->next() );
1343 std::map< const SMDS_MeshNode*,int >::iterator it = anEnforcedNodeToHybridIdMap.find(node);
1344 if (it == anEnforcedNodeToHybridIdMap.end())
1346 it = anExistingEnforcedNodeToHybridIdMap.find(node);
1347 if (it == anEnforcedNodeToHybridIdMap.end())
1348 throw "Node not found";
1350 ntri[index] = it->second;
1353 MGInput->GmfSetLin(idx, GmfTriangles, ntri[0], ntri[1], ntri[2], enforcedTag);
1354 aFaceGroupByHybridId[k] = theEnforcedTriangles.find(elem)->second;
1355 usedEnforcedTriangles++;
1361 if (usedEnforcedTriangles)
1363 MGInput->GmfSetKwd(idx, GmfRequiredTriangles, usedEnforcedTriangles);
1364 for (int enfID=1;enfID<=usedEnforcedTriangles;enfID++)
1365 MGInput->GmfSetLin(idx, GmfRequiredTriangles, anElemSetTri.size()+enfID);
1368 if (anElemSetQuad.size())
1370 MGInput->GmfSetKwd(idx, GmfQuadrilaterals, anElemSetQuad.size());
1372 for(elemSetIt = anElemSetQuad.begin() ; elemSetIt != anElemSetQuad.end() ; ++elemSetIt,++k)
1374 elem = (*elemSetIt);
1375 theFaceByHybridId.push_back( elem );
1376 nodeIt = elem->nodesIterator();
1378 for ( int j = 0; j < 4; ++j )
1381 const SMDS_MeshNode* node = castToNode( nodeIt->next() );
1382 std::map< const SMDS_MeshNode*,int >::iterator it = aNodeToHybridIdMap.find(node);
1383 if (it == aNodeToHybridIdMap.end())
1384 throw "Node not found";
1385 nquad[index] = it->second;
1388 MGInput->GmfSetLin(idx, GmfQuadrilaterals, nquad[0], nquad[1], nquad[2], nquad[3],
1389 /*tag=*/elem->getshapeId() );
1390 // _CEA_cbo what is it for???
1391 //aFaceGroupByHybridId[k] = "";
1395 MGInput->GmfCloseMesh(idx);
1396 chmodUserOnly(theMeshFileName);
1399 MGInput->GmfCloseMesh(idxRequired);
1400 chmodUserOnly(theRequiredFileName);
1404 MGInput->GmfCloseMesh(idxSol);
1405 chmodUserOnly(theSolFileName);
1411 //=============================================================================
1413 *Here we are going to use the HYBRID mesher with geometry
1415 //=============================================================================
1417 bool HYBRIDPlugin_HYBRID::Compute(SMESH_Mesh& theMesh,
1418 const TopoDS_Shape& theShape)
1422 // a unique working file name
1423 // to avoid access to the same files by eg different users
1424 _genericName = HYBRIDPlugin_Hypothesis::GetFileName(_hyp);
1425 std::string aGenericName = _genericName;
1426 std::string aGenericNameRequired = aGenericName + "_required";
1428 std::string aLogFileName = aGenericName + ".log"; // log
1429 std::string aResultFileName;
1431 std::string aGMFFileName, aRequiredVerticesFileName, aSolFileName, aResSolFileName;
1432 aGMFFileName = aGenericName + ".mesh"; // GMF mesh file
1433 aResultFileName = aGenericName + "Vol.mesh"; // GMF mesh file
1434 aResSolFileName = aGenericName + "Vol.sol"; // GMF mesh file
1435 aRequiredVerticesFileName = aGenericNameRequired + ".mesh"; // GMF required vertices mesh file
1436 aSolFileName = aGenericNameRequired + ".sol"; // GMF solution file
1438 std::map <int,int> aNodeId2NodeIndexMap, aSmdsToHybridIdMap, anEnforcedNodeIdToHybridIdMap;
1439 std::map <int, int> nodeID2nodeIndexMap;
1440 std::map<std::vector<double>, std::string> enfVerticesWithGroup;
1441 HYBRIDPlugin_Hypothesis::THYBRIDEnforcedVertexCoordsValues coordsSizeMap = HYBRIDPlugin_Hypothesis::GetEnforcedVerticesCoordsSize(_hyp);
1442 HYBRIDPlugin_Hypothesis::TIDSortedNodeGroupMap enforcedNodes = HYBRIDPlugin_Hypothesis::GetEnforcedNodes(_hyp);
1443 HYBRIDPlugin_Hypothesis::TIDSortedElemGroupMap enforcedEdges = HYBRIDPlugin_Hypothesis::GetEnforcedEdges(_hyp);
1444 HYBRIDPlugin_Hypothesis::TIDSortedElemGroupMap enforcedTriangles = HYBRIDPlugin_Hypothesis::GetEnforcedTriangles(_hyp);
1445 HYBRIDPlugin_Hypothesis::TID2SizeMap nodeIDToSizeMap = HYBRIDPlugin_Hypothesis::GetNodeIDToSizeMap(_hyp);
1447 HYBRIDPlugin_Hypothesis::THYBRIDEnforcedVertexList enfVertices = HYBRIDPlugin_Hypothesis::GetEnforcedVertices(_hyp);
1448 HYBRIDPlugin_Hypothesis::THYBRIDEnforcedVertexList::const_iterator enfVerIt = enfVertices.begin();
1449 std::vector<double> coords;
1451 for ( ; enfVerIt != enfVertices.end() ; ++enfVerIt)
1453 HYBRIDPlugin_Hypothesis::THYBRIDEnforcedVertex* enfVertex = (*enfVerIt);
1454 if (enfVertex->coords.size()) {
1455 coordsSizeMap.insert(std::make_pair(enfVertex->coords,enfVertex->size));
1456 enfVerticesWithGroup.insert(std::make_pair(enfVertex->coords,enfVertex->groupName));
1459 TopoDS_Shape GeomShape = entryToShape(enfVertex->geomEntry);
1460 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
1462 if (it.Value().ShapeType() == TopAbs_VERTEX){
1463 gp_Pnt aPnt = BRep_Tool::Pnt(TopoDS::Vertex(it.Value()));
1464 coords.push_back(aPnt.X());
1465 coords.push_back(aPnt.Y());
1466 coords.push_back(aPnt.Z());
1467 if (coordsSizeMap.find(coords) == coordsSizeMap.end()) {
1468 coordsSizeMap.insert(std::make_pair(coords,enfVertex->size));
1469 enfVerticesWithGroup.insert(std::make_pair(coords,enfVertex->groupName));
1475 int nbEnforcedVertices = coordsSizeMap.size();
1476 int nbEnforcedNodes = enforcedNodes.size();
1479 (nbEnforcedNodes <= 1) ? tmpStr = "node" : "nodes";
1480 std::cout << nbEnforcedNodes << " enforced " << tmpStr << " from hypo" << std::endl;
1481 (nbEnforcedVertices <= 1) ? tmpStr = "vertex" : "vertices";
1482 std::cout << nbEnforcedVertices << " enforced " << tmpStr << " from hypo" << std::endl;
1484 SMESH_MesherHelper helper( theMesh );
1485 helper.SetSubShape( theShape );
1487 std::vector <const SMDS_MeshNode*> aNodeByHybridId, anEnforcedNodeByHybridId;
1488 std::vector <const SMDS_MeshElement*> aFaceByHybridId;
1489 std::map<const SMDS_MeshNode*,int> aNodeToHybridIdMap;
1490 std::vector<std::string> aNodeGroupByHybridId, anEdgeGroupByHybridId, aFaceGroupByHybridId;
1492 SMESH_ProxyMesh::Ptr proxyMesh( new SMESH_ProxyMesh( theMesh ));
1494 MG_HYBRID_API mgHybrid( _computeCanceled, _progress );
1496 Ok = writeGMFFile(&mgHybrid,
1497 aGMFFileName.c_str(),
1498 aRequiredVerticesFileName.c_str(),
1499 aSolFileName.c_str(),
1501 aNodeByHybridId, aFaceByHybridId, aNodeToHybridIdMap,
1502 aNodeGroupByHybridId, anEdgeGroupByHybridId, aFaceGroupByHybridId,
1503 enforcedNodes, enforcedEdges, enforcedTriangles, /*enforcedQuadrangles,*/
1504 enfVerticesWithGroup, coordsSizeMap);
1506 // Write aSmdsToHybridIdMap to temp file
1507 std::string aSmdsToHybridIdMapFileName;
1508 aSmdsToHybridIdMapFileName = aGenericName + ".ids"; // ids relation
1509 std::ofstream aIdsFile ( aSmdsToHybridIdMapFileName , ios::out);
1510 Ok = aIdsFile.rdbuf()->is_open();
1512 INFOS( "Can't write into " << aSmdsToHybridIdMapFileName);
1513 return error(SMESH_Comment("Can't write into ") << aSmdsToHybridIdMapFileName);
1515 INFOS( "Writing ids relation into " << aSmdsToHybridIdMapFileName);
1516 aIdsFile << "Smds Hybrid" << std::endl;
1517 std::map <int,int>::const_iterator myit;
1518 for (myit=aSmdsToHybridIdMap.begin() ; myit != aSmdsToHybridIdMap.end() ; ++myit) {
1519 aIdsFile << myit->first << " " << myit->second << std::endl;
1523 chmodUserOnly(aSmdsToHybridIdMapFileName.c_str());
1526 if ( !_keepFiles ) {
1527 removeFile( aGMFFileName );
1528 removeFile( aRequiredVerticesFileName );
1529 removeFile( aSolFileName );
1530 removeFile( aSmdsToHybridIdMapFileName );
1532 return error(COMPERR_BAD_INPUT_MESH);
1534 removeFile( aResultFileName ); // needed for boundary recovery module usage
1536 // -----------------
1537 // run hybrid mesher
1538 // -----------------
1540 std::string cmd = HYBRIDPlugin_Hypothesis::CommandToRun( _hyp );
1542 if ( mgHybrid.IsExecutable() )
1544 cmd += " --in " + aGMFFileName;
1545 cmd += " --out " + aResultFileName;
1547 std::cout << std::endl;
1548 std::cout << "Hybrid execution with geometry..." << std::endl;
1550 if ( !_logInStandardOutput )
1552 mgHybrid.SetLogFile( aLogFileName );
1553 if ( mgHybrid.IsExecutable() )
1554 cmd += " 1>" + aLogFileName; // dump into file
1555 std::cout << " 1> " << aLogFileName;
1557 std::cout << std::endl;
1559 _computeCanceled = false;
1562 Ok = mgHybrid.Compute( cmd, errStr ); // run
1564 if ( _logInStandardOutput && mgHybrid.IsLibrary() )
1565 std::cout << std::endl << mgHybrid.GetLog() << std::endl;
1567 chmodUserOnly(aLogFileName.c_str());
1569 std::cout << "End of Hybrid execution !" << std::endl;
1571 if ( mgHybrid.IsExecutable() )
1573 chmodUserOnly(aResultFileName.c_str());
1574 chmodUserOnly(aResSolFileName.c_str());
1581 // Mapping the result file
1583 HYBRIDPlugin_Hypothesis::TSetStrings groupsToRemove = HYBRIDPlugin_Hypothesis::GetGroupsToRemove(_hyp);
1585 _hyp ? _hyp->GetToMeshHoles(true) : HYBRIDPlugin_Hypothesis::DefaultMeshHoles();
1586 const bool toMakeGroupsOfDomains = HYBRIDPlugin_Hypothesis::GetToMakeGroupsOfDomains( _hyp );
1588 helper.IsQuadraticSubMesh( theShape );
1589 helper.SetElementsOnShape( false );
1591 Ok = readGMFFile(&mgHybrid, aResultFileName.c_str(),
1593 &helper, aNodeByHybridId, aFaceByHybridId, aNodeToHybridIdMap,
1594 aNodeGroupByHybridId, anEdgeGroupByHybridId, aFaceGroupByHybridId,
1595 groupsToRemove, toMakeGroupsOfDomains, toMeshHoles);
1597 removeEmptyGroupsOfDomains( helper.GetMesh(), !toMakeGroupsOfDomains );
1601 // ---------------------
1602 // remove working files
1603 // ---------------------
1607 if ( _removeLogOnSuccess )
1608 removeFile( aLogFileName );
1610 else if ( mgHybrid.HasLog() )
1612 // get problem description from the log file
1613 _Ghs2smdsConvertor conv( aNodeByHybridId );
1614 storeErrorDescription( _logInStandardOutput ? 0 : aLogFileName.c_str(),
1615 mgHybrid.GetLog(), conv );
1617 else if ( !errStr.empty() )
1619 // the log file is empty
1620 removeFile( aLogFileName );
1621 INFOS( "HYBRID Error, " << errStr );
1622 error(COMPERR_ALGO_FAILED, errStr );
1625 if ( !_keepFiles ) {
1626 if (! Ok && _computeCanceled)
1627 removeFile( aLogFileName );
1628 removeFile( aGMFFileName );
1629 removeFile( aRequiredVerticesFileName );
1630 removeFile( aSolFileName );
1631 removeFile( aResSolFileName );
1632 removeFile( aResultFileName );
1633 removeFile( aSmdsToHybridIdMapFileName );
1635 if ( mgHybrid.IsExecutable() )
1637 std::cout << "<" << aResultFileName << "> HYBRID output file ";
1639 std::cout << "not ";
1640 std::cout << "treated !" << std::endl;
1641 std::cout << std::endl;
1645 std::cout << "MG-HYBRID " << ( Ok ? "succeeded" : "failed") << std::endl;
1651 //=============================================================================
1653 *Here we are going to use the HYBRID mesher w/o geometry
1655 //=============================================================================
1656 bool HYBRIDPlugin_HYBRID::Compute(SMESH_Mesh& theMesh,
1657 SMESH_MesherHelper* theHelper)
1659 theHelper->IsQuadraticSubMesh( theHelper->GetSubShape() );
1661 // a unique working file name
1662 // to avoid access to the same files by eg different users
1663 _genericName = HYBRIDPlugin_Hypothesis::GetFileName(_hyp);
1664 std::string aGenericName((char*) _genericName.c_str() );
1665 std::string aGenericNameRequired = aGenericName + "_required";
1667 std::string aLogFileName = aGenericName + ".log"; // log
1668 std::string aResultFileName;
1671 std::string aGMFFileName, aRequiredVerticesFileName, aSolFileName, aResSolFileName;
1672 aGMFFileName = aGenericName + ".mesh"; // GMF mesh file
1673 aResultFileName = aGenericName + "Vol.mesh"; // GMF mesh file
1674 aResSolFileName = aGenericName + "Vol.sol"; // GMF mesh file
1675 aRequiredVerticesFileName = aGenericNameRequired + ".mesh"; // GMF required vertices mesh file
1676 aSolFileName = aGenericNameRequired + ".sol"; // GMF solution file
1678 std::map <int, int> nodeID2nodeIndexMap;
1679 std::map<std::vector<double>, std::string> enfVerticesWithGroup;
1680 HYBRIDPlugin_Hypothesis::THYBRIDEnforcedVertexCoordsValues coordsSizeMap;
1681 TopoDS_Shape GeomShape;
1682 std::vector<double> coords;
1684 HYBRIDPlugin_Hypothesis::THYBRIDEnforcedVertex* enfVertex;
1686 HYBRIDPlugin_Hypothesis::THYBRIDEnforcedVertexList enfVertices = HYBRIDPlugin_Hypothesis::GetEnforcedVertices(_hyp);
1687 HYBRIDPlugin_Hypothesis::THYBRIDEnforcedVertexList::const_iterator enfVerIt = enfVertices.begin();
1689 for ( ; enfVerIt != enfVertices.end() ; ++enfVerIt)
1691 enfVertex = (*enfVerIt);
1692 if (enfVertex->coords.size()) {
1693 coordsSizeMap.insert(std::make_pair(enfVertex->coords,enfVertex->size));
1694 enfVerticesWithGroup.insert(std::make_pair(enfVertex->coords,enfVertex->groupName));
1697 GeomShape = entryToShape(enfVertex->geomEntry);
1698 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()){
1700 if (it.Value().ShapeType() == TopAbs_VERTEX){
1701 aPnt = BRep_Tool::Pnt(TopoDS::Vertex(it.Value()));
1702 coords.push_back(aPnt.X());
1703 coords.push_back(aPnt.Y());
1704 coords.push_back(aPnt.Z());
1705 if (coordsSizeMap.find(coords) == coordsSizeMap.end()) {
1706 coordsSizeMap.insert(std::make_pair(coords,enfVertex->size));
1707 enfVerticesWithGroup.insert(std::make_pair(coords,enfVertex->groupName));
1714 HYBRIDPlugin_Hypothesis::TIDSortedNodeGroupMap enforcedNodes = HYBRIDPlugin_Hypothesis::GetEnforcedNodes(_hyp);
1715 HYBRIDPlugin_Hypothesis::TIDSortedElemGroupMap enforcedEdges = HYBRIDPlugin_Hypothesis::GetEnforcedEdges(_hyp);
1716 HYBRIDPlugin_Hypothesis::TIDSortedElemGroupMap enforcedTriangles = HYBRIDPlugin_Hypothesis::GetEnforcedTriangles(_hyp);
1717 HYBRIDPlugin_Hypothesis::TID2SizeMap nodeIDToSizeMap = HYBRIDPlugin_Hypothesis::GetNodeIDToSizeMap(_hyp);
1721 int nbEnforcedVertices = coordsSizeMap.size();
1722 int nbEnforcedNodes = enforcedNodes.size();
1723 (nbEnforcedNodes <= 1) ? tmpStr = "node" : tmpStr = "nodes";
1724 std::cout << nbEnforcedNodes << " enforced " << tmpStr << " from hypo" << std::endl;
1725 (nbEnforcedVertices <= 1) ? tmpStr = "vertex" : tmpStr = "vertices";
1726 std::cout << nbEnforcedVertices << " enforced " << tmpStr << " from hypo" << std::endl;
1728 std::vector <const SMDS_MeshNode*> aNodeByHybridId, anEnforcedNodeByHybridId;
1729 std::vector <const SMDS_MeshElement*> aFaceByHybridId;
1730 std::map<const SMDS_MeshNode*,int> aNodeToHybridIdMap;
1731 std::vector<std::string> aNodeGroupByHybridId, anEdgeGroupByHybridId, aFaceGroupByHybridId;
1733 SMESH_ProxyMesh::Ptr proxyMesh( new SMESH_ProxyMesh( theMesh ));
1735 MG_HYBRID_API mgHybrid( _computeCanceled, _progress );
1737 Ok = writeGMFFile(&mgHybrid,
1738 aGMFFileName.c_str(),
1739 aRequiredVerticesFileName.c_str(), aSolFileName.c_str(),
1740 *proxyMesh, *theHelper,
1741 aNodeByHybridId, aFaceByHybridId, aNodeToHybridIdMap,
1742 aNodeGroupByHybridId, anEdgeGroupByHybridId, aFaceGroupByHybridId,
1743 enforcedNodes, enforcedEdges, enforcedTriangles,
1744 enfVerticesWithGroup, coordsSizeMap);
1746 // -----------------
1747 // run hybrid mesher
1748 // -----------------
1750 std::string cmd = HYBRIDPlugin_Hypothesis::CommandToRun( _hyp );
1752 if ( mgHybrid.IsExecutable() )
1754 cmd += " --in " + aGMFFileName;
1755 cmd += " --out " + aResultFileName;
1757 if ( !_logInStandardOutput )
1759 cmd += " 1> " + aLogFileName; // dump into file
1760 mgHybrid.SetLogFile( aLogFileName );
1762 std::cout << std::endl;
1763 std::cout << "Hybrid execution w/o geometry..." << std::endl;
1764 std::cout << cmd << std::endl;
1766 _computeCanceled = false;
1769 Ok = mgHybrid.Compute( cmd, errStr ); // run
1771 if ( _logInStandardOutput && mgHybrid.IsLibrary() )
1772 std::cout << std::endl << mgHybrid.GetLog() << std::endl;
1774 chmodUserOnly(aLogFileName.c_str());
1776 std::cout << "End of Hybrid execution !" << std::endl;
1778 if (mgHybrid.IsExecutable())
1780 chmodUserOnly(aResultFileName.c_str());
1781 chmodUserOnly(aResSolFileName.c_str());
1787 HYBRIDPlugin_Hypothesis::TSetStrings groupsToRemove = HYBRIDPlugin_Hypothesis::GetGroupsToRemove(_hyp);
1788 const bool toMakeGroupsOfDomains = HYBRIDPlugin_Hypothesis::GetToMakeGroupsOfDomains( _hyp );
1790 Ok = Ok && readGMFFile(&mgHybrid,
1791 aResultFileName.c_str(),
1793 theHelper, aNodeByHybridId, aFaceByHybridId, aNodeToHybridIdMap,
1794 aNodeGroupByHybridId, anEdgeGroupByHybridId, aFaceGroupByHybridId,
1795 groupsToRemove, toMakeGroupsOfDomains);
1797 updateMeshGroups(theHelper->GetMesh(), groupsToRemove);
1798 removeEmptyGroupsOfDomains( theHelper->GetMesh(), !toMakeGroupsOfDomains );
1801 HYBRIDPlugin_Hypothesis* that = (HYBRIDPlugin_Hypothesis*)this->_hyp;
1803 that->ClearGroupsToRemove();
1805 // ---------------------
1806 // remove working files
1807 // ---------------------
1811 if ( _removeLogOnSuccess )
1812 removeFile( aLogFileName );
1814 else if ( mgHybrid.HasLog() )
1816 // get problem description from the log file
1817 _Ghs2smdsConvertor conv( aNodeByHybridId );
1818 storeErrorDescription( _logInStandardOutput ? 0 : aLogFileName.c_str(),
1819 mgHybrid.GetLog(), conv );
1822 // the log file is empty
1823 removeFile( aLogFileName );
1824 INFOS( "HYBRID Error, command '" << cmd << "' failed" );
1825 error(COMPERR_ALGO_FAILED, "hybrid: command not found" );
1830 if (! Ok && _computeCanceled)
1831 removeFile( aLogFileName );
1832 removeFile( aGMFFileName );
1833 removeFile( aResultFileName );
1834 removeFile( aRequiredVerticesFileName );
1835 removeFile( aSolFileName );
1836 removeFile( aResSolFileName );
1841 void HYBRIDPlugin_HYBRID::CancelCompute()
1843 _computeCanceled = true;
1844 #if !defined(WIN32) && !defined(__APPLE__)
1845 std::string cmd = "ps xo pid,args | grep " + _genericName;
1846 //cmd += " | grep -e \"^ *[0-9]\\+ \\+" + HYBRIDPlugin_Hypothesis::GetExeName() + "\"";
1847 cmd += " | awk '{print $1}' | xargs kill -9 > /dev/null 2>&1";
1848 system( cmd.c_str() );
1852 //================================================================================
1854 * \brief Provide human readable text by error code reported by hybrid
1856 //================================================================================
1858 static const char* translateError(const int errNum)
1862 return "error distene 0";
1864 return "error distene 1";
1866 return "unknown distene error";
1869 //================================================================================
1871 * \brief Retrieve from a string given number of integers
1873 //================================================================================
1875 static char* getIds( char* ptr, int nbIds, std::vector<int>& ids )
1878 ids.reserve( nbIds );
1881 while ( !isdigit( *ptr )) ++ptr;
1882 if ( ptr[-1] == '-' ) --ptr;
1883 ids.push_back( strtol( ptr, &ptr, 10 ));
1889 //================================================================================
1891 * \brief Retrieve problem description form a log file
1892 * \retval bool - always false
1894 //================================================================================
1896 bool HYBRIDPlugin_HYBRID::storeErrorDescription(const char* logFile,
1897 const std::string& log,
1898 const _Ghs2smdsConvertor & toSmdsConvertor )
1900 if(_computeCanceled)
1901 return error(SMESH_Comment("interruption initiated by user"));
1903 char* ptr = const_cast<char*>( log.c_str() );
1904 char* buf = ptr, * bufEnd = ptr + log.size();
1906 SMESH_Comment errDescription;
1908 enum { NODE = 1, EDGE, TRIA, VOL, SKIP_ID = 1 };
1910 // look for MeshGems version
1911 // Since "MG-TETRA -- MeshGems 1.1-3 (January, 2013)" error codes change.
1912 // To discriminate old codes from new ones we add 1000000 to the new codes.
1913 // This way value of the new codes is same as absolute value of codes printed
1914 // in the log after "MGMESSAGE" string.
1915 int versionAddition = 0;
1918 while ( ++verPtr < bufEnd )
1920 if ( strncmp( verPtr, "MG-TETRA -- MeshGems ", 21 ) != 0 )
1922 if ( strcmp( verPtr, "MG-TETRA -- MeshGems 1.1-3 " ) >= 0 )
1923 versionAddition = 1000000;
1929 // look for errors "ERR #"
1931 std::set<std::string> foundErrorStr; // to avoid reporting same error several times
1932 std::set<int> elemErrorNums; // not to report different types of errors with bad elements
1933 while ( ++ptr < bufEnd )
1935 if ( strncmp( ptr, "ERR ", 4 ) != 0 )
1938 std::list<const SMDS_MeshElement*> badElems;
1939 std::vector<int> nodeIds;
1943 int errNum = strtol(ptr, &ptr, 10) + versionAddition;
1944 // we treat errors enumerated in [SALOME platform 0019316] issue
1945 // and all errors from a new (Release 1.1) MeshGems User Manual
1947 case 0015: // The face number (numfac) with vertices (f 1, f 2, f 3) has a null vertex.
1948 case 1005620 : // a too bad quality face is detected. This face is considered degenerated.
1949 ptr = getIds(ptr, SKIP_ID, nodeIds);
1950 ptr = getIds(ptr, TRIA, nodeIds);
1951 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1953 case 1005621 : // a too bad quality face is detected. This face is degenerated.
1954 // hence the is degenerated it is invisible, add its edges in addition
1955 ptr = getIds(ptr, SKIP_ID, nodeIds);
1956 ptr = getIds(ptr, TRIA, nodeIds);
1957 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1959 std::vector<int> edgeNodes( nodeIds.begin(), --nodeIds.end() ); // 01
1960 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1961 edgeNodes[1] = nodeIds[2]; // 02
1962 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
1963 edgeNodes[0] = nodeIds[1]; // 12
1966 case 1000: // Face (f 1, f 2, f 3) appears more than once in the input surface mesh.
1968 case 1002: // Face (f 1, f 2, f 3) has a vertex negative or null
1969 case 3019: // Constrained face (f 1, f 2, f 3) cannot be enforced
1970 case 1002211: // a face has a vertex negative or null.
1971 case 1005200 : // a surface mesh appears more than once in the input surface mesh.
1972 case 1008423 : // a constrained face cannot be enforced (regeneration phase failed).
1973 ptr = getIds(ptr, TRIA, nodeIds);
1974 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1976 case 1001: // Edge (e1, e2) appears more than once in the input surface mesh
1977 case 3009: // Constrained edge (e1, e2) cannot be enforced (warning).
1978 // ERR 3109 : EDGE 5 6 UNIQUE
1979 case 3109: // Edge (e1, e2) is unique (i.e., bounds a hole in the surface)
1980 case 1005210 : // an edge appears more than once in the input surface mesh.
1981 case 1005820 : // an edge is unique (i.e., bounds a hole in the surface).
1982 case 1008441 : // a constrained edge cannot be enforced.
1983 ptr = getIds(ptr, EDGE, nodeIds);
1984 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1986 case 2004: // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
1987 case 2014: // at least two points whose distance is dist, i.e., considered as coincident
1988 case 2103: // Vertex v1 and vertex v2 are too close to one another or coincident (warning).
1989 // ERR 2103 : 16 WITH 3
1990 case 1005105 : // two vertices are too close to one another or coincident.
1991 case 1005107: // Two vertices are too close to one another or coincident.
1992 ptr = getIds(ptr, NODE, nodeIds);
1993 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1994 ptr = getIds(ptr, NODE, nodeIds);
1995 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
1997 case 2012: // Vertex v1 cannot be inserted (warning).
1998 case 1005106 : // a vertex cannot be inserted.
1999 ptr = getIds(ptr, NODE, nodeIds);
2000 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2002 case 3103: // The surface edge (e1, e2) intersects another surface edge (e3, e4)
2003 case 1005110 : // two surface edges are intersecting.
2004 // ERR 3103 : 1 2 WITH 7 3
2005 ptr = getIds(ptr, EDGE, nodeIds);
2006 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2007 ptr = getIds(ptr, EDGE, nodeIds);
2008 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2010 case 3104: // The surface edge (e1, e2) intersects the surface face (f 1, f 2, f 3)
2011 // ERR 3104 : 9 10 WITH 1 2 3
2012 case 3106: // One surface edge (say e1, e2) intersects a surface face (f 1, f 2, f 3)
2013 case 1005120 : // a surface edge intersects a surface face.
2014 ptr = getIds(ptr, EDGE, nodeIds);
2015 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2016 ptr = getIds(ptr, TRIA, nodeIds);
2017 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2019 case 3105: // One boundary point (say p1) lies within a surface face (f 1, f 2, f 3)
2020 // ERR 3105 : 8 IN 2 3 5
2021 case 1005150 : // a boundary point lies within a surface face.
2022 ptr = getIds(ptr, NODE, nodeIds);
2023 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2024 ptr = getIds(ptr, TRIA, nodeIds);
2025 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2027 case 3107: // One boundary point (say p1) lies within a surface edge (e1, e2) (stop).
2028 // ERR 3107 : 2 IN 4 1
2029 case 1005160 : // a boundary point lies within a surface edge.
2030 ptr = getIds(ptr, NODE, nodeIds);
2031 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2032 ptr = getIds(ptr, EDGE, nodeIds);
2033 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2035 case 9000: // ERR 9000
2036 // ELEMENT 261 WITH VERTICES : 7 396 -8 242
2037 // VOLUME : -1.11325045E+11 W.R.T. EPSILON 0.
2038 // A too small volume element is detected. Are reported the index of the element,
2039 // its four vertex indices, its volume and the tolerance threshold value
2040 ptr = getIds(ptr, SKIP_ID, nodeIds);
2041 ptr = getIds(ptr, VOL, nodeIds);
2042 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2043 // even if all nodes found, volume it most probably invisible,
2044 // add its faces to demonstrate it anyhow
2046 std::vector<int> faceNodes( nodeIds.begin(), --nodeIds.end() ); // 012
2047 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
2048 faceNodes[2] = nodeIds[3]; // 013
2049 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
2050 faceNodes[1] = nodeIds[2]; // 023
2051 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
2052 faceNodes[0] = nodeIds[1]; // 123
2053 badElems.push_back( toSmdsConvertor.getElement(faceNodes));
2056 case 9001: // ERR 9001
2057 // %% NUMBER OF NEGATIVE VOLUME TETS : 1
2058 // %% THE LARGEST NEGATIVE TET : 1.75376581E+11
2059 // %% NUMBER OF NULL VOLUME TETS : 0
2060 // There exists at least a null or negative volume element
2063 // There exist n null or negative volume elements
2066 // A too small volume element is detected
2069 // A too bad quality face is detected. This face is considered degenerated,
2070 // its index, its three vertex indices together with its quality value are reported
2071 break; // same as next
2072 case 9112: // ERR 9112
2073 // FACE 2 WITH VERTICES : 4 2 5
2074 // SMALL INRADIUS : 0.
2075 // A too bad quality face is detected. This face is degenerated,
2076 // its index, its three vertex indices together with its inradius are reported
2077 ptr = getIds(ptr, SKIP_ID, nodeIds);
2078 ptr = getIds(ptr, TRIA, nodeIds);
2079 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2080 // add triangle edges as it most probably has zero area and hence invisible
2082 std::vector<int> edgeNodes(2);
2083 edgeNodes[0] = nodeIds[0]; edgeNodes[1] = nodeIds[1]; // 0-1
2084 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
2085 edgeNodes[1] = nodeIds[2]; // 0-2
2086 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
2087 edgeNodes[0] = nodeIds[1]; // 1-2
2088 badElems.push_back( toSmdsConvertor.getElement(edgeNodes));
2091 case 1005103 : // the vertices of an element are too close to one another or coincident.
2092 ptr = getIds(ptr, TRIA, nodeIds);
2093 if ( nodeIds.back() == 0 ) // index of the third vertex of the element (0 for an edge)
2094 nodeIds.resize( EDGE );
2095 badElems.push_back( toSmdsConvertor.getElement(nodeIds));
2099 bool isNewError = foundErrorStr.insert( std::string( errBeg, ptr )).second;
2101 continue; // not to report same error several times
2103 // const SMDS_MeshElement* nullElem = 0;
2104 // bool allElemsOk = ( find( badElems.begin(), badElems.end(), nullElem) == badElems.end());
2106 // if ( allElemsOk && !badElems.empty() && !elemErrorNums.empty() ) {
2107 // bool oneMoreErrorType = elemErrorNums.insert( errNum ).second;
2108 // if ( oneMoreErrorType )
2109 // continue; // not to report different types of errors with bad elements
2112 // store bad elements
2113 //if ( allElemsOk ) {
2114 std::list<const SMDS_MeshElement*>::iterator elem = badElems.begin();
2115 for ( ; elem != badElems.end(); ++elem )
2116 addBadInputElement( *elem );
2120 std::string text = translateError( errNum );
2121 if ( errDescription.find( text ) == text.npos ) {
2122 if ( !errDescription.empty() )
2123 errDescription << "\n";
2124 errDescription << text;
2129 if ( errDescription.empty() ) { // no errors found
2130 char msgLic1[] = "connection to server failed";
2131 char msgLic2[] = " Dlim ";
2132 char msgLic3[] = "license is not valid";
2133 if ( std::search( &buf[0], bufEnd, msgLic1, msgLic1 + strlen(msgLic1)) != bufEnd ||
2134 std::search( &buf[0], bufEnd, msgLic2, msgLic2 + strlen(msgLic2)) != bufEnd )
2135 errDescription << "Network license problem.";
2136 else if ( std::search( &buf[0], bufEnd, msgLic3, msgLic3 + strlen(msgLic3)) != bufEnd )
2137 errDescription << "License is not valid.";
2140 char msg2[] = "SEGMENTATION FAULT";
2141 if ( std::search( &buf[0], bufEnd, msg2, msg2 + strlen(msg2)) != bufEnd )
2142 errDescription << "hybrid: SEGMENTATION FAULT. ";
2146 if ( logFile && logFile[0] )
2148 if ( errDescription.empty() )
2149 errDescription << "See " << logFile << " for problem description";
2151 errDescription << "\nSee " << logFile << " for more information";
2153 return error( errDescription );
2156 //================================================================================
2158 * \brief Creates _Ghs2smdsConvertor
2160 //================================================================================
2162 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const std::map <int,const SMDS_MeshNode*> & ghs2NodeMap)
2163 :_ghs2NodeMap( & ghs2NodeMap ), _nodeByGhsId( 0 )
2167 //================================================================================
2169 * \brief Creates _Ghs2smdsConvertor
2171 //================================================================================
2173 _Ghs2smdsConvertor::_Ghs2smdsConvertor( const std::vector <const SMDS_MeshNode*> & nodeByGhsId)
2174 : _ghs2NodeMap( 0 ), _nodeByGhsId( &nodeByGhsId )
2178 //================================================================================
2180 * \brief Return SMDS element by ids of HYBRID nodes
2182 //================================================================================
2184 const SMDS_MeshElement* _Ghs2smdsConvertor::getElement(const std::vector<int>& ghsNodes) const
2186 size_t nbNodes = ghsNodes.size();
2187 std::vector<const SMDS_MeshNode*> nodes( nbNodes, 0 );
2188 for ( size_t i = 0; i < nbNodes; ++i ) {
2189 int ghsNode = ghsNodes[ i ];
2190 if ( _ghs2NodeMap ) {
2191 std::map <int,const SMDS_MeshNode*>::const_iterator in = _ghs2NodeMap->find( ghsNode);
2192 if ( in == _ghs2NodeMap->end() )
2194 nodes[ i ] = in->second;
2197 if ( ghsNode < 1 || ghsNode > (int)_nodeByGhsId->size() )
2199 nodes[ i ] = (*_nodeByGhsId)[ ghsNode-1 ];
2205 if ( nbNodes == 2 ) {
2206 const SMDS_MeshElement* edge= SMDS_Mesh::FindEdge( nodes[0], nodes[1] );
2208 edge = new SMDS_LinearEdge( nodes[0], nodes[1] );
2211 if ( nbNodes == 3 ) {
2212 const SMDS_MeshElement* face = SMDS_Mesh::FindFace( nodes );
2214 face = new SMDS_FaceOfNodes( nodes[0], nodes[1], nodes[2] );
2218 return new SMDS_VolumeOfNodes( nodes[0], nodes[1], nodes[2], nodes[3] );
2224 //=============================================================================
2228 //=============================================================================
2229 bool HYBRIDPlugin_HYBRID::Evaluate(SMESH_Mesh& aMesh,
2230 const TopoDS_Shape& aShape,
2231 MapShapeNbElems& aResMap)
2233 int nbtri = 0, nbqua = 0;
2234 double fullArea = 0.0;
2235 for (TopExp_Explorer exp(aShape, TopAbs_FACE); exp.More(); exp.Next()) {
2236 TopoDS_Face F = TopoDS::Face( exp.Current() );
2237 SMESH_subMesh *sm = aMesh.GetSubMesh(F);
2238 MapShapeNbElemsItr anIt = aResMap.find(sm);
2239 if( anIt==aResMap.end() ) {
2240 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
2241 smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,
2242 "Submesh can not be evaluated",this));
2245 std::vector<smIdType> aVec = (*anIt).second;
2246 nbtri += std::max(aVec[SMDSEntity_Triangle],aVec[SMDSEntity_Quad_Triangle]);
2247 nbqua += std::max(aVec[SMDSEntity_Quadrangle],aVec[SMDSEntity_Quad_Quadrangle]);
2249 BRepGProp::SurfaceProperties(F,G);
2250 double anArea = G.Mass();
2254 // collect info from edges
2255 int nb0d_e = 0, nb1d_e = 0;
2256 bool IsQuadratic = false;
2257 bool IsFirst = true;
2258 TopTools_MapOfShape tmpMap;
2259 for (TopExp_Explorer exp(aShape, TopAbs_EDGE); exp.More(); exp.Next()) {
2260 TopoDS_Edge E = TopoDS::Edge(exp.Current());
2261 if( tmpMap.Contains(E) )
2264 SMESH_subMesh *aSubMesh = aMesh.GetSubMesh(exp.Current());
2265 MapShapeNbElemsItr anIt = aResMap.find(aSubMesh);
2266 std::vector<smIdType> aVec = (*anIt).second;
2267 nb0d_e += aVec[SMDSEntity_Node];
2268 nb1d_e += std::max(aVec[SMDSEntity_Edge],aVec[SMDSEntity_Quad_Edge]);
2270 IsQuadratic = (aVec[SMDSEntity_Quad_Edge] > aVec[SMDSEntity_Edge]);
2276 double ELen = sqrt(2.* ( fullArea/(nbtri+nbqua*2) ) / sqrt(3.0) );
2279 BRepGProp::VolumeProperties(aShape,G);
2280 double aVolume = G.Mass();
2281 double tetrVol = 0.1179*ELen*ELen*ELen;
2282 double CoeffQuality = 0.9;
2283 int nbVols = int(aVolume/tetrVol/CoeffQuality);
2284 int nb1d_f = (nbtri*3 + nbqua*4 - nb1d_e) / 2;
2285 int nb1d_in = (int) ( nbVols*6 - nb1d_e - nb1d_f ) / 5;
2286 std::vector<smIdType> aVec(SMDSEntity_Last);
2287 for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i]=0;
2289 aVec[SMDSEntity_Node] = nb1d_in/6 + 1 + nb1d_in;
2290 aVec[SMDSEntity_Quad_Tetra] = nbVols - nbqua*2;
2291 aVec[SMDSEntity_Quad_Pyramid] = nbqua;
2294 aVec[SMDSEntity_Node] = nb1d_in/6 + 1;
2295 aVec[SMDSEntity_Tetra] = nbVols - nbqua*2;
2296 aVec[SMDSEntity_Pyramid] = nbqua;
2298 SMESH_subMesh *sm = aMesh.GetSubMesh(aShape);
2299 aResMap.insert(std::make_pair(sm,aVec));
2304 bool HYBRIDPlugin_HYBRID::importGMFMesh(const char* theGMFFileName, SMESH_Mesh& theMesh)
2306 SMESH_ComputeErrorPtr err = theMesh.GMFToMesh( theGMFFileName, /*makeRequiredGroups =*/ true );
2308 theMesh.GetMeshDS()->Modified();
2310 return ( !err || err->IsOK());
2315 //================================================================================
2317 * \brief Sub-mesh event listener setting enforced elements as soon as an enforced
2320 struct _EnforcedMeshRestorer : public SMESH_subMeshEventListener
2322 _EnforcedMeshRestorer():
2323 SMESH_subMeshEventListener( /*isDeletable = */true, Name() )
2326 //================================================================================
2328 * \brief Returns an ID of listener
2330 static const char* Name() { return "HYBRIDPlugin_HYBRID::_EnforcedMeshRestorer"; }
2332 //================================================================================
2334 * \brief Treat events of the subMesh
2336 void ProcessEvent(const int event,
2337 const int eventType,
2338 SMESH_subMesh* /*subMesh*/,
2339 SMESH_subMeshEventListenerData* data,
2340 const SMESH_Hypothesis* /*hyp*/)
2342 if ( SMESH_subMesh::SUBMESH_LOADED == event &&
2343 SMESH_subMesh::COMPUTE_EVENT == eventType &&
2345 !data->mySubMeshes.empty() )
2347 // An enforced mesh (subMesh->_father) has been loaded from hdf file
2348 if ( HYBRIDPlugin_Hypothesis* hyp = GetGHSHypothesis( data->mySubMeshes.front() ))
2349 hyp->RestoreEnfElemsByMeshes();
2352 //================================================================================
2354 * \brief Returns HYBRIDPlugin_Hypothesis used to compute a subMesh
2356 static HYBRIDPlugin_Hypothesis* GetGHSHypothesis( SMESH_subMesh* subMesh )
2358 SMESH_HypoFilter ghsHypFilter( SMESH_HypoFilter::HasName( "HYBRID_Parameters" ));
2359 return (HYBRIDPlugin_Hypothesis* )
2360 subMesh->GetFather()->GetHypothesis( subMesh->GetSubShape(),
2362 /*visitAncestors=*/true);
2366 //================================================================================
2368 * \brief Sub-mesh event listener removing empty groups created due to "To make
2369 * groups of domains".
2371 struct _GroupsOfDomainsRemover : public SMESH_subMeshEventListener
2373 _GroupsOfDomainsRemover():
2374 SMESH_subMeshEventListener( /*isDeletable = */true,
2375 "HYBRIDPlugin_HYBRID::_GroupsOfDomainsRemover" ) {}
2377 * \brief Treat events of the subMesh
2379 void ProcessEvent(const int /*event*/,
2380 const int eventType,
2381 SMESH_subMesh* subMesh,
2382 SMESH_subMeshEventListenerData* /*data*/,
2383 const SMESH_Hypothesis* /*hyp*/)
2385 if (SMESH_subMesh::ALGO_EVENT == eventType &&
2386 !subMesh->GetAlgo() )
2388 removeEmptyGroupsOfDomains( subMesh->GetFather(), /*notEmptyAsWell=*/true );
2394 //================================================================================
2396 * \brief Set an event listener to set enforced elements as soon as an enforced
2399 //================================================================================
2401 void HYBRIDPlugin_HYBRID::SubmeshRestored(SMESH_subMesh* subMesh)
2403 if ( HYBRIDPlugin_Hypothesis* hyp = _EnforcedMeshRestorer::GetGHSHypothesis( subMesh ))
2405 HYBRIDPlugin_Hypothesis::THYBRIDEnforcedMeshList enfMeshes = hyp->_GetEnforcedMeshes();
2406 HYBRIDPlugin_Hypothesis::THYBRIDEnforcedMeshList::iterator it = enfMeshes.begin();
2407 for(;it != enfMeshes.end();++it) {
2408 HYBRIDPlugin_Hypothesis::THYBRIDEnforcedMesh* enfMesh = *it;
2409 if ( SMESH_Mesh* mesh = GetMeshByPersistentID( enfMesh->persistID ))
2411 SMESH_subMesh* smToListen = mesh->GetSubMesh( mesh->GetShapeToMesh() );
2412 // a listener set to smToListen will care of hypothesis stored in SMESH_EventListenerData
2413 subMesh->SetEventListener( new _EnforcedMeshRestorer(),
2414 SMESH_subMeshEventListenerData::MakeData( subMesh ),
2421 //================================================================================
2423 * \brief Sets an event listener removing empty groups created due to "To make
2424 * groups of domains".
2425 * \param subMesh - submesh where algo is set
2427 * This method is called when a submesh gets HYP_OK algo_state.
2428 * After being set, event listener is notified on each event of a submesh.
2430 //================================================================================
2432 void HYBRIDPlugin_HYBRID::SetEventListener(SMESH_subMesh* subMesh)
2434 subMesh->SetEventListener( new _GroupsOfDomainsRemover(), 0, subMesh );