void copyMesh(GEdge *from, GEdge *to, int direction) { Range<double> u_bounds = from->parBounds(0); double u_min = u_bounds.low(); double u_max = u_bounds.high(); Range<double> to_u_bounds = to->parBounds(0); double to_u_min = to_u_bounds.low(); double to_u_max = to_u_bounds.high(); for(unsigned int i = 0; i < from->mesh_vertices.size(); i++){ int index = (direction < 0) ? (from->mesh_vertices.size() - 1 - i) : i; MVertex *v = from->mesh_vertices[index]; double u; v->getParameter(0, u); double newu = (direction > 0) ? (u-u_min+to_u_min) : (u_max-u+to_u_min); GPoint gp = to->point(newu); MEdgeVertex *vv = new MEdgeVertex(gp.x(), gp.y(), gp.z(), to, newu); to->mesh_vertices.push_back(vv); to->correspondingVertices[vv] = v; } for(unsigned int i = 0; i < to->mesh_vertices.size() + 1; i++){ MVertex *v0 = (i == 0) ? to->getBeginVertex()->mesh_vertices[0] : to->mesh_vertices[i - 1]; MVertex *v1 = (i == to->mesh_vertices.size()) ? to->getEndVertex()->mesh_vertices[0] : to->mesh_vertices[i]; to->lines.push_back(new MLine(v0, v1)); } }
void backgroundMesh2D::updateSizes() { DoubleStorageType::iterator itv = sizeField.begin(); for ( ; itv != sizeField.end(); ++itv) { SPoint2 p; MVertex *v = _2Dto3D[itv->first]; double lc; if (v->onWhat()->dim() == 0) { lc = sizeFactor * BGM_MeshSize(v->onWhat(), 0,0,v->x(),v->y(),v->z()); } else if (v->onWhat()->dim() == 1) { double u; v->getParameter(0, u); lc = sizeFactor * BGM_MeshSize(v->onWhat(), u, 0, v->x(), v->y(), v->z()); } else { GFace *face = dynamic_cast<GFace*>(gf); if(!face) { Msg::Error("Entity is not a face in background mesh"); return; } reparamMeshVertexOnFace(v, face, p); lc = sizeFactor * BGM_MeshSize(face, p.x(), p.y(), v->x(), v->y(), v->z()); } // printf("2D -- %g %g 3D -- %g %g\n",p.x(),p.y(),v->x(),v->y()); itv->second = min(lc,itv->second); itv->second = max(itv->second, sizeFactor * CTX::instance()->mesh.lcMin); itv->second = min(itv->second, sizeFactor * CTX::instance()->mesh.lcMax); } // do not allow large variations in the size field // (Int. J. Numer. Meth. Engng. 43, 1143-1165 (1998) MESH GRADATION // CONTROL, BOROUCHAKI, HECHT, FREY) std::set<MEdge,Less_Edge> edges; for (unsigned int i = 0; i < getNumMeshElements(); i++) { for (int j = 0; j < getElement(i)->getNumEdges(); j++) { edges.insert(getElement(i)->getEdge(j)); } } const double _beta = 1.3; for (int i=0; i<0; i++) { std::set<MEdge,Less_Edge>::iterator it = edges.begin(); for ( ; it != edges.end(); ++it) { MVertex *v0 = it->getVertex(0); MVertex *v1 = it->getVertex(1); MVertex *V0 = _2Dto3D[v0]; MVertex *V1 = _2Dto3D[v1]; DoubleStorageType::iterator s0 = sizeField.find(V0); DoubleStorageType::iterator s1 = sizeField.find(V1); if (s0->second < s1->second)s1->second = min(s1->second,_beta*s0->second); else s0->second = min(s0->second,_beta*s1->second); } } }
void GEdge::relocateMeshVertices() { for(unsigned int i = 0; i < mesh_vertices.size(); i++){ MVertex *v = mesh_vertices[i]; double u0 = 0; if(v->getParameter(0, u0)){ GPoint p = point(u0); v->x() = p.x(); v->y() = p.y(); v->z() = p.z(); } } }
void Mesh::getGEntityPositions(std::vector<SPoint3> &xyz, std::vector<SPoint3> &uvw) { xyz.resize(nVert()); uvw.resize(nFV()); for (int iV = 0; iV < nVert(); iV++) xyz[iV] = SPoint3(_vert[iV]->x(),_vert[iV]->y(),_vert[iV]->z()); for (int iFV = 0; iFV < nFV(); iFV++){ MVertex *v = _freeVert[iFV]; if (v->onWhat()->dim() == 1){ double t; v->getParameter(0,t); uvw[iFV] = SPoint3(t,0,0); } if (v->onWhat()->dim() == 2){ double uu,vv; v->getParameter(0,uu); v->getParameter(1,vv); uvw[iFV] = SPoint3(uu,vv,0); } } }
bool GEdge::computeDistanceFromMeshToGeometry (double &d2, double &dmax) { d2 = 0.0; dmax = 0.0; if (geomType() == Line) return true; if (!lines.size())return false; IntPt *pts; int npts; lines[0]->getIntegrationPoints(2*lines[0]->getPolynomialOrder(), &npts, &pts); for (unsigned int i = 0; i < lines.size(); i++){ MLine *l = lines[i]; double t[256]; for (int j=0; j< l->getNumVertices();j++){ MVertex *v = l->getVertex(j); if (v->onWhat() == getBeginVertex()){ t[j] = getLowerBound(); } else if (v->onWhat() == getEndVertex()){ t[j] = getUpperBound(); } else { v->getParameter(0,t[j]); } } for (int j=0;j<npts;j++){ SPoint3 p; l->pnt(pts[j].pt[0],0,0,p); double tinit = l->interpolate(t,pts[j].pt[0],0,0); GPoint pc = closestPoint(p, tinit); if (!pc.succeeded())continue; double dsq = (pc.x()-p.x())*(pc.x()-p.x()) + (pc.y()-p.y())*(pc.y()-p.y()) + (pc.z()-p.z())*(pc.z()-p.z()); d2 += pts[i].weight * fabs(l->getJacobianDeterminant(pts[j].pt[0],0,0)) * dsq; dmax = std::max(dmax,sqrt(dsq)); } } d2 = sqrt(d2); return true; }
bool OptHOM::addBndObjGrad(double factor, double &Obj, alglib::real_1d_array &gradObj) { // set the mesh to its present position std::vector<SPoint3> xyz,uvw; mesh.getGEntityPositions(xyz,uvw); mesh.updateGEntityPositions(); //could be better (e.g. store the model in the Mesh:: datastrucure) GModel *gm = GModel::current(); // for all model edges, compute the error between the geometry and the mesh maxDistCAD = 0.0; double distCAD = 0.0; for (GModel::eiter it = gm->firstEdge(); it != gm->lastEdge(); ++it){ // do not do straight lines if ((*it)->geomType() == GEntity::Line)continue; // look at all mesh lines std::vector<bool> doWeCompute((*it)->lines.size()); for (unsigned int i=0;i<(*it)->lines.size(); i++){ doWeCompute[i] = false; for (unsigned int j=0;j<(*it)->lines[i]->getNumVertices(); j++){ int index = mesh.getFreeVertexStartIndex((*it)->lines[i]->getVertex(j)); if (index >=0){ doWeCompute[i] = true; continue; } } } std::vector<double> dist((*it)->lines.size()); for (unsigned int i=0;i<(*it)->lines.size(); i++){ if (doWeCompute[i]){ // compute the distance from the geometry to the mesh dist[i] = MLineGEdgeDistance ( (*it)->lines[i] , *it ); maxDistCAD = std::max(maxDistCAD,dist[i]); distCAD += dist [i] * factor; } } // be clever to compute the derivative : iterate on all // Distance = \sum_{lines} Distance (line, GEdge) // For a high order vertex, we compute the derivative only by // recomputing the distance to one only line const double eps = 1.e-6; for (unsigned int i=0;i<(*it)->lines.size(); i++){ if (doWeCompute[i]){ for (int j=2 ; j<(*it)->lines[i]->getNumVertices() ; j++){ MVertex *v = (*it)->lines[i]->getVertex(j); int index = mesh.getFreeVertexStartIndex(v); // printf("%d %d (%d %d)\n",v->getNum(),index,v->onWhat()->tag(),v->onWhat()->dim()); if (index >= 0){ double t; v->getParameter(0,t); SPoint3 pp (v->x(),v->y(),v->z()); GPoint gp = (*it)->point(t+eps); v->setParameter(0,t+eps); v->setXYZ(gp.x(),gp.y(),gp.z()); double dist2 = MLineGEdgeDistance ( (*it)->lines[i] , *it ); double deriv = (dist2 - dist[i])/eps; v->setXYZ(pp.x(),pp.y(),pp.z()); v->setParameter(0,t); // printf("%g %g %g\n",dist[i],dist2, MLineGEdgeDistance ( (*it)->lines[i] , *it )); // get the index of the vertex gradObj[index] += deriv * factor; } } } // printf("done\n"); // For a low order vertex classified on the GEdge, we recompute // two distances for the two MLines connected to the vertex for (unsigned int i=0;i<(*it)->lines.size()-1; i++){ MVertex *v = (*it)->lines[i]->getVertex(1); int index = mesh.getFreeVertexStartIndex(v); if (index >= 0){ double t; v->getParameter(0,t); SPoint3 pp (v->x(),v->y(),v->z()); GPoint gp = (*it)->point(t+eps); v->setParameter(0,t+eps); v->setXYZ(gp.x(),gp.y(),gp.z()); MLine *l1 = (*it)->lines[i]; MLine *l2 = (*it)->lines[i+1]; // printf("%d %d -- %d %d\n",l1->getVertex(0)->getNum(),l1->getVertex(1)->getNum(),l2->getVertex(0)->getNum(),l2->getVertex(1)->getNum()); double deriv = (MLineGEdgeDistance ( l1 , *it ) - dist[i]) /eps + (MLineGEdgeDistance ( l2 , *it ) - dist[i+1])/eps; v->setXYZ(pp.x(),pp.y(),pp.z()); v->setParameter(0,t); gradObj[index] += deriv * factor; } } } } // printf("computing distance : 1D part %12.5E\n",distCAD); // now the 3D part ! std::vector<std::vector<SVector3> > gsfT; computeGradSFAtNodes ( (*gm->firstFace())->triangles[0],gsfT); std::map<MVertex*,SVector3> normalsToCAD; for(GModel::fiter it = gm->firstFace(); it != gm->lastFace(); ++it){ // do not do plane surfaces if ((*it)->geomType() == GEntity::Plane)continue; std::map<MTriangle*,double> dist; std::vector<bool> doWeCompute((*it)->triangles.size()); for (unsigned int i=0;i<(*it)->triangles.size(); i++){ doWeCompute[i] = false; for (unsigned int j=0;j<(*it)->triangles[i]->getNumVertices(); j++){ int index = mesh.getFreeVertexStartIndex((*it)->triangles[i]->getVertex(j)); if (index >=0){ doWeCompute[i] = true; } } if (doWeCompute[i]){ for (unsigned int j=0;j<(*it)->triangles[i]->getNumVertices(); j++){ MVertex *v = (*it)->triangles[i]->getVertex(j); if (normalsToCAD.find(v) == normalsToCAD.end()){ SPoint2 p_cad; reparamMeshVertexOnFace(v, *it, p_cad); SVector3 tg_cad = (*it)->normal(p_cad); tg_cad.normalize(); normalsToCAD[v] = tg_cad; } } } } for (unsigned int i=0;i<(*it)->triangles.size(); i++){ // compute the distance from the geometry to the mesh if(doWeCompute[i]){ const double d = MFaceGFaceDistanceOld((*it)->triangles[i], *it, &gsfT, &normalsToCAD); dist[(*it)->triangles[i]] = d; maxDistCAD = std::max(maxDistCAD,d); distCAD += d * factor; } } // be clever again to compute the derivatives const double eps = 1.e-6; for (unsigned int i=0;i<(*it)->triangles.size(); i++){ if(doWeCompute[i]){ for (unsigned int j=0;j<(*it)->triangles[i]->getNumVertices(); j++){ // for (; itm !=v2t.end(); ++itm){ MVertex *v = (*it)->triangles[i]->getVertex(j); if(v->onWhat()->dim() == 1){ int index = mesh.getFreeVertexStartIndex(v); if (index >= 0){ MTriangle *t = (*it)->triangles[i]; GEdge *ge = v->onWhat()->cast2Edge(); double t_; v->getParameter(0,t_); SPoint3 pp (v->x(),v->y(),v->z()); GPoint gp = ge->point(t_+eps); v->setParameter(0,t_+eps); v->setXYZ(gp.x(),gp.y(),gp.z()); const double distT = dist[t]; double deriv = (MFaceGFaceDistanceOld(t, *it, &gsfT, &normalsToCAD) - distT) /eps; v->setXYZ(pp.x(),pp.y(),pp.z()); v->setParameter(0,t_); gradObj[index] += deriv * factor; } } if(v->onWhat() == *it){ int index = mesh.getFreeVertexStartIndex(v); if (index >= 0){ MTriangle *t = (*it)->triangles[i]; double uu,vv; v->getParameter(0,uu); v->getParameter(1,vv); SPoint3 pp (v->x(),v->y(),v->z()); const double distT = dist[t]; GPoint gp = (*it)->point(uu+eps,vv); v->setParameter(0,uu+eps); v->setXYZ(gp.x(),gp.y(),gp.z()); double deriv = (MFaceGFaceDistanceOld(t, *it, &gsfT, &normalsToCAD) - distT) /eps; v->setXYZ(pp.x(),pp.y(),pp.z()); v->setParameter(0,uu); gradObj[index] += deriv * factor; gp = (*it)->point(uu,vv+eps); v->setParameter(1,vv+eps); v->setXYZ(gp.x(),gp.y(),gp.z()); deriv = (MFaceGFaceDistanceOld(t, *it, &gsfT, &normalsToCAD) - distT) /eps; v->setXYZ(pp.x(),pp.y(),pp.z()); v->setParameter(1,vv); gradObj[index+1] += deriv * factor; } } } } } } mesh.updateGEntityPositions(xyz,uvw); Obj +=distCAD; // printf("computing distance : 2D part %12.5E\n",distCAD); // printf("%22.15E\n",distCAD); return true; }