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 Unexpect aCatch(SalomeException);
139 SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
140 aMesh.GetSubMesh(aShape);
142 SMESH_MesherHelper helper(aMesh);
145 _quadraticMesh = myTool->IsQuadraticSubMesh(aShape);
147 FaceQuadStruct *quad = CheckNbEdges( aMesh, aShape );
148 std::auto_ptr<FaceQuadStruct> quadDeleter( quad ); // to delete quad at exit from Compute()
152 if(myQuadranglePreference) {
153 int n1 = quad->side[0]->NbPoints();
154 int n2 = quad->side[1]->NbPoints();
155 int n3 = quad->side[2]->NbPoints();
156 int n4 = quad->side[3]->NbPoints();
157 int nfull = n1+n2+n3+n4;
160 if( nfull==ntmp && ( (n1!=n3) || (n2!=n4) ) ) {
161 // special path for using only quandrangle faces
162 bool ok = ComputeQuadPref(aMesh, aShape, quad);
167 // set normalized grid on unit square in parametric domain
169 if (!SetNormalizedGrid(aMesh, aShape, quad))
172 // --- compute 3D values on points, store points & quadrangles
174 int nbdown = quad->side[0]->NbPoints();
175 int nbup = quad->side[2]->NbPoints();
177 int nbright = quad->side[1]->NbPoints();
178 int nbleft = quad->side[3]->NbPoints();
180 int nbhoriz = Min(nbdown, nbup);
181 int nbvertic = Min(nbright, nbleft);
183 const TopoDS_Face& F = TopoDS::Face(aShape);
184 Handle(Geom_Surface) S = BRep_Tool::Surface(F);
186 // internal mesh nodes
187 int i, j, geomFaceID = meshDS->ShapeToIndex( F );
188 for (i = 1; i < nbhoriz - 1; i++) {
189 for (j = 1; j < nbvertic - 1; j++) {
190 int ij = j * nbhoriz + i;
191 double u = quad->uv_grid[ij].u;
192 double v = quad->uv_grid[ij].v;
193 gp_Pnt P = S->Value(u, v);
194 SMDS_MeshNode * node = meshDS->AddNode(P.X(), P.Y(), P.Z());
195 meshDS->SetNodeOnFace(node, geomFaceID, u, v);
196 quad->uv_grid[ij].node = node;
203 // --.--.--.--.--.-- nbvertic
209 // ---.----.----.--- 0
210 // 0 > > > > > > > > nbhoriz
216 int iup = nbhoriz - 1;
217 if (quad->isEdgeOut[3]) { ilow++; } else { if (quad->isEdgeOut[1]) iup--; }
220 int jup = nbvertic - 1;
221 if (quad->isEdgeOut[0]) { jlow++; } else { if (quad->isEdgeOut[2]) jup--; }
223 // regular quadrangles
224 for (i = ilow; i < iup; i++) {
225 for (j = jlow; j < jup; j++) {
226 const SMDS_MeshNode *a, *b, *c, *d;
227 a = quad->uv_grid[j * nbhoriz + i].node;
228 b = quad->uv_grid[j * nbhoriz + i + 1].node;
229 c = quad->uv_grid[(j + 1) * nbhoriz + i + 1].node;
230 d = quad->uv_grid[(j + 1) * nbhoriz + i].node;
231 SMDS_MeshFace* face = myTool->AddFace(a, b, c, d);
232 meshDS->SetMeshElementOnShape(face, geomFaceID);
236 const vector<UVPtStruct>& uv_e0 = quad->side[0]->GetUVPtStruct(true,0 );
237 const vector<UVPtStruct>& uv_e1 = quad->side[1]->GetUVPtStruct(false,1);
238 const vector<UVPtStruct>& uv_e2 = quad->side[2]->GetUVPtStruct(true,1 );
239 const vector<UVPtStruct>& uv_e3 = quad->side[3]->GetUVPtStruct(false,0);
241 double eps = Precision::Confusion();
243 // Boundary quadrangles
245 if (quad->isEdgeOut[0]) {
248 // |___|___|___|___|___|___|
250 // |___|___|___|___|___|___|
252 // |___|___|___|___|___|___| __ first row of the regular grid
253 // . . . . . . . . . __ down edge nodes
255 // >->->->->->->->->->->->-> -- direction of processing
257 int g = 0; // number of last processed node in the regular grid
259 // number of last node of the down edge to be processed
260 int stop = nbdown - 1;
261 // if right edge is out, we will stop at a node, previous to the last one
262 if (quad->isEdgeOut[1]) stop--;
264 // for each node of the down edge find nearest node
265 // in the first row of the regular grid and link them
266 for (i = 0; i < stop; i++) {
267 const SMDS_MeshNode *a, *b, *c, *d;
269 b = uv_e0[i + 1].node;
270 gp_Pnt pb (b->X(), b->Y(), b->Z());
272 // find node c in the regular grid, which will be linked with node b
275 // right bound reached, link with the rightmost node
277 c = quad->uv_grid[nbhoriz + iup].node;
280 // find in the grid node c, nearest to the b
281 double mind = RealLast();
282 for (int k = g; k <= iup; k++) {
284 const SMDS_MeshNode *nk;
285 if (k < ilow) // this can be, if left edge is out
286 nk = uv_e3[1].node; // get node from the left edge
288 nk = quad->uv_grid[nbhoriz + k].node; // get one of middle nodes
290 gp_Pnt pnk (nk->X(), nk->Y(), nk->Z());
291 double dist = pb.Distance(pnk);
292 if (dist < mind - eps) {
302 if (near == g) { // make triangle
303 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c);
304 SMDS_MeshFace* face = myTool->AddFace(a, b, c);
305 meshDS->SetMeshElementOnShape(face, geomFaceID);
307 else { // make quadrangle
311 d = quad->uv_grid[nbhoriz + near - 1].node;
312 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c, d);
313 SMDS_MeshFace* face = myTool->AddFace(a, b, c, d);
314 meshDS->SetMeshElementOnShape(face, geomFaceID);
316 // if node d is not at position g - make additional triangles
318 for (int k = near - 1; k > g; k--) {
319 c = quad->uv_grid[nbhoriz + k].node;
323 d = quad->uv_grid[nbhoriz + k - 1].node;
324 //SMDS_MeshFace* face = meshDS->AddFace(a, c, d);
325 SMDS_MeshFace* face = myTool->AddFace(a, c, d);
326 meshDS->SetMeshElementOnShape(face, geomFaceID);
333 if (quad->isEdgeOut[2]) {
336 // <-<-<-<-<-<-<-<-<-<-<-<-< -- direction of processing
338 // . . . . . . . . . __ up edge nodes
339 // ___ ___ ___ ___ ___ ___ __ first row of the regular grid
341 // |___|___|___|___|___|___|
343 // |___|___|___|___|___|___|
346 int g = nbhoriz - 1; // last processed node in the regular grid
349 // if left edge is out, we will stop at a second node
350 if (quad->isEdgeOut[3]) stop++;
352 // for each node of the up edge find nearest node
353 // in the first row of the regular grid and link them
354 for (i = nbup - 1; i > stop; i--) {
355 const SMDS_MeshNode *a, *b, *c, *d;
357 b = uv_e2[i - 1].node;
358 gp_Pnt pb (b->X(), b->Y(), b->Z());
360 // find node c in the grid, which will be linked with node b
362 if (i == stop + 1) { // left bound reached, link with the leftmost node
363 c = quad->uv_grid[nbhoriz*(nbvertic - 2) + ilow].node;
366 // find node c in the grid, nearest to the b
367 double mind = RealLast();
368 for (int k = g; k >= ilow; k--) {
369 const SMDS_MeshNode *nk;
371 nk = uv_e1[nbright - 2].node;
373 nk = quad->uv_grid[nbhoriz*(nbvertic - 2) + k].node;
374 gp_Pnt pnk (nk->X(), nk->Y(), nk->Z());
375 double dist = pb.Distance(pnk);
376 if (dist < mind - eps) {
386 if (near == g) { // make triangle
387 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c);
388 SMDS_MeshFace* face = myTool->AddFace(a, b, c);
389 meshDS->SetMeshElementOnShape(face, geomFaceID);
391 else { // make quadrangle
393 d = uv_e1[nbright - 2].node;
395 d = quad->uv_grid[nbhoriz*(nbvertic - 2) + near + 1].node;
396 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c, d);
397 SMDS_MeshFace* face = myTool->AddFace(a, b, c, d);
398 meshDS->SetMeshElementOnShape(face, geomFaceID);
400 if (near + 1 < g) { // if d not is at g - make additional triangles
401 for (int k = near + 1; k < g; k++) {
402 c = quad->uv_grid[nbhoriz*(nbvertic - 2) + k].node;
404 d = uv_e1[nbright - 2].node;
406 d = quad->uv_grid[nbhoriz*(nbvertic - 2) + k + 1].node;
407 //SMDS_MeshFace* face = meshDS->AddFace(a, c, d);
408 SMDS_MeshFace* face = myTool->AddFace(a, c, d);
409 meshDS->SetMeshElementOnShape(face, geomFaceID);
418 // right or left boundary quadrangles
419 if (quad->isEdgeOut[1]) {
420 // MESSAGE("right edge is out");
421 int g = 0; // last processed node in the grid
422 int stop = nbright - 1;
423 if (quad->isEdgeOut[2]) stop--;
424 for (i = 0; i < stop; i++) {
425 const SMDS_MeshNode *a, *b, *c, *d;
427 b = uv_e1[i + 1].node;
428 gp_Pnt pb (b->X(), b->Y(), b->Z());
430 // find node c in the grid, nearest to the b
432 if (i == stop - 1) { // up bondary reached
433 c = quad->uv_grid[nbhoriz*(jup + 1) - 2].node;
436 double mind = RealLast();
437 for (int k = g; k <= jup; k++) {
438 const SMDS_MeshNode *nk;
440 nk = uv_e0[nbdown - 2].node;
442 nk = quad->uv_grid[nbhoriz*(k + 1) - 2].node;
443 gp_Pnt pnk (nk->X(), nk->Y(), nk->Z());
444 double dist = pb.Distance(pnk);
445 if (dist < mind - eps) {
455 if (near == g) { // make triangle
456 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c);
457 SMDS_MeshFace* face = myTool->AddFace(a, b, c);
458 meshDS->SetMeshElementOnShape(face, geomFaceID);
460 else { // make quadrangle
462 d = uv_e0[nbdown - 2].node;
464 d = quad->uv_grid[nbhoriz*near - 2].node;
465 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c, d);
466 SMDS_MeshFace* face = myTool->AddFace(a, b, c, d);
467 meshDS->SetMeshElementOnShape(face, geomFaceID);
469 if (near - 1 > g) { // if d not is at g - make additional triangles
470 for (int k = near - 1; k > g; k--) {
471 c = quad->uv_grid[nbhoriz*(k + 1) - 2].node;
473 d = uv_e0[nbdown - 2].node;
475 d = quad->uv_grid[nbhoriz*k - 2].node;
476 //SMDS_MeshFace* face = meshDS->AddFace(a, c, d);
477 SMDS_MeshFace* face = myTool->AddFace(a, c, d);
478 meshDS->SetMeshElementOnShape(face, geomFaceID);
485 if (quad->isEdgeOut[3]) {
486 // MESSAGE("left edge is out");
487 int g = nbvertic - 1; // last processed node in the grid
489 if (quad->isEdgeOut[0]) stop++;
490 for (i = nbleft - 1; i > stop; i--) {
491 const SMDS_MeshNode *a, *b, *c, *d;
493 b = uv_e3[i - 1].node;
494 gp_Pnt pb (b->X(), b->Y(), b->Z());
496 // find node c in the grid, nearest to the b
498 if (i == stop + 1) { // down bondary reached
499 c = quad->uv_grid[nbhoriz*jlow + 1].node;
502 double mind = RealLast();
503 for (int k = g; k >= jlow; k--) {
504 const SMDS_MeshNode *nk;
508 nk = quad->uv_grid[nbhoriz*k + 1].node;
509 gp_Pnt pnk (nk->X(), nk->Y(), nk->Z());
510 double dist = pb.Distance(pnk);
511 if (dist < mind - eps) {
521 if (near == g) { // make triangle
522 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c);
523 SMDS_MeshFace* face = myTool->AddFace(a, b, c);
524 meshDS->SetMeshElementOnShape(face, geomFaceID);
526 else { // make quadrangle
530 d = quad->uv_grid[nbhoriz*(near + 1) + 1].node;
531 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c, d);
532 SMDS_MeshFace* face = myTool->AddFace(a, b, c, d);
533 meshDS->SetMeshElementOnShape(face, geomFaceID);
535 if (near + 1 < g) { // if d not is at g - make additional triangles
536 for (int k = near + 1; k < g; k++) {
537 c = quad->uv_grid[nbhoriz*k + 1].node;
541 d = quad->uv_grid[nbhoriz*(k + 1) + 1].node;
542 //SMDS_MeshFace* face = meshDS->AddFace(a, c, d);
543 SMDS_MeshFace* face = myTool->AddFace(a, c, d);
544 meshDS->SetMeshElementOnShape(face, geomFaceID);
557 //=============================================================================
561 //=============================================================================
563 FaceQuadStruct* StdMeshers_Quadrangle_2D::CheckNbEdges(SMESH_Mesh & aMesh,
564 const TopoDS_Shape & aShape)
565 throw(SALOME_Exception)
567 Unexpect aCatch(SalomeException);
569 const TopoDS_Face & F = TopoDS::Face(aShape);
570 const bool ignoreMediumNodes = _quadraticMesh;
572 // verify 1 wire only, with 4 edges
574 list< TopoDS_Edge > edges;
575 list< int > nbEdgesInWire;
576 int nbWire = SMESH_Block::GetOrderedEdges (F, V, edges, nbEdgesInWire);
578 error(COMPERR_BAD_SHAPE, TComm("Wrong number of wires: ") << nbWire);
581 FaceQuadStruct* quad = new FaceQuadStruct;
583 quad->side.reserve(nbEdgesInWire.front());
586 list< TopoDS_Edge >::iterator edgeIt = edges.begin();
587 if ( nbEdgesInWire.front() == 4 ) { // exactly 4 edges
588 for ( ; edgeIt != edges.end(); ++edgeIt, nbSides++ )
589 quad->side.push_back( new StdMeshers_FaceSide(F, *edgeIt, &aMesh,
590 nbSides<TOP_SIDE, ignoreMediumNodes));
592 else if ( nbEdgesInWire.front() > 4 ) { // more than 4 edges - try to unite some
593 list< TopoDS_Edge > sideEdges;
594 while ( !edges.empty()) {
596 sideEdges.splice( sideEdges.end(), edges, edges.begin()); // edges.front() -> sideEdges.end()
597 bool sameSide = true;
598 while ( !edges.empty() && sameSide ) {
599 GeomAbs_Shape cont = SMESH_Algo::Continuity( sideEdges.back(), edges.front() );
600 sameSide = ( cont >= GeomAbs_G1 );
602 sideEdges.splice( sideEdges.end(), edges, edges.begin());
604 if ( nbSides == 0 ) { // go backward from the first edge
606 while ( !edges.empty() && sameSide ) {
607 GeomAbs_Shape cont = SMESH_Algo::Continuity( sideEdges.front(), edges.back() );
608 sameSide = ( cont >= GeomAbs_G1 );
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 )) << " ";
629 error(COMPERR_BAD_SHAPE, TComm("Face must have 4 side but not ") << nbSides);
637 //=============================================================================
641 //=============================================================================
643 FaceQuadStruct *StdMeshers_Quadrangle_2D::CheckAnd2Dcompute
645 const TopoDS_Shape & aShape,
646 const bool CreateQuadratic) throw(SALOME_Exception)
648 Unexpect aCatch(SalomeException);
650 _quadraticMesh = CreateQuadratic;
652 FaceQuadStruct *quad = CheckNbEdges(aMesh, aShape);
656 // set normalized grid on unit square in parametric domain
657 SetNormalizedGrid(aMesh, aShape, quad);
662 //=============================================================================
666 //=============================================================================
668 faceQuadStruct::~faceQuadStruct()
670 for (int i = 0; i < side.size(); i++) {
671 if (side[i]) delete side[i];
673 if (uv_grid) delete [] uv_grid;
677 inline const vector<UVPtStruct>& GetUVPtStructIn(FaceQuadStruct* quad, int i, int nbSeg)
679 bool isXConst = ( i == BOTTOM_SIDE || i == TOP_SIDE );
680 double constValue = ( i == BOTTOM_SIDE || i == LEFT_SIDE ) ? 0 : 1;
683 quad->side[i]->SimulateUVPtStruct(nbSeg,isXConst,constValue) :
684 quad->side[i]->GetUVPtStruct(isXConst,constValue);
688 //=============================================================================
692 //=============================================================================
694 bool StdMeshers_Quadrangle_2D::SetNormalizedGrid (SMESH_Mesh & aMesh,
695 const TopoDS_Shape& aShape,
696 FaceQuadStruct* & quad) throw (SALOME_Exception)
698 Unexpect aCatch(SalomeException);
699 // Algorithme décrit dans "Génération automatique de maillages"
700 // P.L. GEORGE, MASSON, § 6.4.1 p. 84-85
701 // traitement dans le domaine paramétrique 2d u,v
702 // transport - projection sur le carré unité
704 // MESSAGE("StdMeshers_Quadrangle_2D::SetNormalizedGrid");
705 // const TopoDS_Face& F = TopoDS::Face(aShape);
707 // 1 --- find orientation of the 4 edges, by test on extrema
710 // |<----north-2-------^ a3 -------------> a2
712 // west-3 east-1 =right | |
716 // v----south-0--------> a0 -------------> a1
721 // 3 --- 2D normalized values on unit square [0..1][0..1]
723 int nbhoriz = Min(quad->side[0]->NbPoints(), quad->side[2]->NbPoints());
724 int nbvertic = Min(quad->side[1]->NbPoints(), quad->side[3]->NbPoints());
726 quad->isEdgeOut[0] = (quad->side[0]->NbPoints() > quad->side[2]->NbPoints());
727 quad->isEdgeOut[1] = (quad->side[1]->NbPoints() > quad->side[3]->NbPoints());
728 quad->isEdgeOut[2] = (quad->side[2]->NbPoints() > quad->side[0]->NbPoints());
729 quad->isEdgeOut[3] = (quad->side[3]->NbPoints() > quad->side[1]->NbPoints());
731 UVPtStruct *uv_grid = quad->uv_grid = new UVPtStruct[nbvertic * nbhoriz];
733 const vector<UVPtStruct>& uv_e0 = GetUVPtStructIn( quad, 0, nbhoriz - 1 );
734 const vector<UVPtStruct>& uv_e1 = GetUVPtStructIn( quad, 1, nbvertic - 1 );
735 const vector<UVPtStruct>& uv_e2 = GetUVPtStructIn( quad, 2, nbhoriz - 1 );
736 const vector<UVPtStruct>& uv_e3 = GetUVPtStructIn( quad, 3, nbvertic - 1 );
738 if ( uv_e0.empty() || uv_e1.empty() || uv_e2.empty() || uv_e3.empty() )
739 return error(dfltErr(), "Can't find nodes on sides");
741 // nodes Id on "in" edges
742 if (! quad->isEdgeOut[0]) {
744 for (int i = 0; i < nbhoriz; i++) { // down
745 int ij = j * nbhoriz + i;
746 uv_grid[ij].node = uv_e0[i].node;
749 if (! quad->isEdgeOut[1]) {
751 for (int j = 0; j < nbvertic; j++) { // right
752 int ij = j * nbhoriz + i;
753 uv_grid[ij].node = uv_e1[j].node;
756 if (! quad->isEdgeOut[2]) {
757 int j = nbvertic - 1;
758 for (int i = 0; i < nbhoriz; i++) { // up
759 int ij = j * nbhoriz + i;
760 uv_grid[ij].node = uv_e2[i].node;
763 if (! quad->isEdgeOut[3]) {
765 for (int j = 0; j < nbvertic; j++) { // left
766 int ij = j * nbhoriz + i;
767 uv_grid[ij].node = uv_e3[j].node;
771 // normalized 2d values on grid
772 for (int i = 0; i < nbhoriz; i++)
774 for (int j = 0; j < nbvertic; j++)
776 int ij = j * nbhoriz + i;
777 // --- droite i cste : x = x0 + y(x1-x0)
778 double x0 = uv_e0[i].normParam; // bas - sud
779 double x1 = uv_e2[i].normParam; // haut - nord
780 // --- droite j cste : y = y0 + x(y1-y0)
781 double y0 = uv_e3[j].normParam; // gauche-ouest
782 double y1 = uv_e1[j].normParam; // droite - est
783 // --- intersection : x=x0+(y0+x(y1-y0))(x1-x0)
784 double x = (x0 + y0 * (x1 - x0)) / (1 - (y1 - y0) * (x1 - x0));
785 double y = y0 + x * (y1 - y0);
788 //MESSAGE("-xy-01 "<<x0<<" "<<x1<<" "<<y0<<" "<<y1);
789 //MESSAGE("-xy-norm "<<i<<" "<<j<<" "<<x<<" "<<y);
793 // 4 --- projection on 2d domain (u,v)
794 gp_UV a0( uv_e0.front().u, uv_e0.front().v );
795 gp_UV a1( uv_e0.back().u, uv_e0.back().v );
796 gp_UV a2( uv_e2.back().u, uv_e2.back().v );
797 gp_UV a3( uv_e2.front().u, uv_e2.front().v );
799 for (int i = 0; i < nbhoriz; i++)
801 for (int j = 0; j < nbvertic; j++)
803 int ij = j * nbhoriz + i;
804 double x = uv_grid[ij].x;
805 double y = uv_grid[ij].y;
806 double param_0 = uv_e0[0].normParam + x * (uv_e0.back().normParam - uv_e0[0].normParam); // sud
807 double param_2 = uv_e2[0].normParam + x * (uv_e2.back().normParam - uv_e2[0].normParam); // nord
808 double param_1 = uv_e1[0].normParam + y * (uv_e1.back().normParam - uv_e1[0].normParam); // est
809 double param_3 = uv_e3[0].normParam + y * (uv_e3.back().normParam - uv_e3[0].normParam); // ouest
811 //MESSAGE("params "<<param_0<<" "<<param_1<<" "<<param_2<<" "<<param_3);
812 gp_UV p0 = quad->side[0]->Value2d(param_0).XY();
813 gp_UV p1 = quad->side[1]->Value2d(param_1).XY();
814 gp_UV p2 = quad->side[2]->Value2d(param_2).XY();
815 gp_UV p3 = quad->side[3]->Value2d(param_3).XY();
817 gp_UV uv = (1 - y) * p0 + x * p1 + y * p2 + (1 - x) * p3;
818 uv -= (1 - x) * (1 - y) * a0 + x * (1 - y) * a1 + x * y * a2 + (1 - x) * y * a3;
820 uv_grid[ij].u = uv.X();
821 uv_grid[ij].v = uv.Y();
827 //=======================================================================
828 //function : ShiftQuad
829 //purpose : auxilary function for ComputeQuadPref
830 //=======================================================================
832 static void ShiftQuad(FaceQuadStruct* quad, const int num, bool)
834 StdMeshers_FaceSide* side[4] = { quad->side[0], quad->side[1], quad->side[2], quad->side[3] };
835 for (int i = BOTTOM_SIDE; i < NB_SIDES; ++i ) {
836 int id = ( i + num ) % NB_SIDES;
837 bool wasForward = ( i < TOP_SIDE );
838 bool newForward = ( id < TOP_SIDE );
839 if ( wasForward != newForward )
840 side[ i ]->Reverse();
841 quad->side[ id ] = side[ i ];
845 //=======================================================================
847 //purpose : auxilary function for ComputeQuadPref
848 //=======================================================================
850 static gp_UV CalcUV(double x0, double x1, double y0, double y1,
851 FaceQuadStruct* quad,
852 const gp_UV& a0, const gp_UV& a1,
853 const gp_UV& a2, const gp_UV& a3)
855 const vector<UVPtStruct>& uv_eb = quad->side[0]->GetUVPtStruct(true,0 );
856 const vector<UVPtStruct>& uv_er = quad->side[1]->GetUVPtStruct(false,1);
857 const vector<UVPtStruct>& uv_et = quad->side[2]->GetUVPtStruct(true,1 );
858 const vector<UVPtStruct>& uv_el = quad->side[3]->GetUVPtStruct(false,0);
860 double x = (x0 + y0 * (x1 - x0)) / (1 - (y1 - y0) * (x1 - x0));
861 double y = y0 + x * (y1 - y0);
863 double param_b = uv_eb[0].normParam + x * (uv_eb.back().normParam - uv_eb[0].normParam);
864 double param_t = uv_et[0].normParam + x * (uv_et.back().normParam - uv_et[0].normParam);
865 double param_r = uv_er[0].normParam + y * (uv_er.back().normParam - uv_er[0].normParam);
866 double param_l = uv_el[0].normParam + y * (uv_el.back().normParam - uv_el[0].normParam);
868 gp_UV p0 = quad->side[BOTTOM_SIDE]->Value2d(param_b).XY();
869 gp_UV p1 = quad->side[RIGHT_SIDE ]->Value2d(param_r).XY();
870 gp_UV p2 = quad->side[TOP_SIDE ]->Value2d(param_t).XY();
871 gp_UV p3 = quad->side[LEFT_SIDE ]->Value2d(param_l).XY();
873 gp_UV uv = p0 * (1 - y) + p1 * x + p2 * y + p3 * (1 - x);
875 uv -= (1 - x) * (1 - y) * a0 + x * (1 - y) * a1 + x * y * a2 + (1 - x) * y * a3;
880 //=======================================================================
882 * Create only quandrangle faces
884 //=======================================================================
886 bool StdMeshers_Quadrangle_2D::ComputeQuadPref (SMESH_Mesh & aMesh,
887 const TopoDS_Shape& aShape,
888 FaceQuadStruct* quad)
889 throw (SALOME_Exception)
891 Unexpect aCatch(SalomeException);
893 SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
894 const TopoDS_Face& F = TopoDS::Face(aShape);
895 Handle(Geom_Surface) S = BRep_Tool::Surface(F);
896 const TopoDS_Wire& W = BRepTools::OuterWire(F);
898 if(W.Orientation()==TopAbs_FORWARD)
900 //if(WisF) cout<<"W is FORWARD"<<endl;
901 //else cout<<"W is REVERSED"<<endl;
902 bool FisF = (F.Orientation()==TopAbs_FORWARD);
903 if(!FisF) WisF = !WisF;
904 int i,j,geomFaceID = meshDS->ShapeToIndex( F );
906 int nb = quad->side[0]->NbPoints();
907 int nr = quad->side[1]->NbPoints();
908 int nt = quad->side[2]->NbPoints();
909 int nl = quad->side[3]->NbPoints();
915 // it is a base case => not shift quad but me be replacement is need
916 ShiftQuad(quad,0,WisF);
919 // we have to shift quad on 2
920 ShiftQuad(quad,2,WisF);
925 // we have to shift quad on 1
926 ShiftQuad(quad,1,WisF);
929 // we have to shift quad on 3
930 ShiftQuad(quad,3,WisF);
934 nb = quad->side[0]->NbPoints();
935 nr = quad->side[1]->NbPoints();
936 nt = quad->side[2]->NbPoints();
937 nl = quad->side[3]->NbPoints();
940 int nbh = Max(nb,nt);
941 int nbv = Max(nr,nl);
945 // orientation of face and 3 main domain for future faces
951 // left | | | | rigth
967 const vector<UVPtStruct>& uv_eb = quad->side[0]->GetUVPtStruct(true,0 );
968 const vector<UVPtStruct>& uv_er = quad->side[1]->GetUVPtStruct(false,1);
969 const vector<UVPtStruct>& uv_et = quad->side[2]->GetUVPtStruct(true,1 );
970 const vector<UVPtStruct>& uv_el = quad->side[3]->GetUVPtStruct(false,0);
972 // arrays for normalized params
973 //cout<<"Dump B:"<<endl;
974 TColStd_SequenceOfReal npb, npr, npt, npl;
975 for(i=0; i<nb; i++) {
976 npb.Append(uv_eb[i].normParam);
977 //cout<<"i="<<i<<" par="<<uv_eb[i].normParam<<" npar="<<uv_eb[i].normParam;
978 //const SMDS_MeshNode* N = uv_eb[i].node;
979 //cout<<" node("<<N->X()<<","<<N->Y()<<","<<N->Z()<<")"<<endl;
981 for(i=0; i<nr; i++) {
982 npr.Append(uv_er[i].normParam);
984 for(i=0; i<nt; i++) {
985 npt.Append(uv_et[i].normParam);
987 for(i=0; i<nl; i++) {
988 npl.Append(uv_el[i].normParam);
991 // add some params to right and left after the first param
994 double dpr = (npr.Value(2) - npr.Value(1))/(dr+1);
995 for(i=1; i<=dr; i++) {
996 npr.InsertAfter(1,npr.Value(2)-dpr);
1000 dpr = (npl.Value(2) - npl.Value(1))/(dl+1);
1001 for(i=1; i<=dl; i++) {
1002 npl.InsertAfter(1,npl.Value(2)-dpr);
1005 //for(i=1; i<=npb.Length(); i++) {
1006 // cout<<" "<<npb.Value(i);
1010 gp_XY a0( uv_eb.front().u, uv_eb.front().v );
1011 gp_XY a1( uv_eb.back().u, uv_eb.back().v );
1012 gp_XY a2( uv_et.back().u, uv_et.back().v );
1013 gp_XY a3( uv_et.front().u, uv_et.front().v );
1014 //cout<<" a0("<<a0.X()<<","<<a0.Y()<<")"<<" a1("<<a1.X()<<","<<a1.Y()<<")"
1015 // <<" a2("<<a2.X()<<","<<a2.Y()<<")"<<" a3("<<a3.X()<<","<<a3.Y()<<")"<<endl;
1017 int nnn = Min(nr,nl);
1018 // auxilary sequence of XY for creation nodes
1019 // in the bottom part of central domain
1020 // it's length must be == nbv-nnn-1
1021 TColgp_SequenceOfXY UVL;
1022 TColgp_SequenceOfXY UVR;
1024 // step1: create faces for left domain
1025 StdMeshers_Array2OfNode NodesL(1,dl+1,1,nl);
1027 for(j=1; j<=nl; j++)
1028 NodesL.SetValue(1,j,uv_el[j-1].node);
1031 for(i=1; i<=dl; i++)
1032 NodesL.SetValue(i+1,nl,uv_et[i].node);
1033 // create and add needed nodes
1034 TColgp_SequenceOfXY UVtmp;
1035 for(i=1; i<=dl; i++) {
1036 double x0 = npt.Value(i+1);
1039 double y0 = npl.Value(i+1);
1040 double y1 = npr.Value(i+1);
1041 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1042 gp_Pnt P = S->Value(UV.X(),UV.Y());
1043 SMDS_MeshNode * N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1044 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1045 NodesL.SetValue(i+1,1,N);
1046 if(UVL.Length()<nbv-nnn-1) UVL.Append(UV);
1048 for(j=2; j<nl; j++) {
1049 double y0 = npl.Value(dl+j);
1050 double y1 = npr.Value(dl+j);
1051 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1052 gp_Pnt P = S->Value(UV.X(),UV.Y());
1053 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1054 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1055 NodesL.SetValue(i+1,j,N);
1056 if( i==dl ) UVtmp.Append(UV);
1059 for(i=1; i<=UVtmp.Length() && UVL.Length()<nbv-nnn-1; i++) {
1060 UVL.Append(UVtmp.Value(i));
1062 //cout<<"Dump NodesL:"<<endl;
1063 //for(i=1; i<=dl+1; i++) {
1065 // for(j=1; j<=nl; j++) {
1066 // cout<<" ("<<NodesL.Value(i,j)->X()<<","<<NodesL.Value(i,j)->Y()<<","<<NodesL.Value(i,j)->Z()<<")";
1071 for(i=1; i<=dl; i++) {
1072 for(j=1; j<nl; j++) {
1075 myTool->AddFace(NodesL.Value(i,j), NodesL.Value(i+1,j),
1076 NodesL.Value(i+1,j+1), NodesL.Value(i,j+1));
1077 meshDS->SetMeshElementOnShape(F, geomFaceID);
1081 myTool->AddFace(NodesL.Value(i,j), NodesL.Value(i,j+1),
1082 NodesL.Value(i+1,j+1), NodesL.Value(i+1,j));
1083 meshDS->SetMeshElementOnShape(F, geomFaceID);
1089 // fill UVL using c2d
1090 for(i=1; i<npl.Length() && UVL.Length()<nbv-nnn-1; i++) {
1091 UVL.Append( gp_UV ( uv_el[i].u, uv_el[i].v ));
1095 // step2: create faces for right domain
1096 StdMeshers_Array2OfNode NodesR(1,dr+1,1,nr);
1098 for(j=1; j<=nr; j++)
1099 NodesR.SetValue(1,j,uv_er[nr-j].node);
1102 for(i=1; i<=dr; i++)
1103 NodesR.SetValue(i+1,1,uv_et[nt-1-i].node);
1104 // create and add needed nodes
1105 TColgp_SequenceOfXY UVtmp;
1106 for(i=1; i<=dr; i++) {
1107 double x0 = npt.Value(nt-i);
1110 double y0 = npl.Value(i+1);
1111 double y1 = npr.Value(i+1);
1112 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1113 gp_Pnt P = S->Value(UV.X(),UV.Y());
1114 SMDS_MeshNode * N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1115 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1116 NodesR.SetValue(i+1,nr,N);
1117 if(UVR.Length()<nbv-nnn-1) UVR.Append(UV);
1119 for(j=2; j<nr; j++) {
1120 double y0 = npl.Value(nbv-j+1);
1121 double y1 = npr.Value(nbv-j+1);
1122 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1123 gp_Pnt P = S->Value(UV.X(),UV.Y());
1124 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1125 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1126 NodesR.SetValue(i+1,j,N);
1127 if( i==dr ) UVtmp.Prepend(UV);
1130 for(i=1; i<=UVtmp.Length() && UVR.Length()<nbv-nnn-1; i++) {
1131 UVR.Append(UVtmp.Value(i));
1134 for(i=1; i<=dr; i++) {
1135 for(j=1; j<nr; j++) {
1138 myTool->AddFace(NodesR.Value(i,j), NodesR.Value(i+1,j),
1139 NodesR.Value(i+1,j+1), NodesR.Value(i,j+1));
1140 meshDS->SetMeshElementOnShape(F, geomFaceID);
1144 myTool->AddFace(NodesR.Value(i,j), NodesR.Value(i,j+1),
1145 NodesR.Value(i+1,j+1), NodesR.Value(i+1,j));
1146 meshDS->SetMeshElementOnShape(F, geomFaceID);
1152 // fill UVR using c2d
1153 for(i=1; i<npr.Length() && UVR.Length()<nbv-nnn-1; i++) {
1154 UVR.Append( gp_UV( uv_er[i].u, uv_er[i].v ));
1158 // step3: create faces for central domain
1159 StdMeshers_Array2OfNode NodesC(1,nb,1,nbv);
1160 // add first string using NodesL
1161 for(i=1; i<=dl+1; i++)
1162 NodesC.SetValue(1,i,NodesL(i,1));
1163 for(i=2; i<=nl; i++)
1164 NodesC.SetValue(1,dl+i,NodesL(dl+1,i));
1165 // add last string using NodesR
1166 for(i=1; i<=dr+1; i++)
1167 NodesC.SetValue(nb,i,NodesR(i,nr));
1169 NodesC.SetValue(nb,dr+i+1,NodesR(dr+1,nr-i));
1170 // add top nodes (last columns)
1171 for(i=dl+2; i<nbh-dr; i++)
1172 NodesC.SetValue(i-dl,nbv,uv_et[i-1].node);
1173 // add bottom nodes (first columns)
1175 NodesC.SetValue(i,1,uv_eb[i-1].node);
1177 // create and add needed nodes
1178 // add linear layers
1179 for(i=2; i<nb; i++) {
1180 double x0 = npt.Value(dl+i);
1182 for(j=1; j<nnn; j++) {
1183 double y0 = npl.Value(nbv-nnn+j);
1184 double y1 = npr.Value(nbv-nnn+j);
1185 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1186 gp_Pnt P = S->Value(UV.X(),UV.Y());
1187 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1188 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1189 NodesC.SetValue(i,nbv-nnn+j,N);
1192 // add diagonal layers
1193 //cout<<"UVL.Length()="<<UVL.Length()<<" UVR.Length()="<<UVR.Length()<<endl;
1194 //cout<<"Dump UVL:"<<endl;
1195 //for(i=1; i<=UVL.Length(); i++) {
1196 // cout<<" ("<<UVL.Value(i).X()<<","<<UVL.Value(i).Y()<<")";
1199 for(i=1; i<nbv-nnn; i++) {
1200 double du = UVR.Value(i).X() - UVL.Value(i).X();
1201 double dv = UVR.Value(i).Y() - UVL.Value(i).Y();
1202 for(j=2; j<nb; j++) {
1203 double u = UVL.Value(i).X() + du*npb.Value(j);
1204 double v = UVL.Value(i).Y() + dv*npb.Value(j);
1205 gp_Pnt P = S->Value(u,v);
1206 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1207 meshDS->SetNodeOnFace(N, geomFaceID, u, v);
1208 NodesC.SetValue(j,i+1,N);
1212 for(i=1; i<nb; i++) {
1213 for(j=1; j<nbv; j++) {
1216 myTool->AddFace(NodesC.Value(i,j), NodesC.Value(i+1,j),
1217 NodesC.Value(i+1,j+1), NodesC.Value(i,j+1));
1218 meshDS->SetMeshElementOnShape(F, geomFaceID);
1222 myTool->AddFace(NodesC.Value(i,j), NodesC.Value(i,j+1),
1223 NodesC.Value(i+1,j+1), NodesC.Value(i+1,j));
1224 meshDS->SetMeshElementOnShape(F, geomFaceID);