// SMESH SMESH : implementaion of SMESH idl descriptions // // Copyright (C) 2003 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN, // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS // // This library is free software; you can redistribute it and/or // modify it under the terms of the GNU Lesser General Public // License as published by the Free Software Foundation; either // version 2.1 of the License. // // This library is distributed in the hope that it will be useful, // but WITHOUT ANY WARRANTY; without even the implied warranty of // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU // Lesser General Public License for more details. // // You should have received a copy of the GNU Lesser General Public // License along with this library; if not, write to the Free Software // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA // // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com // // // // File : StdMeshers_Quadrangle_2D.cxx // Moved here from SMESH_Quadrangle_2D.cxx // Author : Paul RASCLE, EDF // Module : SMESH // $Header$ #include "StdMeshers_Quadrangle_2D.hxx" #include "StdMeshers_FaceSide.hxx" #include "SMESH_Gen.hxx" #include "SMESH_Mesh.hxx" #include "SMESH_subMesh.hxx" #include "SMESH_MesherHelper.hxx" #include "SMESH_Block.hxx" #include "SMESH_Comment.hxx" #include "SMDS_MeshElement.hxx" #include "SMDS_MeshNode.hxx" #include "SMDS_EdgePosition.hxx" #include "SMDS_FacePosition.hxx" #include #include #include #include #include #include #include #include #include #include #include "utilities.h" #include "Utils_ExceptHandlers.hxx" #ifndef StdMeshers_Array2OfNode_HeaderFile #define StdMeshers_Array2OfNode_HeaderFile typedef const SMDS_MeshNode* SMDS_MeshNodePtr; DEFINE_BASECOLLECTION (StdMeshers_BaseCollectionNodePtr, SMDS_MeshNodePtr) DEFINE_ARRAY2(StdMeshers_Array2OfNode, StdMeshers_BaseCollectionNodePtr, SMDS_MeshNodePtr) #endif using namespace std; typedef gp_XY gp_UV; typedef SMESH_Comment TComm; //============================================================================= /*! * */ //============================================================================= StdMeshers_Quadrangle_2D::StdMeshers_Quadrangle_2D (int hypId, int studyId, SMESH_Gen* gen) : SMESH_2D_Algo(hypId, studyId, gen) { MESSAGE("StdMeshers_Quadrangle_2D::StdMeshers_Quadrangle_2D"); _name = "Quadrangle_2D"; _shapeType = (1 << TopAbs_FACE); _compatibleHypothesis.push_back("QuadranglePreference"); myTool = 0; } //============================================================================= /*! * */ //============================================================================= StdMeshers_Quadrangle_2D::~StdMeshers_Quadrangle_2D() { MESSAGE("StdMeshers_Quadrangle_2D::~StdMeshers_Quadrangle_2D"); } //============================================================================= /*! * */ //============================================================================= bool StdMeshers_Quadrangle_2D::CheckHypothesis (SMESH_Mesh& aMesh, const TopoDS_Shape& aShape, SMESH_Hypothesis::Hypothesis_Status& aStatus) { bool isOk = true; aStatus = SMESH_Hypothesis::HYP_OK; // there is only one compatible Hypothesis so far const list &hyps = GetUsedHypothesis(aMesh, aShape, false); myQuadranglePreference = hyps.size() > 0; return isOk; } //============================================================================= /*! * */ //============================================================================= bool StdMeshers_Quadrangle_2D::Compute (SMESH_Mesh& aMesh, const TopoDS_Shape& aShape)// throw (SALOME_Exception) { // PAL14921. Enable catching std::bad_alloc and Standard_OutOfMemory outside //Unexpect aCatchSalomeException); SMESHDS_Mesh * meshDS = aMesh.GetMeshDS(); aMesh.GetSubMesh(aShape); SMESH_MesherHelper helper(aMesh); myTool = &helper; _quadraticMesh = myTool->IsQuadraticSubMesh(aShape); FaceQuadStruct *quad = CheckNbEdges( aMesh, aShape ); std::auto_ptr quadDeleter( quad ); // to delete quad at exit from Compute() if (!quad) return false; if(myQuadranglePreference) { int n1 = quad->side[0]->NbPoints(); int n2 = quad->side[1]->NbPoints(); int n3 = quad->side[2]->NbPoints(); int n4 = quad->side[3]->NbPoints(); int nfull = n1+n2+n3+n4; int ntmp = nfull/2; ntmp = ntmp*2; if( nfull==ntmp && ( (n1!=n3) || (n2!=n4) ) ) { // special path for using only quandrangle faces bool ok = ComputeQuadPref(aMesh, aShape, quad); return ok; } } // set normalized grid on unit square in parametric domain if (!SetNormalizedGrid(aMesh, aShape, quad)) return false; // --- compute 3D values on points, store points & quadrangles int nbdown = quad->side[0]->NbPoints(); int nbup = quad->side[2]->NbPoints(); int nbright = quad->side[1]->NbPoints(); int nbleft = quad->side[3]->NbPoints(); int nbhoriz = Min(nbdown, nbup); int nbvertic = Min(nbright, nbleft); const TopoDS_Face& F = TopoDS::Face(aShape); Handle(Geom_Surface) S = BRep_Tool::Surface(F); // internal mesh nodes int i, j, geomFaceID = meshDS->ShapeToIndex( F ); for (i = 1; i < nbhoriz - 1; i++) { for (j = 1; j < nbvertic - 1; j++) { int ij = j * nbhoriz + i; double u = quad->uv_grid[ij].u; double v = quad->uv_grid[ij].v; gp_Pnt P = S->Value(u, v); SMDS_MeshNode * node = meshDS->AddNode(P.X(), P.Y(), P.Z()); meshDS->SetNodeOnFace(node, geomFaceID, u, v); quad->uv_grid[ij].node = node; } } // mesh faces // [2] // --.--.--.--.--.-- nbvertic // | | ^ // | | ^ // [3] | | ^ j [1] // | | ^ // | | ^ // ---.----.----.--- 0 // 0 > > > > > > > > nbhoriz // i // [0] i = 0; int ilow = 0; int iup = nbhoriz - 1; if (quad->isEdgeOut[3]) { ilow++; } else { if (quad->isEdgeOut[1]) iup--; } int jlow = 0; int jup = nbvertic - 1; if (quad->isEdgeOut[0]) { jlow++; } else { if (quad->isEdgeOut[2]) jup--; } // regular quadrangles for (i = ilow; i < iup; i++) { for (j = jlow; j < jup; j++) { const SMDS_MeshNode *a, *b, *c, *d; a = quad->uv_grid[j * nbhoriz + i].node; b = quad->uv_grid[j * nbhoriz + i + 1].node; c = quad->uv_grid[(j + 1) * nbhoriz + i + 1].node; d = quad->uv_grid[(j + 1) * nbhoriz + i].node; SMDS_MeshFace* face = myTool->AddFace(a, b, c, d); meshDS->SetMeshElementOnShape(face, geomFaceID); } } const vector& uv_e0 = quad->side[0]->GetUVPtStruct(true,0 ); const vector& uv_e1 = quad->side[1]->GetUVPtStruct(false,1); const vector& uv_e2 = quad->side[2]->GetUVPtStruct(true,1 ); const vector& uv_e3 = quad->side[3]->GetUVPtStruct(false,0); if ( uv_e0.empty() || uv_e1.empty() || uv_e2.empty() || uv_e3.empty() ) return error( COMPERR_BAD_INPUT_MESH ); double eps = Precision::Confusion(); // Boundary quadrangles if (quad->isEdgeOut[0]) { // Down edge is out // // |___|___|___|___|___|___| // | | | | | | | // |___|___|___|___|___|___| // | | | | | | | // |___|___|___|___|___|___| __ first row of the regular grid // . . . . . . . . . __ down edge nodes // // >->->->->->->->->->->->-> -- direction of processing int g = 0; // number of last processed node in the regular grid // number of last node of the down edge to be processed int stop = nbdown - 1; // if right edge is out, we will stop at a node, previous to the last one if (quad->isEdgeOut[1]) stop--; // for each node of the down edge find nearest node // in the first row of the regular grid and link them for (i = 0; i < stop; i++) { const SMDS_MeshNode *a, *b, *c, *d; a = uv_e0[i].node; b = uv_e0[i + 1].node; gp_Pnt pb (b->X(), b->Y(), b->Z()); // find node c in the regular grid, which will be linked with node b int near = g; if (i == stop - 1) { // right bound reached, link with the rightmost node near = iup; c = quad->uv_grid[nbhoriz + iup].node; } else { // find in the grid node c, nearest to the b double mind = RealLast(); for (int k = g; k <= iup; k++) { const SMDS_MeshNode *nk; if (k < ilow) // this can be, if left edge is out nk = uv_e3[1].node; // get node from the left edge else nk = quad->uv_grid[nbhoriz + k].node; // get one of middle nodes gp_Pnt pnk (nk->X(), nk->Y(), nk->Z()); double dist = pb.Distance(pnk); if (dist < mind - eps) { c = nk; near = k; mind = dist; } else { break; } } } if (near == g) { // make triangle //SMDS_MeshFace* face = meshDS->AddFace(a, b, c); SMDS_MeshFace* face = myTool->AddFace(a, b, c); meshDS->SetMeshElementOnShape(face, geomFaceID); } else { // make quadrangle if (near - 1 < ilow) d = uv_e3[1].node; else d = quad->uv_grid[nbhoriz + near - 1].node; //SMDS_MeshFace* face = meshDS->AddFace(a, b, c, d); SMDS_MeshFace* face = myTool->AddFace(a, b, c, d); meshDS->SetMeshElementOnShape(face, geomFaceID); // if node d is not at position g - make additional triangles if (near - 1 > g) { for (int k = near - 1; k > g; k--) { c = quad->uv_grid[nbhoriz + k].node; if (k - 1 < ilow) d = uv_e3[1].node; else d = quad->uv_grid[nbhoriz + k - 1].node; //SMDS_MeshFace* face = meshDS->AddFace(a, c, d); SMDS_MeshFace* face = myTool->AddFace(a, c, d); meshDS->SetMeshElementOnShape(face, geomFaceID); } } g = near; } } } else { if (quad->isEdgeOut[2]) { // Up edge is out // // <-<-<-<-<-<-<-<-<-<-<-<-< -- direction of processing // // . . . . . . . . . __ up edge nodes // ___ ___ ___ ___ ___ ___ __ first row of the regular grid // | | | | | | | // |___|___|___|___|___|___| // | | | | | | | // |___|___|___|___|___|___| // | | | | | | | int g = nbhoriz - 1; // last processed node in the regular grid int stop = 0; // if left edge is out, we will stop at a second node if (quad->isEdgeOut[3]) stop++; // for each node of the up edge find nearest node // in the first row of the regular grid and link them for (i = nbup - 1; i > stop; i--) { const SMDS_MeshNode *a, *b, *c, *d; a = uv_e2[i].node; b = uv_e2[i - 1].node; gp_Pnt pb (b->X(), b->Y(), b->Z()); // find node c in the grid, which will be linked with node b int near = g; if (i == stop + 1) { // left bound reached, link with the leftmost node c = quad->uv_grid[nbhoriz*(nbvertic - 2) + ilow].node; near = ilow; } else { // find node c in the grid, nearest to the b double mind = RealLast(); for (int k = g; k >= ilow; k--) { const SMDS_MeshNode *nk; if (k > iup) nk = uv_e1[nbright - 2].node; else nk = quad->uv_grid[nbhoriz*(nbvertic - 2) + k].node; gp_Pnt pnk (nk->X(), nk->Y(), nk->Z()); double dist = pb.Distance(pnk); if (dist < mind - eps) { c = nk; near = k; mind = dist; } else { break; } } } if (near == g) { // make triangle //SMDS_MeshFace* face = meshDS->AddFace(a, b, c); SMDS_MeshFace* face = myTool->AddFace(a, b, c); meshDS->SetMeshElementOnShape(face, geomFaceID); } else { // make quadrangle if (near + 1 > iup) d = uv_e1[nbright - 2].node; else d = quad->uv_grid[nbhoriz*(nbvertic - 2) + near + 1].node; //SMDS_MeshFace* face = meshDS->AddFace(a, b, c, d); SMDS_MeshFace* face = myTool->AddFace(a, b, c, d); meshDS->SetMeshElementOnShape(face, geomFaceID); if (near + 1 < g) { // if d not is at g - make additional triangles for (int k = near + 1; k < g; k++) { c = quad->uv_grid[nbhoriz*(nbvertic - 2) + k].node; if (k + 1 > iup) d = uv_e1[nbright - 2].node; else d = quad->uv_grid[nbhoriz*(nbvertic - 2) + k + 1].node; //SMDS_MeshFace* face = meshDS->AddFace(a, c, d); SMDS_MeshFace* face = myTool->AddFace(a, c, d); meshDS->SetMeshElementOnShape(face, geomFaceID); } } g = near; } } } } // right or left boundary quadrangles if (quad->isEdgeOut[1]) { // MESSAGE("right edge is out"); int g = 0; // last processed node in the grid int stop = nbright - 1; if (quad->isEdgeOut[2]) stop--; for (i = 0; i < stop; i++) { const SMDS_MeshNode *a, *b, *c, *d; a = uv_e1[i].node; b = uv_e1[i + 1].node; gp_Pnt pb (b->X(), b->Y(), b->Z()); // find node c in the grid, nearest to the b int near = g; if (i == stop - 1) { // up bondary reached c = quad->uv_grid[nbhoriz*(jup + 1) - 2].node; near = jup; } else { double mind = RealLast(); for (int k = g; k <= jup; k++) { const SMDS_MeshNode *nk; if (k < jlow) nk = uv_e0[nbdown - 2].node; else nk = quad->uv_grid[nbhoriz*(k + 1) - 2].node; gp_Pnt pnk (nk->X(), nk->Y(), nk->Z()); double dist = pb.Distance(pnk); if (dist < mind - eps) { c = nk; near = k; mind = dist; } else { break; } } } if (near == g) { // make triangle //SMDS_MeshFace* face = meshDS->AddFace(a, b, c); SMDS_MeshFace* face = myTool->AddFace(a, b, c); meshDS->SetMeshElementOnShape(face, geomFaceID); } else { // make quadrangle if (near - 1 < jlow) d = uv_e0[nbdown - 2].node; else d = quad->uv_grid[nbhoriz*near - 2].node; //SMDS_MeshFace* face = meshDS->AddFace(a, b, c, d); SMDS_MeshFace* face = myTool->AddFace(a, b, c, d); meshDS->SetMeshElementOnShape(face, geomFaceID); if (near - 1 > g) { // if d not is at g - make additional triangles for (int k = near - 1; k > g; k--) { c = quad->uv_grid[nbhoriz*(k + 1) - 2].node; if (k - 1 < jlow) d = uv_e0[nbdown - 2].node; else d = quad->uv_grid[nbhoriz*k - 2].node; //SMDS_MeshFace* face = meshDS->AddFace(a, c, d); SMDS_MeshFace* face = myTool->AddFace(a, c, d); meshDS->SetMeshElementOnShape(face, geomFaceID); } } g = near; } } } else { if (quad->isEdgeOut[3]) { // MESSAGE("left edge is out"); int g = nbvertic - 1; // last processed node in the grid int stop = 0; if (quad->isEdgeOut[0]) stop++; for (i = nbleft - 1; i > stop; i--) { const SMDS_MeshNode *a, *b, *c, *d; a = uv_e3[i].node; b = uv_e3[i - 1].node; gp_Pnt pb (b->X(), b->Y(), b->Z()); // find node c in the grid, nearest to the b int near = g; if (i == stop + 1) { // down bondary reached c = quad->uv_grid[nbhoriz*jlow + 1].node; near = jlow; } else { double mind = RealLast(); for (int k = g; k >= jlow; k--) { const SMDS_MeshNode *nk; if (k > jup) nk = uv_e2[1].node; else nk = quad->uv_grid[nbhoriz*k + 1].node; gp_Pnt pnk (nk->X(), nk->Y(), nk->Z()); double dist = pb.Distance(pnk); if (dist < mind - eps) { c = nk; near = k; mind = dist; } else { break; } } } if (near == g) { // make triangle //SMDS_MeshFace* face = meshDS->AddFace(a, b, c); SMDS_MeshFace* face = myTool->AddFace(a, b, c); meshDS->SetMeshElementOnShape(face, geomFaceID); } else { // make quadrangle if (near + 1 > jup) d = uv_e2[1].node; else d = quad->uv_grid[nbhoriz*(near + 1) + 1].node; //SMDS_MeshFace* face = meshDS->AddFace(a, b, c, d); SMDS_MeshFace* face = myTool->AddFace(a, b, c, d); meshDS->SetMeshElementOnShape(face, geomFaceID); if (near + 1 < g) { // if d not is at g - make additional triangles for (int k = near + 1; k < g; k++) { c = quad->uv_grid[nbhoriz*k + 1].node; if (k + 1 > jup) d = uv_e2[1].node; else d = quad->uv_grid[nbhoriz*(k + 1) + 1].node; //SMDS_MeshFace* face = meshDS->AddFace(a, c, d); SMDS_MeshFace* face = myTool->AddFace(a, c, d); meshDS->SetMeshElementOnShape(face, geomFaceID); } } g = near; } } } } bool isOk = true; return isOk; } //============================================================================= /*! * */ //============================================================================= FaceQuadStruct* StdMeshers_Quadrangle_2D::CheckNbEdges(SMESH_Mesh & aMesh, const TopoDS_Shape & aShape) //throw(SALOME_Exception) { const TopoDS_Face & F = TopoDS::Face(aShape); const bool ignoreMediumNodes = _quadraticMesh; // verify 1 wire only, with 4 edges TopoDS_Vertex V; list< TopoDS_Edge > edges; list< int > nbEdgesInWire; int nbWire = SMESH_Block::GetOrderedEdges (F, V, edges, nbEdgesInWire); if (nbWire != 1) { error(COMPERR_BAD_SHAPE, TComm("Wrong number of wires: ") << nbWire); return 0; } FaceQuadStruct* quad = new FaceQuadStruct; quad->uv_grid = 0; quad->side.reserve(nbEdgesInWire.front()); int nbSides = 0; list< TopoDS_Edge >::iterator edgeIt = edges.begin(); if ( nbEdgesInWire.front() == 4 ) { // exactly 4 edges for ( ; edgeIt != edges.end(); ++edgeIt, nbSides++ ) quad->side.push_back( new StdMeshers_FaceSide(F, *edgeIt, &aMesh, nbSides 4 ) { // more than 4 edges - try to unite some list< TopoDS_Edge > sideEdges; while ( !edges.empty()) { sideEdges.clear(); sideEdges.splice( sideEdges.end(), edges, edges.begin()); // edges.front() -> sideEdges.end() bool sameSide = true; while ( !edges.empty() && sameSide ) { sameSide = SMESH_Algo::IsContinuous( sideEdges.back(), edges.front() ); if ( sameSide ) sideEdges.splice( sideEdges.end(), edges, edges.begin()); } if ( nbSides == 0 ) { // go backward from the first edge sameSide = true; while ( !edges.empty() && sameSide ) { sameSide = SMESH_Algo::IsContinuous( sideEdges.front(), edges.back() ); if ( sameSide ) sideEdges.splice( sideEdges.begin(), edges, --edges.end()); } } quad->side.push_back( new StdMeshers_FaceSide(F, sideEdges, &aMesh, nbSidesside[i]->NbEdges(); ++e ) cout << myTool->GetMeshDS()->ShapeToIndex( quad->side[i]->Edge( e )) << " "; cout << ")"; } cout << endl; #endif if ( !nbSides ) nbSides = nbEdgesInWire.front(); error(COMPERR_BAD_SHAPE, TComm("Face must have 4 sides but not ") << nbSides); delete quad; quad = 0; } return quad; } //============================================================================= /*! * CheckAnd2Dcompute */ //============================================================================= FaceQuadStruct *StdMeshers_Quadrangle_2D::CheckAnd2Dcompute (SMESH_Mesh & aMesh, const TopoDS_Shape & aShape, const bool CreateQuadratic) //throw(SALOME_Exception) { _quadraticMesh = CreateQuadratic; FaceQuadStruct *quad = CheckNbEdges(aMesh, aShape); if(!quad) return 0; // set normalized grid on unit square in parametric domain bool stat = SetNormalizedGrid(aMesh, aShape, quad); if(!stat) { if(!quad) delete quad; quad = 0; } return quad; } //============================================================================= /*! * */ //============================================================================= faceQuadStruct::~faceQuadStruct() { for (int i = 0; i < side.size(); i++) { if (side[i]) delete side[i]; } if (uv_grid) delete [] uv_grid; } namespace { inline const vector& GetUVPtStructIn(FaceQuadStruct* quad, int i, int nbSeg) { bool isXConst = ( i == BOTTOM_SIDE || i == TOP_SIDE ); double constValue = ( i == BOTTOM_SIDE || i == LEFT_SIDE ) ? 0 : 1; return quad->isEdgeOut[i] ? quad->side[i]->SimulateUVPtStruct(nbSeg,isXConst,constValue) : quad->side[i]->GetUVPtStruct(isXConst,constValue); } } //============================================================================= /*! * */ //============================================================================= bool StdMeshers_Quadrangle_2D::SetNormalizedGrid (SMESH_Mesh & aMesh, const TopoDS_Shape& aShape, FaceQuadStruct* & quad) //throw (SALOME_Exception) { // Algorithme décrit dans "Génération automatique de maillages" // P.L. GEORGE, MASSON, § 6.4.1 p. 84-85 // traitement dans le domaine paramétrique 2d u,v // transport - projection sur le carré unité // MESSAGE("StdMeshers_Quadrangle_2D::SetNormalizedGrid"); // const TopoDS_Face& F = TopoDS::Face(aShape); // 1 --- find orientation of the 4 edges, by test on extrema // max min 0 x1 1 // |<----north-2-------^ a3 -------------> a2 // | | ^1 1^ // west-3 east-1 =right | | // | | ==> | | // y0 | | y1 | | // | | |0 0| // v----south-0--------> a0 -------------> a1 // min max 0 x0 1 // =down // // 3 --- 2D normalized values on unit square [0..1][0..1] int nbhoriz = Min(quad->side[0]->NbPoints(), quad->side[2]->NbPoints()); int nbvertic = Min(quad->side[1]->NbPoints(), quad->side[3]->NbPoints()); quad->isEdgeOut[0] = (quad->side[0]->NbPoints() > quad->side[2]->NbPoints()); quad->isEdgeOut[1] = (quad->side[1]->NbPoints() > quad->side[3]->NbPoints()); quad->isEdgeOut[2] = (quad->side[2]->NbPoints() > quad->side[0]->NbPoints()); quad->isEdgeOut[3] = (quad->side[3]->NbPoints() > quad->side[1]->NbPoints()); UVPtStruct *uv_grid = quad->uv_grid = new UVPtStruct[nbvertic * nbhoriz]; const vector& uv_e0 = GetUVPtStructIn( quad, 0, nbhoriz - 1 ); const vector& uv_e1 = GetUVPtStructIn( quad, 1, nbvertic - 1 ); const vector& uv_e2 = GetUVPtStructIn( quad, 2, nbhoriz - 1 ); const vector& uv_e3 = GetUVPtStructIn( quad, 3, nbvertic - 1 ); if ( uv_e0.empty() || uv_e1.empty() || uv_e2.empty() || uv_e3.empty() ) //return error( "Can't find nodes on sides"); return error( COMPERR_BAD_INPUT_MESH ); // nodes Id on "in" edges if (! quad->isEdgeOut[0]) { int j = 0; for (int i = 0; i < nbhoriz; i++) { // down int ij = j * nbhoriz + i; uv_grid[ij].node = uv_e0[i].node; } } if (! quad->isEdgeOut[1]) { int i = nbhoriz - 1; for (int j = 0; j < nbvertic; j++) { // right int ij = j * nbhoriz + i; uv_grid[ij].node = uv_e1[j].node; } } if (! quad->isEdgeOut[2]) { int j = nbvertic - 1; for (int i = 0; i < nbhoriz; i++) { // up int ij = j * nbhoriz + i; uv_grid[ij].node = uv_e2[i].node; } } if (! quad->isEdgeOut[3]) { int i = 0; for (int j = 0; j < nbvertic; j++) { // left int ij = j * nbhoriz + i; uv_grid[ij].node = uv_e3[j].node; } } // normalized 2d values on grid for (int i = 0; i < nbhoriz; i++) { for (int j = 0; j < nbvertic; j++) { int ij = j * nbhoriz + i; // --- droite i cste : x = x0 + y(x1-x0) double x0 = uv_e0[i].normParam; // bas - sud double x1 = uv_e2[i].normParam; // haut - nord // --- droite j cste : y = y0 + x(y1-y0) double y0 = uv_e3[j].normParam; // gauche-ouest double y1 = uv_e1[j].normParam; // droite - est // --- intersection : x=x0+(y0+x(y1-y0))(x1-x0) double x = (x0 + y0 * (x1 - x0)) / (1 - (y1 - y0) * (x1 - x0)); double y = y0 + x * (y1 - y0); uv_grid[ij].x = x; uv_grid[ij].y = y; //MESSAGE("-xy-01 "<side[0]->Value2d(param_0).XY(); gp_UV p1 = quad->side[1]->Value2d(param_1).XY(); gp_UV p2 = quad->side[2]->Value2d(param_2).XY(); gp_UV p3 = quad->side[3]->Value2d(param_3).XY(); gp_UV uv = (1 - y) * p0 + x * p1 + y * p2 + (1 - x) * p3; uv -= (1 - x) * (1 - y) * a0 + x * (1 - y) * a1 + x * y * a2 + (1 - x) * y * a3; uv_grid[ij].u = uv.X(); uv_grid[ij].v = uv.Y(); } } return true; } //======================================================================= //function : ShiftQuad //purpose : auxilary function for ComputeQuadPref //======================================================================= static void ShiftQuad(FaceQuadStruct* quad, const int num, bool) { StdMeshers_FaceSide* side[4] = { quad->side[0], quad->side[1], quad->side[2], quad->side[3] }; for (int i = BOTTOM_SIDE; i < NB_SIDES; ++i ) { int id = ( i + num ) % NB_SIDES; bool wasForward = ( i < TOP_SIDE ); bool newForward = ( id < TOP_SIDE ); if ( wasForward != newForward ) side[ i ]->Reverse(); quad->side[ id ] = side[ i ]; } } //======================================================================= //function : CalcUV //purpose : auxilary function for ComputeQuadPref //======================================================================= static gp_UV CalcUV(double x0, double x1, double y0, double y1, FaceQuadStruct* quad, const gp_UV& a0, const gp_UV& a1, const gp_UV& a2, const gp_UV& a3) { const vector& uv_eb = quad->side[0]->GetUVPtStruct(true,0 ); const vector& uv_er = quad->side[1]->GetUVPtStruct(false,1); const vector& uv_et = quad->side[2]->GetUVPtStruct(true,1 ); const vector& uv_el = quad->side[3]->GetUVPtStruct(false,0); double x = (x0 + y0 * (x1 - x0)) / (1 - (y1 - y0) * (x1 - x0)); double y = y0 + x * (y1 - y0); double param_b = uv_eb[0].normParam + x * (uv_eb.back().normParam - uv_eb[0].normParam); double param_t = uv_et[0].normParam + x * (uv_et.back().normParam - uv_et[0].normParam); double param_r = uv_er[0].normParam + y * (uv_er.back().normParam - uv_er[0].normParam); double param_l = uv_el[0].normParam + y * (uv_el.back().normParam - uv_el[0].normParam); gp_UV p0 = quad->side[BOTTOM_SIDE]->Value2d(param_b).XY(); gp_UV p1 = quad->side[RIGHT_SIDE ]->Value2d(param_r).XY(); gp_UV p2 = quad->side[TOP_SIDE ]->Value2d(param_t).XY(); gp_UV p3 = quad->side[LEFT_SIDE ]->Value2d(param_l).XY(); gp_UV uv = p0 * (1 - y) + p1 * x + p2 * y + p3 * (1 - x); uv -= (1 - x) * (1 - y) * a0 + x * (1 - y) * a1 + x * y * a2 + (1 - x) * y * a3; return uv; } //======================================================================= /*! * Create only quandrangle faces */ //======================================================================= bool StdMeshers_Quadrangle_2D::ComputeQuadPref (SMESH_Mesh & aMesh, const TopoDS_Shape& aShape, FaceQuadStruct* quad) { SMESHDS_Mesh * meshDS = aMesh.GetMeshDS(); const TopoDS_Face& F = TopoDS::Face(aShape); Handle(Geom_Surface) S = BRep_Tool::Surface(F); // const TopoDS_Wire& W = BRepTools::OuterWire(F); bool WisF = true; // if(W.Orientation()==TopAbs_FORWARD) // WisF = true; //if(WisF) cout<<"W is FORWARD"<ShapeToIndex( F ); int nb = quad->side[0]->NbPoints(); int nr = quad->side[1]->NbPoints(); int nt = quad->side[2]->NbPoints(); int nl = quad->side[3]->NbPoints(); int dh = abs(nb-nt); int dv = abs(nr-nl); if( dh>=dv ) { if( nt>nb ) { // it is a base case => not shift quad but me be replacement is need ShiftQuad(quad,0,WisF); } else { // we have to shift quad on 2 ShiftQuad(quad,2,WisF); } } else { if( nr>nl ) { // we have to shift quad on 1 ShiftQuad(quad,1,WisF); } else { // we have to shift quad on 3 ShiftQuad(quad,3,WisF); } } nb = quad->side[0]->NbPoints(); nr = quad->side[1]->NbPoints(); nt = quad->side[2]->NbPoints(); nl = quad->side[3]->NbPoints(); dh = abs(nb-nt); dv = abs(nr-nl); int nbh = Max(nb,nt); int nbv = Max(nr,nl); int addh = 0; int addv = 0; // orientation of face and 3 main domain for future faces // 0 top 1 // 1------------1 // | | | | // | | | | // | L | | R | // left | | | | rigth // | / \ | // | / C \ | // |/ \| // 0------------0 // 0 bottom 1 if(dh>dv) { addv = (dh-dv)/2; nbv = nbv + addv; } else { // dv>=dh addh = (dv-dh)/2; nbh = nbh + addh; } const vector& uv_eb = quad->side[0]->GetUVPtStruct(true,0 ); const vector& uv_er = quad->side[1]->GetUVPtStruct(false,1); const vector& uv_et = quad->side[2]->GetUVPtStruct(true,1 ); const vector& uv_el = quad->side[3]->GetUVPtStruct(false,0); // arrays for normalized params //cout<<"Dump B:"<X()<<","<Y()<<","<Z()<<")"<0) { // add top nodes for(i=1; i<=dl; i++) NodesL.SetValue(i+1,nl,uv_et[i].node); // create and add needed nodes TColgp_SequenceOfXY UVtmp; for(i=1; i<=dl; i++) { double x0 = npt.Value(i+1); double x1 = x0; // diagonal node double y0 = npl.Value(i+1); double y1 = npr.Value(i+1); gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3); gp_Pnt P = S->Value(UV.X(),UV.Y()); SMDS_MeshNode * N = meshDS->AddNode(P.X(), P.Y(), P.Z()); meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y()); NodesL.SetValue(i+1,1,N); if(UVL.Length()Value(UV.X(),UV.Y()); SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z()); meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y()); NodesL.SetValue(i+1,j,N); if( i==dl ) UVtmp.Append(UV); } } for(i=1; i<=UVtmp.Length() && UVL.Length()X()<<","<Y()<<","<Z()<<")"; // } // cout<AddFace(NodesL.Value(i,j), NodesL.Value(i+1,j), NodesL.Value(i+1,j+1), NodesL.Value(i,j+1)); meshDS->SetMeshElementOnShape(F, geomFaceID); } else { SMDS_MeshFace* F = myTool->AddFace(NodesL.Value(i,j), NodesL.Value(i,j+1), NodesL.Value(i+1,j+1), NodesL.Value(i+1,j)); meshDS->SetMeshElementOnShape(F, geomFaceID); } } } } else { // fill UVL using c2d for(i=1; i0) { // add top nodes for(i=1; i<=dr; i++) NodesR.SetValue(i+1,1,uv_et[nt-1-i].node); // create and add needed nodes TColgp_SequenceOfXY UVtmp; for(i=1; i<=dr; i++) { double x0 = npt.Value(nt-i); double x1 = x0; // diagonal node double y0 = npl.Value(i+1); double y1 = npr.Value(i+1); gp_UV UV = CalcUV(x0, x1, y0, y1, quad, a0, a1, a2, a3); gp_Pnt P = S->Value(UV.X(),UV.Y()); SMDS_MeshNode * N = meshDS->AddNode(P.X(), P.Y(), P.Z()); meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y()); NodesR.SetValue(i+1,nr,N); if(UVR.Length()Value(UV.X(),UV.Y()); SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z()); meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y()); NodesR.SetValue(i+1,j,N); if( i==dr ) UVtmp.Prepend(UV); } } for(i=1; i<=UVtmp.Length() && UVR.Length()AddFace(NodesR.Value(i,j), NodesR.Value(i+1,j), NodesR.Value(i+1,j+1), NodesR.Value(i,j+1)); meshDS->SetMeshElementOnShape(F, geomFaceID); } else { SMDS_MeshFace* F = myTool->AddFace(NodesR.Value(i,j), NodesR.Value(i,j+1), NodesR.Value(i+1,j+1), NodesR.Value(i+1,j)); meshDS->SetMeshElementOnShape(F, geomFaceID); } } } } else { // fill UVR using c2d for(i=1; iValue(UV.X(),UV.Y()); SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z()); meshDS->SetNodeOnFace(N, geomFaceID, UV.X(), UV.Y()); NodesC.SetValue(i,nbv-nnn+j,N); } } // add diagonal layers //cout<<"UVL.Length()="<Value(u,v); SMDS_MeshNode* N = meshDS->AddNode(P.X(), P.Y(), P.Z()); meshDS->SetNodeOnFace(N, geomFaceID, u, v); NodesC.SetValue(j,i+1,N); } } // create faces for(i=1; iAddFace(NodesC.Value(i,j), NodesC.Value(i+1,j), NodesC.Value(i+1,j+1), NodesC.Value(i,j+1)); meshDS->SetMeshElementOnShape(F, geomFaceID); } else { SMDS_MeshFace* F = myTool->AddFace(NodesC.Value(i,j), NodesC.Value(i,j+1), NodesC.Value(i+1,j+1), NodesC.Value(i+1,j)); meshDS->SetMeshElementOnShape(F, geomFaceID); } } } bool isOk = true; return isOk; }