1 // Copyright (C) 2007-2011 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
19 // SMESH SMESH : implementaion of SMESH idl descriptions
20 // File : StdMeshers_RadialQuadrangle_1D2D.cxx
22 // Created : Fri Oct 20 11:37:07 2006
23 // Author : Edward AGAPOV (eap)
25 #include "StdMeshers_RadialQuadrangle_1D2D.hxx"
27 #include "StdMeshers_NumberOfLayers.hxx"
28 #include "StdMeshers_LayerDistribution.hxx"
29 #include "StdMeshers_Regular_1D.hxx"
30 #include "StdMeshers_NumberOfSegments.hxx"
32 #include "SMDS_MeshNode.hxx"
33 #include "SMESHDS_SubMesh.hxx"
34 #include "SMESH_Gen.hxx"
35 #include "SMESH_HypoFilter.hxx"
36 #include "SMESH_Mesh.hxx"
37 #include "SMESH_MesherHelper.hxx"
38 #include "SMESH_subMesh.hxx"
39 #include "SMESH_subMeshEventListener.hxx"
41 #include "utilities.h"
43 #include <BRepAdaptor_Curve.hxx>
44 #include <BRepBuilderAPI_MakeEdge.hxx>
45 #include <BRep_Tool.hxx>
46 #include <GeomAPI_ProjectPointOnSurf.hxx>
47 #include <Geom_Circle.hxx>
48 #include <Geom_Line.hxx>
49 #include <Geom_TrimmedCurve.hxx>
50 #include <TColgp_SequenceOfPnt.hxx>
51 #include <TColgp_SequenceOfPnt2d.hxx>
53 #include <TopExp_Explorer.hxx>
54 #include <TopTools_ListIteratorOfListOfShape.hxx>
60 #define RETURN_BAD_RESULT(msg) { MESSAGE(")-: Error: " << msg); return false; }
61 #define gpXYZ(n) gp_XYZ(n->X(),n->Y(),n->Z())
64 //=======================================================================
65 //function : StdMeshers_RadialQuadrangle_1D2D
67 //=======================================================================
69 StdMeshers_RadialQuadrangle_1D2D::StdMeshers_RadialQuadrangle_1D2D(int hypId,
72 :SMESH_2D_Algo(hypId, studyId, gen)
74 _name = "RadialQuadrangle_1D2D";
75 _shapeType = (1 << TopAbs_FACE); // 1 bit per shape type
77 _compatibleHypothesis.push_back("LayerDistribution2D");
78 _compatibleHypothesis.push_back("NumberOfLayers2D");
80 myDistributionHypo = 0;
81 _requireDescretBoundary = false;
82 _supportSubmeshes = true;
86 //================================================================================
90 //================================================================================
92 StdMeshers_RadialQuadrangle_1D2D::~StdMeshers_RadialQuadrangle_1D2D()
96 //=======================================================================
97 //function : CheckHypothesis
99 //=======================================================================
101 bool StdMeshers_RadialQuadrangle_1D2D::CheckHypothesis
103 const TopoDS_Shape& aShape,
104 SMESH_Hypothesis::Hypothesis_Status& aStatus)
108 myDistributionHypo = 0;
110 list <const SMESHDS_Hypothesis * >::const_iterator itl;
112 const list <const SMESHDS_Hypothesis * >&hyps = GetUsedHypothesis(aMesh, aShape);
113 if ( hyps.size() == 0 ) {
114 aStatus = SMESH_Hypothesis::HYP_OK;
115 return true; // can work with no hypothesis
118 if ( hyps.size() > 1 ) {
119 aStatus = SMESH_Hypothesis::HYP_ALREADY_EXIST;
123 const SMESHDS_Hypothesis *theHyp = hyps.front();
125 string hypName = theHyp->GetName();
127 if (hypName == "NumberOfLayers2D") {
128 myNbLayerHypo = static_cast<const StdMeshers_NumberOfLayers *>(theHyp);
129 aStatus = SMESH_Hypothesis::HYP_OK;
132 if (hypName == "LayerDistribution2D") {
133 myDistributionHypo = static_cast<const StdMeshers_LayerDistribution *>(theHyp);
134 aStatus = SMESH_Hypothesis::HYP_OK;
137 aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
143 // ------------------------------------------------------------------------------
145 * \brief Listener used to mark edges meshed by StdMeshers_RadialQuadrangle_1D2D
147 class TEdgeMarker : public SMESH_subMeshEventListener
149 TEdgeMarker(): SMESH_subMeshEventListener(/*isDeletable=*/false) {}
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 = -1000;
299 map < int, SMESH_1D_Algo * > & algoMap = aMesh.GetGen()->_map1D_Algo;
300 map < int, SMESH_1D_Algo * >::iterator id_algo = algoMap.find( myID );
301 if ( id_algo == algoMap.end() )
302 return new TNodeDistributor( myID, 0, aMesh.GetGen() );
303 return static_cast< TNodeDistributor* >( id_algo->second );
305 // -----------------------------------------------------------------------------
306 //! Computes distribution of nodes on a straight line ending at pIn and pOut
307 bool Compute( vector< double > & positions,
311 const SMESH_Hypothesis* hyp1d)
313 if ( !hyp1d ) return error( "Invalid LayerDistribution hypothesis");
315 double len = pIn.Distance( pOut );
316 if ( len <= DBL_MIN ) return error("Too close points of inner and outer shells");
319 myUsedHyps.push_back( hyp1d );
321 TopoDS_Edge edge = BRepBuilderAPI_MakeEdge( pIn, pOut );
322 SMESH_Hypothesis::Hypothesis_Status aStatus;
323 if ( !StdMeshers_Regular_1D::CheckHypothesis( aMesh, edge, aStatus ))
324 return error( "StdMeshers_Regular_1D::CheckHypothesis() failed "
325 "with LayerDistribution hypothesis");
327 BRepAdaptor_Curve C3D(edge);
328 double f = C3D.FirstParameter(), l = C3D.LastParameter();
329 list< double > params;
330 if ( !StdMeshers_Regular_1D::computeInternalParameters( aMesh, C3D, len, f, l, params, false ))
331 return error("StdMeshers_Regular_1D failed to compute layers distribution");
334 positions.reserve( params.size() );
335 for (list<double>::iterator itU = params.begin(); itU != params.end(); itU++)
336 positions.push_back( *itU / len );
339 // -----------------------------------------------------------------------------
340 //! Make mesh on an adge using assigned 1d hyp or defaut nb of segments
341 bool ComputeCircularEdge(SMESH_Mesh& aMesh,
342 const TopoDS_Edge& anEdge)
344 _gen->Compute( aMesh, anEdge);
345 SMESH_subMesh *sm = aMesh.GetSubMesh(anEdge);
346 if ( sm->GetComputeState() != SMESH_subMesh::COMPUTE_OK)
348 // find any 1d hyp assigned (there can be a hyp w/o algo)
349 myUsedHyps = SMESH_Algo::GetUsedHypothesis(aMesh, anEdge, /*ignoreAux=*/true);
350 Hypothesis_Status aStatus;
351 if ( !StdMeshers_Regular_1D::CheckHypothesis( aMesh, anEdge, aStatus ))
353 // no valid 1d hyp assigned, use default nb of segments
354 _hypType = NB_SEGMENTS;
355 _ivalue[ DISTR_TYPE_IND ] = StdMeshers_NumberOfSegments::DT_Regular;
356 _ivalue[ NB_SEGMENTS_IND ] = _gen->GetDefaultNbSegments();
358 return StdMeshers_Regular_1D::Compute( aMesh, anEdge );
362 // -----------------------------------------------------------------------------
363 //! Make mesh on an adge using assigned 1d hyp or defaut nb of segments
364 bool EvaluateCircularEdge(SMESH_Mesh& aMesh,
365 const TopoDS_Edge& anEdge,
366 MapShapeNbElems& aResMap)
368 _gen->Evaluate( aMesh, anEdge, aResMap );
369 if ( aResMap.count( aMesh.GetSubMesh( anEdge )))
372 // find any 1d hyp assigned
373 myUsedHyps = SMESH_Algo::GetUsedHypothesis(aMesh, anEdge, /*ignoreAux=*/true);
374 Hypothesis_Status aStatus;
375 if ( !StdMeshers_Regular_1D::CheckHypothesis( aMesh, anEdge, aStatus ))
377 // no valid 1d hyp assigned, use default nb of segments
378 _hypType = NB_SEGMENTS;
379 _ivalue[ DISTR_TYPE_IND ] = StdMeshers_NumberOfSegments::DT_Regular;
380 _ivalue[ NB_SEGMENTS_IND ] = _gen->GetDefaultNbSegments();
382 return StdMeshers_Regular_1D::Evaluate( aMesh, anEdge, aResMap );
385 // -----------------------------------------------------------------------------
386 TNodeDistributor( int hypId, int studyId, SMESH_Gen* gen)
387 : StdMeshers_Regular_1D( hypId, studyId, gen)
390 // -----------------------------------------------------------------------------
391 virtual const list <const SMESHDS_Hypothesis *> &
392 GetUsedHypothesis(SMESH_Mesh &, const TopoDS_Shape &, const bool)
396 // -----------------------------------------------------------------------------
400 //=======================================================================
402 * \brief Allow algo to do something after persistent restoration
403 * \param subMesh - restored submesh
405 * call markEdgeAsComputedByMe()
407 //=======================================================================
409 void StdMeshers_RadialQuadrangle_1D2D::SubmeshRestored(SMESH_subMesh* faceSubMesh)
411 if ( !faceSubMesh->IsEmpty() )
413 TopoDS_Edge CircEdge, LinEdge1, LinEdge2;
414 analyseFace( faceSubMesh->GetSubShape(), CircEdge, LinEdge1, LinEdge2 );
415 if ( !CircEdge.IsNull() ) markEdgeAsComputedByMe( CircEdge, faceSubMesh );
416 if ( !LinEdge1.IsNull() ) markEdgeAsComputedByMe( LinEdge1, faceSubMesh );
417 if ( !LinEdge2.IsNull() ) markEdgeAsComputedByMe( LinEdge2, faceSubMesh );
421 //=======================================================================
424 //=======================================================================
426 bool StdMeshers_RadialQuadrangle_1D2D::Compute(SMESH_Mesh& aMesh,
427 const TopoDS_Shape& aShape)
430 SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
432 myHelper = new SMESH_MesherHelper( aMesh );
433 myHelper->IsQuadraticSubMesh( aShape );
434 // to delete helper at exit from Compute()
435 auto_ptr<SMESH_MesherHelper> helperDeleter( myHelper );
437 TNodeDistributor* algo1d = TNodeDistributor::GetDistributor(aMesh);
439 TopoDS_Edge CircEdge, LinEdge1, LinEdge2;
440 int nbe = analyseFace( aShape, CircEdge, LinEdge1, LinEdge2 );
441 Handle(Geom_Circle) aCirc = Handle(Geom_Circle)::DownCast( getCurve( CircEdge ));
442 if( nbe>3 || nbe < 1 || aCirc.IsNull() )
443 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)");
446 // points for rotation
447 TColgp_SequenceOfPnt Points;
448 // angles for rotation
449 TColStd_SequenceOfReal Angles;
450 // Nodes1 and Nodes2 - nodes along radiuses
451 // CNodes - nodes on circle edge
452 vector< const SMDS_MeshNode* > Nodes1, Nodes2, CNodes;
454 // parameters edge nodes on face
455 TColgp_SequenceOfPnt2d Pnts2d1;
458 int faceID = meshDS->ShapeToIndex(aShape);
459 TopoDS_Face F = TopoDS::Face(aShape);
460 Handle(Geom_Surface) S = BRep_Tool::Surface(F);
465 if (!algo1d->ComputeCircularEdge( aMesh, CircEdge ))
466 return error( algo1d->GetComputeError() );
467 map< double, const SMDS_MeshNode* > theNodes;
468 if ( !GetSortedNodesOnEdge(aMesh.GetMeshDS(),CircEdge,true,theNodes))
469 return error("Circular edge is incorrectly meshed");
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 exp.Init( CircEdge, TopAbs_VERTEX );
493 TopoDS_Vertex V1 = TopoDS::Vertex( exp.Current() );
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 bool linEdge1Computed, linEdge2Computed;
701 if ( !computeLayerPositions(P0,P1,LinEdge1,&linEdge1Computed))
704 Nodes1.resize( myLayerPositions.size()+1 );
705 Nodes2.resize( myLayerPositions.size()+1 );
707 // check that both linear edges have same hypotheses
708 if ( !computeLayerPositions(P0,P1,LinEdge2, &linEdge2Computed))
710 if ( Nodes1.size() != myLayerPositions.size()+1 )
711 return error("Different hypotheses apply to radial edges");
713 exp.Init( LinEdge1, TopAbs_VERTEX );
714 TopoDS_Vertex V1 = TopoDS::Vertex( exp.Current() );
716 TopoDS_Vertex V2 = TopoDS::Vertex( exp.Current() );
717 gp_Pnt PE1 = BRep_Tool::Pnt(V1);
718 gp_Pnt PE2 = BRep_Tool::Pnt(V2);
719 if( ( P1.Distance(PE1) > Precision::Confusion() ) &&
720 ( P1.Distance(PE2) > Precision::Confusion() ) )
722 std::swap( LinEdge1, LinEdge2 );
723 std::swap( linEdge1Computed, linEdge2Computed );
725 TopoDS_Vertex VC = V2;
726 if( ( P1.Distance(PE1) > Precision::Confusion() ) &&
727 ( P2.Distance(PE1) > Precision::Confusion() ) )
729 int vertID = meshDS->ShapeToIndex(VC);
732 if ( linEdge1Computed )
734 if (!GetSortedNodesOnEdge(aMesh.GetMeshDS(),LinEdge1,true,theNodes))
735 return error("Invalid mesh on a straight edge");
737 bool nodesFromP0ToP1 = ( theNodes.rbegin()->second == NF );
738 NC = const_cast<SMDS_MeshNode*>
739 ( nodesFromP0ToP1 ? theNodes.begin()->second : theNodes.rbegin()->second );
740 int i = 0, ir = Nodes1.size()-1;
741 int * pi = nodesFromP0ToP1 ? &i : &ir;
742 itn = theNodes.begin();
743 if ( nodesFromP0ToP1 ) ++itn;
744 for ( ; i < Nodes1.size(); ++i, --ir, ++itn )
746 Nodes1[*pi] = itn->second;
748 for ( i = 0; i < Nodes1.size()-1; ++i )
750 Points.Append( gpXYZ( Nodes1[i]));
751 Pnts2d1.Append( myHelper->GetNodeUV( F, Nodes1[i]));
756 int edgeID = meshDS->ShapeToIndex(LinEdge1);
759 Handle(Geom_Curve) Crv = BRep_Tool::Curve(LinEdge1,fp,lp);
760 gp_Pnt Ptmp = Crv->Value(fp);
762 if( P1.Distance(Ptmp) > Precision::Confusion() )
764 // get UV points for edge
766 BRep_Tool::UVPoints( LinEdge1, TopoDS::Face(aShape), PF, PL );
769 V2d = gp_Vec2d(PF,PL);
773 V2d = gp_Vec2d(PL,PF);
776 NC = const_cast<SMDS_MeshNode*>( VertexNode( VC, meshDS ));
779 NC = meshDS->AddNode(P0.X(), P0.Y(), P0.Z());
780 meshDS->SetNodeOnVertex(NC, vertID);
784 for(; i<myLayerPositions.size(); i++) {
785 gp_Pnt P( P0.X() + aVec.X()*myLayerPositions[i],
786 P0.Y() + aVec.Y()*myLayerPositions[i],
787 P0.Z() + aVec.Z()*myLayerPositions[i] );
789 SMDS_MeshNode * node = meshDS->AddNode(P.X(), P.Y(), P.Z());
793 param = fp + dp*(1-myLayerPositions[i]);
795 param = fp + dp*myLayerPositions[i];
796 meshDS->SetNodeOnEdge(node, edgeID, param);
797 // parameters on face
798 gp_Pnt2d P2d( PC.X() + V2d.X()*myLayerPositions[i],
799 PC.Y() + V2d.Y()*myLayerPositions[i] );
802 Nodes1[ myLayerPositions.size() ] = NF;
803 // create 1D elements on edge
804 SMDS_MeshEdge* ME = myHelper->AddEdge( NC, Nodes1[0] );
805 if(ME) meshDS->SetMeshElementOnShape(ME, edgeID);
806 for(i=1; i<Nodes1.size(); i++) {
807 ME = myHelper->AddEdge( Nodes1[i-1], Nodes1[i] );
808 if(ME) meshDS->SetMeshElementOnShape(ME, edgeID);
810 if (nbe==2 && LinEdge1.Orientation() == TopAbs_INTERNAL )
813 markEdgeAsComputedByMe( LinEdge1, aMesh.GetSubMesh( F ));
816 if ( linEdge2Computed )
818 if (!GetSortedNodesOnEdge(aMesh.GetMeshDS(),LinEdge2,true,theNodes))
819 return error("Invalid mesh on a straight edge");
821 bool nodesFromP0ToP2 = ( theNodes.rbegin()->second == NL );
822 int i = 0, ir = Nodes1.size()-1;
823 int * pi = nodesFromP0ToP2 ? &i : &ir;
824 itn = theNodes.begin();
825 if ( nodesFromP0ToP2 ) ++itn;
826 for ( ; i < Nodes2.size(); ++i, --ir, ++itn )
827 Nodes2[*pi] = itn->second;
831 int edgeID = meshDS->ShapeToIndex(LinEdge2);
832 gp_Vec aVec = gp_Vec(P0,P2);
834 Handle(Geom_Curve) Crv = BRep_Tool::Curve(LinEdge2,fp,lp);
835 gp_Pnt Ptmp = Crv->Value(fp);
837 if( P2.Distance(Ptmp) > Precision::Confusion() )
839 // get UV points for edge
841 BRep_Tool::UVPoints( LinEdge2, TopoDS::Face(aShape), PF, PL );
844 V2d = gp_Vec2d(PF,PL);
848 V2d = gp_Vec2d(PL,PF);
852 for(int i=0; i<myLayerPositions.size(); i++) {
853 gp_Pnt P( P0.X() + aVec.X()*myLayerPositions[i],
854 P0.Y() + aVec.Y()*myLayerPositions[i],
855 P0.Z() + aVec.Z()*myLayerPositions[i] );
856 SMDS_MeshNode * node = meshDS->AddNode(P.X(), P.Y(), P.Z());
860 param = fp + dp*(1-myLayerPositions[i]);
862 param = fp + dp*myLayerPositions[i];
863 meshDS->SetNodeOnEdge(node, edgeID, param);
864 // parameters on face
865 gp_Pnt2d P2d( PC.X() + V2d.X()*myLayerPositions[i],
866 PC.Y() + V2d.Y()*myLayerPositions[i] );
868 Nodes2[ myLayerPositions.size() ] = NL;
869 // create 1D elements on edge
870 SMDS_MeshEdge* ME = myHelper->AddEdge( NC, Nodes2[0] );
871 if(ME) meshDS->SetMeshElementOnShape(ME, edgeID);
872 for(int i=1; i<Nodes2.size(); i++) {
873 ME = myHelper->AddEdge( Nodes2[i-1], Nodes2[i] );
874 if(ME) meshDS->SetMeshElementOnShape(ME, edgeID);
877 markEdgeAsComputedByMe( LinEdge2, aMesh.GetSubMesh( F ));
879 markEdgeAsComputedByMe( CircEdge, aMesh.GetSubMesh( F ));
882 bool IsForward = ( CircEdge.Orientation()==TopAbs_FORWARD );
884 // create nodes and mesh elements on face
885 // find axis of rotation
886 gp_Pnt P2 = gp_Pnt( CNodes[1]->X(), CNodes[1]->Y(), CNodes[1]->Z() );
889 gp_Vec Axis = Vec1.Crossed(Vec2);
892 //cout<<"Angles.Length() = "<<Angles.Length()<<" Points.Length() = "<<Points.Length()<<endl;
893 //cout<<"Nodes1.size() = "<<Nodes1.size()<<" Pnts2d1.Length() = "<<Pnts2d1.Length()<<endl;
894 for(; i<Angles.Length(); i++) {
895 vector< const SMDS_MeshNode* > tmpNodes;
896 tmpNodes.reserve(Nodes1.size());
898 gp_Ax1 theAxis(P0,gp_Dir(Axis));
899 aTrsf.SetRotation( theAxis, Angles.Value(i) );
901 aTrsf2d.SetRotation( PC, Angles.Value(i) );
904 for(; j<=Points.Length(); j++) {
906 Points.Value(j).Coord( cx, cy, cz );
907 aTrsf.Transforms( cx, cy, cz );
908 SMDS_MeshNode* node = myHelper->AddNode( cx, cy, cz );
909 // find parameters on face
910 Pnts2d1.Value(j).Coord( cx, cy );
911 aTrsf2d.Transforms( cx, cy );
913 meshDS->SetNodeOnFace( node, faceID, cx, cy );
914 tmpNodes[j-1] = node;
917 tmpNodes[Points.Length()] = CNodes[i];
919 for(j=0; j<Nodes1.size()-1; j++) {
922 MF = myHelper->AddFace( tmpNodes[j], Nodes1[j],
923 Nodes1[j+1], tmpNodes[j+1] );
925 MF = myHelper->AddFace( tmpNodes[j], tmpNodes[j+1],
926 Nodes1[j+1], Nodes1[j] );
927 if(MF) meshDS->SetMeshElementOnShape(MF, faceID);
932 MF = myHelper->AddFace( NC, Nodes1[0], tmpNodes[0] );
934 MF = myHelper->AddFace( NC, tmpNodes[0], Nodes1[0] );
935 if(MF) meshDS->SetMeshElementOnShape(MF, faceID);
936 for(j=0; j<Nodes1.size(); j++) {
937 Nodes1[j] = tmpNodes[j];
942 for(i=0; i<Nodes1.size()-1; i++) {
945 MF = myHelper->AddFace( Nodes2[i], Nodes1[i],
946 Nodes1[i+1], Nodes2[i+1] );
948 MF = myHelper->AddFace( Nodes2[i], Nodes2[i+1],
949 Nodes1[i+1], Nodes1[i] );
950 if(MF) meshDS->SetMeshElementOnShape(MF, faceID);
955 MF = myHelper->AddFace( NC, Nodes1[0], Nodes2[0] );
957 MF = myHelper->AddFace( NC, Nodes2[0], Nodes1[0] );
958 if(MF) meshDS->SetMeshElementOnShape(MF, faceID);
963 //================================================================================
965 * \brief Compute positions of nodes on the radial edge
966 * \retval bool - is a success
968 //================================================================================
970 bool StdMeshers_RadialQuadrangle_1D2D::computeLayerPositions(const gp_Pnt& p1,
972 const TopoDS_Edge& linEdge,
973 bool* linEdgeComputed)
975 // First, try to compute positions of layers
977 myLayerPositions.clear();
979 SMESH_Mesh * mesh = myHelper->GetMesh();
981 const SMESH_Hypothesis* hyp1D = myDistributionHypo ? myDistributionHypo->GetLayerDistribution() : 0;
982 int nbLayers = myNbLayerHypo ? myNbLayerHypo->GetNumberOfLayers() : 0;
984 if ( !hyp1D && !nbLayers )
986 // No own algo hypotheses assigned, so first try to find any 1D hypothesis.
988 TopoDS_Shape edge = linEdge;
989 if ( edge.IsNull() && !myHelper->GetSubShape().IsNull())
990 for ( TopExp_Explorer e(myHelper->GetSubShape(), TopAbs_EDGE); e.More(); e.Next())
992 if ( !edge.IsNull() )
994 // find a hyp usable by TNodeDistributor
995 SMESH_HypoFilter hypKind;
996 TNodeDistributor::GetDistributor(*mesh)->InitCompatibleHypoFilter(hypKind,/*ignoreAux=*/1);
997 hyp1D = mesh->GetHypothesis( edge, hypKind, /*fromAncestors=*/true);
1000 if ( hyp1D ) // try to compute with hyp1D
1002 if ( !TNodeDistributor::GetDistributor(*mesh)->Compute( myLayerPositions,p1,p2,*mesh,hyp1D )) {
1003 if ( myDistributionHypo ) { // bad hyp assigned
1004 return error( TNodeDistributor::GetDistributor(*mesh)->GetComputeError() );
1007 // bad hyp found, its Ok, lets try with default nb of segnents
1012 if ( myLayerPositions.empty() ) // try to use nb of layers
1015 nbLayers = _gen->GetDefaultNbSegments();
1019 myLayerPositions.resize( nbLayers - 1 );
1020 for ( int z = 1; z < nbLayers; ++z )
1021 myLayerPositions[ z - 1 ] = double( z )/ double( nbLayers );
1025 // Second, check presence of a mesh built by other algo on linEdge
1026 // and mesh conformity to my hypothesis
1028 bool meshComputed = (!linEdge.IsNull() && !mesh->GetSubMesh(linEdge)->IsEmpty() );
1029 if ( linEdgeComputed ) *linEdgeComputed = meshComputed;
1033 vector< double > nodeParams;
1034 GetNodeParamOnEdge( mesh->GetMeshDS(), linEdge, nodeParams );
1036 // nb of present nodes must be different in cases of 1 and 2 straight edges
1038 TopoDS_Vertex VV[2];
1039 TopExp::Vertices( linEdge, VV[0], VV[1]);
1040 const gp_Pnt* points[] = { &p1, &p2 };
1041 gp_Pnt vPoints[] = { BRep_Tool::Pnt(VV[0]), BRep_Tool::Pnt(VV[1]) };
1042 const double tol[] = { BRep_Tool::Tolerance(VV[0]), BRep_Tool::Tolerance(VV[1]) };
1043 bool pointsAreOnVertices = true;
1044 for ( int iP = 0; iP < 2 && pointsAreOnVertices; ++iP )
1045 pointsAreOnVertices = ( points[iP]->Distance( vPoints[0] ) < tol[0] ||
1046 points[iP]->Distance( vPoints[1] ) < tol[1] );
1048 int nbNodes = nodeParams.size() - 2; // 2 straight edges
1049 if ( !pointsAreOnVertices )
1050 nbNodes = ( nodeParams.size() - 3 ) / 2; // 1 straight edge
1052 if ( myLayerPositions.empty() )
1054 myLayerPositions.resize( nbNodes );
1056 else if ( myDistributionHypo || myNbLayerHypo )
1058 // linEdge is computed by other algo. Check if there is a meshed face
1059 // using nodes on linEdge
1060 bool nodesAreUsed = false;
1061 TopTools_ListIteratorOfListOfShape ancestIt = mesh->GetAncestors( linEdge );
1062 for ( ; ancestIt.More() && !nodesAreUsed; ancestIt.Next() )
1063 if ( ancestIt.Value().ShapeType() == TopAbs_FACE )
1064 nodesAreUsed = (!mesh->GetSubMesh( ancestIt.Value() )->IsEmpty());
1065 if ( !nodesAreUsed ) {
1067 mesh->GetSubMesh( linEdge )->ComputeStateEngine( SMESH_subMesh::CLEAN );
1068 if ( linEdgeComputed ) *linEdgeComputed = false;
1072 if ( myLayerPositions.size() != nbNodes )
1073 return error("Radial edge is meshed by other algorithm");
1078 return !myLayerPositions.empty();
1082 //=======================================================================
1083 //function : Evaluate
1085 //=======================================================================
1087 bool StdMeshers_RadialQuadrangle_1D2D::Evaluate(SMESH_Mesh& aMesh,
1088 const TopoDS_Shape& aShape,
1089 MapShapeNbElems& aResMap)
1091 if( aShape.ShapeType() != TopAbs_FACE ) {
1094 SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
1095 if( aResMap.count(sm) )
1098 vector<int>& aResVec =
1099 aResMap.insert( make_pair(sm, vector<int>(SMDSEntity_Last,0))).first->second;
1101 myHelper = new SMESH_MesherHelper( aMesh );
1102 myHelper->SetSubShape( aShape );
1103 auto_ptr<SMESH_MesherHelper> helperDeleter( myHelper );
1105 TNodeDistributor* algo1d = TNodeDistributor::GetDistributor(aMesh);
1107 TopoDS_Edge CircEdge, LinEdge1, LinEdge2;
1108 int nbe = analyseFace( aShape, CircEdge, LinEdge1, LinEdge2 );
1109 if( nbe>3 || nbe < 1 || CircEdge.IsNull() )
1112 Handle(Geom_Circle) aCirc = Handle(Geom_Circle)::DownCast( getCurve( CircEdge ));
1113 if( aCirc.IsNull() )
1114 return error(COMPERR_BAD_SHAPE);
1116 gp_Pnt P0 = aCirc->Location();
1117 gp_Pnt P1 = aCirc->Value(0.);
1118 computeLayerPositions( P0, P1, LinEdge1 );
1120 int nb0d=0, nb2d_tria=0, nb2d_quad=0;
1121 bool isQuadratic = false, ok = true;
1124 // C1 must be a circle
1125 ok = algo1d->EvaluateCircularEdge( aMesh, CircEdge, aResMap );
1127 const vector<int>& aVec = aResMap[aMesh.GetSubMesh(CircEdge)];
1128 isQuadratic = aVec[SMDSEntity_Quad_Edge]>aVec[SMDSEntity_Edge];
1131 nb0d = (aVec[SMDSEntity_Node]+1) * myLayerPositions.size();
1132 // radial medium nodes
1133 nb0d += (aVec[SMDSEntity_Node]+1) * (myLayerPositions.size()+1);
1134 // other medium nodes
1135 nb0d += (aVec[SMDSEntity_Node]+1) * myLayerPositions.size();
1138 nb0d = (aVec[SMDSEntity_Node]+1) * myLayerPositions.size();
1140 nb2d_tria = aVec[SMDSEntity_Node] + 1;
1144 else if(nbe==2 && LinEdge1.Orientation() != TopAbs_INTERNAL)
1146 // one curve must be a half of circle and other curve must be
1147 // a segment of line
1149 Handle(Geom_Circle) aCirc = Handle(Geom_Circle)::DownCast( getCurve( CircEdge, &fp, &lp ));
1150 if( fabs(fabs(lp-fp)-M_PI) > Precision::Confusion() ) {
1151 // not half of circle
1152 return error(COMPERR_BAD_SHAPE);
1154 Handle(Geom_Line) aLine = Handle(Geom_Line)::DownCast( getCurve( LinEdge1 ));
1155 if( aLine.IsNull() ) {
1156 // other curve not line
1157 return error(COMPERR_BAD_SHAPE);
1159 ok = !aResMap.count( aMesh.GetSubMesh(LinEdge1) );
1161 const vector<int>& aVec = aResMap[ aMesh.GetSubMesh(LinEdge1) ];
1162 ok = ( aVec[SMDSEntity_Node] == myLayerPositions.size() );
1165 ok = algo1d->EvaluateCircularEdge( aMesh, CircEdge, aResMap );
1168 const vector<int>& aVec = aResMap[ aMesh.GetSubMesh(CircEdge) ];
1169 isQuadratic = aVec[SMDSEntity_Quad_Edge] > aVec[SMDSEntity_Edge];
1172 nb0d = aVec[SMDSEntity_Node] * myLayerPositions.size();
1173 // radial medium nodes
1174 nb0d += aVec[SMDSEntity_Node] * (myLayerPositions.size()+1);
1175 // other medium nodes
1176 nb0d += (aVec[SMDSEntity_Node]+1) * myLayerPositions.size();
1179 nb0d = aVec[SMDSEntity_Node] * myLayerPositions.size();
1181 nb2d_tria = aVec[SMDSEntity_Node] + 1;
1182 nb2d_quad = nb2d_tria * myLayerPositions.size();
1183 // add evaluation for edges
1184 vector<int> aResVec(SMDSEntity_Last,0);
1186 aResVec[SMDSEntity_Node] = 4*myLayerPositions.size() + 3;
1187 aResVec[SMDSEntity_Quad_Edge] = 2*myLayerPositions.size() + 2;
1190 aResVec[SMDSEntity_Node] = 2*myLayerPositions.size() + 1;
1191 aResVec[SMDSEntity_Edge] = 2*myLayerPositions.size() + 2;
1193 aResMap[ aMesh.GetSubMesh(LinEdge1) ] = aResVec;
1196 else // nbe==3 or ( nbe==2 && linEdge is INTERNAL )
1198 if (nbe==2 && LinEdge1.Orientation() == TopAbs_INTERNAL )
1199 LinEdge2 = LinEdge1;
1201 // one curve must be a part of circle and other curves must be
1203 Handle(Geom_Line) aLine1 = Handle(Geom_Line)::DownCast( getCurve( LinEdge1 ));
1204 Handle(Geom_Line) aLine2 = Handle(Geom_Line)::DownCast( getCurve( LinEdge2 ));
1205 if( aLine1.IsNull() || aLine2.IsNull() ) {
1206 // other curve not line
1207 return error(COMPERR_BAD_SHAPE);
1209 int nbLayers = myLayerPositions.size();
1210 computeLayerPositions( P0, P1, LinEdge2 );
1211 if ( nbLayers != myLayerPositions.size() )
1212 return error("Different hypotheses apply to radial edges");
1214 bool ok = !aResMap.count( aMesh.GetSubMesh(LinEdge1));
1216 if ( myDistributionHypo || myNbLayerHypo )
1217 ok = true; // override other 1d hyps
1219 const vector<int>& aVec = aResMap[ aMesh.GetSubMesh(LinEdge1) ];
1220 ok = ( aVec[SMDSEntity_Node] == myLayerPositions.size() );
1223 if( ok && aResMap.count( aMesh.GetSubMesh(LinEdge2) )) {
1224 if ( myDistributionHypo || myNbLayerHypo )
1225 ok = true; // override other 1d hyps
1227 const vector<int>& aVec = aResMap[ aMesh.GetSubMesh(LinEdge2) ];
1228 ok = ( aVec[SMDSEntity_Node] == myLayerPositions.size() );
1232 ok = algo1d->EvaluateCircularEdge( aMesh, CircEdge, aResMap );
1235 const vector<int>& aVec = aResMap[ aMesh.GetSubMesh(CircEdge) ];
1236 isQuadratic = aVec[SMDSEntity_Quad_Edge]>aVec[SMDSEntity_Edge];
1239 nb0d = aVec[SMDSEntity_Node] * myLayerPositions.size();
1240 // radial medium nodes
1241 nb0d += aVec[SMDSEntity_Node] * (myLayerPositions.size()+1);
1242 // other medium nodes
1243 nb0d += (aVec[SMDSEntity_Node]+1) * myLayerPositions.size();
1246 nb0d = aVec[SMDSEntity_Node] * myLayerPositions.size();
1248 nb2d_tria = aVec[SMDSEntity_Node] + 1;
1249 nb2d_quad = nb2d_tria * myLayerPositions.size();
1250 // add evaluation for edges
1251 vector<int> aResVec(SMDSEntity_Last, 0);
1253 aResVec[SMDSEntity_Node] = 2*myLayerPositions.size() + 1;
1254 aResVec[SMDSEntity_Quad_Edge] = myLayerPositions.size() + 1;
1257 aResVec[SMDSEntity_Node] = myLayerPositions.size();
1258 aResVec[SMDSEntity_Edge] = myLayerPositions.size() + 1;
1260 sm = aMesh.GetSubMesh(LinEdge1);
1261 aResMap[sm] = aResVec;
1262 sm = aMesh.GetSubMesh(LinEdge2);
1263 aResMap[sm] = aResVec;
1270 aResVec[SMDSEntity_Quad_Triangle] = nb2d_tria;
1271 aResVec[SMDSEntity_Quad_Quadrangle] = nb2d_quad;
1274 aResVec[SMDSEntity_Triangle] = nb2d_tria;
1275 aResVec[SMDSEntity_Quadrangle] = nb2d_quad;
1281 sm = aMesh.GetSubMesh(aShape);
1282 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1283 smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,
1284 "Submesh can not be evaluated",this));