#include <Precision.hxx>
#include <gp_Pnt2d.hxx>
#include <TColStd_ListIteratorOfListOfInteger.hxx>
+#include <TColStd_SequenceOfReal.hxx>
+#include <TColgp_SequenceOfXY.hxx>
#include "utilities.h"
#include "Utils_ExceptHandlers.hxx"
+#ifndef StdMeshers_Array2OfNode_HeaderFile
+#define StdMeshers_Array2OfNode_HeaderFile
+typedef const SMDS_MeshNode* SMDS_MeshNodePtr;
+#include <NCollection_DefineArray2.hxx>
+DEFINE_BASECOLLECTION (StdMeshers_BaseCollectionNodePtr, SMDS_MeshNodePtr)
+DEFINE_ARRAY2(StdMeshers_Array2OfNode,
+ StdMeshers_BaseCollectionNodePtr, SMDS_MeshNodePtr)
+#endif
+
//=============================================================================
/*!
SMESHDS_Mesh * meshDS = aMesh.GetMeshDS();
aMesh.GetSubMesh(aShape);
- FaceQuadStruct *quad = CheckAnd2Dcompute(aMesh, aShape);
+ //FaceQuadStruct *quad = CheckAnd2Dcompute(aMesh, aShape);
+ FaceQuadStruct* quad = CheckNbEdges(aMesh, aShape);
+
+ if (!quad)
+ return false;
+
+ if(myQuadranglePreference) {
+ int n1 = quad->nbPts[0];
+ int n2 = quad->nbPts[1];
+ int n3 = quad->nbPts[2];
+ int n4 = quad->nbPts[3];
+ 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
+ return ComputeQuadPref(aMesh, aShape, quad);
+ }
+ }
+
+ // set normalized grid on unit square in parametric domain
+ SetNormalizedGrid(aMesh, aShape, quad);
if (!quad)
return false;
return isOk;
}
+
//=============================================================================
/*!
*
*/
//=============================================================================
-FaceQuadStruct *StdMeshers_Quadrangle_2D::CheckAnd2Dcompute
- (SMESH_Mesh & aMesh, const TopoDS_Shape & aShape) throw(SALOME_Exception)
+FaceQuadStruct* StdMeshers_Quadrangle_2D::CheckNbEdges(SMESH_Mesh & aMesh,
+ const TopoDS_Shape & aShape)
+ throw(SALOME_Exception)
{
Unexpect aCatch(SalomeException);
// verify 1 wire only, with 4 edges
- if (NumberOfWires(F) != 1)
- {
+ if (NumberOfWires(F) != 1) {
INFOS("only 1 wire by face (quadrangles)");
return 0;
}
const TopoDS_Wire& W = BRepTools::OuterWire(F);
BRepTools_WireExplorer wexp (W, F);
- FaceQuadStruct *quad = new FaceQuadStruct;
+ FaceQuadStruct* quad = new FaceQuadStruct;
for (int i = 0; i < 4; i++)
quad->uv_edges[i] = 0;
quad->uv_grid = 0;
int nbEdges = 0;
- for (wexp.Init(W, F); wexp.More(); wexp.Next())
- {
+ for (wexp.Init(W, F); wexp.More(); wexp.Next()) {
const TopoDS_Edge& E = wexp.Current();
int nb = aMesh.GetSubMesh(E)->GetSubMeshDS()->NbNodes();
- if (nbEdges < 4)
- {
+ if (nbEdges < 4) {
quad->edge[nbEdges] = E;
quad->nbPts[nbEdges] = nb + 2; // internal points + 2 extrema
}
nbEdges++;
}
- if (nbEdges != 4)
- {
+ if (nbEdges != 4) {
INFOS("face must have 4 edges /quadrangles");
QuadDelete(quad);
return 0;
}
- // set normalized grid on unit square in parametric domain
+ return quad;
+}
+
+
+//=============================================================================
+/*!
+ *
+ */
+//=============================================================================
+
+FaceQuadStruct *StdMeshers_Quadrangle_2D::CheckAnd2Dcompute
+ (SMESH_Mesh & aMesh, const TopoDS_Shape & aShape) throw(SALOME_Exception)
+{
+ Unexpect aCatch(SalomeException);
+
+ FaceQuadStruct *quad = CheckNbEdges(aMesh, aShape);
+
+ if(!quad) return 0;
- SetNormalizedGrid(aMesh, F, quad);
+ // set normalized grid on unit square in parametric domain
+ SetNormalizedGrid(aMesh, aShape, quad);
return quad;
}
}
}
+
+//=======================================================================
+//function : ShiftQuad
+//purpose : auxilary function for ComputeQuadPref
+//=======================================================================
+static void ShiftQuad(FaceQuadStruct* quad, const int num, bool WisF)
+{
+ if(num>3) return;
+ int i;
+ for(i=1; i<=num; i++) {
+ int nbPts3 = quad->nbPts[0];
+ quad->nbPts[0] = quad->nbPts[1];
+ quad->nbPts[1] = quad->nbPts[2];
+ quad->nbPts[2] = quad->nbPts[3];
+ quad->nbPts[3] = nbPts3;
+ TopoDS_Edge edge3 = quad->edge[0];
+ quad->edge[0] = quad->edge[1];
+ quad->edge[1] = quad->edge[2];
+ quad->edge[2] = quad->edge[3];
+ quad->edge[3] = edge3;
+ double first3 = quad->first[0];
+ quad->first[0] = quad->first[1];
+ quad->first[1] = quad->first[2];
+ quad->first[2] = quad->first[3];
+ quad->first[3] = first3;
+ double last3 = quad->last[0];
+ quad->last[0] = quad->last[1];
+ quad->last[1] = quad->last[2];
+ quad->last[2] = quad->last[3];
+ quad->last[3] = last3;
+ bool isEdgeForward3 = quad->isEdgeForward[0];
+ quad->isEdgeForward[0] = quad->isEdgeForward[1];
+ quad->isEdgeForward[1] = quad->isEdgeForward[2];
+ quad->isEdgeForward[2] = quad->isEdgeForward[3];
+ quad->isEdgeForward[3] = isEdgeForward3;
+ bool isEdgeOut3 = quad->isEdgeOut[0];
+ quad->isEdgeOut[0] = quad->isEdgeOut[1];
+ quad->isEdgeOut[1] = quad->isEdgeOut[2];
+ quad->isEdgeOut[2] = quad->isEdgeOut[3];
+ quad->isEdgeOut[3] = isEdgeOut3;
+ UVPtStruct* uv_edges3 = quad->uv_edges[0];
+ quad->uv_edges[0] = quad->uv_edges[1];
+ quad->uv_edges[1] = quad->uv_edges[2];
+ quad->uv_edges[2] = quad->uv_edges[3];
+ quad->uv_edges[3] = uv_edges3;
+ }
+ if(!WisF) {
+ // replacement left and right edges
+ int nbPts3 = quad->nbPts[1];
+ quad->nbPts[1] = quad->nbPts[3];
+ quad->nbPts[3] = nbPts3;
+ TopoDS_Edge edge3 = quad->edge[1];
+ quad->edge[1] = quad->edge[3];
+ quad->edge[3] = edge3;
+ double first3 = quad->first[1];
+ quad->first[1] = quad->first[3];
+ quad->first[3] = first3;
+ double last3 = quad->last[1];
+ quad->last[1] = quad->last[2];
+ quad->last[3] = last3;
+ bool isEdgeForward3 = quad->isEdgeForward[1];
+ quad->isEdgeForward[1] = quad->isEdgeForward[3];
+ quad->isEdgeForward[3] = isEdgeForward3;
+ bool isEdgeOut3 = quad->isEdgeOut[1];
+ quad->isEdgeOut[1] = quad->isEdgeOut[3];
+ quad->isEdgeOut[3] = isEdgeOut3;
+ UVPtStruct* uv_edges3 = quad->uv_edges[1];
+ quad->uv_edges[1] = quad->uv_edges[3];
+ quad->uv_edges[3] = uv_edges3;
+ }
+}
+
+
+//=======================================================================
+//function : CalcUV
+//purpose : auxilary function for ComputeQuadPref
+//=======================================================================
+static gp_XY CalcUV(double x0, double x1, double y0, double y1,
+ FaceQuadStruct* quad,
+ const gp_Pnt2d& a0, const gp_Pnt2d& a1,
+ const gp_Pnt2d& a2, const gp_Pnt2d& a3,
+ const Handle(Geom2d_Curve)& c2db,
+ const Handle(Geom2d_Curve)& c2dr,
+ const Handle(Geom2d_Curve)& c2dt,
+ const Handle(Geom2d_Curve)& c2dl)
+{
+ int nb = quad->nbPts[0];
+ int nr = quad->nbPts[1];
+ int nt = quad->nbPts[2];
+ int nl = quad->nbPts[3];
+
+ UVPtStruct* uv_eb = quad->uv_edges[0];
+ UVPtStruct* uv_er = quad->uv_edges[1];
+ UVPtStruct* uv_et = quad->uv_edges[2];
+ UVPtStruct* uv_el = quad->uv_edges[3];
+
+ double x = (x0 + y0 * (x1 - x0)) / (1 - (y1 - y0) * (x1 - x0));
+ double y = y0 + x * (y1 - y0);
+
+ double param_b = uv_eb[0].param + x * (uv_eb[nb-1].param - uv_eb[0].param);
+ double param_t = uv_et[0].param + x * (uv_et[nt-1].param - uv_et[0].param);
+ double param_r = uv_er[0].param + y * (uv_er[nr-1].param - uv_er[0].param);
+ double param_l = uv_el[0].param + y * (uv_el[nl-1].param - uv_el[0].param);
+
+ gp_Pnt2d p0 = c2db->Value(param_b);
+ gp_Pnt2d p1 = c2dr->Value(param_r);
+ gp_Pnt2d p2 = c2dt->Value(param_t);
+ gp_Pnt2d p3 = c2dl->Value(param_l);
+
+ double u = (1 - y) * p0.X() + x * p1.X() + y * p2.X() + (1 - x) * p3.X();
+ double v = (1 - y) * p0.Y() + x * p1.Y() + y * p2.Y() + (1 - x) * p3.Y();
+
+ u -= (1 - x) * (1 - y) * a0.X() + x * (1 - y) * a1.X() +
+ x * y * a2.X() + (1 - x) * y * a3.X();
+ v -= (1 - x) * (1 - y) * a0.Y() + x * (1 - y) * a1.Y() +
+ x * y * a2.Y() + (1 - x) * y * a3.Y();
+
+ //cout<<"x0="<<x0<<" x1="<<x1<<" y0="<<y0<<" y1="<<y1<<endl;
+ //cout<<"x="<<x<<" y="<<y<<endl;
+ //cout<<"param_b="<<param_b<<" param_t="<<param_t<<" param_r="<<param_r<<" param_l="<<param_l<<endl;
+ //cout<<"u="<<u<<" v="<<v<<endl;
+
+ return gp_XY(u,v);
+}
+
+
+//=======================================================================
+//function : ComputeQuadPref
+//purpose :
+//=======================================================================
+/*!
+ * Special function for creation only quandrangle faces
+ */
+bool StdMeshers_Quadrangle_2D::ComputeQuadPref
+ (SMESH_Mesh & aMesh,
+ const TopoDS_Shape& aShape,
+ FaceQuadStruct* quad) throw (SALOME_Exception)
+{
+ Unexpect aCatch(SalomeException);
+
+ 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 = false;
+ if(W.Orientation()==TopAbs_FORWARD)
+ WisF = true;
+ //if(WisF) cout<<"W is FORWARD"<<endl;
+ //else cout<<"W is REVERSED"<<endl;
+ bool FisF = (F.Orientation()==TopAbs_FORWARD);
+ if(!FisF) WisF = !WisF;
+ int i,j,geomFaceID = meshDS->ShapeToIndex( F );
+
+ int nb = quad->nbPts[0];
+ int nr = quad->nbPts[1];
+ int nt = quad->nbPts[2];
+ int nl = quad->nbPts[3];
+ 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 3
+ ShiftQuad(quad,3,WisF);
+ }
+ else {
+ // we have to shift quad on 1
+ ShiftQuad(quad,1,WisF);
+ }
+ }
+
+ nb = quad->nbPts[0];
+ nr = quad->nbPts[1];
+ nt = quad->nbPts[2];
+ nl = quad->nbPts[3];
+ 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;
+ }
+
+ Handle(Geom2d_Curve) c2d[4];
+ for(i=0; i<4; i++) {
+ c2d[i] = BRep_Tool::CurveOnSurface(quad->edge[i], F,
+ quad->first[i], quad->last[i]);
+ }
+
+ bool loadOk = true;
+ for(i=0; i<2; i++) {
+ quad->uv_edges[i] = LoadEdgePoints2(aMesh, F, quad->edge[i], false);
+ if(!quad->uv_edges[i]) loadOk = false;
+ }
+ for(i=2; i<4; i++) {
+ quad->uv_edges[i] = LoadEdgePoints2(aMesh, F, quad->edge[i], true);
+ if (!quad->uv_edges[i]) loadOk = false;
+ }
+ if (!loadOk) {
+ INFOS("StdMeshers_Quadrangle_2D::ComputeQuadPref - LoadEdgePoints failed");
+ QuadDelete( quad );
+ quad = 0;
+ return false;
+ }
+
+ UVPtStruct* uv_eb = quad->uv_edges[0];
+ UVPtStruct* uv_er = quad->uv_edges[1];
+ UVPtStruct* uv_et = quad->uv_edges[2];
+ UVPtStruct* uv_el = quad->uv_edges[3];
+
+ // arrays for normalized params
+ //cout<<"Dump B:"<<endl;
+ TColStd_SequenceOfReal npb, npr, npt, npl;
+ for(i=0; i<nb; i++) {
+ npb.Append(uv_eb[i].normParam);
+ //cout<<"i="<<i<<" par="<<uv_eb[i].param<<" npar="<<uv_eb[i].normParam;
+ //const SMDS_MeshNode* N = uv_eb[i].node;
+ //cout<<" node("<<N->X()<<","<<N->Y()<<","<<N->Z()<<")"<<endl;
+ }
+ for(i=0; i<nr; i++) {
+ npr.Append(uv_er[i].normParam);
+ }
+ for(i=0; i<nt; i++) {
+ npt.Append(uv_et[i].normParam);
+ }
+ for(i=0; i<nl; i++) {
+ npl.Append(uv_el[i].normParam);
+ }
+
+ // we have to add few values of params to right and left
+ // insert them after first param
+ // insert to right
+ int dr = nbv - nr;
+ double dpr = (npr.Value(2) - npr.Value(1))/(dr+1);
+ for(i=1; i<=dr; i++) {
+ npr.InsertAfter(1,npr.Value(2)-dpr);
+ }
+ // insert to left
+ int dl = nbv - nl;
+ dpr = (npl.Value(2) - npl.Value(1))/(dl+1);
+ for(i=1; i<=dl; i++) {
+ npl.InsertAfter(1,npl.Value(2)-dpr);
+ }
+ //cout<<"npb:";
+ //for(i=1; i<=npb.Length(); i++) {
+ // cout<<" "<<npb.Value(i);
+ //}
+ //cout<<endl;
+
+ gp_Pnt2d a[4];
+ c2d[0]->D0(uv_eb[0].param,a[0]);
+ c2d[0]->D0(uv_eb[nb-1].param,a[1]);
+ c2d[2]->D0(uv_et[nt-1].param,a[2]);
+ c2d[2]->D0(uv_et[0].param,a[3]);
+ //cout<<" a[0]("<<a[0].X()<<","<<a[0].Y()<<")"<<" a[1]("<<a[1].X()<<","<<a[1].Y()<<")"
+ // <<" a[2]("<<a[2].X()<<","<<a[2].Y()<<")"<<" a[3]("<<a[3].X()<<","<<a[3].Y()<<")"<<endl;
+
+ int nnn = Min(nr,nl);
+ // auxilary sequence of XY for creation nodes
+ // in the bottom part of central domain
+ // it's length must be == nbv-nnn-1
+ TColgp_SequenceOfXY UVL;
+ TColgp_SequenceOfXY UVR;
+
+ // step1: create faces for left domain
+ StdMeshers_Array2OfNode NodesL(1,dl+1,1,nl);
+ // add left nodes
+ for(j=1; j<=nl; j++)
+ NodesL.SetValue(1,j,uv_el[j-1].node);
+ if(dl>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_XY UV = CalcUV(x0, x1, y0, y1, quad, a[0], a[1], a[2], a[3],
+ c2d[0], c2d[1], c2d[2], c2d[3]);
+ 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()<nbv-nnn-1) UVL.Append(UV);
+ // internal nodes
+ for(j=2; j<nl; j++) {
+ double y0 = npl.Value(dl+j);
+ double y1 = npr.Value(dl+j);
+ gp_XY UV = CalcUV(x0, x1, y0, y1, quad, a[0], a[1], a[2], a[3],
+ c2d[0], c2d[1], c2d[2], c2d[3]);
+ 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,j,N);
+ if( i==dl ) UVtmp.Append(UV);
+ }
+ }
+ for(i=1; i<=UVtmp.Length() && UVL.Length()<nbv-nnn-1; i++) {
+ UVL.Append(UVtmp.Value(i));
+ }
+ //cout<<"Dump NodesL:"<<endl;
+ //for(i=1; i<=dl+1; i++) {
+ // cout<<"i="<<i;
+ // for(j=1; j<=nl; j++) {
+ // cout<<" ("<<NodesL.Value(i,j)->X()<<","<<NodesL.Value(i,j)->Y()<<","<<NodesL.Value(i,j)->Z()<<")";
+ // }
+ // cout<<endl;
+ //}
+ // create faces
+ for(i=1; i<=dl; i++) {
+ for(j=1; j<nl; j++) {
+ if(WisF) {
+ SMDS_MeshFace* F =
+ meshDS->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 =
+ meshDS->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; i<npl.Length() && UVL.Length()<nbv-nnn-1; i++) {
+ gp_Pnt2d p2d;
+ c2d[3]->D0(uv_el[i].param,p2d);
+ UVL.Append(p2d.XY());
+ }
+ }
+
+ // step2: create faces for right domain
+ StdMeshers_Array2OfNode NodesR(1,dr+1,1,nr);
+ // add right nodes
+ for(j=1; j<=nr; j++)
+ NodesR.SetValue(1,j,uv_er[nr-j].node);
+ if(dr>0) {
+ // 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_XY UV = CalcUV(x0, x1, y0, y1, quad, a[0], a[1], a[2], a[3],
+ c2d[0], c2d[1], c2d[2], c2d[3]);
+ 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()<nbv-nnn-1) UVR.Append(UV);
+ // internal nodes
+ for(j=2; j<nr; j++) {
+ double y0 = npl.Value(nbv-j+1);
+ double y1 = npr.Value(nbv-j+1);
+ gp_XY UV = CalcUV(x0, x1, y0, y1, quad, a[0], a[1], a[2], a[3],
+ c2d[0], c2d[1], c2d[2], c2d[3]);
+ 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,j,N);
+ if( i==dr ) UVtmp.Prepend(UV);
+ }
+ }
+ for(i=1; i<=UVtmp.Length() && UVR.Length()<nbv-nnn-1; i++) {
+ UVR.Append(UVtmp.Value(i));
+ }
+ // create faces
+ for(i=1; i<=dr; i++) {
+ for(j=1; j<nr; j++) {
+ if(WisF) {
+ SMDS_MeshFace* F =
+ meshDS->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 =
+ meshDS->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; i<npr.Length() && UVR.Length()<nbv-nnn-1; i++) {
+ gp_Pnt2d p2d;
+ c2d[1]->D0(uv_er[i].param,p2d);
+ UVR.Append(p2d.XY());
+ }
+ }
+
+ // step3: create faces for central domain
+ StdMeshers_Array2OfNode NodesC(1,nb,1,nbv);
+ // add first string using NodesL
+ for(i=1; i<=dl+1; i++)
+ NodesC.SetValue(1,i,NodesL(i,1));
+ for(i=2; i<=nl; i++)
+ NodesC.SetValue(1,dl+i,NodesL(dl+1,i));
+ // add last string using NodesR
+ for(i=1; i<=dr+1; i++)
+ NodesC.SetValue(nb,i,NodesR(i,nr));
+ for(i=1; i<nr; i++)
+ NodesC.SetValue(nb,dr+i+1,NodesR(dr+1,nr-i));
+ // add top nodes (last columns)
+ for(i=dl+2; i<nbh-dr; i++)
+ NodesC.SetValue(i-dl,nbv,uv_et[i-1].node);
+ // add bottom nodes (first columns)
+ for(i=2; i<nb; i++) {
+ NodesC.SetValue(i,1,uv_eb[i-1].node);
+ gp_Pnt2d p2d;
+ c2d[0]->D0(uv_eb[i-1].param,p2d);
+ }
+ // create and add needed nodes
+ // add linear layers
+ for(i=2; i<nb; i++) {
+ double x0 = npt.Value(dl+i);
+ double x1 = x0;
+ for(j=1; j<nnn; j++) {
+ double y0 = npl.Value(nbv-nnn+j);
+ double y1 = npr.Value(nbv-nnn+j);
+ gp_XY UV = CalcUV(x0, x1, y0, y1, quad, a[0], a[1], a[2], a[3],
+ c2d[0], c2d[1], c2d[2], c2d[3]);
+ 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());
+ NodesC.SetValue(i,nbv-nnn+j,N);
+ }
+ }
+ // add diagonal layers
+ //cout<<"UVL.Length()="<<UVL.Length()<<" UVR.Length()="<<UVR.Length()<<endl;
+ //cout<<"Dump UVL:"<<endl;
+ //for(i=1; i<=UVL.Length(); i++) {
+ // cout<<" ("<<UVL.Value(i).X()<<","<<UVL.Value(i).Y()<<")";
+ //}
+ //cout<<endl;
+ for(i=1; i<nbv-nnn; i++) {
+ double du = UVR.Value(i).X() - UVL.Value(i).X();
+ double dv = UVR.Value(i).Y() - UVL.Value(i).Y();
+ for(j=2; j<nb; j++) {
+ double u = UVL.Value(i).X() + du*npb.Value(j);
+ double v = UVL.Value(i).Y() + dv*npb.Value(j);
+ gp_Pnt P = S->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; i<nb; i++) {
+ for(j=1; j<nbv; j++) {
+ if(WisF) {
+ SMDS_MeshFace* F =
+ meshDS->AddFace(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 =
+ meshDS->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);
+ }
+ }
+ }
+
+ QuadDelete(quad);
+ bool isOk = true;
+ return isOk;
+}
+
+
+//=============================================================================
+/*!
+ * LoadEdgePoints2
+ */
+//=============================================================================
+UVPtStruct* StdMeshers_Quadrangle_2D::LoadEdgePoints2 (SMESH_Mesh & aMesh,
+ const TopoDS_Face& F,
+ const TopoDS_Edge& E,
+ bool IsReverse)
+{
+ //MESSAGE("StdMeshers_Quadrangle_2D::LoadEdgePoints");
+ // --- IDNodes of first and last Vertex
+ TopoDS_Vertex VFirst, VLast;
+ TopExp::Vertices(E, VFirst, VLast); // corresponds to f and l
+
+ ASSERT(!VFirst.IsNull());
+ SMDS_NodeIteratorPtr lid = aMesh.GetSubMesh(VFirst)->GetSubMeshDS()->GetNodes();
+ if (!lid->more()) {
+ MESSAGE ( "NO NODE BUILT ON VERTEX" );
+ return 0;
+ }
+ const SMDS_MeshNode* idFirst = lid->next();
+
+ ASSERT(!VLast.IsNull());
+ lid = aMesh.GetSubMesh(VLast)->GetSubMeshDS()->GetNodes();
+ if (!lid->more()) {
+ MESSAGE ( "NO NODE BUILT ON VERTEX" );
+ return 0;
+ }
+ const SMDS_MeshNode* idLast = lid->next();
+
+ // --- edge internal IDNodes (relies on good order storage, not checked)
+
+ map<double, const SMDS_MeshNode *> params;
+ SMDS_NodeIteratorPtr ite = aMesh.GetSubMesh(E)->GetSubMeshDS()->GetNodes();
+
+ while(ite->more()) {
+ const SMDS_MeshNode* node = ite->next();
+ const SMDS_EdgePosition* epos =
+ static_cast<const SMDS_EdgePosition*>(node->GetPosition().get());
+ double param = epos->GetUParameter();
+ params[param] = node;
+ }
+
+ int nbPoints = aMesh.GetSubMesh(E)->GetSubMeshDS()->NbNodes();
+ if (nbPoints != params.size()) {
+ MESSAGE( "BAD NODE ON EDGE POSITIONS" );
+ return 0;
+ }
+ UVPtStruct* uvslf = new UVPtStruct[nbPoints + 2];
+
+ double f, l;
+ Handle(Geom2d_Curve) C2d = BRep_Tool::CurveOnSurface(E, F, f, l);
+
+ const TopoDS_Wire& W = BRepTools::OuterWire(F);
+ bool FisF = (F.Orientation()==TopAbs_FORWARD);
+ bool WisF = (W.Orientation()==TopAbs_FORWARD);
+ bool isForward = (E.Orientation()==TopAbs_FORWARD);
+ //if(isForward) cout<<"E is FORWARD"<<endl;
+ //else cout<<"E is REVERSED"<<endl;
+ if(!WisF) isForward = !isForward;
+ if(!FisF) isForward = !isForward;
+ //bool isForward = !(E.Orientation()==TopAbs_FORWARD);
+ if(IsReverse) isForward = !isForward;
+ double paramin = 0;
+ double paramax = 0;
+ if (isForward) {
+ paramin = f;
+ paramax = l;
+ gp_Pnt2d p = C2d->Value(f); // first point = Vertex Forward
+ uvslf[0].x = p.X();
+ uvslf[0].y = p.Y();
+ uvslf[0].param = f;
+ uvslf[0].node = idFirst;
+ //MESSAGE("__ f "<<f<<" "<<uvslf[0].x <<" "<<uvslf[0].y);
+ map < double, const SMDS_MeshNode* >::iterator itp = params.begin();
+ for (int i = 1; i <= nbPoints; i++) { // nbPoints internal
+ double param = (*itp).first;
+ gp_Pnt2d p = C2d->Value(param);
+ uvslf[i].x = p.X();
+ uvslf[i].y = p.Y();
+ uvslf[i].param = param;
+ uvslf[i].node = (*itp).second;
+ //MESSAGE("__ "<<i<<" "<<param<<" "<<uvslf[i].x <<" "<<uvslf[i].y);
+ itp++;
+ }
+ p = C2d->Value(l); // last point = Vertex Reversed
+ uvslf[nbPoints + 1].x = p.X();
+ uvslf[nbPoints + 1].y = p.Y();
+ uvslf[nbPoints + 1].param = l;
+ uvslf[nbPoints + 1].node = idLast;
+ //MESSAGE("__ l "<<l<<" "<<uvslf[nbPoints+1].x <<" "<<uvslf[nbPoints+1].y);
+ }
+ else {
+ paramin = l;
+ paramax = f;
+ gp_Pnt2d p = C2d->Value(l); // first point = Vertex Reversed
+ uvslf[0].x = p.X();
+ uvslf[0].y = p.Y();
+ uvslf[0].param = l;
+ uvslf[0].node = idLast;
+ //MESSAGE("__ l "<<l<<" "<<uvslf[0].x <<" "<<uvslf[0].y);
+ map < double, const SMDS_MeshNode* >::reverse_iterator itp = params.rbegin();
+ for (int j = nbPoints; j >= 1; j--) { // nbPoints internal
+ double param = (*itp).first;
+ int i = nbPoints + 1 - j;
+ gp_Pnt2d p = C2d->Value(param);
+ uvslf[i].x = p.X();
+ uvslf[i].y = p.Y();
+ uvslf[i].param = param;
+ uvslf[i].node = (*itp).second;
+ //MESSAGE("__ "<<i<<" "<<param<<" "<<uvslf[i].x <<" "<<uvslf[i].y);
+ itp++;
+ }
+ p = C2d->Value(f); // last point = Vertex Forward
+ uvslf[nbPoints + 1].x = p.X();
+ uvslf[nbPoints + 1].y = p.Y();
+ uvslf[nbPoints + 1].param = f;
+ uvslf[nbPoints + 1].node = idFirst;
+ //MESSAGE("__ f "<<f<<" "<<uvslf[nbPoints+1].x <<" "<<uvslf[nbPoints+1].y);
+ }
+
+ ASSERT(paramin != paramax);
+ for (int i = 0; i < nbPoints + 2; i++) {
+ uvslf[i].normParam = (uvslf[i].param - paramin) / (paramax - paramin);
+ }
+
+ return uvslf;
+}
+
+
//=============================================================================
/*!
* LoadEdgePoints