// Recursively computes the B-spline basis function at interval u. // Handles division by zero as zero and stops recursing in a branch when multiplication by zero would occur. float BsplineBasis(int j, int d, float u, float * t) { if (d == 1) if (t[j] <= u && u < t[j+1]) return 1.0f; else return 0.0f; else return ((t[j+d-1] == t[j]) ? 0.0f : (u-t[j]) / (t[j+d-1]-t[j]) * BsplineBasis(j, d-1, u, t)) + ((t[j+d] == t[j+1]) ? 0.0f : (t[j+d]-u)/(t[j+d]-t[j+1]) * BsplineBasis(j+1, d-1, u, t)); }
// Computes a point on a N-th degree (d - 1) B-spline curve at interval u. // Returns the point through resx and resy. void computeNOrderBsplinelf(int d, struct list * l, float u, float * t, float * resx, float * resy) { int j = 0; *resx = 0.0f; *resy = 0.0f; j = 0; while (l) { *resx += l->x * BsplineBasis(j, d, u, t); *resy += l->y * BsplineBasis(j, d, u, t); l = l->next; j++; } }
//=========================================================================== void SplineCurve::makeKnotStartRegular() //=========================================================================== { // Testing whether knotstart is already d+1-regular. if (basis_.begin()[0] < basis_.begin()[order() - 1]) { double tpar = basis_.startparam(); int ti = order() - 1; // Index of last occurence of tpar (in other_curve). int mt = 1; // Multiplicity of tpar. while ((basis_.begin()[ti - mt] == tpar) && (mt < order())) ++mt; std::vector<double> new_knots; for (int i = 0; i < order() - mt; ++i) new_knots.push_back(tpar); insertKnot(new_knots); coefs_.erase(coefs_begin(), coefs_begin() + (order() - mt) * dim_); basis_ = BsplineBasis(order(), basis_.begin() + order() - mt, basis_.end()); } }
//=========================================================================== void SplineCurve::appendCurve(ParamCurve* other_curve, int continuity, double& dist, bool repar) //=========================================================================== { SplineCurve* other_cv = dynamic_cast<SplineCurve*>(other_curve); ALWAYS_ERROR_IF(other_cv == 0, "Given an empty curve or not a SplineCurve."); ALWAYS_ERROR_IF(dim_ != other_cv->dimension(), "The curves must lie in the same space."); ALWAYS_ERROR_IF(continuity < -1 || continuity + 1 > order(), "Specified continuity not attainable."); #ifdef DEBUG // DEBUG OUTPUT std::ofstream of0("basis0.g2"); writeStandardHeader(of0); write(of0); other_cv->writeStandardHeader(of0); other_cv->write(of0); #endif // Making sure the curves have the same order. Raise if necessary. int diff_order = order() - other_cv->order(); if (diff_order > 0) other_cv->raiseOrder(diff_order); else if (diff_order < 0) raiseOrder(abs(diff_order)); // Make sure that we have k-regularity at meeting ends. makeKnotEndRegular(); other_cv->makeKnotStartRegular(); // Ensure that either none of the curves or both are rational if (rational_ && !other_cv->rational()) other_cv->representAsRational(); if (!rational_ && other_cv->rational()) representAsRational(); #ifdef DEBUG // DEBUG OUTPUT std::ofstream of1("basis1.g2"); writeStandardHeader(of1); write(of1); other_cv->writeStandardHeader(of1); other_cv->write(of1); #endif if (rational_) { // Set end weight to 1 setBdWeight(1.0, false); other_cv->setBdWeight(1.0, true); } // Reparametrization (translatation and mult.) of other_cv->basis().knots_ . if (repar && continuity > 0) { if (rational_) { // The weights corresponding to the two coefficients close to the // joints should be equal equalBdWeights(false); other_cv->equalBdWeights(true); } } #ifdef DEBUG // DEBUG OUTPUT std::ofstream of1_2("basis1_2.g2"); writeStandardHeader(of1_2); write(of1_2); other_cv->writeStandardHeader(of1_2); other_cv->write(of1_2); #endif #ifdef DEBUG // DEBUG OUTPUT std::ofstream of2("basis2.dat"); basis_.write(of2); of2 << std::endl; other_cv->basis_.write(of2); of2 << std::endl; #endif if (continuity > -1 && rational_) { // Make sure that the rational curves have equal weights in the end int k2 = (numCoefs() - 1) * (dim_+1); double frac = rcoefs_[k2+dim_]/other_cv->rcoefs_[dim_]; int kn = other_cv->numCoefs(); for (int j=0; j<kn*(dim_+1); ++j) other_cv->rcoefs_[j] *= frac; } #ifdef DEBUG // DEBUG OUTPUT std::ofstream of3("basis3.dat"); basis_.write(of3); of3 << std::endl; other_cv->basis_.write(of3); of3 << std::endl; #endif if (repar && continuity > 0) { double sum1 = 0; double sum2 = 0; if (rational_) { int k2 = (numCoefs() - 1) * (dim_+1); int k1 = (numCoefs() - 2) * (dim_+1); for (int j = 0; j < dim_; ++j) { double t0 = (rcoefs_[k2 + j] - rcoefs_[k1 + j])*rcoefs_[k2 + dim_] - rcoefs_[k2 + j]*(rcoefs_[k2 + dim_] - rcoefs_[k1 + dim_]); sum1 += t0*t0; } sum1 = sqrt(sum1)/(rcoefs_[k2+dim_]*rcoefs_[k2+dim_]); k2 = dim_+1; k1 = 0; for (int j = 0; j < dim_; ++j) { double t0 = (other_cv->rcoefs_[k2 + j] - other_cv->rcoefs_[k1 + j])*other_cv->rcoefs_[k1 + dim_] - other_cv->rcoefs_[k1 + j]*(other_cv->rcoefs_[k2 + dim_] - other_cv->rcoefs_[k1 + dim_]); sum2 += t0*t0; } sum2 = sqrt(sum2)/(other_cv->rcoefs_[k1+dim_]*other_cv->rcoefs_[k1+dim_]); } else { for (int j = 0; j < dim_; ++j) { sum1 += (coefs_[(numCoefs() - 1) * dim_ + j] - coefs_[(numCoefs() - 2) * dim_ + j]) * (coefs_[(numCoefs() - 1) * dim_ + j] - coefs_[(numCoefs() - 2) * dim_ + j]); sum2 += (other_cv->coefs_[dim_ + j] - other_cv->coefs_[j]) * (other_cv->coefs_[dim_ + j] - other_cv->coefs_[j]); } sum1 = sqrt(sum1); sum2 = sqrt(sum2); } #ifdef DEBUG // DEBUG OUTPUT std::ofstream of4("basis4.dat"); basis_.write(of4); of4 << std::endl; other_cv->basis_.write(of4); of4 << std::endl; #endif if (sum1 > 1.0e-14) { // @@sbr We should have a universal noise-tolerance. double del1 = basis_.begin()[numCoefs()] - basis_.begin()[numCoefs() - 1]; double del2 = other_cv->basis_.begin()[order()] - other_cv->basis_.begin()[order() - 1]; double k = sum2*del1/(sum1*del2); other_cv->basis_.rescale(endparam(), endparam() + k * (other_cv->basis_.begin() [other_cv->numCoefs() + order() - 1] - other_cv->basis_.startparam())); } else { MESSAGE("Curve seems to be degenerated in end pt!"); } } else { other_cv->basis_.rescale(endparam(), endparam() + (other_cv->basis_.begin() [other_cv->numCoefs() + order() - 1] - other_cv->basis_.startparam())); } // Join the curve-segments (i.e. set endpoints equal), given that... if (continuity != -1) { for (int j = 0; j < dim_; ++j) { other_cv->coefs_[j] = coefs_[(numCoefs() - 1)*dim_ + j] = (coefs_[(numCoefs() - 1)*dim_ + j] + other_cv->coefs_[j])/2; } } double tpar = basis_.endparam(); int ti = numCoefs() + order() - 1; // Index of last occurence of tpar. #ifdef DEBUG // DEBUG OUTPUT std::ofstream of5("basis5.dat"); basis_.write(of5); of5 << std::endl; other_cv->basis_.write(of5); of5 << std::endl; #endif // Add other_cv's coefs. if (rational_) { // Ensure identity of the weight in the joint other_cv->setBdWeight(rcoefs_[rcoefs_.size()-1], true); rcoefs_.insert(rcoefs_end(), other_cv->rcoefs_begin(), other_cv->rcoefs_end()); } else coefs_.insert(coefs_end(), other_cv->coefs_begin(), other_cv->coefs_end()); #ifdef DEBUG // DEBUG OUTPUT std::ofstream of("basis.dat"); basis_.write(of); of << std::endl; other_cv->basis_.write(of); of << std::endl; #endif // Make an updated basis_ . std::vector<double> new_knotvector; std::copy(basis_.begin(), basis_.end(), std::back_inserter(new_knotvector)); std::copy(other_cv->basis_.begin() + order(), other_cv->basis_.end(), std::back_inserter(new_knotvector)); basis_ = BsplineBasis(order(), new_knotvector.begin(), new_knotvector.end()); if (rational_) updateCoefsFromRcoefs(); SplineCurve orig_curve = *this; // Save curve for later estimates. // Obtain wanted smoothness. int i; try { for (i = 0; i < continuity + 1; ++i) removeKnot(tpar); } catch (...) { // Leave the knots } // Estimate distance between curve and smoothed curve: // Raise (copy of) smoothed curve to original spline space // and calculate max distance between corresponding spline-coefs. SplineCurve raised_smooth_curve = *this; std::vector<double> knots; for (i = 0; i < continuity + 1; ++i) knots.push_back(tpar); raised_smooth_curve.insertKnot(knots); double sum, root_sum; dist = 0; for (i = std::max(0, ti - (continuity + 1) - order()); i < ti - order() + continuity + 1; ++i) { sum = 0; for (int j = 0; j < dim_; ++j) sum += (orig_curve.coefs_[i * dim_ + j] - raised_smooth_curve.coefs_[i * dim_ + j]) * (orig_curve.coefs_[i * dim_ + j] - raised_smooth_curve.coefs_[i * dim_ + j]); // to avoid use of the max function, which is likely to cause // trouble with the Microsoft Visual C++ Compiler, the following // two lines are added, and the third one is commented out. root_sum = sqrt(sum); dist = dist > root_sum ? dist : root_sum; //dist = std::max(dist, sqrt(sum)); } }
//=========================================================================== void SplineInterpolator::interpolate(int num_points, int dimension, const double* param_start, const double* data_start, std::vector<double>& coefs) //=========================================================================== { // Check that we have reasonable conditions set and a good param sequence ALWAYS_ERROR_IF(ctype_ == None, "No end conditions set."); for (int i = 1; i < num_points; ++i) { ALWAYS_ERROR_IF(param_start[i] <= param_start[i-1], "Parameter sequence must be strictly increasing."); } ALWAYS_ERROR_IF(ctype_ == Hermite && (dimension != start_tangent_->dimension() || dimension != end_tangent_->dimension()), "In Hermite interpolation, the end tangents must have the " "same dimension as the interpolation data."); // First we make a knot vector and define the spline space int additional_coefs = (ctype_ == Free) ? 0 : (((ctype_ == NaturalAtStart && end_tangent_.get() == 0) || (ctype_ == NaturalAtEnd && start_tangent_.get() == 0)) ? 1 : 2); int num_coefs = num_points + additional_coefs; int order = std::min(4, num_coefs); ALWAYS_ERROR_IF(num_coefs < 2,"Insufficient number of points."); std::vector<double> knots; knots.reserve(num_coefs + order); knots.insert(knots.end(), order, param_start[0]); if (ctype_ == Free) { knots.insert(knots.end(), param_start + 2, param_start + num_points - 2); } else if (ctype_ == NaturalAtStart) { if (end_tangent_.get() != 0) knots.insert(knots.end(), param_start + 1, param_start + num_points - 1); else knots.insert(knots.end(), param_start + 1, param_start + num_points - 2); } else if (ctype_ == NaturalAtEnd) { if (start_tangent_.get() != 0) knots.insert(knots.end(), param_start + 1, param_start + num_points - 1); else knots.insert(knots.end(), param_start + 2, param_start + num_points - 1); } else { // ctype_ == Natural or Hermite knots.insert(knots.end(), param_start + 1, param_start + num_points - 1); } knots.insert(knots.end(), order, param_start[num_points-1]); basis_ = BsplineBasis(num_coefs, order, &knots[0]); // Create the interpolation matrix. // The first and last row (equation) depends on the boundary // conditions (for Hermite and Natural conditions) or are // nonexisting (for Free conditions, the matrix has two fewer // rows/equations). // The interpolation matrix is a tridiagonal matrix in the Free // and Hermite cases, and a tridiagonal matrix plus two elements // in the Natural case. // This code dates from the time when this code was dependent on the NEWMAT // library. This is no longer the case, and the code has been substituted with // the code below. However, it is kept here for reference reasons. //---------------------NEWMAT dependent----------------------------- // #if 0 // cerr << "NEWMAT DEPENDENT" << endl; // BandMatrix A(num_coefs, 2, 2); // A = 0; // ColumnVector c(num_coefs); // ColumnVector b(num_coefs); // // Boundary conditions // if (ctype_ == Hermite) { // // Derivative conditions // double tmp[8]; // basis_.computeBasisValues(param_start[0], tmp, 1); // A.element(0, 0) = tmp[1]; // derivative of first B-spline // A.element(0, 1) = tmp[3]; // derivative of second B-spline // basis_.computeBasisValues(param_start[num_points-1], tmp, 1); // A.element(num_coefs - 1, num_coefs - 2) = tmp[5]; // A.element(num_coefs - 1, num_coefs - 1) = tmp[7]; // // Boundary element conditions // A.element(1, 0) = 1.0; // A.element(num_coefs - 2, num_coefs - 1) = 1.0; // } else if (ctype_ == Natural) { // // Derivative conditions // double tmp[12]; // basis_.computeBasisValues(param_start[0], tmp, 2); // A.element(0, 0) = tmp[2]; // second derivative of first B-spline // A.element(0, 1) = tmp[5]; // A.element(0, 2) = tmp[8]; // basis_.computeBasisValues(param_start[num_points-1], tmp, 2); // A.element(num_coefs - 1, num_coefs - 3) = tmp[5]; // A.element(num_coefs - 1, num_coefs - 2) = tmp[8]; // A.element(num_coefs - 1, num_coefs - 1) = tmp[11]; // // Boundary element conditions // A.element(1, 0) = 1.0; // A.element(num_coefs - 2, num_coefs - 1) = 1.0; // } else if (ctype_ == NaturalAtStart) { // // Derivative condition // double tmp[12]; // basis_.computeBasisValues(param_start[0], tmp, 2); // A.element(0, 0) = tmp[2]; // second derivative of first B-spline // A.element(0, 1) = tmp[5]; // A.element(0, 2) = tmp[8]; // if (end_tangent_.get() != 0) { // double tmp[8]; // basis_.computeBasisValues(param_start[num_points-1], tmp, 1); // A.element(num_coefs - 1, num_coefs - 2) = tmp[5]; // A.element(num_coefs - 1, num_coefs - 1) = tmp[7]; // A.element(num_coefs - 2, num_coefs - 1) = 1.0; // } else { // A.element(num_coefs - 1, num_coefs - 1) = 1.0; // } // // Boundary element conditions // A.element(1, 0) = 1.0; // } else if (ctype_ == NaturalAtEnd) { // // Derivative condition // double tmp[12]; // basis_.computeBasisValues(param_start[num_points-1], tmp, 2); // A.element(num_coefs - 1, num_coefs - 3) = tmp[5]; // A.element(num_coefs - 1, num_coefs - 2) = tmp[8]; // A.element(num_coefs - 1, num_coefs - 1) = tmp[11]; // if (start_tangent_.get() != 0) { // basis_.computeBasisValues(param_start[0], tmp, 1); // A.element(0, 0) = tmp[1]; // derivative of first B-spline // A.element(0, 1) = tmp[3]; // derivative of second B-spline // A.element(1, 0) = 1.0; // } else { // A.element(0, 0) = 1.0; // } // // Boundary element conditions // A.element(num_coefs - 2, num_coefs - 1) = 1.0; // } else if (ctype_ == Free) { // // Boundary element conditions // A.element(0, 0) = 1.0; // A.element(num_coefs - 1, num_coefs - 1) = 1.0; // } else { // THROW("Unknown boundary condition type: " << ctype_); // } // // Interior conditions // int rowoffset = ((ctype_ == Free) || // ((ctype_ == NaturalAtEnd) && (start_tangent_.get() == 0)) ? 1 : 2); // int j; // for (j = 0; j < num_points-2; ++j) { // double tmp[4]; // basis_.computeBasisValues(param_start[j+1], tmp, 0); // int column = 1 + (basis_.lastKnotInterval() - 4); // A.element(j + rowoffset, column) = tmp[0]; // A.element(j + rowoffset, column + 1) = tmp[1]; // A.element(j + rowoffset, column + 2) = tmp[2]; // A.element(j + rowoffset, column + 3) = tmp[3]; // } // // Now we are ready to repeatedly solve Ac = b for # = dimension different // // right-hand-sides b. // BandLUMatrix ALUfact = A; // Computes LU factorization // ALWAYS_ERROR_IF(ALUfact.IsSingular(), // "Matrix is singular! This should never happen!"); // coefs.resize(dimension*num_coefs); // for (int dd = 0; dd < dimension; ++dd) { // // Make the b vector // int offset = (ctype_ == Free || // ((ctype_ == NaturalAtEnd) && start_tangent_.get() == 0) ? 0 : 1); // if (ctype_ == Hermite) { // b.element(0) = (*start_tangent_)[dd]; // b.element(num_coefs - 1) = (*end_tangent_)[dd]; // } else if (ctype_ == Natural) { // b.element(0) = 0; // b.element(num_coefs - 1) = 0; // } else if (ctype_ == NaturalAtStart) { // b.element(0) = 0; // if (end_tangent_.get() != 0) // b.element(num_coefs - 1) = (*end_tangent_)[dd]; // } else if (ctype_ == NaturalAtEnd) { // b.element(num_coefs - 1) = 0; // if (start_tangent_.get() != 0) // b.element(0) = (*start_tangent_)[dd]; // } // for (j = 0; j < num_points; ++j) { // b.element(j+offset) = data_start[j*dimension + dd]; // } // // Solve // c = ALUfact.i() * b; // // Copy results // for (j = 0; j < num_coefs; ++j) { // coefs[j*dimension + dd] = c.element(j); // } // } // -------------------NEWMAT INDEPENDENT------------------------------ //#else vector<vector<double> > A(num_coefs, vector<double>(num_coefs, 0)); vector<vector<double> > b(num_coefs, vector<double>(dimension)); double tmp[12]; // boundary conditions switch (ctype_) { case Hermite: basis_.computeBasisValues(param_start[0], tmp, 1); A[0][0] = tmp[1]; // derivative of first B-spline A[0][1] = tmp[3]; // derivative of second B-spline basis_.computeBasisValues(param_start[num_points-1], tmp, 1); A[num_coefs - 1][num_coefs - 2] = tmp[5]; A[num_coefs - 1][num_coefs - 1] = tmp[7]; // Boundary element conditions A[1][0] = 1.0; A[num_coefs - 2][num_coefs - 1] = 1.0; break; case Natural: // Derivative conditions basis_.computeBasisValues(param_start[0], tmp, 2); A[0][0] = tmp[2]; // second derivative of first B-spline A[0][1] = tmp[5]; A[0][2] = tmp[8]; basis_.computeBasisValues(param_start[num_points-1], tmp, 2); A[num_coefs - 1][num_coefs - 3] = tmp[5]; A[num_coefs - 1][num_coefs - 2] = tmp[8]; A[num_coefs - 1][num_coefs - 1] = tmp[11]; // Boundary element conditions A[1][0] = 1.0; A[num_coefs - 2][num_coefs - 1] = 1.0; break; case NaturalAtStart: basis_.computeBasisValues(param_start[0], tmp, 2); A[0][0] = tmp[2]; // second derivative of first B-spline A[0][1] = tmp[5]; A[0][2] = tmp[8]; if (end_tangent_.get() != 0) { double tmp[8]; basis_.computeBasisValues(param_start[num_points-1], tmp, 1); A[num_coefs - 1][num_coefs - 2] = tmp[5]; A[num_coefs - 1][num_coefs - 1] = tmp[7]; A[num_coefs - 2][num_coefs - 1] = 1.0; } else { A[num_coefs - 1][num_coefs - 1] = 1.0; } // Boundary element conditions A[1][0] = 1.0; break; case NaturalAtEnd: basis_.computeBasisValues(param_start[num_points-1], tmp, 2); A[num_coefs - 1][num_coefs - 3] = tmp[5]; A[num_coefs - 1][num_coefs - 2] = tmp[8]; A[num_coefs - 1][num_coefs - 1] = tmp[11]; if (start_tangent_.get() != 0) { basis_.computeBasisValues(param_start[0], tmp, 1); A[0][0] = tmp[1]; // derivative of first B-spline A[0][1] = tmp[3]; // derivative of second B-spline A[1][0] = 1.0; } else { A[0][0] = 1.0; } // Boundary element conditions A[num_coefs - 2][num_coefs - 1] = 1.0; break; case Free: // Boundary element conditions A[0][0] = 1.0; A[num_coefs - 1][num_coefs - 1] = 1.0; break; default: THROW("Unknown boundary condition type." << ctype_); } // interior conditions int rowoffset = ((ctype_ == Free) || ((ctype_ == NaturalAtEnd) && (start_tangent_.get() == 0)) ? 1 : 2); int j; for (j = 0; j < num_points-2; ++j) { basis_.computeBasisValues(param_start[j+1], tmp, 0); int column = 1 + (basis_.lastKnotInterval() - 4); A[j + rowoffset][column] = tmp[0]; A[j + rowoffset][column + 1] = tmp[1]; A[j + rowoffset][column + 2] = tmp[2]; A[j + rowoffset][column + 3] = tmp[3]; } // make the b vectors boundary condition int offset = (ctype_ == Free || ((ctype_ == NaturalAtEnd) && start_tangent_.get() == 0) ? 0 : 1); switch(ctype_) { case Hermite: copy(start_tangent_->begin(), start_tangent_->end(), b[0].begin()); copy(end_tangent_->begin(), end_tangent_->end(), b[num_coefs-1].begin()); break; case Natural: b[0] = b[num_coefs-1] = vector<double>(dimension, 0); break; case NaturalAtStart: b[0] = vector<double>(dimension, 0); if (end_tangent_.get() != 0) copy(end_tangent_->begin(), end_tangent_->end(), b[num_coefs-1].begin()); break; case NaturalAtEnd: b[num_coefs-1] = vector<double>(dimension, 0); if (start_tangent_.get() != 0) copy(start_tangent_->begin(), start_tangent_->end(), b[0].begin()); break; default: // do nothing break; } // fill in interior of the b vectors for (j = 0; j < num_points; ++j) { copy(&(data_start[j * dimension]), &(data_start[(j+1) * dimension]), b[j + offset].begin()); } // computing the unknown vector A c = b. b is overwritten by this unknown vector LUsolveSystem(A, num_coefs, &b[0]); // writing result to coefficients coefs.resize(dimension*num_coefs); for (int j = 0; j < num_coefs; ++j) { copy(b[j].begin(), b[j].end(), &coefs[j*dimension]); } //#endif }
//=========================================================================== void SplineInterpolator::makeBasis(const std::vector<double>& params, const std::vector<int>& tangent_index, int order) //=========================================================================== { int tsize = (int)tangent_index.size(); int nmb_points = (int)params.size(); DEBUG_ERROR_IF(nmb_points + tsize < order, "Mismatch between interpolation conditions and order."); int i; int ti = 0; // Variable used to iterate through the tangent_index. int di; int k2 = order / 2; // Half the order. int kstop; // Control variable of the loop. int kpar; // Values in params are not repeated; derivatives are given by tangent_index. int nmb_int_cond = nmb_points + tsize; double startparam = params[0]; double endparam = params[nmb_points - 1]; // We remember start and end multiplicities. int start_mult = (tsize != 0 && tangent_index[0] == 0) ? 2 : 1; int end_mult = (tsize != 0 && tangent_index[tsize - 1] == nmb_points - 1) ? 2 : 1; std::vector<double> knots(order + nmb_int_cond); int ki = 0; // Variabel used when setting values in the knot vector; // We make the knot vector k-regular in the startparam. for (ki = 0; ki < order; ++ki) knots[ki] = startparam; double dummy1, dummy2; if (order % 2 == 0) { // Even order: place the internal knots at the parameter values. dummy1 = 0.5 * (params[1] + startparam); dummy2 = 0.5 * (params[nmb_points - 2] + endparam); if (dummy1 == dummy2) { dummy1 = (params[1] + startparam + startparam)/3.0; dummy2 = (params[nmb_points - 2] + endparam + endparam)/3.0; } if (start_mult > 1) ++ti; for (i = 0; i < start_mult - k2; ++i, ++ki) knots[ki] = dummy1; // #ifdef _MSC_VER // int kss = (k2-start_mult > 0) ? k2-start_mult : 0; // int kse = (k2-end_mult > 0) ? k2-end_mult : 0; // #else int kss = std::max(0, k2 - start_mult); int kse = std::max(0, k2 - end_mult); // #endif for (kpar = 1 + kss, kstop = nmb_points - 1 - kse; kpar < kstop; kpar++, ++ki) { knots[ki] = params[kpar]; di = ki; if (tsize != 0) while (tangent_index[ti] == kpar) { knots[++ki] = knots[di]; ++ti; } } for (i = 0; i < end_mult - k2; ++i, ++ki) knots[ki] = dummy2; } else { double delta; // The order is odd: place internal knots between parameter values. if (start_mult > 1) ++ti; if (start_mult - k2 > 0) { delta = (params[1] - startparam) / (start_mult - k2 + 1); dummy1 = startparam + delta; for (i = 0; i < start_mult - k2; ++i, dummy1 += delta, ++ki) knots[ki] = dummy1; } // #ifdef _MSC_VER // int kss = (k2-start_mult > 0) ? k2-start_mult : 0; // int kse = (k2-end_mult > 0) ? k2-end_mult : 0; // #else int kss = std::max(0, k2 - start_mult); int kse = std::max(0, k2 - end_mult); // #endif for (kpar = start_mult + kss, kstop = nmb_points - end_mult - kse; kpar < kstop; kpar++, ++ki) { knots[ki] = 0.5 *(params[kpar] + params[kpar + 1]); di = ki; if (tsize != 0) while (tangent_index[ti] == kpar) { knots[++ki] = knots[di]; ++ti; } } if (end_mult - k2 > 0) { // delta = (endparam - params[nmb_points -end_mult - 1]) / delta = (endparam - params[nmb_points - 2]) / (end_mult - k2 + 1); // dummy2 = params[nmb_points -end_mult - 1] +delta; dummy2 = params[nmb_points - 2] +delta; for (i = 0; i < end_mult - k2; ++i, dummy2 += delta, ++ki) knots[ki] = dummy2; } } for (ki = 0; ki < order; ki++) knots[nmb_int_cond+ki] = endparam; basis_ = BsplineBasis(order, knots.begin(), knots.end()); basis_set_ = true; }