Skip to content
Snippets Groups Projects
Commit a30e5248 authored by Emilie Marchandise's avatar Emilie Marchandise
Browse files

Rbf works fine with eps=0.5*distMin, delta=0.3*distMin, epsUV and deltaUV for...

Rbf works fine with eps=0.5*distMin, delta=0.3*distMin, epsUV and deltaUV for the Rbf interpolatin in the newton
parent 5b67c181
No related branches found
No related tags found
No related merge requests found
...@@ -665,7 +665,7 @@ bool GFaceCompound::parametrize() const ...@@ -665,7 +665,7 @@ bool GFaceCompound::parametrize() const
int variableEps = 0; int variableEps = 0;
int radFunInd = 1; //MQ RBF int radFunInd = 1; //MQ RBF
double delta = 0.33*getDistMin(); double delta = 0.33*getDistMin();
double epsilon = 0.15*(allNodes.size()/150.0)/delta; //max(2.5, 0.5*(nbNodes/150.0)/dist_min); double epsilon = 0.5/getDistMin(); //0.15*(allNodes.size()/150.0)/delta; //max(2.5, 0.5*(nbNodes/150.0)/dist_min);
double radius= 3.*getSizeH()/sqrt(allNodes.size()); double radius= 3.*getSizeH()/sqrt(allNodes.size());
Msg::Info("*****************************************"); Msg::Info("*****************************************");
...@@ -685,7 +685,6 @@ bool GFaceCompound::parametrize() const ...@@ -685,7 +685,6 @@ bool GFaceCompound::parametrize() const
_rbf->solveHarmonicMap(Oper, _ordered, _coords, coordinates); _rbf->solveHarmonicMap(Oper, _ordered, _coords, coordinates);
printStuff(); printStuff();
exit(1);
} }
...@@ -781,16 +780,16 @@ double GFaceCompound::getSizeH() const ...@@ -781,16 +780,16 @@ double GFaceCompound::getSizeH() const
double GFaceCompound::getDistMin() const double GFaceCompound::getDistMin() const
{ {
double dist_min = 1.e6; double dist_min = 1.e6;
double tol = 1.e-8;
for(std::set<MVertex *>::iterator itv = allNodes.begin(); itv !=allNodes.end() ; itv++){ for(std::set<MVertex *>::iterator itv = allNodes.begin(); itv !=allNodes.end() ; itv++){
for(std::set<MVertex *>::iterator itv2 = allNodes.begin(); itv2 !=allNodes.end() ; itv2++){ for(std::set<MVertex *>::iterator itv2 = allNodes.begin(); itv2 !=allNodes.end() ; itv2++){
MVertex *v1 = *itv; MVertex *v1 = *itv;
MVertex *v2 = *itv2; MVertex *v2 = *itv2;
double dist = sqrt((v1->x()-v2->x())*(v1->x()-v2->x())+(v1->y()-v2->y())*(v1->y()-v2->y()) double dist = sqrt((v1->x()-v2->x())*(v1->x()-v2->x())+(v1->y()-v2->y())*(v1->y()-v2->y())
+(v1->z()-v2->z())*(v1->z()-v2->z())); +(v1->z()-v2->z())*(v1->z()-v2->z()));
if (dist<dist_min && dist != 0.0) dist_min = dist; if (dist<dist_min && dist > tol) dist_min = dist;
} }
} }
return dist_min; return dist_min;
} }
double GFaceCompound::getSizeBB(const std::list<GEdge* > &elist) const double GFaceCompound::getSizeBB(const std::list<GEdge* > &elist) const
...@@ -1446,11 +1445,15 @@ GPoint GFaceCompound::point(double par1, double par2) const ...@@ -1446,11 +1445,15 @@ GPoint GFaceCompound::point(double par1, double par2) const
if(!oct) parametrize(); if(!oct) parametrize();
if (_mapping == RBF){ if (_mapping == RBF){
if (fabs(par1) > 1 || fabs(par2) > 1){
GPoint gp (3,3,0,this);
gp.setNoSuccess();
return gp;
}
double x, y, z; double x, y, z;
SVector3 dXdu, dXdv; SVector3 dXdu, dXdv;
_rbf->UVStoXYZ(par1, par2,x,y,z, dXdu, dXdv); bool conv = _rbf->UVStoXYZ(par1, par2,x,y,z, dXdu, dXdv);
//printf("XYZ =%g %g %g \n", x,y,z); if (!conv) printf("UV=%g %g \n", par1, par2);
//exit(1);
return GPoint(x,y,z); return GPoint(x,y,z);
} }
...@@ -1543,7 +1546,7 @@ Pair<SVector3,SVector3> GFaceCompound::firstDer(const SPoint2 &param) const ...@@ -1543,7 +1546,7 @@ Pair<SVector3,SVector3> GFaceCompound::firstDer(const SPoint2 &param) const
if (_mapping == RBF){ if (_mapping == RBF){
double x, y, z; double x, y, z;
SVector3 dXdu, dXdv ; SVector3 dXdu, dXdv ;
_rbf->UVStoXYZ(param.x(), param.y(), x,y,z, dXdu, dXdv); bool conv = _rbf->UVStoXYZ(param.x(), param.y(), x,y,z, dXdu, dXdv);
return Pair<SVector3, SVector3>(dXdu,dXdv); return Pair<SVector3, SVector3>(dXdu,dXdv);
} }
...@@ -1892,6 +1895,8 @@ double GFaceCompound::checkAspectRatio() const ...@@ -1892,6 +1895,8 @@ double GFaceCompound::checkAspectRatio() const
void GFaceCompound::coherencePatches() const void GFaceCompound::coherencePatches() const
{ {
if (_mapping == RBF) return;
Msg::Info("Re-orient all %d compound patches normals coherently", _compound.size()); Msg::Info("Re-orient all %d compound patches normals coherently", _compound.size());
std::map<MEdge, std::set<MElement*>, Less_Edge > edge2elems; std::map<MEdge, std::set<MElement*>, Less_Edge > edge2elems;
...@@ -2033,7 +2038,7 @@ int GFaceCompound::genusGeom() const ...@@ -2033,7 +2038,7 @@ int GFaceCompound::genusGeom() const
void GFaceCompound::printStuff(int iNewton) const void GFaceCompound::printStuff(int iNewton) const
{ {
if( !CTX::instance()->mesh.saveAll) return; //if( !CTX::instance()->mesh.saveAll) return;
std::list<GFace*>::const_iterator it = _compound.begin(); std::list<GFace*>::const_iterator it = _compound.begin();
......
...@@ -45,7 +45,7 @@ static int SphereInEle(void *a, double*c){ ...@@ -45,7 +45,7 @@ static int SphereInEle(void *a, double*c){
GRbf::GRbf (double eps, double del, double rad, int variableEps, int rbfFun, GRbf::GRbf (double eps, double del, double rad, int variableEps, int rbfFun,
std::map<MVertex*, SVector3> _normals, std::set<MVertex *> allNodes) std::map<MVertex*, SVector3> _normals, std::set<MVertex *> allNodes)
: ep_scalar(eps), delta(del), radius (rad), variableShapeParam(variableEps), : epsilonXYZ(eps), epsilonUV(eps), delta(del), radius (rad), variableShapeParam(variableEps),
radialFunctionIndex (rbfFun), XYZkdtree(0) radialFunctionIndex (rbfFun), XYZkdtree(0)
{ {
...@@ -231,19 +231,22 @@ double GRbf::evalRadialFnDer (int p, double dx, double dy, double dz, double ep) ...@@ -231,19 +231,22 @@ double GRbf::evalRadialFnDer (int p, double dx, double dy, double dz, double ep)
fullMatrix<double> GRbf::generateRbfMat(int p, fullMatrix<double> GRbf::generateRbfMat(int p,
const fullMatrix<double> &nodes1, const fullMatrix<double> &nodes1,
const fullMatrix<double> &nodes2) { const fullMatrix<double> &nodes2,
int inUV) {
int m = nodes2.size1(); int m = nodes2.size1();
int n = nodes1.size1(); int n = nodes1.size1();
fullMatrix<double> rbfMat(m,n); fullMatrix<double> rbfMat(m,n);
fullVector<double > epsilon; //fullVector<double > epsilon;
computeEpsilon(nodes1, epsilon); //computeEpsilon(nodes1, epsilon, inUN);
double eps = epsilonXYZ;
if (inUV) eps = epsilonUV;
for (int i = 0; i < m; i++) { for (int i = 0; i < m; i++) {
for (int j = 0; j < n; j++) { for (int j = 0; j < n; j++) {
double dx = -nodes2(i,0)+nodes1(j,0); double dx = -nodes2(i,0)+nodes1(j,0);
double dy = -nodes2(i,1)+nodes1(j,1); double dy = -nodes2(i,1)+nodes1(j,1);
double dz = -nodes2(i,2)+nodes1(j,2); double dz = -nodes2(i,2)+nodes1(j,2);
rbfMat(i,j) = evalRadialFnDer(p,dx,dy,dz,epsilon(j)); rbfMat(i,j) = evalRadialFnDer(p,dx,dy,dz,eps);
} }
} }
...@@ -254,28 +257,28 @@ fullMatrix<double> GRbf::generateRbfMat(int p, ...@@ -254,28 +257,28 @@ fullMatrix<double> GRbf::generateRbfMat(int p,
void GRbf::RbfOp(int p, void GRbf::RbfOp(int p,
const fullMatrix<double> &cntrs, const fullMatrix<double> &cntrs,
const fullMatrix<double> &nodes, const fullMatrix<double> &nodes,
fullMatrix<double> &D) { fullMatrix<double> &D, int inUV) {
fullMatrix<double> rbfInvA, rbfMatB; fullMatrix<double> rbfInvA, rbfMatB;
D.resize(nodes.size1(), cntrs.size1()); D.resize(nodes.size1(), cntrs.size1());
if (isLocal){ if (isLocal){
rbfInvA = generateRbfMat(0,cntrs,cntrs); rbfInvA = generateRbfMat(0,cntrs,cntrs, inUV);
rbfInvA.invertInPlace(); rbfInvA.invertInPlace();
} }
else{ else{
if (cntrs.size1() == nbNodes ) if (cntrs.size1() == nbNodes && !inUV)
rbfInvA = matAInv; rbfInvA = matAInv;
else if (cntrs.size1() == 3*nbNodes ) else if (cntrs.size1() == 3*nbNodes && !inUV)
rbfInvA = matAInv_nn; rbfInvA = matAInv_nn;
else{ else{
rbfInvA = generateRbfMat(0,cntrs,cntrs); rbfInvA = generateRbfMat(0,cntrs,cntrs, inUV);
rbfInvA.invertInPlace(); rbfInvA.invertInPlace();
} }
} }
rbfMatB = generateRbfMat(p,cntrs,nodes); rbfMatB = generateRbfMat(p,cntrs,nodes, inUV);
D.gemm(rbfMatB, rbfInvA, 1.0, 0.0); D.gemm(rbfMatB, rbfInvA, 1.0, 0.0);
} }
...@@ -284,19 +287,22 @@ void GRbf::evalRbfDer(int p, ...@@ -284,19 +287,22 @@ void GRbf::evalRbfDer(int p,
const fullMatrix<double> &cntrs, const fullMatrix<double> &cntrs,
const fullMatrix<double> &nodes, const fullMatrix<double> &nodes,
const fullMatrix<double> &fValues, const fullMatrix<double> &fValues,
fullMatrix<double> &fApprox) { fullMatrix<double> &fApprox,
int inUV) {
fApprox.resize(nodes.size1(),fValues.size2()); fApprox.resize(nodes.size1(),fValues.size2());
fullMatrix<double> D; fullMatrix<double> D;
RbfOp(p,cntrs,nodes,D); RbfOp(p,cntrs,nodes,D, inUV);
fApprox.gemm(D,fValues, 1.0, 0.0); fApprox.gemm(D,fValues, 1.0, 0.0);
} }
void GRbf::computeEpsilon(const fullMatrix<double> &cntrs, fullVector<double> &epsilon){ void GRbf::computeEpsilon(const fullMatrix<double> &cntrs, fullVector<double> &epsilon,
int inUV){
epsilon.resize(nbNodes*3); epsilon.resize(nbNodes*3);
epsilon.setAll(ep_scalar); if (inUV) epsilon.setAll(epsilonUV);
epsilon.setAll(epsilonXYZ);
if (variableShapeParam) { if (variableShapeParam) {
...@@ -603,20 +609,12 @@ void GRbf::solveHarmonicMap(fullMatrix<double> Oper, ...@@ -603,20 +609,12 @@ void GRbf::solveHarmonicMap(fullMatrix<double> Oper,
int nb = Oper.size2(); int nb = Oper.size2();
UV.resize(nb,2); UV.resize(nb,2);
fullMatrix<double> dirichletBC(ordered.size(),2);
for(unsigned int i = 0; i < ordered.size(); i++){
MVertex *v = ordered[i];
std::map<MVertex*, int>::iterator itm = _mapV.find(v);
const double theta = 2 * M_PI * coords[i];
dirichletBC(i,0) = itm->second;
dirichletBC(i,1) = theta;
}
fullMatrix<double> F(nb,2); fullMatrix<double> F(nb,2);
F.scale(0.0); F.scale(0.0);
for (int i=0; i < ordered.size(); i++){ for (int i=0; i < ordered.size(); i++){
int iFix = (int)dirichletBC(i,0); std::map<MVertex*, int>::iterator itm = _mapV.find(ordered[i]);
double theta = dirichletBC(i,1); double theta = 2 * M_PI * coords[i];
int iFix = itm->second;
for (int j=0;j<nb ; ++j) Oper(iFix,j) = 0.0; for (int j=0;j<nb ; ++j) Oper(iFix,j) = 0.0;
Oper(iFix,iFix) = 1.0; Oper(iFix,iFix) = 1.0;
F(iFix,0) = cos(theta); F(iFix,0) = cos(theta);
...@@ -626,17 +624,26 @@ void GRbf::solveHarmonicMap(fullMatrix<double> Oper, ...@@ -626,17 +624,26 @@ void GRbf::solveHarmonicMap(fullMatrix<double> Oper,
Oper.invertInPlace(); Oper.invertInPlace();
Oper.mult(F,UV); Oper.mult(F,UV);
//ANN UVtree + dist_min
double dist_min = 1.e6;
#if defined (HAVE_ANN) #if defined (HAVE_ANN)
UVnodes = annAllocPts(nbNodes, 3); UVnodes = annAllocPts(nbNodes, 3);
for(int i = 0; i < nbNodes; i++){ for(int i = 0; i < nbNodes; i++){
UVnodes[i][0] = UV(i,0); UVnodes[i][0] = UV(i,0);
UVnodes[i][1] = UV(i,1); UVnodes[i][1] = UV(i,1);
UVnodes[i][2] = 0.0; UVnodes[i][2] = 0.0;
for(int j = i+1; j < nbNodes; j++){
double dist = sqrt((UV(i,0)-UV(j,0))*(UV(i,0)-UV(j,0))+
(UV(i,1)-UV(j,1))*(UV(i,1)-UV(j,1)));
if (dist<dist_min) dist_min = dist;
}
} }
UVkdtree = new ANNkd_tree(UVnodes, nbNodes, 3); UVkdtree = new ANNkd_tree(UVnodes, nbNodes, 3);
index = new ANNidx[num_neighbours]; index = new ANNidx[num_neighbours];
dist = new ANNdist[num_neighbours]; dist = new ANNdist[num_neighbours];
#endif #endif
epsilonUV = 0.5/dist_min;
deltaUV = 0.6*dist_min;
//fill rbf_param //fill rbf_param
std::map<MVertex*, int>::iterator itm = _mapV.begin(); std::map<MVertex*, int>::iterator itm = _mapV.begin();
...@@ -728,7 +735,7 @@ void GRbf::UVStoXYZ_global(const double u_eval, const double v_eval, ...@@ -728,7 +735,7 @@ void GRbf::UVStoXYZ_global(const double u_eval, const double v_eval,
} }
void GRbf::UVStoXYZ(const double u_eval, const double v_eval, bool GRbf::UVStoXYZ(const double u_eval, const double v_eval,
double &XX, double &YY, double &ZZ, double &XX, double &YY, double &ZZ,
SVector3 &dXdu, SVector3& dXdv){ SVector3 &dXdu, SVector3& dXdv){
...@@ -736,6 +743,7 @@ void GRbf::UVStoXYZ(const double u_eval, const double v_eval, ...@@ -736,6 +743,7 @@ void GRbf::UVStoXYZ(const double u_eval, const double v_eval,
//Thus in total, we're working with '3*num_neighbours' nodes //Thus in total, we're working with '3*num_neighbours' nodes
//Say that the vector 'index' gives us the indices of the closest points //Say that the vector 'index' gives us the indices of the closest points
//TO FILL IN //TO FILL IN
bool converged = true;
#if defined (HAVE_ANN) #if defined (HAVE_ANN)
double uvw[3] = { u_eval, v_eval, 0.0 }; double uvw[3] = { u_eval, v_eval, 0.0 };
...@@ -753,15 +761,16 @@ void GRbf::UVStoXYZ(const double u_eval, const double v_eval, ...@@ -753,15 +761,16 @@ void GRbf::UVStoXYZ(const double u_eval, const double v_eval,
//THE LOCAL UVS //THE LOCAL UVS
u_vec(i,0) = UV(index[i],0); u_vec(i,0) = UV(index[i],0);
u_vec(i,1) = UV(index[i],1); u_vec(i,1) = UV(index[i],1);
u_vec(i,2) = surfInterp(index[i],0); u_vec(i,2) = 0.0;
u_vec(i+num_neighbours,0) = UV(index[i]+nbNodes,0); u_vec(i+num_neighbours,0) = UV(index[i]+nbNodes,0);
u_vec(i+num_neighbours,1) = UV(index[i]+nbNodes,1); u_vec(i+num_neighbours,1) = UV(index[i]+nbNodes,1);
u_vec(i+num_neighbours,2) = surfInterp(index[i]+nbNodes,0); u_vec(i+num_neighbours,2) = surfInterp(index[i]+nbNodes,0)*deltaUV;
u_vec(i+2*num_neighbours,0) = UV(index[i]+2*nbNodes,0); u_vec(i+2*num_neighbours,0) = UV(index[i]+2*nbNodes,0);
u_vec(i+2*num_neighbours,1) = UV(index[i]+2*nbNodes,1); u_vec(i+2*num_neighbours,1) = UV(index[i]+2*nbNodes,1);
u_vec(i+2*num_neighbours,2) = surfInterp(index[i]+2*nbNodes,0); u_vec(i+2*num_neighbours,2) = surfInterp(index[i]+2*nbNodes,0)*deltaUV;
//THE LOCAL XYZ //THE LOCAL XYZ
xyz_local(i,0) = extendedX(index[i],0); xyz_local(i,0) = extendedX(index[i],0);
...@@ -784,7 +793,13 @@ void GRbf::UVStoXYZ(const double u_eval, const double v_eval, ...@@ -784,7 +793,13 @@ void GRbf::UVStoXYZ(const double u_eval, const double v_eval,
u_vec_eval(0,2) = 0.0; u_vec_eval(0,2) = 0.0;
//we will use a local interpolation to find the corresponding XYZ point to (u_eval,v_eval). //we will use a local interpolation to find the corresponding XYZ point to (u_eval,v_eval).
evalRbfDer(0, u_vec, u_vec_eval,xyz_local, nodes_eval); evalRbfDer(0, u_vec, u_vec_eval,xyz_local, nodes_eval, 1);
//printf("nodes eval =%g %g %g \n", nodes_eval(0,0), nodes_eval(0,1), nodes_eval(0,2));
//exit(1);
// nodes_eval(0,0) = extendedX(index[0],0);
// nodes_eval(0,1) = extendedX(index[0],1);
// nodes_eval(0,2) = extendedX(index[0],2);
u_temp(0,0) = u_eval; u_temp(0,0) = u_eval;
u_temp(0,1) = v_eval; u_temp(0,1) = v_eval;
...@@ -793,7 +808,7 @@ void GRbf::UVStoXYZ(const double u_eval, const double v_eval, ...@@ -793,7 +808,7 @@ void GRbf::UVStoXYZ(const double u_eval, const double v_eval,
int incr = 0; int incr = 0;
double norm_s = 0.0; double norm_s = 0.0;
fullMatrix<double> Jac(3,3); fullMatrix<double> Jac(3,3);
while(norm_s <5 && incr < 5){ while(norm_s <5 && incr < 10){
// Find the entries of the m Jacobians // Find the entries of the m Jacobians
evalRbfDer(1,xyz_local,nodes_eval,u_vec,ux); evalRbfDer(1,xyz_local,nodes_eval,u_vec,ux);
...@@ -822,14 +837,16 @@ void GRbf::UVStoXYZ(const double u_eval, const double v_eval, ...@@ -822,14 +837,16 @@ void GRbf::UVStoXYZ(const double u_eval, const double v_eval,
incr++; incr++;
} }
if (norm_s < 5 ) if (norm_s < 5 ){
printf("Newton not converged for point (uv)=(%g,%g -> norm_s =%g )\n", u_eval, v_eval, norm_s); printf("Newton not converged for point (uv)=(%g,%g -> norm_s =%g ) XYZ =%g %g %g \n", u_eval, v_eval, norm_s, nodes_eval(0,0), nodes_eval(0,1), nodes_eval(0,2));
converged = false;
}
XX = nodes_eval(0,0); XX = nodes_eval(0,0);
YY = nodes_eval(0,1); YY = nodes_eval(0,1);
ZZ = nodes_eval(0,2); ZZ = nodes_eval(0,2);
dXdu = SVector3(Jac(0,0), Jac(1,0), Jac(2,0)); dXdu = SVector3(Jac(0,0), Jac(1,0), Jac(2,0));
dXdv = SVector3(Jac(0,1), Jac(1,1), Jac(2,1)); dXdv = SVector3(Jac(0,1), Jac(1,1), Jac(2,1));
return converged;
} }
...@@ -32,8 +32,10 @@ class GRbf { ...@@ -32,8 +32,10 @@ class GRbf {
int nn; int nn;
int num_neighbours; int num_neighbours;
double ep_scalar; // Shape parameter double epsilonXYZ; // Shape parameter
double epsilonUV; // Shape parameter
double delta; //offset level set double delta; //offset level set
double deltaUV; //offset level set
double radius; double radius;
int variableShapeParam; // 1 if one chooses epsilon to vary spatially, 0 if one chooses it to be constant int variableShapeParam; // 1 if one chooses epsilon to vary spatially, 0 if one chooses it to be constant
int radialFunctionIndex; // Index corresponding to the radial function used (0 - GA,1 - MQ, ... ) int radialFunctionIndex; // Index corresponding to the radial function used (0 - GA,1 - MQ, ... )
...@@ -79,22 +81,23 @@ class GRbf { ...@@ -79,22 +81,23 @@ class GRbf {
//(p)th derivative of the radial function w.r.t. the (q)th variable //(p)th derivative of the radial function w.r.t. the (q)th variable
fullMatrix<double> generateRbfMat(int p, fullMatrix<double> generateRbfMat(int p,
const fullMatrix<double> &nodes1, const fullMatrix<double> &nodes1,
const fullMatrix<double> &nodes2); const fullMatrix<double> &nodes2,
int inUV=0);
// Computes the interpolation(p==0) or the derivative (p!=0) operator(mxn) (n:number of centers, m: number of evaluation nodes) // Computes the interpolation(p==0) or the derivative (p!=0) operator(mxn) (n:number of centers, m: number of evaluation nodes)
void RbfOp(int p, // (p)th derivatives void RbfOp(int p, // (p)th derivatives
const fullMatrix<double> &cntrs, const fullMatrix<double> &cntrs,
const fullMatrix<double> &nodes, const fullMatrix<double> &nodes,
fullMatrix<double> &D); fullMatrix<double> &D, int inUV=0);
// Computes the interpolant(p==0) or the derivative (p!=0) of the function values entered and evaluates it at the new nodes // Computes the interpolant(p==0) or the derivative (p!=0) of the function values entered and evaluates it at the new nodes
void evalRbfDer(int p, // (p)th derivatives void evalRbfDer(int p, // (p)th derivatives
const fullMatrix<double> &cntrs, const fullMatrix<double> &cntrs,
const fullMatrix<double> &nodes, const fullMatrix<double> &nodes,
const fullMatrix<double> &fValues, const fullMatrix<double> &fValues,
fullMatrix<double> &fApprox); fullMatrix<double> &fApprox, int inUV=0);
void computeEpsilon(const fullMatrix<double> &cntrs, fullVector<double> &epsilon); void computeEpsilon(const fullMatrix<double> &cntrs, fullVector<double> &epsilon, int inUV=0);
// Finds surface differentiation matrix using the LOCAL projection method // Finds surface differentiation matrix using the LOCAL projection method
void RbfLapSurface_local_projection(const fullMatrix<double> &cntrs, void RbfLapSurface_local_projection(const fullMatrix<double> &cntrs,
...@@ -129,12 +132,11 @@ class GRbf { ...@@ -129,12 +132,11 @@ class GRbf {
const fullMatrix<double> &node, const fullMatrix<double> &node,
fullMatrix<double> &curvature); fullMatrix<double> &curvature);
virtual void UVStoXYZ_global(const double u_eval, void UVStoXYZ_global(const double u_eval, const double v_eval,
const double v_eval,
double &XX, double &YY, double &ZZ, double &XX, double &YY, double &ZZ,
SVector3 &dXdu, SVector3& dxdv); SVector3 &dXdu, SVector3& dxdv);
virtual void UVStoXYZ(const double u_eval,
const double v_eval, bool UVStoXYZ(const double u_eval, const double v_eval,
double &XX, double &YY, double &ZZ, double &XX, double &YY, double &ZZ,
SVector3 &dXdu, SVector3& dxdv); SVector3 &dXdu, SVector3& dxdv);
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment