1 // Copyright (C) 2007-2012 CEA/DEN, EDF R&D, OPEN CASCADE
3 // Copyright (C) 2003-2007 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
23 // File : StdMeshers_Quadrangle_2D.cxx
24 // Author : Paul RASCLE, EDF
27 #include "StdMeshers_Quadrangle_2D.hxx"
29 #include "StdMeshers_FaceSide.hxx"
31 #include "StdMeshers_QuadrangleParams.hxx"
33 #include "SMESH_Gen.hxx"
34 #include "SMESH_Mesh.hxx"
35 #include "SMESH_subMesh.hxx"
36 #include "SMESH_MesherHelper.hxx"
37 #include "SMESH_Block.hxx"
38 #include "SMESH_Comment.hxx"
40 #include "SMDS_MeshElement.hxx"
41 #include "SMDS_MeshNode.hxx"
42 #include "SMDS_EdgePosition.hxx"
43 #include "SMDS_FacePosition.hxx"
45 #include <BRep_Tool.hxx>
46 #include <Geom_Surface.hxx>
47 #include <NCollection_DefineArray2.hxx>
48 #include <Precision.hxx>
49 #include <TColStd_SequenceOfReal.hxx>
50 #include <TColStd_SequenceOfInteger.hxx>
51 #include <TColgp_SequenceOfXY.hxx>
53 #include <TopExp_Explorer.hxx>
54 #include <TopTools_ListIteratorOfListOfShape.hxx>
55 #include <TopTools_MapOfShape.hxx>
58 #include "utilities.h"
59 #include "Utils_ExceptHandlers.hxx"
61 #ifndef StdMeshers_Array2OfNode_HeaderFile
62 #define StdMeshers_Array2OfNode_HeaderFile
63 typedef const SMDS_MeshNode* SMDS_MeshNodePtr;
64 DEFINE_BASECOLLECTION (StdMeshers_BaseCollectionNodePtr, SMDS_MeshNodePtr)
65 DEFINE_ARRAY2(StdMeshers_Array2OfNode,
66 StdMeshers_BaseCollectionNodePtr, SMDS_MeshNodePtr)
72 typedef SMESH_Comment TComm;
74 //=============================================================================
78 //=============================================================================
80 StdMeshers_Quadrangle_2D::StdMeshers_Quadrangle_2D (int hypId, int studyId,
82 : SMESH_2D_Algo(hypId, studyId, gen)
84 MESSAGE("StdMeshers_Quadrangle_2D::StdMeshers_Quadrangle_2D");
85 _name = "Quadrangle_2D";
86 _shapeType = (1 << TopAbs_FACE);
87 _compatibleHypothesis.push_back("QuadrangleParams");
88 _compatibleHypothesis.push_back("QuadranglePreference");
89 _compatibleHypothesis.push_back("TrianglePreference");
93 //=============================================================================
97 //=============================================================================
99 StdMeshers_Quadrangle_2D::~StdMeshers_Quadrangle_2D()
101 MESSAGE("StdMeshers_Quadrangle_2D::~StdMeshers_Quadrangle_2D");
104 //=============================================================================
108 //=============================================================================
110 bool StdMeshers_Quadrangle_2D::CheckHypothesis
112 const TopoDS_Shape& aShape,
113 SMESH_Hypothesis::Hypothesis_Status& aStatus)
116 aStatus = SMESH_Hypothesis::HYP_OK;
118 const list <const SMESHDS_Hypothesis * >& hyps =
119 GetUsedHypothesis(aMesh, aShape, false);
120 const SMESHDS_Hypothesis * aHyp = 0;
123 myQuadType = QUAD_STANDARD;
124 myQuadranglePreference = false;
125 myTrianglePreference = false;
127 bool isFirstParams = true;
129 // First assigned hypothesis (if any) is processed now
130 if (hyps.size() > 0) {
132 if (strcmp("QuadrangleParams", aHyp->GetName()) == 0) {
133 const StdMeshers_QuadrangleParams* aHyp1 =
134 (const StdMeshers_QuadrangleParams*)aHyp;
135 myTriaVertexID = aHyp1->GetTriaVertex();
136 myQuadType = aHyp1->GetQuadType();
137 if (myQuadType == QUAD_QUADRANGLE_PREF ||
138 myQuadType == QUAD_QUADRANGLE_PREF_REVERSED)
139 myQuadranglePreference = true;
140 else if (myQuadType == QUAD_TRIANGLE_PREF)
141 myTrianglePreference = true;
143 else if (strcmp("QuadranglePreference", aHyp->GetName()) == 0) {
144 isFirstParams = false;
145 myQuadranglePreference = true;
147 else if (strcmp("TrianglePreference", aHyp->GetName()) == 0){
148 isFirstParams = false;
149 myTrianglePreference = true;
152 isFirstParams = false;
156 // Second(last) assigned hypothesis (if any) is processed now
157 if (hyps.size() > 1) {
160 if (strcmp("QuadranglePreference", aHyp->GetName()) == 0) {
161 myQuadranglePreference = true;
162 myTrianglePreference = false;
163 myQuadType = QUAD_STANDARD;
165 else if (strcmp("TrianglePreference", aHyp->GetName()) == 0){
166 myQuadranglePreference = false;
167 myTrianglePreference = true;
168 myQuadType = QUAD_STANDARD;
172 const StdMeshers_QuadrangleParams* aHyp2 =
173 (const StdMeshers_QuadrangleParams*)aHyp;
174 myTriaVertexID = aHyp2->GetTriaVertex();
176 if (!myQuadranglePreference && !myTrianglePreference) { // priority of hypos
177 myQuadType = aHyp2->GetQuadType();
178 if (myQuadType == QUAD_QUADRANGLE_PREF ||
179 myQuadType == QUAD_QUADRANGLE_PREF_REVERSED)
180 myQuadranglePreference = true;
181 else if (myQuadType == QUAD_TRIANGLE_PREF)
182 myTrianglePreference = true;
190 //=============================================================================
194 //=============================================================================
196 bool StdMeshers_Quadrangle_2D::Compute (SMESH_Mesh& aMesh,
197 const TopoDS_Shape& aShape)// throw (SALOME_Exception)
199 // PAL14921. Enable catching std::bad_alloc and Standard_OutOfMemory outside
200 //Unexpect aCatchSalomeException);
202 SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
203 aMesh.GetSubMesh(aShape);
205 SMESH_MesherHelper helper (aMesh);
208 _quadraticMesh = myHelper->IsQuadraticSubMesh(aShape);
209 myNeedSmooth = false;
211 FaceQuadStruct *quad = CheckNbEdges(aMesh, aShape);
212 std::auto_ptr<FaceQuadStruct> quadDeleter (quad); // to delete quad at exit from Compute()
216 if (myQuadranglePreference) {
217 int n1 = quad->side[0]->NbPoints();
218 int n2 = quad->side[1]->NbPoints();
219 int n3 = quad->side[2]->NbPoints();
220 int n4 = quad->side[3]->NbPoints();
221 int nfull = n1+n2+n3+n4;
224 if (nfull == ntmp && ((n1 != n3) || (n2 != n4))) {
225 // special path for using only quandrangle faces
226 bool ok = ComputeQuadPref(aMesh, aShape, quad);
227 if ( ok && myNeedSmooth )
232 else if (myQuadType == QUAD_REDUCED) {
233 int n1 = quad->side[0]->NbPoints();
234 int n2 = quad->side[1]->NbPoints();
235 int n3 = quad->side[2]->NbPoints();
236 int n4 = quad->side[3]->NbPoints();
239 int n13tmp = n13/2; n13tmp = n13tmp*2;
240 int n24tmp = n24/2; n24tmp = n24tmp*2;
241 if ((n1 == n3 && n2 != n4 && n24tmp == n24) ||
242 (n2 == n4 && n1 != n3 && n13tmp == n13)) {
243 bool ok = ComputeReduced(aMesh, aShape, quad);
244 if ( ok && myNeedSmooth )
250 // set normalized grid on unit square in parametric domain
252 if (!SetNormalizedGrid(aMesh, aShape, quad))
255 // --- compute 3D values on points, store points & quadrangles
257 int nbdown = quad->side[0]->NbPoints();
258 int nbup = quad->side[2]->NbPoints();
260 int nbright = quad->side[1]->NbPoints();
261 int nbleft = quad->side[3]->NbPoints();
263 int nbhoriz = Min(nbdown, nbup);
264 int nbvertic = Min(nbright, nbleft);
266 const TopoDS_Face& F = TopoDS::Face(aShape);
267 Handle(Geom_Surface) S = BRep_Tool::Surface(F);
269 // internal mesh nodes
270 int i, j, geomFaceID = meshDS->ShapeToIndex(F);
271 for (i = 1; i < nbhoriz - 1; i++) {
272 for (j = 1; j < nbvertic - 1; j++) {
273 int ij = j * nbhoriz + i;
274 double u = quad->uv_grid[ij].u;
275 double v = quad->uv_grid[ij].v;
276 gp_Pnt P = S->Value(u, v);
277 SMDS_MeshNode * node = meshDS->AddNode(P.X(), P.Y(), P.Z());
278 meshDS->SetNodeOnFace(node, geomFaceID, u, v);
279 quad->uv_grid[ij].node = node;
286 // --.--.--.--.--.-- nbvertic
292 // ---.----.----.--- 0
293 // 0 > > > > > > > > nbhoriz
299 int iup = nbhoriz - 1;
300 if (quad->isEdgeOut[3]) { ilow++; } else { if (quad->isEdgeOut[1]) iup--; }
303 int jup = nbvertic - 1;
304 if (quad->isEdgeOut[0]) { jlow++; } else { if (quad->isEdgeOut[2]) jup--; }
306 // regular quadrangles
307 for (i = ilow; i < iup; i++) {
308 for (j = jlow; j < jup; j++) {
309 const SMDS_MeshNode *a, *b, *c, *d;
310 a = quad->uv_grid[j * nbhoriz + i].node;
311 b = quad->uv_grid[j * nbhoriz + i + 1].node;
312 c = quad->uv_grid[(j + 1) * nbhoriz + i + 1].node;
313 d = quad->uv_grid[(j + 1) * nbhoriz + i].node;
314 SMDS_MeshFace* face = myHelper->AddFace(a, b, c, d);
316 meshDS->SetMeshElementOnShape(face, geomFaceID);
321 const vector<UVPtStruct>& uv_e0 = quad->side[0]->GetUVPtStruct(true,0);
322 const vector<UVPtStruct>& uv_e1 = quad->side[1]->GetUVPtStruct(false,1);
323 const vector<UVPtStruct>& uv_e2 = quad->side[2]->GetUVPtStruct(true,1);
324 const vector<UVPtStruct>& uv_e3 = quad->side[3]->GetUVPtStruct(false,0);
326 if (uv_e0.empty() || uv_e1.empty() || uv_e2.empty() || uv_e3.empty())
327 return error(COMPERR_BAD_INPUT_MESH);
329 double eps = Precision::Confusion();
331 // Boundary quadrangles
333 if (quad->isEdgeOut[0]) {
336 // |___|___|___|___|___|___|
338 // |___|___|___|___|___|___|
340 // |___|___|___|___|___|___| __ first row of the regular grid
341 // . . . . . . . . . __ down edge nodes
343 // >->->->->->->->->->->->-> -- direction of processing
345 int g = 0; // number of last processed node in the regular grid
347 // number of last node of the down edge to be processed
348 int stop = nbdown - 1;
349 // if right edge is out, we will stop at a node, previous to the last one
350 if (quad->isEdgeOut[1]) stop--;
352 // for each node of the down edge find nearest node
353 // in the first row of the regular grid and link them
354 for (i = 0; i < stop; i++) {
355 const SMDS_MeshNode *a, *b, *c, *d;
357 b = uv_e0[i + 1].node;
358 gp_Pnt pb (b->X(), b->Y(), b->Z());
360 // find node c in the regular grid, which will be linked with node b
363 // right bound reached, link with the rightmost node
365 c = quad->uv_grid[nbhoriz + iup].node;
368 // find in the grid node c, nearest to the b
369 double mind = RealLast();
370 for (int k = g; k <= iup; k++) {
372 const SMDS_MeshNode *nk;
373 if (k < ilow) // this can be, if left edge is out
374 nk = uv_e3[1].node; // get node from the left edge
376 nk = quad->uv_grid[nbhoriz + k].node; // get one of middle nodes
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 = myHelper->AddFace(a, b, c);
392 if (face) meshDS->SetMeshElementOnShape(face, geomFaceID);
394 else { // make quadrangle
398 d = quad->uv_grid[nbhoriz + near - 1].node;
399 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c, d);
401 if (!myTrianglePreference){
402 SMDS_MeshFace* face = myHelper->AddFace(a, b, c, d);
403 if (face) meshDS->SetMeshElementOnShape(face, geomFaceID);
406 SplitQuad(meshDS, geomFaceID, a, b, c, d);
409 // if node d is not at position g - make additional triangles
411 for (int k = near - 1; k > g; k--) {
412 c = quad->uv_grid[nbhoriz + k].node;
416 d = quad->uv_grid[nbhoriz + k - 1].node;
417 SMDS_MeshFace* face = myHelper->AddFace(a, c, d);
418 if (face) meshDS->SetMeshElementOnShape(face, geomFaceID);
425 if (quad->isEdgeOut[2]) {
428 // <-<-<-<-<-<-<-<-<-<-<-<-< -- direction of processing
430 // . . . . . . . . . __ up edge nodes
431 // ___ ___ ___ ___ ___ ___ __ first row of the regular grid
433 // |___|___|___|___|___|___|
435 // |___|___|___|___|___|___|
438 int g = nbhoriz - 1; // last processed node in the regular grid
441 // if left edge is out, we will stop at a second node
442 if (quad->isEdgeOut[3]) stop++;
444 // for each node of the up edge find nearest node
445 // in the first row of the regular grid and link them
446 for (i = nbup - 1; i > stop; i--) {
447 const SMDS_MeshNode *a, *b, *c, *d;
449 b = uv_e2[i - 1].node;
450 gp_Pnt pb (b->X(), b->Y(), b->Z());
452 // find node c in the grid, which will be linked with node b
454 if (i == stop + 1) { // left bound reached, link with the leftmost node
455 c = quad->uv_grid[nbhoriz*(nbvertic - 2) + ilow].node;
458 // find node c in the grid, nearest to the b
459 double mind = RealLast();
460 for (int k = g; k >= ilow; k--) {
461 const SMDS_MeshNode *nk;
463 nk = uv_e1[nbright - 2].node;
465 nk = quad->uv_grid[nbhoriz*(nbvertic - 2) + k].node;
466 gp_Pnt pnk (nk->X(), nk->Y(), nk->Z());
467 double dist = pb.Distance(pnk);
468 if (dist < mind - eps) {
478 if (near == g) { // make triangle
479 SMDS_MeshFace* face = myHelper->AddFace(a, b, c);
480 if (face) meshDS->SetMeshElementOnShape(face, geomFaceID);
482 else { // make quadrangle
484 d = uv_e1[nbright - 2].node;
486 d = quad->uv_grid[nbhoriz*(nbvertic - 2) + near + 1].node;
487 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c, d);
488 if (!myTrianglePreference){
489 SMDS_MeshFace* face = myHelper->AddFace(a, b, c, d);
490 if (face) meshDS->SetMeshElementOnShape(face, geomFaceID);
493 SplitQuad(meshDS, geomFaceID, a, b, c, d);
496 if (near + 1 < g) { // if d not is at g - make additional triangles
497 for (int k = near + 1; k < g; k++) {
498 c = quad->uv_grid[nbhoriz*(nbvertic - 2) + k].node;
500 d = uv_e1[nbright - 2].node;
502 d = quad->uv_grid[nbhoriz*(nbvertic - 2) + k + 1].node;
503 SMDS_MeshFace* face = myHelper->AddFace(a, c, d);
504 if (face) meshDS->SetMeshElementOnShape(face, geomFaceID);
513 // right or left boundary quadrangles
514 if (quad->isEdgeOut[1]) {
515 // MESSAGE("right edge is out");
516 int g = 0; // last processed node in the grid
517 int stop = nbright - 1;
518 if (quad->isEdgeOut[2]) stop--;
519 for (i = 0; i < stop; i++) {
520 const SMDS_MeshNode *a, *b, *c, *d;
522 b = uv_e1[i + 1].node;
523 gp_Pnt pb (b->X(), b->Y(), b->Z());
525 // find node c in the grid, nearest to the b
527 if (i == stop - 1) { // up bondary reached
528 c = quad->uv_grid[nbhoriz*(jup + 1) - 2].node;
531 double mind = RealLast();
532 for (int k = g; k <= jup; k++) {
533 const SMDS_MeshNode *nk;
535 nk = uv_e0[nbdown - 2].node;
537 nk = quad->uv_grid[nbhoriz*(k + 1) - 2].node;
538 gp_Pnt pnk (nk->X(), nk->Y(), nk->Z());
539 double dist = pb.Distance(pnk);
540 if (dist < mind - eps) {
550 if (near == g) { // make triangle
551 SMDS_MeshFace* face = myHelper->AddFace(a, b, c);
552 if (face) meshDS->SetMeshElementOnShape(face, geomFaceID);
554 else { // make quadrangle
556 d = uv_e0[nbdown - 2].node;
558 d = quad->uv_grid[nbhoriz*near - 2].node;
559 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c, d);
561 if (!myTrianglePreference){
562 SMDS_MeshFace* face = myHelper->AddFace(a, b, c, d);
563 if (face) meshDS->SetMeshElementOnShape(face, geomFaceID);
566 SplitQuad(meshDS, geomFaceID, a, b, c, d);
569 if (near - 1 > g) { // if d not is at g - make additional triangles
570 for (int k = near - 1; k > g; k--) {
571 c = quad->uv_grid[nbhoriz*(k + 1) - 2].node;
573 d = uv_e0[nbdown - 2].node;
575 d = quad->uv_grid[nbhoriz*k - 2].node;
576 SMDS_MeshFace* face = myHelper->AddFace(a, c, d);
577 if (face) meshDS->SetMeshElementOnShape(face, geomFaceID);
584 if (quad->isEdgeOut[3]) {
585 // MESSAGE("left edge is out");
586 int g = nbvertic - 1; // last processed node in the grid
588 if (quad->isEdgeOut[0]) stop++;
589 for (i = nbleft - 1; i > stop; i--) {
590 const SMDS_MeshNode *a, *b, *c, *d;
592 b = uv_e3[i - 1].node;
593 gp_Pnt pb (b->X(), b->Y(), b->Z());
595 // find node c in the grid, nearest to the b
597 if (i == stop + 1) { // down bondary reached
598 c = quad->uv_grid[nbhoriz*jlow + 1].node;
601 double mind = RealLast();
602 for (int k = g; k >= jlow; k--) {
603 const SMDS_MeshNode *nk;
607 nk = quad->uv_grid[nbhoriz*k + 1].node;
608 gp_Pnt pnk (nk->X(), nk->Y(), nk->Z());
609 double dist = pb.Distance(pnk);
610 if (dist < mind - eps) {
620 if (near == g) { // make triangle
621 SMDS_MeshFace* face = myHelper->AddFace(a, b, c);
622 if (face) meshDS->SetMeshElementOnShape(face, geomFaceID);
624 else { // make quadrangle
628 d = quad->uv_grid[nbhoriz*(near + 1) + 1].node;
629 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c, d);
630 if (!myTrianglePreference){
631 SMDS_MeshFace* face = myHelper->AddFace(a, b, c, d);
632 if (face) meshDS->SetMeshElementOnShape(face, geomFaceID);
635 SplitQuad(meshDS, geomFaceID, a, b, c, d);
638 if (near + 1 < g) { // if d not is at g - make additional triangles
639 for (int k = near + 1; k < g; k++) {
640 c = quad->uv_grid[nbhoriz*k + 1].node;
644 d = quad->uv_grid[nbhoriz*(k + 1) + 1].node;
645 SMDS_MeshFace* face = myHelper->AddFace(a, c, d);
646 if (face) meshDS->SetMeshElementOnShape(face, geomFaceID);
663 //=============================================================================
667 //=============================================================================
669 bool StdMeshers_Quadrangle_2D::Evaluate(SMESH_Mesh& aMesh,
670 const TopoDS_Shape& aShape,
671 MapShapeNbElems& aResMap)
674 aMesh.GetSubMesh(aShape);
676 std::vector<int> aNbNodes(4);
677 bool IsQuadratic = false;
678 if (!CheckNbEdgesForEvaluate(aMesh, aShape, aResMap, aNbNodes, IsQuadratic)) {
679 std::vector<int> aResVec(SMDSEntity_Last);
680 for (int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
681 SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
682 aResMap.insert(std::make_pair(sm,aResVec));
683 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
684 smError.reset(new SMESH_ComputeError(COMPERR_ALGO_FAILED,"Submesh can not be evaluated",this));
688 if (myQuadranglePreference) {
689 int n1 = aNbNodes[0];
690 int n2 = aNbNodes[1];
691 int n3 = aNbNodes[2];
692 int n4 = aNbNodes[3];
693 int nfull = n1+n2+n3+n4;
696 if (nfull==ntmp && ((n1!=n3) || (n2!=n4))) {
697 // special path for using only quandrangle faces
698 return EvaluateQuadPref(aMesh, aShape, aNbNodes, aResMap, IsQuadratic);
703 int nbdown = aNbNodes[0];
704 int nbup = aNbNodes[2];
706 int nbright = aNbNodes[1];
707 int nbleft = aNbNodes[3];
709 int nbhoriz = Min(nbdown, nbup);
710 int nbvertic = Min(nbright, nbleft);
712 int dh = Max(nbdown, nbup) - nbhoriz;
713 int dv = Max(nbright, nbleft) - nbvertic;
720 int nbNodes = (nbhoriz-2)*(nbvertic-2);
721 //int nbFaces3 = dh + dv + kdh*(nbvertic-1)*2 + kdv*(nbhoriz-1)*2;
722 int nbFaces3 = dh + dv;
723 //if (kdh==1 && kdv==1) nbFaces3 -= 2;
724 //if (dh>0 && dv>0) nbFaces3 -= 2;
725 //int nbFaces4 = (nbhoriz-1-kdh)*(nbvertic-1-kdv);
726 int nbFaces4 = (nbhoriz-1)*(nbvertic-1);
728 std::vector<int> aVec(SMDSEntity_Last);
729 for (int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i] = 0;
731 aVec[SMDSEntity_Quad_Triangle] = nbFaces3;
732 aVec[SMDSEntity_Quad_Quadrangle] = nbFaces4;
733 int nbbndedges = nbdown + nbup + nbright + nbleft -4;
734 int nbintedges = (nbFaces4*4 + nbFaces3*3 - nbbndedges) / 2;
735 aVec[SMDSEntity_Node] = nbNodes + nbintedges;
736 if (aNbNodes.size()==5) {
737 aVec[SMDSEntity_Quad_Triangle] = nbFaces3 + aNbNodes[3] -1;
738 aVec[SMDSEntity_Quad_Quadrangle] = nbFaces4 - aNbNodes[3] +1;
742 aVec[SMDSEntity_Node] = nbNodes;
743 aVec[SMDSEntity_Triangle] = nbFaces3;
744 aVec[SMDSEntity_Quadrangle] = nbFaces4;
745 if (aNbNodes.size()==5) {
746 aVec[SMDSEntity_Triangle] = nbFaces3 + aNbNodes[3] - 1;
747 aVec[SMDSEntity_Quadrangle] = nbFaces4 - aNbNodes[3] + 1;
750 SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
751 aResMap.insert(std::make_pair(sm,aVec));
757 //================================================================================
759 * \brief Return true if only two given edges meat at their common vertex
761 //================================================================================
763 static bool twoEdgesMeatAtVertex(const TopoDS_Edge& e1,
764 const TopoDS_Edge& e2,
768 if (!TopExp::CommonVertex(e1, e2, v))
770 TopTools_ListIteratorOfListOfShape ancestIt(mesh.GetAncestors(v));
771 for (; ancestIt.More() ; ancestIt.Next())
772 if (ancestIt.Value().ShapeType() == TopAbs_EDGE)
773 if (!e1.IsSame(ancestIt.Value()) && !e2.IsSame(ancestIt.Value()))
778 //=============================================================================
782 //=============================================================================
784 FaceQuadStruct* StdMeshers_Quadrangle_2D::CheckNbEdges(SMESH_Mesh & aMesh,
785 const TopoDS_Shape & aShape)
786 //throw(SALOME_Exception)
788 TopoDS_Face F = TopoDS::Face(aShape);
789 if ( F.Orientation() >= TopAbs_INTERNAL ) F.Orientation( TopAbs_FORWARD );
790 const bool ignoreMediumNodes = _quadraticMesh;
792 // verify 1 wire only, with 4 edges
794 list< TopoDS_Edge > edges;
795 list< int > nbEdgesInWire;
796 int nbWire = SMESH_Block::GetOrderedEdges (F, V, edges, nbEdgesInWire);
798 error(COMPERR_BAD_SHAPE, TComm("Wrong number of wires: ") << nbWire);
801 FaceQuadStruct* quad = new FaceQuadStruct;
803 quad->side.reserve(nbEdgesInWire.front());
807 list< TopoDS_Edge >::iterator edgeIt = edges.begin();
808 if (nbEdgesInWire.front() == 3) // exactly 3 edges
810 SMESH_Comment comment;
811 SMESHDS_Mesh* meshDS = aMesh.GetMeshDS();
812 if (myTriaVertexID == -1)
814 comment << "No Base vertex parameter provided for a trilateral geometrical face";
818 TopoDS_Vertex V = TopoDS::Vertex(meshDS->IndexToShape(myTriaVertexID));
820 TopoDS_Edge E1,E2,E3;
821 for (; edgeIt != edges.end(); ++edgeIt) {
822 TopoDS_Edge E = *edgeIt;
823 TopoDS_Vertex VF, VL;
824 TopExp::Vertices(E, VF, VL, true);
827 else if (VL.IsSame(V))
832 if (!E1.IsNull() && !E2.IsNull() && !E3.IsNull())
834 quad->side.push_back(new StdMeshers_FaceSide(F, E1, &aMesh, true, ignoreMediumNodes));
835 quad->side.push_back(new StdMeshers_FaceSide(F, E2, &aMesh, true, ignoreMediumNodes));
836 quad->side.push_back(new StdMeshers_FaceSide(F, E3, &aMesh, false,ignoreMediumNodes));
837 const vector<UVPtStruct>& UVPSleft = quad->side[0]->GetUVPtStruct(true,0);
838 /* vector<UVPtStruct>& UVPStop = */quad->side[1]->GetUVPtStruct(false,1);
839 /* vector<UVPtStruct>& UVPSright = */quad->side[2]->GetUVPtStruct(true,1);
840 const SMDS_MeshNode* aNode = UVPSleft[0].node;
841 gp_Pnt2d aPnt2d(UVPSleft[0].u, UVPSleft[0].v);
842 quad->side.push_back(new StdMeshers_FaceSide(aNode, aPnt2d, quad->side[1]));
846 comment << "Invalid Base vertex parameter: " << myTriaVertexID << " is not among [";
847 TopTools_MapOfShape vMap;
848 for (TopExp_Explorer v(aShape, TopAbs_VERTEX); v.More(); v.Next())
849 if (vMap.Add(v.Current()))
850 comment << meshDS->ShapeToIndex(v.Current()) << (vMap.Extent()==3 ? "]" : ", ");
856 else if (nbEdgesInWire.front() == 4) // exactly 4 edges
858 for (; edgeIt != edges.end(); ++edgeIt, nbSides++)
859 quad->side.push_back(new StdMeshers_FaceSide(F, *edgeIt, &aMesh,
860 nbSides<TOP_SIDE, ignoreMediumNodes));
862 else if (nbEdgesInWire.front() > 4) // more than 4 edges - try to unite some
864 list< TopoDS_Edge > sideEdges;
865 vector< int > degenSides;
866 while (!edges.empty()) {
868 sideEdges.splice(sideEdges.end(), edges, edges.begin()); // edges.front() -> sideEdges.end()
869 bool sameSide = true;
870 while (!edges.empty() && sameSide) {
871 sameSide = SMESH_Algo::IsContinuous(sideEdges.back(), edges.front());
873 sideEdges.splice(sideEdges.end(), edges, edges.begin());
875 if (nbSides == 0) { // go backward from the first edge
877 while (!edges.empty() && sameSide) {
878 sameSide = SMESH_Algo::IsContinuous(sideEdges.front(), edges.back());
880 sideEdges.splice(sideEdges.begin(), edges, --edges.end());
883 if ( sideEdges.size() == 1 && BRep_Tool::Degenerated( sideEdges.front() ))
884 degenSides.push_back( nbSides );
886 quad->side.push_back(new StdMeshers_FaceSide(F, sideEdges, &aMesh,
887 nbSides<TOP_SIDE, ignoreMediumNodes));
890 if ( !degenSides.empty() && nbSides - degenSides.size() == 4 )
893 for ( unsigned i = TOP_SIDE; i < quad->side.size(); ++i )
894 quad->side[i]->Reverse();
896 for ( int i = degenSides.size()-1; i > -1; --i )
898 StdMeshers_FaceSide* degenSide = quad->side[ degenSides[ i ]];
900 quad->side.erase( quad->side.begin() + degenSides[ i ] );
902 for ( unsigned i = TOP_SIDE; i < quad->side.size(); ++i )
903 quad->side[i]->Reverse();
905 nbSides -= degenSides.size();
907 // issue 20222. Try to unite only edges shared by two same faces
909 // delete found sides
910 { FaceQuadStruct cleaner(*quad); }
912 quad->side.reserve(nbEdgesInWire.front());
915 SMESH_Block::GetOrderedEdges (F, V, edges, nbEdgesInWire);
916 while (!edges.empty()) {
918 sideEdges.splice(sideEdges.end(), edges, edges.begin());
919 bool sameSide = true;
920 while (!edges.empty() && sameSide) {
922 SMESH_Algo::IsContinuous(sideEdges.back(), edges.front()) &&
923 twoEdgesMeatAtVertex(sideEdges.back(), edges.front(), aMesh);
925 sideEdges.splice(sideEdges.end(), edges, edges.begin());
927 if (nbSides == 0) { // go backward from the first edge
929 while (!edges.empty() && sameSide) {
931 SMESH_Algo::IsContinuous(sideEdges.front(), edges.back()) &&
932 twoEdgesMeatAtVertex(sideEdges.front(), edges.back(), aMesh);
934 sideEdges.splice(sideEdges.begin(), edges, --edges.end());
937 quad->side.push_back(new StdMeshers_FaceSide(F, sideEdges, &aMesh,
938 nbSides<TOP_SIDE, ignoreMediumNodes));
945 MESSAGE ("StdMeshers_Quadrangle_2D. Edge IDs of " << nbSides << " sides:\n");
946 for (int i = 0; i < nbSides; ++i) {
948 for (int e = 0; e < quad->side[i]->NbEdges(); ++e)
949 MESSAGE (myHelper->GetMeshDS()->ShapeToIndex(quad->side[i]->Edge(e)) << " ");
955 nbSides = nbEdgesInWire.front();
956 error(COMPERR_BAD_SHAPE, TComm("Face must have 4 sides but not ") << nbSides);
965 //=============================================================================
969 //=============================================================================
971 bool StdMeshers_Quadrangle_2D::CheckNbEdgesForEvaluate(SMESH_Mesh& aMesh,
972 const TopoDS_Shape & aShape,
973 MapShapeNbElems& aResMap,
974 std::vector<int>& aNbNodes,
978 const TopoDS_Face & F = TopoDS::Face(aShape);
980 // verify 1 wire only, with 4 edges
982 list< TopoDS_Edge > edges;
983 list< int > nbEdgesInWire;
984 int nbWire = SMESH_Block::GetOrderedEdges (F, V, edges, nbEdgesInWire);
992 list< TopoDS_Edge >::iterator edgeIt = edges.begin();
993 SMESH_subMesh * sm = aMesh.GetSubMesh(*edgeIt);
994 MapShapeNbElemsItr anIt = aResMap.find(sm);
995 if (anIt==aResMap.end()) {
998 std::vector<int> aVec = (*anIt).second;
999 IsQuadratic = (aVec[SMDSEntity_Quad_Edge] > aVec[SMDSEntity_Edge]);
1000 if (nbEdgesInWire.front() == 3) { // exactly 3 edges
1001 if (myTriaVertexID>0) {
1002 SMESHDS_Mesh* meshDS = aMesh.GetMeshDS();
1003 TopoDS_Vertex V = TopoDS::Vertex(meshDS->IndexToShape(myTriaVertexID));
1005 TopoDS_Edge E1,E2,E3;
1006 for (; edgeIt != edges.end(); ++edgeIt) {
1007 TopoDS_Edge E = TopoDS::Edge(*edgeIt);
1008 TopoDS_Vertex VF, VL;
1009 TopExp::Vertices(E, VF, VL, true);
1012 else if (VL.IsSame(V))
1017 SMESH_subMesh * sm = aMesh.GetSubMesh(E1);
1018 MapShapeNbElemsItr anIt = aResMap.find(sm);
1019 if (anIt==aResMap.end()) return false;
1020 std::vector<int> aVec = (*anIt).second;
1022 aNbNodes[0] = (aVec[SMDSEntity_Node]-1)/2 + 2;
1024 aNbNodes[0] = aVec[SMDSEntity_Node] + 2;
1025 sm = aMesh.GetSubMesh(E2);
1026 anIt = aResMap.find(sm);
1027 if (anIt==aResMap.end()) return false;
1028 aVec = (*anIt).second;
1030 aNbNodes[1] = (aVec[SMDSEntity_Node]-1)/2 + 2;
1032 aNbNodes[1] = aVec[SMDSEntity_Node] + 2;
1033 sm = aMesh.GetSubMesh(E3);
1034 anIt = aResMap.find(sm);
1035 if (anIt==aResMap.end()) return false;
1036 aVec = (*anIt).second;
1038 aNbNodes[2] = (aVec[SMDSEntity_Node]-1)/2 + 2;
1040 aNbNodes[2] = aVec[SMDSEntity_Node] + 2;
1041 aNbNodes[3] = aNbNodes[1];
1047 if (nbEdgesInWire.front() == 4) { // exactly 4 edges
1048 for (; edgeIt != edges.end(); edgeIt++) {
1049 SMESH_subMesh * sm = aMesh.GetSubMesh(*edgeIt);
1050 MapShapeNbElemsItr anIt = aResMap.find(sm);
1051 if (anIt==aResMap.end()) {
1054 std::vector<int> aVec = (*anIt).second;
1056 aNbNodes[nbSides] = (aVec[SMDSEntity_Node]-1)/2 + 2;
1058 aNbNodes[nbSides] = aVec[SMDSEntity_Node] + 2;
1062 else if (nbEdgesInWire.front() > 4) { // more than 4 edges - try to unite some
1063 list< TopoDS_Edge > sideEdges;
1064 while (!edges.empty()) {
1066 sideEdges.splice(sideEdges.end(), edges, edges.begin()); // edges.front() -> sideEdges.end()
1067 bool sameSide = true;
1068 while (!edges.empty() && sameSide) {
1069 sameSide = SMESH_Algo::IsContinuous(sideEdges.back(), edges.front());
1071 sideEdges.splice(sideEdges.end(), edges, edges.begin());
1073 if (nbSides == 0) { // go backward from the first edge
1075 while (!edges.empty() && sameSide) {
1076 sameSide = SMESH_Algo::IsContinuous(sideEdges.front(), edges.back());
1078 sideEdges.splice(sideEdges.begin(), edges, --edges.end());
1081 list<TopoDS_Edge>::iterator ite = sideEdges.begin();
1082 aNbNodes[nbSides] = 1;
1083 for (; ite!=sideEdges.end(); ite++) {
1084 SMESH_subMesh * sm = aMesh.GetSubMesh(*ite);
1085 MapShapeNbElemsItr anIt = aResMap.find(sm);
1086 if (anIt==aResMap.end()) {
1089 std::vector<int> aVec = (*anIt).second;
1091 aNbNodes[nbSides] += (aVec[SMDSEntity_Node]-1)/2 + 1;
1093 aNbNodes[nbSides] += aVec[SMDSEntity_Node] + 1;
1097 // issue 20222. Try to unite only edges shared by two same faces
1100 SMESH_Block::GetOrderedEdges (F, V, edges, nbEdgesInWire);
1101 while (!edges.empty()) {
1103 sideEdges.splice(sideEdges.end(), edges, edges.begin());
1104 bool sameSide = true;
1105 while (!edges.empty() && sameSide) {
1107 SMESH_Algo::IsContinuous(sideEdges.back(), edges.front()) &&
1108 twoEdgesMeatAtVertex(sideEdges.back(), edges.front(), aMesh);
1110 sideEdges.splice(sideEdges.end(), edges, edges.begin());
1112 if (nbSides == 0) { // go backward from the first edge
1114 while (!edges.empty() && sameSide) {
1116 SMESH_Algo::IsContinuous(sideEdges.front(), edges.back()) &&
1117 twoEdgesMeatAtVertex(sideEdges.front(), edges.back(), aMesh);
1119 sideEdges.splice(sideEdges.begin(), edges, --edges.end());
1122 list<TopoDS_Edge>::iterator ite = sideEdges.begin();
1123 aNbNodes[nbSides] = 1;
1124 for (; ite!=sideEdges.end(); ite++) {
1125 SMESH_subMesh * sm = aMesh.GetSubMesh(*ite);
1126 MapShapeNbElemsItr anIt = aResMap.find(sm);
1127 if (anIt==aResMap.end()) {
1130 std::vector<int> aVec = (*anIt).second;
1132 aNbNodes[nbSides] += (aVec[SMDSEntity_Node]-1)/2 + 1;
1134 aNbNodes[nbSides] += aVec[SMDSEntity_Node] + 1;
1142 nbSides = nbEdgesInWire.front();
1143 error(COMPERR_BAD_SHAPE, TComm("Face must have 4 sides but not ") << nbSides);
1151 //=============================================================================
1155 //=============================================================================
1157 FaceQuadStruct *StdMeshers_Quadrangle_2D::CheckAnd2Dcompute
1158 (SMESH_Mesh & aMesh,
1159 const TopoDS_Shape & aShape,
1160 const bool CreateQuadratic) //throw(SALOME_Exception)
1162 _quadraticMesh = CreateQuadratic;
1164 FaceQuadStruct *quad = CheckNbEdges(aMesh, aShape);
1166 if (!quad) return 0;
1168 // set normalized grid on unit square in parametric domain
1169 bool stat = SetNormalizedGrid(aMesh, aShape, quad);
1171 if (quad) delete quad;
1178 //=============================================================================
1182 //=============================================================================
1184 faceQuadStruct::~faceQuadStruct()
1186 for (int i = 0; i < side.size(); i++) {
1187 if (side[i]) delete side[i];
1189 if (uv_grid) delete [] uv_grid;
1193 inline const vector<UVPtStruct>& GetUVPtStructIn(FaceQuadStruct* quad, int i, int nbSeg)
1195 bool isXConst = (i == BOTTOM_SIDE || i == TOP_SIDE);
1196 double constValue = (i == BOTTOM_SIDE || i == LEFT_SIDE) ? 0 : 1;
1198 quad->isEdgeOut[i] ?
1199 quad->side[i]->SimulateUVPtStruct(nbSeg,isXConst,constValue) :
1200 quad->side[i]->GetUVPtStruct(isXConst,constValue);
1202 inline gp_UV CalcUV(double x, double y,
1203 const gp_UV& a0,const gp_UV& a1,const gp_UV& a2,const gp_UV& a3,
1204 const gp_UV& p0,const gp_UV& p1,const gp_UV& p2,const gp_UV& p3)
1207 ((1 - y) * p0 + x * p1 + y * p2 + (1 - x) * p3 ) -
1208 ((1 - x) * (1 - y) * a0 + x * (1 - y) * a1 + x * y * a2 + (1 - x) * y * a3);
1212 //=============================================================================
1216 //=============================================================================
1218 bool StdMeshers_Quadrangle_2D::SetNormalizedGrid (SMESH_Mesh & aMesh,
1219 const TopoDS_Shape& aShape,
1220 FaceQuadStruct* & quad) //throw (SALOME_Exception)
1222 // Algorithme décrit dans "Génération automatique de maillages"
1223 // P.L. GEORGE, MASSON, § 6.4.1 p. 84-85
1224 // traitement dans le domaine paramétrique 2d u,v
1225 // transport - projection sur le carré unité
1227 // MESSAGE("StdMeshers_Quadrangle_2D::SetNormalizedGrid");
1228 // const TopoDS_Face& F = TopoDS::Face(aShape);
1230 // 1 --- find orientation of the 4 edges, by test on extrema
1233 // |<----north-2-------^ a3 -------------> a2
1235 // west-3 east-1 =right | |
1239 // v----south-0--------> a0 -------------> a1
1244 // 3 --- 2D normalized values on unit square [0..1][0..1]
1246 int nbhoriz = Min(quad->side[0]->NbPoints(), quad->side[2]->NbPoints());
1247 int nbvertic = Min(quad->side[1]->NbPoints(), quad->side[3]->NbPoints());
1249 quad->isEdgeOut[0] = (quad->side[0]->NbPoints() > quad->side[2]->NbPoints());
1250 quad->isEdgeOut[1] = (quad->side[1]->NbPoints() > quad->side[3]->NbPoints());
1251 quad->isEdgeOut[2] = (quad->side[2]->NbPoints() > quad->side[0]->NbPoints());
1252 quad->isEdgeOut[3] = (quad->side[3]->NbPoints() > quad->side[1]->NbPoints());
1254 UVPtStruct *uv_grid = quad->uv_grid = new UVPtStruct[nbvertic * nbhoriz];
1256 const vector<UVPtStruct>& uv_e0 = GetUVPtStructIn(quad, 0, nbhoriz - 1);
1257 const vector<UVPtStruct>& uv_e1 = GetUVPtStructIn(quad, 1, nbvertic - 1);
1258 const vector<UVPtStruct>& uv_e2 = GetUVPtStructIn(quad, 2, nbhoriz - 1);
1259 const vector<UVPtStruct>& uv_e3 = GetUVPtStructIn(quad, 3, nbvertic - 1);
1261 if (uv_e0.empty() || uv_e1.empty() || uv_e2.empty() || uv_e3.empty())
1262 //return error("Can't find nodes on sides");
1263 return error(COMPERR_BAD_INPUT_MESH);
1266 UpdateDegenUV( quad );
1268 // nodes Id on "in" edges
1269 if (! quad->isEdgeOut[0]) {
1271 for (int i = 0; i < nbhoriz; i++) { // down
1272 int ij = j * nbhoriz + i;
1273 uv_grid[ij].node = uv_e0[i].node;
1276 if (! quad->isEdgeOut[1]) {
1277 int i = nbhoriz - 1;
1278 for (int j = 0; j < nbvertic; j++) { // right
1279 int ij = j * nbhoriz + i;
1280 uv_grid[ij].node = uv_e1[j].node;
1283 if (! quad->isEdgeOut[2]) {
1284 int j = nbvertic - 1;
1285 for (int i = 0; i < nbhoriz; i++) { // up
1286 int ij = j * nbhoriz + i;
1287 uv_grid[ij].node = uv_e2[i].node;
1290 if (! quad->isEdgeOut[3]) {
1292 for (int j = 0; j < nbvertic; j++) { // left
1293 int ij = j * nbhoriz + i;
1294 uv_grid[ij].node = uv_e3[j].node;
1298 // normalized 2d values on grid
1299 for (int i = 0; i < nbhoriz; i++) {
1300 for (int j = 0; j < nbvertic; j++) {
1301 int ij = j * nbhoriz + i;
1302 // --- droite i cste : x = x0 + y(x1-x0)
1303 double x0 = uv_e0[i].normParam; // bas - sud
1304 double x1 = uv_e2[i].normParam; // haut - nord
1305 // --- droite j cste : y = y0 + x(y1-y0)
1306 double y0 = uv_e3[j].normParam; // gauche-ouest
1307 double y1 = uv_e1[j].normParam; // droite - est
1308 // --- intersection : x=x0+(y0+x(y1-y0))(x1-x0)
1309 double x = (x0 + y0 * (x1 - x0)) / (1 - (y1 - y0) * (x1 - x0));
1310 double y = y0 + x * (y1 - y0);
1313 //MESSAGE("-xy-01 "<<x0<<" "<<x1<<" "<<y0<<" "<<y1);
1314 //MESSAGE("-xy-norm "<<i<<" "<<j<<" "<<x<<" "<<y);
1318 // 4 --- projection on 2d domain (u,v)
1319 gp_UV a0(uv_e0.front().u, uv_e0.front().v);
1320 gp_UV a1(uv_e0.back().u, uv_e0.back().v);
1321 gp_UV a2(uv_e2.back().u, uv_e2.back().v);
1322 gp_UV a3(uv_e2.front().u, uv_e2.front().v);
1324 for (int i = 0; i < nbhoriz; i++) {
1325 for (int j = 0; j < nbvertic; j++) {
1326 int ij = j * nbhoriz + i;
1327 double x = uv_grid[ij].x;
1328 double y = uv_grid[ij].y;
1329 double param_0 = uv_e0[0].normParam + x * (uv_e0.back().normParam - uv_e0[0].normParam); // sud
1330 double param_2 = uv_e2[0].normParam + x * (uv_e2.back().normParam - uv_e2[0].normParam); // nord
1331 double param_1 = uv_e1[0].normParam + y * (uv_e1.back().normParam - uv_e1[0].normParam); // est
1332 double param_3 = uv_e3[0].normParam + y * (uv_e3.back().normParam - uv_e3[0].normParam); // ouest
1334 //MESSAGE("params "<<param_0<<" "<<param_1<<" "<<param_2<<" "<<param_3);
1335 gp_UV p0 = quad->side[0]->Value2d(param_0).XY();
1336 gp_UV p1 = quad->side[1]->Value2d(param_1).XY();
1337 gp_UV p2 = quad->side[2]->Value2d(param_2).XY();
1338 gp_UV p3 = quad->side[3]->Value2d(param_3).XY();
1340 gp_UV uv = CalcUV(x,y, a0,a1,a2,a3, p0,p1,p2,p3);
1342 uv_grid[ij].u = uv.X();
1343 uv_grid[ij].v = uv.Y();
1349 //=======================================================================
1350 //function : ShiftQuad
1351 //purpose : auxilary function for ComputeQuadPref
1352 //=======================================================================
1354 static void ShiftQuad(FaceQuadStruct* quad, const int num, bool)
1356 StdMeshers_FaceSide* side[4] = { quad->side[0], quad->side[1], quad->side[2], quad->side[3] };
1357 for (int i = BOTTOM_SIDE; i < NB_SIDES; ++i) {
1358 int id = (i + num) % NB_SIDES;
1359 bool wasForward = (i < TOP_SIDE);
1360 bool newForward = (id < TOP_SIDE);
1361 if (wasForward != newForward)
1362 side[ i ]->Reverse();
1363 quad->side[ id ] = side[ i ];
1367 //=======================================================================
1369 //purpose : auxilary function for ComputeQuadPref
1370 //=======================================================================
1372 static gp_UV CalcUV(double x0, double x1, double y0, double y1,
1373 FaceQuadStruct* quad,
1374 const gp_UV& a0, const gp_UV& a1,
1375 const gp_UV& a2, const gp_UV& a3)
1377 const vector<UVPtStruct>& uv_eb = quad->side[0]->GetUVPtStruct(true,0);
1378 const vector<UVPtStruct>& uv_er = quad->side[1]->GetUVPtStruct(false,1);
1379 const vector<UVPtStruct>& uv_et = quad->side[2]->GetUVPtStruct(true,1);
1380 const vector<UVPtStruct>& uv_el = quad->side[3]->GetUVPtStruct(false,0);
1382 double x = (x0 + y0 * (x1 - x0)) / (1 - (y1 - y0) * (x1 - x0));
1383 double y = y0 + x * (y1 - y0);
1385 double param_b = uv_eb[0].normParam + x * (uv_eb.back().normParam - uv_eb[0].normParam);
1386 double param_t = uv_et[0].normParam + x * (uv_et.back().normParam - uv_et[0].normParam);
1387 double param_r = uv_er[0].normParam + y * (uv_er.back().normParam - uv_er[0].normParam);
1388 double param_l = uv_el[0].normParam + y * (uv_el.back().normParam - uv_el[0].normParam);
1390 gp_UV p0 = quad->side[BOTTOM_SIDE]->Value2d(param_b).XY();
1391 gp_UV p1 = quad->side[RIGHT_SIDE ]->Value2d(param_r).XY();
1392 gp_UV p2 = quad->side[TOP_SIDE ]->Value2d(param_t).XY();
1393 gp_UV p3 = quad->side[LEFT_SIDE ]->Value2d(param_l).XY();
1395 gp_UV uv = CalcUV(x,y, a0,a1,a2,a3, p0,p1,p2,p3);
1400 //=======================================================================
1401 //function : CalcUV2
1402 //purpose : auxilary function for ComputeQuadPref
1403 //=======================================================================
1405 static gp_UV CalcUV2(double x, double y,
1406 FaceQuadStruct* quad,
1407 const gp_UV& a0, const gp_UV& a1,
1408 const gp_UV& a2, const gp_UV& a3)
1410 gp_UV p0 = quad->side[BOTTOM_SIDE]->Value2d(x).XY();
1411 gp_UV p1 = quad->side[RIGHT_SIDE ]->Value2d(y).XY();
1412 gp_UV p2 = quad->side[TOP_SIDE ]->Value2d(x).XY();
1413 gp_UV p3 = quad->side[LEFT_SIDE ]->Value2d(y).XY();
1415 gp_UV uv = CalcUV(x,y, a0,a1,a2,a3, p0,p1,p2,p3);
1421 //=======================================================================
1423 * Create only quandrangle faces
1425 //=======================================================================
1427 bool StdMeshers_Quadrangle_2D::ComputeQuadPref (SMESH_Mesh & aMesh,
1428 const TopoDS_Shape& aShape,
1429 FaceQuadStruct* quad)
1431 // Auxilary key in order to keep old variant
1432 // of meshing after implementation new variant
1433 // for bug 0016220 from Mantis.
1434 bool OldVersion = false;
1435 if (myQuadType == QUAD_QUADRANGLE_PREF_REVERSED)
1438 SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
1439 const TopoDS_Face& F = TopoDS::Face(aShape);
1440 Handle(Geom_Surface) S = BRep_Tool::Surface(F);
1442 int i,j,geomFaceID = meshDS->ShapeToIndex(F);
1444 int nb = quad->side[0]->NbPoints();
1445 int nr = quad->side[1]->NbPoints();
1446 int nt = quad->side[2]->NbPoints();
1447 int nl = quad->side[3]->NbPoints();
1448 int dh = abs(nb-nt);
1449 int dv = abs(nr-nl);
1453 // it is a base case => not shift quad but me be replacement is need
1454 ShiftQuad(quad,0,WisF);
1457 // we have to shift quad on 2
1458 ShiftQuad(quad,2,WisF);
1463 // we have to shift quad on 1
1464 ShiftQuad(quad,1,WisF);
1467 // we have to shift quad on 3
1468 ShiftQuad(quad,3,WisF);
1472 nb = quad->side[0]->NbPoints();
1473 nr = quad->side[1]->NbPoints();
1474 nt = quad->side[2]->NbPoints();
1475 nl = quad->side[3]->NbPoints();
1478 int nbh = Max(nb,nt);
1479 int nbv = Max(nr,nl);
1483 // ----------- Old version ---------------
1484 // orientation of face and 3 main domain for future faces
1490 // left | | | | rigth
1497 // ----------- New version ---------------
1498 // orientation of face and 3 main domain for future faces
1504 // left |/________\| rigth
1520 const vector<UVPtStruct>& uv_eb = quad->side[0]->GetUVPtStruct(true,0);
1521 const vector<UVPtStruct>& uv_er = quad->side[1]->GetUVPtStruct(false,1);
1522 const vector<UVPtStruct>& uv_et = quad->side[2]->GetUVPtStruct(true,1);
1523 const vector<UVPtStruct>& uv_el = quad->side[3]->GetUVPtStruct(false,0);
1525 if (uv_eb.size() != nb || uv_er.size() != nr || uv_et.size() != nt || uv_el.size() != nl)
1526 return error(COMPERR_BAD_INPUT_MESH);
1529 UpdateDegenUV( quad );
1531 // arrays for normalized params
1532 //cout<<"Dump B:"<<endl;
1533 TColStd_SequenceOfReal npb, npr, npt, npl;
1534 for (i=0; i<nb; i++) {
1535 npb.Append(uv_eb[i].normParam);
1536 //cout<<"i="<<i<<" par="<<uv_eb[i].normParam<<" npar="<<uv_eb[i].normParam;
1537 //const SMDS_MeshNode* N = uv_eb[i].node;
1538 //cout<<" node("<<N->X()<<","<<N->Y()<<","<<N->Z()<<")"<<endl;
1540 for (i=0; i<nr; i++) {
1541 npr.Append(uv_er[i].normParam);
1543 for (i=0; i<nt; i++) {
1544 npt.Append(uv_et[i].normParam);
1546 for (i=0; i<nl; i++) {
1547 npl.Append(uv_el[i].normParam);
1552 // add some params to right and left after the first param
1555 double dpr = (npr.Value(2) - npr.Value(1))/(dr+1);
1556 for (i=1; i<=dr; i++) {
1557 npr.InsertAfter(1,npr.Value(2)-dpr);
1561 dpr = (npl.Value(2) - npl.Value(1))/(dl+1);
1562 for (i=1; i<=dl; i++) {
1563 npl.InsertAfter(1,npl.Value(2)-dpr);
1567 //for (i=1; i<=npb.Length(); i++) {
1568 // cout<<" "<<npb.Value(i);
1572 gp_XY a0(uv_eb.front().u, uv_eb.front().v);
1573 gp_XY a1(uv_eb.back().u, uv_eb.back().v);
1574 gp_XY a2(uv_et.back().u, uv_et.back().v);
1575 gp_XY a3(uv_et.front().u, uv_et.front().v);
1576 //cout<<" a0("<<a0.X()<<","<<a0.Y()<<")"<<" a1("<<a1.X()<<","<<a1.Y()<<")"
1577 // <<" a2("<<a2.X()<<","<<a2.Y()<<")"<<" a3("<<a3.X()<<","<<a3.Y()<<")"<<endl;
1579 int nnn = Min(nr,nl);
1580 // auxilary sequence of XY for creation nodes
1581 // in the bottom part of central domain
1582 // Length of UVL and UVR must be == nbv-nnn
1583 TColgp_SequenceOfXY UVL, UVR, UVT;
1586 // step1: create faces for left domain
1587 StdMeshers_Array2OfNode NodesL(1,dl+1,1,nl);
1589 for (j=1; j<=nl; j++)
1590 NodesL.SetValue(1,j,uv_el[j-1].node);
1593 for (i=1; i<=dl; i++)
1594 NodesL.SetValue(i+1,nl,uv_et[i].node);
1595 // create and add needed nodes
1596 TColgp_SequenceOfXY UVtmp;
1597 for (i=1; i<=dl; i++) {
1598 double x0 = npt.Value(i+1);
1601 double y0 = npl.Value(i+1);
1602 double y1 = npr.Value(i+1);
1603 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1604 gp_Pnt P = S->Value(UV.X(),UV.Y());
1605 SMDS_MeshNode * N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1606 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1607 NodesL.SetValue(i+1,1,N);
1608 if (UVL.Length()<nbv-nnn) UVL.Append(UV);
1610 for (j=2; j<nl; j++) {
1611 double y0 = npl.Value(dl+j);
1612 double y1 = npr.Value(dl+j);
1613 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1614 gp_Pnt P = S->Value(UV.X(),UV.Y());
1615 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1616 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1617 NodesL.SetValue(i+1,j,N);
1618 if (i==dl) UVtmp.Append(UV);
1621 for (i=1; i<=UVtmp.Length() && UVL.Length()<nbv-nnn; i++) {
1622 UVL.Append(UVtmp.Value(i));
1624 //cout<<"Dump NodesL:"<<endl;
1625 //for (i=1; i<=dl+1; i++) {
1627 // for (j=1; j<=nl; j++) {
1628 // cout<<" ("<<NodesL.Value(i,j)->X()<<","<<NodesL.Value(i,j)->Y()<<","<<NodesL.Value(i,j)->Z()<<")";
1633 for (i=1; i<=dl; i++) {
1634 for (j=1; j<nl; j++) {
1637 myHelper->AddFace(NodesL.Value(i,j), NodesL.Value(i+1,j),
1638 NodesL.Value(i+1,j+1), NodesL.Value(i,j+1));
1639 if (F) meshDS->SetMeshElementOnShape(F, geomFaceID);
1643 myHelper->AddFace(NodesL.Value(i,j), NodesL.Value(i,j+1),
1644 NodesL.Value(i+1,j+1), NodesL.Value(i+1,j));
1645 if (F) meshDS->SetMeshElementOnShape(F, geomFaceID);
1651 // fill UVL using c2d
1652 for (i=1; i<npl.Length() && UVL.Length()<nbv-nnn; i++) {
1653 UVL.Append(gp_UV (uv_el[i].u, uv_el[i].v));
1657 // step2: create faces for right domain
1658 StdMeshers_Array2OfNode NodesR(1,dr+1,1,nr);
1660 for (j=1; j<=nr; j++)
1661 NodesR.SetValue(1,j,uv_er[nr-j].node);
1664 for (i=1; i<=dr; i++)
1665 NodesR.SetValue(i+1,1,uv_et[nt-1-i].node);
1666 // create and add needed nodes
1667 TColgp_SequenceOfXY UVtmp;
1668 for (i=1; i<=dr; i++) {
1669 double x0 = npt.Value(nt-i);
1672 double y0 = npl.Value(i+1);
1673 double y1 = npr.Value(i+1);
1674 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1675 gp_Pnt P = S->Value(UV.X(),UV.Y());
1676 SMDS_MeshNode * N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1677 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1678 NodesR.SetValue(i+1,nr,N);
1679 if (UVR.Length()<nbv-nnn) UVR.Append(UV);
1681 for (j=2; j<nr; j++) {
1682 double y0 = npl.Value(nbv-j+1);
1683 double y1 = npr.Value(nbv-j+1);
1684 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1685 gp_Pnt P = S->Value(UV.X(),UV.Y());
1686 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1687 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1688 NodesR.SetValue(i+1,j,N);
1689 if (i==dr) UVtmp.Prepend(UV);
1692 for (i=1; i<=UVtmp.Length() && UVR.Length()<nbv-nnn; i++) {
1693 UVR.Append(UVtmp.Value(i));
1696 for (i=1; i<=dr; i++) {
1697 for (j=1; j<nr; j++) {
1700 myHelper->AddFace(NodesR.Value(i,j), NodesR.Value(i+1,j),
1701 NodesR.Value(i+1,j+1), NodesR.Value(i,j+1));
1702 if (F) meshDS->SetMeshElementOnShape(F, geomFaceID);
1706 myHelper->AddFace(NodesR.Value(i,j), NodesR.Value(i,j+1),
1707 NodesR.Value(i+1,j+1), NodesR.Value(i+1,j));
1708 if (F) meshDS->SetMeshElementOnShape(F, geomFaceID);
1714 // fill UVR using c2d
1715 for (i=1; i<npr.Length() && UVR.Length()<nbv-nnn; i++) {
1716 UVR.Append(gp_UV(uv_er[i].u, uv_er[i].v));
1720 // step3: create faces for central domain
1721 StdMeshers_Array2OfNode NodesC(1,nb,1,nbv);
1722 // add first line using NodesL
1723 for (i=1; i<=dl+1; i++)
1724 NodesC.SetValue(1,i,NodesL(i,1));
1725 for (i=2; i<=nl; i++)
1726 NodesC.SetValue(1,dl+i,NodesL(dl+1,i));
1727 // add last line using NodesR
1728 for (i=1; i<=dr+1; i++)
1729 NodesC.SetValue(nb,i,NodesR(i,nr));
1730 for (i=1; i<nr; i++)
1731 NodesC.SetValue(nb,dr+i+1,NodesR(dr+1,nr-i));
1732 // add top nodes (last columns)
1733 for (i=dl+2; i<nbh-dr; i++)
1734 NodesC.SetValue(i-dl,nbv,uv_et[i-1].node);
1735 // add bottom nodes (first columns)
1736 for (i=2; i<nb; i++)
1737 NodesC.SetValue(i,1,uv_eb[i-1].node);
1739 // create and add needed nodes
1740 // add linear layers
1741 for (i=2; i<nb; i++) {
1742 double x0 = npt.Value(dl+i);
1744 for (j=1; j<nnn; j++) {
1745 double y0 = npl.Value(nbv-nnn+j);
1746 double y1 = npr.Value(nbv-nnn+j);
1747 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1748 gp_Pnt P = S->Value(UV.X(),UV.Y());
1749 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1750 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1751 NodesC.SetValue(i,nbv-nnn+j,N);
1756 // add diagonal layers
1757 //cout<<"UVL.Length()="<<UVL.Length()<<" UVR.Length()="<<UVR.Length()<<endl;
1758 //cout<<"Dump UVL:"<<endl;
1759 //for (i=1; i<=UVL.Length(); i++) {
1760 // cout<<" ("<<UVL.Value(i).X()<<","<<UVL.Value(i).Y()<<")";
1763 gp_UV A2 = UVR.Value(nbv-nnn);
1764 gp_UV A3 = UVL.Value(nbv-nnn);
1765 for (i=1; i<nbv-nnn; i++) {
1766 gp_UV p1 = UVR.Value(i);
1767 gp_UV p3 = UVL.Value(i);
1768 double y = i / double(nbv-nnn);
1769 for (j=2; j<nb; j++) {
1770 double x = npb.Value(j);
1771 gp_UV p0( uv_eb[j-1].u, uv_eb[j-1].v );
1772 gp_UV p2 = UVT.Value( j-1 );
1773 gp_UV UV = CalcUV(x, y, a0, a1, A2, A3, p0,p1,p2,p3 );
1774 gp_Pnt P = S->Value(UV.X(),UV.Y());
1775 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1776 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(),UV.Y());
1777 NodesC.SetValue(j,i+1,N);
1781 for (i=1; i<nb; i++) {
1782 for (j=1; j<nbv; j++) {
1785 myHelper->AddFace(NodesC.Value(i,j), NodesC.Value(i+1,j),
1786 NodesC.Value(i+1,j+1), NodesC.Value(i,j+1));
1787 if (F) meshDS->SetMeshElementOnShape(F, geomFaceID);
1791 myHelper->AddFace(NodesC.Value(i,j), NodesC.Value(i,j+1),
1792 NodesC.Value(i+1,j+1), NodesC.Value(i+1,j));
1793 if (F) meshDS->SetMeshElementOnShape(F, geomFaceID);
1799 else { // New version (!OldVersion)
1800 // step1: create faces for bottom rectangle domain
1801 StdMeshers_Array2OfNode NodesBRD(1,nb,1,nnn-1);
1802 // fill UVL and UVR using c2d
1803 for (j=0; j<nb; j++) {
1804 NodesBRD.SetValue(j+1,1,uv_eb[j].node);
1806 for (i=1; i<nnn-1; i++) {
1807 NodesBRD.SetValue(1,i+1,uv_el[i].node);
1808 NodesBRD.SetValue(nb,i+1,uv_er[i].node);
1809 for (j=2; j<nb; j++) {
1810 double x = npb.Value(j);
1811 double y = (1-x) * npl.Value(i+1) + x * npr.Value(i+1);
1812 gp_UV UV = CalcUV2(x, y, quad, a0, a1, a2, a3);
1813 gp_Pnt P = S->Value(UV.X(),UV.Y());
1814 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1815 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(),UV.Y());
1816 NodesBRD.SetValue(j,i+1,N);
1819 for (j=1; j<nnn-1; j++) {
1820 for (i=1; i<nb; i++) {
1823 myHelper->AddFace(NodesBRD.Value(i,j), NodesBRD.Value(i+1,j),
1824 NodesBRD.Value(i+1,j+1), NodesBRD.Value(i,j+1));
1825 if (F) meshDS->SetMeshElementOnShape(F, geomFaceID);
1829 myHelper->AddFace(NodesBRD.Value(i,j), NodesBRD.Value(i,j+1),
1830 NodesBRD.Value(i+1,j+1), NodesBRD.Value(i+1,j));
1831 if (F) meshDS->SetMeshElementOnShape(F, geomFaceID);
1835 int drl = abs(nr-nl);
1836 // create faces for region C
1837 StdMeshers_Array2OfNode NodesC(1,nb,1,drl+1+addv);
1838 // add nodes from previous region
1839 for (j=1; j<=nb; j++) {
1840 NodesC.SetValue(j,1,NodesBRD.Value(j,nnn-1));
1842 if ((drl+addv) > 0) {
1847 TColgp_SequenceOfXY UVtmp;
1848 double drparam = npr.Value(nr) - npr.Value(nnn-1);
1849 double dlparam = npl.Value(nnn) - npl.Value(nnn-1);
1851 for (i=1; i<=drl; i++) {
1852 // add existed nodes from right edge
1853 NodesC.SetValue(nb,i+1,uv_er[nnn+i-2].node);
1854 //double dtparam = npt.Value(i+1);
1855 y1 = npr.Value(nnn+i-1); // param on right edge
1856 double dpar = (y1 - npr.Value(nnn-1))/drparam;
1857 y0 = npl.Value(nnn-1) + dpar*dlparam; // param on left edge
1858 double dy = y1 - y0;
1859 for (j=1; j<nb; j++) {
1860 double x = npt.Value(i+1) + npb.Value(j)*(1-npt.Value(i+1));
1861 double y = y0 + dy*x;
1862 gp_UV UV = CalcUV2(x, y, quad, a0, a1, a2, a3);
1863 gp_Pnt P = S->Value(UV.X(),UV.Y());
1864 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1865 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1866 NodesC.SetValue(j,i+1,N);
1869 double dy0 = (1-y0)/(addv+1);
1870 double dy1 = (1-y1)/(addv+1);
1871 for (i=1; i<=addv; i++) {
1872 double yy0 = y0 + dy0*i;
1873 double yy1 = y1 + dy1*i;
1874 double dyy = yy1 - yy0;
1875 for (j=1; j<=nb; j++) {
1876 double x = npt.Value(i+1+drl) +
1877 npb.Value(j) * (npt.Value(nt-i) - npt.Value(i+1+drl));
1878 double y = yy0 + dyy*x;
1879 gp_UV UV = CalcUV2(x, y, quad, a0, a1, a2, a3);
1880 gp_Pnt P = S->Value(UV.X(),UV.Y());
1881 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1882 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1883 NodesC.SetValue(j,i+drl+1,N);
1890 TColgp_SequenceOfXY UVtmp;
1891 double dlparam = npl.Value(nl) - npl.Value(nnn-1);
1892 double drparam = npr.Value(nnn) - npr.Value(nnn-1);
1893 double y0 = npl.Value(nnn-1);
1894 double y1 = npr.Value(nnn-1);
1895 for (i=1; i<=drl; i++) {
1896 // add existed nodes from right edge
1897 NodesC.SetValue(1,i+1,uv_el[nnn+i-2].node);
1898 y0 = npl.Value(nnn+i-1); // param on left edge
1899 double dpar = (y0 - npl.Value(nnn-1))/dlparam;
1900 y1 = npr.Value(nnn-1) + dpar*drparam; // param on right edge
1901 double dy = y1 - y0;
1902 for (j=2; j<=nb; j++) {
1903 double x = npb.Value(j)*npt.Value(nt-i);
1904 double y = y0 + dy*x;
1905 gp_UV UV = CalcUV2(x, y, quad, a0, a1, a2, a3);
1906 gp_Pnt P = S->Value(UV.X(),UV.Y());
1907 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1908 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1909 NodesC.SetValue(j,i+1,N);
1912 double dy0 = (1-y0)/(addv+1);
1913 double dy1 = (1-y1)/(addv+1);
1914 for (i=1; i<=addv; i++) {
1915 double yy0 = y0 + dy0*i;
1916 double yy1 = y1 + dy1*i;
1917 double dyy = yy1 - yy0;
1918 for (j=1; j<=nb; j++) {
1919 double x = npt.Value(i+1) +
1920 npb.Value(j) * (npt.Value(nt-i-drl) - npt.Value(i+1));
1921 double y = yy0 + dyy*x;
1922 gp_UV UV = CalcUV2(x, y, quad, a0, a1, a2, a3);
1923 gp_Pnt P = S->Value(UV.X(),UV.Y());
1924 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1925 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1926 NodesC.SetValue(j,i+drl+1,N);
1931 for (j=1; j<=drl+addv; j++) {
1932 for (i=1; i<nb; i++) {
1935 myHelper->AddFace(NodesC.Value(i,j), NodesC.Value(i+1,j),
1936 NodesC.Value(i+1,j+1), NodesC.Value(i,j+1));
1937 if (F) meshDS->SetMeshElementOnShape(F, geomFaceID);
1941 myHelper->AddFace(NodesC.Value(i,j), NodesC.Value(i,j+1),
1942 NodesC.Value(i+1,j+1), NodesC.Value(i+1,j));
1943 if (F) meshDS->SetMeshElementOnShape(F, geomFaceID);
1948 StdMeshers_Array2OfNode NodesLast(1,nt,1,2);
1949 for (i=1; i<=nt; i++) {
1950 NodesLast.SetValue(i,2,uv_et[i-1].node);
1953 for (i=n1; i<drl+addv+1; i++) {
1955 NodesLast.SetValue(nnn,1,NodesC.Value(1,i));
1957 for (i=1; i<=nb; i++) {
1959 NodesLast.SetValue(nnn,1,NodesC.Value(i,drl+addv+1));
1961 for (i=drl+addv; i>=n2; i--) {
1963 NodesLast.SetValue(nnn,1,NodesC.Value(nb,i));
1965 for (i=1; i<nt; i++) {
1968 myHelper->AddFace(NodesLast.Value(i,1), NodesLast.Value(i+1,1),
1969 NodesLast.Value(i+1,2), NodesLast.Value(i,2));
1970 if (F) meshDS->SetMeshElementOnShape(F, geomFaceID);
1974 myHelper->AddFace(NodesLast.Value(i,1), NodesLast.Value(i,2),
1975 NodesLast.Value(i+1,2), NodesLast.Value(i+1,2));
1976 if (F) meshDS->SetMeshElementOnShape(F, geomFaceID);
1979 } // if ((drl+addv) > 0)
1981 } // end new version implementation
1988 //=======================================================================
1990 * Evaluate only quandrangle faces
1992 //=======================================================================
1994 bool StdMeshers_Quadrangle_2D::EvaluateQuadPref(SMESH_Mesh & aMesh,
1995 const TopoDS_Shape& aShape,
1996 std::vector<int>& aNbNodes,
1997 MapShapeNbElems& aResMap,
2000 // Auxilary key in order to keep old variant
2001 // of meshing after implementation new variant
2002 // for bug 0016220 from Mantis.
2003 bool OldVersion = false;
2004 if (myQuadType == QUAD_QUADRANGLE_PREF_REVERSED)
2007 const TopoDS_Face& F = TopoDS::Face(aShape);
2008 Handle(Geom_Surface) S = BRep_Tool::Surface(F);
2010 int nb = aNbNodes[0];
2011 int nr = aNbNodes[1];
2012 int nt = aNbNodes[2];
2013 int nl = aNbNodes[3];
2014 int dh = abs(nb-nt);
2015 int dv = abs(nr-nl);
2019 // it is a base case => not shift
2022 // we have to shift on 2
2031 // we have to shift quad on 1
2038 // we have to shift quad on 3
2048 int nbh = Max(nb,nt);
2049 int nbv = Max(nr,nl);
2064 // add some params to right and left after the first param
2071 int nnn = Min(nr,nl);
2076 // step1: create faces for left domain
2078 nbNodes += dl*(nl-1);
2079 nbFaces += dl*(nl-1);
2081 // step2: create faces for right domain
2083 nbNodes += dr*(nr-1);
2084 nbFaces += dr*(nr-1);
2086 // step3: create faces for central domain
2087 nbNodes += (nb-2)*(nnn-1) + (nbv-nnn-1)*(nb-2);
2088 nbFaces += (nb-1)*(nbv-1);
2090 else { // New version (!OldVersion)
2091 nbNodes += (nnn-2)*(nb-2);
2092 nbFaces += (nnn-2)*(nb-1);
2093 int drl = abs(nr-nl);
2094 nbNodes += drl*(nb-1) + addv*nb;
2095 nbFaces += (drl+addv)*(nb-1) + (nt-1);
2096 } // end new version implementation
2098 std::vector<int> aVec(SMDSEntity_Last);
2099 for (int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aVec[i] = 0;
2101 aVec[SMDSEntity_Quad_Quadrangle] = nbFaces;
2102 aVec[SMDSEntity_Node] = nbNodes + nbFaces*4;
2103 if (aNbNodes.size()==5) {
2104 aVec[SMDSEntity_Quad_Triangle] = aNbNodes[3] - 1;
2105 aVec[SMDSEntity_Quad_Quadrangle] = nbFaces - aNbNodes[3] + 1;
2109 aVec[SMDSEntity_Node] = nbNodes;
2110 aVec[SMDSEntity_Quadrangle] = nbFaces;
2111 if (aNbNodes.size()==5) {
2112 aVec[SMDSEntity_Triangle] = aNbNodes[3] - 1;
2113 aVec[SMDSEntity_Quadrangle] = nbFaces - aNbNodes[3] + 1;
2116 SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
2117 aResMap.insert(std::make_pair(sm,aVec));
2123 //=============================================================================
2124 /*! Split quadrangle in to 2 triangles by smallest diagonal
2127 //=============================================================================
2128 void StdMeshers_Quadrangle_2D::SplitQuad(SMESHDS_Mesh *theMeshDS,
2130 const SMDS_MeshNode* theNode1,
2131 const SMDS_MeshNode* theNode2,
2132 const SMDS_MeshNode* theNode3,
2133 const SMDS_MeshNode* theNode4)
2135 gp_Pnt a(theNode1->X(),theNode1->Y(),theNode1->Z());
2136 gp_Pnt b(theNode2->X(),theNode2->Y(),theNode2->Z());
2137 gp_Pnt c(theNode3->X(),theNode3->Y(),theNode3->Z());
2138 gp_Pnt d(theNode4->X(),theNode4->Y(),theNode4->Z());
2139 SMDS_MeshFace* face;
2140 if (a.Distance(c) > b.Distance(d)){
2141 face = myHelper->AddFace(theNode2, theNode4 , theNode1);
2142 if (face) theMeshDS->SetMeshElementOnShape(face, theFaceID);
2143 face = myHelper->AddFace(theNode2, theNode3, theNode4);
2144 if (face) theMeshDS->SetMeshElementOnShape(face, theFaceID);
2148 face = myHelper->AddFace(theNode1, theNode2 ,theNode3);
2149 if (face) theMeshDS->SetMeshElementOnShape(face, theFaceID);
2150 face = myHelper->AddFace(theNode1, theNode3, theNode4);
2151 if (face) theMeshDS->SetMeshElementOnShape(face, theFaceID);
2157 enum uvPos { UV_A0, UV_A1, UV_A2, UV_A3, UV_B, UV_R, UV_T, UV_L, UV_SIZE };
2159 inline SMDS_MeshNode* makeNode( UVPtStruct & uvPt,
2161 FaceQuadStruct* quad,
2163 SMESH_MesherHelper* helper,
2164 Handle(Geom_Surface) S)
2166 const vector<UVPtStruct>& uv_eb = quad->side[BOTTOM_SIDE]->GetUVPtStruct();
2167 const vector<UVPtStruct>& uv_et = quad->side[TOP_SIDE ]->GetUVPtStruct();
2168 double rBot = ( uv_eb.size() - 1 ) * uvPt.normParam;
2169 double rTop = ( uv_et.size() - 1 ) * uvPt.normParam;
2170 int iBot = int( rBot );
2171 int iTop = int( rTop );
2172 double xBot = uv_eb[ iBot ].normParam + ( rBot - iBot ) * ( uv_eb[ iBot+1 ].normParam - uv_eb[ iBot ].normParam );
2173 double xTop = uv_et[ iTop ].normParam + ( rTop - iTop ) * ( uv_et[ iTop+1 ].normParam - uv_et[ iTop ].normParam );
2174 double x = xBot + y * ( xTop - xBot );
2176 gp_UV uv = CalcUV(/*x,y=*/x, y,
2177 /*a0,...=*/UVs[UV_A0], UVs[UV_A1], UVs[UV_A2], UVs[UV_A3],
2178 /*p0=*/quad->side[BOTTOM_SIDE]->Value2d( x ).XY(),
2180 /*p2=*/quad->side[TOP_SIDE ]->Value2d( x ).XY(),
2181 /*p3=*/UVs[ UV_L ]);
2182 gp_Pnt P = S->Value( uv.X(), uv.Y() );
2185 return helper->AddNode(P.X(), P.Y(), P.Z(), 0, uv.X(), uv.Y() );
2188 void reduce42( const vector<UVPtStruct>& curr_base,
2189 vector<UVPtStruct>& next_base,
2191 int & next_base_len,
2192 FaceQuadStruct* quad,
2195 SMESH_MesherHelper* helper,
2196 Handle(Geom_Surface)& S)
2198 // add one "HH": nodes a,b,c,d,e and faces 1,2,3,4,5,6
2200 // .-----a-----b i + 1
2211 const SMDS_MeshNode*& Na = next_base[ ++next_base_len ].node;
2213 Na = makeNode( next_base[ next_base_len ], y, quad, UVs, helper, S );
2216 const SMDS_MeshNode*& Nb = next_base[ ++next_base_len ].node;
2218 Nb = makeNode( next_base[ next_base_len ], y, quad, UVs, helper, S );
2221 double u = (curr_base[j + 2].u + next_base[next_base_len - 2].u) / 2.0;
2222 double v = (curr_base[j + 2].v + next_base[next_base_len - 2].v) / 2.0;
2223 gp_Pnt P = S->Value(u,v);
2224 SMDS_MeshNode* Nc = helper->AddNode(P.X(), P.Y(), P.Z(), 0, u, v);
2227 u = (curr_base[j + 2].u + next_base[next_base_len - 1].u) / 2.0;
2228 v = (curr_base[j + 2].v + next_base[next_base_len - 1].v) / 2.0;
2230 SMDS_MeshNode* Nd = helper->AddNode(P.X(), P.Y(), P.Z(), 0, u, v);
2233 u = (curr_base[j + 2].u + next_base[next_base_len].u) / 2.0;
2234 v = (curr_base[j + 2].v + next_base[next_base_len].v) / 2.0;
2236 SMDS_MeshNode* Ne = helper->AddNode(P.X(), P.Y(), P.Z(), 0, u, v);
2239 helper->AddFace(curr_base[j + 0].node,
2240 curr_base[j + 1].node, Nc,
2241 next_base[next_base_len - 2].node);
2243 helper->AddFace(curr_base[j + 1].node,
2244 curr_base[j + 2].node, Nd, Nc);
2246 helper->AddFace(curr_base[j + 2].node,
2247 curr_base[j + 3].node, Ne, Nd);
2249 helper->AddFace(curr_base[j + 3].node,
2250 curr_base[j + 4].node, Nb, Ne);
2252 helper->AddFace(Nc, Nd, Na, next_base[next_base_len - 2].node);
2254 helper->AddFace(Nd, Ne, Nb, Na);
2257 void reduce31( const vector<UVPtStruct>& curr_base,
2258 vector<UVPtStruct>& next_base,
2260 int & next_base_len,
2261 FaceQuadStruct* quad,
2264 SMESH_MesherHelper* helper,
2265 Handle(Geom_Surface)& S)
2267 // add one "H": nodes b,c,e and faces 1,2,4,5
2269 // .---------b i + 1
2280 const SMDS_MeshNode*& Nb = next_base[ ++next_base_len ].node;
2282 Nb = makeNode( next_base[ next_base_len ], y, quad, UVs, helper, S );
2285 double u1 = (curr_base[ j ].u + next_base[ next_base_len-1 ].u ) / 2.0;
2286 double u2 = (curr_base[ j+3 ].u + next_base[ next_base_len ].u ) / 2.0;
2287 double u3 = (u2 - u1) / 3.0;
2289 double v1 = (curr_base[ j ].v + next_base[ next_base_len-1 ].v ) / 2.0;
2290 double v2 = (curr_base[ j+3 ].v + next_base[ next_base_len ].v ) / 2.0;
2291 double v3 = (v2 - v1) / 3.0;
2295 gp_Pnt P = S->Value(u,v);
2296 SMDS_MeshNode* Nc = helper->AddNode( P.X(), P.Y(), P.Z(), 0, u, v );
2301 SMDS_MeshNode* Ne = helper->AddNode( P.X(), P.Y(), P.Z(), 0, u, v );
2305 helper->AddFace( curr_base[ j + 0 ].node,
2306 curr_base[ j + 1 ].node,
2308 next_base[ next_base_len - 1 ].node);
2310 helper->AddFace( curr_base[ j + 1 ].node,
2311 curr_base[ j + 2 ].node, Ne, Nc);
2313 helper->AddFace( curr_base[ j + 2 ].node,
2314 curr_base[ j + 3 ].node, Nb, Ne);
2316 helper->AddFace(Nc, Ne, Nb,
2317 next_base[ next_base_len - 1 ].node);
2320 typedef void (* PReduceFunction) ( const vector<UVPtStruct>& curr_base,
2321 vector<UVPtStruct>& next_base,
2323 int & next_base_len,
2324 FaceQuadStruct* quad,
2327 SMESH_MesherHelper* helper,
2328 Handle(Geom_Surface)& S);
2332 //=======================================================================
2334 * Implementation of Reduced algorithm (meshing with quadrangles only)
2336 //=======================================================================
2337 bool StdMeshers_Quadrangle_2D::ComputeReduced (SMESH_Mesh & aMesh,
2338 const TopoDS_Shape& aShape,
2339 FaceQuadStruct* quad)
2341 SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
2342 const TopoDS_Face& F = TopoDS::Face(aShape);
2343 Handle(Geom_Surface) S = BRep_Tool::Surface(F);
2344 int i,j,geomFaceID = meshDS->ShapeToIndex(F);
2346 int nb = quad->side[0]->NbPoints();
2347 int nr = quad->side[1]->NbPoints();
2348 int nt = quad->side[2]->NbPoints();
2349 int nl = quad->side[3]->NbPoints();
2351 // Simple Reduce 10->8->6->4 (3 steps) Multiple Reduce 10->4 (1 step)
2353 // .-----.-----.-----.-----. .-----.-----.-----.-----.
2354 // | / \ | / \ | | / \ | / \ |
2355 // | / .--.--. \ | | / \ | / \ |
2356 // | / / | \ \ | | / .----.----. \ |
2357 // .---.---.---.---.---.---. | / / \ | / \ \ |
2358 // | / / \ | / \ \ | | / / \ | / \ \ |
2359 // | / / .-.-. \ \ | | / / .---.---. \ \ |
2360 // | / / / | \ \ \ | | / / / \ | / \ \ \ |
2361 // .--.--.--.--.--.--.--.--. | / / / \ | / \ \ \ |
2362 // | / / / \ | / \ \ \ | | / / / .-.-. \ \ \ |
2363 // | / / / .-.-. \ \ \ | | / / / / | \ \ \ \ |
2364 // | / / / / | \ \ \ \ | | / / / / | \ \ \ \ |
2365 // .-.-.-.--.--.--.--.-.-.-. .-.-.-.--.--.--.--.-.-.-.
2367 bool MultipleReduce = false;
2379 else if (nb == nt) {
2380 nr1 = nb; // and == nt
2394 // number of rows and columns
2395 int nrows = nr1 - 1;
2396 int ncol_top = nt1 - 1;
2397 int ncol_bot = nb1 - 1;
2398 // number of rows needed to reduce ncol_bot to ncol_top using simple 3->1 "tree" (see below)
2399 int nrows_tree31 = int( log( (double)(ncol_bot / ncol_top) ) / log((double) 3 )); // = log x base 3
2400 if ( nrows < nrows_tree31 )
2401 MultipleReduce = true;
2404 if (MultipleReduce) { // == ComputeQuadPref QUAD_QUADRANGLE_PREF_REVERSED
2405 //==================================================
2406 int dh = abs(nb-nt);
2407 int dv = abs(nr-nl);
2411 // it is a base case => not shift quad but may be replacement is need
2412 ShiftQuad(quad,0,true);
2415 // we have to shift quad on 2
2416 ShiftQuad(quad,2,true);
2421 // we have to shift quad on 1
2422 ShiftQuad(quad,1,true);
2425 // we have to shift quad on 3
2426 ShiftQuad(quad,3,true);
2430 nb = quad->side[0]->NbPoints();
2431 nr = quad->side[1]->NbPoints();
2432 nt = quad->side[2]->NbPoints();
2433 nl = quad->side[3]->NbPoints();
2436 int nbh = Max(nb,nt);
2437 int nbv = Max(nr,nl);
2450 const vector<UVPtStruct>& uv_eb = quad->side[0]->GetUVPtStruct(true,0);
2451 const vector<UVPtStruct>& uv_er = quad->side[1]->GetUVPtStruct(false,1);
2452 const vector<UVPtStruct>& uv_et = quad->side[2]->GetUVPtStruct(true,1);
2453 const vector<UVPtStruct>& uv_el = quad->side[3]->GetUVPtStruct(false,0);
2455 if (uv_eb.size() != nb || uv_er.size() != nr || uv_et.size() != nt || uv_el.size() != nl)
2456 return error(COMPERR_BAD_INPUT_MESH);
2459 UpdateDegenUV( quad );
2461 // arrays for normalized params
2462 TColStd_SequenceOfReal npb, npr, npt, npl;
2463 for (j = 0; j < nb; j++) {
2464 npb.Append(uv_eb[j].normParam);
2466 for (i = 0; i < nr; i++) {
2467 npr.Append(uv_er[i].normParam);
2469 for (j = 0; j < nt; j++) {
2470 npt.Append(uv_et[j].normParam);
2472 for (i = 0; i < nl; i++) {
2473 npl.Append(uv_el[i].normParam);
2477 // orientation of face and 3 main domain for future faces
2483 // left | | | | rigth
2490 // add some params to right and left after the first param
2493 double dpr = (npr.Value(2) - npr.Value(1))/(dr+1);
2494 for (i=1; i<=dr; i++) {
2495 npr.InsertAfter(1,npr.Value(2)-dpr);
2499 dpr = (npl.Value(2) - npl.Value(1))/(dl+1);
2500 for (i=1; i<=dl; i++) {
2501 npl.InsertAfter(1,npl.Value(2)-dpr);
2504 gp_XY a0 (uv_eb.front().u, uv_eb.front().v);
2505 gp_XY a1 (uv_eb.back().u, uv_eb.back().v);
2506 gp_XY a2 (uv_et.back().u, uv_et.back().v);
2507 gp_XY a3 (uv_et.front().u, uv_et.front().v);
2509 int nnn = Min(nr,nl);
2510 // auxilary sequence of XY for creation of nodes
2511 // in the bottom part of central domain
2512 // it's length must be == nbv-nnn-1
2513 TColgp_SequenceOfXY UVL;
2514 TColgp_SequenceOfXY UVR;
2515 //==================================================
2517 // step1: create faces for left domain
2518 StdMeshers_Array2OfNode NodesL(1,dl+1,1,nl);
2520 for (j=1; j<=nl; j++)
2521 NodesL.SetValue(1,j,uv_el[j-1].node);
2524 for (i=1; i<=dl; i++)
2525 NodesL.SetValue(i+1,nl,uv_et[i].node);
2526 // create and add needed nodes
2527 TColgp_SequenceOfXY UVtmp;
2528 for (i=1; i<=dl; i++) {
2529 double x0 = npt.Value(i+1);
2532 double y0 = npl.Value(i+1);
2533 double y1 = npr.Value(i+1);
2534 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
2535 gp_Pnt P = S->Value(UV.X(),UV.Y());
2536 SMDS_MeshNode * N = meshDS->AddNode(P.X(), P.Y(), P.Z());
2537 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
2538 NodesL.SetValue(i+1,1,N);
2539 if (UVL.Length()<nbv-nnn-1) UVL.Append(UV);
2541 for (j=2; j<nl; j++) {
2542 double y0 = npl.Value(dl+j);
2543 double y1 = npr.Value(dl+j);
2544 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
2545 gp_Pnt P = S->Value(UV.X(),UV.Y());
2546 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
2547 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
2548 NodesL.SetValue(i+1,j,N);
2549 if (i==dl) UVtmp.Append(UV);
2552 for (i=1; i<=UVtmp.Length() && UVL.Length()<nbv-nnn-1; i++) {
2553 UVL.Append(UVtmp.Value(i));
2556 for (i=1; i<=dl; i++) {
2557 for (j=1; j<nl; j++) {
2559 myHelper->AddFace(NodesL.Value(i,j), NodesL.Value(i+1,j),
2560 NodesL.Value(i+1,j+1), NodesL.Value(i,j+1));
2561 if (F) meshDS->SetMeshElementOnShape(F, geomFaceID);
2566 // fill UVL using c2d
2567 for (i=1; i<npl.Length() && UVL.Length()<nbv-nnn-1; i++) {
2568 UVL.Append(gp_UV (uv_el[i].u, uv_el[i].v));
2572 // step2: create faces for right domain
2573 StdMeshers_Array2OfNode NodesR(1,dr+1,1,nr);
2575 for (j=1; j<=nr; j++)
2576 NodesR.SetValue(1,j,uv_er[nr-j].node);
2579 for (i=1; i<=dr; i++)
2580 NodesR.SetValue(i+1,1,uv_et[nt-1-i].node);
2581 // create and add needed nodes
2582 TColgp_SequenceOfXY UVtmp;
2583 for (i=1; i<=dr; i++) {
2584 double x0 = npt.Value(nt-i);
2587 double y0 = npl.Value(i+1);
2588 double y1 = npr.Value(i+1);
2589 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
2590 gp_Pnt P = S->Value(UV.X(),UV.Y());
2591 SMDS_MeshNode * N = meshDS->AddNode(P.X(), P.Y(), P.Z());
2592 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
2593 NodesR.SetValue(i+1,nr,N);
2594 if (UVR.Length()<nbv-nnn-1) UVR.Append(UV);
2596 for (j=2; j<nr; j++) {
2597 double y0 = npl.Value(nbv-j+1);
2598 double y1 = npr.Value(nbv-j+1);
2599 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
2600 gp_Pnt P = S->Value(UV.X(),UV.Y());
2601 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
2602 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
2603 NodesR.SetValue(i+1,j,N);
2604 if (i==dr) UVtmp.Prepend(UV);
2607 for (i=1; i<=UVtmp.Length() && UVR.Length()<nbv-nnn-1; i++) {
2608 UVR.Append(UVtmp.Value(i));
2611 for (i=1; i<=dr; i++) {
2612 for (j=1; j<nr; j++) {
2614 myHelper->AddFace(NodesR.Value(i,j), NodesR.Value(i+1,j),
2615 NodesR.Value(i+1,j+1), NodesR.Value(i,j+1));
2616 if (F) meshDS->SetMeshElementOnShape(F, geomFaceID);
2621 // fill UVR using c2d
2622 for (i=1; i<npr.Length() && UVR.Length()<nbv-nnn-1; i++) {
2623 UVR.Append(gp_UV(uv_er[i].u, uv_er[i].v));
2627 // step3: create faces for central domain
2628 StdMeshers_Array2OfNode NodesC(1,nb,1,nbv);
2629 // add first line using NodesL
2630 for (i=1; i<=dl+1; i++)
2631 NodesC.SetValue(1,i,NodesL(i,1));
2632 for (i=2; i<=nl; i++)
2633 NodesC.SetValue(1,dl+i,NodesL(dl+1,i));
2634 // add last line using NodesR
2635 for (i=1; i<=dr+1; i++)
2636 NodesC.SetValue(nb,i,NodesR(i,nr));
2637 for (i=1; i<nr; i++)
2638 NodesC.SetValue(nb,dr+i+1,NodesR(dr+1,nr-i));
2639 // add top nodes (last columns)
2640 for (i=dl+2; i<nbh-dr; i++)
2641 NodesC.SetValue(i-dl,nbv,uv_et[i-1].node);
2642 // add bottom nodes (first columns)
2643 for (i=2; i<nb; i++)
2644 NodesC.SetValue(i,1,uv_eb[i-1].node);
2646 // create and add needed nodes
2647 // add linear layers
2648 for (i=2; i<nb; i++) {
2649 double x0 = npt.Value(dl+i);
2651 for (j=1; j<nnn; j++) {
2652 double y0 = npl.Value(nbv-nnn+j);
2653 double y1 = npr.Value(nbv-nnn+j);
2654 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
2655 gp_Pnt P = S->Value(UV.X(),UV.Y());
2656 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
2657 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
2658 NodesC.SetValue(i,nbv-nnn+j,N);
2661 // add diagonal layers
2662 for (i=1; i<nbv-nnn; i++) {
2663 double du = UVR.Value(i).X() - UVL.Value(i).X();
2664 double dv = UVR.Value(i).Y() - UVL.Value(i).Y();
2665 for (j=2; j<nb; j++) {
2666 double u = UVL.Value(i).X() + du*npb.Value(j);
2667 double v = UVL.Value(i).Y() + dv*npb.Value(j);
2668 gp_Pnt P = S->Value(u,v);
2669 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
2670 meshDS->SetNodeOnFace(N, geomFaceID, u, v);
2671 NodesC.SetValue(j,i+1,N);
2675 for (i=1; i<nb; i++) {
2676 for (j=1; j<nbv; j++) {
2678 myHelper->AddFace(NodesC.Value(i,j), NodesC.Value(i+1,j),
2679 NodesC.Value(i+1,j+1), NodesC.Value(i,j+1));
2680 if (F) meshDS->SetMeshElementOnShape(F, geomFaceID);
2684 } // end Multiple Reduce implementation
2685 else { // Simple Reduce (!MultipleReduce)
2686 //=========================================================
2689 // it is a base case => not shift quad
2690 //ShiftQuad(quad,0,true);
2693 // we have to shift quad on 2
2694 ShiftQuad(quad,2,true);
2699 // we have to shift quad on 1
2700 ShiftQuad(quad,1,true);
2703 // we have to shift quad on 3
2704 ShiftQuad(quad,3,true);
2708 nb = quad->side[0]->NbPoints();
2709 nr = quad->side[1]->NbPoints();
2710 nt = quad->side[2]->NbPoints();
2711 nl = quad->side[3]->NbPoints();
2713 // number of rows and columns
2714 int nrows = nr - 1; // and also == nl - 1
2715 int ncol_top = nt - 1;
2716 int ncol_bot = nb - 1;
2717 int npair_top = ncol_top / 2;
2718 // maximum number of bottom elements for "linear" simple reduce 4->2
2719 int max_lin42 = ncol_top + npair_top * 2 * nrows;
2720 // maximum number of bottom elements for "linear" simple reduce 3->1
2721 int max_lin31 = ncol_top + ncol_top * 2 * nrows;
2722 // maximum number of bottom elements for "tree" simple reduce 4->2
2724 // number of rows needed to reduce ncol_bot to ncol_top using simple 4->2 "tree"
2725 int nrows_tree42 = int( log( (double)(ncol_bot / ncol_top) )/log((double)2) ); // needed to avoid overflow at pow(2) while computing max_tree42
2726 if (nrows_tree42 < nrows) {
2727 max_tree42 = npair_top * pow(2.0, nrows + 1);
2728 if ( ncol_top > npair_top * 2 ) {
2729 int delta = ncol_bot - max_tree42;
2730 for (int irow = 1; irow < nrows; irow++) {
2731 int nfour = delta / 4;
2734 if (delta <= (ncol_top - npair_top * 2))
2735 max_tree42 = ncol_bot;
2738 // maximum number of bottom elements for "tree" simple reduce 3->1
2739 //int max_tree31 = ncol_top * pow(3.0, nrows);
2740 bool is_lin_31 = false;
2741 bool is_lin_42 = false;
2742 bool is_tree_31 = false;
2743 bool is_tree_42 = false;
2744 int max_lin = max_lin42;
2745 if (ncol_bot > max_lin42) {
2746 if (ncol_bot <= max_lin31) {
2748 max_lin = max_lin31;
2752 // if ncol_bot is a 3*n or not 2*n
2753 if ((ncol_bot/3)*3 == ncol_bot || (ncol_bot/2)*2 != ncol_bot) {
2755 max_lin = max_lin31;
2761 if (ncol_bot > max_lin) { // not "linear"
2762 is_tree_31 = (ncol_bot > max_tree42);
2763 if (ncol_bot <= max_tree42) {
2764 if ((ncol_bot/3)*3 == ncol_bot || (ncol_bot/2)*2 != ncol_bot) {
2773 const vector<UVPtStruct>& uv_eb = quad->side[0]->GetUVPtStruct(true,0);
2774 const vector<UVPtStruct>& uv_er = quad->side[1]->GetUVPtStruct(false,1);
2775 const vector<UVPtStruct>& uv_et = quad->side[2]->GetUVPtStruct(true,1);
2776 const vector<UVPtStruct>& uv_el = quad->side[3]->GetUVPtStruct(false,0);
2778 if (uv_eb.size() != nb || uv_er.size() != nr || uv_et.size() != nt || uv_el.size() != nl)
2779 return error(COMPERR_BAD_INPUT_MESH);
2781 myHelper->SetElementsOnShape( true );
2783 gp_UV uv[ UV_SIZE ];
2784 uv[ UV_A0 ].SetCoord( uv_eb.front().u, uv_eb.front().v);
2785 uv[ UV_A1 ].SetCoord( uv_eb.back().u, uv_eb.back().v );
2786 uv[ UV_A2 ].SetCoord( uv_et.back().u, uv_et.back().v );
2787 uv[ UV_A3 ].SetCoord( uv_et.front().u, uv_et.front().v);
2789 vector<UVPtStruct> curr_base = uv_eb, next_base;
2791 UVPtStruct nullUVPtStruct; nullUVPtStruct.node = 0;
2793 int curr_base_len = nb;
2794 int next_base_len = 0;
2797 { // ------------------------------------------------------------------
2798 // New algorithm implemented by request of IPAL22856
2799 // "2D quadrangle mesher of reduced type works wrong"
2800 // http://bugtracker.opencascade.com/show_bug.cgi?id=22856
2802 // the algorithm is following: all reduces are centred in horizontal
2803 // direction and are distributed among all rows
2805 if (ncol_bot > max_tree42) {
2809 if ((ncol_top/3)*3 == ncol_top ) {
2817 const int col_top_size = is_lin_42 ? 2 : 1;
2818 const int col_base_size = is_lin_42 ? 4 : 3;
2820 // Compute nb of "columns" (like in "linear" simple reducing) in all rows
2822 vector<int> nb_col_by_row;
2824 int delta_all = nb - nt;
2825 int delta_one_col = nrows * 2;
2826 int nb_col = delta_all / delta_one_col;
2827 int remainder = delta_all - nb_col * delta_one_col;
2828 if (remainder > 0) {
2831 if ( nb_col * col_top_size >= nt ) // == "tree" reducing situation
2833 // top row is full (all elements reduced), add "columns" one by one
2834 // in rows below until all bottom elements are reduced
2835 nb_col = ( nt - 1 ) / col_top_size;
2836 nb_col_by_row.resize( nrows, nb_col );
2837 int nbrows_not_full = nrows - 1;
2838 int cur_top_size = nt - 1;
2839 remainder = delta_all - nb_col * delta_one_col;
2840 while ( remainder > 0 )
2842 delta_one_col = nbrows_not_full * 2;
2843 int nb_col_add = remainder / delta_one_col;
2844 cur_top_size += 2 * nb_col_by_row[ nbrows_not_full ];
2845 int nb_col_free = cur_top_size / col_top_size - nb_col_by_row[ nbrows_not_full-1 ];
2846 if ( nb_col_add > nb_col_free )
2847 nb_col_add = nb_col_free;
2848 for ( int irow = 0; irow < nbrows_not_full; ++irow )
2849 nb_col_by_row[ irow ] += nb_col_add;
2851 remainder -= nb_col_add * delta_one_col;
2854 else // == "linear" reducing situation
2856 nb_col_by_row.resize( nrows, nb_col );
2858 for ( int irow = remainder / 2; irow < nrows; ++irow )
2859 nb_col_by_row[ irow ]--;
2864 PReduceFunction reduceFunction = & ( is_lin_42 ? reduce42 : reduce31 );
2866 const int reduce_grp_size = is_lin_42 ? 4 : 3;
2868 for (i = 1; i < nr; i++) // layer by layer
2870 nb_col = nb_col_by_row[ i-1 ];
2871 int nb_next = curr_base_len - nb_col * 2;
2872 if (nb_next < nt) nb_next = nt;
2874 const double y = uv_el[ i ].normParam;
2876 if ( i + 1 == nr ) // top
2883 next_base.resize( nb_next, nullUVPtStruct );
2884 next_base.front() = uv_el[i];
2885 next_base.back() = uv_er[i];
2887 // compute normalized param u
2888 double du = 1. / ( nb_next - 1 );
2889 next_base[0].normParam = 0.;
2890 for ( j = 1; j < nb_next; ++j )
2891 next_base[j].normParam = next_base[j-1].normParam + du;
2893 uv[ UV_L ].SetCoord( next_base.front().u, next_base.front().v );
2894 uv[ UV_R ].SetCoord( next_base.back().u, next_base.back().v );
2896 int free_left = ( curr_base_len - 1 - nb_col * col_base_size ) / 2;
2897 int free_middle = curr_base_len - 1 - nb_col * col_base_size - 2 * free_left;
2899 // not reduced left elements
2900 for (j = 0; j < free_left; j++)
2903 const SMDS_MeshNode*& Nf = next_base[++next_base_len].node;
2905 Nf = makeNode( next_base[ next_base_len ], y, quad, uv, myHelper, S );
2907 myHelper->AddFace(curr_base[ j ].node,
2908 curr_base[ j+1 ].node,
2910 next_base[ next_base_len-1 ].node);
2913 for (int icol = 1; icol <= nb_col; icol++)
2916 reduceFunction( curr_base, next_base, j, next_base_len, quad, uv, y, myHelper, S );
2918 j += reduce_grp_size;
2920 // elements in the middle of "columns" added for symmetry
2921 if ( free_middle > 0 && ( nb_col % 2 == 0 ) && icol == nb_col / 2 )
2923 for (int imiddle = 1; imiddle <= free_middle; imiddle++) {
2924 // f (i + 1, j + imiddle)
2925 const SMDS_MeshNode*& Nf = next_base[++next_base_len].node;
2927 Nf = makeNode( next_base[ next_base_len ], y, quad, uv, myHelper, S );
2929 myHelper->AddFace(curr_base[ j-1+imiddle ].node,
2930 curr_base[ j +imiddle ].node,
2932 next_base[ next_base_len-1 ].node);
2938 // not reduced right elements
2939 for (; j < curr_base_len-1; j++) {
2941 const SMDS_MeshNode*& Nf = next_base[++next_base_len].node;
2943 Nf = makeNode( next_base[ next_base_len ], y, quad, uv, myHelper, S );
2945 myHelper->AddFace(curr_base[ j ].node,
2946 curr_base[ j+1 ].node,
2948 next_base[ next_base_len-1 ].node);
2951 curr_base_len = next_base_len + 1;
2953 curr_base.swap( next_base );
2957 else if ( is_tree_42 || is_tree_31 )
2959 // "tree" simple reduce "42": 2->4->8->16->32->...
2961 // .-------------------------------.-------------------------------. nr
2963 // | \ .---------------.---------------. / |
2965 // .---------------.---------------.---------------.---------------.
2966 // | \ | / | \ | / |
2967 // | \ .-------.-------. / | \ .-------.-------. / |
2968 // | | | | | | | | |
2969 // .-------.-------.-------.-------.-------.-------.-------.-------. i
2970 // |\ | /|\ | /|\ | /|\ | /|
2971 // | \.---.---./ | \.---.---./ | \.---.---./ | \.---.---./ |
2972 // | | | | | | | | | | | | | | | | |
2973 // .---.---.---.---.---.---.---.---.---.---.---.---.---.---.---.---.
2974 // |\ | /|\ | /|\ | /|\ | /|\ | /|\ | /|\ | /|\ | /|
2975 // | .-.-. | .-.-. | .-.-. | .-.-. | .-.-. | .-.-. | .-.-. | .-.-. |
2976 // | | | | | | | | | | | | | | | | | | | | | | | | | | | | | | | | |
2977 // .-.-.-.-.-.-.-.-.-.-.-.-.-.-.-.-.-.-.-.-.-.-.-.-.-.-.-.-.-.-.-.-. 1
2980 // "tree" simple reduce "31": 1->3->9->27->...
2982 // .-----------------------------------------------------. nr
2984 // | .-----------------. |
2986 // .-----------------.-----------------.-----------------.
2987 // | \ / | \ / | \ / |
2988 // | .-----. | .-----. | .-----. | i
2989 // | | | | | | | | | |
2990 // .-----.-----.-----.-----.-----.-----.-----.-----.-----.
2991 // |\ /|\ /|\ /|\ /|\ /|\ /|\ /|\ /|\ /|
2992 // | .-. | .-. | .-. | .-. | .-. | .-. | .-. | .-. | .-. |
2993 // | | | | | | | | | | | | | | | | | | | | | | | | | | | |
2994 // .-.-.-.-.-.-.-.-.-.-.-.-.-.-.-.-.-.-.-.-.-.-.-.-.-.-.-. 1
2997 PReduceFunction reduceFunction = & ( is_tree_42 ? reduce42 : reduce31 );
2999 const int reduce_grp_size = is_tree_42 ? 4 : 3;
3001 for (i = 1; i < nr; i++) // layer by layer
3003 // to stop reducing, if number of nodes reaches nt
3004 int delta = curr_base_len - nt;
3006 // to calculate normalized parameter, we must know number of points in next layer
3007 int nb_reduce_groups = (curr_base_len - 1) / reduce_grp_size;
3008 int nb_next = nb_reduce_groups * (reduce_grp_size-2) + (curr_base_len - nb_reduce_groups*reduce_grp_size);
3009 if (nb_next < nt) nb_next = nt;
3011 const double y = uv_el[ i ].normParam;
3013 if ( i + 1 == nr ) // top
3020 next_base.resize( nb_next, nullUVPtStruct );
3021 next_base.front() = uv_el[i];
3022 next_base.back() = uv_er[i];
3024 // compute normalized param u
3025 double du = 1. / ( nb_next - 1 );
3026 next_base[0].normParam = 0.;
3027 for ( j = 1; j < nb_next; ++j )
3028 next_base[j].normParam = next_base[j-1].normParam + du;
3030 uv[ UV_L ].SetCoord( next_base.front().u, next_base.front().v );
3031 uv[ UV_R ].SetCoord( next_base.back().u, next_base.back().v );
3033 for (j = 0; j+reduce_grp_size < curr_base_len && delta > 0; j+=reduce_grp_size, delta-=2)
3035 reduceFunction( curr_base, next_base, j, next_base_len, quad, uv, y, myHelper, S );
3038 // not reduced side elements (if any)
3039 for (; j < curr_base_len-1; j++)
3042 const SMDS_MeshNode*& Nf = next_base[++next_base_len].node;
3044 Nf = makeNode( next_base[ next_base_len ], y, quad, uv, myHelper, S );
3046 myHelper->AddFace(curr_base[ j ].node,
3047 curr_base[ j+1 ].node,
3049 next_base[ next_base_len-1 ].node);
3051 curr_base_len = next_base_len + 1;
3053 curr_base.swap( next_base );
3055 } // end "tree" simple reduce
3057 else if ( is_lin_42 || is_lin_31 ) {
3058 // "linear" simple reduce "31": 2->6->10->14
3060 // .-----------------------------.-----------------------------. nr
3062 // | .---------. | .---------. |
3064 // .---------.---------.---------.---------.---------.---------.
3065 // | / \ / \ | / \ / \ |
3066 // | / .-----. \ | / .-----. \ | i
3067 // | / | | \ | / | | \ |
3068 // .-----.-----.-----.-----.-----.-----.-----.-----.-----.-----.
3069 // | / / \ / \ \ | / / \ / \ \ |
3070 // | / / .-. \ \ | / / .-. \ \ |
3071 // | / / / \ \ \ | / / / \ \ \ |
3072 // .--.----.---.-----.---.-----.-.--.----.---.-----.---.-----.-. 1
3075 // "linear" simple reduce "42": 4->8->12->16
3077 // .---------------.---------------.---------------.---------------. nr
3078 // | \ | / | \ | / |
3079 // | \ .-------.-------. / | \ .-------.-------. / |
3080 // | | | | | | | | |
3081 // .-------.-------.-------.-------.-------.-------.-------.-------.
3082 // | / \ | / \ | / \ | / \ |
3083 // | / \.----.----./ \ | / \.----.----./ \ | i
3084 // | / | | | \ | / | | | \ |
3085 // .-----.----.----.----.----.-----.-----.----.----.----.----.-----.
3086 // | / / \ | / \ \ | / / \ | / \ \ |
3087 // | / / .-.-. \ \ | / / .-.-. \ \ |
3088 // | / / / | \ \ \ | / / / | \ \ \ |
3089 // .---.---.---.---.---.---.---.---.---.---.---.---.---.---.---.---. 1
3092 // nt = 5, nb = 7, nr = 4
3093 //int delta_all = 2;
3094 //int delta_one_col = 6;
3096 //int remainder = 2;
3097 //if (remainder > 0) nb_col++;
3099 //int free_left = 1;
3101 //int free_middle = 4;
3103 int delta_all = nb - nt;
3104 int delta_one_col = (nr - 1) * 2;
3105 int nb_col = delta_all / delta_one_col;
3106 int remainder = delta_all - nb_col * delta_one_col;
3107 if (remainder > 0) {
3110 const int col_top_size = is_lin_42 ? 2 : 1;
3111 int free_left = ((nt - 1) - nb_col * col_top_size) / 2;
3112 free_left += nr - 2;
3113 int free_middle = (nr - 2) * 2;
3114 if (remainder > 0 && nb_col == 1) {
3115 int nb_rows_short_col = remainder / 2;
3116 int nb_rows_thrown = (nr - 1) - nb_rows_short_col;
3117 free_left -= nb_rows_thrown;
3120 // nt = 5, nb = 17, nr = 4
3121 //int delta_all = 12;
3122 //int delta_one_col = 6;
3124 //int remainder = 0;
3125 //int free_left = 2;
3126 //int free_middle = 4;
3128 PReduceFunction reduceFunction = & ( is_lin_42 ? reduce42 : reduce31 );
3130 const int reduce_grp_size = is_lin_42 ? 4 : 3;
3132 for (i = 1; i < nr; i++, free_middle -= 2, free_left -= 1) // layer by layer
3134 // to calculate normalized parameter, we must know number of points in next layer
3135 int nb_next = curr_base_len - nb_col * 2;
3136 if (remainder > 0 && i > remainder / 2)
3137 // take into account short "column"
3139 if (nb_next < nt) nb_next = nt;
3141 const double y = uv_el[ i ].normParam;
3143 if ( i + 1 == nr ) // top
3150 next_base.resize( nb_next, nullUVPtStruct );
3151 next_base.front() = uv_el[i];
3152 next_base.back() = uv_er[i];
3154 // compute normalized param u
3155 double du = 1. / ( nb_next - 1 );
3156 next_base[0].normParam = 0.;
3157 for ( j = 1; j < nb_next; ++j )
3158 next_base[j].normParam = next_base[j-1].normParam + du;
3160 uv[ UV_L ].SetCoord( next_base.front().u, next_base.front().v );
3161 uv[ UV_R ].SetCoord( next_base.back().u, next_base.back().v );
3163 // not reduced left elements
3164 for (j = 0; j < free_left; j++)
3167 const SMDS_MeshNode*& Nf = next_base[++next_base_len].node;
3169 Nf = makeNode( next_base[ next_base_len ], y, quad, uv, myHelper, S );
3171 myHelper->AddFace(curr_base[ j ].node,
3172 curr_base[ j+1 ].node,
3174 next_base[ next_base_len-1 ].node);
3177 for (int icol = 1; icol <= nb_col; icol++) {
3179 if (remainder > 0 && icol == nb_col && i > remainder / 2)
3180 // stop short "column"
3184 reduceFunction( curr_base, next_base, j, next_base_len, quad, uv, y, myHelper, S );
3186 j += reduce_grp_size;
3188 // not reduced middle elements
3189 if (icol < nb_col) {
3190 if (remainder > 0 && icol == nb_col - 1 && i > remainder / 2)
3191 // pass middle elements before stopped short "column"
3194 int free_add = free_middle;
3195 if (remainder > 0 && icol == nb_col - 1)
3196 // next "column" is short
3197 free_add -= (nr - 1) - (remainder / 2);
3199 for (int imiddle = 1; imiddle <= free_add; imiddle++) {
3200 // f (i + 1, j + imiddle)
3201 const SMDS_MeshNode*& Nf = next_base[++next_base_len].node;
3203 Nf = makeNode( next_base[ next_base_len ], y, quad, uv, myHelper, S );
3205 myHelper->AddFace(curr_base[ j-1+imiddle ].node,
3206 curr_base[ j +imiddle ].node,
3208 next_base[ next_base_len-1 ].node);
3214 // not reduced right elements
3215 for (; j < curr_base_len-1; j++) {
3217 const SMDS_MeshNode*& Nf = next_base[++next_base_len].node;
3219 Nf = makeNode( next_base[ next_base_len ], y, quad, uv, myHelper, S );
3221 myHelper->AddFace(curr_base[ j ].node,
3222 curr_base[ j+1 ].node,
3224 next_base[ next_base_len-1 ].node);
3227 curr_base_len = next_base_len + 1;
3229 curr_base.swap( next_base );
3232 } // end "linear" simple reduce
3237 } // end Simple Reduce implementation
3243 //================================================================================
3244 namespace // data for smoothing
3247 // --------------------------------------------------------------------------------
3249 * \brief Structure used to check validity of node position after smoothing.
3250 * It holds two nodes connected to a smoothed node and belonging to
3257 TTriangle( TSmoothNode* n1=0, TSmoothNode* n2=0 ): _n1(n1), _n2(n2) {}
3259 inline bool IsForward( gp_UV uv ) const;
3261 // --------------------------------------------------------------------------------
3263 * \brief Data of a smoothed node
3268 vector< TTriangle > _triangles; // if empty, then node is not movable
3270 // --------------------------------------------------------------------------------
3271 inline bool TTriangle::IsForward( gp_UV uv ) const
3273 gp_Vec2d v1( uv, _n1->_uv ), v2( uv, _n2->_uv );
3279 //================================================================================
3281 * \brief Set UV of nodes on degenerated VERTEXes in the middle of degenerated EDGE
3283 * WARNING: this method must be called AFTER retrieving UVPtStruct's from quad
3285 //================================================================================
3287 void StdMeshers_Quadrangle_2D::UpdateDegenUV(FaceQuadStruct* quad)
3289 for ( unsigned i = 0; i < quad->side.size(); ++i )
3291 StdMeshers_FaceSide* side = quad->side[i];
3292 const vector<UVPtStruct>& uvVec = side->GetUVPtStruct();
3294 // find which end of the side is on degenerated shape
3296 if ( myHelper->IsDegenShape( uvVec[0].node->getshapeId() ))
3298 else if ( myHelper->IsDegenShape( uvVec.back().node->getshapeId() ))
3299 degenInd = uvVec.size() - 1;
3303 // find another side sharing the degenerated shape
3304 bool isPrev = ( degenInd == 0 );
3305 if ( i >= TOP_SIDE )
3307 int i2 = ( isPrev ? ( i + 3 ) : ( i + 1 )) % 4;
3308 StdMeshers_FaceSide* side2 = quad->side[ i2 ];
3309 const vector<UVPtStruct>& uvVec2 = side2->GetUVPtStruct();
3311 if ( uvVec[ degenInd ].node == uvVec2[0].node )
3313 else if ( uvVec[ degenInd ].node == uvVec2.back().node )
3314 degenInd2 = uvVec2.size() - 1;
3316 throw SALOME_Exception( LOCALIZED( "Logical error" ));
3318 // move UV in the middle
3319 uvPtStruct& uv1 = const_cast<uvPtStruct&>( uvVec [ degenInd ]);
3320 uvPtStruct& uv2 = const_cast<uvPtStruct&>( uvVec2[ degenInd2 ]);
3321 uv1.u = uv2.u = 0.5 * ( uv1.u + uv2.u );
3322 uv1.v = uv2.v = 0.5 * ( uv1.v + uv2.v );
3326 //================================================================================
3328 * \brief Perform smoothing of 2D elements on a FACE with ignored degenerated EDGE
3330 //================================================================================
3332 void StdMeshers_Quadrangle_2D::Smooth (FaceQuadStruct* quad)
3334 if ( !myNeedSmooth ) return;
3336 // Get nodes to smooth
3338 typedef map< const SMDS_MeshNode*, TSmoothNode, TIDCompare > TNo2SmooNoMap;
3339 TNo2SmooNoMap smooNoMap;
3341 const TopoDS_Face& geomFace = TopoDS::Face( myHelper->GetSubShape() );
3342 SMESHDS_Mesh* meshDS = myHelper->GetMeshDS();
3343 SMESHDS_SubMesh* fSubMesh = meshDS->MeshElements( geomFace );
3344 SMDS_NodeIteratorPtr nIt = fSubMesh->GetNodes();
3345 while ( nIt->more() ) // loop on nodes bound to a FACE
3347 const SMDS_MeshNode* node = nIt->next();
3348 TSmoothNode & sNode = smooNoMap[ node ];
3349 sNode._uv = myHelper->GetNodeUV( geomFace, node );
3351 // set sNode._triangles
3352 SMDS_ElemIteratorPtr fIt = node->GetInverseElementIterator( SMDSAbs_Face );
3353 while ( fIt->more() )
3355 const SMDS_MeshElement* face = fIt->next();
3356 const int nbN = face->NbCornerNodes();
3357 const int nInd = face->GetNodeIndex( node );
3358 const int prevInd = myHelper->WrapIndex( nInd - 1, nbN );
3359 const int nextInd = myHelper->WrapIndex( nInd + 1, nbN );
3360 const SMDS_MeshNode* prevNode = face->GetNode( prevInd );
3361 const SMDS_MeshNode* nextNode = face->GetNode( nextInd );
3362 sNode._triangles.push_back( TTriangle( & smooNoMap[ prevNode ],
3363 & smooNoMap[ nextNode ]));
3366 // set _uv of smooth nodes on FACE boundary
3367 for ( unsigned i = 0; i < quad->side.size(); ++i )
3369 const vector<UVPtStruct>& uvVec = quad->side[i]->GetUVPtStruct();
3370 for ( unsigned j = 0; j < uvVec.size(); ++j )
3372 TSmoothNode & sNode = smooNoMap[ uvVec[j].node ];
3373 sNode._uv.SetCoord( uvVec[j].u, uvVec[j].v );
3377 // define refernce orientation in 2D
3378 TNo2SmooNoMap::iterator n2sn = smooNoMap.begin();
3379 for ( ; n2sn != smooNoMap.end(); ++n2sn )
3380 if ( !n2sn->second._triangles.empty() )
3382 if ( n2sn == smooNoMap.end() ) return;
3383 const TSmoothNode & sampleNode = n2sn->second;
3384 const bool refForward = ( sampleNode._triangles[0].IsForward( sampleNode._uv ));
3388 for ( int iLoop = 0; iLoop < 5; ++iLoop )
3390 for ( n2sn = smooNoMap.begin(); n2sn != smooNoMap.end(); ++n2sn )
3392 TSmoothNode& sNode = n2sn->second;
3393 if ( sNode._triangles.empty() )
3394 continue; // not movable node
3398 for ( unsigned i = 0; i < sNode._triangles.size(); ++i )
3399 newUV += sNode._triangles[i]._n1->_uv;
3400 newUV /= sNode._triangles.size();
3402 // check validity of the newUV
3403 bool isValid = true;
3404 for ( unsigned i = 0; i < sNode._triangles.size() && isValid; ++i )
3405 isValid = ( sNode._triangles[i].IsForward( newUV ) == refForward );
3412 // Set new XYZ to the smoothed nodes
3414 Handle(Geom_Surface) surface = BRep_Tool::Surface( geomFace );
3416 for ( n2sn = smooNoMap.begin(); n2sn != smooNoMap.end(); ++n2sn )
3418 TSmoothNode& sNode = n2sn->second;
3419 if ( sNode._triangles.empty() )
3420 continue; // not movable node
3422 SMDS_MeshNode* node = const_cast< SMDS_MeshNode*>( n2sn->first );
3423 gp_Pnt xyz = surface->Value( sNode._uv.X(), sNode._uv.Y() );
3424 meshDS->MoveNode( node, xyz.X(), xyz.Y(), xyz.Z() );
3427 node->SetPosition( SMDS_PositionPtr( new SMDS_FacePosition( sNode._uv.X(), sNode._uv.Y() )));
3430 // Move medium nodes in quadratic mesh
3431 if ( _quadraticMesh )
3433 const TLinkNodeMap& links = myHelper->GetTLinkNodeMap();
3434 TLinkNodeMap::const_iterator linkIt = links.begin();
3435 for ( ; linkIt != links.end(); ++linkIt )
3437 const SMESH_TLink& link = linkIt->first;
3438 SMDS_MeshNode* node = const_cast< SMDS_MeshNode*>( linkIt->second );
3440 if ( node->getshapeId() != myHelper->GetSubShapeID() )
3441 continue; // medium node is on EDGE or VERTEX
3443 gp_XY uv1 = myHelper->GetNodeUV( geomFace, link.node1(), node );
3444 gp_XY uv2 = myHelper->GetNodeUV( geomFace, link.node2(), node );
3446 gp_XY uv = myHelper->GetMiddleUV( surface, uv1, uv2 );
3447 node->SetPosition( SMDS_PositionPtr( new SMDS_FacePosition( uv.X(), uv.Y() )));
3449 gp_Pnt xyz = surface->Value( uv.X(), uv.Y() );
3450 meshDS->MoveNode( node, xyz.X(), xyz.Y(), xyz.Z() );