1 // Copyright (C) 2007-2012 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.
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
23 // Created : Fri Oct 20 11:37:07 2006
24 // Author : Edward AGAPOV (eap)
26 #include "StdMeshers_RadialQuadrangle_1D2D.hxx"
28 #include "StdMeshers_NumberOfLayers.hxx"
29 #include "StdMeshers_LayerDistribution.hxx"
30 #include "StdMeshers_Regular_1D.hxx"
31 #include "StdMeshers_NumberOfSegments.hxx"
33 #include "SMDS_MeshNode.hxx"
34 #include "SMESHDS_SubMesh.hxx"
35 #include "SMESH_Gen.hxx"
36 #include "SMESH_HypoFilter.hxx"
37 #include "SMESH_Mesh.hxx"
38 #include "SMESH_MesherHelper.hxx"
39 #include "SMESH_subMesh.hxx"
40 #include "SMESH_subMeshEventListener.hxx"
42 #include "utilities.h"
44 #include <BRepAdaptor_Curve.hxx>
45 #include <BRepBuilderAPI_MakeEdge.hxx>
46 #include <BRep_Tool.hxx>
47 #include <GeomAPI_ProjectPointOnSurf.hxx>
48 #include <Geom_Circle.hxx>
49 #include <Geom_Line.hxx>
50 #include <Geom_TrimmedCurve.hxx>
51 #include <TColgp_SequenceOfPnt.hxx>
52 #include <TColgp_SequenceOfPnt2d.hxx>
54 #include <TopExp_Explorer.hxx>
55 #include <TopTools_ListIteratorOfListOfShape.hxx>
61 #define RETURN_BAD_RESULT(msg) { MESSAGE(")-: Error: " << msg); return false; }
62 #define gpXYZ(n) gp_XYZ(n->X(),n->Y(),n->Z())
65 //=======================================================================
66 //function : StdMeshers_RadialQuadrangle_1D2D
68 //=======================================================================
70 StdMeshers_RadialQuadrangle_1D2D::StdMeshers_RadialQuadrangle_1D2D(int hypId,
73 :SMESH_2D_Algo(hypId, studyId, gen)
75 _name = "RadialQuadrangle_1D2D";
76 _shapeType = (1 << TopAbs_FACE); // 1 bit per shape type
78 _compatibleHypothesis.push_back("LayerDistribution2D");
79 _compatibleHypothesis.push_back("NumberOfLayers2D");
81 myDistributionHypo = 0;
82 _requireDiscreteBoundary = false;
83 _supportSubmeshes = true;
87 //================================================================================
91 //================================================================================
93 StdMeshers_RadialQuadrangle_1D2D::~StdMeshers_RadialQuadrangle_1D2D()
97 //=======================================================================
98 //function : CheckHypothesis
100 //=======================================================================
102 bool StdMeshers_RadialQuadrangle_1D2D::CheckHypothesis
104 const TopoDS_Shape& aShape,
105 SMESH_Hypothesis::Hypothesis_Status& aStatus)
109 myDistributionHypo = 0;
111 list <const SMESHDS_Hypothesis * >::const_iterator itl;
113 const list <const SMESHDS_Hypothesis * >&hyps = GetUsedHypothesis(aMesh, aShape);
114 if ( hyps.size() == 0 ) {
115 aStatus = SMESH_Hypothesis::HYP_OK;
116 return true; // can work with no hypothesis
119 if ( hyps.size() > 1 ) {
120 aStatus = SMESH_Hypothesis::HYP_ALREADY_EXIST;
124 const SMESHDS_Hypothesis *theHyp = hyps.front();
126 string hypName = theHyp->GetName();
128 if (hypName == "NumberOfLayers2D") {
129 myNbLayerHypo = static_cast<const StdMeshers_NumberOfLayers *>(theHyp);
130 aStatus = SMESH_Hypothesis::HYP_OK;
133 if (hypName == "LayerDistribution2D") {
134 myDistributionHypo = static_cast<const StdMeshers_LayerDistribution *>(theHyp);
135 aStatus = SMESH_Hypothesis::HYP_OK;
138 aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
144 // ------------------------------------------------------------------------------
146 * \brief Listener used to mark edges meshed by StdMeshers_RadialQuadrangle_1D2D
148 class TEdgeMarker : public SMESH_subMeshEventListener
150 TEdgeMarker(): SMESH_subMeshEventListener(/*isDeletable=*/false,
151 "StdMeshers_RadialQuadrangle_1D2D::TEdgeMarker") {}
153 //!< Return static listener
154 static SMESH_subMeshEventListener* getListener()
156 static TEdgeMarker theEdgeMarker;
157 return &theEdgeMarker;
159 //! Clear face sumbesh if something happens on edges
160 void ProcessEvent(const int event,
162 SMESH_subMesh* edgeSubMesh,
163 EventListenerData* data,
164 const SMESH_Hypothesis* /*hyp*/)
166 if ( data && !data->mySubMeshes.empty() && eventType == SMESH_subMesh::ALGO_EVENT)
168 ASSERT( data->mySubMeshes.front() != edgeSubMesh );
169 SMESH_subMesh* faceSubMesh = data->mySubMeshes.front();
170 faceSubMesh->ComputeStateEngine( SMESH_subMesh::CLEAN );
175 // ------------------------------------------------------------------------------
177 * \brief Mark an edge as computed by StdMeshers_RadialQuadrangle_1D2D
179 void markEdgeAsComputedByMe(const TopoDS_Edge& edge, SMESH_subMesh* faceSubMesh)
181 if ( SMESH_subMesh* edgeSM = faceSubMesh->GetFather()->GetSubMeshContaining( edge ))
183 if ( !edgeSM->GetEventListenerData( TEdgeMarker::getListener() ))
184 faceSubMesh->SetEventListener( TEdgeMarker::getListener(),
185 SMESH_subMeshEventListenerData::MakeData(faceSubMesh),
189 // ------------------------------------------------------------------------------
191 * \brief Return true if a radial edge was meshed with StdMeshers_RadialQuadrangle_1D2D with
192 * the same radial distribution
194 // bool isEdgeCompatiballyMeshed(const TopoDS_Edge& edge, SMESH_subMesh* faceSubMesh)
196 // if ( SMESH_subMesh* edgeSM = faceSubMesh->GetFather()->GetSubMeshContaining( edge ))
198 // if ( SMESH_subMeshEventListenerData* otherFaceData =
199 // edgeSM->GetEventListenerData( TEdgeMarker::getListener() ))
201 // // compare hypothesis aplied to two disk faces sharing radial edges
202 // SMESH_Mesh& mesh = *faceSubMesh->GetFather();
203 // SMESH_Algo* radialQuadAlgo = mesh.GetGen()->GetAlgo(mesh, faceSubMesh->GetSubShape() );
204 // SMESH_subMesh* otherFaceSubMesh = otherFaceData->mySubMeshes.front();
205 // list <const SMESHDS_Hypothesis *> hyps1 =
206 // radialQuadAlgo->GetUsedHypothesis( mesh, faceSubMesh->GetSubShape());
207 // list <const SMESHDS_Hypothesis *> hyps2 =
208 // radialQuadAlgo->GetUsedHypothesis( mesh, otherFaceSubMesh->GetSubShape());
209 // if( hyps1.empty() && hyps2.empty() )
210 // return true; // defaul hyps
211 // if ( hyps1.size() != hyps2.size() )
213 // return *hyps1.front() == *hyps2.front();
219 //================================================================================
221 * \brief Return base curve of the edge and extremum parameters
223 //================================================================================
225 Handle(Geom_Curve) getCurve(const TopoDS_Edge& edge, double* f=0, double* l=0)
227 Handle(Geom_Curve) C;
228 if ( !edge.IsNull() )
230 double first = 0., last = 0.;
231 C = BRep_Tool::Curve(edge, first, last);
234 Handle(Geom_TrimmedCurve) tc = Handle(Geom_TrimmedCurve)::DownCast(C);
235 while( !tc.IsNull() ) {
236 C = tc->BasisCurve();
237 tc = Handle(Geom_TrimmedCurve)::DownCast(C);
246 //================================================================================
248 * \brief Return edges of the face
249 * \retval int - nb of edges
251 //================================================================================
253 int analyseFace(const TopoDS_Shape& face,
254 TopoDS_Edge& CircEdge,
255 TopoDS_Edge& LinEdge1,
256 TopoDS_Edge& LinEdge2)
258 CircEdge.Nullify(); LinEdge1.Nullify(); LinEdge2.Nullify();
261 for ( TopExp_Explorer exp( face, TopAbs_EDGE ); exp.More(); exp.Next(), ++nbe )
263 const TopoDS_Edge& E = TopoDS::Edge( exp.Current() );
265 Handle(Geom_Curve) C = getCurve(E,&f,&l);
268 if ( C->IsKind( STANDARD_TYPE(Geom_Circle)))
270 if ( CircEdge.IsNull() )
275 else if ( LinEdge1.IsNull() )
284 //================================================================================
285 //================================================================================
287 * \brief Class computing layers distribution using data of
288 * StdMeshers_LayerDistribution hypothesis
290 //================================================================================
291 //================================================================================
293 class TNodeDistributor: public StdMeshers_Regular_1D
295 list <const SMESHDS_Hypothesis *> myUsedHyps;
297 // -----------------------------------------------------------------------------
298 static TNodeDistributor* GetDistributor(SMESH_Mesh& aMesh)
300 const int myID = -1000;
301 map < int, SMESH_1D_Algo * > & algoMap = aMesh.GetGen()->_map1D_Algo;
302 map < int, SMESH_1D_Algo * >::iterator id_algo = algoMap.find( myID );
303 if ( id_algo == algoMap.end() )
304 return new TNodeDistributor( myID, 0, aMesh.GetGen() );
305 return static_cast< TNodeDistributor* >( id_algo->second );
307 // -----------------------------------------------------------------------------
308 //! Computes distribution of nodes on a straight line ending at pIn and pOut
309 bool Compute( vector< double > & positions,
313 const SMESH_Hypothesis* hyp1d)
315 if ( !hyp1d ) return error( "Invalid LayerDistribution hypothesis");
317 double len = pIn.Distance( pOut );
318 if ( len <= DBL_MIN ) return error("Too close points of inner and outer shells");
321 myUsedHyps.push_back( hyp1d );
323 TopoDS_Edge edge = BRepBuilderAPI_MakeEdge( pIn, pOut );
324 SMESH_Hypothesis::Hypothesis_Status aStatus;
325 if ( !StdMeshers_Regular_1D::CheckHypothesis( aMesh, edge, aStatus ))
326 return error( "StdMeshers_Regular_1D::CheckHypothesis() failed "
327 "with LayerDistribution hypothesis");
329 BRepAdaptor_Curve C3D(edge);
330 double f = C3D.FirstParameter(), l = C3D.LastParameter();
331 list< double > params;
332 if ( !StdMeshers_Regular_1D::computeInternalParameters( aMesh, C3D, len, f, l, params, false ))
333 return error("StdMeshers_Regular_1D failed to compute layers distribution");
336 positions.reserve( params.size() );
337 for (list<double>::iterator itU = params.begin(); itU != params.end(); itU++)
338 positions.push_back( *itU / len );
341 // -----------------------------------------------------------------------------
342 //! Make mesh on an adge using assigned 1d hyp or defaut nb of segments
343 bool ComputeCircularEdge(SMESH_Mesh& aMesh,
344 const TopoDS_Edge& anEdge)
346 _gen->Compute( aMesh, anEdge);
347 SMESH_subMesh *sm = aMesh.GetSubMesh(anEdge);
348 if ( sm->GetComputeState() != SMESH_subMesh::COMPUTE_OK)
350 // find any 1d hyp assigned (there can be a hyp w/o algo)
351 myUsedHyps = SMESH_Algo::GetUsedHypothesis(aMesh, anEdge, /*ignoreAux=*/true);
352 Hypothesis_Status aStatus;
353 if ( !StdMeshers_Regular_1D::CheckHypothesis( aMesh, anEdge, aStatus ))
355 // no valid 1d hyp assigned, use default nb of segments
356 _hypType = NB_SEGMENTS;
357 _ivalue[ DISTR_TYPE_IND ] = StdMeshers_NumberOfSegments::DT_Regular;
358 _ivalue[ NB_SEGMENTS_IND ] = _gen->GetDefaultNbSegments();
360 return StdMeshers_Regular_1D::Compute( aMesh, anEdge );
364 // -----------------------------------------------------------------------------
365 //! Make mesh on an adge using assigned 1d hyp or defaut nb of segments
366 bool EvaluateCircularEdge(SMESH_Mesh& aMesh,
367 const TopoDS_Edge& anEdge,
368 MapShapeNbElems& aResMap)
370 _gen->Evaluate( aMesh, anEdge, aResMap );
371 if ( aResMap.count( aMesh.GetSubMesh( anEdge )))
374 // find any 1d hyp assigned
375 myUsedHyps = SMESH_Algo::GetUsedHypothesis(aMesh, anEdge, /*ignoreAux=*/true);
376 Hypothesis_Status aStatus;
377 if ( !StdMeshers_Regular_1D::CheckHypothesis( aMesh, anEdge, aStatus ))
379 // no valid 1d hyp assigned, use default nb of segments
380 _hypType = NB_SEGMENTS;
381 _ivalue[ DISTR_TYPE_IND ] = StdMeshers_NumberOfSegments::DT_Regular;
382 _ivalue[ NB_SEGMENTS_IND ] = _gen->GetDefaultNbSegments();
384 return StdMeshers_Regular_1D::Evaluate( aMesh, anEdge, aResMap );
387 // -----------------------------------------------------------------------------
388 TNodeDistributor( int hypId, int studyId, SMESH_Gen* gen)
389 : StdMeshers_Regular_1D( hypId, studyId, gen)
392 // -----------------------------------------------------------------------------
393 virtual const list <const SMESHDS_Hypothesis *> &
394 GetUsedHypothesis(SMESH_Mesh &, const TopoDS_Shape &, const bool)
398 // -----------------------------------------------------------------------------
402 //=======================================================================
404 * \brief Allow algo to do something after persistent restoration
405 * \param subMesh - restored submesh
407 * call markEdgeAsComputedByMe()
409 //=======================================================================
411 void StdMeshers_RadialQuadrangle_1D2D::SubmeshRestored(SMESH_subMesh* faceSubMesh)
413 if ( !faceSubMesh->IsEmpty() )
415 TopoDS_Edge CircEdge, LinEdge1, LinEdge2;
416 analyseFace( faceSubMesh->GetSubShape(), CircEdge, LinEdge1, LinEdge2 );
417 if ( !CircEdge.IsNull() ) markEdgeAsComputedByMe( CircEdge, faceSubMesh );
418 if ( !LinEdge1.IsNull() ) markEdgeAsComputedByMe( LinEdge1, faceSubMesh );
419 if ( !LinEdge2.IsNull() ) markEdgeAsComputedByMe( LinEdge2, faceSubMesh );
423 //=======================================================================
426 //=======================================================================
428 bool StdMeshers_RadialQuadrangle_1D2D::Compute(SMESH_Mesh& aMesh,
429 const TopoDS_Shape& aShape)
431 SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
433 myHelper = new SMESH_MesherHelper( aMesh );
434 myHelper->IsQuadraticSubMesh( aShape );
435 // to delete helper at exit from Compute()
436 auto_ptr<SMESH_MesherHelper> helperDeleter( myHelper );
438 TNodeDistributor* algo1d = TNodeDistributor::GetDistributor(aMesh);
440 TopoDS_Edge CircEdge, LinEdge1, LinEdge2;
441 int nbe = analyseFace( aShape, CircEdge, LinEdge1, LinEdge2 );
442 Handle(Geom_Circle) aCirc = Handle(Geom_Circle)::DownCast( getCurve( CircEdge ));
443 if( nbe>3 || nbe < 1 || aCirc.IsNull() )
444 return error("The face must be a full circle or a part of circle (i.e. the number of edges is less or equal to 3 and one of them is a circle curve)");
447 // points for rotation
448 TColgp_SequenceOfPnt Points;
449 // angles for rotation
450 TColStd_SequenceOfReal Angles;
451 // Nodes1 and Nodes2 - nodes along radiuses
452 // CNodes - nodes on circle edge
453 vector< const SMDS_MeshNode* > Nodes1, Nodes2, CNodes;
455 // parameters edge nodes on face
456 TColgp_SequenceOfPnt2d Pnts2d1;
459 int faceID = meshDS->ShapeToIndex(aShape);
460 TopoDS_Face F = TopoDS::Face(aShape);
461 Handle(Geom_Surface) S = BRep_Tool::Surface(F);
466 if (!algo1d->ComputeCircularEdge( aMesh, CircEdge ))
467 return error( algo1d->GetComputeError() );
468 map< double, const SMDS_MeshNode* > theNodes;
469 if ( !GetSortedNodesOnEdge(aMesh.GetMeshDS(),CircEdge,true,theNodes))
470 return error("Circular edge is incorrectly meshed");
473 map< double, const SMDS_MeshNode* >::iterator itn = theNodes.begin();
474 const SMDS_MeshNode* NF = (*itn).second;
475 CNodes.push_back( (*itn).second );
476 double fang = (*itn).first;
477 if ( itn != theNodes.end() ) {
479 for(; itn != theNodes.end(); itn++ ) {
480 CNodes.push_back( (*itn).second );
481 double ang = (*itn).first - fang;
482 if( ang>M_PI ) ang = ang - 2.*M_PI;
483 if( ang<-M_PI ) ang = ang + 2.*M_PI;
484 Angles.Append( ang );
487 P1 = gp_Pnt( NF->X(), NF->Y(), NF->Z() );
488 P0 = aCirc->Location();
490 if ( !computeLayerPositions(P0,P1))
493 TopoDS_Vertex V1 = myHelper->IthVertex(0, CircEdge );
494 gp_Pnt2d p2dV = BRep_Tool::Parameters( V1, TopoDS::Face(aShape) );
496 NC = meshDS->AddNode(P0.X(), P0.Y(), P0.Z());
497 GeomAPI_ProjectPointOnSurf PPS(P0,S);
499 PPS.Parameters(1,U0,V0);
500 meshDS->SetNodeOnFace(NC, faceID, U0, V0);
501 PC = gp_Pnt2d(U0,V0);
504 gp_Vec2d aVec2d(PC,p2dV);
505 Nodes1.resize( myLayerPositions.size()+1 );
506 Nodes2.resize( myLayerPositions.size()+1 );
508 for(; i<myLayerPositions.size(); i++) {
509 gp_Pnt P( P0.X() + aVec.X()*myLayerPositions[i],
510 P0.Y() + aVec.Y()*myLayerPositions[i],
511 P0.Z() + aVec.Z()*myLayerPositions[i] );
513 SMDS_MeshNode * node = meshDS->AddNode(P.X(), P.Y(), P.Z());
516 double U = PC.X() + aVec2d.X()*myLayerPositions[i];
517 double V = PC.Y() + aVec2d.Y()*myLayerPositions[i];
518 meshDS->SetNodeOnFace( node, faceID, U, V );
519 Pnts2d1.Append(gp_Pnt2d(U,V));
521 Nodes1[Nodes1.size()-1] = NF;
522 Nodes2[Nodes1.size()-1] = NF;
524 else if(nbe==2 && LinEdge1.Orientation() != TopAbs_INTERNAL )
526 // one curve must be a half of circle and other curve must be
529 Handle(Geom_Circle) aCirc = Handle(Geom_Circle)::DownCast( getCurve( CircEdge, &fp, &lp ));
530 if( fabs(fabs(lp-fp)-M_PI) > Precision::Confusion() ) {
531 // not half of circle
532 return error(COMPERR_BAD_SHAPE);
534 Handle(Geom_Line) aLine = Handle(Geom_Line)::DownCast( getCurve( LinEdge1 ));
535 if( aLine.IsNull() ) {
536 // other curve not line
537 return error(COMPERR_BAD_SHAPE);
540 if ( !algo1d->ComputeCircularEdge( aMesh, CircEdge ))
541 return error( algo1d->GetComputeError() );
542 map< double, const SMDS_MeshNode* > theNodes;
543 if ( !GetSortedNodesOnEdge(aMesh.GetMeshDS(),CircEdge,true,theNodes) )
544 return error("Circular edge is incorrectly meshed");
546 map< double, const SMDS_MeshNode* >::iterator itn = theNodes.begin();
548 CNodes.push_back( itn->second );
549 double fang = (*itn).first;
551 for(; itn != theNodes.end(); itn++ ) {
552 CNodes.push_back( (*itn).second );
553 double ang = (*itn).first - fang;
554 if( ang>M_PI ) ang = ang - 2.*M_PI;
555 if( ang<-M_PI ) ang = ang + 2.*M_PI;
556 Angles.Append( ang );
558 const SMDS_MeshNode* NF = theNodes.begin()->second;
559 const SMDS_MeshNode* NL = theNodes.rbegin()->second;
560 P1 = gp_Pnt( NF->X(), NF->Y(), NF->Z() );
561 gp_Pnt P2( NL->X(), NL->Y(), NL->Z() );
562 P0 = aCirc->Location();
564 bool linEdgeComputed;
565 if ( !computeLayerPositions(P0,P1,LinEdge1,&linEdgeComputed))
568 if ( linEdgeComputed )
570 if (!GetSortedNodesOnEdge(aMesh.GetMeshDS(),LinEdge1,true,theNodes))
571 return error("Invalid mesh on a straight edge");
573 Nodes1.resize( myLayerPositions.size()+1 );
574 Nodes2.resize( myLayerPositions.size()+1 );
575 vector< const SMDS_MeshNode* > *pNodes1 = &Nodes1, *pNodes2 = &Nodes2;
576 bool nodesFromP0ToP1 = ( theNodes.rbegin()->second == NF );
577 if ( !nodesFromP0ToP1 ) std::swap( pNodes1, pNodes2 );
579 map< double, const SMDS_MeshNode* >::reverse_iterator ritn = theNodes.rbegin();
580 itn = theNodes.begin();
581 for ( int i = Nodes1.size()-1; i > -1; ++itn, ++ritn, --i )
583 (*pNodes1)[i] = ritn->second;
584 (*pNodes2)[i] = itn->second;
585 Points.Prepend( gpXYZ( Nodes1[i]));
586 Pnts2d1.Prepend( myHelper->GetNodeUV( F, Nodes1[i]));
588 NC = const_cast<SMDS_MeshNode*>( itn->second );
589 Points.Remove( Nodes1.size() );
594 int edgeID = meshDS->ShapeToIndex(LinEdge1);
596 Handle(Geom_Curve) Crv = BRep_Tool::Curve(LinEdge1,fp,lp);
600 if( P1.Distance(Ptmp) > Precision::Confusion() )
602 // get UV points for edge
604 BRep_Tool::UVPoints( LinEdge1, TopoDS::Face(aShape), PF, PL );
605 PC = gp_Pnt2d( (PF.X()+PL.X())/2, (PF.Y()+PL.Y())/2 );
607 if(ori) V2d = gp_Vec2d(PC,PF);
608 else V2d = gp_Vec2d(PC,PL);
610 double cp = (fp+lp)/2;
611 double dp2 = (lp-fp)/2;
612 NC = meshDS->AddNode(P0.X(), P0.Y(), P0.Z());
613 meshDS->SetNodeOnEdge(NC, edgeID, cp);
614 Nodes1.resize( myLayerPositions.size()+1 );
615 Nodes2.resize( myLayerPositions.size()+1 );
617 for(; i<myLayerPositions.size(); i++) {
618 gp_Pnt P( P0.X() + aVec.X()*myLayerPositions[i],
619 P0.Y() + aVec.Y()*myLayerPositions[i],
620 P0.Z() + aVec.Z()*myLayerPositions[i] );
622 SMDS_MeshNode * node = meshDS->AddNode(P.X(), P.Y(), P.Z());
626 param = fp + dp2*(1-myLayerPositions[i]);
628 param = cp + dp2*myLayerPositions[i];
629 meshDS->SetNodeOnEdge(node, edgeID, param);
630 P = gp_Pnt( P0.X() - aVec.X()*myLayerPositions[i],
631 P0.Y() - aVec.Y()*myLayerPositions[i],
632 P0.Z() - aVec.Z()*myLayerPositions[i] );
633 node = meshDS->AddNode(P.X(), P.Y(), P.Z());
636 param = fp + dp2*(1-myLayerPositions[i]);
638 param = cp + dp2*myLayerPositions[i];
639 meshDS->SetNodeOnEdge(node, edgeID, param);
640 // parameters on face
641 gp_Pnt2d P2d( PC.X() + V2d.X()*myLayerPositions[i],
642 PC.Y() + V2d.Y()*myLayerPositions[i] );
645 Nodes1[ myLayerPositions.size() ] = NF;
646 Nodes2[ myLayerPositions.size() ] = NL;
647 // create 1D elements on edge
648 vector< const SMDS_MeshNode* > tmpNodes;
649 tmpNodes.resize(2*Nodes1.size()+1);
650 for(i=0; i<Nodes2.size(); i++)
651 tmpNodes[Nodes2.size()-i-1] = Nodes2[i];
652 tmpNodes[Nodes2.size()] = NC;
653 for(i=0; i<Nodes1.size(); i++)
654 tmpNodes[Nodes2.size()+1+i] = Nodes1[i];
655 for(i=1; i<tmpNodes.size(); i++) {
656 SMDS_MeshEdge* ME = myHelper->AddEdge( tmpNodes[i-1], tmpNodes[i] );
657 if(ME) meshDS->SetMeshElementOnShape(ME, edgeID);
659 markEdgeAsComputedByMe( LinEdge1, aMesh.GetSubMesh( F ));
662 else // nbe==3 or ( nbe==2 && linEdge is INTERNAL )
664 if (nbe==2 && LinEdge1.Orientation() == TopAbs_INTERNAL )
667 // one curve must be a part of circle and other curves must be
670 Handle(Geom_Circle) aCirc = Handle(Geom_Circle)::DownCast( getCurve( CircEdge ));
671 Handle(Geom_Line) aLine1 = Handle(Geom_Line)::DownCast( getCurve( LinEdge1 ));
672 Handle(Geom_Line) aLine2 = Handle(Geom_Line)::DownCast( getCurve( LinEdge2 ));
673 if( aCirc.IsNull() || aLine1.IsNull() || aLine2.IsNull() )
674 return error(COMPERR_BAD_SHAPE);
676 if ( !algo1d->ComputeCircularEdge( aMesh, CircEdge ))
677 return error( algo1d->GetComputeError() );
678 map< double, const SMDS_MeshNode* > theNodes;
679 if ( !GetSortedNodesOnEdge(aMesh.GetMeshDS(),CircEdge,true,theNodes))
680 return error("Circular edge is incorrectly meshed");
682 const SMDS_MeshNode* NF = theNodes.begin()->second;
683 const SMDS_MeshNode* NL = theNodes.rbegin()->second;
685 CNodes.push_back( NF );
686 map< double, const SMDS_MeshNode* >::iterator itn = theNodes.begin();
687 double fang = (*itn).first;
689 for(; itn != theNodes.end(); itn++ ) {
690 CNodes.push_back( (*itn).second );
691 double ang = (*itn).first - fang;
692 if( ang>M_PI ) ang = ang - 2.*M_PI;
693 if( ang<-M_PI ) ang = ang + 2.*M_PI;
694 Angles.Append( ang );
696 P1 = gp_Pnt( NF->X(), NF->Y(), NF->Z() );
697 gp_Pnt P2( NL->X(), NL->Y(), NL->Z() );
698 P0 = aCirc->Location();
700 // make P1 belong to LinEdge1
701 TopoDS_Vertex V1 = myHelper->IthVertex( 0, LinEdge1 );
702 TopoDS_Vertex V2 = myHelper->IthVertex( 1, LinEdge1 );
703 gp_Pnt PE1 = BRep_Tool::Pnt(V1);
704 gp_Pnt PE2 = BRep_Tool::Pnt(V2);
705 if( ( P1.Distance(PE1) > Precision::Confusion() ) &&
706 ( P1.Distance(PE2) > Precision::Confusion() ) )
707 std::swap( LinEdge1, LinEdge2 );
709 bool linEdge1Computed, linEdge2Computed;
710 if ( !computeLayerPositions(P0,P1,LinEdge1,&linEdge1Computed))
713 Nodes1.resize( myLayerPositions.size()+1 );
714 Nodes2.resize( myLayerPositions.size()+1 );
716 // check that both linear edges have same hypotheses
717 if ( !computeLayerPositions(P0,P2,LinEdge2, &linEdge2Computed))
719 if ( Nodes1.size() != myLayerPositions.size()+1 )
720 return error("Different hypotheses apply to radial edges");
722 // find the central vertex
723 TopoDS_Vertex VC = V2;
724 if( ( P1.Distance(PE1) > Precision::Confusion() ) &&
725 ( P2.Distance(PE1) > Precision::Confusion() ) )
727 int vertID = meshDS->ShapeToIndex(VC);
730 if ( linEdge1Computed )
732 if (!GetSortedNodesOnEdge(aMesh.GetMeshDS(),LinEdge1,true,theNodes))
733 return error("Invalid mesh on a straight edge");
735 bool nodesFromP0ToP1 = ( theNodes.rbegin()->second == NF );
736 NC = const_cast<SMDS_MeshNode*>
737 ( nodesFromP0ToP1 ? theNodes.begin()->second : theNodes.rbegin()->second );
738 int i = 0, ir = Nodes1.size()-1;
739 int * pi = nodesFromP0ToP1 ? &i : &ir;
740 itn = theNodes.begin();
741 if ( nodesFromP0ToP1 ) ++itn;
742 for ( ; i < Nodes1.size(); ++i, --ir, ++itn )
744 Nodes1[*pi] = itn->second;
746 for ( i = 0; i < Nodes1.size()-1; ++i )
748 Points.Append( gpXYZ( Nodes1[i]));
749 Pnts2d1.Append( myHelper->GetNodeUV( F, Nodes1[i]));
754 int edgeID = meshDS->ShapeToIndex(LinEdge1);
757 Handle(Geom_Curve) Crv = BRep_Tool::Curve(LinEdge1,fp,lp);
758 gp_Pnt Ptmp = Crv->Value(fp);
760 if( P1.Distance(Ptmp) > Precision::Confusion() )
762 // get UV points for edge
764 BRep_Tool::UVPoints( LinEdge1, TopoDS::Face(aShape), PF, PL );
767 V2d = gp_Vec2d(PF,PL);
771 V2d = gp_Vec2d(PL,PF);
774 NC = const_cast<SMDS_MeshNode*>( VertexNode( VC, meshDS ));
777 NC = meshDS->AddNode(P0.X(), P0.Y(), P0.Z());
778 meshDS->SetNodeOnVertex(NC, vertID);
782 for(; i<myLayerPositions.size(); i++) {
783 gp_Pnt P( P0.X() + aVec.X()*myLayerPositions[i],
784 P0.Y() + aVec.Y()*myLayerPositions[i],
785 P0.Z() + aVec.Z()*myLayerPositions[i] );
787 SMDS_MeshNode * node = meshDS->AddNode(P.X(), P.Y(), P.Z());
791 param = fp + dp*(1-myLayerPositions[i]);
793 param = fp + dp*myLayerPositions[i];
794 meshDS->SetNodeOnEdge(node, edgeID, param);
795 // parameters on face
796 gp_Pnt2d P2d( PC.X() + V2d.X()*myLayerPositions[i],
797 PC.Y() + V2d.Y()*myLayerPositions[i] );
800 Nodes1[ myLayerPositions.size() ] = NF;
801 // create 1D elements on edge
802 SMDS_MeshEdge* ME = myHelper->AddEdge( NC, Nodes1[0] );
803 if(ME) meshDS->SetMeshElementOnShape(ME, edgeID);
804 for(i=1; i<Nodes1.size(); i++) {
805 ME = myHelper->AddEdge( Nodes1[i-1], Nodes1[i] );
806 if(ME) meshDS->SetMeshElementOnShape(ME, edgeID);
808 if (nbe==2 && LinEdge1.Orientation() == TopAbs_INTERNAL )
811 markEdgeAsComputedByMe( LinEdge1, aMesh.GetSubMesh( F ));
814 if ( linEdge2Computed )
816 if (!GetSortedNodesOnEdge(aMesh.GetMeshDS(),LinEdge2,true,theNodes))
817 return error("Invalid mesh on a straight edge");
819 bool nodesFromP0ToP2 = ( theNodes.rbegin()->second == NL );
820 int i = 0, ir = Nodes1.size()-1;
821 int * pi = nodesFromP0ToP2 ? &i : &ir;
822 itn = theNodes.begin();
823 if ( nodesFromP0ToP2 ) ++itn;
824 for ( ; i < Nodes2.size(); ++i, --ir, ++itn )
825 Nodes2[*pi] = itn->second;
829 int edgeID = meshDS->ShapeToIndex(LinEdge2);
830 gp_Vec aVec = gp_Vec(P0,P2);
832 Handle(Geom_Curve) Crv = BRep_Tool::Curve(LinEdge2,fp,lp);
833 gp_Pnt Ptmp = Crv->Value(fp);
835 if( P2.Distance(Ptmp) > Precision::Confusion() )
837 // get UV points for edge
839 BRep_Tool::UVPoints( LinEdge2, TopoDS::Face(aShape), PF, PL );
842 V2d = gp_Vec2d(PF,PL);
846 V2d = gp_Vec2d(PL,PF);
850 for(int i=0; i<myLayerPositions.size(); i++) {
851 gp_Pnt P( P0.X() + aVec.X()*myLayerPositions[i],
852 P0.Y() + aVec.Y()*myLayerPositions[i],
853 P0.Z() + aVec.Z()*myLayerPositions[i] );
854 SMDS_MeshNode * node = meshDS->AddNode(P.X(), P.Y(), P.Z());
858 param = fp + dp*(1-myLayerPositions[i]);
860 param = fp + dp*myLayerPositions[i];
861 meshDS->SetNodeOnEdge(node, edgeID, param);
862 // parameters on face
863 gp_Pnt2d P2d( PC.X() + V2d.X()*myLayerPositions[i],
864 PC.Y() + V2d.Y()*myLayerPositions[i] );
866 Nodes2[ myLayerPositions.size() ] = NL;
867 // create 1D elements on edge
868 SMDS_MeshEdge* ME = myHelper->AddEdge( NC, Nodes2[0] );
869 if(ME) meshDS->SetMeshElementOnShape(ME, edgeID);
870 for(int i=1; i<Nodes2.size(); i++) {
871 ME = myHelper->AddEdge( Nodes2[i-1], Nodes2[i] );
872 if(ME) meshDS->SetMeshElementOnShape(ME, edgeID);
875 markEdgeAsComputedByMe( LinEdge2, aMesh.GetSubMesh( F ));
877 markEdgeAsComputedByMe( CircEdge, aMesh.GetSubMesh( F ));
880 bool IsForward = ( CircEdge.Orientation()==TopAbs_FORWARD );
882 // create nodes and mesh elements on face
883 // find axis of rotation
884 gp_Pnt P2 = gp_Pnt( CNodes[1]->X(), CNodes[1]->Y(), CNodes[1]->Z() );
887 gp_Vec Axis = Vec1.Crossed(Vec2);
890 //cout<<"Angles.Length() = "<<Angles.Length()<<" Points.Length() = "<<Points.Length()<<endl;
891 //cout<<"Nodes1.size() = "<<Nodes1.size()<<" Pnts2d1.Length() = "<<Pnts2d1.Length()<<endl;
892 for(; i<Angles.Length(); i++) {
893 vector< const SMDS_MeshNode* > tmpNodes;
894 tmpNodes.reserve(Nodes1.size());
896 gp_Ax1 theAxis(P0,gp_Dir(Axis));
897 aTrsf.SetRotation( theAxis, Angles.Value(i) );
899 aTrsf2d.SetRotation( PC, Angles.Value(i) );
902 for(; j<=Points.Length(); j++) {
904 Points.Value(j).Coord( cx, cy, cz );
905 aTrsf.Transforms( cx, cy, cz );
906 SMDS_MeshNode* node = myHelper->AddNode( cx, cy, cz );
907 // find parameters on face
908 Pnts2d1.Value(j).Coord( cx, cy );
909 aTrsf2d.Transforms( cx, cy );
911 meshDS->SetNodeOnFace( node, faceID, cx, cy );
912 tmpNodes[j-1] = node;
915 tmpNodes[Points.Length()] = CNodes[i];
917 for(j=0; j<Nodes1.size()-1; j++) {
920 MF = myHelper->AddFace( tmpNodes[j], Nodes1[j],
921 Nodes1[j+1], tmpNodes[j+1] );
923 MF = myHelper->AddFace( tmpNodes[j], tmpNodes[j+1],
924 Nodes1[j+1], Nodes1[j] );
925 if(MF) meshDS->SetMeshElementOnShape(MF, faceID);
930 MF = myHelper->AddFace( NC, Nodes1[0], tmpNodes[0] );
932 MF = myHelper->AddFace( NC, tmpNodes[0], Nodes1[0] );
933 if(MF) meshDS->SetMeshElementOnShape(MF, faceID);
934 for(j=0; j<Nodes1.size(); j++) {
935 Nodes1[j] = tmpNodes[j];
940 for(i=0; i<Nodes1.size()-1; i++) {
943 MF = myHelper->AddFace( Nodes2[i], Nodes1[i],
944 Nodes1[i+1], Nodes2[i+1] );
946 MF = myHelper->AddFace( Nodes2[i], Nodes2[i+1],
947 Nodes1[i+1], Nodes1[i] );
948 if(MF) meshDS->SetMeshElementOnShape(MF, faceID);
953 MF = myHelper->AddFace( NC, Nodes1[0], Nodes2[0] );
955 MF = myHelper->AddFace( NC, Nodes2[0], Nodes1[0] );
956 if(MF) meshDS->SetMeshElementOnShape(MF, faceID);
961 //================================================================================
963 * \brief Compute positions of nodes on the radial edge
964 * \retval bool - is a success
966 //================================================================================
968 bool StdMeshers_RadialQuadrangle_1D2D::computeLayerPositions(const gp_Pnt& p1,
970 const TopoDS_Edge& linEdge,
971 bool* linEdgeComputed)
973 // First, try to compute positions of layers
975 myLayerPositions.clear();
977 SMESH_Mesh * mesh = myHelper->GetMesh();
979 const SMESH_Hypothesis* hyp1D = myDistributionHypo ? myDistributionHypo->GetLayerDistribution() : 0;
980 int nbLayers = myNbLayerHypo ? myNbLayerHypo->GetNumberOfLayers() : 0;
982 if ( !hyp1D && !nbLayers )
984 // No own algo hypotheses assigned, so first try to find any 1D hypothesis.
986 TopoDS_Shape edge = linEdge;
987 if ( edge.IsNull() && !myHelper->GetSubShape().IsNull())
988 for ( TopExp_Explorer e(myHelper->GetSubShape(), TopAbs_EDGE); e.More(); e.Next())
990 if ( !edge.IsNull() )
992 // find a hyp usable by TNodeDistributor
993 SMESH_HypoFilter hypKind;
994 TNodeDistributor::GetDistributor(*mesh)->InitCompatibleHypoFilter(hypKind,/*ignoreAux=*/1);
995 hyp1D = mesh->GetHypothesis( edge, hypKind, /*fromAncestors=*/true);
998 if ( hyp1D ) // try to compute with hyp1D
1000 if ( !TNodeDistributor::GetDistributor(*mesh)->Compute( myLayerPositions,p1,p2,*mesh,hyp1D )) {
1001 if ( myDistributionHypo ) { // bad hyp assigned
1002 return error( TNodeDistributor::GetDistributor(*mesh)->GetComputeError() );
1005 // bad hyp found, its Ok, lets try with default nb of segnents
1010 if ( myLayerPositions.empty() ) // try to use nb of layers
1013 nbLayers = _gen->GetDefaultNbSegments();
1017 myLayerPositions.resize( nbLayers - 1 );
1018 for ( int z = 1; z < nbLayers; ++z )
1019 myLayerPositions[ z - 1 ] = double( z )/ double( nbLayers );
1023 // Second, check presence of a mesh built by other algo on linEdge
1024 // and mesh conformity to my hypothesis
1026 bool meshComputed = (!linEdge.IsNull() && !mesh->GetSubMesh(linEdge)->IsEmpty() );
1027 if ( linEdgeComputed ) *linEdgeComputed = meshComputed;
1031 vector< double > nodeParams;
1032 GetNodeParamOnEdge( mesh->GetMeshDS(), linEdge, nodeParams );
1034 // nb of present nodes must be different in cases of 1 and 2 straight edges
1036 TopoDS_Vertex VV[2];
1037 TopExp::Vertices( linEdge, VV[0], VV[1]);
1038 const gp_Pnt* points[] = { &p1, &p2 };
1039 gp_Pnt vPoints[] = { BRep_Tool::Pnt(VV[0]), BRep_Tool::Pnt(VV[1]) };
1040 const double tol[] = { BRep_Tool::Tolerance(VV[0]), BRep_Tool::Tolerance(VV[1]) };
1041 bool pointsAreOnVertices = true;
1042 for ( int iP = 0; iP < 2 && pointsAreOnVertices; ++iP )
1043 pointsAreOnVertices = ( points[iP]->Distance( vPoints[0] ) < tol[0] ||
1044 points[iP]->Distance( vPoints[1] ) < tol[1] );
1046 int nbNodes = nodeParams.size() - 2; // 2 straight edges
1047 if ( !pointsAreOnVertices )
1048 nbNodes = ( nodeParams.size() - 3 ) / 2; // 1 straight edge
1050 if ( myLayerPositions.empty() )
1052 myLayerPositions.resize( nbNodes );
1054 else if ( myDistributionHypo || myNbLayerHypo )
1056 // linEdge is computed by other algo. Check if there is a meshed face
1057 // using nodes on linEdge
1058 bool nodesAreUsed = false;
1059 TopTools_ListIteratorOfListOfShape ancestIt = mesh->GetAncestors( linEdge );
1060 for ( ; ancestIt.More() && !nodesAreUsed; ancestIt.Next() )
1061 if ( ancestIt.Value().ShapeType() == TopAbs_FACE )
1062 nodesAreUsed = (!mesh->GetSubMesh( ancestIt.Value() )->IsEmpty());
1063 if ( !nodesAreUsed ) {
1065 mesh->GetSubMesh( linEdge )->ComputeStateEngine( SMESH_subMesh::CLEAN );
1066 if ( linEdgeComputed ) *linEdgeComputed = false;
1070 if ( myLayerPositions.size() != nbNodes )
1071 return error("Radial edge is meshed by other algorithm");
1076 return !myLayerPositions.empty();
1080 //=======================================================================
1081 //function : Evaluate
1083 //=======================================================================
1085 bool StdMeshers_RadialQuadrangle_1D2D::Evaluate(SMESH_Mesh& aMesh,
1086 const TopoDS_Shape& aShape,
1087 MapShapeNbElems& aResMap)
1089 if( aShape.ShapeType() != TopAbs_FACE ) {
1092 SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
1093 if( aResMap.count(sm) )
1096 vector<int>& aResVec =
1097 aResMap.insert( make_pair(sm, vector<int>(SMDSEntity_Last,0))).first->second;
1099 myHelper = new SMESH_MesherHelper( aMesh );
1100 myHelper->SetSubShape( aShape );
1101 auto_ptr<SMESH_MesherHelper> helperDeleter( myHelper );
1103 TNodeDistributor* algo1d = TNodeDistributor::GetDistributor(aMesh);
1105 TopoDS_Edge CircEdge, LinEdge1, LinEdge2;
1106 int nbe = analyseFace( aShape, CircEdge, LinEdge1, LinEdge2 );
1107 if( nbe>3 || nbe < 1 || CircEdge.IsNull() )
1110 Handle(Geom_Circle) aCirc = Handle(Geom_Circle)::DownCast( getCurve( CircEdge ));
1111 if( aCirc.IsNull() )
1112 return error(COMPERR_BAD_SHAPE);
1114 gp_Pnt P0 = aCirc->Location();
1115 gp_Pnt P1 = aCirc->Value(0.);
1116 computeLayerPositions( P0, P1, LinEdge1 );
1118 int nb0d=0, nb2d_tria=0, nb2d_quad=0;
1119 bool isQuadratic = false, ok = true;
1122 // C1 must be a circle
1123 ok = algo1d->EvaluateCircularEdge( aMesh, CircEdge, aResMap );
1125 const vector<int>& aVec = aResMap[aMesh.GetSubMesh(CircEdge)];
1126 isQuadratic = aVec[SMDSEntity_Quad_Edge]>aVec[SMDSEntity_Edge];
1129 nb0d = (aVec[SMDSEntity_Node]+1) * myLayerPositions.size();
1130 // radial medium nodes
1131 nb0d += (aVec[SMDSEntity_Node]+1) * (myLayerPositions.size()+1);
1132 // other medium nodes
1133 nb0d += (aVec[SMDSEntity_Node]+1) * myLayerPositions.size();
1136 nb0d = (aVec[SMDSEntity_Node]+1) * myLayerPositions.size();
1138 nb2d_tria = aVec[SMDSEntity_Node] + 1;
1142 else if(nbe==2 && LinEdge1.Orientation() != TopAbs_INTERNAL)
1144 // one curve must be a half of circle and other curve must be
1145 // a segment of line
1147 Handle(Geom_Circle) aCirc = Handle(Geom_Circle)::DownCast( getCurve( CircEdge, &fp, &lp ));
1148 if( fabs(fabs(lp-fp)-M_PI) > Precision::Confusion() ) {
1149 // not half of circle
1150 return error(COMPERR_BAD_SHAPE);
1152 Handle(Geom_Line) aLine = Handle(Geom_Line)::DownCast( getCurve( LinEdge1 ));
1153 if( aLine.IsNull() ) {
1154 // other curve not line
1155 return error(COMPERR_BAD_SHAPE);
1157 ok = !aResMap.count( aMesh.GetSubMesh(LinEdge1) );
1159 const vector<int>& aVec = aResMap[ aMesh.GetSubMesh(LinEdge1) ];
1160 ok = ( aVec[SMDSEntity_Node] == myLayerPositions.size() );
1163 ok = algo1d->EvaluateCircularEdge( aMesh, CircEdge, aResMap );
1166 const vector<int>& aVec = aResMap[ aMesh.GetSubMesh(CircEdge) ];
1167 isQuadratic = aVec[SMDSEntity_Quad_Edge] > aVec[SMDSEntity_Edge];
1170 nb0d = aVec[SMDSEntity_Node] * myLayerPositions.size();
1171 // radial medium nodes
1172 nb0d += aVec[SMDSEntity_Node] * (myLayerPositions.size()+1);
1173 // other medium nodes
1174 nb0d += (aVec[SMDSEntity_Node]+1) * myLayerPositions.size();
1177 nb0d = aVec[SMDSEntity_Node] * myLayerPositions.size();
1179 nb2d_tria = aVec[SMDSEntity_Node] + 1;
1180 nb2d_quad = nb2d_tria * myLayerPositions.size();
1181 // add evaluation for edges
1182 vector<int> aResVec(SMDSEntity_Last,0);
1184 aResVec[SMDSEntity_Node] = 4*myLayerPositions.size() + 3;
1185 aResVec[SMDSEntity_Quad_Edge] = 2*myLayerPositions.size() + 2;
1188 aResVec[SMDSEntity_Node] = 2*myLayerPositions.size() + 1;
1189 aResVec[SMDSEntity_Edge] = 2*myLayerPositions.size() + 2;
1191 aResMap[ aMesh.GetSubMesh(LinEdge1) ] = aResVec;
1194 else // nbe==3 or ( nbe==2 && linEdge is INTERNAL )
1196 if (nbe==2 && LinEdge1.Orientation() == TopAbs_INTERNAL )
1197 LinEdge2 = LinEdge1;
1199 // one curve must be a part of circle and other curves must be
1201 Handle(Geom_Line) aLine1 = Handle(Geom_Line)::DownCast( getCurve( LinEdge1 ));
1202 Handle(Geom_Line) aLine2 = Handle(Geom_Line)::DownCast( getCurve( LinEdge2 ));
1203 if( aLine1.IsNull() || aLine2.IsNull() ) {
1204 // other curve not line
1205 return error(COMPERR_BAD_SHAPE);
1207 int nbLayers = myLayerPositions.size();
1208 computeLayerPositions( P0, P1, LinEdge2 );
1209 if ( nbLayers != myLayerPositions.size() )
1210 return error("Different hypotheses apply to radial edges");
1212 bool ok = !aResMap.count( aMesh.GetSubMesh(LinEdge1));
1214 if ( myDistributionHypo || myNbLayerHypo )
1215 ok = true; // override other 1d hyps
1217 const vector<int>& aVec = aResMap[ aMesh.GetSubMesh(LinEdge1) ];
1218 ok = ( aVec[SMDSEntity_Node] == myLayerPositions.size() );
1221 if( ok && aResMap.count( aMesh.GetSubMesh(LinEdge2) )) {
1222 if ( myDistributionHypo || myNbLayerHypo )
1223 ok = true; // override other 1d hyps
1225 const vector<int>& aVec = aResMap[ aMesh.GetSubMesh(LinEdge2) ];
1226 ok = ( aVec[SMDSEntity_Node] == myLayerPositions.size() );
1230 ok = algo1d->EvaluateCircularEdge( aMesh, CircEdge, aResMap );
1233 const vector<int>& aVec = aResMap[ aMesh.GetSubMesh(CircEdge) ];
1234 isQuadratic = aVec[SMDSEntity_Quad_Edge]>aVec[SMDSEntity_Edge];
1237 nb0d = aVec[SMDSEntity_Node] * myLayerPositions.size();
1238 // radial medium nodes
1239 nb0d += aVec[SMDSEntity_Node] * (myLayerPositions.size()+1);
1240 // other medium nodes
1241 nb0d += (aVec[SMDSEntity_Node]+1) * myLayerPositions.size();
1244 nb0d = aVec[SMDSEntity_Node] * myLayerPositions.size();
1246 nb2d_tria = aVec[SMDSEntity_Node] + 1;
1247 nb2d_quad = nb2d_tria * myLayerPositions.size();
1248 // add evaluation for edges
1249 vector<int> aResVec(SMDSEntity_Last, 0);
1251 aResVec[SMDSEntity_Node] = 2*myLayerPositions.size() + 1;
1252 aResVec[SMDSEntity_Quad_Edge] = myLayerPositions.size() + 1;
1255 aResVec[SMDSEntity_Node] = myLayerPositions.size();
1256 aResVec[SMDSEntity_Edge] = myLayerPositions.size() + 1;
1258 sm = aMesh.GetSubMesh(LinEdge1);
1259 aResMap[sm] = aResVec;
1260 sm = aMesh.GetSubMesh(LinEdge2);
1261 aResMap[sm] = aResVec;
1268 aResVec[SMDSEntity_Quad_Triangle] = nb2d_tria;
1269 aResVec[SMDSEntity_Quad_Quadrangle] = nb2d_quad;
1272 aResVec[SMDSEntity_Triangle] = nb2d_tria;
1273 aResVec[SMDSEntity_Quadrangle] = nb2d_quad;
1279 sm = aMesh.GetSubMesh(aShape);
1280 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1281 smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,
1282 "Submesh can not be evaluated",this));