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"
40 #include "SMDS_MeshElement.hxx"
41 #include "SMDS_MeshNode.hxx"
42 #include "SMDS_EdgePosition.hxx"
43 #include "SMDS_FacePosition.hxx"
45 #include <BRepAdaptor_Curve.hxx>
46 #include <BRep_Tool.hxx>
47 #include <BRepLProp.hxx>
48 #include <BRepTools.hxx>
49 #include <BRepTools_WireExplorer.hxx>
50 #include <Geom_Surface.hxx>
51 #include <Geom_Curve.hxx>
52 #include <Geom2d_Curve.hxx>
53 #include <GeomAdaptor_Curve.hxx>
54 #include <GCPnts_UniformAbscissa.hxx>
56 #include <Precision.hxx>
57 #include <gp_Pnt2d.hxx>
58 #include <TColStd_ListIteratorOfListOfInteger.hxx>
59 #include <TColStd_SequenceOfReal.hxx>
60 #include <TColgp_SequenceOfXY.hxx>
61 #include <NCollection_DefineArray2.hxx>
63 #include "utilities.h"
64 #include "Utils_ExceptHandlers.hxx"
66 #ifndef StdMeshers_Array2OfNode_HeaderFile
67 #define StdMeshers_Array2OfNode_HeaderFile
68 typedef const SMDS_MeshNode* SMDS_MeshNodePtr;
69 DEFINE_BASECOLLECTION (StdMeshers_BaseCollectionNodePtr, SMDS_MeshNodePtr)
70 DEFINE_ARRAY2(StdMeshers_Array2OfNode,
71 StdMeshers_BaseCollectionNodePtr, SMDS_MeshNodePtr)
78 //=============================================================================
82 //=============================================================================
84 StdMeshers_Quadrangle_2D::StdMeshers_Quadrangle_2D (int hypId, int studyId, SMESH_Gen* gen)
85 : SMESH_2D_Algo(hypId, studyId, gen)
87 MESSAGE("StdMeshers_Quadrangle_2D::StdMeshers_Quadrangle_2D");
88 _name = "Quadrangle_2D";
89 _shapeType = (1 << TopAbs_FACE);
90 _compatibleHypothesis.push_back("QuadranglePreference");
94 //=============================================================================
98 //=============================================================================
100 StdMeshers_Quadrangle_2D::~StdMeshers_Quadrangle_2D()
102 MESSAGE("StdMeshers_Quadrangle_2D::~StdMeshers_Quadrangle_2D");
105 //=============================================================================
109 //=============================================================================
111 bool StdMeshers_Quadrangle_2D::CheckHypothesis
113 const TopoDS_Shape& aShape,
114 SMESH_Hypothesis::Hypothesis_Status& aStatus)
117 aStatus = SMESH_Hypothesis::HYP_OK;
119 // there is only one compatible Hypothesis so far
120 const list <const SMESHDS_Hypothesis * >&hyps = GetUsedHypothesis(aMesh, aShape);
121 myQuadranglePreference = hyps.size() > 0;
126 //=============================================================================
130 //=============================================================================
132 bool StdMeshers_Quadrangle_2D::Compute (SMESH_Mesh& aMesh,
133 const TopoDS_Shape& aShape) throw (SALOME_Exception)
135 Unexpect aCatch(SalomeException);
137 SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
138 aMesh.GetSubMesh(aShape);
140 myTool = new SMESH_MesherHelper(aMesh);
141 std::auto_ptr<SMESH_MesherHelper> helperDeleter( myTool );// to delete helper at exit from Compute()
143 _quadraticMesh = myTool->IsQuadraticSubMesh(aShape);
145 FaceQuadStruct *quad = CheckNbEdges( aMesh, aShape );
146 std::auto_ptr<FaceQuadStruct> quadDeleter( quad ); // to delete quad at exit from Compute()
150 if(myQuadranglePreference) {
151 int n1 = quad->side[0]->NbPoints();
152 int n2 = quad->side[1]->NbPoints();
153 int n3 = quad->side[2]->NbPoints();
154 int n4 = quad->side[3]->NbPoints();
155 int nfull = n1+n2+n3+n4;
158 if( nfull==ntmp && ( (n1!=n3) || (n2!=n4) ) ) {
159 // special path for using only quandrangle faces
160 bool ok = ComputeQuadPref(aMesh, aShape, quad);
165 // set normalized grid on unit square in parametric domain
167 if (!SetNormalizedGrid(aMesh, aShape, quad))
170 // --- compute 3D values on points, store points & quadrangles
172 int nbdown = quad->side[0]->NbPoints();
173 int nbup = quad->side[2]->NbPoints();
175 int nbright = quad->side[1]->NbPoints();
176 int nbleft = quad->side[3]->NbPoints();
178 int nbhoriz = Min(nbdown, nbup);
179 int nbvertic = Min(nbright, nbleft);
181 const TopoDS_Face& F = TopoDS::Face(aShape);
182 Handle(Geom_Surface) S = BRep_Tool::Surface(F);
184 // internal mesh nodes
185 int i, j, geomFaceID = meshDS->ShapeToIndex( F );
186 for (i = 1; i < nbhoriz - 1; i++) {
187 for (j = 1; j < nbvertic - 1; j++) {
188 int ij = j * nbhoriz + i;
189 double u = quad->uv_grid[ij].u;
190 double v = quad->uv_grid[ij].v;
191 gp_Pnt P = S->Value(u, v);
192 SMDS_MeshNode * node = meshDS->AddNode(P.X(), P.Y(), P.Z());
193 meshDS->SetNodeOnFace(node, geomFaceID, u, v);
194 quad->uv_grid[ij].node = node;
201 // --.--.--.--.--.-- nbvertic
207 // ---.----.----.--- 0
208 // 0 > > > > > > > > nbhoriz
214 int iup = nbhoriz - 1;
215 if (quad->isEdgeOut[3]) { ilow++; } else { if (quad->isEdgeOut[1]) iup--; }
218 int jup = nbvertic - 1;
219 if (quad->isEdgeOut[0]) { jlow++; } else { if (quad->isEdgeOut[2]) jup--; }
221 // regular quadrangles
222 for (i = ilow; i < iup; i++) {
223 for (j = jlow; j < jup; j++) {
224 const SMDS_MeshNode *a, *b, *c, *d;
225 a = quad->uv_grid[j * nbhoriz + i].node;
226 b = quad->uv_grid[j * nbhoriz + i + 1].node;
227 c = quad->uv_grid[(j + 1) * nbhoriz + i + 1].node;
228 d = quad->uv_grid[(j + 1) * nbhoriz + i].node;
229 SMDS_MeshFace* face = myTool->AddFace(a, b, c, d);
230 meshDS->SetMeshElementOnShape(face, geomFaceID);
234 const vector<UVPtStruct>& uv_e0 = quad->side[0]->GetUVPtStruct(true,0 );
235 const vector<UVPtStruct>& uv_e1 = quad->side[1]->GetUVPtStruct(false,1);
236 const vector<UVPtStruct>& uv_e2 = quad->side[2]->GetUVPtStruct(true,1 );
237 const vector<UVPtStruct>& uv_e3 = quad->side[3]->GetUVPtStruct(false,0);
239 double eps = Precision::Confusion();
241 // Boundary quadrangles
243 if (quad->isEdgeOut[0]) {
246 // |___|___|___|___|___|___|
248 // |___|___|___|___|___|___|
250 // |___|___|___|___|___|___| __ first row of the regular grid
251 // . . . . . . . . . __ down edge nodes
253 // >->->->->->->->->->->->-> -- direction of processing
255 int g = 0; // number of last processed node in the regular grid
257 // number of last node of the down edge to be processed
258 int stop = nbdown - 1;
259 // if right edge is out, we will stop at a node, previous to the last one
260 if (quad->isEdgeOut[1]) stop--;
262 // for each node of the down edge find nearest node
263 // in the first row of the regular grid and link them
264 for (i = 0; i < stop; i++) {
265 const SMDS_MeshNode *a, *b, *c, *d;
267 b = uv_e0[i + 1].node;
268 gp_Pnt pb (b->X(), b->Y(), b->Z());
270 // find node c in the regular grid, which will be linked with node b
273 // right bound reached, link with the rightmost node
275 c = quad->uv_grid[nbhoriz + iup].node;
278 // find in the grid node c, nearest to the b
279 double mind = RealLast();
280 for (int k = g; k <= iup; k++) {
282 const SMDS_MeshNode *nk;
283 if (k < ilow) // this can be, if left edge is out
284 nk = uv_e3[1].node; // get node from the left edge
286 nk = quad->uv_grid[nbhoriz + k].node; // get one of middle nodes
288 gp_Pnt pnk (nk->X(), nk->Y(), nk->Z());
289 double dist = pb.Distance(pnk);
290 if (dist < mind - eps) {
300 if (near == g) { // make triangle
301 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c);
302 SMDS_MeshFace* face = myTool->AddFace(a, b, c);
303 meshDS->SetMeshElementOnShape(face, geomFaceID);
305 else { // make quadrangle
309 d = quad->uv_grid[nbhoriz + near - 1].node;
310 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c, d);
311 SMDS_MeshFace* face = myTool->AddFace(a, b, c, d);
312 meshDS->SetMeshElementOnShape(face, geomFaceID);
314 // if node d is not at position g - make additional triangles
316 for (int k = near - 1; k > g; k--) {
317 c = quad->uv_grid[nbhoriz + k].node;
321 d = quad->uv_grid[nbhoriz + k - 1].node;
322 //SMDS_MeshFace* face = meshDS->AddFace(a, c, d);
323 SMDS_MeshFace* face = myTool->AddFace(a, c, d);
324 meshDS->SetMeshElementOnShape(face, geomFaceID);
331 if (quad->isEdgeOut[2]) {
334 // <-<-<-<-<-<-<-<-<-<-<-<-< -- direction of processing
336 // . . . . . . . . . __ up edge nodes
337 // ___ ___ ___ ___ ___ ___ __ first row of the regular grid
339 // |___|___|___|___|___|___|
341 // |___|___|___|___|___|___|
344 int g = nbhoriz - 1; // last processed node in the regular grid
347 // if left edge is out, we will stop at a second node
348 if (quad->isEdgeOut[3]) stop++;
350 // for each node of the up edge find nearest node
351 // in the first row of the regular grid and link them
352 for (i = nbup - 1; i > stop; i--) {
353 const SMDS_MeshNode *a, *b, *c, *d;
355 b = uv_e2[i - 1].node;
356 gp_Pnt pb (b->X(), b->Y(), b->Z());
358 // find node c in the grid, which will be linked with node b
360 if (i == stop + 1) { // left bound reached, link with the leftmost node
361 c = quad->uv_grid[nbhoriz*(nbvertic - 2) + ilow].node;
364 // find node c in the grid, nearest to the b
365 double mind = RealLast();
366 for (int k = g; k >= ilow; k--) {
367 const SMDS_MeshNode *nk;
369 nk = uv_e1[nbright - 2].node;
371 nk = quad->uv_grid[nbhoriz*(nbvertic - 2) + k].node;
372 gp_Pnt pnk (nk->X(), nk->Y(), nk->Z());
373 double dist = pb.Distance(pnk);
374 if (dist < mind - eps) {
384 if (near == g) { // make triangle
385 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c);
386 SMDS_MeshFace* face = myTool->AddFace(a, b, c);
387 meshDS->SetMeshElementOnShape(face, geomFaceID);
389 else { // make quadrangle
391 d = uv_e1[nbright - 2].node;
393 d = quad->uv_grid[nbhoriz*(nbvertic - 2) + near + 1].node;
394 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c, d);
395 SMDS_MeshFace* face = myTool->AddFace(a, b, c, d);
396 meshDS->SetMeshElementOnShape(face, geomFaceID);
398 if (near + 1 < g) { // if d not is at g - make additional triangles
399 for (int k = near + 1; k < g; k++) {
400 c = quad->uv_grid[nbhoriz*(nbvertic - 2) + k].node;
402 d = uv_e1[nbright - 2].node;
404 d = quad->uv_grid[nbhoriz*(nbvertic - 2) + k + 1].node;
405 //SMDS_MeshFace* face = meshDS->AddFace(a, c, d);
406 SMDS_MeshFace* face = myTool->AddFace(a, c, d);
407 meshDS->SetMeshElementOnShape(face, geomFaceID);
416 // right or left boundary quadrangles
417 if (quad->isEdgeOut[1]) {
418 // MESSAGE("right edge is out");
419 int g = 0; // last processed node in the grid
420 int stop = nbright - 1;
421 if (quad->isEdgeOut[2]) stop--;
422 for (i = 0; i < stop; i++) {
423 const SMDS_MeshNode *a, *b, *c, *d;
425 b = uv_e1[i + 1].node;
426 gp_Pnt pb (b->X(), b->Y(), b->Z());
428 // find node c in the grid, nearest to the b
430 if (i == stop - 1) { // up bondary reached
431 c = quad->uv_grid[nbhoriz*(jup + 1) - 2].node;
434 double mind = RealLast();
435 for (int k = g; k <= jup; k++) {
436 const SMDS_MeshNode *nk;
438 nk = uv_e0[nbdown - 2].node;
440 nk = quad->uv_grid[nbhoriz*(k + 1) - 2].node;
441 gp_Pnt pnk (nk->X(), nk->Y(), nk->Z());
442 double dist = pb.Distance(pnk);
443 if (dist < mind - eps) {
453 if (near == g) { // make triangle
454 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c);
455 SMDS_MeshFace* face = myTool->AddFace(a, b, c);
456 meshDS->SetMeshElementOnShape(face, geomFaceID);
458 else { // make quadrangle
460 d = uv_e0[nbdown - 2].node;
462 d = quad->uv_grid[nbhoriz*near - 2].node;
463 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c, d);
464 SMDS_MeshFace* face = myTool->AddFace(a, b, c, d);
465 meshDS->SetMeshElementOnShape(face, geomFaceID);
467 if (near - 1 > g) { // if d not is at g - make additional triangles
468 for (int k = near - 1; k > g; k--) {
469 c = quad->uv_grid[nbhoriz*(k + 1) - 2].node;
471 d = uv_e0[nbdown - 2].node;
473 d = quad->uv_grid[nbhoriz*k - 2].node;
474 //SMDS_MeshFace* face = meshDS->AddFace(a, c, d);
475 SMDS_MeshFace* face = myTool->AddFace(a, c, d);
476 meshDS->SetMeshElementOnShape(face, geomFaceID);
483 if (quad->isEdgeOut[3]) {
484 // MESSAGE("left edge is out");
485 int g = nbvertic - 1; // last processed node in the grid
487 if (quad->isEdgeOut[0]) stop++;
488 for (i = nbleft - 1; i > stop; i--) {
489 const SMDS_MeshNode *a, *b, *c, *d;
491 b = uv_e3[i - 1].node;
492 gp_Pnt pb (b->X(), b->Y(), b->Z());
494 // find node c in the grid, nearest to the b
496 if (i == stop + 1) { // down bondary reached
497 c = quad->uv_grid[nbhoriz*jlow + 1].node;
500 double mind = RealLast();
501 for (int k = g; k >= jlow; k--) {
502 const SMDS_MeshNode *nk;
506 nk = quad->uv_grid[nbhoriz*k + 1].node;
507 gp_Pnt pnk (nk->X(), nk->Y(), nk->Z());
508 double dist = pb.Distance(pnk);
509 if (dist < mind - eps) {
519 if (near == g) { // make triangle
520 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c);
521 SMDS_MeshFace* face = myTool->AddFace(a, b, c);
522 meshDS->SetMeshElementOnShape(face, geomFaceID);
524 else { // make quadrangle
528 d = quad->uv_grid[nbhoriz*(near + 1) + 1].node;
529 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c, d);
530 SMDS_MeshFace* face = myTool->AddFace(a, b, c, d);
531 meshDS->SetMeshElementOnShape(face, geomFaceID);
533 if (near + 1 < g) { // if d not is at g - make additional triangles
534 for (int k = near + 1; k < g; k++) {
535 c = quad->uv_grid[nbhoriz*k + 1].node;
539 d = quad->uv_grid[nbhoriz*(k + 1) + 1].node;
540 //SMDS_MeshFace* face = meshDS->AddFace(a, c, d);
541 SMDS_MeshFace* face = myTool->AddFace(a, c, d);
542 meshDS->SetMeshElementOnShape(face, geomFaceID);
555 //=============================================================================
559 //=============================================================================
561 FaceQuadStruct* StdMeshers_Quadrangle_2D::CheckNbEdges(SMESH_Mesh & aMesh,
562 const TopoDS_Shape & aShape)
563 throw(SALOME_Exception)
565 Unexpect aCatch(SalomeException);
567 const TopoDS_Face & F = TopoDS::Face(aShape);
568 const bool ignoreMediumNodes = _quadraticMesh;
570 // verify 1 wire only, with 4 edges
572 list< TopoDS_Edge > edges;
573 list< int > nbEdgesInWire;
574 int nbWire = SMESH_Block::GetOrderedEdges (F, V, edges, nbEdgesInWire);
576 INFOS("only 1 wire by face (quadrangles)");
579 FaceQuadStruct* quad = new FaceQuadStruct;
581 for ( int i = 0; i < NB_SIDES; ++i )
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[nbSides] = 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 ( edgeIt != edges.end()) {
595 sideEdges.push_back( *edgeIt++ );
596 bool sameSide = true;
597 while ( edgeIt != edges.end() && sameSide ) {
598 GeomAbs_Shape cont = SMESH_Algo::Continuity( sideEdges.back(), *edgeIt );
599 sameSide = ( cont >= GeomAbs_C1 );
601 sideEdges.push_back( *edgeIt++ );
603 quad->side[nbSides] = new StdMeshers_FaceSide(F, sideEdges, &aMesh,
604 nbSides<TOP_SIDE, ignoreMediumNodes);
609 INFOS("face must have 4 edges / quadrangle");
617 //=============================================================================
621 //=============================================================================
623 FaceQuadStruct *StdMeshers_Quadrangle_2D::CheckAnd2Dcompute
625 const TopoDS_Shape & aShape,
626 const bool CreateQuadratic) throw(SALOME_Exception)
628 Unexpect aCatch(SalomeException);
630 _quadraticMesh = CreateQuadratic;
632 FaceQuadStruct *quad = CheckNbEdges(aMesh, aShape);
636 // set normalized grid on unit square in parametric domain
637 SetNormalizedGrid(aMesh, aShape, quad);
642 //=============================================================================
646 //=============================================================================
648 faceQuadStruct::~faceQuadStruct()
650 for (int i = 0; i < 4; i++) {
651 if (side[i]) delete side[i];
652 //if (uv_edges[i]) delete [] uv_edges[i];
654 if (uv_grid) delete [] uv_grid;
658 inline const vector<UVPtStruct>& GetUVPtStructIn(FaceQuadStruct* quad, int i, int nbSeg)
660 bool isXConst = ( i == BOTTOM_SIDE || i == TOP_SIDE );
661 double constValue = ( i == BOTTOM_SIDE || i == LEFT_SIDE ) ? 0 : 1;
664 quad->side[i]->SimulateUVPtStruct(nbSeg,isXConst,constValue) :
665 quad->side[i]->GetUVPtStruct(isXConst,constValue);
669 //=============================================================================
673 //=============================================================================
675 bool StdMeshers_Quadrangle_2D::SetNormalizedGrid (SMESH_Mesh & aMesh,
676 const TopoDS_Shape& aShape,
677 FaceQuadStruct* & quad) throw (SALOME_Exception)
679 Unexpect aCatch(SalomeException);
680 // Algorithme décrit dans "Génération automatique de maillages"
681 // P.L. GEORGE, MASSON, § 6.4.1 p. 84-85
682 // traitement dans le domaine paramétrique 2d u,v
683 // transport - projection sur le carré unité
685 // MESSAGE("StdMeshers_Quadrangle_2D::SetNormalizedGrid");
686 // const TopoDS_Face& F = TopoDS::Face(aShape);
688 // 1 --- find orientation of the 4 edges, by test on extrema
691 // |<----north-2-------^ a3 -------------> a2
693 // west-3 east-1 =right | |
697 // v----south-0--------> a0 -------------> a1
702 // 3 --- 2D normalized values on unit square [0..1][0..1]
704 int nbhoriz = Min(quad->side[0]->NbPoints(), quad->side[2]->NbPoints());
705 int nbvertic = Min(quad->side[1]->NbPoints(), quad->side[3]->NbPoints());
707 quad->isEdgeOut[0] = (quad->side[0]->NbPoints() > quad->side[2]->NbPoints());
708 quad->isEdgeOut[1] = (quad->side[1]->NbPoints() > quad->side[3]->NbPoints());
709 quad->isEdgeOut[2] = (quad->side[2]->NbPoints() > quad->side[0]->NbPoints());
710 quad->isEdgeOut[3] = (quad->side[3]->NbPoints() > quad->side[1]->NbPoints());
712 UVPtStruct *uv_grid = quad->uv_grid = new UVPtStruct[nbvertic * nbhoriz];
714 const vector<UVPtStruct>& uv_e0 = GetUVPtStructIn( quad, 0, nbhoriz - 1 );
715 const vector<UVPtStruct>& uv_e1 = GetUVPtStructIn( quad, 1, nbvertic - 1 );
716 const vector<UVPtStruct>& uv_e2 = GetUVPtStructIn( quad, 2, nbhoriz - 1 );
717 const vector<UVPtStruct>& uv_e3 = GetUVPtStructIn( quad, 3, nbvertic - 1 );
719 if ( uv_e0.empty() || uv_e1.empty() || uv_e2.empty() || uv_e3.empty() ) {
720 MESSAGE("Nodes from edges not loaded");
724 // nodes Id on "in" edges
725 if (! quad->isEdgeOut[0]) {
727 for (int i = 0; i < nbhoriz; i++) { // down
728 int ij = j * nbhoriz + i;
729 uv_grid[ij].node = uv_e0[i].node;
732 if (! quad->isEdgeOut[1]) {
734 for (int j = 0; j < nbvertic; j++) { // right
735 int ij = j * nbhoriz + i;
736 uv_grid[ij].node = uv_e1[j].node;
739 if (! quad->isEdgeOut[2]) {
740 int j = nbvertic - 1;
741 for (int i = 0; i < nbhoriz; i++) { // up
742 int ij = j * nbhoriz + i;
743 uv_grid[ij].node = uv_e2[i].node;
746 if (! quad->isEdgeOut[3]) {
748 for (int j = 0; j < nbvertic; j++) { // left
749 int ij = j * nbhoriz + i;
750 uv_grid[ij].node = uv_e3[j].node;
754 // normalized 2d values on grid
755 for (int i = 0; i < nbhoriz; i++)
757 for (int j = 0; j < nbvertic; j++)
759 int ij = j * nbhoriz + i;
760 // --- droite i cste : x = x0 + y(x1-x0)
761 double x0 = uv_e0[i].normParam; // bas - sud
762 double x1 = uv_e2[i].normParam; // haut - nord
763 // --- droite j cste : y = y0 + x(y1-y0)
764 double y0 = uv_e3[j].normParam; // gauche-ouest
765 double y1 = uv_e1[j].normParam; // droite - est
766 // --- intersection : x=x0+(y0+x(y1-y0))(x1-x0)
767 double x = (x0 + y0 * (x1 - x0)) / (1 - (y1 - y0) * (x1 - x0));
768 double y = y0 + x * (y1 - y0);
771 //MESSAGE("-xy-01 "<<x0<<" "<<x1<<" "<<y0<<" "<<y1);
772 //MESSAGE("-xy-norm "<<i<<" "<<j<<" "<<x<<" "<<y);
776 // 4 --- projection on 2d domain (u,v)
777 gp_UV a0( uv_e0.front().u, uv_e0.front().v );
778 gp_UV a1( uv_e0.back().u, uv_e0.back().v );
779 gp_UV a2( uv_e2.back().u, uv_e2.back().v );
780 gp_UV a3( uv_e2.front().u, uv_e2.front().v );
782 for (int i = 0; i < nbhoriz; i++)
784 for (int j = 0; j < nbvertic; j++)
786 int ij = j * nbhoriz + i;
787 double x = uv_grid[ij].x;
788 double y = uv_grid[ij].y;
789 double param_0 = uv_e0[0].normParam + x * (uv_e0.back().normParam - uv_e0[0].normParam); // sud
790 double param_2 = uv_e2[0].normParam + x * (uv_e2.back().normParam - uv_e2[0].normParam); // nord
791 double param_1 = uv_e1[0].normParam + y * (uv_e1.back().normParam - uv_e1[0].normParam); // est
792 double param_3 = uv_e3[0].normParam + y * (uv_e3.back().normParam - uv_e3[0].normParam); // ouest
794 //MESSAGE("params "<<param_0<<" "<<param_1<<" "<<param_2<<" "<<param_3);
795 gp_UV p0 = quad->side[0]->Value2d(param_0).XY();
796 gp_UV p1 = quad->side[1]->Value2d(param_1).XY();
797 gp_UV p2 = quad->side[2]->Value2d(param_2).XY();
798 gp_UV p3 = quad->side[3]->Value2d(param_3).XY();
800 gp_UV uv = (1 - y) * p0 + x * p1 + y * p2 + (1 - x) * p3;
801 uv -= (1 - x) * (1 - y) * a0 + x * (1 - y) * a1 + x * y * a2 + (1 - x) * y * a3;
803 uv_grid[ij].u = uv.X();
804 uv_grid[ij].v = uv.Y();
810 //=======================================================================
811 //function : ShiftQuad
812 //purpose : auxilary function for ComputeQuadPref
813 //=======================================================================
815 static void ShiftQuad(FaceQuadStruct* quad, const int num, bool)
817 StdMeshers_FaceSide* side[4] = { quad->side[0], quad->side[1], quad->side[2], quad->side[3] };
818 for (int i = BOTTOM_SIDE; i < NB_SIDES; ++i ) {
819 int id = ( i + num ) % NB_SIDES;
820 bool wasForward = ( i < TOP_SIDE );
821 bool newForward = ( id < TOP_SIDE );
822 if ( wasForward != newForward )
823 side[ i ]->Reverse();
824 quad->side[ id ] = side[ i ];
828 //=======================================================================
830 //purpose : auxilary function for ComputeQuadPref
831 //=======================================================================
833 static gp_UV CalcUV(double x0, double x1, double y0, double y1,
834 FaceQuadStruct* quad,
835 const gp_UV& a0, const gp_UV& a1,
836 const gp_UV& a2, const gp_UV& a3)
838 const vector<UVPtStruct>& uv_eb = quad->side[0]->GetUVPtStruct(true,0 );
839 const vector<UVPtStruct>& uv_er = quad->side[1]->GetUVPtStruct(false,1);
840 const vector<UVPtStruct>& uv_et = quad->side[2]->GetUVPtStruct(true,1 );
841 const vector<UVPtStruct>& uv_el = quad->side[3]->GetUVPtStruct(false,0);
843 double x = (x0 + y0 * (x1 - x0)) / (1 - (y1 - y0) * (x1 - x0));
844 double y = y0 + x * (y1 - y0);
846 double param_b = uv_eb[0].normParam + x * (uv_eb.back().normParam - uv_eb[0].normParam);
847 double param_t = uv_et[0].normParam + x * (uv_et.back().normParam - uv_et[0].normParam);
848 double param_r = uv_er[0].normParam + y * (uv_er.back().normParam - uv_er[0].normParam);
849 double param_l = uv_el[0].normParam + y * (uv_el.back().normParam - uv_el[0].normParam);
851 gp_UV p0 = quad->side[BOTTOM_SIDE]->Value2d(param_b).XY();
852 gp_UV p1 = quad->side[RIGHT_SIDE ]->Value2d(param_r).XY();
853 gp_UV p2 = quad->side[TOP_SIDE ]->Value2d(param_t).XY();
854 gp_UV p3 = quad->side[LEFT_SIDE ]->Value2d(param_l).XY();
856 gp_UV uv = p0 * (1 - y) + p1 * x + p2 * y + p3 * (1 - x);
858 uv -= (1 - x) * (1 - y) * a0 + x * (1 - y) * a1 + x * y * a2 + (1 - x) * y * a3;
860 //cout<<"x0="<<x0<<" x1="<<x1<<" y0="<<y0<<" y1="<<y1<<endl;
861 //cout<<"x="<<x<<" y="<<y<<endl;
862 //cout<<"param_b="<<param_b<<" param_t="<<param_t<<" param_r="<<param_r<<" param_l="<<param_l<<endl;
863 //cout<<"u="<<u<<" v="<<v<<endl;
868 //=======================================================================
869 //function : ComputeQuadPref
871 //=======================================================================
873 * Special function for creation only quandrangle faces
875 bool StdMeshers_Quadrangle_2D::ComputeQuadPref
877 const TopoDS_Shape& aShape,
878 FaceQuadStruct* quad) throw (SALOME_Exception)
880 Unexpect aCatch(SalomeException);
882 SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
883 const TopoDS_Face& F = TopoDS::Face(aShape);
884 Handle(Geom_Surface) S = BRep_Tool::Surface(F);
885 const TopoDS_Wire& W = BRepTools::OuterWire(F);
887 if(W.Orientation()==TopAbs_FORWARD)
889 //if(WisF) cout<<"W is FORWARD"<<endl;
890 //else cout<<"W is REVERSED"<<endl;
891 bool FisF = (F.Orientation()==TopAbs_FORWARD);
892 if(!FisF) WisF = !WisF;
893 int i,j,geomFaceID = meshDS->ShapeToIndex( F );
895 int nb = quad->side[0]->NbPoints();
896 int nr = quad->side[1]->NbPoints();
897 int nt = quad->side[2]->NbPoints();
898 int nl = quad->side[3]->NbPoints();
904 // it is a base case => not shift quad but me be replacement is need
905 ShiftQuad(quad,0,WisF);
908 // we have to shift quad on 2
909 ShiftQuad(quad,2,WisF);
914 // we have to shift quad on 1
915 ShiftQuad(quad,1,WisF);
918 // we have to shift quad on 3
919 ShiftQuad(quad,3,WisF);
923 nb = quad->side[0]->NbPoints();
924 nr = quad->side[1]->NbPoints();
925 nt = quad->side[2]->NbPoints();
926 nl = quad->side[3]->NbPoints();
929 int nbh = Max(nb,nt);
930 int nbv = Max(nr,nl);
934 // orientation of face and 3 main domain for future faces
940 // left | | | | rigth
956 const vector<UVPtStruct>& uv_eb = quad->side[0]->GetUVPtStruct(true,0 );
957 const vector<UVPtStruct>& uv_er = quad->side[1]->GetUVPtStruct(false,1);
958 const vector<UVPtStruct>& uv_et = quad->side[2]->GetUVPtStruct(true,1 );
959 const vector<UVPtStruct>& uv_el = quad->side[3]->GetUVPtStruct(false,0);
961 // arrays for normalized params
962 //cout<<"Dump B:"<<endl;
963 TColStd_SequenceOfReal npb, npr, npt, npl;
964 for(i=0; i<nb; i++) {
965 npb.Append(uv_eb[i].normParam);
966 //cout<<"i="<<i<<" par="<<uv_eb[i].normParam<<" npar="<<uv_eb[i].normParam;
967 //const SMDS_MeshNode* N = uv_eb[i].node;
968 //cout<<" node("<<N->X()<<","<<N->Y()<<","<<N->Z()<<")"<<endl;
970 for(i=0; i<nr; i++) {
971 npr.Append(uv_er[i].normParam);
973 for(i=0; i<nt; i++) {
974 npt.Append(uv_et[i].normParam);
976 for(i=0; i<nl; i++) {
977 npl.Append(uv_el[i].normParam);
980 // add some params to right and left after the first param
983 double dpr = (npr.Value(2) - npr.Value(1))/(dr+1);
984 for(i=1; i<=dr; i++) {
985 npr.InsertAfter(1,npr.Value(2)-dpr);
989 dpr = (npl.Value(2) - npl.Value(1))/(dl+1);
990 for(i=1; i<=dl; i++) {
991 npl.InsertAfter(1,npl.Value(2)-dpr);
994 //for(i=1; i<=npb.Length(); i++) {
995 // cout<<" "<<npb.Value(i);
999 gp_XY a0( uv_eb.front().u, uv_eb.front().v );
1000 gp_XY a1( uv_eb.back().u, uv_eb.back().v );
1001 gp_XY a2( uv_et.back().u, uv_et.back().v );
1002 gp_XY a3( uv_et.front().u, uv_et.front().v );
1003 //cout<<" a0("<<a0.X()<<","<<a0.Y()<<")"<<" a1("<<a1.X()<<","<<a1.Y()<<")"
1004 // <<" a2("<<a2.X()<<","<<a2.Y()<<")"<<" a3("<<a3.X()<<","<<a3.Y()<<")"<<endl;
1006 int nnn = Min(nr,nl);
1007 // auxilary sequence of XY for creation nodes
1008 // in the bottom part of central domain
1009 // it's length must be == nbv-nnn-1
1010 TColgp_SequenceOfXY UVL;
1011 TColgp_SequenceOfXY UVR;
1013 // step1: create faces for left domain
1014 StdMeshers_Array2OfNode NodesL(1,dl+1,1,nl);
1016 for(j=1; j<=nl; j++)
1017 NodesL.SetValue(1,j,uv_el[j-1].node);
1020 for(i=1; i<=dl; i++)
1021 NodesL.SetValue(i+1,nl,uv_et[i].node);
1022 // create and add needed nodes
1023 TColgp_SequenceOfXY UVtmp;
1024 for(i=1; i<=dl; i++) {
1025 double x0 = npt.Value(i+1);
1028 double y0 = npl.Value(i+1);
1029 double y1 = npr.Value(i+1);
1030 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1031 gp_Pnt P = S->Value(UV.X(),UV.Y());
1032 SMDS_MeshNode * N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1033 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1034 NodesL.SetValue(i+1,1,N);
1035 if(UVL.Length()<nbv-nnn-1) UVL.Append(UV);
1037 for(j=2; j<nl; j++) {
1038 double y0 = npl.Value(dl+j);
1039 double y1 = npr.Value(dl+j);
1040 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1041 gp_Pnt P = S->Value(UV.X(),UV.Y());
1042 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1043 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1044 NodesL.SetValue(i+1,j,N);
1045 if( i==dl ) UVtmp.Append(UV);
1048 for(i=1; i<=UVtmp.Length() && UVL.Length()<nbv-nnn-1; i++) {
1049 UVL.Append(UVtmp.Value(i));
1051 //cout<<"Dump NodesL:"<<endl;
1052 //for(i=1; i<=dl+1; i++) {
1054 // for(j=1; j<=nl; j++) {
1055 // cout<<" ("<<NodesL.Value(i,j)->X()<<","<<NodesL.Value(i,j)->Y()<<","<<NodesL.Value(i,j)->Z()<<")";
1060 for(i=1; i<=dl; i++) {
1061 for(j=1; j<nl; j++) {
1064 myTool->AddFace(NodesL.Value(i,j), NodesL.Value(i+1,j),
1065 NodesL.Value(i+1,j+1), NodesL.Value(i,j+1));
1066 meshDS->SetMeshElementOnShape(F, geomFaceID);
1070 myTool->AddFace(NodesL.Value(i,j), NodesL.Value(i,j+1),
1071 NodesL.Value(i+1,j+1), NodesL.Value(i+1,j));
1072 meshDS->SetMeshElementOnShape(F, geomFaceID);
1078 // fill UVL using c2d
1079 for(i=1; i<npl.Length() && UVL.Length()<nbv-nnn-1; i++) {
1080 UVL.Append( gp_UV ( uv_el[i].u, uv_el[i].v ));
1084 // step2: create faces for right domain
1085 StdMeshers_Array2OfNode NodesR(1,dr+1,1,nr);
1087 for(j=1; j<=nr; j++)
1088 NodesR.SetValue(1,j,uv_er[nr-j].node);
1091 for(i=1; i<=dr; i++)
1092 NodesR.SetValue(i+1,1,uv_et[nt-1-i].node);
1093 // create and add needed nodes
1094 TColgp_SequenceOfXY UVtmp;
1095 for(i=1; i<=dr; i++) {
1096 double x0 = npt.Value(nt-i);
1099 double y0 = npl.Value(i+1);
1100 double y1 = npr.Value(i+1);
1101 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1102 gp_Pnt P = S->Value(UV.X(),UV.Y());
1103 SMDS_MeshNode * N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1104 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1105 NodesR.SetValue(i+1,nr,N);
1106 if(UVR.Length()<nbv-nnn-1) UVR.Append(UV);
1108 for(j=2; j<nr; j++) {
1109 double y0 = npl.Value(nbv-j+1);
1110 double y1 = npr.Value(nbv-j+1);
1111 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1112 gp_Pnt P = S->Value(UV.X(),UV.Y());
1113 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1114 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1115 NodesR.SetValue(i+1,j,N);
1116 if( i==dr ) UVtmp.Prepend(UV);
1119 for(i=1; i<=UVtmp.Length() && UVR.Length()<nbv-nnn-1; i++) {
1120 UVR.Append(UVtmp.Value(i));
1123 for(i=1; i<=dr; i++) {
1124 for(j=1; j<nr; j++) {
1127 myTool->AddFace(NodesR.Value(i,j), NodesR.Value(i+1,j),
1128 NodesR.Value(i+1,j+1), NodesR.Value(i,j+1));
1129 meshDS->SetMeshElementOnShape(F, geomFaceID);
1133 myTool->AddFace(NodesR.Value(i,j), NodesR.Value(i,j+1),
1134 NodesR.Value(i+1,j+1), NodesR.Value(i+1,j));
1135 meshDS->SetMeshElementOnShape(F, geomFaceID);
1141 // fill UVR using c2d
1142 for(i=1; i<npr.Length() && UVR.Length()<nbv-nnn-1; i++) {
1143 UVR.Append( gp_UV( uv_er[i].u, uv_er[i].v ));
1147 // step3: create faces for central domain
1148 StdMeshers_Array2OfNode NodesC(1,nb,1,nbv);
1149 // add first string using NodesL
1150 for(i=1; i<=dl+1; i++)
1151 NodesC.SetValue(1,i,NodesL(i,1));
1152 for(i=2; i<=nl; i++)
1153 NodesC.SetValue(1,dl+i,NodesL(dl+1,i));
1154 // add last string using NodesR
1155 for(i=1; i<=dr+1; i++)
1156 NodesC.SetValue(nb,i,NodesR(i,nr));
1158 NodesC.SetValue(nb,dr+i+1,NodesR(dr+1,nr-i));
1159 // add top nodes (last columns)
1160 for(i=dl+2; i<nbh-dr; i++)
1161 NodesC.SetValue(i-dl,nbv,uv_et[i-1].node);
1162 // add bottom nodes (first columns)
1164 NodesC.SetValue(i,1,uv_eb[i-1].node);
1166 // create and add needed nodes
1167 // add linear layers
1168 for(i=2; i<nb; i++) {
1169 double x0 = npt.Value(dl+i);
1171 for(j=1; j<nnn; j++) {
1172 double y0 = npl.Value(nbv-nnn+j);
1173 double y1 = npr.Value(nbv-nnn+j);
1174 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1175 gp_Pnt P = S->Value(UV.X(),UV.Y());
1176 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1177 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1178 NodesC.SetValue(i,nbv-nnn+j,N);
1181 // add diagonal layers
1182 //cout<<"UVL.Length()="<<UVL.Length()<<" UVR.Length()="<<UVR.Length()<<endl;
1183 //cout<<"Dump UVL:"<<endl;
1184 //for(i=1; i<=UVL.Length(); i++) {
1185 // cout<<" ("<<UVL.Value(i).X()<<","<<UVL.Value(i).Y()<<")";
1188 for(i=1; i<nbv-nnn; i++) {
1189 double du = UVR.Value(i).X() - UVL.Value(i).X();
1190 double dv = UVR.Value(i).Y() - UVL.Value(i).Y();
1191 for(j=2; j<nb; j++) {
1192 double u = UVL.Value(i).X() + du*npb.Value(j);
1193 double v = UVL.Value(i).Y() + dv*npb.Value(j);
1194 gp_Pnt P = S->Value(u,v);
1195 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1196 meshDS->SetNodeOnFace(N, geomFaceID, u, v);
1197 NodesC.SetValue(j,i+1,N);
1201 for(i=1; i<nb; i++) {
1202 for(j=1; j<nbv; j++) {
1205 myTool->AddFace(NodesC.Value(i,j), NodesC.Value(i+1,j),
1206 NodesC.Value(i+1,j+1), NodesC.Value(i,j+1));
1207 meshDS->SetMeshElementOnShape(F, geomFaceID);
1211 myTool->AddFace(NodesC.Value(i,j), NodesC.Value(i,j+1),
1212 NodesC.Value(i+1,j+1), NodesC.Value(i+1,j));
1213 meshDS->SetMeshElementOnShape(F, geomFaceID);