1 // Copyright (C) 2007-2015 CEA/DEN, EDF R&D, OPEN CASCADE
3 // Copyright (C) 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
6 // This library is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License, or (at your option) any later version.
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14 // Lesser General Public License for more details.
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
22 // File : StdMeshers_QuadFromMedialAxis_1D2D.cxx
23 // Created : Wed Jun 3 17:33:45 2015
24 // Author : Edward AGAPOV (eap)
26 #include "StdMeshers_QuadFromMedialAxis_1D2D.hxx"
28 #include "SMESH_Block.hxx"
29 #include "SMESH_Gen.hxx"
30 #include "SMESH_MAT2d.hxx"
31 #include "SMESH_Mesh.hxx"
32 #include "SMESH_MeshEditor.hxx"
33 #include "SMESH_MesherHelper.hxx"
34 #include "SMESH_ProxyMesh.hxx"
35 #include "SMESH_subMesh.hxx"
36 #include "SMESH_subMeshEventListener.hxx"
37 #include "StdMeshers_FaceSide.hxx"
38 #include "StdMeshers_LayerDistribution.hxx"
39 #include "StdMeshers_NumberOfLayers.hxx"
40 #include "StdMeshers_Regular_1D.hxx"
41 #include "StdMeshers_ViscousLayers2D.hxx"
43 #include <BRepAdaptor_Curve.hxx>
44 #include <BRepBuilderAPI_MakeEdge.hxx>
45 #include <BRepTools.hxx>
46 #include <BRep_Tool.hxx>
47 #include <GeomAPI_Interpolate.hxx>
48 #include <Geom_Surface.hxx>
49 #include <Precision.hxx>
50 #include <TColgp_HArray1OfPnt.hxx>
52 #include <TopExp_Explorer.hxx>
53 #include <TopLoc_Location.hxx>
54 #include <TopTools_MapOfShape.hxx>
56 #include <TopoDS_Edge.hxx>
57 #include <TopoDS_Face.hxx>
58 #include <TopoDS_Vertex.hxx>
64 //================================================================================
68 class StdMeshers_QuadFromMedialAxis_1D2D::Algo1D : public StdMeshers_Regular_1D
71 Algo1D(int studyId, SMESH_Gen* gen):
72 StdMeshers_Regular_1D( gen->GetANewId(), studyId, gen )
75 void SetSegmentLength( double len )
77 SMESH_Algo::_usedHypList.clear();
78 _value[ BEG_LENGTH_IND ] = len;
79 _value[ PRECISION_IND ] = 1e-7;
80 _hypType = LOCAL_LENGTH;
82 void SetRadialDistribution( const SMESHDS_Hypothesis* hyp )
84 SMESH_Algo::_usedHypList.clear();
88 if ( const StdMeshers_NumberOfLayers* nl =
89 dynamic_cast< const StdMeshers_NumberOfLayers* >( hyp ))
91 _ivalue[ NB_SEGMENTS_IND ] = nl->GetNumberOfLayers();
92 _ivalue[ DISTR_TYPE_IND ] = 0;
93 _hypType = NB_SEGMENTS;
95 if ( const StdMeshers_LayerDistribution* ld =
96 dynamic_cast< const StdMeshers_LayerDistribution* >( hyp ))
98 if ( SMESH_Hypothesis* h = ld->GetLayerDistribution() )
100 SMESH_Algo::_usedHypList.clear();
101 SMESH_Algo::_usedHypList.push_back( h );
105 void ComputeDistribution(SMESH_MesherHelper& theHelper,
106 const gp_Pnt& thePnt1,
107 const gp_Pnt& thePnt2,
108 list< double >& theParams)
110 SMESH_Mesh& mesh = *theHelper.GetMesh();
111 TopoDS_Edge edge = BRepBuilderAPI_MakeEdge( thePnt1, thePnt2 );
113 SMESH_Hypothesis::Hypothesis_Status aStatus;
114 CheckHypothesis( mesh, edge, aStatus );
117 BRepAdaptor_Curve C3D(edge);
118 double f = C3D.FirstParameter(), l = C3D.LastParameter(), len = thePnt1.Distance( thePnt2 );
119 if ( !StdMeshers_Regular_1D::computeInternalParameters( mesh, C3D, len, f, l, theParams, false))
121 for ( size_t i = 1; i < 15; ++i )
122 theParams.push_back( i/15 );
126 for (list<double>::iterator itU = theParams.begin(); itU != theParams.end(); ++itU )
130 virtual const list <const SMESHDS_Hypothesis *> &
131 GetUsedHypothesis(SMESH_Mesh &, const TopoDS_Shape &, const bool)
133 return SMESH_Algo::_usedHypList;
135 virtual bool CheckHypothesis(SMESH_Mesh& aMesh,
136 const TopoDS_Shape& aShape,
137 SMESH_Hypothesis::Hypothesis_Status& aStatus)
139 if ( !SMESH_Algo::_usedHypList.empty() )
140 return StdMeshers_Regular_1D::CheckHypothesis( aMesh, aShape, aStatus );
145 //================================================================================
147 * \brief Constructor sets algo features
149 //================================================================================
151 StdMeshers_QuadFromMedialAxis_1D2D::StdMeshers_QuadFromMedialAxis_1D2D(int hypId,
154 : StdMeshers_Quadrangle_2D(hypId, studyId, gen),
157 _name = "QuadFromMedialAxis_1D2D";
158 _shapeType = (1 << TopAbs_FACE);
159 _onlyUnaryInput = true; // FACE by FACE so far
160 _requireDiscreteBoundary = false; // make 1D by myself
161 _supportSubmeshes = true; // make 1D by myself
162 _neededLowerHyps[ 1 ] = true; // suppress warning on hiding a global 1D algo
163 _neededLowerHyps[ 2 ] = true; // suppress warning on hiding a global 2D algo
164 _compatibleHypothesis.clear();
165 _compatibleHypothesis.push_back("ViscousLayers2D");
166 _compatibleHypothesis.push_back("LayerDistribution2D");
167 _compatibleHypothesis.push_back("NumberOfLayers2D");
170 //================================================================================
174 //================================================================================
176 StdMeshers_QuadFromMedialAxis_1D2D::~StdMeshers_QuadFromMedialAxis_1D2D()
182 //================================================================================
184 * \brief Check if needed hypotheses are present
186 //================================================================================
188 bool StdMeshers_QuadFromMedialAxis_1D2D::CheckHypothesis(SMESH_Mesh& aMesh,
189 const TopoDS_Shape& aShape,
190 Hypothesis_Status& aStatus)
194 // get one main optional hypothesis
195 const list <const SMESHDS_Hypothesis * >& hyps = GetUsedHypothesis(aMesh, aShape);
196 _hyp2D = hyps.empty() ? 0 : hyps.front();
198 return true; // does not require hypothesis
203 typedef map< const SMDS_MeshNode*, list< const SMDS_MeshNode* > > TMergeMap;
205 //================================================================================
207 * \brief Sinuous face
211 FaceQuadStruct::Ptr _quad;
212 vector< TopoDS_Edge > _edges;
213 vector< TopoDS_Edge > _sinuSide[2], _shortSide[2];
214 vector< TopoDS_Edge > _sinuEdges;
215 vector< Handle(Geom_Curve) > _sinuCurves;
217 list< int > _nbEdgesInWire;
218 TMergeMap _nodesToMerge;
220 SinuousFace( const TopoDS_Face& f ): _quad( new FaceQuadStruct )
222 list< TopoDS_Edge > edges;
223 _nbWires = SMESH_Block::GetOrderedEdges (f, edges, _nbEdgesInWire);
224 _edges.assign( edges.begin(), edges.end() );
226 _quad->side.resize( 4 );
229 const TopoDS_Face& Face() const { return _quad->face; }
230 bool IsRing() const { return _shortSide[0].empty() && !_sinuSide[0].empty(); }
233 //================================================================================
235 * \brief Temporary mesh
237 struct TmpMesh : public SMESH_Mesh
241 _myMeshDS = new SMESHDS_Mesh(/*id=*/0, /*isEmbeddedMode=*/true);
245 //================================================================================
247 * \brief Event listener which removes mesh from EDGEs when 2D hyps change
249 struct EdgeCleaner : public SMESH_subMeshEventListener
253 SMESH_subMeshEventListener( /*isDeletable=*/true,
254 "StdMeshers_QuadFromMedialAxis_1D2D::EdgeCleaner")
258 virtual void ProcessEvent(const int event,
260 SMESH_subMesh* faceSubMesh,
261 SMESH_subMeshEventListenerData* data,
262 const SMESH_Hypothesis* hyp)
264 if ( eventType == SMESH_subMesh::ALGO_EVENT )
266 _prevAlgoEvent = event;
269 // SMESH_subMesh::COMPUTE_EVENT
270 if ( _prevAlgoEvent == SMESH_subMesh::REMOVE_HYP ||
271 _prevAlgoEvent == SMESH_subMesh::REMOVE_ALGO ||
272 _prevAlgoEvent == SMESH_subMesh::MODIF_HYP )
274 SMESH_subMeshIteratorPtr smIt = faceSubMesh->getDependsOnIterator(/*includeSelf=*/false);
275 while ( smIt->more() )
276 smIt->next()->ComputeStateEngine( SMESH_subMesh::CLEAN );
282 //================================================================================
284 * \brief Return a member of a std::pair
286 //================================================================================
288 template< typename T >
289 T& get( std::pair< T, T >& thePair, bool is2nd )
291 return is2nd ? thePair.second : thePair.first;
294 //================================================================================
296 * \brief Select two EDGEs from a map, either mapped to least values or to max values
298 //================================================================================
300 // template< class TVal2EdgesMap >
301 // void getTwo( bool least,
302 // TVal2EdgesMap& map,
303 // vector<TopoDS_Edge>& twoEdges,
304 // vector<TopoDS_Edge>& otherEdges)
307 // otherEdges.clear();
310 // TVal2EdgesMap::iterator i = map.begin();
311 // twoEdges.push_back( i->second );
312 // twoEdges.push_back( ++i->second );
313 // for ( ; i != map.end(); ++i )
314 // otherEdges.push_back( i->second );
318 // TVal2EdgesMap::reverse_iterator i = map.rbegin();
319 // twoEdges.push_back( i->second );
320 // twoEdges.push_back( ++i->second );
321 // for ( ; i != map.rend(); ++i )
322 // otherEdges.push_back( i->second );
325 // if ( TopExp::CommonVertex( twoEdges[0], twoEdges[1], v ))
327 // twoEdges.clear(); // two EDGEs must not be connected
328 // otherEdges.clear();
332 //================================================================================
334 * \brief Finds out a minimal segment length given EDGEs will be divided into.
335 * This length is further used to discretize the Medial Axis
337 //================================================================================
339 double getMinSegLen(SMESH_MesherHelper& theHelper,
340 const vector<TopoDS_Edge>& theEdges)
343 SMESH_Mesh* mesh = theHelper.GetMesh();
345 vector< SMESH_Algo* > algos( theEdges.size() );
346 for ( size_t i = 0; i < theEdges.size(); ++i )
348 SMESH_subMesh* sm = mesh->GetSubMesh( theEdges[i] );
349 algos[i] = sm->GetAlgo();
352 int nbSegDflt = mesh->GetGen() ? mesh->GetGen()->GetDefaultNbSegments() : 15;
353 double minSegLen = Precision::Infinite();
355 for ( size_t i = 0; i < theEdges.size(); ++i )
357 SMESH_subMesh* sm = mesh->GetSubMesh( theEdges[i] );
358 if ( SMESH_Algo::IsStraight( theEdges[i], /*degenResult=*/true ))
361 size_t iOpp = ( theEdges.size() == 4 ? (i+2)%4 : i );
362 SMESH_Algo* algo = sm->GetAlgo();
363 if ( !algo ) algo = algos[ iOpp ];
365 SMESH_Hypothesis::Hypothesis_Status status = SMESH_Hypothesis::HYP_MISSING;
368 if ( !algo->CheckHypothesis( *mesh, theEdges[i], status ))
369 algo->CheckHypothesis( *mesh, theEdges[iOpp], status );
372 if ( status != SMESH_Hypothesis::HYP_OK )
374 minSegLen = Min( minSegLen, SMESH_Algo::EdgeLength( theEdges[i] ) / nbSegDflt );
379 tmpMesh.ShapeToMesh( TopoDS_Shape());
380 tmpMesh.ShapeToMesh( theEdges[i] );
382 if ( !mesh->GetGen() ) continue; // tmp mesh
383 mesh->GetGen()->Compute( tmpMesh, theEdges[i], true, true ); // make nodes on VERTEXes
384 if ( !algo->Compute( tmpMesh, theEdges[i] ))
390 SMDS_EdgeIteratorPtr segIt = tmpMesh.GetMeshDS()->edgesIterator();
391 while ( segIt->more() )
393 const SMDS_MeshElement* seg = segIt->next();
394 double len = SMESH_TNodeXYZ( seg->GetNode(0) ).Distance( seg->GetNode(1) );
395 minSegLen = Min( minSegLen, len );
399 if ( Precision::IsInfinite( minSegLen ))
400 minSegLen = mesh->GetShapeDiagonalSize() / nbSegDflt;
405 //================================================================================
407 * \brief Returns EDGEs located between two VERTEXes at which given MA branches end
408 * \param [in] br1 - one MA branch
409 * \param [in] br2 - one more MA branch
410 * \param [in] allEdges - all EDGEs of a FACE
411 * \param [out] shortEdges - the found EDGEs
412 * \return bool - is OK or not
414 //================================================================================
416 bool getConnectedEdges( const SMESH_MAT2d::Branch* br1,
417 const SMESH_MAT2d::Branch* br2,
418 const vector<TopoDS_Edge>& allEdges,
419 vector<TopoDS_Edge>& shortEdges)
421 vector< size_t > edgeIDs[4];
422 br1->getGeomEdges( edgeIDs[0], edgeIDs[1] );
423 br2->getGeomEdges( edgeIDs[2], edgeIDs[3] );
425 // EDGEs returned by a Branch form a connected chain with a VERTEX where
426 // the Branch ends at the chain middle. One of end EDGEs of the chain is common
427 // with either end EDGE of the chain of the other Branch, or the chains are connected
428 // at a common VERTEX;
430 // Get indices of end EDGEs of the branches
431 bool vAtStart1 = ( br1->getEnd(0)->_type == SMESH_MAT2d::BE_ON_VERTEX );
432 bool vAtStart2 = ( br2->getEnd(0)->_type == SMESH_MAT2d::BE_ON_VERTEX );
434 vAtStart1 ? edgeIDs[0].back() : edgeIDs[0][0],
435 vAtStart1 ? edgeIDs[1].back() : edgeIDs[1][0],
436 vAtStart2 ? edgeIDs[2].back() : edgeIDs[2][0],
437 vAtStart2 ? edgeIDs[3].back() : edgeIDs[3][0]
440 set< size_t > connectedIDs;
441 TopoDS_Vertex vCommon;
442 // look for the same EDGEs
443 for ( int i = 0; i < 2; ++i )
444 for ( int j = 2; j < 4; ++j )
445 if ( iEnd[i] == iEnd[j] )
447 connectedIDs.insert( edgeIDs[i].begin(), edgeIDs[i].end() );
448 connectedIDs.insert( edgeIDs[j].begin(), edgeIDs[j].end() );
451 if ( connectedIDs.empty() )
452 // look for connected EDGEs
453 for ( int i = 0; i < 2; ++i )
454 for ( int j = 2; j < 4; ++j )
455 if ( TopExp::CommonVertex( allEdges[ iEnd[i]], allEdges[ iEnd[j]], vCommon ))
457 connectedIDs.insert( edgeIDs[i].begin(), edgeIDs[i].end() );
458 connectedIDs.insert( edgeIDs[j].begin(), edgeIDs[j].end() );
461 if ( connectedIDs.empty() || // nothing
462 allEdges.size() - connectedIDs.size() < 2 ) // too many
465 // set shortEdges in the order as in allEdges
466 if ( connectedIDs.count( 0 ) &&
467 connectedIDs.count( allEdges.size()-1 ))
469 size_t iE = allEdges.size()-1;
470 while ( connectedIDs.count( iE-1 ))
472 for ( size_t i = 0; i < connectedIDs.size(); ++i )
474 shortEdges.push_back( allEdges[ iE ]);
475 iE = ( iE + 1 ) % allEdges.size();
480 set< size_t >::iterator i = connectedIDs.begin();
481 for ( ; i != connectedIDs.end(); ++i )
482 shortEdges.push_back( allEdges[ *i ]);
487 //================================================================================
489 * \brief Find EDGEs to discretize using projection from MA
490 * \param [in,out] theSinuFace - the FACE to be meshed
491 * \return bool - OK or not
493 * It separates all EDGEs into four sides of a quadrangle connected in the order:
494 * theSinuEdges[0], theShortEdges[0], theSinuEdges[1], theShortEdges[1]
496 //================================================================================
498 bool getSinuousEdges( SMESH_MesherHelper& theHelper,
499 SinuousFace& theSinuFace)
501 vector<TopoDS_Edge> * theSinuEdges = & theSinuFace._sinuSide [0];
502 vector<TopoDS_Edge> * theShortEdges = & theSinuFace._shortSide[0];
503 theSinuEdges[0].clear();
504 theSinuEdges[1].clear();
505 theShortEdges[0].clear();
506 theShortEdges[1].clear();
508 vector<TopoDS_Edge> & allEdges = theSinuFace._edges;
509 const size_t nbEdges = allEdges.size();
510 if ( nbEdges < 4 && theSinuFace._nbWires == 1 )
513 if ( theSinuFace._nbWires == 2 ) // ring
515 size_t nbOutEdges = theSinuFace._nbEdgesInWire.front();
516 theSinuEdges[0].assign ( allEdges.begin(), allEdges.begin() + nbOutEdges );
517 theSinuEdges[1].assign ( allEdges.begin() + nbOutEdges, allEdges.end() );
518 theSinuFace._sinuEdges = allEdges;
521 if ( theSinuFace._nbWires > 2 )
524 // create MedialAxis to find short edges by analyzing MA branches
525 double minSegLen = getMinSegLen( theHelper, allEdges );
526 SMESH_MAT2d::MedialAxis ma( theSinuFace.Face(), allEdges, minSegLen * 3 );
528 // in an initial request case, theFace represents a part of a river with almost parallel banks
529 // so there should be two branch points
530 using SMESH_MAT2d::BranchEnd;
531 using SMESH_MAT2d::Branch;
532 const vector< const BranchEnd* >& braPoints = ma.getBranchPoints();
533 if ( braPoints.size() < 2 )
535 TopTools_MapOfShape shortMap;
536 size_t nbBranchPoints = 0;
537 for ( size_t i = 0; i < braPoints.size(); ++i )
539 vector< const Branch* > vertBranches; // branches with an end on VERTEX
540 for ( size_t ib = 0; ib < braPoints[i]->_branches.size(); ++ib )
542 const Branch* branch = braPoints[i]->_branches[ ib ];
543 if ( branch->hasEndOfType( SMESH_MAT2d::BE_ON_VERTEX ))
544 vertBranches.push_back( branch );
546 if ( vertBranches.size() != 2 || braPoints[i]->_branches.size() != 3)
549 // get common EDGEs of two branches
550 if ( !getConnectedEdges( vertBranches[0], vertBranches[1],
551 allEdges, theShortEdges[ nbBranchPoints > 0 ] ))
554 for ( size_t iS = 0; iS < theShortEdges[ nbBranchPoints ].size(); ++iS )
555 shortMap.Add( theShortEdges[ nbBranchPoints ][ iS ]);
560 if ( nbBranchPoints != 2 )
563 // add to theSinuEdges all edges that are not theShortEdges
564 vector< vector<TopoDS_Edge> > sinuEdges(1);
565 TopoDS_Vertex vCommon;
566 for ( size_t i = 0; i < allEdges.size(); ++i )
568 if ( !shortMap.Contains( allEdges[i] ))
570 if ( !sinuEdges.back().empty() )
571 if ( !TopExp::CommonVertex( sinuEdges.back().back(), allEdges[ i ], vCommon ))
572 sinuEdges.resize( sinuEdges.size() + 1 );
574 sinuEdges.back().push_back( allEdges[i] );
577 if ( sinuEdges.size() == 3 )
579 if ( !TopExp::CommonVertex( sinuEdges.back().back(), sinuEdges[0][0], vCommon ))
581 vector<TopoDS_Edge>& last = sinuEdges.back();
582 last.insert( last.end(), sinuEdges[0].begin(), sinuEdges[0].end() );
583 sinuEdges[0].swap( last );
584 sinuEdges.resize( 2 );
586 if ( sinuEdges.size() != 2 )
589 theSinuEdges[0].swap( sinuEdges[0] );
590 theSinuEdges[1].swap( sinuEdges[1] );
592 if ( !TopExp::CommonVertex( theSinuEdges[0].back(), theShortEdges[0][0], vCommon ) ||
593 !vCommon.IsSame( theHelper.IthVertex( 1, theSinuEdges[0].back() )))
594 theShortEdges[0].swap( theShortEdges[1] );
596 theSinuFace._sinuEdges = theSinuEdges[0];
597 theSinuFace._sinuEdges.insert( theSinuFace._sinuEdges.end(),
598 theSinuEdges[1].begin(), theSinuEdges[1].end() );
600 return ( theShortEdges[0].size() > 0 && theShortEdges[1].size() > 0 &&
601 theSinuEdges [0].size() > 0 && theSinuEdges [1].size() > 0 );
603 // the sinuous EDGEs can be composite and C0 continuous,
604 // therefor we use a complex criterion to find TWO short non-sinuous EDGEs
605 // and the rest EDGEs will be treated as sinuous.
606 // A short edge should have the following features:
609 // c) with convex corners at ends
610 // d) far from the other short EDGE
612 // vector< double > isStraightEdge( nbEdges, 0 ); // criterion value
614 // // a0) evaluate continuity
615 // const double contiWgt = 0.5; // weight of continuity in the criterion
616 // multimap< int, TopoDS_Edge > continuity;
617 // for ( size_t i = 0; i < nbEdges; ++I )
619 // BRepAdaptor_Curve curve( allEdges[i] );
620 // GeomAbs_Shape C = GeomAbs_CN;
622 // C = curve.Continuity(); // C0, G1, C1, G2, C2, C3, CN
623 // catch ( Standard_Failure ) {}
624 // continuity.insert( make_pair( C, allEdges[i] ));
625 // isStraight[i] += double( C ) / double( CN ) * contiWgt;
628 // // try to choose by continuity
629 // int mostStraight = (int) continuity.rbegin()->first;
630 // int lessStraight = (int) continuity.begin()->first;
631 // if ( mostStraight != lessStraight )
633 // int nbStraight = continuity.count( mostStraight );
634 // if ( nbStraight == 2 )
636 // getTwo( /*least=*/false, continuity, theShortEdges, theSinuEdges );
638 // else if ( nbStraight == 3 && nbEdges == 4 )
640 // theSinuEdges.push_back( continuity.begin()->second );
641 // vector<TopoDS_Edge>::iterator it =
642 // std::find( allEdges.begin(), allEdges.end(), theSinuEdges[0] );
643 // int i = std::distance( allEdges.begin(), it );
644 // theSinuEdges .push_back( allEdges[( i+2 )%4 ]);
645 // theShortEdges.push_back( allEdges[( i+1 )%4 ]);
646 // theShortEdges.push_back( allEdges[( i+3 )%4 ]);
648 // if ( theShortEdges.size() == 2 )
652 // // a) curvature; evaluate aspect ratio
654 // const double curvWgt = 0.5;
655 // for ( size_t i = 0; i < nbEdges; ++I )
657 // BRepAdaptor_Curve curve( allEdges[i] );
658 // double curvature = 1;
659 // if ( !curve.IsClosed() )
661 // const double f = curve.FirstParameter(), l = curve.LastParameter();
662 // gp_Pnt pf = curve.Value( f ), pl = curve.Value( l );
663 // gp_Lin line( pf, pl.XYZ() - pf.XYZ() );
664 // double distMax = 0;
665 // for ( double u = f; u < l; u += (l-f)/30. )
666 // distMax = Max( distMax, line.SquareDistance( curve.Value( u )));
667 // curvature = Sqrt( distMax ) / ( pf.Distance( pl ));
669 // isStraight[i] += curvWgt / ( curvature + 1e-20 );
674 // const double lenWgt = 0.5;
675 // for ( size_t i = 0; i < nbEdges; ++I )
677 // double length = SMESH_Algo::Length( allEdges[i] );
679 // isStraight[i] += lenWgt / length;
682 // // c) with convex corners at ends
684 // const double cornerWgt = 0.25;
685 // for ( size_t i = 0; i < nbEdges; ++I )
687 // double convex = 0;
688 // int iPrev = SMESH_MesherHelper::WrapIndex( int(i)-1, nbEdges );
689 // int iNext = SMESH_MesherHelper::WrapIndex( int(i)+1, nbEdges );
690 // TopoDS_Vertex v = helper.IthVertex( 0, allEdges[i] );
691 // double angle = SMESH_MesherHelper::GetAngle( allEdges[iPrev], allEdges[i], theFace, v );
692 // if ( angle < M_PI ) // [-PI; PI]
693 // convex += ( angle + M_PI ) / M_PI / M_PI;
694 // v = helper.IthVertex( 1, allEdges[i] );
695 // angle = SMESH_MesherHelper::GetAngle( allEdges[iNext], allEdges[i], theFace, v );
696 // if ( angle < M_PI ) // [-PI; PI]
697 // convex += ( angle + M_PI ) / M_PI / M_PI;
698 // isStraight[i] += cornerWgt * convex;
703 //================================================================================
705 * \brief Creates an EDGE from a sole branch of MA
707 //================================================================================
709 TopoDS_Edge makeEdgeFromMA( SMESH_MesherHelper& theHelper,
710 const SMESH_MAT2d::MedialAxis& theMA,
711 const double theMinSegLen)
713 if ( theMA.nbBranches() != 1 )
714 return TopoDS_Edge();
717 theMA.getPoints( theMA.getBranch(0), uv );
719 return TopoDS_Edge();
721 TopoDS_Face face = TopoDS::Face( theHelper.GetSubShape() );
722 Handle(Geom_Surface) surface = BRep_Tool::Surface( face );
724 vector< gp_Pnt > pnt;
725 pnt.reserve( uv.size() * 2 );
726 pnt.push_back( surface->Value( uv[0].X(), uv[0].Y() ));
727 for ( size_t i = 1; i < uv.size(); ++i )
729 gp_Pnt p = surface->Value( uv[i].X(), uv[i].Y() );
730 int nbDiv = int( p.Distance( pnt.back() ) / theMinSegLen );
731 for ( int iD = 1; iD < nbDiv; ++iD )
733 double R = iD / double( nbDiv );
734 gp_XY uvR = uv[i-1] * (1 - R) + uv[i] * R;
735 pnt.push_back( surface->Value( uvR.X(), uvR.Y() ));
740 // cout << "from salome.geom import geomBuilder" << endl;
741 // cout << "geompy = geomBuilder.New(salome.myStudy)" << endl;
742 Handle(TColgp_HArray1OfPnt) points = new TColgp_HArray1OfPnt(1, pnt.size());
743 for ( size_t i = 0; i < pnt.size(); ++i )
746 points->SetValue( i+1, p );
747 // cout << "geompy.MakeVertex( "<< p.X()<<", " << p.Y()<<", " << p.Z()
748 // <<" theName = 'p_" << i << "')" << endl;
751 GeomAPI_Interpolate interpol( points, /*isClosed=*/false, gp::Resolution());
753 if ( !interpol.IsDone())
754 return TopoDS_Edge();
756 TopoDS_Edge branchEdge = BRepBuilderAPI_MakeEdge(interpol.Curve());
760 //================================================================================
762 * \brief Returns a type of shape, to which a hypothesis used to mesh a given edge is assigned
764 //================================================================================
766 TopAbs_ShapeEnum getHypShape( SMESH_Mesh* mesh, const TopoDS_Shape& edge )
768 TopAbs_ShapeEnum shapeType = TopAbs_SHAPE;
770 SMESH_subMesh* sm = mesh->GetSubMesh( edge );
771 SMESH_Algo* algo = sm->GetAlgo();
772 if ( !algo ) return shapeType;
774 const list <const SMESHDS_Hypothesis *> & hyps =
775 algo->GetUsedHypothesis( *mesh, edge, /*ignoreAuxiliary=*/true );
776 if ( hyps.empty() ) return shapeType;
778 TopoDS_Shape shapeOfHyp =
779 SMESH_MesherHelper::GetShapeOfHypothesis( hyps.front(), edge, mesh);
781 return SMESH_MesherHelper::GetGroupType( shapeOfHyp, /*woCompound=*/true);
784 //================================================================================
786 * \brief Discretize a sole branch of MA an returns parameters of divisions on MA
788 //================================================================================
790 bool divideMA( SMESH_MesherHelper& theHelper,
791 const SMESH_MAT2d::MedialAxis& theMA,
792 const SinuousFace& theSinuFace,
793 SMESH_Algo* the1dAlgo,
794 const double theMinSegLen,
795 vector<double>& theMAParams )
797 // Check if all EDGEs of one size are meshed, then MA discretization is not needed
798 SMESH_Mesh* mesh = theHelper.GetMesh();
799 size_t nbComputedEdges[2] = { 0, 0 };
800 for ( size_t iS = 0; iS < 2; ++iS )
801 for ( size_t i = 0; i < theSinuFace._sinuSide[iS].size(); ++i )
803 const TopoDS_Edge& sinuEdge = theSinuFace._sinuSide[iS][i];
804 SMESH_subMesh* sm = mesh->GetSubMesh( sinuEdge );
805 bool isComputed = ( !sm->IsEmpty() );
808 TopAbs_ShapeEnum shape = getHypShape( mesh, sinuEdge );
809 if ( shape == TopAbs_SHAPE || shape <= TopAbs_FACE )
811 // EDGE computed using global hypothesis -> clear it
812 bool hasComputedFace = false;
813 PShapeIteratorPtr faceIt = theHelper.GetAncestors( sinuEdge, *mesh, TopAbs_FACE );
814 while ( const TopoDS_Shape* face = faceIt->next() )
815 if (( !face->IsSame( theSinuFace.Face() )) &&
816 ( hasComputedFace = !mesh->GetSubMesh( *face )->IsEmpty() ))
818 if ( !hasComputedFace )
820 sm->ComputeStateEngine( SMESH_subMesh::CLEAN );
825 nbComputedEdges[ iS ] += isComputed;
827 if ( nbComputedEdges[0] == theSinuFace._sinuSide[0].size() ||
828 nbComputedEdges[1] == theSinuFace._sinuSide[1].size() )
829 return true; // discretization is not needed
832 TopoDS_Edge branchEdge = makeEdgeFromMA( theHelper, theMA, theMinSegLen );
833 if ( branchEdge.IsNull() )
836 // const char* file = "/misc/dn25/salome/eap/salome/misc/tmp/MAedge.brep";
837 // BRepTools::Write( branchEdge, file);
838 // cout << "Write " << file << endl;
841 // Find 1D algo to mesh branchEdge
843 // look for a most local 1D hyp assigned to the FACE
844 int mostSimpleShape = -1, maxShape = TopAbs_EDGE;
846 for ( size_t i = 0; i < theSinuFace._sinuEdges.size(); ++i )
848 TopAbs_ShapeEnum shapeType = getHypShape( mesh, theSinuFace._sinuEdges[i] );
849 if ( mostSimpleShape < shapeType && shapeType < maxShape )
851 edge = theSinuFace._sinuEdges[i];
852 mostSimpleShape = shapeType;
856 SMESH_Algo* algo = the1dAlgo;
857 if ( mostSimpleShape > -1 )
859 algo = mesh->GetSubMesh( edge )->GetAlgo();
860 SMESH_Hypothesis::Hypothesis_Status status;
861 if ( !algo->CheckHypothesis( *mesh, edge, status ))
866 tmpMesh.ShapeToMesh( branchEdge );
868 mesh->GetGen()->Compute( tmpMesh, branchEdge, true, true ); // make nodes on VERTEXes
869 if ( !algo->Compute( tmpMesh, branchEdge ))
875 return SMESH_Algo::GetNodeParamOnEdge( tmpMesh.GetMeshDS(), branchEdge, theMAParams );
878 //================================================================================
880 * \brief Select division parameters on MA and make them coincide at ends with
881 * projections of VERTEXes to MA for a given pair of opposite EDGEs
882 * \param [in] theEdgePairInd - index of the EDGE pair
883 * \param [in] theDivPoints - the BranchPoint's dividing MA into parts each
884 * corresponding to a unique pair of opposite EDGEs
885 * \param [in] theMAParams - the MA division parameters
886 * \param [out] theSelectedMAParams - the selected MA parameters
887 * \return bool - is OK
889 //================================================================================
891 bool getParamsForEdgePair( const size_t theEdgePairInd,
892 const vector< SMESH_MAT2d::BranchPoint >& theDivPoints,
893 const vector<double>& theMAParams,
894 vector<double>& theSelectedMAParams)
896 if ( theDivPoints.empty() )
898 theSelectedMAParams = theMAParams;
901 if ( theEdgePairInd > theDivPoints.size() || theMAParams.empty() )
904 // find a range of params to copy
908 if ( theEdgePairInd > 0 )
910 const SMESH_MAT2d::BranchPoint& bp = theDivPoints[ theEdgePairInd-1 ];
911 bp._branch->getParameter( bp, par1 );
912 while ( theMAParams[ iPar1 ] < par1 ) ++iPar1;
913 if ( par1 - theMAParams[ iPar1-1 ] < theMAParams[ iPar1 ] - par1 )
918 size_t iPar2 = theMAParams.size() - 1;
919 if ( theEdgePairInd < theDivPoints.size() )
921 const SMESH_MAT2d::BranchPoint& bp = theDivPoints[ theEdgePairInd ];
922 bp._branch->getParameter( bp, par2 );
924 while ( theMAParams[ iPar2 ] < par2 ) ++iPar2;
925 if ( par2 - theMAParams[ iPar2-1 ] < theMAParams[ iPar2 ] - par2 )
929 theSelectedMAParams.assign( theMAParams.begin() + iPar1,
930 theMAParams.begin() + iPar2 + 1 );
932 // adjust theSelectedMAParams to fit between par1 and par2
934 double d = par1 - theSelectedMAParams[0];
935 double f = ( par2 - par1 ) / ( theSelectedMAParams.back() - theSelectedMAParams[0] );
937 for ( size_t i = 0; i < theSelectedMAParams.size(); ++i )
939 theSelectedMAParams[i] += d;
940 theSelectedMAParams[i] = par1 + ( theSelectedMAParams[i] - par1 ) * f;
946 //--------------------------------------------------------------------------------
947 // node or node parameter on EDGE
950 const SMDS_MeshNode* _node;
952 int _edgeInd; // index in theSinuEdges vector
954 NodePoint(): _node(0), _u(0), _edgeInd(-1) {}
955 NodePoint(const SMDS_MeshNode* n, double u, size_t iEdge ): _node(n), _u(u), _edgeInd(iEdge) {}
956 NodePoint(double u, size_t iEdge) : _node(0), _u(u), _edgeInd(iEdge) {}
957 NodePoint(const SMESH_MAT2d::BoundaryPoint& p) : _node(0), _u(p._param), _edgeInd(p._edgeIndex) {}
958 gp_Pnt Point(const vector< Handle(Geom_Curve) >& curves) const
960 return _node ? SMESH_TNodeXYZ(_node) : curves[ _edgeInd ]->Value( _u );
963 typedef multimap< double, pair< NodePoint, NodePoint > > TMAPar2NPoints;
965 //================================================================================
967 * \brief Finds a VERTEX corresponding to a point on EDGE, which is also filled
968 * with a node on the VERTEX, present or created
969 * \param [in,out] theNodePnt - the node position on the EDGE
970 * \param [in] theSinuEdges - the sinuous EDGEs
971 * \param [in] theMeshDS - the mesh
972 * \return bool - true if the \a theBndPnt is on VERTEX
974 //================================================================================
976 bool findVertexAndNode( NodePoint& theNodePnt,
977 const vector<TopoDS_Edge>& theSinuEdges,
978 SMESHDS_Mesh* theMeshDS = 0,
979 size_t theEdgeIndPrev = 0,
980 size_t theEdgeIndNext = 0)
982 if ( theNodePnt._edgeInd >= theSinuEdges.size() )
986 BRep_Tool::Range( theSinuEdges[ theNodePnt._edgeInd ], f,l );
987 const double tol = 1e-3 * ( l - f );
990 if ( Abs( f - theNodePnt._u ) < tol )
991 V = SMESH_MesherHelper::IthVertex( 0, theSinuEdges[ theNodePnt._edgeInd ], /*CumOri=*/false);
992 else if ( Abs( l - theNodePnt._u ) < tol )
993 V = SMESH_MesherHelper::IthVertex( 1, theSinuEdges[ theNodePnt._edgeInd ], /*CumOri=*/false);
994 else if ( theEdgeIndPrev != theEdgeIndNext )
995 TopExp::CommonVertex( theSinuEdges[theEdgeIndPrev], theSinuEdges[theEdgeIndNext], V );
997 if ( !V.IsNull() && theMeshDS )
999 theNodePnt._node = SMESH_Algo::VertexNode( V, theMeshDS );
1000 if ( !theNodePnt._node )
1002 gp_Pnt p = BRep_Tool::Pnt( V );
1003 theNodePnt._node = theMeshDS->AddNode( p.X(), p.Y(), p.Z() );
1004 theMeshDS->SetNodeOnVertex( theNodePnt._node, V );
1010 //================================================================================
1012 * \brief Add to the map of NodePoint's those on VERTEXes
1013 * \param [in,out] theHelper - the helper
1014 * \param [in] theMA - Medial Axis
1015 * \param [in] theMinSegLen - minimal segment length
1016 * \param [in] theDivPoints - projections of VERTEXes to MA
1017 * \param [in] theSinuEdges - the sinuous EDGEs
1018 * \param [in] theSideEdgeIDs - indices of sinuous EDGEs per side
1019 * \param [in] theIsEdgeComputed - is sinuous EGDE is meshed
1020 * \param [in,out] thePointsOnE - the map to fill
1021 * \param [out] theNodes2Merge - the map of nodes to merge
1023 //================================================================================
1025 bool projectVertices( SMESH_MesherHelper& theHelper,
1026 const SMESH_MAT2d::MedialAxis& theMA,
1027 const vector< SMESH_MAT2d::BranchPoint >& theDivPoints,
1028 const vector< std::size_t > & theEdgeIDs1,
1029 const vector< std::size_t > & theEdgeIDs2,
1030 const vector< bool >& theIsEdgeComputed,
1031 TMAPar2NPoints & thePointsOnE,
1032 SinuousFace& theSinuFace)
1034 SMESHDS_Mesh* meshDS = theHelper.GetMeshDS();
1035 const vector<TopoDS_Edge>& theSinuEdges = theSinuFace._sinuEdges;
1036 const vector< Handle(Geom_Curve) >& theCurves = theSinuFace._sinuCurves;
1039 SMESH_MAT2d::BoundaryPoint bp[2];
1040 const SMESH_MAT2d::Branch& branch = *theMA.getBranch(0);
1042 // add to thePointsOnE NodePoint's of ends of theSinuEdges
1043 if ( !branch.getBoundaryPoints( 0., bp[0], bp[1] ) ||
1044 !theMA.getBoundary().moveToClosestEdgeEnd( bp[0] )) return false;
1045 if ( !theSinuFace.IsRing() &&
1046 !theMA.getBoundary().moveToClosestEdgeEnd( bp[1] )) return false;
1047 NodePoint np0( bp[0] ), np1( bp[1] );
1048 findVertexAndNode( np0, theSinuEdges, meshDS );
1049 findVertexAndNode( np1, theSinuEdges, meshDS );
1050 thePointsOnE.insert( make_pair( -0.1, make_pair( np0, np1 )));
1052 if ( !theSinuFace.IsRing() )
1054 if ( !branch.getBoundaryPoints( 1., bp[0], bp[1] ) ||
1055 !theMA.getBoundary().moveToClosestEdgeEnd( bp[0] ) ||
1056 !theMA.getBoundary().moveToClosestEdgeEnd( bp[1] )) return false;
1057 np0 = bp[0]; np1 = bp[1];
1058 findVertexAndNode( np0, theSinuEdges, meshDS );
1059 findVertexAndNode( np1, theSinuEdges, meshDS );
1060 thePointsOnE.insert( make_pair( 1.1, make_pair( np0, np1)));
1062 // project theDivPoints
1064 if ( theDivPoints.empty() )
1067 for ( size_t i = 0; i < theDivPoints.size(); ++i )
1069 if ( !branch.getParameter( theDivPoints[i], uMA ))
1071 if ( !branch.getBoundaryPoints( theDivPoints[i], bp[0], bp[1] ))
1078 bool isVertex[2] = {
1079 findVertexAndNode( np[0], theSinuEdges, meshDS, theEdgeIDs1[i], theEdgeIDs1[i+1] ),
1080 findVertexAndNode( np[1], theSinuEdges, meshDS, theEdgeIDs2[i], theEdgeIDs2[i+1] )
1083 TMAPar2NPoints::iterator u2NP =
1084 thePointsOnE.insert( make_pair( uMA, make_pair( np[0], np[1])));//.first;
1086 if ( !isVertex[0] && !isVertex[1] ) return false; // error
1087 if ( isVertex[0] && isVertex[1] )
1089 const size_t iVert = isVertex[0] ? 0 : 1;
1090 const size_t iNode = 1 - iVert;
1092 bool isOppComputed = theIsEdgeComputed[ np[ iNode ]._edgeInd ];
1093 if ( !isOppComputed )
1096 // a VERTEX is projected on a meshed EDGE; there are two options:
1097 // 1) a projected point is joined with a closet node if a strip between this and neighbor
1098 // projection is WIDE enough; joining is done by creating a node coincident with the
1099 // existing node which will be merged together after all;
1100 // 2) a neighbor projection is merged with this one if it is TOO CLOSE; a node of deleted
1101 // projection is set to the BoundaryPoint of this projection
1103 // evaluate distance to neighbor projections
1104 const double rShort = 0.2;
1105 bool isShortPrev[2], isShortNext[2];
1106 TMAPar2NPoints::iterator u2NPPrev = u2NP, u2NPNext = u2NP;
1107 --u2NPPrev; ++u2NPNext;
1108 // bool hasPrev = ( u2NP != thePointsOnE.begin() );
1109 // bool hasNext = ( u2NPNext != thePointsOnE.end() );
1110 // if ( !hasPrev ) u2NPPrev = u2NP0;
1111 // if ( !hasNext ) u2NPNext = u2NP1;
1112 for ( int iS = 0; iS < 2; ++iS ) // side with Vertex and side with Nodes
1114 NodePoint np = get( u2NP->second, iS );
1115 NodePoint npPrev = get( u2NPPrev->second, iS );
1116 NodePoint npNext = get( u2NPNext->second, iS );
1117 gp_Pnt p = np .Point( theCurves );
1118 gp_Pnt pPrev = npPrev.Point( theCurves );
1119 gp_Pnt pNext = npNext.Point( theCurves );
1120 double distPrev = p.Distance( pPrev );
1121 double distNext = p.Distance( pNext );
1122 double r = distPrev / ( distPrev + distNext );
1123 isShortPrev[iS] = ( r < rShort );
1124 isShortNext[iS] = (( 1 - r ) > ( 1 - rShort ));
1126 // if ( !hasPrev ) isShortPrev[0] = isShortPrev[1] = false;
1127 // if ( !hasNext ) isShortNext[0] = isShortNext[1] = false;
1129 TMAPar2NPoints::iterator u2NPClose;
1131 if (( isShortPrev[0] && isShortPrev[1] ) || // option 2) -> remove a too close projection
1132 ( isShortNext[0] && isShortNext[1] ))
1134 u2NPClose = isShortPrev[0] ? u2NPPrev : u2NPNext;
1135 NodePoint& npProj = get( u2NP->second, iNode ); // NP of VERTEX projection
1136 NodePoint npCloseN = get( u2NPClose->second, iNode ); // NP close to npProj
1137 NodePoint npCloseV = get( u2NPClose->second, iVert ); // NP close to VERTEX
1138 if ( !npCloseV._node )
1141 thePointsOnE.erase( isShortPrev[0] ? u2NPPrev : u2NPNext );
1146 // can't remove the neighbor projection as it is also from VERTEX, -> option 1)
1149 // else: option 1) - wide enough -> "duplicate" existing node
1151 u2NPClose = isShortPrev[ iNode ] ? u2NPPrev : u2NPNext;
1152 NodePoint& npProj = get( u2NP->second, iNode ); // NP of VERTEX projection
1153 NodePoint& npCloseN = get( u2NPClose->second, iNode ); // NP close to npProj
1156 //npProj._edgeInd = npCloseN._edgeInd;
1157 // npProj._u = npCloseN._u + 1e-3 * Abs( get( u2NPPrev->second, iNode )._u -
1158 // get( u2NPNext->second, iNode )._u );
1159 // gp_Pnt p = npProj.Point( theCurves );
1160 // npProj._node = meshDS->AddNode( p.X(), p.Y(), p.Z() );
1161 // meshDS->SetNodeOnEdge( npProj._node, theSinuEdges[ npProj._edgeInd ], npProj._u );
1163 //theNodes2Merge[ npCloseN._node ].push_back( npProj._node );
1169 double getUOnEdgeByPoint( const size_t iEdge,
1170 const NodePoint* point,
1171 SinuousFace& sinuFace )
1173 if ( point->_edgeInd == iEdge )
1176 TopoDS_Vertex V0 = TopExp::FirstVertex( sinuFace._sinuEdges[ iEdge ]);
1177 TopoDS_Vertex V1 = TopExp::LastVertex ( sinuFace._sinuEdges[ iEdge ]);
1178 gp_Pnt p0 = BRep_Tool::Pnt( V0 );
1179 gp_Pnt p1 = BRep_Tool::Pnt( V1 );
1180 gp_Pnt p = point->Point( sinuFace._sinuCurves );
1183 BRep_Tool::Range( sinuFace._sinuEdges[ iEdge ], f,l );
1184 return p.SquareDistance( p0 ) < p.SquareDistance( p1 ) ? f : l;
1187 //================================================================================
1189 * \brief Move coincident nodes to make node params on EDGE unique
1190 * \param [in] theHelper - the helper
1191 * \param [in] thePointsOnE - nodes on two opposite river sides
1192 * \param [in] theSinuFace - the sinuous FACE
1193 * \param [out] theNodes2Merge - the map of nodes to merge
1195 //================================================================================
1197 void separateNodes( SMESH_MesherHelper& theHelper,
1198 const SMESH_MAT2d::MedialAxis& theMA,
1199 TMAPar2NPoints & thePointsOnE,
1200 SinuousFace& theSinuFace )
1202 if ( thePointsOnE.size() < 2 )
1205 SMESHDS_Mesh* meshDS = theHelper.GetMeshDS();
1206 const vector<TopoDS_Edge>& theSinuEdges = theSinuFace._sinuEdges;
1207 const vector< Handle(Geom_Curve) >& curves = theSinuFace._sinuCurves;
1209 SMESH_MAT2d::BoundaryPoint bp[2];
1210 const SMESH_MAT2d::Branch& branch = *theMA.getBranch(0);
1212 typedef TMAPar2NPoints::iterator TIterator;
1214 for ( int iSide = 0; iSide < 2; ++iSide ) // loop on two sinuous sides
1216 // get a tolerance to compare points
1217 double tol = Precision::Confusion();
1218 for ( size_t i = 0; i < theSinuFace._sinuSide[ iSide ].size(); ++i )
1219 tol = Max( tol , BRep_Tool::Tolerance( theSinuFace._sinuSide[ iSide ][ i ]));
1221 // find coincident points
1222 TIterator u2NP = thePointsOnE.begin();
1223 vector< TIterator > sameU2NP( 1, u2NP++ );
1224 while ( u2NP != thePointsOnE.end() )
1226 for ( ; u2NP != thePointsOnE.end(); ++u2NP )
1228 NodePoint& np1 = get( sameU2NP.back()->second, iSide );
1229 NodePoint& np2 = get( u2NP ->second, iSide );
1231 if (( !np1._node || !np2._node ) &&
1232 ( np1.Point( curves ).SquareDistance( np2.Point( curves )) < tol*tol ))
1234 sameU2NP.push_back( u2NP );
1236 else if ( sameU2NP.size() == 1 )
1238 sameU2NP[ 0 ] = u2NP;
1246 if ( sameU2NP.size() > 1 )
1248 // find an existing node on VERTEX among sameU2NP and get underlying EDGEs
1249 const SMDS_MeshNode* existingNode = 0;
1250 set< int > edgeInds;
1252 for ( size_t i = 0; i < sameU2NP.size(); ++i )
1254 np = &get( sameU2NP[i]->second, iSide );
1256 if ( !existingNode || np->_node->GetPosition()->GetDim() == 0 )
1257 existingNode = np->_node;
1258 edgeInds.insert( np->_edgeInd );
1260 list< const SMDS_MeshNode* >& mergeNodes = theSinuFace._nodesToMerge[ existingNode ];
1262 TIterator u2NPprev = sameU2NP.front();
1263 TIterator u2NPnext = sameU2NP.back() ;
1264 if ( u2NPprev->first > 0. ) --u2NPprev;
1265 if ( u2NPnext->first < 1. ) ++u2NPprev;
1267 set< int >::iterator edgeID = edgeInds.begin();
1268 for ( ; edgeID != edgeInds.end(); ++edgeID )
1270 // get U range on iEdge within which the equal points will be distributed
1272 np = &get( u2NPprev->second, iSide );
1273 u0 = getUOnEdgeByPoint( *edgeID, np, theSinuFace );
1275 np = &get( u2NPnext->second, iSide );
1276 u1 = getUOnEdgeByPoint( *edgeID, np, theSinuFace );
1278 // distribute points and create nodes
1279 double du = ( u1 - u0 ) / ( sameU2NP.size() + 1 );
1281 for ( size_t i = 0; i < sameU2NP.size(); ++i )
1283 np = &get( sameU2NP[i]->second, iSide );
1284 if ( !np->_node && *edgeID == np->_edgeInd )
1288 gp_Pnt p = np->Point( curves );
1289 np->_node = meshDS->AddNode( p.X(), p.Y(), p.Z() );
1290 meshDS->SetNodeOnEdge( np->_node, theSinuEdges[ *edgeID ], np->_u );
1291 //mergeNodes.push_back( np->_node );
1296 sameU2NP.resize( 1 );
1297 u2NP = ++sameU2NP.back();
1298 sameU2NP[ 0 ] = u2NP;
1300 } // if ( sameU2NP.size() > 1 )
1301 } // while ( u2NP != thePointsOnE.end() )
1302 } // for ( int iSide = 0; iSide < 2; ++iSide )
1305 } // separateNodes()
1307 //================================================================================
1309 * \brief Setup sides of SinuousFace::_quad
1310 * \param [in] theHelper - helper
1311 * \param [in] thePointsOnEdges - NodePoint's on sinuous sides
1312 * \param [in,out] theSinuFace - the FACE
1313 * \param [in] the1dAlgo - algorithm to use for radial discretization of a ring FACE
1314 * \return bool - is OK
1316 //================================================================================
1318 bool setQuadSides(SMESH_MesherHelper& theHelper,
1319 const TMAPar2NPoints& thePointsOnEdges,
1320 SinuousFace& theFace,
1321 SMESH_Algo* the1dAlgo)
1323 SMESH_Mesh* mesh = theHelper.GetMesh();
1324 const TopoDS_Face& face = theFace._quad->face;
1325 SMESH_ProxyMesh::Ptr proxyMesh = StdMeshers_ViscousLayers2D::Compute( *mesh, face );
1329 list< TopoDS_Edge > side[4];
1330 side[0].insert( side[0].end(), theFace._shortSide[0].begin(), theFace._shortSide[0].end() );
1331 side[1].insert( side[1].end(), theFace._sinuSide[1].begin(), theFace._sinuSide[1].end() );
1332 side[2].insert( side[2].end(), theFace._shortSide[1].begin(), theFace._shortSide[1].end() );
1333 side[3].insert( side[3].end(), theFace._sinuSide[0].begin(), theFace._sinuSide[0].end() );
1335 for ( int i = 0; i < 4; ++i )
1337 theFace._quad->side[i] = StdMeshers_FaceSide::New( face, side[i], mesh, i < QUAD_TOP_SIDE,
1338 /*skipMediumNodes=*/true, proxyMesh );
1341 if ( theFace.IsRing() )
1343 // --------------------------------------
1344 // Discretize a ring in radial direction
1345 // --------------------------------------
1347 if ( thePointsOnEdges.size() < 4 )
1350 // find most distant opposite nodes
1351 double maxDist = 0, dist;
1352 TMAPar2NPoints::const_iterator u2NPdist, u2NP = thePointsOnEdges.begin();
1353 for ( ; u2NP != thePointsOnEdges.end(); ++u2NP )
1355 SMESH_TNodeXYZ xyz( u2NP->second.first._node );
1356 dist = xyz.SquareDistance( u2NP->second.second._node );
1357 if ( dist > maxDist )
1363 // compute distribution of radial nodes
1364 list< double > params; // normalized params
1365 static_cast< StdMeshers_QuadFromMedialAxis_1D2D::Algo1D* >
1366 ( the1dAlgo )->ComputeDistribution( theHelper,
1367 SMESH_TNodeXYZ( u2NPdist->second.first._node ),
1368 SMESH_TNodeXYZ( u2NPdist->second.second._node ),
1371 // add a radial quad side
1372 u2NP = thePointsOnEdges.begin();
1373 const SMDS_MeshNode* nOut = u2NP->second.first._node;
1374 const SMDS_MeshNode* nIn = u2NP->second.second._node;
1375 nOut = proxyMesh->GetProxyNode( nOut );
1376 nIn = proxyMesh->GetProxyNode( nIn );
1377 gp_XY uvOut = theHelper.GetNodeUV( face, nOut );
1378 gp_XY uvIn = theHelper.GetNodeUV( face, nIn );
1379 Handle(Geom_Surface) surface = BRep_Tool::Surface( face );
1380 UVPtStructVec uvsNew; UVPtStruct uvPt;
1384 uvsNew.push_back( uvPt );
1385 for (list<double>::iterator itU = params.begin(); itU != params.end(); ++itU )
1387 gp_XY uv = ( 1 - *itU ) * uvOut + *itU * uvIn;
1388 gp_Pnt p = surface->Value( uv.X(), uv.Y() );
1389 uvPt.node = theHelper.AddNode( p.X(), p.Y(), p.Z(), /*id=*/0, uv.X(), uv.Y() );
1392 uvsNew.push_back( uvPt );
1397 uvsNew.push_back( uvPt );
1399 theFace._quad->side[ 0 ] = StdMeshers_FaceSide::New( uvsNew );
1400 theFace._quad->side[ 2 ] = theFace._quad->side[ 0 ];
1402 // rotate the IN side if opposite nodes of IN and OUT sides don't match
1403 const SMDS_MeshNode * nIn0 = theFace._quad->side[ 1 ].First().node;
1406 nIn = proxyMesh->GetProxyNode( nIn );
1407 const UVPtStructVec& uvsIn = theFace._quad->side[ 1 ].GetUVPtStruct(); // _sinuSide[1]
1408 size_t i; // find UVPtStruct holding nIn
1409 for ( i = 0; i < uvsIn.size(); ++i )
1410 if ( nIn == uvsIn[i].node )
1412 if ( i == uvsIn.size() )
1415 // create a new IN quad side
1417 uvsNew.reserve( uvsIn.size() );
1418 uvsNew.insert( uvsNew.end(), uvsIn.begin() + i, uvsIn.end() );
1419 uvsNew.insert( uvsNew.end(), uvsIn.begin() + 1, uvsIn.begin() + i + 1);
1420 theFace._quad->side[ 1 ] = StdMeshers_FaceSide::New( uvsNew );
1422 } // if ( theShortEdges[0].empty() )
1428 //================================================================================
1430 * \brief Divide the sinuous EDGEs by projecting the division point of Medial
1432 * \param [in] theHelper - the helper
1433 * \param [in] theMinSegLen - minimal segment length
1434 * \param [in] theMA - the Medial Axis
1435 * \param [in] theMAParams - parameters of division points of \a theMA
1436 * \param [in] theSinuEdges - the EDGEs to make nodes on
1437 * \param [in] theSinuSide0Size - the number of EDGEs in the 1st sinuous side
1438 * \param [in] the1dAlgo - algorithm to use for radial discretization of a ring FACE
1439 * \return bool - is OK or not
1441 //================================================================================
1443 bool computeSinuEdges( SMESH_MesherHelper& theHelper,
1444 double /*theMinSegLen*/,
1445 SMESH_MAT2d::MedialAxis& theMA,
1446 vector<double>& theMAParams,
1447 SinuousFace& theSinuFace,
1448 SMESH_Algo* the1dAlgo)
1450 if ( theMA.nbBranches() != 1 )
1453 // normalize theMAParams
1454 for ( size_t i = 0; i < theMAParams.size(); ++i )
1455 theMAParams[i] /= theMAParams.back();
1458 SMESH_Mesh* mesh = theHelper.GetMesh();
1459 SMESHDS_Mesh* meshDS = theHelper.GetMeshDS();
1462 // get data of sinuous EDGEs and remove unnecessary nodes
1463 const vector< TopoDS_Edge >& theSinuEdges = theSinuFace._sinuEdges;
1464 vector< Handle(Geom_Curve) >& curves = theSinuFace._sinuCurves;
1465 vector< int > edgeIDs ( theSinuEdges.size() ); // IDs in the main shape
1466 vector< bool > isComputed( theSinuEdges.size() );
1467 curves.resize( theSinuEdges.size(), 0 );
1468 for ( size_t i = 0; i < theSinuEdges.size(); ++i )
1470 curves[i] = BRep_Tool::Curve( theSinuEdges[i], f,l );
1473 SMESH_subMesh* sm = mesh->GetSubMesh( theSinuEdges[i] );
1474 edgeIDs [i] = sm->GetId();
1475 isComputed[i] = ( !sm->IsEmpty() );
1478 const SMESH_MAT2d::Branch& branch = *theMA.getBranch(0);
1479 SMESH_MAT2d::BoundaryPoint bp[2];
1481 vector< std::size_t > edgeIDs1, edgeIDs2; // indices in theSinuEdges
1482 vector< SMESH_MAT2d::BranchPoint > divPoints;
1483 branch.getOppositeGeomEdges( edgeIDs1, edgeIDs2, divPoints );
1484 for ( size_t i = 0; i < edgeIDs1.size(); ++i )
1485 if ( isComputed[ edgeIDs1[i]] &&
1486 isComputed[ edgeIDs2[i]] )
1488 int nbNodes1 = meshDS->MeshElements(edgeIDs[ edgeIDs1[i]] )->NbNodes();
1489 int nbNodes2 = meshDS->MeshElements(edgeIDs[ edgeIDs2[i]] )->NbNodes();
1490 if ( nbNodes1 != nbNodes2 )
1493 ( edgeIDs1[i-1] == edgeIDs1[i] ||
1494 edgeIDs2[i-1] == edgeIDs2[i] ))
1496 if (( i+1 < edgeIDs1.size() ) &&
1497 ( edgeIDs1[i+1] == edgeIDs1[i] ||
1498 edgeIDs2[i+1] == edgeIDs2[i] ))
1502 // map param on MA to parameters of nodes on a pair of theSinuEdges
1503 TMAPar2NPoints pointsOnE;
1504 vector<double> maParams;
1506 // compute params of nodes on EDGEs by projecting division points from MA
1508 for ( size_t iEdgePair = 0; iEdgePair < edgeIDs1.size(); ++iEdgePair )
1509 // loop on pairs of opposite EDGEs
1511 // --------------------------------------------------------------------------------
1512 if ( isComputed[ edgeIDs1[ iEdgePair ]] != // one EDGE is meshed
1513 isComputed[ edgeIDs2[ iEdgePair ]])
1515 // "projection" from one side to the other
1517 size_t iEdgeComputed = edgeIDs1[iEdgePair], iSideComputed = 0;
1518 if ( !isComputed[ iEdgeComputed ])
1519 ++iSideComputed, iEdgeComputed = edgeIDs2[iEdgePair];
1521 map< double, const SMDS_MeshNode* > nodeParams; // params of existing nodes
1522 if ( !SMESH_Algo::GetSortedNodesOnEdge( meshDS, theSinuEdges[ iEdgeComputed ], /*skipMedium=*/true, nodeParams ))
1525 SMESH_MAT2d::BoundaryPoint& bndPnt = bp[ 1-iSideComputed ];
1526 SMESH_MAT2d::BranchPoint brp;
1527 NodePoint npN, npB; // NodePoint's initialized by node and BoundaryPoint
1528 NodePoint& np0 = iSideComputed ? npB : npN;
1529 NodePoint& np1 = iSideComputed ? npN : npB;
1531 double maParam1st, maParamLast, maParam;
1532 if ( !theMA.getBoundary().getBranchPoint( iEdgeComputed, nodeParams.begin()->first, brp ))
1534 branch.getParameter( brp, maParam1st );
1535 if ( !theMA.getBoundary().getBranchPoint( iEdgeComputed, nodeParams.rbegin()->first, brp ))
1537 branch.getParameter( brp, maParamLast );
1539 map< double, const SMDS_MeshNode* >::iterator u2n = nodeParams.begin(), u2nEnd = nodeParams.end();
1540 TMAPar2NPoints::iterator end = pointsOnE.end(), pos = end;
1541 TMAPar2NPoints::iterator & hint = (maParamLast > maParam1st) ? end : pos;
1542 for ( ++u2n, --u2nEnd; u2n != u2nEnd; ++u2n )
1544 // point on EDGE (u2n) --> MA point (brp)
1545 if ( !theMA.getBoundary().getBranchPoint( iEdgeComputed, u2n->first, brp ))
1547 // MA point --> points on 2 EDGEs (bp)
1548 if ( !branch.getBoundaryPoints( brp, bp[0], bp[1] ) ||
1549 !branch.getParameter( brp, maParam ))
1552 npN = NodePoint( u2n->second, u2n->first, iEdgeComputed );
1553 npB = NodePoint( bndPnt );
1554 pos = pointsOnE.insert( hint, make_pair( maParam, make_pair( np0, np1 )));
1557 // move iEdgePair forward;
1558 // find divPoints most close to max MA param
1559 if ( edgeIDs1.size() > 1 )
1561 maParamLast = pointsOnE.rbegin()->first;
1563 double minDist = 1.;
1564 for ( ; iEdgePair < edgeIDs1.size()-1; ++iEdgePair )
1566 branch.getParameter( divPoints[iEdgePair], maParam );
1567 double d = Abs( maParamLast - maParam );
1569 minDist = d, iClosest = iEdgePair;
1573 if ( Abs( maParamLast - 1. ) < minDist )
1574 break; // the last pair treated
1576 iEdgePair = iClosest;
1579 // --------------------------------------------------------------------------------
1580 else if ( !isComputed[ edgeIDs1[ iEdgePair ]] && // none of EDGEs is meshed
1581 !isComputed[ edgeIDs2[ iEdgePair ]])
1583 // "projection" from MA
1585 if ( !getParamsForEdgePair( iEdgePair, divPoints, theMAParams, maParams ))
1588 for ( size_t i = 1; i < maParams.size()-1; ++i )
1590 if ( !branch.getBoundaryPoints( maParams[i], bp[0], bp[1] ))
1593 pointsOnE.insert( pointsOnE.end(), make_pair( maParams[i], make_pair( NodePoint(bp[0]),
1594 NodePoint(bp[1]))));
1597 // --------------------------------------------------------------------------------
1598 else if ( isComputed[ edgeIDs1[ iEdgePair ]] && // equally meshed EDGES
1599 isComputed[ edgeIDs2[ iEdgePair ]])
1601 // add existing nodes
1603 size_t iE0 = edgeIDs1[ iEdgePair ];
1604 size_t iE1 = edgeIDs2[ iEdgePair ];
1605 map< double, const SMDS_MeshNode* > nodeParams[2]; // params of existing nodes
1606 if ( !SMESH_Algo::GetSortedNodesOnEdge( meshDS, theSinuEdges[ iE0 ],
1607 /*skipMedium=*/false, nodeParams[0] ) ||
1608 !SMESH_Algo::GetSortedNodesOnEdge( meshDS, theSinuEdges[ iE1 ],
1609 /*skipMedium=*/false, nodeParams[1] ) ||
1610 nodeParams[0].size() != nodeParams[1].size() )
1613 if ( nodeParams[0].size() <= 2 )
1614 continue; // nodes on VERTEXes only
1616 bool reverse = ( theSinuEdges[0].Orientation() == theSinuEdges[1].Orientation() );
1618 SMESH_MAT2d::BranchPoint brp;
1619 std::pair< NodePoint, NodePoint > npPair;
1621 map< double, const SMDS_MeshNode* >::iterator
1622 u2n0F = ++nodeParams[0].begin(),
1623 u2n1F = ++nodeParams[1].begin();
1624 map< double, const SMDS_MeshNode* >::reverse_iterator
1625 u2n1R = ++nodeParams[1].rbegin();
1626 for ( ; u2n0F != nodeParams[0].end(); ++u2n0F )
1628 if ( !theMA.getBoundary().getBranchPoint( iE0, u2n0F->first, brp ) ||
1629 !branch.getParameter( brp, maParam ))
1632 npPair.first = NodePoint( u2n0F->second, u2n0F->first, iE0 );
1635 npPair.second = NodePoint( u2n1R->second, u2n1R->first, iE1 );
1640 npPair.second = NodePoint( u2n1F->second, u2n1F->first, iE1 );
1643 pointsOnE.insert( make_pair( maParam, npPair ));
1646 } // loop on pairs of opposite EDGEs
1648 if ( !projectVertices( theHelper, theMA, divPoints, edgeIDs1, edgeIDs2,
1649 isComputed, pointsOnE, theSinuFace ))
1652 separateNodes( theHelper, theMA, pointsOnE, theSinuFace );
1655 TMAPar2NPoints::iterator u2np = pointsOnE.begin();
1656 for ( ; u2np != pointsOnE.end(); ++u2np )
1658 NodePoint* np[2] = { & u2np->second.first, & u2np->second.second };
1659 for ( int iSide = 0; iSide < 2; ++iSide )
1661 if ( np[ iSide ]->_node ) continue;
1662 size_t iEdge = np[ iSide ]->_edgeInd;
1663 double u = np[ iSide ]->_u;
1664 gp_Pnt p = curves[ iEdge ]->Value( u );
1665 np[ iSide ]->_node = meshDS->AddNode( p.X(), p.Y(), p.Z() );
1666 meshDS->SetNodeOnEdge( np[ iSide ]->_node, edgeIDs[ iEdge ], u );
1670 // create mesh segments on EDGEs
1671 theHelper.SetElementsOnShape( false );
1672 TopoDS_Face face = TopoDS::Face( theHelper.GetSubShape() );
1673 for ( size_t i = 0; i < theSinuEdges.size(); ++i )
1675 SMESH_subMesh* sm = mesh->GetSubMesh( theSinuEdges[i] );
1676 if ( sm->GetSubMeshDS() && sm->GetSubMeshDS()->NbElements() > 0 )
1679 StdMeshers_FaceSide side( face, theSinuEdges[i], mesh,
1680 /*isFwd=*/true, /*skipMediumNodes=*/true );
1681 vector<const SMDS_MeshNode*> nodes = side.GetOrderedNodes();
1682 for ( size_t in = 1; in < nodes.size(); ++in )
1684 const SMDS_MeshElement* seg = theHelper.AddEdge( nodes[in-1], nodes[in], 0, false );
1685 meshDS->SetMeshElementOnShape( seg, edgeIDs[ i ] );
1689 // update sub-meshes on VERTEXes
1690 for ( size_t i = 0; i < theSinuEdges.size(); ++i )
1692 mesh->GetSubMesh( theHelper.IthVertex( 0, theSinuEdges[i] ))
1693 ->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
1694 mesh->GetSubMesh( theHelper.IthVertex( 1, theSinuEdges[i] ))
1695 ->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
1698 // Setup sides of a quadrangle
1699 if ( !setQuadSides( theHelper, pointsOnE, theSinuFace, the1dAlgo ))
1705 //================================================================================
1707 * \brief Mesh short EDGEs
1709 //================================================================================
1711 bool computeShortEdges( SMESH_MesherHelper& theHelper,
1712 const vector<TopoDS_Edge>& theShortEdges,
1713 SMESH_Algo* the1dAlgo,
1714 const bool theHasRadialHyp,
1715 const bool theIs2nd)
1717 SMESH_Hypothesis::Hypothesis_Status aStatus;
1718 for ( size_t i = 0; i < theShortEdges.size(); ++i )
1720 if ( !theHasRadialHyp )
1722 theHelper.GetGen()->Compute( *theHelper.GetMesh(), theShortEdges[i], true, true );
1724 SMESH_subMesh* sm = theHelper.GetMesh()->GetSubMesh(theShortEdges[i] );
1725 if ( sm->IsEmpty() )
1727 // use 2D hyp or minSegLen
1730 SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(/*includeSelf=*/false);
1731 while ( smIt->more() )
1732 smIt->next()->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1735 the1dAlgo->CheckHypothesis( *theHelper.GetMesh(), theShortEdges[i], aStatus );
1736 if ( !the1dAlgo->Compute( *theHelper.GetMesh(), theShortEdges[i] ))
1742 sm->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
1743 if ( sm->IsEmpty() )
1750 inline double area( const UVPtStruct& p1, const UVPtStruct& p2, const UVPtStruct& p3 )
1752 gp_XY v1 = p2.UV() - p1.UV();
1753 gp_XY v2 = p3.UV() - p1.UV();
1757 bool ellipticSmooth( FaceQuadStruct::Ptr quad, int nbLoops )
1760 if ( quad->uv_grid.empty() )
1763 int nbhoriz = quad->iSize;
1764 int nbvertic = quad->jSize;
1766 const double dksi = 0.5, deta = 0.5;
1767 const double dksi2 = dksi*dksi, deta2 = deta*deta;
1768 double err = 0., g11, g22, g12;
1771 FaceQuadStruct& q = *quad;
1774 double refArea = area( q.UVPt(0,0), q.UVPt(1,0), q.UVPt(1,1) );
1776 for ( int iLoop = 0; iLoop < nbLoops; ++iLoop )
1779 for ( int i = 1; i < nbhoriz - 1; i++ )
1780 for ( int j = 1; j < nbvertic - 1; j++ )
1782 g11 = ( (q.U(i,j+1) - q.U(i,j-1))*(q.U(i,j+1) - q.U(i,j-1))/dksi2 +
1783 (q.V(i,j+1) - q.V(i,j-1))*(q.V(i,j+1) - q.V(i,j-1))/deta2 )/4;
1785 g22 = ( (q.U(i+1,j) - q.U(i-1,j))*(q.U(i+1,j) - q.U(i-1,j))/dksi2 +
1786 (q.V(i+1,j) - q.V(i-1,j))*(q.V(i+1,j) - q.V(i-1,j))/deta2 )/4;
1788 g12 = ( (q.U(i+1,j) - q.U(i-1,j))*(q.U(i,j+1) - q.U(i,j-1))/dksi2 +
1789 (q.V(i+1,j) - q.V(i-1,j))*(q.V(i,j+1) - q.V(i,j-1))/deta2 )/(4*dksi*deta);
1791 pNew.u = dksi2/(2*(g11+g22)) * (g11*(q.U(i+1,j) + q.U(i-1,j))/dksi2 +
1792 g22*(q.U(i,j+1) + q.U(i,j-1))/dksi2
1793 - 0.5*g12*q.U(i+1,j+1) + 0.5*g12*q.U(i-1,j+1) +
1794 - 0.5*g12*q.U(i-1,j-1) + 0.5*g12*q.U(i+1,j-1));
1796 pNew.v = deta2/(2*(g11+g22)) * (g11*(q.V(i+1,j) + q.V(i-1,j))/deta2 +
1797 g22*(q.V(i,j+1) + q.V(i,j-1))/deta2
1798 - 0.5*g12*q.V(i+1,j+1) + 0.5*g12*q.V(i-1,j+1) +
1799 - 0.5*g12*q.V(i-1,j-1) + 0.5*g12*q.V(i+1,j-1));
1801 // if (( refArea * area( q.UVPt(i-1,j-1), q.UVPt(i,j-1), pNew ) > 0 ) &&
1802 // ( refArea * area( q.UVPt(i+1,j-1), q.UVPt(i+1,j), pNew ) > 0 ) &&
1803 // ( refArea * area( q.UVPt(i+1,j+1), q.UVPt(i,j+1), pNew ) > 0 ) &&
1804 // ( refArea * area( q.UVPt(i-1,j), q.UVPt(i-1,j-1), pNew ) > 0 ))
1806 err += sqrt(( q.U(i,j) - pNew.u ) * ( q.U(i,j) - pNew.u ) +
1807 ( q.V(i,j) - pNew.v ) * ( q.V(i,j) - pNew.v ));
1811 // else if ( ++nbErr < 10 )
1813 // cout << i << ", " << j << endl;
1815 // << "[ " << q.U(i-1,j-1) << ", " <<q.U(i,j-1) << ", " << q.U(i+1,j-1) << " ],"
1816 // << "[ " << q.U(i-1,j-0) << ", " <<q.U(i,j-0) << ", " << q.U(i+1,j-0) << " ],"
1817 // << "[ " << q.U(i-1,j+1) << ", " <<q.U(i,j+1) << ", " << q.U(i+1,j+1) << " ]]" << endl;
1819 // << "[ " << q.V(i-1,j-1) << ", " <<q.V(i,j-1) << ", " << q.V(i+1,j-1) << " ],"
1820 // << "[ " << q.V(i-1,j-0) << ", " <<q.V(i,j-0) << ", " << q.V(i+1,j-0) << " ],"
1821 // << "[ " << q.V(i-1,j+1) << ", " <<q.V(i,j+1) << ", " << q.V(i+1,j+1) << " ]]" << endl<<endl;
1825 if ( err / ( nbhoriz - 2 ) / ( nbvertic - 2 ) < 1e-6 )
1828 //cout << " ERR " << err / ( nbhoriz - 2 ) / ( nbvertic - 2 ) << endl;
1833 //================================================================================
1835 * \brief Remove temporary node
1837 //================================================================================
1839 void mergeNodes( SMESH_MesherHelper& theHelper,
1840 SinuousFace& theSinuFace )
1842 SMESH_MeshEditor editor( theHelper.GetMesh() );
1843 SMESH_MeshEditor::TListOfListOfNodes nodesGroups;
1845 TMergeMap::iterator n2nn = theSinuFace._nodesToMerge.begin();
1846 for ( ; n2nn != theSinuFace._nodesToMerge.end(); ++n2nn )
1848 if ( !n2nn->first ) continue;
1849 nodesGroups.push_back( list< const SMDS_MeshNode* >() );
1850 list< const SMDS_MeshNode* > & group = nodesGroups.back();
1852 group.push_back( n2nn->first );
1853 group.splice( group.end(), n2nn->second );
1855 editor.MergeNodes( nodesGroups );
1860 //================================================================================
1862 * \brief Sets event listener which removes mesh from EDGEs when 2D hyps change
1864 //================================================================================
1866 void StdMeshers_QuadFromMedialAxis_1D2D::SetEventListener(SMESH_subMesh* faceSubMesh)
1868 faceSubMesh->SetEventListener( new EdgeCleaner, 0, faceSubMesh );
1871 //================================================================================
1873 * \brief Create quadrangle elements
1874 * \param [in] theHelper - the helper
1875 * \param [in] theFace - the face to mesh
1876 * \param [in] theSinuEdges - the sinuous EDGEs
1877 * \param [in] theShortEdges - the short EDGEs
1878 * \return bool - is OK or not
1880 //================================================================================
1882 bool StdMeshers_QuadFromMedialAxis_1D2D::computeQuads( SMESH_MesherHelper& theHelper,
1883 FaceQuadStruct::Ptr theQuad)
1885 StdMeshers_Quadrangle_2D::myHelper = &theHelper;
1886 StdMeshers_Quadrangle_2D::myNeedSmooth = false;
1887 StdMeshers_Quadrangle_2D::myCheckOri = false;
1888 StdMeshers_Quadrangle_2D::myQuadList.clear();
1890 int nbNodesShort0 = theQuad->side[0].NbPoints();
1891 int nbNodesShort1 = theQuad->side[2].NbPoints();
1893 // compute UV of internal points
1894 myQuadList.push_back( theQuad );
1895 if ( !StdMeshers_Quadrangle_2D::setNormalizedGrid( theQuad ))
1898 // elliptic smooth of internal points to get boundary cell normal to the boundary
1899 bool isRing = theQuad->side[0].grid->Edge(0).IsNull();
1901 ellipticSmooth( theQuad, 1 );
1903 // create quadrangles
1905 if ( nbNodesShort0 == nbNodesShort1 )
1906 ok = StdMeshers_Quadrangle_2D::computeQuadDominant( *theHelper.GetMesh(),
1907 theQuad->face, theQuad );
1909 ok = StdMeshers_Quadrangle_2D::computeTriangles( *theHelper.GetMesh(),
1910 theQuad->face, theQuad );
1912 StdMeshers_Quadrangle_2D::myHelper = 0;
1917 //================================================================================
1919 * \brief Generate quadrangle mesh
1921 //================================================================================
1923 bool StdMeshers_QuadFromMedialAxis_1D2D::Compute(SMESH_Mesh& theMesh,
1924 const TopoDS_Shape& theShape)
1926 SMESH_MesherHelper helper( theMesh );
1927 helper.SetSubShape( theShape );
1929 TopoDS_Face F = TopoDS::Face( theShape );
1930 if ( F.Orientation() >= TopAbs_INTERNAL ) F.Orientation( TopAbs_FORWARD );
1932 SinuousFace sinuFace( F );
1936 if ( getSinuousEdges( helper, sinuFace ))
1940 double minSegLen = getMinSegLen( helper, sinuFace._sinuEdges );
1941 SMESH_MAT2d::MedialAxis ma( F, sinuFace._sinuEdges, minSegLen, /*ignoreCorners=*/true );
1944 _regular1D = new Algo1D( _studyId, _gen );
1945 _regular1D->SetSegmentLength( minSegLen );
1947 vector<double> maParams;
1948 if ( ! divideMA( helper, ma, sinuFace, _regular1D, minSegLen, maParams ))
1949 return error(COMPERR_BAD_SHAPE);
1953 _regular1D->SetRadialDistribution( _hyp2D );
1955 if ( !computeShortEdges( helper, sinuFace._shortSide[0], _regular1D, _hyp2D, 0 ) ||
1956 !computeShortEdges( helper, sinuFace._shortSide[1], _regular1D, _hyp2D, 1 ))
1957 return error("Failed to mesh short edges");
1961 if ( !computeSinuEdges( helper, minSegLen, ma, maParams, sinuFace, _regular1D ))
1962 return error("Failed to mesh sinuous edges");
1966 bool ok = computeQuads( helper, sinuFace._quad );
1969 mergeNodes( helper, sinuFace );
1976 return error(COMPERR_BAD_SHAPE, "Not implemented so far");
1979 //================================================================================
1981 * \brief Predict nb of elements
1983 //================================================================================
1985 bool StdMeshers_QuadFromMedialAxis_1D2D::Evaluate(SMESH_Mesh & theMesh,
1986 const TopoDS_Shape & theShape,
1987 MapShapeNbElems& theResMap)
1989 return StdMeshers_Quadrangle_2D::Evaluate(theMesh,theShape,theResMap);
1992 //================================================================================
1994 * \brief Return true if the algorithm can mesh this shape
1995 * \param [in] aShape - shape to check
1996 * \param [in] toCheckAll - if true, this check returns OK if all shapes are OK,
1997 * else, returns OK if at least one shape is OK
1999 //================================================================================
2001 bool StdMeshers_QuadFromMedialAxis_1D2D::IsApplicable( const TopoDS_Shape & aShape,
2005 SMESH_MesherHelper helper( tmpMesh );
2007 int nbFoundFaces = 0;
2008 for (TopExp_Explorer exp( aShape, TopAbs_FACE ); exp.More(); exp.Next(), ++nbFoundFaces )
2010 const TopoDS_Face& face = TopoDS::Face( exp.Current() );
2011 SinuousFace sinuFace( face );
2012 bool isApplicable = getSinuousEdges( helper, sinuFace );
2014 if ( toCheckAll && !isApplicable ) return false;
2015 if ( !toCheckAll && isApplicable ) return true;
2017 return ( toCheckAll && nbFoundFaces != 0 );