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, false);
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 if ( uv_e0.empty() || uv_e1.empty() || uv_e2.empty() || uv_e3.empty() )
243 return error( COMPERR_BAD_INPUT_MESH );
245 double eps = Precision::Confusion();
247 // Boundary quadrangles
249 if (quad->isEdgeOut[0]) {
252 // |___|___|___|___|___|___|
254 // |___|___|___|___|___|___|
256 // |___|___|___|___|___|___| __ first row of the regular grid
257 // . . . . . . . . . __ down edge nodes
259 // >->->->->->->->->->->->-> -- direction of processing
261 int g = 0; // number of last processed node in the regular grid
263 // number of last node of the down edge to be processed
264 int stop = nbdown - 1;
265 // if right edge is out, we will stop at a node, previous to the last one
266 if (quad->isEdgeOut[1]) stop--;
268 // for each node of the down edge find nearest node
269 // in the first row of the regular grid and link them
270 for (i = 0; i < stop; i++) {
271 const SMDS_MeshNode *a, *b, *c, *d;
273 b = uv_e0[i + 1].node;
274 gp_Pnt pb (b->X(), b->Y(), b->Z());
276 // find node c in the regular grid, which will be linked with node b
279 // right bound reached, link with the rightmost node
281 c = quad->uv_grid[nbhoriz + iup].node;
284 // find in the grid node c, nearest to the b
285 double mind = RealLast();
286 for (int k = g; k <= iup; k++) {
288 const SMDS_MeshNode *nk;
289 if (k < ilow) // this can be, if left edge is out
290 nk = uv_e3[1].node; // get node from the left edge
292 nk = quad->uv_grid[nbhoriz + k].node; // get one of middle nodes
294 gp_Pnt pnk (nk->X(), nk->Y(), nk->Z());
295 double dist = pb.Distance(pnk);
296 if (dist < mind - eps) {
306 if (near == g) { // make triangle
307 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c);
308 SMDS_MeshFace* face = myTool->AddFace(a, b, c);
309 meshDS->SetMeshElementOnShape(face, geomFaceID);
311 else { // make quadrangle
315 d = quad->uv_grid[nbhoriz + near - 1].node;
316 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c, d);
317 SMDS_MeshFace* face = myTool->AddFace(a, b, c, d);
318 meshDS->SetMeshElementOnShape(face, geomFaceID);
320 // if node d is not at position g - make additional triangles
322 for (int k = near - 1; k > g; k--) {
323 c = quad->uv_grid[nbhoriz + k].node;
327 d = quad->uv_grid[nbhoriz + k - 1].node;
328 //SMDS_MeshFace* face = meshDS->AddFace(a, c, d);
329 SMDS_MeshFace* face = myTool->AddFace(a, c, d);
330 meshDS->SetMeshElementOnShape(face, geomFaceID);
337 if (quad->isEdgeOut[2]) {
340 // <-<-<-<-<-<-<-<-<-<-<-<-< -- direction of processing
342 // . . . . . . . . . __ up edge nodes
343 // ___ ___ ___ ___ ___ ___ __ first row of the regular grid
345 // |___|___|___|___|___|___|
347 // |___|___|___|___|___|___|
350 int g = nbhoriz - 1; // last processed node in the regular grid
353 // if left edge is out, we will stop at a second node
354 if (quad->isEdgeOut[3]) stop++;
356 // for each node of the up edge find nearest node
357 // in the first row of the regular grid and link them
358 for (i = nbup - 1; i > stop; i--) {
359 const SMDS_MeshNode *a, *b, *c, *d;
361 b = uv_e2[i - 1].node;
362 gp_Pnt pb (b->X(), b->Y(), b->Z());
364 // find node c in the grid, which will be linked with node b
366 if (i == stop + 1) { // left bound reached, link with the leftmost node
367 c = quad->uv_grid[nbhoriz*(nbvertic - 2) + ilow].node;
370 // find node c in the grid, nearest to the b
371 double mind = RealLast();
372 for (int k = g; k >= ilow; k--) {
373 const SMDS_MeshNode *nk;
375 nk = uv_e1[nbright - 2].node;
377 nk = quad->uv_grid[nbhoriz*(nbvertic - 2) + k].node;
378 gp_Pnt pnk (nk->X(), nk->Y(), nk->Z());
379 double dist = pb.Distance(pnk);
380 if (dist < mind - eps) {
390 if (near == g) { // make triangle
391 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c);
392 SMDS_MeshFace* face = myTool->AddFace(a, b, c);
393 meshDS->SetMeshElementOnShape(face, geomFaceID);
395 else { // make quadrangle
397 d = uv_e1[nbright - 2].node;
399 d = quad->uv_grid[nbhoriz*(nbvertic - 2) + near + 1].node;
400 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c, d);
401 SMDS_MeshFace* face = myTool->AddFace(a, b, c, d);
402 meshDS->SetMeshElementOnShape(face, geomFaceID);
404 if (near + 1 < g) { // if d not is at g - make additional triangles
405 for (int k = near + 1; k < g; k++) {
406 c = quad->uv_grid[nbhoriz*(nbvertic - 2) + k].node;
408 d = uv_e1[nbright - 2].node;
410 d = quad->uv_grid[nbhoriz*(nbvertic - 2) + k + 1].node;
411 //SMDS_MeshFace* face = meshDS->AddFace(a, c, d);
412 SMDS_MeshFace* face = myTool->AddFace(a, c, d);
413 meshDS->SetMeshElementOnShape(face, geomFaceID);
422 // right or left boundary quadrangles
423 if (quad->isEdgeOut[1]) {
424 // MESSAGE("right edge is out");
425 int g = 0; // last processed node in the grid
426 int stop = nbright - 1;
427 if (quad->isEdgeOut[2]) stop--;
428 for (i = 0; i < stop; i++) {
429 const SMDS_MeshNode *a, *b, *c, *d;
431 b = uv_e1[i + 1].node;
432 gp_Pnt pb (b->X(), b->Y(), b->Z());
434 // find node c in the grid, nearest to the b
436 if (i == stop - 1) { // up bondary reached
437 c = quad->uv_grid[nbhoriz*(jup + 1) - 2].node;
440 double mind = RealLast();
441 for (int k = g; k <= jup; k++) {
442 const SMDS_MeshNode *nk;
444 nk = uv_e0[nbdown - 2].node;
446 nk = quad->uv_grid[nbhoriz*(k + 1) - 2].node;
447 gp_Pnt pnk (nk->X(), nk->Y(), nk->Z());
448 double dist = pb.Distance(pnk);
449 if (dist < mind - eps) {
459 if (near == g) { // make triangle
460 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c);
461 SMDS_MeshFace* face = myTool->AddFace(a, b, c);
462 meshDS->SetMeshElementOnShape(face, geomFaceID);
464 else { // make quadrangle
466 d = uv_e0[nbdown - 2].node;
468 d = quad->uv_grid[nbhoriz*near - 2].node;
469 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c, d);
470 SMDS_MeshFace* face = myTool->AddFace(a, b, c, d);
471 meshDS->SetMeshElementOnShape(face, geomFaceID);
473 if (near - 1 > g) { // if d not is at g - make additional triangles
474 for (int k = near - 1; k > g; k--) {
475 c = quad->uv_grid[nbhoriz*(k + 1) - 2].node;
477 d = uv_e0[nbdown - 2].node;
479 d = quad->uv_grid[nbhoriz*k - 2].node;
480 //SMDS_MeshFace* face = meshDS->AddFace(a, c, d);
481 SMDS_MeshFace* face = myTool->AddFace(a, c, d);
482 meshDS->SetMeshElementOnShape(face, geomFaceID);
489 if (quad->isEdgeOut[3]) {
490 // MESSAGE("left edge is out");
491 int g = nbvertic - 1; // last processed node in the grid
493 if (quad->isEdgeOut[0]) stop++;
494 for (i = nbleft - 1; i > stop; i--) {
495 const SMDS_MeshNode *a, *b, *c, *d;
497 b = uv_e3[i - 1].node;
498 gp_Pnt pb (b->X(), b->Y(), b->Z());
500 // find node c in the grid, nearest to the b
502 if (i == stop + 1) { // down bondary reached
503 c = quad->uv_grid[nbhoriz*jlow + 1].node;
506 double mind = RealLast();
507 for (int k = g; k >= jlow; k--) {
508 const SMDS_MeshNode *nk;
512 nk = quad->uv_grid[nbhoriz*k + 1].node;
513 gp_Pnt pnk (nk->X(), nk->Y(), nk->Z());
514 double dist = pb.Distance(pnk);
515 if (dist < mind - eps) {
525 if (near == g) { // make triangle
526 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c);
527 SMDS_MeshFace* face = myTool->AddFace(a, b, c);
528 meshDS->SetMeshElementOnShape(face, geomFaceID);
530 else { // make quadrangle
534 d = quad->uv_grid[nbhoriz*(near + 1) + 1].node;
535 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c, d);
536 SMDS_MeshFace* face = myTool->AddFace(a, b, c, d);
537 meshDS->SetMeshElementOnShape(face, geomFaceID);
539 if (near + 1 < g) { // if d not is at g - make additional triangles
540 for (int k = near + 1; k < g; k++) {
541 c = quad->uv_grid[nbhoriz*k + 1].node;
545 d = quad->uv_grid[nbhoriz*(k + 1) + 1].node;
546 //SMDS_MeshFace* face = meshDS->AddFace(a, c, d);
547 SMDS_MeshFace* face = myTool->AddFace(a, c, d);
548 meshDS->SetMeshElementOnShape(face, geomFaceID);
561 //=============================================================================
565 //=============================================================================
567 FaceQuadStruct* StdMeshers_Quadrangle_2D::CheckNbEdges(SMESH_Mesh & aMesh,
568 const TopoDS_Shape & aShape)
569 //throw(SALOME_Exception)
571 const TopoDS_Face & F = TopoDS::Face(aShape);
572 const bool ignoreMediumNodes = _quadraticMesh;
574 // verify 1 wire only, with 4 edges
576 list< TopoDS_Edge > edges;
577 list< int > nbEdgesInWire;
578 int nbWire = SMESH_Block::GetOrderedEdges (F, V, edges, nbEdgesInWire);
580 error(COMPERR_BAD_SHAPE, TComm("Wrong number of wires: ") << nbWire);
583 FaceQuadStruct* quad = new FaceQuadStruct;
585 quad->side.reserve(nbEdgesInWire.front());
588 list< TopoDS_Edge >::iterator edgeIt = edges.begin();
589 if ( nbEdgesInWire.front() == 4 ) { // exactly 4 edges
590 for ( ; edgeIt != edges.end(); ++edgeIt, nbSides++ )
591 quad->side.push_back( new StdMeshers_FaceSide(F, *edgeIt, &aMesh,
592 nbSides<TOP_SIDE, ignoreMediumNodes));
594 else if ( nbEdgesInWire.front() > 4 ) { // more than 4 edges - try to unite some
595 list< TopoDS_Edge > sideEdges;
596 while ( !edges.empty()) {
598 sideEdges.splice( sideEdges.end(), edges, edges.begin()); // edges.front() -> sideEdges.end()
599 bool sameSide = true;
600 while ( !edges.empty() && sameSide ) {
601 sameSide = SMESH_Algo::IsContinuous( sideEdges.back(), edges.front() );
603 sideEdges.splice( sideEdges.end(), edges, edges.begin());
605 if ( nbSides == 0 ) { // go backward from the first edge
607 while ( !edges.empty() && sameSide ) {
608 sameSide = SMESH_Algo::IsContinuous( sideEdges.front(), edges.back() );
610 sideEdges.splice( sideEdges.begin(), edges, --edges.end());
613 quad->side.push_back( new StdMeshers_FaceSide(F, sideEdges, &aMesh,
614 nbSides<TOP_SIDE, ignoreMediumNodes));
620 cout << endl << "StdMeshers_Quadrangle_2D. Edge IDs of " << nbSides << " sides:";
621 for ( int i = 0; i < nbSides; ++i ) {
623 for ( int e = 0; e < quad->side[i]->NbEdges(); ++e )
624 cout << myTool->GetMeshDS()->ShapeToIndex( quad->side[i]->Edge( e )) << " ";
630 nbSides = nbEdgesInWire.front();
631 error(COMPERR_BAD_SHAPE, TComm("Face must have 4 sides but not ") << nbSides);
639 //=============================================================================
643 //=============================================================================
645 FaceQuadStruct *StdMeshers_Quadrangle_2D::CheckAnd2Dcompute
647 const TopoDS_Shape & aShape,
648 const bool CreateQuadratic) //throw(SALOME_Exception)
650 _quadraticMesh = CreateQuadratic;
652 FaceQuadStruct *quad = CheckNbEdges(aMesh, aShape);
656 // set normalized grid on unit square in parametric domain
657 bool stat = SetNormalizedGrid(aMesh, aShape, quad);
667 //=============================================================================
671 //=============================================================================
673 faceQuadStruct::~faceQuadStruct()
675 for (int i = 0; i < side.size(); i++) {
676 if (side[i]) delete side[i];
678 if (uv_grid) delete [] uv_grid;
682 inline const vector<UVPtStruct>& GetUVPtStructIn(FaceQuadStruct* quad, int i, int nbSeg)
684 bool isXConst = ( i == BOTTOM_SIDE || i == TOP_SIDE );
685 double constValue = ( i == BOTTOM_SIDE || i == LEFT_SIDE ) ? 0 : 1;
688 quad->side[i]->SimulateUVPtStruct(nbSeg,isXConst,constValue) :
689 quad->side[i]->GetUVPtStruct(isXConst,constValue);
693 //=============================================================================
697 //=============================================================================
699 bool StdMeshers_Quadrangle_2D::SetNormalizedGrid (SMESH_Mesh & aMesh,
700 const TopoDS_Shape& aShape,
701 FaceQuadStruct* & quad) //throw (SALOME_Exception)
703 // Algorithme décrit dans "Génération automatique de maillages"
704 // P.L. GEORGE, MASSON, § 6.4.1 p. 84-85
705 // traitement dans le domaine paramétrique 2d u,v
706 // transport - projection sur le carré unité
708 // MESSAGE("StdMeshers_Quadrangle_2D::SetNormalizedGrid");
709 // const TopoDS_Face& F = TopoDS::Face(aShape);
711 // 1 --- find orientation of the 4 edges, by test on extrema
714 // |<----north-2-------^ a3 -------------> a2
716 // west-3 east-1 =right | |
720 // v----south-0--------> a0 -------------> a1
725 // 3 --- 2D normalized values on unit square [0..1][0..1]
727 int nbhoriz = Min(quad->side[0]->NbPoints(), quad->side[2]->NbPoints());
728 int nbvertic = Min(quad->side[1]->NbPoints(), quad->side[3]->NbPoints());
730 quad->isEdgeOut[0] = (quad->side[0]->NbPoints() > quad->side[2]->NbPoints());
731 quad->isEdgeOut[1] = (quad->side[1]->NbPoints() > quad->side[3]->NbPoints());
732 quad->isEdgeOut[2] = (quad->side[2]->NbPoints() > quad->side[0]->NbPoints());
733 quad->isEdgeOut[3] = (quad->side[3]->NbPoints() > quad->side[1]->NbPoints());
735 UVPtStruct *uv_grid = quad->uv_grid = new UVPtStruct[nbvertic * nbhoriz];
737 const vector<UVPtStruct>& uv_e0 = GetUVPtStructIn( quad, 0, nbhoriz - 1 );
738 const vector<UVPtStruct>& uv_e1 = GetUVPtStructIn( quad, 1, nbvertic - 1 );
739 const vector<UVPtStruct>& uv_e2 = GetUVPtStructIn( quad, 2, nbhoriz - 1 );
740 const vector<UVPtStruct>& uv_e3 = GetUVPtStructIn( quad, 3, nbvertic - 1 );
742 if ( uv_e0.empty() || uv_e1.empty() || uv_e2.empty() || uv_e3.empty() )
743 //return error( "Can't find nodes on sides");
744 return error( COMPERR_BAD_INPUT_MESH );
746 // nodes Id on "in" edges
747 if (! quad->isEdgeOut[0]) {
749 for (int i = 0; i < nbhoriz; i++) { // down
750 int ij = j * nbhoriz + i;
751 uv_grid[ij].node = uv_e0[i].node;
754 if (! quad->isEdgeOut[1]) {
756 for (int j = 0; j < nbvertic; j++) { // right
757 int ij = j * nbhoriz + i;
758 uv_grid[ij].node = uv_e1[j].node;
761 if (! quad->isEdgeOut[2]) {
762 int j = nbvertic - 1;
763 for (int i = 0; i < nbhoriz; i++) { // up
764 int ij = j * nbhoriz + i;
765 uv_grid[ij].node = uv_e2[i].node;
768 if (! quad->isEdgeOut[3]) {
770 for (int j = 0; j < nbvertic; j++) { // left
771 int ij = j * nbhoriz + i;
772 uv_grid[ij].node = uv_e3[j].node;
776 // normalized 2d values on grid
777 for (int i = 0; i < nbhoriz; i++)
779 for (int j = 0; j < nbvertic; j++)
781 int ij = j * nbhoriz + i;
782 // --- droite i cste : x = x0 + y(x1-x0)
783 double x0 = uv_e0[i].normParam; // bas - sud
784 double x1 = uv_e2[i].normParam; // haut - nord
785 // --- droite j cste : y = y0 + x(y1-y0)
786 double y0 = uv_e3[j].normParam; // gauche-ouest
787 double y1 = uv_e1[j].normParam; // droite - est
788 // --- intersection : x=x0+(y0+x(y1-y0))(x1-x0)
789 double x = (x0 + y0 * (x1 - x0)) / (1 - (y1 - y0) * (x1 - x0));
790 double y = y0 + x * (y1 - y0);
793 //MESSAGE("-xy-01 "<<x0<<" "<<x1<<" "<<y0<<" "<<y1);
794 //MESSAGE("-xy-norm "<<i<<" "<<j<<" "<<x<<" "<<y);
798 // 4 --- projection on 2d domain (u,v)
799 gp_UV a0( uv_e0.front().u, uv_e0.front().v );
800 gp_UV a1( uv_e0.back().u, uv_e0.back().v );
801 gp_UV a2( uv_e2.back().u, uv_e2.back().v );
802 gp_UV a3( uv_e2.front().u, uv_e2.front().v );
804 for (int i = 0; i < nbhoriz; i++)
806 for (int j = 0; j < nbvertic; j++)
808 int ij = j * nbhoriz + i;
809 double x = uv_grid[ij].x;
810 double y = uv_grid[ij].y;
811 double param_0 = uv_e0[0].normParam + x * (uv_e0.back().normParam - uv_e0[0].normParam); // sud
812 double param_2 = uv_e2[0].normParam + x * (uv_e2.back().normParam - uv_e2[0].normParam); // nord
813 double param_1 = uv_e1[0].normParam + y * (uv_e1.back().normParam - uv_e1[0].normParam); // est
814 double param_3 = uv_e3[0].normParam + y * (uv_e3.back().normParam - uv_e3[0].normParam); // ouest
816 //MESSAGE("params "<<param_0<<" "<<param_1<<" "<<param_2<<" "<<param_3);
817 gp_UV p0 = quad->side[0]->Value2d(param_0).XY();
818 gp_UV p1 = quad->side[1]->Value2d(param_1).XY();
819 gp_UV p2 = quad->side[2]->Value2d(param_2).XY();
820 gp_UV p3 = quad->side[3]->Value2d(param_3).XY();
822 gp_UV uv = (1 - y) * p0 + x * p1 + y * p2 + (1 - x) * p3;
823 uv -= (1 - x) * (1 - y) * a0 + x * (1 - y) * a1 + x * y * a2 + (1 - x) * y * a3;
825 uv_grid[ij].u = uv.X();
826 uv_grid[ij].v = uv.Y();
832 //=======================================================================
833 //function : ShiftQuad
834 //purpose : auxilary function for ComputeQuadPref
835 //=======================================================================
837 static void ShiftQuad(FaceQuadStruct* quad, const int num, bool)
839 StdMeshers_FaceSide* side[4] = { quad->side[0], quad->side[1], quad->side[2], quad->side[3] };
840 for (int i = BOTTOM_SIDE; i < NB_SIDES; ++i ) {
841 int id = ( i + num ) % NB_SIDES;
842 bool wasForward = ( i < TOP_SIDE );
843 bool newForward = ( id < TOP_SIDE );
844 if ( wasForward != newForward )
845 side[ i ]->Reverse();
846 quad->side[ id ] = side[ i ];
850 //=======================================================================
852 //purpose : auxilary function for ComputeQuadPref
853 //=======================================================================
855 static gp_UV CalcUV(double x0, double x1, double y0, double y1,
856 FaceQuadStruct* quad,
857 const gp_UV& a0, const gp_UV& a1,
858 const gp_UV& a2, const gp_UV& a3)
860 const vector<UVPtStruct>& uv_eb = quad->side[0]->GetUVPtStruct(true,0 );
861 const vector<UVPtStruct>& uv_er = quad->side[1]->GetUVPtStruct(false,1);
862 const vector<UVPtStruct>& uv_et = quad->side[2]->GetUVPtStruct(true,1 );
863 const vector<UVPtStruct>& uv_el = quad->side[3]->GetUVPtStruct(false,0);
865 double x = (x0 + y0 * (x1 - x0)) / (1 - (y1 - y0) * (x1 - x0));
866 double y = y0 + x * (y1 - y0);
868 double param_b = uv_eb[0].normParam + x * (uv_eb.back().normParam - uv_eb[0].normParam);
869 double param_t = uv_et[0].normParam + x * (uv_et.back().normParam - uv_et[0].normParam);
870 double param_r = uv_er[0].normParam + y * (uv_er.back().normParam - uv_er[0].normParam);
871 double param_l = uv_el[0].normParam + y * (uv_el.back().normParam - uv_el[0].normParam);
873 gp_UV p0 = quad->side[BOTTOM_SIDE]->Value2d(param_b).XY();
874 gp_UV p1 = quad->side[RIGHT_SIDE ]->Value2d(param_r).XY();
875 gp_UV p2 = quad->side[TOP_SIDE ]->Value2d(param_t).XY();
876 gp_UV p3 = quad->side[LEFT_SIDE ]->Value2d(param_l).XY();
878 gp_UV uv = p0 * (1 - y) + p1 * x + p2 * y + p3 * (1 - x);
880 uv -= (1 - x) * (1 - y) * a0 + x * (1 - y) * a1 + x * y * a2 + (1 - x) * y * a3;
885 //=======================================================================
887 * Create only quandrangle faces
889 //=======================================================================
891 bool StdMeshers_Quadrangle_2D::ComputeQuadPref (SMESH_Mesh & aMesh,
892 const TopoDS_Shape& aShape,
893 FaceQuadStruct* quad)
895 SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
896 const TopoDS_Face& F = TopoDS::Face(aShape);
897 Handle(Geom_Surface) S = BRep_Tool::Surface(F);
898 // const TopoDS_Wire& W = BRepTools::OuterWire(F);
900 // if(W.Orientation()==TopAbs_FORWARD)
902 //if(WisF) cout<<"W is FORWARD"<<endl;
903 //else cout<<"W is REVERSED"<<endl;
904 // bool FisF = (F.Orientation()==TopAbs_FORWARD);
905 // if(!FisF) WisF = !WisF;
907 int i,j,geomFaceID = meshDS->ShapeToIndex( F );
909 int nb = quad->side[0]->NbPoints();
910 int nr = quad->side[1]->NbPoints();
911 int nt = quad->side[2]->NbPoints();
912 int nl = quad->side[3]->NbPoints();
918 // it is a base case => not shift quad but me be replacement is need
919 ShiftQuad(quad,0,WisF);
922 // we have to shift quad on 2
923 ShiftQuad(quad,2,WisF);
928 // we have to shift quad on 1
929 ShiftQuad(quad,1,WisF);
932 // we have to shift quad on 3
933 ShiftQuad(quad,3,WisF);
937 nb = quad->side[0]->NbPoints();
938 nr = quad->side[1]->NbPoints();
939 nt = quad->side[2]->NbPoints();
940 nl = quad->side[3]->NbPoints();
943 int nbh = Max(nb,nt);
944 int nbv = Max(nr,nl);
948 // orientation of face and 3 main domain for future faces
954 // left | | | | rigth
970 const vector<UVPtStruct>& uv_eb = quad->side[0]->GetUVPtStruct(true,0 );
971 const vector<UVPtStruct>& uv_er = quad->side[1]->GetUVPtStruct(false,1);
972 const vector<UVPtStruct>& uv_et = quad->side[2]->GetUVPtStruct(true,1 );
973 const vector<UVPtStruct>& uv_el = quad->side[3]->GetUVPtStruct(false,0);
975 // arrays for normalized params
976 //cout<<"Dump B:"<<endl;
977 TColStd_SequenceOfReal npb, npr, npt, npl;
978 for(i=0; i<nb; i++) {
979 npb.Append(uv_eb[i].normParam);
980 //cout<<"i="<<i<<" par="<<uv_eb[i].normParam<<" npar="<<uv_eb[i].normParam;
981 //const SMDS_MeshNode* N = uv_eb[i].node;
982 //cout<<" node("<<N->X()<<","<<N->Y()<<","<<N->Z()<<")"<<endl;
984 for(i=0; i<nr; i++) {
985 npr.Append(uv_er[i].normParam);
987 for(i=0; i<nt; i++) {
988 npt.Append(uv_et[i].normParam);
990 for(i=0; i<nl; i++) {
991 npl.Append(uv_el[i].normParam);
994 // add some params to right and left after the first param
997 double dpr = (npr.Value(2) - npr.Value(1))/(dr+1);
998 for(i=1; i<=dr; i++) {
999 npr.InsertAfter(1,npr.Value(2)-dpr);
1003 dpr = (npl.Value(2) - npl.Value(1))/(dl+1);
1004 for(i=1; i<=dl; i++) {
1005 npl.InsertAfter(1,npl.Value(2)-dpr);
1008 //for(i=1; i<=npb.Length(); i++) {
1009 // cout<<" "<<npb.Value(i);
1013 gp_XY a0( uv_eb.front().u, uv_eb.front().v );
1014 gp_XY a1( uv_eb.back().u, uv_eb.back().v );
1015 gp_XY a2( uv_et.back().u, uv_et.back().v );
1016 gp_XY a3( uv_et.front().u, uv_et.front().v );
1017 //cout<<" a0("<<a0.X()<<","<<a0.Y()<<")"<<" a1("<<a1.X()<<","<<a1.Y()<<")"
1018 // <<" a2("<<a2.X()<<","<<a2.Y()<<")"<<" a3("<<a3.X()<<","<<a3.Y()<<")"<<endl;
1020 int nnn = Min(nr,nl);
1021 // auxilary sequence of XY for creation nodes
1022 // in the bottom part of central domain
1023 // it's length must be == nbv-nnn-1
1024 TColgp_SequenceOfXY UVL;
1025 TColgp_SequenceOfXY UVR;
1027 // step1: create faces for left domain
1028 StdMeshers_Array2OfNode NodesL(1,dl+1,1,nl);
1030 for(j=1; j<=nl; j++)
1031 NodesL.SetValue(1,j,uv_el[j-1].node);
1034 for(i=1; i<=dl; i++)
1035 NodesL.SetValue(i+1,nl,uv_et[i].node);
1036 // create and add needed nodes
1037 TColgp_SequenceOfXY UVtmp;
1038 for(i=1; i<=dl; i++) {
1039 double x0 = npt.Value(i+1);
1042 double y0 = npl.Value(i+1);
1043 double y1 = npr.Value(i+1);
1044 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1045 gp_Pnt P = S->Value(UV.X(),UV.Y());
1046 SMDS_MeshNode * N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1047 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1048 NodesL.SetValue(i+1,1,N);
1049 if(UVL.Length()<nbv-nnn-1) UVL.Append(UV);
1051 for(j=2; j<nl; j++) {
1052 double y0 = npl.Value(dl+j);
1053 double y1 = npr.Value(dl+j);
1054 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1055 gp_Pnt P = S->Value(UV.X(),UV.Y());
1056 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1057 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1058 NodesL.SetValue(i+1,j,N);
1059 if( i==dl ) UVtmp.Append(UV);
1062 for(i=1; i<=UVtmp.Length() && UVL.Length()<nbv-nnn-1; i++) {
1063 UVL.Append(UVtmp.Value(i));
1065 //cout<<"Dump NodesL:"<<endl;
1066 //for(i=1; i<=dl+1; i++) {
1068 // for(j=1; j<=nl; j++) {
1069 // cout<<" ("<<NodesL.Value(i,j)->X()<<","<<NodesL.Value(i,j)->Y()<<","<<NodesL.Value(i,j)->Z()<<")";
1074 for(i=1; i<=dl; i++) {
1075 for(j=1; j<nl; j++) {
1078 myTool->AddFace(NodesL.Value(i,j), NodesL.Value(i+1,j),
1079 NodesL.Value(i+1,j+1), NodesL.Value(i,j+1));
1080 meshDS->SetMeshElementOnShape(F, geomFaceID);
1084 myTool->AddFace(NodesL.Value(i,j), NodesL.Value(i,j+1),
1085 NodesL.Value(i+1,j+1), NodesL.Value(i+1,j));
1086 meshDS->SetMeshElementOnShape(F, geomFaceID);
1092 // fill UVL using c2d
1093 for(i=1; i<npl.Length() && UVL.Length()<nbv-nnn-1; i++) {
1094 UVL.Append( gp_UV ( uv_el[i].u, uv_el[i].v ));
1098 // step2: create faces for right domain
1099 StdMeshers_Array2OfNode NodesR(1,dr+1,1,nr);
1101 for(j=1; j<=nr; j++)
1102 NodesR.SetValue(1,j,uv_er[nr-j].node);
1105 for(i=1; i<=dr; i++)
1106 NodesR.SetValue(i+1,1,uv_et[nt-1-i].node);
1107 // create and add needed nodes
1108 TColgp_SequenceOfXY UVtmp;
1109 for(i=1; i<=dr; i++) {
1110 double x0 = npt.Value(nt-i);
1113 double y0 = npl.Value(i+1);
1114 double y1 = npr.Value(i+1);
1115 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1116 gp_Pnt P = S->Value(UV.X(),UV.Y());
1117 SMDS_MeshNode * N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1118 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1119 NodesR.SetValue(i+1,nr,N);
1120 if(UVR.Length()<nbv-nnn-1) UVR.Append(UV);
1122 for(j=2; j<nr; j++) {
1123 double y0 = npl.Value(nbv-j+1);
1124 double y1 = npr.Value(nbv-j+1);
1125 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1126 gp_Pnt P = S->Value(UV.X(),UV.Y());
1127 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1128 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1129 NodesR.SetValue(i+1,j,N);
1130 if( i==dr ) UVtmp.Prepend(UV);
1133 for(i=1; i<=UVtmp.Length() && UVR.Length()<nbv-nnn-1; i++) {
1134 UVR.Append(UVtmp.Value(i));
1137 for(i=1; i<=dr; i++) {
1138 for(j=1; j<nr; j++) {
1141 myTool->AddFace(NodesR.Value(i,j), NodesR.Value(i+1,j),
1142 NodesR.Value(i+1,j+1), NodesR.Value(i,j+1));
1143 meshDS->SetMeshElementOnShape(F, geomFaceID);
1147 myTool->AddFace(NodesR.Value(i,j), NodesR.Value(i,j+1),
1148 NodesR.Value(i+1,j+1), NodesR.Value(i+1,j));
1149 meshDS->SetMeshElementOnShape(F, geomFaceID);
1155 // fill UVR using c2d
1156 for(i=1; i<npr.Length() && UVR.Length()<nbv-nnn-1; i++) {
1157 UVR.Append( gp_UV( uv_er[i].u, uv_er[i].v ));
1161 // step3: create faces for central domain
1162 StdMeshers_Array2OfNode NodesC(1,nb,1,nbv);
1163 // add first string using NodesL
1164 for(i=1; i<=dl+1; i++)
1165 NodesC.SetValue(1,i,NodesL(i,1));
1166 for(i=2; i<=nl; i++)
1167 NodesC.SetValue(1,dl+i,NodesL(dl+1,i));
1168 // add last string using NodesR
1169 for(i=1; i<=dr+1; i++)
1170 NodesC.SetValue(nb,i,NodesR(i,nr));
1172 NodesC.SetValue(nb,dr+i+1,NodesR(dr+1,nr-i));
1173 // add top nodes (last columns)
1174 for(i=dl+2; i<nbh-dr; i++)
1175 NodesC.SetValue(i-dl,nbv,uv_et[i-1].node);
1176 // add bottom nodes (first columns)
1178 NodesC.SetValue(i,1,uv_eb[i-1].node);
1180 // create and add needed nodes
1181 // add linear layers
1182 for(i=2; i<nb; i++) {
1183 double x0 = npt.Value(dl+i);
1185 for(j=1; j<nnn; j++) {
1186 double y0 = npl.Value(nbv-nnn+j);
1187 double y1 = npr.Value(nbv-nnn+j);
1188 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1189 gp_Pnt P = S->Value(UV.X(),UV.Y());
1190 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1191 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1192 NodesC.SetValue(i,nbv-nnn+j,N);
1195 // add diagonal layers
1196 //cout<<"UVL.Length()="<<UVL.Length()<<" UVR.Length()="<<UVR.Length()<<endl;
1197 //cout<<"Dump UVL:"<<endl;
1198 //for(i=1; i<=UVL.Length(); i++) {
1199 // cout<<" ("<<UVL.Value(i).X()<<","<<UVL.Value(i).Y()<<")";
1202 for(i=1; i<nbv-nnn; i++) {
1203 double du = UVR.Value(i).X() - UVL.Value(i).X();
1204 double dv = UVR.Value(i).Y() - UVL.Value(i).Y();
1205 for(j=2; j<nb; j++) {
1206 double u = UVL.Value(i).X() + du*npb.Value(j);
1207 double v = UVL.Value(i).Y() + dv*npb.Value(j);
1208 gp_Pnt P = S->Value(u,v);
1209 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1210 meshDS->SetNodeOnFace(N, geomFaceID, u, v);
1211 NodesC.SetValue(j,i+1,N);
1215 for(i=1; i<nb; i++) {
1216 for(j=1; j<nbv; j++) {
1219 myTool->AddFace(NodesC.Value(i,j), NodesC.Value(i+1,j),
1220 NodesC.Value(i+1,j+1), NodesC.Value(i,j+1));
1221 meshDS->SetMeshElementOnShape(F, geomFaceID);
1225 myTool->AddFace(NodesC.Value(i,j), NodesC.Value(i,j+1),
1226 NodesC.Value(i+1,j+1), NodesC.Value(i+1,j));
1227 meshDS->SetMeshElementOnShape(F, geomFaceID);