1 // Copyright (C) 2007-2008 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
22 // SMESH SMESH : implementaion of SMESH idl descriptions
23 // File : StdMeshers_Quadrangle_2D.cxx
24 // Moved here from SMESH_Quadrangle_2D.cxx
25 // Author : Paul RASCLE, EDF
29 #include "StdMeshers_Quadrangle_2D.hxx"
31 #include "StdMeshers_FaceSide.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 <BRepTools.hxx>
46 #include <BRepTools_WireExplorer.hxx>
47 #include <BRep_Tool.hxx>
48 #include <Geom_Surface.hxx>
49 #include <NCollection_DefineArray2.hxx>
50 #include <Precision.hxx>
51 #include <TColStd_SequenceOfReal.hxx>
52 #include <TColgp_SequenceOfXY.hxx>
56 #include "utilities.h"
57 #include "Utils_ExceptHandlers.hxx"
59 #ifndef StdMeshers_Array2OfNode_HeaderFile
60 #define StdMeshers_Array2OfNode_HeaderFile
61 typedef const SMDS_MeshNode* SMDS_MeshNodePtr;
62 DEFINE_BASECOLLECTION (StdMeshers_BaseCollectionNodePtr, SMDS_MeshNodePtr)
63 DEFINE_ARRAY2(StdMeshers_Array2OfNode,
64 StdMeshers_BaseCollectionNodePtr, SMDS_MeshNodePtr)
70 typedef SMESH_Comment TComm;
72 //=============================================================================
76 //=============================================================================
78 StdMeshers_Quadrangle_2D::StdMeshers_Quadrangle_2D (int hypId, int studyId, SMESH_Gen* gen)
79 : SMESH_2D_Algo(hypId, studyId, gen)
81 MESSAGE("StdMeshers_Quadrangle_2D::StdMeshers_Quadrangle_2D");
82 _name = "Quadrangle_2D";
83 _shapeType = (1 << TopAbs_FACE);
84 _compatibleHypothesis.push_back("QuadranglePreference");
85 _compatibleHypothesis.push_back("TrianglePreference");
89 //=============================================================================
93 //=============================================================================
95 StdMeshers_Quadrangle_2D::~StdMeshers_Quadrangle_2D()
97 MESSAGE("StdMeshers_Quadrangle_2D::~StdMeshers_Quadrangle_2D");
100 //=============================================================================
104 //=============================================================================
106 bool StdMeshers_Quadrangle_2D::CheckHypothesis
108 const TopoDS_Shape& aShape,
109 SMESH_Hypothesis::Hypothesis_Status& aStatus)
112 aStatus = SMESH_Hypothesis::HYP_OK;
115 const list <const SMESHDS_Hypothesis * >&hyps = GetUsedHypothesis(aMesh, aShape, false);
116 const SMESHDS_Hypothesis *theHyp = 0;
119 theHyp = *hyps.begin();
120 if(strcmp("QuadranglePreference", theHyp->GetName()) == 0) {
121 myQuadranglePreference= true;
122 myTrianglePreference= false;
124 else if(strcmp("TrianglePreference", theHyp->GetName()) == 0){
125 myQuadranglePreference= false;
126 myTrianglePreference= true;
130 myQuadranglePreference = false;
131 myTrianglePreference = false;
136 //=============================================================================
140 //=============================================================================
142 bool StdMeshers_Quadrangle_2D::Compute (SMESH_Mesh& aMesh,
143 const TopoDS_Shape& aShape)// throw (SALOME_Exception)
145 // PAL14921. Enable catching std::bad_alloc and Standard_OutOfMemory outside
146 //Unexpect aCatchSalomeException);
148 SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
149 aMesh.GetSubMesh(aShape);
151 SMESH_MesherHelper helper(aMesh);
154 _quadraticMesh = myTool->IsQuadraticSubMesh(aShape);
156 FaceQuadStruct *quad = CheckNbEdges( aMesh, aShape );
157 std::auto_ptr<FaceQuadStruct> quadDeleter( quad ); // to delete quad at exit from Compute()
161 if(myQuadranglePreference) {
162 int n1 = quad->side[0]->NbPoints();
163 int n2 = quad->side[1]->NbPoints();
164 int n3 = quad->side[2]->NbPoints();
165 int n4 = quad->side[3]->NbPoints();
166 int nfull = n1+n2+n3+n4;
169 if( nfull==ntmp && ( (n1!=n3) || (n2!=n4) ) ) {
170 // special path for using only quandrangle faces
171 bool ok = ComputeQuadPref(aMesh, aShape, quad);
176 // set normalized grid on unit square in parametric domain
178 if (!SetNormalizedGrid(aMesh, aShape, quad))
181 // --- compute 3D values on points, store points & quadrangles
183 int nbdown = quad->side[0]->NbPoints();
184 int nbup = quad->side[2]->NbPoints();
186 int nbright = quad->side[1]->NbPoints();
187 int nbleft = quad->side[3]->NbPoints();
189 int nbhoriz = Min(nbdown, nbup);
190 int nbvertic = Min(nbright, nbleft);
192 const TopoDS_Face& F = TopoDS::Face(aShape);
193 Handle(Geom_Surface) S = BRep_Tool::Surface(F);
195 // internal mesh nodes
196 int i, j, geomFaceID = meshDS->ShapeToIndex( F );
197 for (i = 1; i < nbhoriz - 1; i++) {
198 for (j = 1; j < nbvertic - 1; j++) {
199 int ij = j * nbhoriz + i;
200 double u = quad->uv_grid[ij].u;
201 double v = quad->uv_grid[ij].v;
202 gp_Pnt P = S->Value(u, v);
203 SMDS_MeshNode * node = meshDS->AddNode(P.X(), P.Y(), P.Z());
204 meshDS->SetNodeOnFace(node, geomFaceID, u, v);
205 quad->uv_grid[ij].node = node;
212 // --.--.--.--.--.-- nbvertic
218 // ---.----.----.--- 0
219 // 0 > > > > > > > > nbhoriz
225 int iup = nbhoriz - 1;
226 if (quad->isEdgeOut[3]) { ilow++; } else { if (quad->isEdgeOut[1]) iup--; }
229 int jup = nbvertic - 1;
230 if (quad->isEdgeOut[0]) { jlow++; } else { if (quad->isEdgeOut[2]) jup--; }
232 // regular quadrangles
233 for (i = ilow; i < iup; i++) {
234 for (j = jlow; j < jup; j++) {
235 const SMDS_MeshNode *a, *b, *c, *d;
236 a = quad->uv_grid[j * nbhoriz + i].node;
237 b = quad->uv_grid[j * nbhoriz + i + 1].node;
238 c = quad->uv_grid[(j + 1) * nbhoriz + i + 1].node;
239 d = quad->uv_grid[(j + 1) * nbhoriz + i].node;
240 SMDS_MeshFace* face = myTool->AddFace(a, b, c, d);
241 meshDS->SetMeshElementOnShape(face, geomFaceID);
245 const vector<UVPtStruct>& uv_e0 = quad->side[0]->GetUVPtStruct(true,0 );
246 const vector<UVPtStruct>& uv_e1 = quad->side[1]->GetUVPtStruct(false,1);
247 const vector<UVPtStruct>& uv_e2 = quad->side[2]->GetUVPtStruct(true,1 );
248 const vector<UVPtStruct>& uv_e3 = quad->side[3]->GetUVPtStruct(false,0);
250 if ( uv_e0.empty() || uv_e1.empty() || uv_e2.empty() || uv_e3.empty() )
251 return error( COMPERR_BAD_INPUT_MESH );
253 double eps = Precision::Confusion();
255 // Boundary quadrangles
257 if (quad->isEdgeOut[0]) {
260 // |___|___|___|___|___|___|
262 // |___|___|___|___|___|___|
264 // |___|___|___|___|___|___| __ first row of the regular grid
265 // . . . . . . . . . __ down edge nodes
267 // >->->->->->->->->->->->-> -- direction of processing
269 int g = 0; // number of last processed node in the regular grid
271 // number of last node of the down edge to be processed
272 int stop = nbdown - 1;
273 // if right edge is out, we will stop at a node, previous to the last one
274 if (quad->isEdgeOut[1]) stop--;
276 // for each node of the down edge find nearest node
277 // in the first row of the regular grid and link them
278 for (i = 0; i < stop; i++) {
279 const SMDS_MeshNode *a, *b, *c, *d;
281 b = uv_e0[i + 1].node;
282 gp_Pnt pb (b->X(), b->Y(), b->Z());
284 // find node c in the regular grid, which will be linked with node b
287 // right bound reached, link with the rightmost node
289 c = quad->uv_grid[nbhoriz + iup].node;
292 // find in the grid node c, nearest to the b
293 double mind = RealLast();
294 for (int k = g; k <= iup; k++) {
296 const SMDS_MeshNode *nk;
297 if (k < ilow) // this can be, if left edge is out
298 nk = uv_e3[1].node; // get node from the left edge
300 nk = quad->uv_grid[nbhoriz + k].node; // get one of middle nodes
302 gp_Pnt pnk (nk->X(), nk->Y(), nk->Z());
303 double dist = pb.Distance(pnk);
304 if (dist < mind - eps) {
314 if (near == g) { // make triangle
315 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c);
316 SMDS_MeshFace* face = myTool->AddFace(a, b, c);
317 meshDS->SetMeshElementOnShape(face, geomFaceID);
319 else { // make quadrangle
323 d = quad->uv_grid[nbhoriz + near - 1].node;
324 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c, d);
326 if(!myTrianglePreference){
327 SMDS_MeshFace* face = myTool->AddFace(a, b, c, d);
328 meshDS->SetMeshElementOnShape(face, geomFaceID);
331 SplitQuad(meshDS, geomFaceID, a, b, c, d);
334 // if node d is not at position g - make additional triangles
336 for (int k = near - 1; k > g; k--) {
337 c = quad->uv_grid[nbhoriz + k].node;
341 d = quad->uv_grid[nbhoriz + k - 1].node;
342 //SMDS_MeshFace* face = meshDS->AddFace(a, c, d);
343 SMDS_MeshFace* face = myTool->AddFace(a, c, d);
344 meshDS->SetMeshElementOnShape(face, geomFaceID);
351 if (quad->isEdgeOut[2]) {
354 // <-<-<-<-<-<-<-<-<-<-<-<-< -- direction of processing
356 // . . . . . . . . . __ up edge nodes
357 // ___ ___ ___ ___ ___ ___ __ first row of the regular grid
359 // |___|___|___|___|___|___|
361 // |___|___|___|___|___|___|
364 int g = nbhoriz - 1; // last processed node in the regular grid
367 // if left edge is out, we will stop at a second node
368 if (quad->isEdgeOut[3]) stop++;
370 // for each node of the up edge find nearest node
371 // in the first row of the regular grid and link them
372 for (i = nbup - 1; i > stop; i--) {
373 const SMDS_MeshNode *a, *b, *c, *d;
375 b = uv_e2[i - 1].node;
376 gp_Pnt pb (b->X(), b->Y(), b->Z());
378 // find node c in the grid, which will be linked with node b
380 if (i == stop + 1) { // left bound reached, link with the leftmost node
381 c = quad->uv_grid[nbhoriz*(nbvertic - 2) + ilow].node;
384 // find node c in the grid, nearest to the b
385 double mind = RealLast();
386 for (int k = g; k >= ilow; k--) {
387 const SMDS_MeshNode *nk;
389 nk = uv_e1[nbright - 2].node;
391 nk = quad->uv_grid[nbhoriz*(nbvertic - 2) + k].node;
392 gp_Pnt pnk (nk->X(), nk->Y(), nk->Z());
393 double dist = pb.Distance(pnk);
394 if (dist < mind - eps) {
404 if (near == g) { // make triangle
405 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c);
406 SMDS_MeshFace* face = myTool->AddFace(a, b, c);
407 meshDS->SetMeshElementOnShape(face, geomFaceID);
409 else { // make quadrangle
411 d = uv_e1[nbright - 2].node;
413 d = quad->uv_grid[nbhoriz*(nbvertic - 2) + near + 1].node;
414 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c, d);
415 if(!myTrianglePreference){
416 SMDS_MeshFace* face = myTool->AddFace(a, b, c, d);
417 meshDS->SetMeshElementOnShape(face, geomFaceID);
420 SplitQuad(meshDS, geomFaceID, a, b, c, d);
423 if (near + 1 < g) { // if d not is at g - make additional triangles
424 for (int k = near + 1; k < g; k++) {
425 c = quad->uv_grid[nbhoriz*(nbvertic - 2) + k].node;
427 d = uv_e1[nbright - 2].node;
429 d = quad->uv_grid[nbhoriz*(nbvertic - 2) + k + 1].node;
430 //SMDS_MeshFace* face = meshDS->AddFace(a, c, d);
431 SMDS_MeshFace* face = myTool->AddFace(a, c, d);
432 meshDS->SetMeshElementOnShape(face, geomFaceID);
441 // right or left boundary quadrangles
442 if (quad->isEdgeOut[1]) {
443 // MESSAGE("right edge is out");
444 int g = 0; // last processed node in the grid
445 int stop = nbright - 1;
446 if (quad->isEdgeOut[2]) stop--;
447 for (i = 0; i < stop; i++) {
448 const SMDS_MeshNode *a, *b, *c, *d;
450 b = uv_e1[i + 1].node;
451 gp_Pnt pb (b->X(), b->Y(), b->Z());
453 // find node c in the grid, nearest to the b
455 if (i == stop - 1) { // up bondary reached
456 c = quad->uv_grid[nbhoriz*(jup + 1) - 2].node;
459 double mind = RealLast();
460 for (int k = g; k <= jup; k++) {
461 const SMDS_MeshNode *nk;
463 nk = uv_e0[nbdown - 2].node;
465 nk = quad->uv_grid[nbhoriz*(k + 1) - 2].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 = meshDS->AddFace(a, b, c);
480 SMDS_MeshFace* face = myTool->AddFace(a, b, c);
481 meshDS->SetMeshElementOnShape(face, geomFaceID);
483 else { // make quadrangle
485 d = uv_e0[nbdown - 2].node;
487 d = quad->uv_grid[nbhoriz*near - 2].node;
488 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c, d);
490 if(!myTrianglePreference){
491 SMDS_MeshFace* face = myTool->AddFace(a, b, c, d);
492 meshDS->SetMeshElementOnShape(face, geomFaceID);
495 SplitQuad(meshDS, geomFaceID, a, b, c, d);
498 if (near - 1 > g) { // if d not is at g - make additional triangles
499 for (int k = near - 1; k > g; k--) {
500 c = quad->uv_grid[nbhoriz*(k + 1) - 2].node;
502 d = uv_e0[nbdown - 2].node;
504 d = quad->uv_grid[nbhoriz*k - 2].node;
505 //SMDS_MeshFace* face = meshDS->AddFace(a, c, d);
506 SMDS_MeshFace* face = myTool->AddFace(a, c, d);
507 meshDS->SetMeshElementOnShape(face, geomFaceID);
514 if (quad->isEdgeOut[3]) {
515 // MESSAGE("left edge is out");
516 int g = nbvertic - 1; // last processed node in the grid
518 if (quad->isEdgeOut[0]) stop++;
519 for (i = nbleft - 1; i > stop; i--) {
520 const SMDS_MeshNode *a, *b, *c, *d;
522 b = uv_e3[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) { // down bondary reached
528 c = quad->uv_grid[nbhoriz*jlow + 1].node;
531 double mind = RealLast();
532 for (int k = g; k >= jlow; k--) {
533 const SMDS_MeshNode *nk;
537 nk = quad->uv_grid[nbhoriz*k + 1].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 = meshDS->AddFace(a, b, c);
552 SMDS_MeshFace* face = myTool->AddFace(a, b, c);
553 meshDS->SetMeshElementOnShape(face, geomFaceID);
555 else { // make quadrangle
559 d = quad->uv_grid[nbhoriz*(near + 1) + 1].node;
560 //SMDS_MeshFace* face = meshDS->AddFace(a, b, c, d);
561 if(!myTrianglePreference){
562 SMDS_MeshFace* face = myTool->AddFace(a, b, c, d);
563 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].node;
575 d = quad->uv_grid[nbhoriz*(k + 1) + 1].node;
576 //SMDS_MeshFace* face = meshDS->AddFace(a, c, d);
577 SMDS_MeshFace* face = myTool->AddFace(a, c, d);
578 meshDS->SetMeshElementOnShape(face, geomFaceID);
591 //=============================================================================
595 //=============================================================================
597 FaceQuadStruct* StdMeshers_Quadrangle_2D::CheckNbEdges(SMESH_Mesh & aMesh,
598 const TopoDS_Shape & aShape)
599 //throw(SALOME_Exception)
601 const TopoDS_Face & F = TopoDS::Face(aShape);
602 const bool ignoreMediumNodes = _quadraticMesh;
604 // verify 1 wire only, with 4 edges
606 list< TopoDS_Edge > edges;
607 list< int > nbEdgesInWire;
608 int nbWire = SMESH_Block::GetOrderedEdges (F, V, edges, nbEdgesInWire);
610 error(COMPERR_BAD_SHAPE, TComm("Wrong number of wires: ") << nbWire);
613 FaceQuadStruct* quad = new FaceQuadStruct;
615 quad->side.reserve(nbEdgesInWire.front());
618 list< TopoDS_Edge >::iterator edgeIt = edges.begin();
619 if ( nbEdgesInWire.front() == 4 ) { // exactly 4 edges
620 for ( ; edgeIt != edges.end(); ++edgeIt, nbSides++ )
621 quad->side.push_back( new StdMeshers_FaceSide(F, *edgeIt, &aMesh,
622 nbSides<TOP_SIDE, ignoreMediumNodes));
624 else if ( nbEdgesInWire.front() > 4 ) { // more than 4 edges - try to unite some
625 list< TopoDS_Edge > sideEdges;
626 while ( !edges.empty()) {
628 sideEdges.splice( sideEdges.end(), edges, edges.begin()); // edges.front() -> sideEdges.end()
629 bool sameSide = true;
630 while ( !edges.empty() && sameSide ) {
631 sameSide = SMESH_Algo::IsContinuous( sideEdges.back(), edges.front() );
633 sideEdges.splice( sideEdges.end(), edges, edges.begin());
635 if ( nbSides == 0 ) { // go backward from the first edge
637 while ( !edges.empty() && sameSide ) {
638 sameSide = SMESH_Algo::IsContinuous( sideEdges.front(), edges.back() );
640 sideEdges.splice( sideEdges.begin(), edges, --edges.end());
643 quad->side.push_back( new StdMeshers_FaceSide(F, sideEdges, &aMesh,
644 nbSides<TOP_SIDE, ignoreMediumNodes));
650 MESSAGE ( "StdMeshers_Quadrangle_2D. Edge IDs of " << nbSides << " sides:\n" );
651 for ( int i = 0; i < nbSides; ++i ) {
653 for ( int e = 0; e < quad->side[i]->NbEdges(); ++e )
654 MESSAGE ( myTool->GetMeshDS()->ShapeToIndex( quad->side[i]->Edge( e )) << " " );
660 nbSides = nbEdgesInWire.front();
661 error(COMPERR_BAD_SHAPE, TComm("Face must have 4 sides but not ") << nbSides);
669 //=============================================================================
673 //=============================================================================
675 FaceQuadStruct *StdMeshers_Quadrangle_2D::CheckAnd2Dcompute
677 const TopoDS_Shape & aShape,
678 const bool CreateQuadratic) //throw(SALOME_Exception)
680 _quadraticMesh = CreateQuadratic;
682 FaceQuadStruct *quad = CheckNbEdges(aMesh, aShape);
686 // set normalized grid on unit square in parametric domain
687 bool stat = SetNormalizedGrid(aMesh, aShape, quad);
697 //=============================================================================
701 //=============================================================================
703 faceQuadStruct::~faceQuadStruct()
705 for (int i = 0; i < side.size(); i++) {
706 if (side[i]) delete side[i];
708 if (uv_grid) delete [] uv_grid;
712 inline const vector<UVPtStruct>& GetUVPtStructIn(FaceQuadStruct* quad, int i, int nbSeg)
714 bool isXConst = ( i == BOTTOM_SIDE || i == TOP_SIDE );
715 double constValue = ( i == BOTTOM_SIDE || i == LEFT_SIDE ) ? 0 : 1;
718 quad->side[i]->SimulateUVPtStruct(nbSeg,isXConst,constValue) :
719 quad->side[i]->GetUVPtStruct(isXConst,constValue);
723 //=============================================================================
727 //=============================================================================
729 bool StdMeshers_Quadrangle_2D::SetNormalizedGrid (SMESH_Mesh & aMesh,
730 const TopoDS_Shape& aShape,
731 FaceQuadStruct* & quad) //throw (SALOME_Exception)
733 // Algorithme décrit dans "Génération automatique de maillages"
734 // P.L. GEORGE, MASSON, § 6.4.1 p. 84-85
735 // traitement dans le domaine paramétrique 2d u,v
736 // transport - projection sur le carré unité
738 // MESSAGE("StdMeshers_Quadrangle_2D::SetNormalizedGrid");
739 // const TopoDS_Face& F = TopoDS::Face(aShape);
741 // 1 --- find orientation of the 4 edges, by test on extrema
744 // |<----north-2-------^ a3 -------------> a2
746 // west-3 east-1 =right | |
750 // v----south-0--------> a0 -------------> a1
755 // 3 --- 2D normalized values on unit square [0..1][0..1]
757 int nbhoriz = Min(quad->side[0]->NbPoints(), quad->side[2]->NbPoints());
758 int nbvertic = Min(quad->side[1]->NbPoints(), quad->side[3]->NbPoints());
760 quad->isEdgeOut[0] = (quad->side[0]->NbPoints() > quad->side[2]->NbPoints());
761 quad->isEdgeOut[1] = (quad->side[1]->NbPoints() > quad->side[3]->NbPoints());
762 quad->isEdgeOut[2] = (quad->side[2]->NbPoints() > quad->side[0]->NbPoints());
763 quad->isEdgeOut[3] = (quad->side[3]->NbPoints() > quad->side[1]->NbPoints());
765 UVPtStruct *uv_grid = quad->uv_grid = new UVPtStruct[nbvertic * nbhoriz];
767 const vector<UVPtStruct>& uv_e0 = GetUVPtStructIn( quad, 0, nbhoriz - 1 );
768 const vector<UVPtStruct>& uv_e1 = GetUVPtStructIn( quad, 1, nbvertic - 1 );
769 const vector<UVPtStruct>& uv_e2 = GetUVPtStructIn( quad, 2, nbhoriz - 1 );
770 const vector<UVPtStruct>& uv_e3 = GetUVPtStructIn( quad, 3, nbvertic - 1 );
772 if ( uv_e0.empty() || uv_e1.empty() || uv_e2.empty() || uv_e3.empty() )
773 //return error( "Can't find nodes on sides");
774 return error( COMPERR_BAD_INPUT_MESH );
776 // nodes Id on "in" edges
777 if (! quad->isEdgeOut[0]) {
779 for (int i = 0; i < nbhoriz; i++) { // down
780 int ij = j * nbhoriz + i;
781 uv_grid[ij].node = uv_e0[i].node;
784 if (! quad->isEdgeOut[1]) {
786 for (int j = 0; j < nbvertic; j++) { // right
787 int ij = j * nbhoriz + i;
788 uv_grid[ij].node = uv_e1[j].node;
791 if (! quad->isEdgeOut[2]) {
792 int j = nbvertic - 1;
793 for (int i = 0; i < nbhoriz; i++) { // up
794 int ij = j * nbhoriz + i;
795 uv_grid[ij].node = uv_e2[i].node;
798 if (! quad->isEdgeOut[3]) {
800 for (int j = 0; j < nbvertic; j++) { // left
801 int ij = j * nbhoriz + i;
802 uv_grid[ij].node = uv_e3[j].node;
806 // normalized 2d values on grid
807 for (int i = 0; i < nbhoriz; i++)
809 for (int j = 0; j < nbvertic; j++)
811 int ij = j * nbhoriz + i;
812 // --- droite i cste : x = x0 + y(x1-x0)
813 double x0 = uv_e0[i].normParam; // bas - sud
814 double x1 = uv_e2[i].normParam; // haut - nord
815 // --- droite j cste : y = y0 + x(y1-y0)
816 double y0 = uv_e3[j].normParam; // gauche-ouest
817 double y1 = uv_e1[j].normParam; // droite - est
818 // --- intersection : x=x0+(y0+x(y1-y0))(x1-x0)
819 double x = (x0 + y0 * (x1 - x0)) / (1 - (y1 - y0) * (x1 - x0));
820 double y = y0 + x * (y1 - y0);
823 //MESSAGE("-xy-01 "<<x0<<" "<<x1<<" "<<y0<<" "<<y1);
824 //MESSAGE("-xy-norm "<<i<<" "<<j<<" "<<x<<" "<<y);
828 // 4 --- projection on 2d domain (u,v)
829 gp_UV a0( uv_e0.front().u, uv_e0.front().v );
830 gp_UV a1( uv_e0.back().u, uv_e0.back().v );
831 gp_UV a2( uv_e2.back().u, uv_e2.back().v );
832 gp_UV a3( uv_e2.front().u, uv_e2.front().v );
834 for (int i = 0; i < nbhoriz; i++)
836 for (int j = 0; j < nbvertic; j++)
838 int ij = j * nbhoriz + i;
839 double x = uv_grid[ij].x;
840 double y = uv_grid[ij].y;
841 double param_0 = uv_e0[0].normParam + x * (uv_e0.back().normParam - uv_e0[0].normParam); // sud
842 double param_2 = uv_e2[0].normParam + x * (uv_e2.back().normParam - uv_e2[0].normParam); // nord
843 double param_1 = uv_e1[0].normParam + y * (uv_e1.back().normParam - uv_e1[0].normParam); // est
844 double param_3 = uv_e3[0].normParam + y * (uv_e3.back().normParam - uv_e3[0].normParam); // ouest
846 //MESSAGE("params "<<param_0<<" "<<param_1<<" "<<param_2<<" "<<param_3);
847 gp_UV p0 = quad->side[0]->Value2d(param_0).XY();
848 gp_UV p1 = quad->side[1]->Value2d(param_1).XY();
849 gp_UV p2 = quad->side[2]->Value2d(param_2).XY();
850 gp_UV p3 = quad->side[3]->Value2d(param_3).XY();
852 gp_UV uv = (1 - y) * p0 + x * p1 + y * p2 + (1 - x) * p3;
853 uv -= (1 - x) * (1 - y) * a0 + x * (1 - y) * a1 + x * y * a2 + (1 - x) * y * a3;
855 uv_grid[ij].u = uv.X();
856 uv_grid[ij].v = uv.Y();
862 //=======================================================================
863 //function : ShiftQuad
864 //purpose : auxilary function for ComputeQuadPref
865 //=======================================================================
867 static void ShiftQuad(FaceQuadStruct* quad, const int num, bool)
869 StdMeshers_FaceSide* side[4] = { quad->side[0], quad->side[1], quad->side[2], quad->side[3] };
870 for (int i = BOTTOM_SIDE; i < NB_SIDES; ++i ) {
871 int id = ( i + num ) % NB_SIDES;
872 bool wasForward = ( i < TOP_SIDE );
873 bool newForward = ( id < TOP_SIDE );
874 if ( wasForward != newForward )
875 side[ i ]->Reverse();
876 quad->side[ id ] = side[ i ];
880 //=======================================================================
882 //purpose : auxilary function for ComputeQuadPref
883 //=======================================================================
885 static gp_UV CalcUV(double x0, double x1, double y0, double y1,
886 FaceQuadStruct* quad,
887 const gp_UV& a0, const gp_UV& a1,
888 const gp_UV& a2, const gp_UV& a3)
890 const vector<UVPtStruct>& uv_eb = quad->side[0]->GetUVPtStruct(true,0 );
891 const vector<UVPtStruct>& uv_er = quad->side[1]->GetUVPtStruct(false,1);
892 const vector<UVPtStruct>& uv_et = quad->side[2]->GetUVPtStruct(true,1 );
893 const vector<UVPtStruct>& uv_el = quad->side[3]->GetUVPtStruct(false,0);
895 double x = (x0 + y0 * (x1 - x0)) / (1 - (y1 - y0) * (x1 - x0));
896 double y = y0 + x * (y1 - y0);
898 double param_b = uv_eb[0].normParam + x * (uv_eb.back().normParam - uv_eb[0].normParam);
899 double param_t = uv_et[0].normParam + x * (uv_et.back().normParam - uv_et[0].normParam);
900 double param_r = uv_er[0].normParam + y * (uv_er.back().normParam - uv_er[0].normParam);
901 double param_l = uv_el[0].normParam + y * (uv_el.back().normParam - uv_el[0].normParam);
903 gp_UV p0 = quad->side[BOTTOM_SIDE]->Value2d(param_b).XY();
904 gp_UV p1 = quad->side[RIGHT_SIDE ]->Value2d(param_r).XY();
905 gp_UV p2 = quad->side[TOP_SIDE ]->Value2d(param_t).XY();
906 gp_UV p3 = quad->side[LEFT_SIDE ]->Value2d(param_l).XY();
908 gp_UV uv = p0 * (1 - y) + p1 * x + p2 * y + p3 * (1 - x);
910 uv -= (1 - x) * (1 - y) * a0 + x * (1 - y) * a1 + x * y * a2 + (1 - x) * y * a3;
916 //=======================================================================
918 //purpose : auxilary function for ComputeQuadPref
919 //=======================================================================
921 static gp_UV CalcUV2(double x, double y,
922 FaceQuadStruct* quad,
923 const gp_UV& a0, const gp_UV& a1,
924 const gp_UV& a2, const gp_UV& a3)
926 const vector<UVPtStruct>& uv_eb = quad->side[0]->GetUVPtStruct(true,0 );
927 const vector<UVPtStruct>& uv_er = quad->side[1]->GetUVPtStruct(false,1);
928 const vector<UVPtStruct>& uv_et = quad->side[2]->GetUVPtStruct(true,1 );
929 const vector<UVPtStruct>& uv_el = quad->side[3]->GetUVPtStruct(false,0);
931 //double x = (x0 + y0 * (x1 - x0)) / (1 - (y1 - y0) * (x1 - x0));
932 //double y = y0 + x * (y1 - y0);
934 double param_b = uv_eb[0].normParam + x * (uv_eb.back().normParam - uv_eb[0].normParam);
935 double param_t = uv_et[0].normParam + x * (uv_et.back().normParam - uv_et[0].normParam);
936 double param_r = uv_er[0].normParam + y * (uv_er.back().normParam - uv_er[0].normParam);
937 double param_l = uv_el[0].normParam + y * (uv_el.back().normParam - uv_el[0].normParam);
939 gp_UV p0 = quad->side[BOTTOM_SIDE]->Value2d(param_b).XY();
940 gp_UV p1 = quad->side[RIGHT_SIDE ]->Value2d(param_r).XY();
941 gp_UV p2 = quad->side[TOP_SIDE ]->Value2d(param_t).XY();
942 gp_UV p3 = quad->side[LEFT_SIDE ]->Value2d(param_l).XY();
944 gp_UV uv = p0 * (1 - y) + p1 * x + p2 * y + p3 * (1 - x);
946 uv -= (1 - x) * (1 - y) * a0 + x * (1 - y) * a1 + x * y * a2 + (1 - x) * y * a3;
951 //=======================================================================
953 * Create only quandrangle faces
955 //=======================================================================
957 bool StdMeshers_Quadrangle_2D::ComputeQuadPref (SMESH_Mesh & aMesh,
958 const TopoDS_Shape& aShape,
959 FaceQuadStruct* quad)
961 // Auxilary key in order to keep old variant
962 // of meshing after implementation new variant
963 // for bug 0016220 from Mantis.
964 bool OldVersion = false;
966 SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
967 const TopoDS_Face& F = TopoDS::Face(aShape);
968 Handle(Geom_Surface) S = BRep_Tool::Surface(F);
969 // const TopoDS_Wire& W = BRepTools::OuterWire(F);
971 // if(W.Orientation()==TopAbs_FORWARD)
973 //if(WisF) cout<<"W is FORWARD"<<endl;
974 //else cout<<"W is REVERSED"<<endl;
975 // bool FisF = (F.Orientation()==TopAbs_FORWARD);
976 // if(!FisF) WisF = !WisF;
978 int i,j,geomFaceID = meshDS->ShapeToIndex( F );
980 int nb = quad->side[0]->NbPoints();
981 int nr = quad->side[1]->NbPoints();
982 int nt = quad->side[2]->NbPoints();
983 int nl = quad->side[3]->NbPoints();
989 // it is a base case => not shift quad but me be replacement is need
990 ShiftQuad(quad,0,WisF);
993 // we have to shift quad on 2
994 ShiftQuad(quad,2,WisF);
999 // we have to shift quad on 1
1000 ShiftQuad(quad,1,WisF);
1003 // we have to shift quad on 3
1004 ShiftQuad(quad,3,WisF);
1008 nb = quad->side[0]->NbPoints();
1009 nr = quad->side[1]->NbPoints();
1010 nt = quad->side[2]->NbPoints();
1011 nl = quad->side[3]->NbPoints();
1014 int nbh = Max(nb,nt);
1015 int nbv = Max(nr,nl);
1019 // ----------- Old version ---------------
1020 // orientation of face and 3 main domain for future faces
1026 // left | | | | rigth
1033 // ----------- New version ---------------
1034 // orientation of face and 3 main domain for future faces
1040 // left |/________\| rigth
1056 const vector<UVPtStruct>& uv_eb = quad->side[0]->GetUVPtStruct(true,0 );
1057 const vector<UVPtStruct>& uv_er = quad->side[1]->GetUVPtStruct(false,1);
1058 const vector<UVPtStruct>& uv_et = quad->side[2]->GetUVPtStruct(true,1 );
1059 const vector<UVPtStruct>& uv_el = quad->side[3]->GetUVPtStruct(false,0);
1061 // arrays for normalized params
1062 //cout<<"Dump B:"<<endl;
1063 TColStd_SequenceOfReal npb, npr, npt, npl;
1064 for(i=0; i<nb; i++) {
1065 npb.Append(uv_eb[i].normParam);
1066 //cout<<"i="<<i<<" par="<<uv_eb[i].normParam<<" npar="<<uv_eb[i].normParam;
1067 //const SMDS_MeshNode* N = uv_eb[i].node;
1068 //cout<<" node("<<N->X()<<","<<N->Y()<<","<<N->Z()<<")"<<endl;
1070 for(i=0; i<nr; i++) {
1071 npr.Append(uv_er[i].normParam);
1073 for(i=0; i<nt; i++) {
1074 npt.Append(uv_et[i].normParam);
1076 for(i=0; i<nl; i++) {
1077 npl.Append(uv_el[i].normParam);
1082 // add some params to right and left after the first param
1085 double dpr = (npr.Value(2) - npr.Value(1))/(dr+1);
1086 for(i=1; i<=dr; i++) {
1087 npr.InsertAfter(1,npr.Value(2)-dpr);
1091 dpr = (npl.Value(2) - npl.Value(1))/(dl+1);
1092 for(i=1; i<=dl; i++) {
1093 npl.InsertAfter(1,npl.Value(2)-dpr);
1097 //for(i=1; i<=npb.Length(); i++) {
1098 // cout<<" "<<npb.Value(i);
1102 gp_XY a0( uv_eb.front().u, uv_eb.front().v );
1103 gp_XY a1( uv_eb.back().u, uv_eb.back().v );
1104 gp_XY a2( uv_et.back().u, uv_et.back().v );
1105 gp_XY a3( uv_et.front().u, uv_et.front().v );
1106 //cout<<" a0("<<a0.X()<<","<<a0.Y()<<")"<<" a1("<<a1.X()<<","<<a1.Y()<<")"
1107 // <<" a2("<<a2.X()<<","<<a2.Y()<<")"<<" a3("<<a3.X()<<","<<a3.Y()<<")"<<endl;
1109 int nnn = Min(nr,nl);
1110 // auxilary sequence of XY for creation nodes
1111 // in the bottom part of central domain
1112 // it's length must be == nbv-nnn-1
1113 TColgp_SequenceOfXY UVL;
1114 TColgp_SequenceOfXY UVR;
1117 // step1: create faces for left domain
1118 StdMeshers_Array2OfNode NodesL(1,dl+1,1,nl);
1120 for(j=1; j<=nl; j++)
1121 NodesL.SetValue(1,j,uv_el[j-1].node);
1124 for(i=1; i<=dl; i++)
1125 NodesL.SetValue(i+1,nl,uv_et[i].node);
1126 // create and add needed nodes
1127 TColgp_SequenceOfXY UVtmp;
1128 for(i=1; i<=dl; i++) {
1129 double x0 = npt.Value(i+1);
1132 double y0 = npl.Value(i+1);
1133 double y1 = npr.Value(i+1);
1134 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1135 gp_Pnt P = S->Value(UV.X(),UV.Y());
1136 SMDS_MeshNode * N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1137 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1138 NodesL.SetValue(i+1,1,N);
1139 if(UVL.Length()<nbv-nnn-1) UVL.Append(UV);
1141 for(j=2; j<nl; j++) {
1142 double y0 = npl.Value(dl+j);
1143 double y1 = npr.Value(dl+j);
1144 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1145 gp_Pnt P = S->Value(UV.X(),UV.Y());
1146 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1147 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1148 NodesL.SetValue(i+1,j,N);
1149 if( i==dl ) UVtmp.Append(UV);
1152 for(i=1; i<=UVtmp.Length() && UVL.Length()<nbv-nnn-1; i++) {
1153 UVL.Append(UVtmp.Value(i));
1155 //cout<<"Dump NodesL:"<<endl;
1156 //for(i=1; i<=dl+1; i++) {
1158 // for(j=1; j<=nl; j++) {
1159 // cout<<" ("<<NodesL.Value(i,j)->X()<<","<<NodesL.Value(i,j)->Y()<<","<<NodesL.Value(i,j)->Z()<<")";
1164 for(i=1; i<=dl; i++) {
1165 for(j=1; j<nl; j++) {
1168 myTool->AddFace(NodesL.Value(i,j), NodesL.Value(i+1,j),
1169 NodesL.Value(i+1,j+1), NodesL.Value(i,j+1));
1170 meshDS->SetMeshElementOnShape(F, geomFaceID);
1174 myTool->AddFace(NodesL.Value(i,j), NodesL.Value(i,j+1),
1175 NodesL.Value(i+1,j+1), NodesL.Value(i+1,j));
1176 meshDS->SetMeshElementOnShape(F, geomFaceID);
1182 // fill UVL using c2d
1183 for(i=1; i<npl.Length() && UVL.Length()<nbv-nnn-1; i++) {
1184 UVL.Append( gp_UV ( uv_el[i].u, uv_el[i].v ));
1188 // step2: create faces for right domain
1189 StdMeshers_Array2OfNode NodesR(1,dr+1,1,nr);
1191 for(j=1; j<=nr; j++)
1192 NodesR.SetValue(1,j,uv_er[nr-j].node);
1195 for(i=1; i<=dr; i++)
1196 NodesR.SetValue(i+1,1,uv_et[nt-1-i].node);
1197 // create and add needed nodes
1198 TColgp_SequenceOfXY UVtmp;
1199 for(i=1; i<=dr; i++) {
1200 double x0 = npt.Value(nt-i);
1203 double y0 = npl.Value(i+1);
1204 double y1 = npr.Value(i+1);
1205 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1206 gp_Pnt P = S->Value(UV.X(),UV.Y());
1207 SMDS_MeshNode * N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1208 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1209 NodesR.SetValue(i+1,nr,N);
1210 if(UVR.Length()<nbv-nnn-1) UVR.Append(UV);
1212 for(j=2; j<nr; j++) {
1213 double y0 = npl.Value(nbv-j+1);
1214 double y1 = npr.Value(nbv-j+1);
1215 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1216 gp_Pnt P = S->Value(UV.X(),UV.Y());
1217 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1218 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1219 NodesR.SetValue(i+1,j,N);
1220 if( i==dr ) UVtmp.Prepend(UV);
1223 for(i=1; i<=UVtmp.Length() && UVR.Length()<nbv-nnn-1; i++) {
1224 UVR.Append(UVtmp.Value(i));
1227 for(i=1; i<=dr; i++) {
1228 for(j=1; j<nr; j++) {
1231 myTool->AddFace(NodesR.Value(i,j), NodesR.Value(i+1,j),
1232 NodesR.Value(i+1,j+1), NodesR.Value(i,j+1));
1233 meshDS->SetMeshElementOnShape(F, geomFaceID);
1237 myTool->AddFace(NodesR.Value(i,j), NodesR.Value(i,j+1),
1238 NodesR.Value(i+1,j+1), NodesR.Value(i+1,j));
1239 meshDS->SetMeshElementOnShape(F, geomFaceID);
1245 // fill UVR using c2d
1246 for(i=1; i<npr.Length() && UVR.Length()<nbv-nnn-1; i++) {
1247 UVR.Append( gp_UV( uv_er[i].u, uv_er[i].v ));
1251 // step3: create faces for central domain
1252 StdMeshers_Array2OfNode NodesC(1,nb,1,nbv);
1253 // add first string using NodesL
1254 for(i=1; i<=dl+1; i++)
1255 NodesC.SetValue(1,i,NodesL(i,1));
1256 for(i=2; i<=nl; i++)
1257 NodesC.SetValue(1,dl+i,NodesL(dl+1,i));
1258 // add last string using NodesR
1259 for(i=1; i<=dr+1; i++)
1260 NodesC.SetValue(nb,i,NodesR(i,nr));
1262 NodesC.SetValue(nb,dr+i+1,NodesR(dr+1,nr-i));
1263 // add top nodes (last columns)
1264 for(i=dl+2; i<nbh-dr; i++)
1265 NodesC.SetValue(i-dl,nbv,uv_et[i-1].node);
1266 // add bottom nodes (first columns)
1268 NodesC.SetValue(i,1,uv_eb[i-1].node);
1270 // create and add needed nodes
1271 // add linear layers
1272 for(i=2; i<nb; i++) {
1273 double x0 = npt.Value(dl+i);
1275 for(j=1; j<nnn; j++) {
1276 double y0 = npl.Value(nbv-nnn+j);
1277 double y1 = npr.Value(nbv-nnn+j);
1278 gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3);
1279 gp_Pnt P = S->Value(UV.X(),UV.Y());
1280 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1281 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1282 NodesC.SetValue(i,nbv-nnn+j,N);
1285 // add diagonal layers
1286 //cout<<"UVL.Length()="<<UVL.Length()<<" UVR.Length()="<<UVR.Length()<<endl;
1287 //cout<<"Dump UVL:"<<endl;
1288 //for(i=1; i<=UVL.Length(); i++) {
1289 // cout<<" ("<<UVL.Value(i).X()<<","<<UVL.Value(i).Y()<<")";
1292 for(i=1; i<nbv-nnn; i++) {
1293 double du = UVR.Value(i).X() - UVL.Value(i).X();
1294 double dv = UVR.Value(i).Y() - UVL.Value(i).Y();
1295 for(j=2; j<nb; j++) {
1296 double u = UVL.Value(i).X() + du*npb.Value(j);
1297 double v = UVL.Value(i).Y() + dv*npb.Value(j);
1298 gp_Pnt P = S->Value(u,v);
1299 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1300 meshDS->SetNodeOnFace(N, geomFaceID, u, v);
1301 NodesC.SetValue(j,i+1,N);
1305 for(i=1; i<nb; i++) {
1306 for(j=1; j<nbv; j++) {
1309 myTool->AddFace(NodesC.Value(i,j), NodesC.Value(i+1,j),
1310 NodesC.Value(i+1,j+1), NodesC.Value(i,j+1));
1311 meshDS->SetMeshElementOnShape(F, geomFaceID);
1315 myTool->AddFace(NodesC.Value(i,j), NodesC.Value(i,j+1),
1316 NodesC.Value(i+1,j+1), NodesC.Value(i+1,j));
1317 meshDS->SetMeshElementOnShape(F, geomFaceID);
1323 else { // New version (!OldVersion)
1324 // step1: create faces for bottom rectangle domain
1325 StdMeshers_Array2OfNode NodesBRD(1,nb,1,nnn-1);
1326 // fill UVL and UVR using c2d
1327 for(j=0; j<nb; j++) {
1328 NodesBRD.SetValue(j+1,1,uv_eb[j].node);
1330 for(i=1; i<nnn-1; i++) {
1331 NodesBRD.SetValue(1,i+1,uv_el[i].node);
1332 NodesBRD.SetValue(nb,i+1,uv_er[i].node);
1333 double du = uv_er[i].u - uv_el[i].u;
1334 double dv = uv_er[i].v - uv_el[i].v;
1335 for(j=2; j<nb; j++) {
1336 double u = uv_el[i].u + du*npb.Value(j);
1337 double v = uv_el[i].v + dv*npb.Value(j);
1338 gp_Pnt P = S->Value(u,v);
1339 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1340 meshDS->SetNodeOnFace(N, geomFaceID, u, v);
1341 NodesBRD.SetValue(j,i+1,N);
1345 for(j=1; j<nnn-1; j++) {
1346 for(i=1; i<nb; i++) {
1349 myTool->AddFace(NodesBRD.Value(i,j), NodesBRD.Value(i+1,j),
1350 NodesBRD.Value(i+1,j+1), NodesBRD.Value(i,j+1));
1351 meshDS->SetMeshElementOnShape(F, geomFaceID);
1355 myTool->AddFace(NodesBRD.Value(i,j), NodesBRD.Value(i,j+1),
1356 NodesBRD.Value(i+1,j+1), NodesBRD.Value(i+1,j));
1357 meshDS->SetMeshElementOnShape(F, geomFaceID);
1362 int drl = abs(nr-nl);
1363 // create faces for region C
1364 StdMeshers_Array2OfNode NodesC(1,nb,1,drl+1+addv);
1365 // add nodes from previous region
1366 for(j=1; j<=nb; j++) {
1367 NodesC.SetValue(j,1,NodesBRD.Value(j,nnn-1));
1369 if( (drl+addv) > 0 ) {
1374 TColgp_SequenceOfXY UVtmp;
1375 double drparam = npr.Value(nr) - npr.Value(nnn-1);
1376 double dlparam = npl.Value(nnn) - npl.Value(nnn-1);
1378 for(i=1; i<=drl; i++) {
1379 // add existed nodes from right edge
1380 NodesC.SetValue(nb,i+1,uv_er[nnn+i-2].node);
1381 //double dtparam = npt.Value(i+1);
1382 y1 = npr.Value(nnn+i-1); // param on right edge
1383 double dpar = (y1 - npr.Value(nnn-1))/drparam;
1384 y0 = npl.Value(nnn-1) + dpar*dlparam; // param on left edge
1385 double dy = y1 - y0;
1386 for(j=1; j<nb; j++) {
1387 double x = npt.Value(i+1) + npb.Value(j)*(1-npt.Value(i+1));
1388 double y = y0 + dy*x;
1389 gp_UV UV = CalcUV2(x, y, quad, a0, a1, a2, a3);
1390 gp_Pnt P = S->Value(UV.X(),UV.Y());
1391 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1392 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1393 NodesC.SetValue(j,i+1,N);
1396 double dy0 = (1-y0)/(addv+1);
1397 double dy1 = (1-y1)/(addv+1);
1398 for(i=1; i<=addv; i++) {
1399 double yy0 = y0 + dy0*i;
1400 double yy1 = y1 + dy1*i;
1401 double dyy = yy1 - yy0;
1402 for(j=1; j<=nb; j++) {
1403 double x = npt.Value(i+1+drl) +
1404 npb.Value(j) * ( npt.Value(nt-i) - npt.Value(i+1+drl) );
1405 double y = yy0 + dyy*x;
1406 gp_UV UV = CalcUV2(x, y, quad, a0, a1, a2, a3);
1407 gp_Pnt P = S->Value(UV.X(),UV.Y());
1408 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1409 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1410 NodesC.SetValue(j,i+drl+1,N);
1417 TColgp_SequenceOfXY UVtmp;
1418 double dlparam = npl.Value(nl) - npl.Value(nnn-1);
1419 double drparam = npr.Value(nnn) - npr.Value(nnn-1);
1420 double y0 = npl.Value(nnn-1);
1421 double y1 = npr.Value(nnn-1);
1422 for(i=1; i<=drl; i++) {
1423 // add existed nodes from right edge
1424 NodesC.SetValue(1,i+1,uv_el[nnn+i-2].node);
1425 y0 = npl.Value(nnn+i-1); // param on left edge
1426 double dpar = (y0 - npl.Value(nnn-1))/dlparam;
1427 y1 = npr.Value(nnn-1) + dpar*drparam; // param on right edge
1428 double dy = y1 - y0;
1429 for(j=2; j<=nb; j++) {
1430 double x = npb.Value(j)*npt.Value(nt-i);
1431 double y = y0 + dy*x;
1432 gp_UV UV = CalcUV2(x, y, quad, a0, a1, a2, a3);
1433 gp_Pnt P = S->Value(UV.X(),UV.Y());
1434 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1435 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1436 NodesC.SetValue(j,i+1,N);
1439 double dy0 = (1-y0)/(addv+1);
1440 double dy1 = (1-y1)/(addv+1);
1441 for(i=1; i<=addv; i++) {
1442 double yy0 = y0 + dy0*i;
1443 double yy1 = y1 + dy1*i;
1444 double dyy = yy1 - yy0;
1445 for(j=1; j<=nb; j++) {
1446 double x = npt.Value(i+1) +
1447 npb.Value(j) * ( npt.Value(nt-i-drl) - npt.Value(i+1) );
1448 double y = yy0 + dyy*x;
1449 gp_UV UV = CalcUV2(x, y, quad, a0, a1, a2, a3);
1450 gp_Pnt P = S->Value(UV.X(),UV.Y());
1451 SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z());
1452 meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y());
1453 NodesC.SetValue(j,i+drl+1,N);
1458 for(j=1; j<=drl+addv; j++) {
1459 for(i=1; i<nb; i++) {
1462 myTool->AddFace(NodesC.Value(i,j), NodesC.Value(i+1,j),
1463 NodesC.Value(i+1,j+1), NodesC.Value(i,j+1));
1464 meshDS->SetMeshElementOnShape(F, geomFaceID);
1468 myTool->AddFace(NodesC.Value(i,j), NodesC.Value(i,j+1),
1469 NodesC.Value(i+1,j+1), NodesC.Value(i+1,j));
1470 meshDS->SetMeshElementOnShape(F, geomFaceID);
1475 StdMeshers_Array2OfNode NodesLast(1,nt,1,2);
1476 for(i=1; i<=nt; i++) {
1477 NodesLast.SetValue(i,2,uv_et[i-1].node);
1480 for(i=n1; i<drl+addv+1; i++) {
1482 NodesLast.SetValue(nnn,1,NodesC.Value(1,i));
1484 for(i=1; i<=nb; i++) {
1486 NodesLast.SetValue(nnn,1,NodesC.Value(i,drl+addv+1));
1488 for(i=drl+addv; i>=n2; i--) {
1490 NodesLast.SetValue(nnn,1,NodesC.Value(nb,i));
1492 for(i=1; i<nt; i++) {
1495 myTool->AddFace(NodesLast.Value(i,1), NodesLast.Value(i+1,1),
1496 NodesLast.Value(i+1,2), NodesLast.Value(i,2));
1497 meshDS->SetMeshElementOnShape(F, geomFaceID);
1501 myTool->AddFace(NodesLast.Value(i,1), NodesLast.Value(i,2),
1502 NodesLast.Value(i+1,2), NodesLast.Value(i+1,2));
1503 meshDS->SetMeshElementOnShape(F, geomFaceID);
1506 } // if( (drl+addv) > 0 )
1508 } // end new version implementation
1514 //=============================================================================
1515 /*! Split quadrangle in to 2 triangles by smallest diagonal
1518 //=============================================================================
1519 void StdMeshers_Quadrangle_2D::SplitQuad(SMESHDS_Mesh *theMeshDS,
1521 const SMDS_MeshNode* theNode1,
1522 const SMDS_MeshNode* theNode2,
1523 const SMDS_MeshNode* theNode3,
1524 const SMDS_MeshNode* theNode4)
1526 gp_Pnt a(theNode1->X(),theNode1->Y(),theNode1->Z());
1527 gp_Pnt b(theNode2->X(),theNode2->Y(),theNode2->Z());
1528 gp_Pnt c(theNode3->X(),theNode3->Y(),theNode3->Z());
1529 gp_Pnt d(theNode4->X(),theNode4->Y(),theNode4->Z());
1530 SMDS_MeshFace* face;
1531 if(a.Distance(c) > b.Distance(d)){
1532 face = myTool->AddFace(theNode2, theNode4, theNode1);
1533 theMeshDS->SetMeshElementOnShape(face, theFaceID );
1534 face = myTool->AddFace(theNode2, theNode3, theNode4);
1535 theMeshDS->SetMeshElementOnShape(face, theFaceID );
1539 face = myTool->AddFace(theNode1, theNode2 ,theNode3);
1540 theMeshDS->SetMeshElementOnShape(face, theFaceID );
1541 face = myTool->AddFace(theNode1, theNode3, theNode4);
1542 theMeshDS->SetMeshElementOnShape(face, theFaceID );