Salome HOME
Merge branch 'V9_9_BR'
[modules/smesh.git] / src / StdMeshers / StdMeshers_MEFISTO_2D.cxx
1 // Copyright (C) 2007-2022  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 // Copyright (C) 2003-2007  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
5 //
6 // This library is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License, or (at your option) any later version.
10 //
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
14 // Lesser General Public License for more details.
15 //
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
19 //
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 //
22
23 //  SMESH SMESH : implementation of SMESH idl descriptions
24 //  File   : StdMeshers_MEFISTO_2D.cxx
25 //           Moved here from SMESH_MEFISTO_2D.cxx
26 //  Author : Paul RASCLE, EDF
27 //  Module : SMESH
28 //
29 #include "StdMeshers_MEFISTO_2D.hxx"
30
31 #include "SMDS_EdgePosition.hxx"
32 #include "SMDS_MeshElement.hxx"
33 #include "SMDS_MeshNode.hxx"
34 #include "SMESHDS_Mesh.hxx"
35 #include "SMESH_Comment.hxx"
36 #include "SMESH_Gen.hxx"
37 #include "SMESH_Mesh.hxx"
38 #include "SMESH_MesherHelper.hxx"
39 #include "SMESH_subMesh.hxx"
40 #include "StdMeshers_FaceSide.hxx"
41 #include "StdMeshers_LengthFromEdges.hxx"
42 #include "StdMeshers_MaxElementArea.hxx"
43 #include "StdMeshers_ViscousLayers2D.hxx"
44
45 #include "utilities.h"
46
47 #include "Rn.h"
48 #include "aptrte.h"
49
50 #include <BRepGProp.hxx>
51 #include <BRepTools.hxx>
52 #include <BRep_Tool.hxx>
53 #include <GProp_GProps.hxx>
54 #include <Geom2d_Curve.hxx>
55 #include <Geom_Curve.hxx>
56 #include <Geom_Surface.hxx>
57 #include <Precision.hxx>
58 #include <TopExp.hxx>
59 #include <TopExp_Explorer.hxx>
60 #include <TopTools_ListIteratorOfListOfShape.hxx>
61 #include <TopTools_ListOfShape.hxx>
62 #include <TopTools_MapOfShape.hxx>
63 #include <TopoDS.hxx>
64 #include <TopoDS_Edge.hxx>
65 #include <TopoDS_Face.hxx>
66 #include <TopoDS_Iterator.hxx>
67 #include <TopoDS_Wire.hxx>
68 #include <gp_Pnt2d.hxx>
69
70 using namespace std;
71
72 #ifdef _DEBUG_
73 //#define DUMP_POINTS // to print coordinates of MEFISTO input
74 #endif
75
76 //=============================================================================
77 /*!
78  *  
79  */
80 //=============================================================================
81
82 StdMeshers_MEFISTO_2D::StdMeshers_MEFISTO_2D(int hypId, SMESH_Gen * gen):
83   SMESH_2D_Algo(hypId, gen)
84 {
85   _name = "MEFISTO_2D";
86   _shapeType = (1 << TopAbs_FACE);
87   _compatibleHypothesis.push_back("MaxElementArea");
88   _compatibleHypothesis.push_back("LengthFromEdges");
89   _compatibleHypothesis.push_back("ViscousLayers2D");
90
91   _edgeLength = 0;
92   _maxElementArea = 0;
93   _hypMaxElementArea = NULL;
94   _hypLengthFromEdges = NULL;
95   _helper = 0;
96 }
97
98 //=============================================================================
99 /*!
100  *  
101  */
102 //=============================================================================
103
104 StdMeshers_MEFISTO_2D::~StdMeshers_MEFISTO_2D()
105 {
106 }
107
108 //=============================================================================
109 /*!
110  *  
111  */
112 //=============================================================================
113
114 bool StdMeshers_MEFISTO_2D::CheckHypothesis
115                          (SMESH_Mesh&                          aMesh,
116                           const TopoDS_Shape&                  aShape,
117                           SMESH_Hypothesis::Hypothesis_Status& aStatus)
118 {
119   _hypMaxElementArea = NULL;
120   _hypLengthFromEdges = NULL;
121   _edgeLength = 0;
122   _maxElementArea = 0;
123
124   if ( !error( StdMeshers_ViscousLayers2D::CheckHypothesis( aMesh, aShape, aStatus )))
125     return false;
126
127   list <const SMESHDS_Hypothesis * >::const_iterator itl;
128   const SMESHDS_Hypothesis *theHyp;
129
130   const list <const SMESHDS_Hypothesis * >&hyps = GetUsedHypothesis(aMesh, aShape);
131   int nbHyp = hyps.size();
132   if (!nbHyp)
133   {
134     aStatus = SMESH_Hypothesis::HYP_OK; //SMESH_Hypothesis::HYP_MISSING;
135     return true;  // (PAL13464) can work with no hypothesis, LengthFromEdges is default one
136   }
137
138   itl = hyps.begin();
139   theHyp = (*itl); // use only the first hypothesis
140
141   string hypName = theHyp->GetName();
142
143   bool isOk = false;
144
145   if (hypName == "MaxElementArea")
146   {
147     _hypMaxElementArea = static_cast<const StdMeshers_MaxElementArea *>(theHyp);
148     ASSERT(_hypMaxElementArea);
149     _maxElementArea = _hypMaxElementArea->GetMaxArea();
150     isOk = true;
151     aStatus = SMESH_Hypothesis::HYP_OK;
152   }
153
154   else if (hypName == "LengthFromEdges")
155   {
156     _hypLengthFromEdges = static_cast<const StdMeshers_LengthFromEdges *>(theHyp);
157     ASSERT(_hypLengthFromEdges);
158     isOk = true;
159     aStatus = SMESH_Hypothesis::HYP_OK;
160   }
161   else
162     aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
163
164   if (isOk)
165   {
166     isOk = false;
167     if (_maxElementArea > 0)
168     {
169       //_edgeLength = 2 * sqrt(_maxElementArea);        // triangles : minorant
170       _edgeLength = sqrt(2. * _maxElementArea/sqrt(3.0));
171       isOk = true;
172     }
173     else
174       isOk = (_hypLengthFromEdges != NULL);     // **** check mode
175     if (!isOk)
176       aStatus = SMESH_Hypothesis::HYP_BAD_PARAMETER;
177   }
178
179   return isOk;
180 }
181
182 //=============================================================================
183 /*!
184  *  
185  */
186 //=============================================================================
187
188 bool StdMeshers_MEFISTO_2D::Compute(SMESH_Mesh & aMesh, const TopoDS_Shape & aShape)
189 {
190   TopoDS_Face F = TopoDS::Face(aShape.Oriented(TopAbs_FORWARD));
191
192   // helper builds quadratic mesh if necessary
193   SMESH_MesherHelper helper(aMesh);
194   _helper = &helper;
195   _quadraticMesh = _helper->IsQuadraticSubMesh(aShape);
196   const bool skipMediumNodes = _quadraticMesh;
197
198   // build viscous layers if required
199   SMESH_ProxyMesh::Ptr proxyMesh = StdMeshers_ViscousLayers2D::Compute( aMesh, F );
200   if ( !proxyMesh )
201     return false;
202
203   // get all edges of a face
204   TError problem;
205   TWireVector wires =
206     StdMeshers_FaceSide::GetFaceWires( F, aMesh, skipMediumNodes, problem, _helper, proxyMesh );
207   int nbWires = wires.size();
208   if ( problem && !problem->IsOK() ) return error( problem );
209   if ( nbWires == 0 ) return error( "Problem in StdMeshers_FaceSide::GetFaceWires()");
210   if ( wires[0]->NbSegments() < 3 ) // ex: a circle with 2 segments
211     return error(COMPERR_BAD_INPUT_MESH,
212                  SMESH_Comment("Too few segments: ")<<wires[0]->NbSegments());
213
214   // compute average edge length
215   if (!_hypMaxElementArea)
216   {
217     _edgeLength = 0;
218     int nbSegments = 0;
219     for ( int iW = 0; iW < nbWires; ++iW )
220     {
221       StdMeshers_FaceSidePtr wire = wires[ iW ];
222       _edgeLength += wire->Length();
223       nbSegments  += wire->NbSegments();
224     }
225     if ( nbSegments )
226       _edgeLength /= nbSegments;
227   }
228
229   if (/*_hypLengthFromEdges &&*/ _edgeLength < DBL_MIN )
230     _edgeLength = 100;
231
232   Z nblf;                 //nombre de lignes fermees (enveloppe en tete)
233   Z *nudslf = NULL;       //numero du dernier sommet de chaque ligne fermee
234   R2 *uvslf = NULL;       
235   Z nbpti = 0;            //nombre points internes futurs sommets de la triangulation
236   R2 *uvpti = NULL;
237   
238   Z nbst;
239   R2 *uvst = NULL;
240   Z nbt;
241   Z *nust = NULL;
242   Z ierr = 0;
243
244   Z nutysu = 1;           // 1: il existe un fonction areteideale_()
245   // Z  nutysu=0;         // 0: on utilise aretmx
246   R aretmx = _edgeLength; // longueur max aretes future triangulation
247   if ( _hypMaxElementArea )
248     aretmx *= 1.5;
249   
250   nblf = nbWires;
251   
252   nudslf = new Z[1 + nblf];
253   nudslf[0] = 0;
254   int iw = 1;
255   int nbpnt = 0;
256
257   // count nb of input points
258   for ( int iW = 0; iW < nbWires; ++iW )
259   {
260     nbpnt += wires[iW]->NbPoints() - 1;
261     nudslf[iw++] = nbpnt;
262   }
263
264   uvslf = new R2[nudslf[nblf]];
265
266   double scalex, scaley;
267   ComputeScaleOnFace(aMesh, F, scalex, scaley);
268
269   // correspondence mefisto index --> Nodes
270   vector< const SMDS_MeshNode*> mefistoToDS(nbpnt, (const SMDS_MeshNode*)0);
271
272   bool isOk = false;
273
274   // fill input points UV
275   if ( LoadPoints(wires, uvslf, mefistoToDS, scalex, scaley) )
276   {
277     // Compute
278     aptrte(nutysu, aretmx,
279            nblf, nudslf, uvslf, nbpti, uvpti, nbst, uvst, nbt, nust, ierr);
280
281     if (ierr == 0)
282     {
283       StoreResult(nbst, uvst, nbt, nust, mefistoToDS, scalex, scaley);
284       isOk = true;
285     }
286     else
287     {
288       error(ierr,"Error in Triangulation (aptrte())");
289     }
290   }
291   if (nudslf != NULL) delete[]nudslf;
292   if (uvslf != NULL)  delete[]uvslf;
293   if (uvst != NULL)   delete[]uvst;
294   if (nust != NULL)   delete[]nust;
295
296   return isOk;
297 }
298
299
300 //=============================================================================
301 /*!
302  *  
303  */
304 //=============================================================================
305
306 bool StdMeshers_MEFISTO_2D::Evaluate(SMESH_Mesh & aMesh,
307                                      const TopoDS_Shape & aShape,
308                                      MapShapeNbElems& aResMap)
309 {
310   TopoDS_Face F = TopoDS::Face(aShape.Oriented(TopAbs_FORWARD));
311
312   double aLen = 0.0;
313   int NbSeg = 0;
314   bool IsQuadratic = false;
315   bool IsFirst = true;
316   TopExp_Explorer exp(F,TopAbs_EDGE);
317   for(; exp.More(); exp.Next()) {
318     TopoDS_Edge E = TopoDS::Edge(exp.Current());
319     MapShapeNbElemsItr anIt = aResMap.find( aMesh.GetSubMesh(E) );
320     if( anIt == aResMap.end() ) continue;
321     std::vector<smIdType> aVec = (*anIt).second;
322     int nbe = Max(aVec[SMDSEntity_Edge],aVec[SMDSEntity_Quad_Edge]);
323     NbSeg += nbe;
324     if(IsFirst) {
325       IsQuadratic = ( aVec[SMDSEntity_Quad_Edge] > aVec[SMDSEntity_Edge] );
326       IsFirst = false;
327     }
328     double a,b;
329     TopLoc_Location L;
330     Handle(Geom_Curve) C = BRep_Tool::Curve(E,L,a,b);
331     gp_Pnt P1;
332     C->D0(a,P1);
333     double dp = (b-a)/nbe;
334     for(int i=1; i<=nbe; i++) {
335       gp_Pnt P2;
336       C->D0(a+i*dp,P2);
337       aLen += P1.Distance(P2);
338       P1 = P2;
339     }
340   }
341   if(NbSeg<1) {
342     std::vector<smIdType> aResVec(SMDSEntity_Last);
343     for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
344     SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
345     aResMap.insert(std::make_pair(sm,aResVec));
346     SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
347     smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,
348                                           "Submesh can not be evaluated",this));
349     return false;
350   }
351   aLen = aLen/NbSeg; // middle length
352
353   _edgeLength = Precision::Infinite();
354   double tmpLength = Min( _edgeLength, aLen );
355
356   GProp_GProps G;
357   BRepGProp::SurfaceProperties(aShape,G);
358   double anArea = G.Mass();
359
360   int nbFaces = Precision::IsInfinite( tmpLength ) ? 0 :
361     (int)( anArea/(tmpLength*tmpLength*sqrt(3.)/4) );
362   int nbNodes = (int) ( nbFaces*3 - (NbSeg-1)*2 ) / 6;
363
364   std::vector<smIdType> aVec(SMDSEntity_Last);
365   for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i] = 0;
366   if(IsQuadratic) {
367     aVec[SMDSEntity_Quad_Triangle] = nbFaces;
368     aVec[SMDSEntity_Node] = (int)( nbNodes + nbFaces*3 - (NbSeg-1) );
369   }
370   else {
371     aVec[SMDSEntity_Node] = nbNodes;
372     aVec[SMDSEntity_Triangle] = nbFaces;
373   }
374   SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
375   aResMap.insert(std::make_pair(sm,aVec));
376
377   return true;
378 }
379
380
381 //=======================================================================
382 //function : fixOverlappedLinkUV
383 //purpose  : prevent failure due to overlapped adjacent links
384 //=======================================================================
385
386 static bool fixOverlappedLinkUV( R2& uv0, const R2& uv1, const R2& uv2 )
387 {
388   gp_XY v1( uv0.x - uv1.x, uv0.y - uv1.y );
389   gp_XY v2( uv2.x - uv1.x, uv2.y - uv1.y );
390
391   double tol2 = DBL_MIN * DBL_MIN;
392   double sqMod1 = v1.SquareModulus();
393   if ( sqMod1 <= tol2 ) return false;
394   double sqMod2 = v2.SquareModulus();
395   if ( sqMod2 <= tol2 ) return false;
396
397   double dot = v1*v2;
398
399   // check sinus >= 1.e-3
400   const double minSin = 1.e-3;
401   if ( dot > 0 && 1 - dot * dot / ( sqMod1 * sqMod2 ) < minSin * minSin ) {
402     MESSAGE(" ___ FIX UV ____" << uv0.x << " " << uv0.y);
403     v1.SetCoord( -v1.Y(), v1.X() );
404     double delta = sqrt( sqMod1 ) * minSin;
405     if ( v1.X() < 0 )
406       uv0.x -= delta;
407     else
408       uv0.x += delta;
409     if ( v1.Y() < 0 )
410       uv0.y -= delta;
411     else
412       uv0.y += delta;
413 // #ifdef _DEBUG_
414 //     MESSAGE(" -> " << uv0.x << " " << uv0.y << " ");
415 //     MESSAGE("v1( " << v1.X() << " " << v1.Y() << " ) " <<
416 //       "v2( " << v2.X() << " " << v2.Y() << " ) ");
417 //    MESSAGE("SIN: " << sqrt(1 - dot * dot / (sqMod1 * sqMod2)));
418 //     v1.SetCoord( uv0.x - uv1.x, uv0.y - uv1.y );
419 //     v2.SetCoord( uv2.x - uv1.x, uv2.y - uv1.y );
420 //     gp_XY v3( uv2.x - uv0.x, uv2.y - uv0.y );
421 //     sqMod1 = v1.SquareModulus();
422 //     sqMod2 = v2.SquareModulus();
423 //     dot = v1*v2;
424 //     double sin = sqrt(1 - dot * dot / (sqMod1 * sqMod2));
425 //     MESSAGE("NEW SIN: " << sin);
426 // #endif
427     return true;
428   }
429   return false;
430 }
431
432 //=======================================================================
433 //function : fixCommonVertexUV
434 //purpose  : 
435 //=======================================================================
436
437 static bool fixCommonVertexUV (R2 &                 theUV,
438                                const TopoDS_Vertex& theV,
439                                const TopoDS_Face&   theF,
440                                const TopTools_IndexedDataMapOfShapeListOfShape & theVWMap,
441                                SMESH_Mesh &         theMesh,
442                                const double         theScaleX,
443                                const double         theScaleY,
444                                const bool           theCreateQuadratic)
445 {
446   if( !theVWMap.Contains( theV )) return false;
447
448   // check if there is another wire sharing theV
449   const TopTools_ListOfShape& WList = theVWMap.FindFromKey( theV );
450   TopTools_ListIteratorOfListOfShape aWIt;
451   TopTools_MapOfShape aWires;
452   for ( aWIt.Initialize( WList ); aWIt.More(); aWIt.Next() )
453     aWires.Add( aWIt.Value() );
454   if ( aWires.Extent() < 2 ) return false;
455
456   TopoDS_Shape anOuterWire = BRepTools::OuterWire(theF);
457   TopoDS_Shape anInnerWire;
458   for ( aWIt.Initialize( WList ); aWIt.More() && anInnerWire.IsNull(); aWIt.Next() )
459     if ( !anOuterWire.IsSame( aWIt.Value() ))
460       anInnerWire = aWIt.Value();
461
462   TopTools_ListOfShape EList;
463   list< double >       UList;
464
465   // find edges of theW sharing theV
466   // and find 2d normal to them at theV
467   gp_Vec2d N(0.,0.);
468   TopoDS_Iterator itE( anInnerWire );
469   for (  ; itE.More(); itE.Next() )
470   {
471     const TopoDS_Edge& E = TopoDS::Edge( itE.Value() );
472     TopoDS_Iterator itV( E );
473     for ( ; itV.More(); itV.Next() )
474     {
475       const TopoDS_Vertex & V = TopoDS::Vertex( itV.Value() );
476       if ( !V.IsSame( theV ))
477         continue;
478       EList.Append( E );
479       Standard_Real u = BRep_Tool::Parameter( V, E );
480       UList.push_back( u );
481       double f, l;
482       Handle(Geom2d_Curve) C2d = BRep_Tool::CurveOnSurface(E, theF, f, l);
483       gp_Vec2d d1;
484       gp_Pnt2d p;
485       C2d->D1( u, p, d1 );
486       gp_Vec2d n( d1.Y() * theScaleX, -d1.X() * theScaleY);
487       if ( E.Orientation() == TopAbs_REVERSED )
488         n.Reverse();
489       N += n.Normalized();
490     }
491   }
492
493   // define step size by which to move theUV
494
495   gp_Pnt2d nextUV; // uv of next node on edge, most distant of the four
496   gp_Pnt2d thisUV( theUV.x, theUV.y );
497   double maxDist = -DBL_MAX;
498   TopTools_ListIteratorOfListOfShape aEIt (EList);
499   list< double >::iterator aUIt = UList.begin();
500   for ( ; aEIt.More(); aEIt.Next(), aUIt++ )
501   {
502     const TopoDS_Edge& E = TopoDS::Edge( aEIt.Value() );
503     double f, l;
504     Handle(Geom2d_Curve) C2d = BRep_Tool::CurveOnSurface(E, theF, f, l);
505
506     double umin = DBL_MAX, umax = -DBL_MAX;
507     SMDS_NodeIteratorPtr nIt = theMesh.GetSubMesh(E)->GetSubMeshDS()->GetNodes();
508     if ( !nIt->more() ) // no nodes on edge, only on vertices
509     {
510       umin = l;
511       umax = f;
512     }
513     else {
514       while ( nIt->more() ) {
515         const SMDS_MeshNode* node = nIt->next();
516         // check if node is medium
517         if ( theCreateQuadratic && SMESH_MesherHelper::IsMedium( node, SMDSAbs_Edge ))
518           continue;
519         SMDS_EdgePositionPtr epos = node->GetPosition();
520         double u = epos->GetUParameter();
521         if ( u < umin )
522           umin = u;
523         if ( u > umax )
524           umax = u;
525       }
526     }
527     bool isFirstCommon = ( *aUIt == f );
528     gp_Pnt2d uv = C2d->Value( isFirstCommon ? umin : umax );
529     double dist = thisUV.SquareDistance( uv );
530     if ( dist > maxDist ) {
531       maxDist = dist;
532       nextUV  = uv;
533     }
534   }
535   R2 uv0, uv1, uv2;
536   uv0.x = thisUV.X();   uv0.y = thisUV.Y();
537   uv1.x = nextUV.X();   uv1.y = nextUV.Y();
538   uv2.x = thisUV.X();   uv2.y = thisUV.Y();
539
540   uv1.x *= theScaleX;   uv1.y *= theScaleY;
541
542   if ( fixOverlappedLinkUV( uv0, uv1, uv2 ))
543   {
544     double step = thisUV.Distance( gp_Pnt2d( uv0.x, uv0.y ));
545
546     // move theUV along the normal by the step
547
548     N *= step;
549
550     MESSAGE("--fixCommonVertexUV move(" << theUV.x << " " << theUV.x
551             << ") by (" << N.X() << " " << N.Y() << ")" 
552             << endl << "--- MAX DIST " << maxDist);
553
554     theUV.x += N.X();
555     theUV.y += N.Y();
556
557     return true;
558   }
559   return false;
560 }
561
562 //=============================================================================
563 /*!
564  *  
565  */
566 //=============================================================================
567
568 bool StdMeshers_MEFISTO_2D::LoadPoints(TWireVector &                 wires,
569                                        R2 *                          uvslf,
570                                        vector<const SMDS_MeshNode*>& mefistoToDS,
571                                        double                        scalex,
572                                        double                        scaley)
573 {
574   // to avoid passing same uv points for a vertex common to 2 wires
575   TopoDS_Face F;
576   TopTools_IndexedDataMapOfShapeListOfShape VWMap;
577   if ( wires.size() > 1 )
578   {
579     F = TopoDS::Face( _helper->GetSubShape() );
580     TopExp::MapShapesAndAncestors( F, TopAbs_VERTEX, TopAbs_WIRE, VWMap );
581     int nbVertices = 0;
582     for ( size_t iW = 0; iW < wires.size(); ++iW )
583       nbVertices += wires[ iW ]->NbEdges();
584     if ( nbVertices == VWMap.Extent() )
585       VWMap.Clear(); // wires have no common vertices
586   }
587
588   int m = 0;
589
590   for ( size_t iW = 0; iW < wires.size(); ++iW )
591   {
592     const vector<UVPtStruct>& uvPtVec = wires[ iW ]->GetUVPtStruct();
593     if ((int) uvPtVec.size() != wires[ iW ]->NbPoints() ) {
594       return error(COMPERR_BAD_INPUT_MESH,SMESH_Comment("Unexpected nb of points on wire ")
595                    << iW << ": " << uvPtVec.size()<<" != "<<wires[ iW ]->NbPoints()
596                    << ", probably because of invalid node parameters on geom edges");
597     }
598     if ( m + uvPtVec.size()-1 > mefistoToDS.size() ) {
599       MESSAGE("Wrong mefistoToDS.size: "<<mefistoToDS.size()<<" < "<<m + uvPtVec.size()-1);
600       return error("Internal error");
601     }
602
603     list< int > mOnVertex;
604     vector<UVPtStruct>::const_iterator uvPt = uvPtVec.begin();
605     for ( ++uvPt; uvPt != uvPtVec.end(); ++uvPt )
606     {
607       // bind mefisto ID to node
608       mefistoToDS[m] = uvPt->node;
609       // set UV
610       uvslf[m].x = uvPt->u * scalex;
611       uvslf[m].y = uvPt->v * scaley;
612       switch ( uvPt->node->GetPosition()->GetTypeOfPosition())
613       {
614       case SMDS_TOP_VERTEX:
615         mOnVertex.push_back( m );
616         break;
617       case SMDS_TOP_EDGE:
618         // In order to detect degenerated faces easily, we replace
619         // nodes on a degenerated edge by node on the vertex of that edge
620         if ( _helper->IsDegenShape( uvPt->node->getshapeId() ))
621         {
622           int edgeID = uvPt->node->getshapeId();
623           SMESH_subMesh* edgeSM = _helper->GetMesh()->GetSubMeshContaining( edgeID );
624           SMESH_subMeshIteratorPtr smIt = edgeSM->getDependsOnIterator( /*includeSelf=*/0,
625                                                                         /*complexShapeFirst=*/0);
626           if ( smIt->more() )
627           {
628             SMESH_subMesh* vertexSM = smIt->next();
629             SMDS_NodeIteratorPtr nIt = vertexSM->GetSubMeshDS()->GetNodes();
630             if ( nIt->more() )
631               mefistoToDS[m] = nIt->next();
632           }
633         }
634         break;
635       default:;
636       }
637       m++;
638     }
639
640     int mFirst = 0, mLast = 0;
641     if (!mOnVertex.empty()) {
642       mFirst = mOnVertex.front();
643       mLast = m - 1;
644     }
645     list< int >::iterator mIt = mOnVertex.begin();
646     for ( ; mIt != mOnVertex.end(); ++mIt)
647     {
648       int m = *mIt;
649       if ( iW && !VWMap.IsEmpty()) { // except outer wire
650         // avoid passing same uv point for a vertex common to 2 wires
651         int vID = mefistoToDS[m]->getshapeId();
652         TopoDS_Vertex V = TopoDS::Vertex( _helper->GetMeshDS()->IndexToShape( vID ));
653         if ( fixCommonVertexUV( uvslf[m], V, F, VWMap, *_helper->GetMesh(),
654                                 scalex, scaley, _quadraticMesh )) {
655           myNodesOnCommonV.push_back( mefistoToDS[m] );
656           continue;
657         }
658       }
659       // prevent failure on overlapped adjacent links,
660       // check only links ending in vertex nodes
661       int mB = m - 1, mA = m + 1; // indices Before and After
662       if ( mB < mFirst ) mB = mLast;
663       if ( mA > mLast )  mA = mFirst;
664       fixOverlappedLinkUV (uvslf[ mB ], uvslf[ m ], uvslf[ mA ]);
665     }
666   }
667
668 #ifdef DUMP_POINTS
669   cout << "MEFISTO INPUT************" << endl;
670   for ( int i =0; i < m; ++i )
671     cout << i << ": \t" << uvslf[i].x << ", " << uvslf[i].y
672          << " Node " << mefistoToDS[i]->GetID()<< endl;
673 #endif
674
675   return true;
676 }
677
678 //=============================================================================
679 /*!
680  *  
681  */
682 //=============================================================================
683
684 void StdMeshers_MEFISTO_2D::ComputeScaleOnFace(SMESH_Mesh &        /*aMesh*/,
685                                                const TopoDS_Face & aFace,
686                                                double &            scalex,
687                                                double &            scaley)
688 {
689   TopoDS_Wire W = BRepTools::OuterWire(aFace);
690
691   double xmin = 1.e300;         // min & max of face 2D parametric coord.
692   double xmax = -1.e300;
693   double ymin = 1.e300;
694   double ymax = -1.e300;
695   const int nbp = 23;
696   scalex = 1;
697   scaley = 1;
698
699   TopExp_Explorer wexp(W, TopAbs_EDGE);
700   for ( ; wexp.More(); wexp.Next())
701   {
702     const TopoDS_Edge & E = TopoDS::Edge( wexp.Current() );
703     double f, l;
704     Handle(Geom2d_Curve) C2d = BRep_Tool::CurveOnSurface(E, aFace, f, l);
705     if ( C2d.IsNull() ) continue;
706     double du = (l - f) / double (nbp);
707     for (int i = 0; i <= nbp; i++)
708     {
709       double param = f + double (i) * du;
710       gp_Pnt2d p = C2d->Value(param);
711       if (p.X() < xmin)
712         xmin = p.X();
713       if (p.X() > xmax)
714         xmax = p.X();
715       if (p.Y() < ymin)
716         ymin = p.Y();
717       if (p.Y() > ymax)
718         ymax = p.Y();
719     }
720   }
721   double xmoy = (xmax + xmin) / 2.;
722   double ymoy = (ymax + ymin) / 2.;
723   double xsize = xmax - xmin;
724   double ysize = ymax - ymin;
725
726   TopLoc_Location L;
727   Handle(Geom_Surface) S = BRep_Tool::Surface(aFace,L);       // 3D surface
728
729   double length_x = 0;
730   double length_y = 0;
731   gp_Pnt PX0 = S->Value(xmin, ymoy);
732   gp_Pnt PY0 = S->Value(xmoy, ymin);
733   double dx = xsize / double (nbp);
734   double dy = ysize / double (nbp);
735   for (int i = 1; i <= nbp; i++)
736   {
737     double x = xmin + double (i) * dx;
738     gp_Pnt PX = S->Value(x, ymoy);
739     double y = ymin + double (i) * dy;
740     gp_Pnt PY = S->Value(xmoy, y);
741     length_x += PX.Distance(PX0);
742     length_y += PY.Distance(PY0);
743     PX0 = PX;
744     PY0 = PY;
745   }
746   scalex = length_x / xsize;
747   scaley = length_y / ysize;
748   double xyratio = xsize*scalex/(ysize*scaley);
749   const double maxratio = 1.e2;
750   if (xyratio > maxratio) {
751     scaley *= xyratio / maxratio;
752   }
753   else if (xyratio < 1./maxratio) {
754     scalex *= 1 / xyratio / maxratio;
755   }
756 }
757
758 // namespace
759 // {
760 //   bool isDegenTria( const SMDS_MeshNode * nn[3] )
761 //   {
762 //     SMESH_TNodeXYZ p1( nn[0] );
763 //     SMESH_TNodeXYZ p2( nn[1] );
764 //     SMESH_TNodeXYZ p3( nn[2] );
765 //     gp_XYZ vec1 = p2 - p1;
766 //     gp_XYZ vec2 = p3 - p1;
767 //     gp_XYZ cross = vec1 ^ vec2;
768 //     const double eps = 1e-100;
769 //     return ( fabs( cross.X() ) < eps &&
770 //              fabs( cross.Y() ) < eps &&
771 //              fabs( cross.Z() ) < eps );
772 //   }
773 // }
774
775 //=============================================================================
776 /*!
777  *  
778  */
779 //=============================================================================
780
781 void StdMeshers_MEFISTO_2D::StoreResult(Z nbst, R2 * uvst, Z nbt, Z * nust,
782                                         vector< const SMDS_MeshNode*>&mefistoToDS,
783                                         double scalex, double scaley)
784 {
785   _helper->SetElementsOnShape( true );
786
787   TopoDS_Face F = TopoDS::Face( _helper->GetSubShape() );
788   Handle(Geom_Surface) S = BRep_Tool::Surface( F );
789
790   //const size_t nbInputNodes = mefistoToDS.size();
791
792   Z n = mefistoToDS.size(); // nb input points
793   mefistoToDS.resize( nbst );
794   for ( ; n < nbst; n++)
795   {
796     if (!mefistoToDS[n])
797     {
798       double u = uvst[n][0] / scalex;
799       double v = uvst[n][1] / scaley;
800       gp_Pnt P = S->Value(u, v);
801
802       mefistoToDS[n] = _helper->AddNode( P.X(), P.Y(), P.Z(), 0, u, v );
803     }
804   }
805
806   Z m = 0;
807
808   // triangle points must be in trigonometric order if face is Forward
809   // else they must be put clockwise
810
811   int i1 = 1, i2 = 2;
812   if ( F.Orientation() != TopAbs_FORWARD )
813     std::swap( i1, i2 );
814
815   const SMDS_MeshNode * nn[3];
816   for (n = 1; n <= nbt; n++)
817   {
818     // const bool allNodesAreOld = ( nust[m + 0] <= nbInputNodes &&
819     //                               nust[m + 1] <= nbInputNodes &&
820     //                               nust[m + 2] <= nbInputNodes );
821     nn[ 0 ] = mefistoToDS[ nust[m++] - 1 ];
822     nn[ 1 ] = mefistoToDS[ nust[m++] - 1 ];
823     nn[ 2 ] = mefistoToDS[ nust[m++] - 1 ];
824     m++;
825
826     // avoid creating degenetrated faces
827     bool isDegen = ( _helper->HasDegeneratedEdges() &&
828                      ( nn[0] == nn[1] || nn[1] == nn[2] || nn[2] == nn[0] ));
829
830     // It was an attempt to fix a problem of a zero area face whose all nodes
831     // are on one staight EDGE. But omitting this face makes a hole in the mesh :(
832     // if ( !isDegen && allNodesAreOld )
833     //   isDegen = isDegenTria( nn );
834
835     if ( !isDegen )
836       _helper->AddFace( nn[0], nn[i1], nn[i2] );
837   }
838
839   // remove bad elements built on vertices shared by wires
840
841   list<const SMDS_MeshNode*>::iterator itN = myNodesOnCommonV.begin();
842   for ( ; itN != myNodesOnCommonV.end(); itN++ )
843   {
844     const SMDS_MeshNode* node = *itN;
845     SMDS_ElemIteratorPtr invElemIt = node->GetInverseElementIterator();
846     while ( invElemIt->more() )
847     {
848       const SMDS_MeshElement* elem = invElemIt->next();
849       SMDS_ElemIteratorPtr itN = elem->nodesIterator();
850       int nbSame = 0;
851       while ( itN->more() )
852         if ( itN->next() == node)
853           nbSame++;
854       if (nbSame > 1) {
855         MESSAGE( "RM bad element " << elem->GetID());
856         _helper->GetMeshDS()->RemoveElement( elem );
857       }
858     }
859   }
860 }