1 // Copyright (C) 2007-2015 CEA/DEN, EDF R&D, OPEN CASCADE
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
20 // SMESH SMESH : implementaion of SMESH idl descriptions
21 // File : StdMeshers_RadialQuadrangle_1D2D.cxx
24 #include "StdMeshers_RadialQuadrangle_1D2D.hxx"
26 #include "StdMeshers_NumberOfLayers.hxx"
27 #include "StdMeshers_LayerDistribution.hxx"
28 #include "StdMeshers_Regular_1D.hxx"
29 #include "StdMeshers_NumberOfSegments.hxx"
31 #include "SMDS_MeshNode.hxx"
32 #include "SMESHDS_SubMesh.hxx"
33 #include "SMESH_Gen.hxx"
34 #include "SMESH_HypoFilter.hxx"
35 #include "SMESH_Mesh.hxx"
36 #include "SMESH_MesherHelper.hxx"
37 #include "SMESH_subMesh.hxx"
38 #include "SMESH_subMeshEventListener.hxx"
40 #include "utilities.h"
42 #include <BRepAdaptor_Curve.hxx>
43 #include <BRepBuilderAPI_MakeEdge.hxx>
44 #include <BRep_Tool.hxx>
45 #include <GeomAPI_ProjectPointOnSurf.hxx>
46 #include <Geom_Circle.hxx>
47 #include <Geom_Line.hxx>
48 #include <Geom_TrimmedCurve.hxx>
49 #include <TColgp_SequenceOfPnt.hxx>
50 #include <TColgp_SequenceOfPnt2d.hxx>
52 #include <TopExp_Explorer.hxx>
53 #include <TopTools_ListIteratorOfListOfShape.hxx>
59 #define RETURN_BAD_RESULT(msg) { MESSAGE(")-: Error: " << msg); return false; }
60 #define gpXYZ(n) gp_XYZ(n->X(),n->Y(),n->Z())
63 //=======================================================================
64 //function : StdMeshers_RadialQuadrangle_1D2D
66 //=======================================================================
68 StdMeshers_RadialQuadrangle_1D2D::StdMeshers_RadialQuadrangle_1D2D(int hypId,
71 :SMESH_2D_Algo(hypId, studyId, gen)
73 _name = "RadialQuadrangle_1D2D";
74 _shapeType = (1 << TopAbs_FACE); // 1 bit per shape type
76 _compatibleHypothesis.push_back("LayerDistribution2D");
77 _compatibleHypothesis.push_back("NumberOfLayers2D");
79 myDistributionHypo = 0;
80 _requireDiscreteBoundary = false;
81 _supportSubmeshes = true;
85 //================================================================================
89 //================================================================================
91 StdMeshers_RadialQuadrangle_1D2D::~StdMeshers_RadialQuadrangle_1D2D()
95 //=======================================================================
96 //function : CheckHypothesis
98 //=======================================================================
100 bool StdMeshers_RadialQuadrangle_1D2D::CheckHypothesis
102 const TopoDS_Shape& aShape,
103 SMESH_Hypothesis::Hypothesis_Status& aStatus)
107 myDistributionHypo = 0;
109 list <const SMESHDS_Hypothesis * >::const_iterator itl;
111 const list <const SMESHDS_Hypothesis * >&hyps = GetUsedHypothesis(aMesh, aShape);
112 if ( hyps.size() == 0 ) {
113 aStatus = SMESH_Hypothesis::HYP_OK;
114 return true; // can work with no hypothesis
117 if ( hyps.size() > 1 ) {
118 aStatus = SMESH_Hypothesis::HYP_ALREADY_EXIST;
122 const SMESHDS_Hypothesis *theHyp = hyps.front();
124 string hypName = theHyp->GetName();
126 if (hypName == "NumberOfLayers2D") {
127 myNbLayerHypo = static_cast<const StdMeshers_NumberOfLayers *>(theHyp);
128 aStatus = SMESH_Hypothesis::HYP_OK;
131 if (hypName == "LayerDistribution2D") {
132 myDistributionHypo = static_cast<const StdMeshers_LayerDistribution *>(theHyp);
133 aStatus = SMESH_Hypothesis::HYP_OK;
136 aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
142 // ------------------------------------------------------------------------------
144 * \brief Listener used to mark edges meshed by StdMeshers_RadialQuadrangle_1D2D
146 class TEdgeMarker : public SMESH_subMeshEventListener
148 TEdgeMarker(): SMESH_subMeshEventListener(/*isDeletable=*/false,
149 "StdMeshers_RadialQuadrangle_1D2D::TEdgeMarker") {}
151 //!< Return static listener
152 static SMESH_subMeshEventListener* getListener()
154 static TEdgeMarker theEdgeMarker;
155 return &theEdgeMarker;
157 //! Clear face sumbesh if something happens on edges
158 void ProcessEvent(const int event,
160 SMESH_subMesh* edgeSubMesh,
161 EventListenerData* data,
162 const SMESH_Hypothesis* /*hyp*/)
164 if ( data && !data->mySubMeshes.empty() && eventType == SMESH_subMesh::ALGO_EVENT)
166 ASSERT( data->mySubMeshes.front() != edgeSubMesh );
167 SMESH_subMesh* faceSubMesh = data->mySubMeshes.front();
168 faceSubMesh->ComputeStateEngine( SMESH_subMesh::CLEAN );
173 // ------------------------------------------------------------------------------
175 * \brief Mark an edge as computed by StdMeshers_RadialQuadrangle_1D2D
177 void markEdgeAsComputedByMe(const TopoDS_Edge& edge, SMESH_subMesh* faceSubMesh)
179 if ( SMESH_subMesh* edgeSM = faceSubMesh->GetFather()->GetSubMeshContaining( edge ))
181 if ( !edgeSM->GetEventListenerData( TEdgeMarker::getListener() ))
182 faceSubMesh->SetEventListener( TEdgeMarker::getListener(),
183 SMESH_subMeshEventListenerData::MakeData(faceSubMesh),
187 // ------------------------------------------------------------------------------
189 * \brief Return true if a radial edge was meshed with StdMeshers_RadialQuadrangle_1D2D with
190 * the same radial distribution
192 // bool isEdgeCompatiballyMeshed(const TopoDS_Edge& edge, SMESH_subMesh* faceSubMesh)
194 // if ( SMESH_subMesh* edgeSM = faceSubMesh->GetFather()->GetSubMeshContaining( edge ))
196 // if ( SMESH_subMeshEventListenerData* otherFaceData =
197 // edgeSM->GetEventListenerData( TEdgeMarker::getListener() ))
199 // // compare hypothesis aplied to two disk faces sharing radial edges
200 // SMESH_Mesh& mesh = *faceSubMesh->GetFather();
201 // SMESH_Algo* radialQuadAlgo = mesh.GetGen()->GetAlgo(mesh, faceSubMesh->GetSubShape() );
202 // SMESH_subMesh* otherFaceSubMesh = otherFaceData->mySubMeshes.front();
203 // list <const SMESHDS_Hypothesis *> hyps1 =
204 // radialQuadAlgo->GetUsedHypothesis( mesh, faceSubMesh->GetSubShape());
205 // list <const SMESHDS_Hypothesis *> hyps2 =
206 // radialQuadAlgo->GetUsedHypothesis( mesh, otherFaceSubMesh->GetSubShape());
207 // if( hyps1.empty() && hyps2.empty() )
208 // return true; // defaul hyps
209 // if ( hyps1.size() != hyps2.size() )
211 // return *hyps1.front() == *hyps2.front();
217 //================================================================================
219 * \brief Return base curve of the edge and extremum parameters
221 //================================================================================
223 Handle(Geom_Curve) getCurve(const TopoDS_Edge& edge, double* f=0, double* l=0)
225 Handle(Geom_Curve) C;
226 if ( !edge.IsNull() )
228 double first = 0., last = 0.;
229 C = BRep_Tool::Curve(edge, first, last);
232 Handle(Geom_TrimmedCurve) tc = Handle(Geom_TrimmedCurve)::DownCast(C);
233 while( !tc.IsNull() ) {
234 C = tc->BasisCurve();
235 tc = Handle(Geom_TrimmedCurve)::DownCast(C);
244 //================================================================================
246 * \brief Return edges of the face
247 * \retval int - nb of edges
249 //================================================================================
251 int analyseFace(const TopoDS_Shape& face,
252 TopoDS_Edge& CircEdge,
253 TopoDS_Edge& LinEdge1,
254 TopoDS_Edge& LinEdge2)
256 CircEdge.Nullify(); LinEdge1.Nullify(); LinEdge2.Nullify();
259 for ( TopExp_Explorer exp( face, TopAbs_EDGE ); exp.More(); exp.Next(), ++nbe )
261 const TopoDS_Edge& E = TopoDS::Edge( exp.Current() );
263 Handle(Geom_Curve) C = getCurve(E,&f,&l);
266 if ( C->IsKind( STANDARD_TYPE(Geom_Circle)))
268 if ( CircEdge.IsNull() )
273 else if ( LinEdge1.IsNull() )
282 //================================================================================
283 //================================================================================
285 * \brief Class computing layers distribution using data of
286 * StdMeshers_LayerDistribution hypothesis
288 //================================================================================
289 //================================================================================
291 class TNodeDistributor: public StdMeshers_Regular_1D
293 list <const SMESHDS_Hypothesis *> myUsedHyps;
295 // -----------------------------------------------------------------------------
296 static TNodeDistributor* GetDistributor(SMESH_Mesh& aMesh)
298 const int myID = -1001;
299 TNodeDistributor* myHyp = dynamic_cast<TNodeDistributor*>( aMesh.GetHypothesis( myID ));
301 myHyp = new TNodeDistributor( myID, 0, aMesh.GetGen() );
304 // -----------------------------------------------------------------------------
305 //! Computes distribution of nodes on a straight line ending at pIn and pOut
306 bool Compute( vector< double > & positions,
310 const SMESH_Hypothesis* hyp1d)
312 if ( !hyp1d ) return error( "Invalid LayerDistribution hypothesis");
314 double len = pIn.Distance( pOut );
315 if ( len <= DBL_MIN ) return error("Too close points of inner and outer shells");
318 myUsedHyps.push_back( hyp1d );
320 TopoDS_Edge edge = BRepBuilderAPI_MakeEdge( pIn, pOut );
321 SMESH_Hypothesis::Hypothesis_Status aStatus;
322 if ( !StdMeshers_Regular_1D::CheckHypothesis( aMesh, edge, aStatus ))
323 return error( "StdMeshers_Regular_1D::CheckHypothesis() failed "
324 "with LayerDistribution hypothesis");
326 BRepAdaptor_Curve C3D(edge);
327 double f = C3D.FirstParameter(), l = C3D.LastParameter();
328 list< double > params;
329 if ( !StdMeshers_Regular_1D::computeInternalParameters( aMesh, C3D, len, f, l, params, false ))
330 return error("StdMeshers_Regular_1D failed to compute layers distribution");
333 positions.reserve( params.size() );
334 for (list<double>::iterator itU = params.begin(); itU != params.end(); itU++)
335 positions.push_back( *itU / len );
338 // -----------------------------------------------------------------------------
339 //! Make mesh on an adge using assigned 1d hyp or defaut nb of segments
340 bool ComputeCircularEdge(SMESH_Mesh& aMesh,
341 const TopoDS_Edge& anEdge)
343 _gen->Compute( aMesh, anEdge);
344 SMESH_subMesh *sm = aMesh.GetSubMesh(anEdge);
345 if ( sm->GetComputeState() != SMESH_subMesh::COMPUTE_OK)
347 // find any 1d hyp assigned (there can be a hyp w/o algo)
348 myUsedHyps = SMESH_Algo::GetUsedHypothesis(aMesh, anEdge, /*ignoreAux=*/true);
349 Hypothesis_Status aStatus;
350 if ( !StdMeshers_Regular_1D::CheckHypothesis( aMesh, anEdge, aStatus ))
352 // no valid 1d hyp assigned, use default nb of segments
353 _hypType = NB_SEGMENTS;
354 _ivalue[ DISTR_TYPE_IND ] = StdMeshers_NumberOfSegments::DT_Regular;
355 _ivalue[ NB_SEGMENTS_IND ] = _gen->GetDefaultNbSegments();
357 return StdMeshers_Regular_1D::Compute( aMesh, anEdge );
361 // -----------------------------------------------------------------------------
362 //! Make mesh on an adge using assigned 1d hyp or defaut nb of segments
363 bool EvaluateCircularEdge(SMESH_Mesh& aMesh,
364 const TopoDS_Edge& anEdge,
365 MapShapeNbElems& aResMap)
367 _gen->Evaluate( aMesh, anEdge, aResMap );
368 if ( aResMap.count( aMesh.GetSubMesh( anEdge )))
371 // find any 1d hyp assigned
372 myUsedHyps = SMESH_Algo::GetUsedHypothesis(aMesh, anEdge, /*ignoreAux=*/true);
373 Hypothesis_Status aStatus;
374 if ( !StdMeshers_Regular_1D::CheckHypothesis( aMesh, anEdge, aStatus ))
376 // no valid 1d hyp assigned, use default nb of segments
377 _hypType = NB_SEGMENTS;
378 _ivalue[ DISTR_TYPE_IND ] = StdMeshers_NumberOfSegments::DT_Regular;
379 _ivalue[ NB_SEGMENTS_IND ] = _gen->GetDefaultNbSegments();
381 return StdMeshers_Regular_1D::Evaluate( aMesh, anEdge, aResMap );
384 // -----------------------------------------------------------------------------
385 TNodeDistributor( int hypId, int studyId, SMESH_Gen* gen)
386 : StdMeshers_Regular_1D( hypId, studyId, gen)
389 // -----------------------------------------------------------------------------
390 virtual const list <const SMESHDS_Hypothesis *> &
391 GetUsedHypothesis(SMESH_Mesh &, const TopoDS_Shape &, const bool)
395 // -----------------------------------------------------------------------------
399 //=======================================================================
401 * \brief Allow algo to do something after persistent restoration
402 * \param subMesh - restored submesh
404 * call markEdgeAsComputedByMe()
406 //=======================================================================
408 void StdMeshers_RadialQuadrangle_1D2D::SubmeshRestored(SMESH_subMesh* faceSubMesh)
410 if ( !faceSubMesh->IsEmpty() )
412 TopoDS_Edge CircEdge, LinEdge1, LinEdge2;
413 analyseFace( faceSubMesh->GetSubShape(), CircEdge, LinEdge1, LinEdge2 );
414 if ( !CircEdge.IsNull() ) markEdgeAsComputedByMe( CircEdge, faceSubMesh );
415 if ( !LinEdge1.IsNull() ) markEdgeAsComputedByMe( LinEdge1, faceSubMesh );
416 if ( !LinEdge2.IsNull() ) markEdgeAsComputedByMe( LinEdge2, faceSubMesh );
420 //=======================================================================
423 //=======================================================================
425 bool StdMeshers_RadialQuadrangle_1D2D::Compute(SMESH_Mesh& aMesh,
426 const TopoDS_Shape& aShape)
428 SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
430 myHelper = new SMESH_MesherHelper( aMesh );
431 // to delete helper at exit from Compute()
432 SMESHUtils::Deleter<SMESH_MesherHelper> helperDeleter( myHelper );
434 TNodeDistributor* algo1d = TNodeDistributor::GetDistributor(aMesh);
436 TopoDS_Edge CircEdge, LinEdge1, LinEdge2;
437 int nbe = analyseFace( aShape, CircEdge, LinEdge1, LinEdge2 );
438 Handle(Geom_Circle) aCirc = Handle(Geom_Circle)::DownCast( getCurve( CircEdge ));
439 if( nbe > 3 || nbe < 1 || aCirc.IsNull() )
440 return error("The face must be a full circle or a part of circle (i.e. the number "
441 "of edges is less or equal to 3 and one of them is a circle curve)");
444 // points for rotation
445 TColgp_SequenceOfPnt Points;
446 // angles for rotation
447 TColStd_SequenceOfReal Angles;
448 // Nodes1 and Nodes2 - nodes along radiuses
449 // CNodes - nodes on circle edge
450 vector< const SMDS_MeshNode* > Nodes1, Nodes2, CNodes;
452 // parameters edge nodes on face
453 TColgp_SequenceOfPnt2d Pnts2d1;
456 int faceID = meshDS->ShapeToIndex(aShape);
457 TopoDS_Face F = TopoDS::Face(aShape);
458 Handle(Geom_Surface) S = BRep_Tool::Surface(F);
463 if (!algo1d->ComputeCircularEdge( aMesh, CircEdge ))
464 return error( algo1d->GetComputeError() );
465 map< double, const SMDS_MeshNode* > theNodes;
466 if ( !GetSortedNodesOnEdge(aMesh.GetMeshDS(),CircEdge,true,theNodes))
467 return error("Circular edge is incorrectly meshed");
469 myHelper->IsQuadraticSubMesh( aShape );
472 map< double, const SMDS_MeshNode* >::iterator itn = theNodes.begin();
473 const SMDS_MeshNode* NF = (*itn).second;
474 CNodes.push_back( (*itn).second );
475 double fang = (*itn).first;
476 if ( itn != theNodes.end() ) {
478 for(; itn != theNodes.end(); itn++ ) {
479 CNodes.push_back( (*itn).second );
480 double ang = (*itn).first - fang;
481 if( ang>M_PI ) ang = ang - 2.*M_PI;
482 if( ang<-M_PI ) ang = ang + 2.*M_PI;
483 Angles.Append( ang );
486 P1 = gp_Pnt( NF->X(), NF->Y(), NF->Z() );
487 P0 = aCirc->Location();
489 if ( !computeLayerPositions(P0,P1))
492 TopoDS_Vertex V1 = myHelper->IthVertex(0, CircEdge );
493 gp_Pnt2d p2dV = BRep_Tool::Parameters( V1, TopoDS::Face(aShape) );
495 NC = meshDS->AddNode(P0.X(), P0.Y(), P0.Z());
496 GeomAPI_ProjectPointOnSurf PPS(P0,S);
498 PPS.Parameters(1,U0,V0);
499 meshDS->SetNodeOnFace(NC, faceID, U0, V0);
500 PC = gp_Pnt2d(U0,V0);
503 gp_Vec2d aVec2d(PC,p2dV);
504 Nodes1.resize( myLayerPositions.size()+1 );
505 Nodes2.resize( myLayerPositions.size()+1 );
507 for(; i<myLayerPositions.size(); i++) {
508 gp_Pnt P( P0.X() + aVec.X()*myLayerPositions[i],
509 P0.Y() + aVec.Y()*myLayerPositions[i],
510 P0.Z() + aVec.Z()*myLayerPositions[i] );
512 SMDS_MeshNode * node = meshDS->AddNode(P.X(), P.Y(), P.Z());
515 double U = PC.X() + aVec2d.X()*myLayerPositions[i];
516 double V = PC.Y() + aVec2d.Y()*myLayerPositions[i];
517 meshDS->SetNodeOnFace( node, faceID, U, V );
518 Pnts2d1.Append(gp_Pnt2d(U,V));
520 Nodes1[Nodes1.size()-1] = NF;
521 Nodes2[Nodes1.size()-1] = NF;
523 else if(nbe==2 && LinEdge1.Orientation() != TopAbs_INTERNAL )
525 // one curve must be a half of circle and other curve must be
528 Handle(Geom_Circle) aCirc = Handle(Geom_Circle)::DownCast( getCurve( CircEdge, &fp, &lp ));
529 if( fabs(fabs(lp-fp)-M_PI) > Precision::Confusion() ) {
530 // not half of circle
531 return error(COMPERR_BAD_SHAPE);
533 Handle(Geom_Line) aLine = Handle(Geom_Line)::DownCast( getCurve( LinEdge1 ));
534 if( aLine.IsNull() ) {
535 // other curve not line
536 return error(COMPERR_BAD_SHAPE);
539 if ( !algo1d->ComputeCircularEdge( aMesh, CircEdge ))
540 return error( algo1d->GetComputeError() );
541 map< double, const SMDS_MeshNode* > theNodes;
542 if ( !GetSortedNodesOnEdge(aMesh.GetMeshDS(),CircEdge,true,theNodes) )
543 return error("Circular edge is incorrectly meshed");
545 myHelper->IsQuadraticSubMesh( aShape );
547 map< double, const SMDS_MeshNode* >::iterator itn = theNodes.begin();
549 CNodes.push_back( itn->second );
550 double fang = (*itn).first;
552 for(; itn != theNodes.end(); itn++ ) {
553 CNodes.push_back( (*itn).second );
554 double ang = (*itn).first - fang;
555 if( ang>M_PI ) ang = ang - 2.*M_PI;
556 if( ang<-M_PI ) ang = ang + 2.*M_PI;
557 Angles.Append( ang );
559 const SMDS_MeshNode* NF = theNodes.begin()->second;
560 const SMDS_MeshNode* NL = theNodes.rbegin()->second;
561 P1 = gp_Pnt( NF->X(), NF->Y(), NF->Z() );
562 gp_Pnt P2( NL->X(), NL->Y(), NL->Z() );
563 P0 = aCirc->Location();
565 bool linEdgeComputed;
566 if ( !computeLayerPositions(P0,P1,LinEdge1,&linEdgeComputed))
569 if ( linEdgeComputed )
571 if (!GetSortedNodesOnEdge(aMesh.GetMeshDS(),LinEdge1,true,theNodes))
572 return error("Invalid mesh on a straight edge");
574 Nodes1.resize( myLayerPositions.size()+1 );
575 Nodes2.resize( myLayerPositions.size()+1 );
576 vector< const SMDS_MeshNode* > *pNodes1 = &Nodes1, *pNodes2 = &Nodes2;
577 bool nodesFromP0ToP1 = ( theNodes.rbegin()->second == NF );
578 if ( !nodesFromP0ToP1 ) std::swap( pNodes1, pNodes2 );
580 map< double, const SMDS_MeshNode* >::reverse_iterator ritn = theNodes.rbegin();
581 itn = theNodes.begin();
582 for ( int i = Nodes1.size()-1; i > -1; ++itn, ++ritn, --i )
584 (*pNodes1)[i] = ritn->second;
585 (*pNodes2)[i] = itn->second;
586 Points.Prepend( gpXYZ( Nodes1[i]));
587 Pnts2d1.Prepend( myHelper->GetNodeUV( F, Nodes1[i]));
589 NC = const_cast<SMDS_MeshNode*>( itn->second );
590 Points.Remove( Nodes1.size() );
595 int edgeID = meshDS->ShapeToIndex(LinEdge1);
597 Handle(Geom_Curve) Crv = BRep_Tool::Curve(LinEdge1,fp,lp);
601 if( P1.Distance(Ptmp) > Precision::Confusion() )
603 // get UV points for edge
605 BRep_Tool::UVPoints( LinEdge1, TopoDS::Face(aShape), PF, PL );
606 PC = gp_Pnt2d( (PF.X()+PL.X())/2, (PF.Y()+PL.Y())/2 );
608 if(ori) V2d = gp_Vec2d(PC,PF);
609 else V2d = gp_Vec2d(PC,PL);
611 double cp = (fp+lp)/2;
612 double dp2 = (lp-fp)/2;
613 NC = meshDS->AddNode(P0.X(), P0.Y(), P0.Z());
614 meshDS->SetNodeOnEdge(NC, edgeID, cp);
615 Nodes1.resize( myLayerPositions.size()+1 );
616 Nodes2.resize( myLayerPositions.size()+1 );
618 for(; i<myLayerPositions.size(); i++) {
619 gp_Pnt P( P0.X() + aVec.X()*myLayerPositions[i],
620 P0.Y() + aVec.Y()*myLayerPositions[i],
621 P0.Z() + aVec.Z()*myLayerPositions[i] );
623 SMDS_MeshNode * node = meshDS->AddNode(P.X(), P.Y(), P.Z());
627 param = fp + dp2*(1-myLayerPositions[i]);
629 param = cp + dp2*myLayerPositions[i];
630 meshDS->SetNodeOnEdge(node, edgeID, param);
631 P = gp_Pnt( P0.X() - aVec.X()*myLayerPositions[i],
632 P0.Y() - aVec.Y()*myLayerPositions[i],
633 P0.Z() - aVec.Z()*myLayerPositions[i] );
634 node = meshDS->AddNode(P.X(), P.Y(), P.Z());
637 param = fp + dp2*(1-myLayerPositions[i]);
639 param = cp + dp2*myLayerPositions[i];
640 meshDS->SetNodeOnEdge(node, edgeID, param);
641 // parameters on face
642 gp_Pnt2d P2d( PC.X() + V2d.X()*myLayerPositions[i],
643 PC.Y() + V2d.Y()*myLayerPositions[i] );
646 Nodes1[ myLayerPositions.size() ] = NF;
647 Nodes2[ myLayerPositions.size() ] = NL;
648 // create 1D elements on edge
649 vector< const SMDS_MeshNode* > tmpNodes;
650 tmpNodes.resize(2*Nodes1.size()+1);
651 for(i=0; i<Nodes2.size(); i++)
652 tmpNodes[Nodes2.size()-i-1] = Nodes2[i];
653 tmpNodes[Nodes2.size()] = NC;
654 for(i=0; i<Nodes1.size(); i++)
655 tmpNodes[Nodes2.size()+1+i] = Nodes1[i];
656 for(i=1; i<tmpNodes.size(); i++) {
657 SMDS_MeshEdge* ME = myHelper->AddEdge( tmpNodes[i-1], tmpNodes[i] );
658 if(ME) meshDS->SetMeshElementOnShape(ME, edgeID);
660 markEdgeAsComputedByMe( LinEdge1, aMesh.GetSubMesh( F ));
663 else // nbe==3 or ( nbe==2 && linEdge is INTERNAL )
665 if (nbe==2 && LinEdge1.Orientation() == TopAbs_INTERNAL )
668 // one curve must be a part of circle and other curves must be
671 Handle(Geom_Circle) aCirc = Handle(Geom_Circle)::DownCast( getCurve( CircEdge ));
672 Handle(Geom_Line) aLine1 = Handle(Geom_Line)::DownCast( getCurve( LinEdge1 ));
673 Handle(Geom_Line) aLine2 = Handle(Geom_Line)::DownCast( getCurve( LinEdge2 ));
674 if( aCirc.IsNull() || aLine1.IsNull() || aLine2.IsNull() )
675 return error(COMPERR_BAD_SHAPE);
677 if ( !algo1d->ComputeCircularEdge( aMesh, CircEdge ))
678 return error( algo1d->GetComputeError() );
679 map< double, const SMDS_MeshNode* > theNodes;
680 if ( !GetSortedNodesOnEdge(aMesh.GetMeshDS(),CircEdge,true,theNodes))
681 return error("Circular edge is incorrectly meshed");
683 myHelper->IsQuadraticSubMesh( aShape );
685 const SMDS_MeshNode* NF = theNodes.begin()->second;
686 const SMDS_MeshNode* NL = theNodes.rbegin()->second;
688 CNodes.push_back( NF );
689 map< double, const SMDS_MeshNode* >::iterator itn = theNodes.begin();
690 double fang = (*itn).first;
692 for(; itn != theNodes.end(); itn++ ) {
693 CNodes.push_back( (*itn).second );
694 double ang = (*itn).first - fang;
695 if( ang>M_PI ) ang = ang - 2.*M_PI;
696 if( ang<-M_PI ) ang = ang + 2.*M_PI;
697 Angles.Append( ang );
699 P1 = gp_Pnt( NF->X(), NF->Y(), NF->Z() );
700 gp_Pnt P2( NL->X(), NL->Y(), NL->Z() );
701 P0 = aCirc->Location();
703 // make P1 belong to LinEdge1
704 TopoDS_Vertex V1 = myHelper->IthVertex( 0, LinEdge1 );
705 TopoDS_Vertex V2 = myHelper->IthVertex( 1, LinEdge1 );
706 gp_Pnt PE1 = BRep_Tool::Pnt(V1);
707 gp_Pnt PE2 = BRep_Tool::Pnt(V2);
708 if( ( P1.Distance(PE1) > Precision::Confusion() ) &&
709 ( P1.Distance(PE2) > Precision::Confusion() ) )
710 std::swap( LinEdge1, LinEdge2 );
712 bool linEdge1Computed, linEdge2Computed;
713 if ( !computeLayerPositions(P0,P1,LinEdge1,&linEdge1Computed))
716 Nodes1.resize( myLayerPositions.size()+1 );
717 Nodes2.resize( myLayerPositions.size()+1 );
719 // check that both linear edges have same hypotheses
720 if ( !computeLayerPositions(P0,P2,LinEdge2, &linEdge2Computed))
722 if ( Nodes1.size() != myLayerPositions.size()+1 )
723 return error("Different hypotheses apply to radial edges");
725 // find the central vertex
726 TopoDS_Vertex VC = V2;
727 if( ( P1.Distance(PE1) > Precision::Confusion() ) &&
728 ( P2.Distance(PE1) > Precision::Confusion() ) )
730 int vertID = meshDS->ShapeToIndex(VC);
733 if ( linEdge1Computed )
735 if (!GetSortedNodesOnEdge(aMesh.GetMeshDS(),LinEdge1,true,theNodes))
736 return error("Invalid mesh on a straight edge");
738 bool nodesFromP0ToP1 = ( theNodes.rbegin()->second == NF );
739 NC = const_cast<SMDS_MeshNode*>
740 ( nodesFromP0ToP1 ? theNodes.begin()->second : theNodes.rbegin()->second );
741 int i = 0, ir = Nodes1.size()-1;
742 int * pi = nodesFromP0ToP1 ? &i : &ir;
743 itn = theNodes.begin();
744 if ( nodesFromP0ToP1 ) ++itn;
745 for ( ; i < Nodes1.size(); ++i, --ir, ++itn )
747 Nodes1[*pi] = itn->second;
749 for ( i = 0; i < Nodes1.size()-1; ++i )
751 Points.Append( gpXYZ( Nodes1[i]));
752 Pnts2d1.Append( myHelper->GetNodeUV( F, Nodes1[i]));
757 int edgeID = meshDS->ShapeToIndex(LinEdge1);
760 Handle(Geom_Curve) Crv = BRep_Tool::Curve(LinEdge1,fp,lp);
761 gp_Pnt Ptmp = Crv->Value(fp);
763 if( P1.Distance(Ptmp) > Precision::Confusion() )
765 // get UV points for edge
767 BRep_Tool::UVPoints( LinEdge1, TopoDS::Face(aShape), PF, PL );
770 V2d = gp_Vec2d(PF,PL);
774 V2d = gp_Vec2d(PL,PF);
777 NC = const_cast<SMDS_MeshNode*>( VertexNode( VC, meshDS ));
780 NC = meshDS->AddNode(P0.X(), P0.Y(), P0.Z());
781 meshDS->SetNodeOnVertex(NC, vertID);
785 for(; i<myLayerPositions.size(); i++) {
786 gp_Pnt P( P0.X() + aVec.X()*myLayerPositions[i],
787 P0.Y() + aVec.Y()*myLayerPositions[i],
788 P0.Z() + aVec.Z()*myLayerPositions[i] );
790 SMDS_MeshNode * node = meshDS->AddNode(P.X(), P.Y(), P.Z());
794 param = fp + dp*(1-myLayerPositions[i]);
796 param = fp + dp*myLayerPositions[i];
797 meshDS->SetNodeOnEdge(node, edgeID, param);
798 // parameters on face
799 gp_Pnt2d P2d( PC.X() + V2d.X()*myLayerPositions[i],
800 PC.Y() + V2d.Y()*myLayerPositions[i] );
803 Nodes1[ myLayerPositions.size() ] = NF;
804 // create 1D elements on edge
805 SMDS_MeshEdge* ME = myHelper->AddEdge( NC, Nodes1[0] );
806 if(ME) meshDS->SetMeshElementOnShape(ME, edgeID);
807 for(i=1; i<Nodes1.size(); i++) {
808 ME = myHelper->AddEdge( Nodes1[i-1], Nodes1[i] );
809 if(ME) meshDS->SetMeshElementOnShape(ME, edgeID);
811 if (nbe==2 && LinEdge1.Orientation() == TopAbs_INTERNAL )
814 markEdgeAsComputedByMe( LinEdge1, aMesh.GetSubMesh( F ));
817 if ( linEdge2Computed )
819 if (!GetSortedNodesOnEdge(aMesh.GetMeshDS(),LinEdge2,true,theNodes))
820 return error("Invalid mesh on a straight edge");
822 bool nodesFromP0ToP2 = ( theNodes.rbegin()->second == NL );
823 int i = 0, ir = Nodes1.size()-1;
824 int * pi = nodesFromP0ToP2 ? &i : &ir;
825 itn = theNodes.begin();
826 if ( nodesFromP0ToP2 ) ++itn;
827 for ( ; i < Nodes2.size(); ++i, --ir, ++itn )
828 Nodes2[*pi] = itn->second;
832 int edgeID = meshDS->ShapeToIndex(LinEdge2);
833 gp_Vec aVec = gp_Vec(P0,P2);
835 Handle(Geom_Curve) Crv = BRep_Tool::Curve(LinEdge2,fp,lp);
836 gp_Pnt Ptmp = Crv->Value(fp);
838 if( P2.Distance(Ptmp) > Precision::Confusion() )
840 // get UV points for edge
842 BRep_Tool::UVPoints( LinEdge2, TopoDS::Face(aShape), PF, PL );
845 V2d = gp_Vec2d(PF,PL);
849 V2d = gp_Vec2d(PL,PF);
853 for(int i=0; i<myLayerPositions.size(); i++) {
854 gp_Pnt P( P0.X() + aVec.X()*myLayerPositions[i],
855 P0.Y() + aVec.Y()*myLayerPositions[i],
856 P0.Z() + aVec.Z()*myLayerPositions[i] );
857 SMDS_MeshNode * node = meshDS->AddNode(P.X(), P.Y(), P.Z());
861 param = fp + dp*(1-myLayerPositions[i]);
863 param = fp + dp*myLayerPositions[i];
864 meshDS->SetNodeOnEdge(node, edgeID, param);
865 // parameters on face
866 gp_Pnt2d P2d( PC.X() + V2d.X()*myLayerPositions[i],
867 PC.Y() + V2d.Y()*myLayerPositions[i] );
869 Nodes2[ myLayerPositions.size() ] = NL;
870 // create 1D elements on edge
871 SMDS_MeshEdge* ME = myHelper->AddEdge( NC, Nodes2[0] );
872 if(ME) meshDS->SetMeshElementOnShape(ME, edgeID);
873 for(int i=1; i<Nodes2.size(); i++) {
874 ME = myHelper->AddEdge( Nodes2[i-1], Nodes2[i] );
875 if(ME) meshDS->SetMeshElementOnShape(ME, edgeID);
878 markEdgeAsComputedByMe( LinEdge2, aMesh.GetSubMesh( F ));
880 markEdgeAsComputedByMe( CircEdge, aMesh.GetSubMesh( F ));
883 bool IsForward = ( CircEdge.Orientation()==TopAbs_FORWARD );
884 const double angleSign = ( F.Orientation() == TopAbs_REVERSED ? -1.0 : 1.0 );
886 // create nodes and mesh elements on face
887 // find axis of rotation
888 gp_Pnt P2 = gp_Pnt( CNodes[1]->X(), CNodes[1]->Y(), CNodes[1]->Z() );
891 gp_Vec Axis = Vec1.Crossed(Vec2);
894 //cout<<"Angles.Length() = "<<Angles.Length()<<" Points.Length() = "<<Points.Length()<<endl;
895 //cout<<"Nodes1.size() = "<<Nodes1.size()<<" Pnts2d1.Length() = "<<Pnts2d1.Length()<<endl;
896 for(; i<Angles.Length(); i++) {
897 vector< const SMDS_MeshNode* > tmpNodes;
898 tmpNodes.reserve(Nodes1.size());
900 gp_Ax1 theAxis(P0,gp_Dir(Axis));
901 aTrsf.SetRotation( theAxis, Angles.Value(i) );
903 aTrsf2d.SetRotation( PC, Angles.Value(i) * angleSign );
906 for(; j<=Points.Length(); j++) {
908 Points.Value(j).Coord( cx, cy, cz );
909 aTrsf.Transforms( cx, cy, cz );
910 SMDS_MeshNode* node = myHelper->AddNode( cx, cy, cz );
911 // find parameters on face
912 Pnts2d1.Value(j).Coord( cx, cy );
913 aTrsf2d.Transforms( cx, cy );
915 meshDS->SetNodeOnFace( node, faceID, cx, cy );
916 tmpNodes[j-1] = node;
919 tmpNodes[Points.Length()] = CNodes[i];
921 for(j=0; j<Nodes1.size()-1; j++) {
924 MF = myHelper->AddFace( tmpNodes[j], Nodes1[j],
925 Nodes1[j+1], tmpNodes[j+1] );
927 MF = myHelper->AddFace( tmpNodes[j], tmpNodes[j+1],
928 Nodes1[j+1], Nodes1[j] );
929 if(MF) meshDS->SetMeshElementOnShape(MF, faceID);
934 MF = myHelper->AddFace( NC, Nodes1[0], tmpNodes[0] );
936 MF = myHelper->AddFace( NC, tmpNodes[0], Nodes1[0] );
937 if(MF) meshDS->SetMeshElementOnShape(MF, faceID);
938 for(j=0; j<Nodes1.size(); j++) {
939 Nodes1[j] = tmpNodes[j];
944 for(i=0; i<Nodes1.size()-1; i++) {
947 MF = myHelper->AddFace( Nodes2[i], Nodes1[i],
948 Nodes1[i+1], Nodes2[i+1] );
950 MF = myHelper->AddFace( Nodes2[i], Nodes2[i+1],
951 Nodes1[i+1], Nodes1[i] );
952 if(MF) meshDS->SetMeshElementOnShape(MF, faceID);
957 MF = myHelper->AddFace( NC, Nodes1[0], Nodes2[0] );
959 MF = myHelper->AddFace( NC, Nodes2[0], Nodes1[0] );
960 if(MF) meshDS->SetMeshElementOnShape(MF, faceID);
965 //================================================================================
967 * \brief Compute positions of nodes on the radial edge
968 * \retval bool - is a success
970 //================================================================================
972 bool StdMeshers_RadialQuadrangle_1D2D::computeLayerPositions(const gp_Pnt& p1,
974 const TopoDS_Edge& linEdge,
975 bool* linEdgeComputed)
977 // First, try to compute positions of layers
979 myLayerPositions.clear();
981 SMESH_Mesh * mesh = myHelper->GetMesh();
983 const SMESH_Hypothesis* hyp1D = myDistributionHypo ? myDistributionHypo->GetLayerDistribution() : 0;
984 int nbLayers = myNbLayerHypo ? myNbLayerHypo->GetNumberOfLayers() : 0;
986 if ( !hyp1D && !nbLayers )
988 // No own algo hypotheses assigned, so first try to find any 1D hypothesis.
990 TopoDS_Shape edge = linEdge;
991 if ( edge.IsNull() && !myHelper->GetSubShape().IsNull())
992 for ( TopExp_Explorer e(myHelper->GetSubShape(), TopAbs_EDGE); e.More(); e.Next())
994 if ( !edge.IsNull() )
996 // find a hyp usable by TNodeDistributor
997 const SMESH_HypoFilter* hypKind =
998 TNodeDistributor::GetDistributor(*mesh)->GetCompatibleHypoFilter(/*ignoreAux=*/true);
999 hyp1D = mesh->GetHypothesis( edge, *hypKind, /*fromAncestors=*/true);
1002 if ( hyp1D ) // try to compute with hyp1D
1004 if ( !TNodeDistributor::GetDistributor(*mesh)->Compute( myLayerPositions,p1,p2,*mesh,hyp1D )) {
1005 if ( myDistributionHypo ) { // bad hyp assigned
1006 return error( TNodeDistributor::GetDistributor(*mesh)->GetComputeError() );
1009 // bad hyp found, its Ok, lets try with default nb of segnents
1014 if ( myLayerPositions.empty() ) // try to use nb of layers
1017 nbLayers = _gen->GetDefaultNbSegments();
1021 myLayerPositions.resize( nbLayers - 1 );
1022 for ( int z = 1; z < nbLayers; ++z )
1023 myLayerPositions[ z - 1 ] = double( z )/ double( nbLayers );
1027 // Second, check presence of a mesh built by other algo on linEdge
1028 // and mesh conformity to my hypothesis
1030 bool meshComputed = (!linEdge.IsNull() && !mesh->GetSubMesh(linEdge)->IsEmpty() );
1031 if ( linEdgeComputed ) *linEdgeComputed = meshComputed;
1035 vector< double > nodeParams;
1036 GetNodeParamOnEdge( mesh->GetMeshDS(), linEdge, nodeParams );
1038 // nb of present nodes must be different in cases of 1 and 2 straight edges
1040 TopoDS_Vertex VV[2];
1041 TopExp::Vertices( linEdge, VV[0], VV[1]);
1042 const gp_Pnt* points[] = { &p1, &p2 };
1043 gp_Pnt vPoints[] = { BRep_Tool::Pnt(VV[0]), BRep_Tool::Pnt(VV[1]) };
1044 const double tol[] = { BRep_Tool::Tolerance(VV[0]), BRep_Tool::Tolerance(VV[1]) };
1045 bool pointsAreOnVertices = true;
1046 for ( int iP = 0; iP < 2 && pointsAreOnVertices; ++iP )
1047 pointsAreOnVertices = ( points[iP]->Distance( vPoints[0] ) < tol[0] ||
1048 points[iP]->Distance( vPoints[1] ) < tol[1] );
1050 int nbNodes = nodeParams.size() - 2; // 2 straight edges
1051 if ( !pointsAreOnVertices )
1052 nbNodes = ( nodeParams.size() - 3 ) / 2; // 1 straight edge
1054 if ( myLayerPositions.empty() )
1056 myLayerPositions.resize( nbNodes );
1058 else if ( myDistributionHypo || myNbLayerHypo )
1060 // linEdge is computed by other algo. Check if there is a meshed face
1061 // using nodes on linEdge
1062 bool nodesAreUsed = false;
1063 TopTools_ListIteratorOfListOfShape ancestIt = mesh->GetAncestors( linEdge );
1064 for ( ; ancestIt.More() && !nodesAreUsed; ancestIt.Next() )
1065 if ( ancestIt.Value().ShapeType() == TopAbs_FACE )
1066 nodesAreUsed = (!mesh->GetSubMesh( ancestIt.Value() )->IsEmpty());
1067 if ( !nodesAreUsed ) {
1069 mesh->GetSubMesh( linEdge )->ComputeStateEngine( SMESH_subMesh::CLEAN );
1070 if ( linEdgeComputed ) *linEdgeComputed = false;
1074 if ( myLayerPositions.size() != nbNodes )
1075 return error("Radial edge is meshed by other algorithm");
1080 return !myLayerPositions.empty();
1084 //=======================================================================
1085 //function : Evaluate
1087 //=======================================================================
1089 bool StdMeshers_RadialQuadrangle_1D2D::Evaluate(SMESH_Mesh& aMesh,
1090 const TopoDS_Shape& aShape,
1091 MapShapeNbElems& aResMap)
1093 if( aShape.ShapeType() != TopAbs_FACE ) {
1096 SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
1097 if( aResMap.count(sm) )
1100 vector<int>& aResVec =
1101 aResMap.insert( make_pair(sm, vector<int>(SMDSEntity_Last,0))).first->second;
1103 myHelper = new SMESH_MesherHelper( aMesh );
1104 myHelper->SetSubShape( aShape );
1105 auto_ptr<SMESH_MesherHelper> helperDeleter( myHelper );
1107 TNodeDistributor* algo1d = TNodeDistributor::GetDistributor(aMesh);
1109 TopoDS_Edge CircEdge, LinEdge1, LinEdge2;
1110 int nbe = analyseFace( aShape, CircEdge, LinEdge1, LinEdge2 );
1111 if( nbe>3 || nbe < 1 || CircEdge.IsNull() )
1114 Handle(Geom_Circle) aCirc = Handle(Geom_Circle)::DownCast( getCurve( CircEdge ));
1115 if( aCirc.IsNull() )
1116 return error(COMPERR_BAD_SHAPE);
1118 gp_Pnt P0 = aCirc->Location();
1119 gp_Pnt P1 = aCirc->Value(0.);
1120 computeLayerPositions( P0, P1, LinEdge1 );
1122 int nb0d=0, nb2d_tria=0, nb2d_quad=0;
1123 bool isQuadratic = false, ok = true;
1126 // C1 must be a circle
1127 ok = algo1d->EvaluateCircularEdge( aMesh, CircEdge, aResMap );
1129 const vector<int>& aVec = aResMap[aMesh.GetSubMesh(CircEdge)];
1130 isQuadratic = aVec[SMDSEntity_Quad_Edge]>aVec[SMDSEntity_Edge];
1133 nb0d = (aVec[SMDSEntity_Node]+1) * myLayerPositions.size();
1134 // radial medium nodes
1135 nb0d += (aVec[SMDSEntity_Node]+1) * (myLayerPositions.size()+1);
1136 // other medium nodes
1137 nb0d += (aVec[SMDSEntity_Node]+1) * myLayerPositions.size();
1140 nb0d = (aVec[SMDSEntity_Node]+1) * myLayerPositions.size();
1142 nb2d_tria = aVec[SMDSEntity_Node] + 1;
1146 else if(nbe==2 && LinEdge1.Orientation() != TopAbs_INTERNAL)
1148 // one curve must be a half of circle and other curve must be
1149 // a segment of line
1151 Handle(Geom_Circle) aCirc = Handle(Geom_Circle)::DownCast( getCurve( CircEdge, &fp, &lp ));
1152 if( fabs(fabs(lp-fp)-M_PI) > Precision::Confusion() ) {
1153 // not half of circle
1154 return error(COMPERR_BAD_SHAPE);
1156 Handle(Geom_Line) aLine = Handle(Geom_Line)::DownCast( getCurve( LinEdge1 ));
1157 if( aLine.IsNull() ) {
1158 // other curve not line
1159 return error(COMPERR_BAD_SHAPE);
1161 ok = !aResMap.count( aMesh.GetSubMesh(LinEdge1) );
1163 const vector<int>& aVec = aResMap[ aMesh.GetSubMesh(LinEdge1) ];
1164 ok = ( aVec[SMDSEntity_Node] == myLayerPositions.size() );
1167 ok = algo1d->EvaluateCircularEdge( aMesh, CircEdge, aResMap );
1170 const vector<int>& aVec = aResMap[ aMesh.GetSubMesh(CircEdge) ];
1171 isQuadratic = aVec[SMDSEntity_Quad_Edge] > aVec[SMDSEntity_Edge];
1174 nb0d = aVec[SMDSEntity_Node] * myLayerPositions.size();
1175 // radial medium nodes
1176 nb0d += aVec[SMDSEntity_Node] * (myLayerPositions.size()+1);
1177 // other medium nodes
1178 nb0d += (aVec[SMDSEntity_Node]+1) * myLayerPositions.size();
1181 nb0d = aVec[SMDSEntity_Node] * myLayerPositions.size();
1183 nb2d_tria = aVec[SMDSEntity_Node] + 1;
1184 nb2d_quad = nb2d_tria * myLayerPositions.size();
1185 // add evaluation for edges
1186 vector<int> aResVec(SMDSEntity_Last,0);
1188 aResVec[SMDSEntity_Node] = 4*myLayerPositions.size() + 3;
1189 aResVec[SMDSEntity_Quad_Edge] = 2*myLayerPositions.size() + 2;
1192 aResVec[SMDSEntity_Node] = 2*myLayerPositions.size() + 1;
1193 aResVec[SMDSEntity_Edge] = 2*myLayerPositions.size() + 2;
1195 aResMap[ aMesh.GetSubMesh(LinEdge1) ] = aResVec;
1198 else // nbe==3 or ( nbe==2 && linEdge is INTERNAL )
1200 if (nbe==2 && LinEdge1.Orientation() == TopAbs_INTERNAL )
1201 LinEdge2 = LinEdge1;
1203 // one curve must be a part of circle and other curves must be
1205 Handle(Geom_Line) aLine1 = Handle(Geom_Line)::DownCast( getCurve( LinEdge1 ));
1206 Handle(Geom_Line) aLine2 = Handle(Geom_Line)::DownCast( getCurve( LinEdge2 ));
1207 if( aLine1.IsNull() || aLine2.IsNull() ) {
1208 // other curve not line
1209 return error(COMPERR_BAD_SHAPE);
1211 int nbLayers = myLayerPositions.size();
1212 computeLayerPositions( P0, P1, LinEdge2 );
1213 if ( nbLayers != myLayerPositions.size() )
1214 return error("Different hypotheses apply to radial edges");
1216 bool ok = !aResMap.count( aMesh.GetSubMesh(LinEdge1));
1218 if ( myDistributionHypo || myNbLayerHypo )
1219 ok = true; // override other 1d hyps
1221 const vector<int>& aVec = aResMap[ aMesh.GetSubMesh(LinEdge1) ];
1222 ok = ( aVec[SMDSEntity_Node] == myLayerPositions.size() );
1225 if( ok && aResMap.count( aMesh.GetSubMesh(LinEdge2) )) {
1226 if ( myDistributionHypo || myNbLayerHypo )
1227 ok = true; // override other 1d hyps
1229 const vector<int>& aVec = aResMap[ aMesh.GetSubMesh(LinEdge2) ];
1230 ok = ( aVec[SMDSEntity_Node] == myLayerPositions.size() );
1234 ok = algo1d->EvaluateCircularEdge( aMesh, CircEdge, aResMap );
1237 const vector<int>& aVec = aResMap[ aMesh.GetSubMesh(CircEdge) ];
1238 isQuadratic = aVec[SMDSEntity_Quad_Edge]>aVec[SMDSEntity_Edge];
1241 nb0d = aVec[SMDSEntity_Node] * myLayerPositions.size();
1242 // radial medium nodes
1243 nb0d += aVec[SMDSEntity_Node] * (myLayerPositions.size()+1);
1244 // other medium nodes
1245 nb0d += (aVec[SMDSEntity_Node]+1) * myLayerPositions.size();
1248 nb0d = aVec[SMDSEntity_Node] * myLayerPositions.size();
1250 nb2d_tria = aVec[SMDSEntity_Node] + 1;
1251 nb2d_quad = nb2d_tria * myLayerPositions.size();
1252 // add evaluation for edges
1253 vector<int> aResVec(SMDSEntity_Last, 0);
1255 aResVec[SMDSEntity_Node] = 2*myLayerPositions.size() + 1;
1256 aResVec[SMDSEntity_Quad_Edge] = myLayerPositions.size() + 1;
1259 aResVec[SMDSEntity_Node] = myLayerPositions.size();
1260 aResVec[SMDSEntity_Edge] = myLayerPositions.size() + 1;
1262 sm = aMesh.GetSubMesh(LinEdge1);
1263 aResMap[sm] = aResVec;
1264 sm = aMesh.GetSubMesh(LinEdge2);
1265 aResMap[sm] = aResVec;
1272 aResVec[SMDSEntity_Quad_Triangle] = nb2d_tria;
1273 aResVec[SMDSEntity_Quad_Quadrangle] = nb2d_quad;
1276 aResVec[SMDSEntity_Triangle] = nb2d_tria;
1277 aResVec[SMDSEntity_Quadrangle] = nb2d_quad;
1283 sm = aMesh.GetSubMesh(aShape);
1284 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1285 smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,
1286 "Submesh can not be evaluated",this));
1291 //================================================================================
1293 * \brief Return true if applied compute mesh on this shape
1295 //================================================================================
1297 bool StdMeshers_RadialQuadrangle_1D2D::IsApplicable( const TopoDS_Shape & aShape, bool toCheckAll )
1299 int nbFoundFaces = 0;
1300 for (TopExp_Explorer exp( aShape, TopAbs_FACE ); exp.More(); exp.Next(), ++nbFoundFaces ){
1301 TopoDS_Edge CircEdge, LinEdge1, LinEdge2;
1302 int nbe = analyseFace( exp.Current(), CircEdge, LinEdge1, LinEdge2 );
1303 Handle(Geom_Circle) aCirc = Handle(Geom_Circle)::DownCast( getCurve( CircEdge ));
1304 bool ok = ( nbe <= 3 && nbe >= 1 && !aCirc.IsNull() );
1305 if( toCheckAll && !ok ) return false;
1306 if( !toCheckAll && ok ) return true;
1308 if( toCheckAll && nbFoundFaces != 0 ) return true;