1 // SMESH SMESH : implementaion of SMESH idl descriptions
3 // Copyright (C) 2003 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
6 // This library is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License.
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14 // Lesser General Public License for more details.
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
24 // File : StdMeshers_Quadrangle_2D.cxx
25 // Moved here from SMESH_Quadrangle_2D.cxx
26 // Author : Paul RASCLE, EDF
30 #include "StdMeshers_Quadrangle_2D.hxx"
32 #include "StdMeshers_FaceSide.hxx"
34 #include "SMESH_Gen.hxx"
35 #include "SMESH_Mesh.hxx"
36 #include "SMESH_subMesh.hxx"
37 #include "SMESH_MesherHelper.hxx"
38 #include "SMESH_Block.hxx"
39 #include "SMESH_Comment.hxx"
41 #include "SMDS_MeshElement.hxx"
42 #include "SMDS_MeshNode.hxx"
43 #include "SMDS_EdgePosition.hxx"
44 #include "SMDS_FacePosition.hxx"
46 #include <BRepAdaptor_Curve.hxx>
47 #include <BRep_Tool.hxx>
48 #include <BRepLProp.hxx>
49 #include <BRepTools.hxx>
50 #include <BRepTools_WireExplorer.hxx>
51 #include <Geom_Surface.hxx>
52 #include <Geom_Curve.hxx>
53 #include <Geom2d_Curve.hxx>
54 #include <GeomAdaptor_Curve.hxx>
55 #include <GCPnts_UniformAbscissa.hxx>
57 #include <Precision.hxx>
58 #include <gp_Pnt2d.hxx>
59 #include <TColStd_ListIteratorOfListOfInteger.hxx>
60 #include <TColStd_SequenceOfReal.hxx>
61 #include <TColgp_SequenceOfXY.hxx>
62 #include <NCollection_DefineArray2.hxx>
64 #include "utilities.h"
65 #include "Utils_ExceptHandlers.hxx"
67 #ifndef StdMeshers_Array2OfNode_HeaderFile
68 #define StdMeshers_Array2OfNode_HeaderFile
69 typedef const SMDS_MeshNode* SMDS_MeshNodePtr;
70 DEFINE_BASECOLLECTION (StdMeshers_BaseCollectionNodePtr, SMDS_MeshNodePtr)
71 DEFINE_ARRAY2(StdMeshers_Array2OfNode,
72 StdMeshers_BaseCollectionNodePtr, SMDS_MeshNodePtr)
78 typedef SMESH_Comment TComm;
80 //=============================================================================
84 //=============================================================================
86 StdMeshers_Quadrangle_2D::StdMeshers_Quadrangle_2D (int hypId, int studyId, SMESH_Gen* gen)
87 : SMESH_2D_Algo(hypId, studyId, gen)
89 MESSAGE("StdMeshers_Quadrangle_2D::StdMeshers_Quadrangle_2D");
90 _name = "Quadrangle_2D";
91 _shapeType = (1 << TopAbs_FACE);
92 _compatibleHypothesis.push_back("QuadranglePreference");
96 //=============================================================================
100 //=============================================================================
102 StdMeshers_Quadrangle_2D::~StdMeshers_Quadrangle_2D()
104 MESSAGE("StdMeshers_Quadrangle_2D::~StdMeshers_Quadrangle_2D");
107 //=============================================================================
111 //=============================================================================
113 bool StdMeshers_Quadrangle_2D::CheckHypothesis
115 const TopoDS_Shape& aShape,
116 SMESH_Hypothesis::Hypothesis_Status& aStatus)
119 aStatus = SMESH_Hypothesis::HYP_OK;
121 // there is only one compatible Hypothesis so far
122 const list <const SMESHDS_Hypothesis * >&hyps = GetUsedHypothesis(aMesh, aShape);
123 myQuadranglePreference = hyps.size() > 0;
128 //=============================================================================
132 //=============================================================================
134 bool StdMeshers_Quadrangle_2D::Compute (SMESH_Mesh& aMesh,
135 const TopoDS_Shape& aShape)// throw (SALOME_Exception)
137 // PAL14921. Enable catching std::bad_alloc and Standard_OutOfMemory outside
138 //Unexpect aCatchSalomeException);
140 SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
141 aMesh.GetSubMesh(aShape);
143 SMESH_MesherHelper helper(aMesh);
146 _quadraticMesh = myTool->IsQuadraticSubMesh(aShape);
148 FaceQuadStruct *quad = CheckNbEdges( aMesh, aShape );
149 std::auto_ptr<FaceQuadStruct> quadDeleter( quad ); // to delete quad at exit from Compute()
153 if(myQuadranglePreference) {
154 int n1 = quad->side[0]->NbPoints();
155 int n2 = quad->side[1]->NbPoints();
156 int n3 = quad->side[2]->NbPoints();
157 int n4 = quad->side[3]->NbPoints();
158 int nfull = n1+n2+n3+n4;
161 if( nfull==ntmp && ( (n1!=n3) || (n2!=n4) ) ) {
162 // special path for using only quandrangle faces
163 bool ok = ComputeQuadPref(aMesh, aShape, quad);
168 // set normalized grid on unit square in parametric domain
170 if (!SetNormalizedGrid(aMesh, aShape, quad))
173 // --- compute 3D values on points, store points & quadrangles
175 int nbdown = quad->side[0]->NbPoints();
176 int nbup = quad->side[2]->NbPoints();
178 int nbright = quad->side[1]->NbPoints();
179 int nbleft = quad->side[3]->NbPoints();
181 int nbhoriz = Min(nbdown, nbup);
182 int nbvertic = Min(nbright, nbleft);
184 const TopoDS_Face& F = TopoDS::Face(aShape);
185 Handle(Geom_Surface) S = BRep_Tool::Surface(F);
187 // internal mesh nodes
188 int i, j, geomFaceID = meshDS->ShapeToIndex( F );
189 for (i = 1; i < nbhoriz - 1; i++) {
190 for (j = 1; j < nbvertic - 1; j++) {
191 int ij = j * nbhoriz + i;
192 double u = quad->uv_grid[ij].u;
193 double v = quad->uv_grid[ij].v;
194 gp_Pnt P = S->Value(u, v);
195 SMDS_MeshNode * node = meshDS->AddNode(P.X(), P.Y(), P.Z());
196 meshDS->SetNodeOnFace(node, geomFaceID, u, v);
197 quad->uv_grid[ij].node = node;
204 // --.--.--.--.--.-- nbvertic
210 // ---.----.----.--- 0
211 // 0 > > > > > > > > nbhoriz
217 int iup = nbhoriz - 1;
218 if (quad->isEdgeOut[3]) { ilow++; } else { if (quad->isEdgeOut[1]) iup--; }
221 int jup = nbvertic - 1;
222 if (quad->isEdgeOut[0]) { jlow++; } else { if (quad->isEdgeOut[2]) jup--; }
224 // regular quadrangles
225 for (i = ilow; i < iup; i++) {
226 for (j = jlow; j < jup; j++) {
227 const SMDS_MeshNode *a, *b, *c, *d;
228 a = quad->uv_grid[j * nbhoriz + i].node;
229 b = quad->uv_grid[j * nbhoriz + i + 1].node;
230 c = quad->uv_grid[(j + 1) * nbhoriz + i + 1].node;
231 d = quad->uv_grid[(j + 1) * nbhoriz + i].node;
232 SMDS_MeshFace* face = myTool->AddFace(a, b, c, d);
233 meshDS->SetMeshElementOnShape(face, geomFaceID);
237 const vector<UVPtStruct>& uv_e0 = quad->side[0]->GetUVPtStruct(true,0 );
238 const vector<UVPtStruct>& uv_e1 = quad->side[1]->GetUVPtStruct(false,1);
239 const vector<UVPtStruct>& uv_e2 = quad->side[2]->GetUVPtStruct(true,1 );
240 const vector<UVPtStruct>& uv_e3 = quad->side[3]->GetUVPtStruct(false,0);
242 double eps = Precision::Confusion();
244 // Boundary quadrangles
246 if (quad->isEdgeOut[0]) {
249 // |___|___|___|___|___|___|
251 // |___|___|___|___|___|___|
253 // |___|___|___|___|___|___| __ first row of the regular grid
254 // . . . . . . . . . __ down edge nodes
256 // >->->->->->->->->->->->-> -- direction of processing
258 int g = 0; // number of last processed node in the regular grid
260 // number of last node of the down edge to be processed
261 int stop = nbdown - 1;
262 // if right edge is out, we will stop at a node, previous to the last one
263 if (quad->isEdgeOut[1]) stop--;
265 // for each node of the down edge find nearest node
266 // in the first row of the regular grid and link them
267 for (i = 0; i < stop; i++) {
268 const SMDS_MeshNode *a, *b, *c, *d;
270 b = uv_e0[i + 1].node;
271 gp_Pnt pb (b->X(), b->Y(), b->Z());
273 // find node c in the regular grid, which will be linked with node b
276 // right bound reached, link with the rightmost node
278 c = quad->uv_grid[nbhoriz + iup].node;
281 // find in the grid node c, nearest to the b
282 double mind = RealLast();
283 for (int k = g; k <= iup; k++) {
285 const SMDS_MeshNode *nk;
286 if (k < ilow) // this can be, if left edge is out
287 nk = uv_e3[1].node; // get node from the left edge
289 nk = quad->uv_grid[nbhoriz + k].node; // get one of middle nodes
291 gp_Pnt pnk (nk->X(), nk->Y(), nk->Z());
292 double dist = pb.Distance(pnk);
293 if (dist < mind - eps) {
303 if (near == g) { // make triangle
304 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c);
305 SMDS_MeshFace* face = myTool->AddFace(a, b, c);
306 meshDS->SetMeshElementOnShape(face, geomFaceID);
308 else { // make quadrangle
312 d = quad->uv_grid[nbhoriz + near - 1].node;
313 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c, d);
314 SMDS_MeshFace* face = myTool->AddFace(a, b, c, d);
315 meshDS->SetMeshElementOnShape(face, geomFaceID);
317 // if node d is not at position g - make additional triangles
319 for (int k = near - 1; k > g; k--) {
320 c = quad->uv_grid[nbhoriz + k].node;
324 d = quad->uv_grid[nbhoriz + k - 1].node;
325 //SMDS_MeshFace* face = meshDS->AddFace(a, c, d);
326 SMDS_MeshFace* face = myTool->AddFace(a, c, d);
327 meshDS->SetMeshElementOnShape(face, geomFaceID);
334 if (quad->isEdgeOut[2]) {
337 // <-<-<-<-<-<-<-<-<-<-<-<-< -- direction of processing
339 // . . . . . . . . . __ up edge nodes
340 // ___ ___ ___ ___ ___ ___ __ first row of the regular grid
342 // |___|___|___|___|___|___|
344 // |___|___|___|___|___|___|
347 int g = nbhoriz - 1; // last processed node in the regular grid
350 // if left edge is out, we will stop at a second node
351 if (quad->isEdgeOut[3]) stop++;
353 // for each node of the up edge find nearest node
354 // in the first row of the regular grid and link them
355 for (i = nbup - 1; i > stop; i--) {
356 const SMDS_MeshNode *a, *b, *c, *d;
358 b = uv_e2[i - 1].node;
359 gp_Pnt pb (b->X(), b->Y(), b->Z());
361 // find node c in the grid, which will be linked with node b
363 if (i == stop + 1) { // left bound reached, link with the leftmost node
364 c = quad->uv_grid[nbhoriz*(nbvertic - 2) + ilow].node;
367 // find node c in the grid, nearest to the b
368 double mind = RealLast();
369 for (int k = g; k >= ilow; k--) {
370 const SMDS_MeshNode *nk;
372 nk = uv_e1[nbright - 2].node;
374 nk = quad->uv_grid[nbhoriz*(nbvertic - 2) + k].node;
375 gp_Pnt pnk (nk->X(), nk->Y(), nk->Z());
376 double dist = pb.Distance(pnk);
377 if (dist < mind - eps) {
387 if (near == g) { // make triangle
388 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c);
389 SMDS_MeshFace* face = myTool->AddFace(a, b, c);
390 meshDS->SetMeshElementOnShape(face, geomFaceID);
392 else { // make quadrangle
394 d = uv_e1[nbright - 2].node;
396 d = quad->uv_grid[nbhoriz*(nbvertic - 2) + near + 1].node;
397 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c, d);
398 SMDS_MeshFace* face = myTool->AddFace(a, b, c, d);
399 meshDS->SetMeshElementOnShape(face, geomFaceID);
401 if (near + 1 < g) { // if d not is at g - make additional triangles
402 for (int k = near + 1; k < g; k++) {
403 c = quad->uv_grid[nbhoriz*(nbvertic - 2) + k].node;
405 d = uv_e1[nbright - 2].node;
407 d = quad->uv_grid[nbhoriz*(nbvertic - 2) + k + 1].node;
408 //SMDS_MeshFace* face = meshDS->AddFace(a, c, d);
409 SMDS_MeshFace* face = myTool->AddFace(a, c, d);
410 meshDS->SetMeshElementOnShape(face, geomFaceID);
419 // right or left boundary quadrangles
420 if (quad->isEdgeOut[1]) {
421 // MESSAGE("right edge is out");
422 int g = 0; // last processed node in the grid
423 int stop = nbright - 1;
424 if (quad->isEdgeOut[2]) stop--;
425 for (i = 0; i < stop; i++) {
426 const SMDS_MeshNode *a, *b, *c, *d;
428 b = uv_e1[i + 1].node;
429 gp_Pnt pb (b->X(), b->Y(), b->Z());
431 // find node c in the grid, nearest to the b
433 if (i == stop - 1) { // up bondary reached
434 c = quad->uv_grid[nbhoriz*(jup + 1) - 2].node;
437 double mind = RealLast();
438 for (int k = g; k <= jup; k++) {
439 const SMDS_MeshNode *nk;
441 nk = uv_e0[nbdown - 2].node;
443 nk = quad->uv_grid[nbhoriz*(k + 1) - 2].node;
444 gp_Pnt pnk (nk->X(), nk->Y(), nk->Z());
445 double dist = pb.Distance(pnk);
446 if (dist < mind - eps) {
456 if (near == g) { // make triangle
457 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c);
458 SMDS_MeshFace* face = myTool->AddFace(a, b, c);
459 meshDS->SetMeshElementOnShape(face, geomFaceID);
461 else { // make quadrangle
463 d = uv_e0[nbdown - 2].node;
465 d = quad->uv_grid[nbhoriz*near - 2].node;
466 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c, d);
467 SMDS_MeshFace* face = myTool->AddFace(a, b, c, d);
468 meshDS->SetMeshElementOnShape(face, geomFaceID);
470 if (near - 1 > g) { // if d not is at g - make additional triangles
471 for (int k = near - 1; k > g; k--) {
472 c = quad->uv_grid[nbhoriz*(k + 1) - 2].node;
474 d = uv_e0[nbdown - 2].node;
476 d = quad->uv_grid[nbhoriz*k - 2].node;
477 //SMDS_MeshFace* face = meshDS->AddFace(a, c, d);
478 SMDS_MeshFace* face = myTool->AddFace(a, c, d);
479 meshDS->SetMeshElementOnShape(face, geomFaceID);
486 if (quad->isEdgeOut[3]) {
487 // MESSAGE("left edge is out");
488 int g = nbvertic - 1; // last processed node in the grid
490 if (quad->isEdgeOut[0]) stop++;
491 for (i = nbleft - 1; i > stop; i--) {
492 const SMDS_MeshNode *a, *b, *c, *d;
494 b = uv_e3[i - 1].node;
495 gp_Pnt pb (b->X(), b->Y(), b->Z());
497 // find node c in the grid, nearest to the b
499 if (i == stop + 1) { // down bondary reached
500 c = quad->uv_grid[nbhoriz*jlow + 1].node;
503 double mind = RealLast();
504 for (int k = g; k >= jlow; k--) {
505 const SMDS_MeshNode *nk;
509 nk = quad->uv_grid[nbhoriz*k + 1].node;
510 gp_Pnt pnk (nk->X(), nk->Y(), nk->Z());
511 double dist = pb.Distance(pnk);
512 if (dist < mind - eps) {
522 if (near == g) { // make triangle
523 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c);
524 SMDS_MeshFace* face = myTool->AddFace(a, b, c);
525 meshDS->SetMeshElementOnShape(face, geomFaceID);
527 else { // make quadrangle
531 d = quad->uv_grid[nbhoriz*(near + 1) + 1].node;
532 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c, d);
533 SMDS_MeshFace* face = myTool->AddFace(a, b, c, d);
534 meshDS->SetMeshElementOnShape(face, geomFaceID);
536 if (near + 1 < g) { // if d not is at g - make additional triangles
537 for (int k = near + 1; k < g; k++) {
538 c = quad->uv_grid[nbhoriz*k + 1].node;
542 d = quad->uv_grid[nbhoriz*(k + 1) + 1].node;
543 //SMDS_MeshFace* face = meshDS->AddFace(a, c, d);
544 SMDS_MeshFace* face = myTool->AddFace(a, c, d);
545 meshDS->SetMeshElementOnShape(face, geomFaceID);
558 //=============================================================================
562 //=============================================================================
564 FaceQuadStruct* StdMeshers_Quadrangle_2D::CheckNbEdges(SMESH_Mesh & aMesh,
565 const TopoDS_Shape & aShape)
566 //throw(SALOME_Exception)
568 const TopoDS_Face & F = TopoDS::Face(aShape);
569 const bool ignoreMediumNodes = _quadraticMesh;
571 // verify 1 wire only, with 4 edges
573 list< TopoDS_Edge > edges;
574 list< int > nbEdgesInWire;
575 int nbWire = SMESH_Block::GetOrderedEdges (F, V, edges, nbEdgesInWire);
577 error(COMPERR_BAD_SHAPE, TComm("Wrong number of wires: ") << nbWire);
580 FaceQuadStruct* quad = new FaceQuadStruct;
582 quad->side.reserve(nbEdgesInWire.front());
585 list< TopoDS_Edge >::iterator edgeIt = edges.begin();
586 if ( nbEdgesInWire.front() == 4 ) { // exactly 4 edges
587 for ( ; edgeIt != edges.end(); ++edgeIt, nbSides++ )
588 quad->side.push_back( new StdMeshers_FaceSide(F, *edgeIt, &aMesh,
589 nbSides<TOP_SIDE, ignoreMediumNodes));
591 else if ( nbEdgesInWire.front() > 4 ) { // more than 4 edges - try to unite some
592 list< TopoDS_Edge > sideEdges;
593 while ( !edges.empty()) {
595 sideEdges.splice( sideEdges.end(), edges, edges.begin()); // edges.front() -> sideEdges.end()
596 bool sameSide = true;
597 while ( !edges.empty() && sameSide ) {
598 sameSide = SMESH_Algo::IsContinuous( sideEdges.back(), edges.front() );
600 sideEdges.splice( sideEdges.end(), edges, edges.begin());
602 if ( nbSides == 0 ) { // go backward from the first edge
604 while ( !edges.empty() && sameSide ) {
605 sameSide = SMESH_Algo::IsContinuous( sideEdges.front(), edges.back() );
607 sideEdges.splice( sideEdges.begin(), edges, --edges.end());
610 quad->side.push_back( new StdMeshers_FaceSide(F, sideEdges, &aMesh,
611 nbSides<TOP_SIDE, ignoreMediumNodes));
617 cout << endl << "StdMeshers_Quadrangle_2D. Edge IDs of " << nbSides << " sides:";
618 for ( int i = 0; i < nbSides; ++i ) {
620 for ( int e = 0; e < quad->side[i]->NbEdges(); ++e )
621 cout << myTool->GetMeshDS()->ShapeToIndex( quad->side[i]->Edge( e )) << " ";
627 nbSides = nbEdgesInWire.front();
628 error(COMPERR_BAD_SHAPE, TComm("Face must have 4 sides but not ") << nbSides);
636 //=============================================================================
640 //=============================================================================
642 FaceQuadStruct *StdMeshers_Quadrangle_2D::CheckAnd2Dcompute
644 const TopoDS_Shape & aShape,
645 const bool CreateQuadratic) //throw(SALOME_Exception)
647 _quadraticMesh = CreateQuadratic;
649 FaceQuadStruct *quad = CheckNbEdges(aMesh, aShape);
653 // set normalized grid on unit square in parametric domain
654 SetNormalizedGrid(aMesh, aShape, quad);
659 //=============================================================================
663 //=============================================================================
665 faceQuadStruct::~faceQuadStruct()
667 for (int i = 0; i < side.size(); i++) {
668 if (side[i]) delete side[i];
670 if (uv_grid) delete [] uv_grid;
674 inline const vector<UVPtStruct>& GetUVPtStructIn(FaceQuadStruct* quad, int i, int nbSeg)
676 bool isXConst = ( i == BOTTOM_SIDE || i == TOP_SIDE );
677 double constValue = ( i == BOTTOM_SIDE || i == LEFT_SIDE ) ? 0 : 1;
680 quad->side[i]->SimulateUVPtStruct(nbSeg,isXConst,constValue) :
681 quad->side[i]->GetUVPtStruct(isXConst,constValue);
685 //=============================================================================
689 //=============================================================================
691 bool StdMeshers_Quadrangle_2D::SetNormalizedGrid (SMESH_Mesh & aMesh,
692 const TopoDS_Shape& aShape,
693 FaceQuadStruct* & quad) //throw (SALOME_Exception)
695 // Algorithme décrit dans "Génération automatique de maillages"
696 // P.L. GEORGE, MASSON, § 6.4.1 p. 84-85
697 // traitement dans le domaine paramétrique 2d u,v
698 // transport - projection sur le carré unité
700 // MESSAGE("StdMeshers_Quadrangle_2D::SetNormalizedGrid");
701 // const TopoDS_Face& F = TopoDS::Face(aShape);
703 // 1 --- find orientation of the 4 edges, by test on extrema
706 // |<----north-2-------^ a3 -------------> a2
708 // west-3 east-1 =right | |
712 // v----south-0--------> a0 -------------> a1
717 // 3 --- 2D normalized values on unit square [0..1][0..1]
719 int nbhoriz = Min(quad->side[0]->NbPoints(), quad->side[2]->NbPoints());
720 int nbvertic = Min(quad->side[1]->NbPoints(), quad->side[3]->NbPoints());
722 quad->isEdgeOut[0] = (quad->side[0]->NbPoints() > quad->side[2]->NbPoints());
723 quad->isEdgeOut[1] = (quad->side[1]->NbPoints() > quad->side[3]->NbPoints());
724 quad->isEdgeOut[2] = (quad->side[2]->NbPoints() > quad->side[0]->NbPoints());
725 quad->isEdgeOut[3] = (quad->side[3]->NbPoints() > quad->side[1]->NbPoints());
727 UVPtStruct *uv_grid = quad->uv_grid = new UVPtStruct[nbvertic * nbhoriz];
729 const vector<UVPtStruct>& uv_e0 = GetUVPtStructIn( quad, 0, nbhoriz - 1 );
730 const vector<UVPtStruct>& uv_e1 = GetUVPtStructIn( quad, 1, nbvertic - 1 );
731 const vector<UVPtStruct>& uv_e2 = GetUVPtStructIn( quad, 2, nbhoriz - 1 );
732 const vector<UVPtStruct>& uv_e3 = GetUVPtStructIn( quad, 3, nbvertic - 1 );
734 if ( uv_e0.empty() || uv_e1.empty() || uv_e2.empty() || uv_e3.empty() )
735 return error( "Can't find nodes on sides");
737 // nodes Id on "in" edges
738 if (! quad->isEdgeOut[0]) {
740 for (int i = 0; i < nbhoriz; i++) { // down
741 int ij = j * nbhoriz + i;
742 uv_grid[ij].node = uv_e0[i].node;
745 if (! quad->isEdgeOut[1]) {
747 for (int j = 0; j < nbvertic; j++) { // right
748 int ij = j * nbhoriz + i;
749 uv_grid[ij].node = uv_e1[j].node;
752 if (! quad->isEdgeOut[2]) {
753 int j = nbvertic - 1;
754 for (int i = 0; i < nbhoriz; i++) { // up
755 int ij = j * nbhoriz + i;
756 uv_grid[ij].node = uv_e2[i].node;
759 if (! quad->isEdgeOut[3]) {
761 for (int j = 0; j < nbvertic; j++) { // left
762 int ij = j * nbhoriz + i;
763 uv_grid[ij].node = uv_e3[j].node;
767 // normalized 2d values on grid
768 for (int i = 0; i < nbhoriz; i++)
770 for (int j = 0; j < nbvertic; j++)
772 int ij = j * nbhoriz + i;
773 // --- droite i cste : x = x0 + y(x1-x0)
774 double x0 = uv_e0[i].normParam; // bas - sud
775 double x1 = uv_e2[i].normParam; // haut - nord
776 // --- droite j cste : y = y0 + x(y1-y0)
777 double y0 = uv_e3[j].normParam; // gauche-ouest
778 double y1 = uv_e1[j].normParam; // droite - est
779 // --- intersection : x=x0+(y0+x(y1-y0))(x1-x0)
780 double x = (x0 + y0 * (x1 - x0)) / (1 - (y1 - y0) * (x1 - x0));
781 double y = y0 + x * (y1 - y0);
784 //MESSAGE("-xy-01 "<<x0<<" "<<x1<<" "<<y0<<" "<<y1);
785 //MESSAGE("-xy-norm "<<i<<" "<<j<<" "<<x<<" "<<y);
789 // 4 --- projection on 2d domain (u,v)
790 gp_UV a0( uv_e0.front().u, uv_e0.front().v );
791 gp_UV a1( uv_e0.back().u, uv_e0.back().v );
792 gp_UV a2( uv_e2.back().u, uv_e2.back().v );
793 gp_UV a3( uv_e2.front().u, uv_e2.front().v );
795 for (int i = 0; i < nbhoriz; i++)
797 for (int j = 0; j < nbvertic; j++)
799 int ij = j * nbhoriz + i;
800 double x = uv_grid[ij].x;
801 double y = uv_grid[ij].y;
802 double param_0 = uv_e0[0].normParam + x * (uv_e0.back().normParam - uv_e0[0].normParam); // sud
803 double param_2 = uv_e2[0].normParam + x * (uv_e2.back().normParam - uv_e2[0].normParam); // nord
804 double param_1 = uv_e1[0].normParam + y * (uv_e1.back().normParam - uv_e1[0].normParam); // est
805 double param_3 = uv_e3[0].normParam + y * (uv_e3.back().normParam - uv_e3[0].normParam); // ouest
807 //MESSAGE("params "<<param_0<<" "<<param_1<<" "<<param_2<<" "<<param_3);
808 gp_UV p0 = quad->side[0]->Value2d(param_0).XY();
809 gp_UV p1 = quad->side[1]->Value2d(param_1).XY();
810 gp_UV p2 = quad->side[2]->Value2d(param_2).XY();
811 gp_UV p3 = quad->side[3]->Value2d(param_3).XY();
813 gp_UV uv = (1 - y) * p0 + x * p1 + y * p2 + (1 - x) * p3;
814 uv -= (1 - x) * (1 - y) * a0 + x * (1 - y) * a1 + x * y * a2 + (1 - x) * y * a3;
816 uv_grid[ij].u = uv.X();
817 uv_grid[ij].v = uv.Y();
823 //=======================================================================
824 //function : ShiftQuad
825 //purpose : auxilary function for ComputeQuadPref
826 //=======================================================================
828 static void ShiftQuad(FaceQuadStruct* quad, const int num, bool)
830 StdMeshers_FaceSide* side[4] = { quad->side[0], quad->side[1], quad->side[2], quad->side[3] };
831 for (int i = BOTTOM_SIDE; i < NB_SIDES; ++i ) {
832 int id = ( i + num ) % NB_SIDES;
833 bool wasForward = ( i < TOP_SIDE );
834 bool newForward = ( id < TOP_SIDE );
835 if ( wasForward != newForward )
836 side[ i ]->Reverse();
837 quad->side[ id ] = side[ i ];
841 //=======================================================================
843 //purpose : auxilary function for ComputeQuadPref
844 //=======================================================================
846 static gp_UV CalcUV(double x0, double x1, double y0, double y1,
847 FaceQuadStruct* quad,
848 const gp_UV& a0, const gp_UV& a1,
849 const gp_UV& a2, const gp_UV& a3)
851 const vector<UVPtStruct>& uv_eb = quad->side[0]->GetUVPtStruct(true,0 );
852 const vector<UVPtStruct>& uv_er = quad->side[1]->GetUVPtStruct(false,1);
853 const vector<UVPtStruct>& uv_et = quad->side[2]->GetUVPtStruct(true,1 );
854 const vector<UVPtStruct>& uv_el = quad->side[3]->GetUVPtStruct(false,0);
856 double x = (x0 + y0 * (x1 - x0)) / (1 - (y1 - y0) * (x1 - x0));
857 double y = y0 + x * (y1 - y0);
859 double param_b = uv_eb[0].normParam + x * (uv_eb.back().normParam - uv_eb[0].normParam);
860 double param_t = uv_et[0].normParam + x * (uv_et.back().normParam - uv_et[0].normParam);
861 double param_r = uv_er[0].normParam + y * (uv_er.back().normParam - uv_er[0].normParam);
862 double param_l = uv_el[0].normParam + y * (uv_el.back().normParam - uv_el[0].normParam);
864 gp_UV p0 = quad->side[BOTTOM_SIDE]->Value2d(param_b).XY();
865 gp_UV p1 = quad->side[RIGHT_SIDE ]->Value2d(param_r).XY();
866 gp_UV p2 = quad->side[TOP_SIDE ]->Value2d(param_t).XY();
867 gp_UV p3 = quad->side[LEFT_SIDE ]->Value2d(param_l).XY();
869 gp_UV uv = p0 * (1 - y) + p1 * x + p2 * y + p3 * (1 - x);
871 uv -= (1 - x) * (1 - y) * a0 + x * (1 - y) * a1 + x * y * a2 + (1 - x) * y * a3;
876 //=======================================================================
878 * Create only quandrangle faces
880 //=======================================================================
882 bool StdMeshers_Quadrangle_2D::ComputeQuadPref (SMESH_Mesh & aMesh,
883 const TopoDS_Shape& aShape,
884 FaceQuadStruct* quad)
886 SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
887 const TopoDS_Face& F = TopoDS::Face(aShape);
888 Handle(Geom_Surface) S = BRep_Tool::Surface(F);
889 const TopoDS_Wire& W = BRepTools::OuterWire(F);
891 // if(W.Orientation()==TopAbs_FORWARD)
893 //if(WisF) cout<<"W is FORWARD"<<endl;
894 //else cout<<"W is REVERSED"<<endl;
895 // bool FisF = (F.Orientation()==TopAbs_FORWARD);
896 // if(!FisF) WisF = !WisF;
898 int i,j,geomFaceID = meshDS->ShapeToIndex( F );
900 int nb = quad->side[0]->NbPoints();
901 int nr = quad->side[1]->NbPoints();
902 int nt = quad->side[2]->NbPoints();
903 int nl = quad->side[3]->NbPoints();
909 // it is a base case => not shift quad but me be replacement is need
910 ShiftQuad(quad,0,WisF);
913 // we have to shift quad on 2
914 ShiftQuad(quad,2,WisF);
919 // we have to shift quad on 1
920 ShiftQuad(quad,1,WisF);
923 // we have to shift quad on 3
924 ShiftQuad(quad,3,WisF);
928 nb = quad->side[0]->NbPoints();
929 nr = quad->side[1]->NbPoints();
930 nt = quad->side[2]->NbPoints();
931 nl = quad->side[3]->NbPoints();
934 int nbh = Max(nb,nt);
935 int nbv = Max(nr,nl);
939 // orientation of face and 3 main domain for future faces
945 // left | | | | rigth
961 const vector<UVPtStruct>& uv_eb = quad->side[0]->GetUVPtStruct(true,0 );
962 const vector<UVPtStruct>& uv_er = quad->side[1]->GetUVPtStruct(false,1);
963 const vector<UVPtStruct>& uv_et = quad->side[2]->GetUVPtStruct(true,1 );
964 const vector<UVPtStruct>& uv_el = quad->side[3]->GetUVPtStruct(false,0);
966 // arrays for normalized params
967 //cout<<"Dump B:"<<endl;
968 TColStd_SequenceOfReal npb, npr, npt, npl;
969 for(i=0; i<nb; i++) {
970 npb.Append(uv_eb[i].normParam);
971 //cout<<"i="<<i<<" par="<<uv_eb[i].normParam<<" npar="<<uv_eb[i].normParam;
972 //const SMDS_MeshNode* N = uv_eb[i].node;
973 //cout<<" node("<<N->X()<<","<<N->Y()<<","<<N->Z()<<")"<<endl;
975 for(i=0; i<nr; i++) {
976 npr.Append(uv_er[i].normParam);
978 for(i=0; i<nt; i++) {
979 npt.Append(uv_et[i].normParam);
981 for(i=0; i<nl; i++) {
982 npl.Append(uv_el[i].normParam);
985 // add some params to right and left after the first param
988 double dpr = (npr.Value(2) - npr.Value(1))/(dr+1);
989 for(i=1; i<=dr; i++) {
990 npr.InsertAfter(1,npr.Value(2)-dpr);
994 dpr = (npl.Value(2) - npl.Value(1))/(dl+1);
995 for(i=1; i<=dl; i++) {
996 npl.InsertAfter(1,npl.Value(2)-dpr);
999 //for(i=1; i<=npb.Length(); i++) {
1000 // cout<<" "<<npb.Value(i);
1004 gp_XY a0( uv_eb.front().u, uv_eb.front().v );
1005 gp_XY a1( uv_eb.back().u, uv_eb.back().v );
1006 gp_XY a2( uv_et.back().u, uv_et.back().v );
1007 gp_XY a3( uv_et.front().u, uv_et.front().v );
1008 //cout<<" a0("<<a0.X()<<","<<a0.Y()<<")"<<" a1("<<a1.X()<<","<<a1.Y()<<")"
1009 // <<" a2("<<a2.X()<<","<<a2.Y()<<")"<<" a3("<<a3.X()<<","<<a3.Y()<<")"<<endl;
1011 int nnn = Min(nr,nl);
1012 // auxilary sequence of XY for creation nodes
1013 // in the bottom part of central domain
1014 // it's length must be == nbv-nnn-1
1015 TColgp_SequenceOfXY UVL;
1016 TColgp_SequenceOfXY UVR;
1018 // step1: create faces for left domain
1019 StdMeshers_Array2OfNode NodesL(1,dl+1,1,nl);
1021 for(j=1; j<=nl; j++)
1022 NodesL.SetValue(1,j,uv_el[j-1].node);
1025 for(i=1; i<=dl; i++)
1026 NodesL.SetValue(i+1,nl,uv_et[i].node);
1027 // create and add needed nodes
1028 TColgp_SequenceOfXY UVtmp;
1029 for(i=1; i<=dl; i++) {
1030 double x0 = npt.Value(i+1);
1033 double y0 = npl.Value(i+1);
1034 double y1 = npr.Value(i+1);
1035 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1036 gp_Pnt P = S->Value(UV.X(),UV.Y());
1037 SMDS_MeshNode * N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1038 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1039 NodesL.SetValue(i+1,1,N);
1040 if(UVL.Length()<nbv-nnn-1) UVL.Append(UV);
1042 for(j=2; j<nl; j++) {
1043 double y0 = npl.Value(dl+j);
1044 double y1 = npr.Value(dl+j);
1045 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1046 gp_Pnt P = S->Value(UV.X(),UV.Y());
1047 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1048 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1049 NodesL.SetValue(i+1,j,N);
1050 if( i==dl ) UVtmp.Append(UV);
1053 for(i=1; i<=UVtmp.Length() && UVL.Length()<nbv-nnn-1; i++) {
1054 UVL.Append(UVtmp.Value(i));
1056 //cout<<"Dump NodesL:"<<endl;
1057 //for(i=1; i<=dl+1; i++) {
1059 // for(j=1; j<=nl; j++) {
1060 // cout<<" ("<<NodesL.Value(i,j)->X()<<","<<NodesL.Value(i,j)->Y()<<","<<NodesL.Value(i,j)->Z()<<")";
1065 for(i=1; i<=dl; i++) {
1066 for(j=1; j<nl; j++) {
1069 myTool->AddFace(NodesL.Value(i,j), NodesL.Value(i+1,j),
1070 NodesL.Value(i+1,j+1), NodesL.Value(i,j+1));
1071 meshDS->SetMeshElementOnShape(F, geomFaceID);
1075 myTool->AddFace(NodesL.Value(i,j), NodesL.Value(i,j+1),
1076 NodesL.Value(i+1,j+1), NodesL.Value(i+1,j));
1077 meshDS->SetMeshElementOnShape(F, geomFaceID);
1083 // fill UVL using c2d
1084 for(i=1; i<npl.Length() && UVL.Length()<nbv-nnn-1; i++) {
1085 UVL.Append( gp_UV ( uv_el[i].u, uv_el[i].v ));
1089 // step2: create faces for right domain
1090 StdMeshers_Array2OfNode NodesR(1,dr+1,1,nr);
1092 for(j=1; j<=nr; j++)
1093 NodesR.SetValue(1,j,uv_er[nr-j].node);
1096 for(i=1; i<=dr; i++)
1097 NodesR.SetValue(i+1,1,uv_et[nt-1-i].node);
1098 // create and add needed nodes
1099 TColgp_SequenceOfXY UVtmp;
1100 for(i=1; i<=dr; i++) {
1101 double x0 = npt.Value(nt-i);
1104 double y0 = npl.Value(i+1);
1105 double y1 = npr.Value(i+1);
1106 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1107 gp_Pnt P = S->Value(UV.X(),UV.Y());
1108 SMDS_MeshNode * N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1109 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1110 NodesR.SetValue(i+1,nr,N);
1111 if(UVR.Length()<nbv-nnn-1) UVR.Append(UV);
1113 for(j=2; j<nr; j++) {
1114 double y0 = npl.Value(nbv-j+1);
1115 double y1 = npr.Value(nbv-j+1);
1116 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1117 gp_Pnt P = S->Value(UV.X(),UV.Y());
1118 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1119 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1120 NodesR.SetValue(i+1,j,N);
1121 if( i==dr ) UVtmp.Prepend(UV);
1124 for(i=1; i<=UVtmp.Length() && UVR.Length()<nbv-nnn-1; i++) {
1125 UVR.Append(UVtmp.Value(i));
1128 for(i=1; i<=dr; i++) {
1129 for(j=1; j<nr; j++) {
1132 myTool->AddFace(NodesR.Value(i,j), NodesR.Value(i+1,j),
1133 NodesR.Value(i+1,j+1), NodesR.Value(i,j+1));
1134 meshDS->SetMeshElementOnShape(F, geomFaceID);
1138 myTool->AddFace(NodesR.Value(i,j), NodesR.Value(i,j+1),
1139 NodesR.Value(i+1,j+1), NodesR.Value(i+1,j));
1140 meshDS->SetMeshElementOnShape(F, geomFaceID);
1146 // fill UVR using c2d
1147 for(i=1; i<npr.Length() && UVR.Length()<nbv-nnn-1; i++) {
1148 UVR.Append( gp_UV( uv_er[i].u, uv_er[i].v ));
1152 // step3: create faces for central domain
1153 StdMeshers_Array2OfNode NodesC(1,nb,1,nbv);
1154 // add first string using NodesL
1155 for(i=1; i<=dl+1; i++)
1156 NodesC.SetValue(1,i,NodesL(i,1));
1157 for(i=2; i<=nl; i++)
1158 NodesC.SetValue(1,dl+i,NodesL(dl+1,i));
1159 // add last string using NodesR
1160 for(i=1; i<=dr+1; i++)
1161 NodesC.SetValue(nb,i,NodesR(i,nr));
1163 NodesC.SetValue(nb,dr+i+1,NodesR(dr+1,nr-i));
1164 // add top nodes (last columns)
1165 for(i=dl+2; i<nbh-dr; i++)
1166 NodesC.SetValue(i-dl,nbv,uv_et[i-1].node);
1167 // add bottom nodes (first columns)
1169 NodesC.SetValue(i,1,uv_eb[i-1].node);
1171 // create and add needed nodes
1172 // add linear layers
1173 for(i=2; i<nb; i++) {
1174 double x0 = npt.Value(dl+i);
1176 for(j=1; j<nnn; j++) {
1177 double y0 = npl.Value(nbv-nnn+j);
1178 double y1 = npr.Value(nbv-nnn+j);
1179 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1180 gp_Pnt P = S->Value(UV.X(),UV.Y());
1181 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1182 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1183 NodesC.SetValue(i,nbv-nnn+j,N);
1186 // add diagonal layers
1187 //cout<<"UVL.Length()="<<UVL.Length()<<" UVR.Length()="<<UVR.Length()<<endl;
1188 //cout<<"Dump UVL:"<<endl;
1189 //for(i=1; i<=UVL.Length(); i++) {
1190 // cout<<" ("<<UVL.Value(i).X()<<","<<UVL.Value(i).Y()<<")";
1193 for(i=1; i<nbv-nnn; i++) {
1194 double du = UVR.Value(i).X() - UVL.Value(i).X();
1195 double dv = UVR.Value(i).Y() - UVL.Value(i).Y();
1196 for(j=2; j<nb; j++) {
1197 double u = UVL.Value(i).X() + du*npb.Value(j);
1198 double v = UVL.Value(i).Y() + dv*npb.Value(j);
1199 gp_Pnt P = S->Value(u,v);
1200 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1201 meshDS->SetNodeOnFace(N, geomFaceID, u, v);
1202 NodesC.SetValue(j,i+1,N);
1206 for(i=1; i<nb; i++) {
1207 for(j=1; j<nbv; j++) {
1210 myTool->AddFace(NodesC.Value(i,j), NodesC.Value(i+1,j),
1211 NodesC.Value(i+1,j+1), NodesC.Value(i,j+1));
1212 meshDS->SetMeshElementOnShape(F, geomFaceID);
1216 myTool->AddFace(NodesC.Value(i,j), NodesC.Value(i,j+1),
1217 NodesC.Value(i+1,j+1), NodesC.Value(i+1,j));
1218 meshDS->SetMeshElementOnShape(F, geomFaceID);