IPP Software Navigation Tools IPP Links Communication Pan-STARRS Links

Ignore:
Timestamp:
Oct 21, 2008, 4:48:50 PM (18 years ago)
Author:
Paul Price
Message:

Moving what is now BIQUADRATIC interpolation (used to be BICUBE) into a 2D kernel again: preserves Gene's definition.

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/psLib/src/imageops/psImageInterpolate.c

    r20306 r20311  
    77 *  @author Paul Price, IfA
    88 *
    9  *  @version $Revision: 1.23 $ $Name: not supported by cvs2svn $
    10  *  @date $Date: 2008-10-22 02:10:37 $
     9 *  @version $Revision: 1.24 $ $Name: not supported by cvs2svn $
     10 *  @date $Date: 2008-10-22 02:48:50 $
    1111 *
    1212 *  Copyright 2004-2007 Institute for Astronomy, University of Hawaii
     
    5555}
    5656
     57#if 0
    5758/// Generate a biquadratic interpolation kernel
    5859/// This reduces Gene's original made-up 2D kernel to 1D
     
    6768    kernel[2] = x + xx;
    6869}
     70#else
     71/// Generate a biquadratic interpolation kernel
     72static inline void interpolationKernelBiquadratic(
     73    psF32 kernel[kernelSizes[PS_INTERPOLATE_BIQUADRATIC]][kernelSizes[PS_INTERPOLATE_BIQUADRATIC]], // Kernel
     74    float xFrac, float yFrac // Fraction of pixel
     75    )
     76{
     77    double xxFrac = xFrac * xFrac / 6.0;
     78    double yyFrac = yFrac * yFrac / 6.0;
     79    double xyFrac = 0.25 * xFrac * yFrac;
     80    xFrac /= 6.0;
     81    yFrac /= 6.0;
     82    kernel[0][0] = - 1.0/9.0 - xFrac - yFrac + xxFrac + yyFrac + xyFrac;
     83    kernel[0][1] = 2.0/9.0 - yFrac - 2.0 * xxFrac + yyFrac;
     84    kernel[0][2] = - 1.0/9.0 + xFrac - yFrac + xxFrac + yyFrac - xyFrac;
     85    kernel[1][0] = 2.0/9.0 - xFrac + xxFrac - 2.0 * yyFrac;
     86    kernel[1][1] = 5.0/9.0 - 2.0 * xxFrac - 2.0 * yyFrac;
     87    kernel[1][2] = 2.0/9.0 + xFrac + xxFrac - 2.0 * yyFrac;
     88    kernel[2][0] = - 1.0/9.0 - xFrac + yFrac + xxFrac + yyFrac - xyFrac;
     89    kernel[2][1] = 2.0/9.0 + yFrac - 2.0 * xxFrac + yyFrac;
     90    kernel[2][2] = - 1.0/9.0 + xFrac + yFrac + xxFrac + yyFrac + xyFrac;
     91}
     92#endif
     93
    6994
    7095#if 0
     
    135160
    136161
    137 // Generate interpolation kernel
     162// Generate 1D interpolation kernel
    138163static inline void interpolationKernel(psF32 *kernel, // Kernel vector to populate
    139164                                       float frac, // Fraction of pixel
     
    148173      case PS_INTERPOLATE_BILINEAR:
    149174        interpolationKernelBilinear(kernel, frac);
    150         break;
    151       case PS_INTERPOLATE_BIQUADRATIC:
    152         interpolationKernelBiquadratic(kernel, frac);
    153175        break;
    154176      case PS_INTERPOLATE_GAUSS:
     
    160182        interpolationKernelLanczos(kernel, kernelSizes[mode], frac);
    161183        break;
     184      case PS_INTERPOLATE_BIQUADRATIC:  // 2D kernel
    162185      default:
    163186        psAbort("Unsupported interpolation mode: %x", mode);
     
    207230    psImage *kernel2 = NULL;            // Kernel^2
    208231    psVector *sumKernel2 = NULL;        // Sum of kernel^2
    209     if (numKernels > 0) {
    210         int size = kernelSizes[mode];      // Size of kernel
    211 
    212         kernel = psImageAlloc(size, numKernels, PS_TYPE_F32); // Kernel
    213         kernel2 = psImageAlloc(size, numKernels, PS_TYPE_F32); // Kernel^2
    214         sumKernel2 = psVectorAlloc(numKernels, PS_TYPE_F32); // Sum of kernel^2
    215 
    216         for (int i = 0; i < numKernels; i++) {
    217             float frac = i / (float)numKernels;   // Fraction of shift
    218             interpolationKernel(kernel->data.F32[i], frac, mode);
    219 
    220             float sum = 0.0;               // Sum of kernel^2
    221             for (int j = 0; j < size; j++) {
    222                 sum += kernel2->data.F32[i][j] = PS_SQR(kernel->data.F32[i][j]);
     232
     233    switch (mode) {
     234      case PS_INTERPOLATE_BILINEAR:
     235      case PS_INTERPOLATE_GAUSS:
     236      case PS_INTERPOLATE_LANCZOS2:
     237      case PS_INTERPOLATE_LANCZOS3:
     238      case PS_INTERPOLATE_LANCZOS4:
     239        if (numKernels > 0) {
     240            int size = kernelSizes[mode];      // Size of kernel
     241
     242            kernel = psImageAlloc(size, numKernels, PS_TYPE_F32); // Kernel
     243            kernel2 = psImageAlloc(size, numKernels, PS_TYPE_F32); // Kernel^2
     244            sumKernel2 = psVectorAlloc(numKernels, PS_TYPE_F32); // Sum of kernel^2
     245
     246            for (int i = 0; i < numKernels; i++) {
     247                float frac = i / (float)numKernels;   // Fraction of shift
     248                interpolationKernel(kernel->data.F32[i], frac, mode);
     249
     250                float sum = 0.0;               // Sum of kernel^2
     251                for (int j = 0; j < size; j++) {
     252                    sum += kernel2->data.F32[i][j] = PS_SQR(kernel->data.F32[i][j]);
     253                }
     254                sumKernel2->data.F32[i] = sum;
    223255            }
    224             sumKernel2->data.F32[i] = sum;
    225         }
    226     }
     256        }
     257        break;
     258      case PS_INTERPOLATE_BIQUADRATIC:
     259        // 2D kernel, would cost too much memory to pre-calculate
     260        break;
     261      default:
     262        psAbort("Unsupported interpolation mode: %x", mode);
     263    }
     264
    227265    interp->kernel = kernel;
    228266    interp->kernel2 = kernel2;
     
    316354
    317355// Set up the kernel parameters; defines some useful values
    318 #define INTERPOLATE_KERNEL_SETUP(X, Y) \
     356#define INTERPOLATE_SETUP(X, Y) \
    319357    /* Central pixel is the pixel below the point of interest */ \
    320358    int xCentral = floor((X) - 0.5), yCentral = floor((Y) - 0.5); /* Central pixel */ \
    321     float xFrac = (X) - xCentral - 0.5, yFrac = (Y) - yCentral - 0.5; /* Fraction of pixel */
    322 
    323 
    324 // Interpolation engine for (separable) interpolation kernels
    325 static psImageInterpolateStatus interpolateKernel(double *imageValue, double *varianceValue,
    326                                                   psMaskType *maskValue, float x, float y,
    327                                                   const psImageInterpolation *interp)
     359    float xFrac = (X) - xCentral - 0.5, yFrac = (Y) - yCentral - 0.5; /* Fraction of pixel */ \
     360    bool xExact = fabsf(xFrac) < FLT_EPSILON, yExact = fabsf(yFrac) < FLT_EPSILON; /* Are shifts exact? */ \
     361
     362// Determine the range of the kernel; defines some useful values
     363#define INTERPOLATE_RANGE() \
     364    int xLast = image->numCols - 1, yLast = image->numRows - 1; /* Last pixels of image */ \
     365    /* Proper range of kernel on image; not all may be defined */ \
     366    int xStart = xCentral - (size - 1) / 2, xStop = xCentral + size / 2; \
     367    int yStart = yCentral - (size - 1) / 2, yStop = yCentral + size / 2; \
     368    if (xStart > xLast || xStop < 0 || yStart > yLast || yStop < 0) { \
     369        /* Interpolation kernel is entirely off the image */ \
     370        INTERPOLATE_OFF(); /* Returns */ \
     371    } \
     372    int xMin, xMax, yMin, yMax;         /* Minimum and maximum valid pixels on image */ \
     373    int iMin, iMax, jMin, jMax;         /* Minimum and maximum valid pixels on kernel */ \
     374    bool offImage = false;              /* At least one pixel of the kernel is off the image */ \
     375    /* XXX Haven't got single dimension exact shifts implemented properly yet. */ \
     376    /* XXX When it is ready, note that the limit checks ([ij]{Min,Max}) below are wrong: */ \
     377    /* should probably refer to [xy]{Start,Stop}. */ \
     378    xExact = yExact = false; \
     379    if (xExact) { \
     380        if (iMin > 0 || iMax < 1) { /* This is wrong */ \
     381            INTERPOLATE_OFF(); /* Returns */ \
     382        } \
     383        iMin = (size - 1) / 2; \
     384        iMax = iMin + 1; \
     385    } else { \
     386        if (xStart < 0) { \
     387            xMin = 0; \
     388            iMin = -xStart; \
     389            offImage = true; \
     390        } else { \
     391            xMin = xStart; \
     392            iMin = 0; \
     393        } \
     394        if (xStop > xLast) { \
     395            xMax = xLast; \
     396            iMax = xLast - xStart + 1; \
     397            offImage = true; \
     398        } else { \
     399            xMax = xStop; \
     400            iMax = size; \
     401        } \
     402    } \
     403    if (yExact) { \
     404        if (jMin > 0 || jMax < 1) { /* This is wrong */ \
     405            INTERPOLATE_OFF(); /* Returns */ \
     406        } \
     407        jMin = (size - 1) / 2; \
     408        jMax = jMin + 1; \
     409    } else { \
     410        if (yStart < 0) { \
     411            yMin = 0; \
     412            jMin = -yStart; \
     413            offImage = true; \
     414        } else { \
     415            yMin = yStart; \
     416            jMin = 0; \
     417        } \
     418        if (yStop > yLast) { \
     419            yMax = yLast; \
     420            jMax = yLast - yStart + 1; \
     421            offImage = true; \
     422        } else { \
     423            yMax = yStop; \
     424            jMax = size; \
     425        } \
     426    }
     427
     428// Determine the result of the interpolation after all the math has been done
     429#define INTERPOLATE_RESULT() \
     430    psImageInterpolateStatus status = PS_INTERPOLATE_STATUS_ERROR; /* Status of interpolation */ \
     431    *imageValue = sumImage / sumKernel; \
     432    if (wantVariance) { \
     433        *varianceValue = sumVariance / sumKernel2; \
     434    } \
     435    if (sumBad == 0) { \
     436        /* Completely good pixel */ \
     437        status = PS_INTERPOLATE_STATUS_GOOD; \
     438    } else if (sumBad < PS_SQR(interp->poorFrac) * (sumBad + sumKernel2) ) { \
     439        /* Some pixels masked: poor pixel */ \
     440        if (haveMask && maskValue) { \
     441            *maskValue |= interp->poorMask; \
     442        } \
     443        status = PS_INTERPOLATE_STATUS_POOR; \
     444    } else { \
     445        /* Many pixels (or a few important pixels) masked: bad pixel */ \
     446        if (haveMask && maskValue) { \
     447            *maskValue |= interp->badMask; \
     448        } \
     449        status = PS_INTERPOLATE_STATUS_BAD; \
     450    }
     451
     452// Interpolation engine for separable interpolation kernels
     453static psImageInterpolateStatus interpolateSeparable(double *imageValue, double *varianceValue,
     454                                                     psMaskType *maskValue, float x, float y,
     455                                                     const psImageInterpolation *interp)
    328456{
    329457    // Parameters have been checked by psImageInterpolate()
     
    339467    // Kernel basics
    340468    int size = kernelSizes[mode];       // Size of kernel
    341     INTERPOLATE_KERNEL_SETUP(x, y);
    342 
    343     // Extent of the kernel on the image
    344     bool xExact = fabsf(xFrac) < FLT_EPSILON, yExact = fabsf(yFrac) < FLT_EPSILON; // Are shifts exact?
     469    INTERPOLATE_SETUP(x, y);
    345470    if (xExact && yExact) {
    346         // Both shifts are exact
    347471        return interpolateFlat(imageValue, varianceValue, maskValue, x, y, interp);
    348472    }
    349 
    350     int xLast = image->numCols - 1, yLast = image->numRows - 1; // Last pixels of image
    351     int xStart = xCentral - (size - 1) / 2, xStop = xCentral + size / 2; // Start and stop of kernel on image
    352     int yStart = yCentral - (size - 1) / 2, yStop = yCentral + size / 2; // Start and stop of kernel on image
    353     if (xStart > xLast || xStop < 0 || yStart > yLast || yStop < 0) {
    354         // Interpolation kernel is entirely off the image
    355         INTERPOLATE_OFF();              // Returns
    356     }
    357     int xMin, xMax, yMin, yMax;         // Minimum and maximum valid pixels on image
    358     int iMin, iMax, jMin, jMax;         // Minimum and maximum valid pixels on kernel
    359     bool offImage = false;              // At least one pixel of the kernel is off the image
    360 
    361     // XXX Haven't got single dimension exact shifts implemented properly yet
    362     // XXX When it is ready, note that the check below limits
    363     xExact = yExact = false;
    364     if (xExact) {
    365         if (iMin > 0 || iMax < 1) {
    366             INTERPOLATE_OFF();          // Returns
    367         }
    368         iMin = (size - 1) / 2;
    369         iMax = iMin + 1;
    370     } else {
    371         if (xStart < 0) {
    372             xMin = 0;
    373             iMin = -xStart;
    374             offImage = true;
    375         } else {
    376             xMin = xStart;
    377             iMin = 0;
    378         }
    379         if (xStop > xLast) {
    380             xMax = xLast;
    381             iMax = xLast - xStart + 1;
    382             offImage = true;
    383         } else {
    384             xMax = xStop;
    385             iMax = size;
    386         }
    387     }
    388     if (yExact) {
    389         if (jMin > 0 || jMax < 1) {
    390             INTERPOLATE_OFF();          // Returns
    391         }
    392         jMin = (size - 1) / 2;
    393         jMax = jMin + 1;
    394     } else {
    395         if (yStart < 0) {
    396             yMin = 0;
    397             jMin = -yStart;
    398             offImage = true;
    399         } else {
    400             yMin = yStart;
    401             jMin = 0;
    402         }
    403         if (yStop > yLast) {
    404             yMax = yLast;
    405             jMax = yLast - yStart + 1;
    406             offImage = true;
    407         } else {
    408             yMax = yStop;
    409             jMax = size;
    410         }
    411     }
     473    INTERPOLATE_RANGE();
    412474
    413475    // Get the appropriate kernels
     
    510572    }
    511573
    512 #define INTERPOLATE_KERNEL_CASE(TYPE) \
     574#define INTERPOLATE_SEPARATE_CASE(TYPE) \
    513575  case PS_TYPE_##TYPE: { \
    514576      if (wantVariance) { \
     
    630692
    631693    switch (image->type.type) {
     694        INTERPOLATE_SEPARATE_CASE(U8);
     695        INTERPOLATE_SEPARATE_CASE(U16);
     696        INTERPOLATE_SEPARATE_CASE(U32);
     697        INTERPOLATE_SEPARATE_CASE(U64);
     698        INTERPOLATE_SEPARATE_CASE(S8);
     699        INTERPOLATE_SEPARATE_CASE(S16);
     700        INTERPOLATE_SEPARATE_CASE(S32);
     701        INTERPOLATE_SEPARATE_CASE(S64);
     702        INTERPOLATE_SEPARATE_CASE(F32);
     703        INTERPOLATE_SEPARATE_CASE(F64);
     704      default:
     705        psError(PS_ERR_BAD_PARAMETER_TYPE, true, "Unrecognised type for image: %x",
     706                image->type.type);
     707        return PS_INTERPOLATE_STATUS_ERROR;
     708    }
     709
     710    INTERPOLATE_RESULT();
     711
     712    psFree(xKernelNew);
     713    psFree(yKernelNew);
     714    psFree(xKernel2New);
     715    psFree(yKernel2New);
     716
     717    return status;
     718}
     719
     720// Interpolation engine for (separable) interpolation kernels
     721static psImageInterpolateStatus interpolateKernel(double *imageValue, double *varianceValue,
     722                                                  psMaskType *maskValue, float x, float y,
     723                                                  const psImageInterpolation *interp)
     724{
     725    // Parameters have been checked by psImageInterpolate()
     726
     727    psImageInterpolateMode mode = interp->mode; // Mode of interpolation
     728    const psImage *image = interp->image; // Image of interest
     729    const psImage *mask = interp->mask; // Image mask
     730    const psImage *variance = interp->variance; // Image variance
     731    psMaskType maskVal = interp->maskVal; // Value to mask
     732    bool wantVariance = variance && varianceValue; // Does the user want the variance value?
     733    bool haveMask = mask && maskVal; // Does the user want the variance value?
     734
     735    // Kernel basics
     736    int size = kernelSizes[mode];       // Size of kernel
     737    INTERPOLATE_SETUP(x, y);
     738    if (xExact && yExact) {
     739        return interpolateFlat(imageValue, varianceValue, maskValue, x, y, interp);
     740    }
     741    INTERPOLATE_RANGE();
     742
     743    // Get the appropriate kernels
     744    psF32 kernel[size][size];
     745    psAssert(mode == PS_INTERPOLATE_BIQUADRATIC, "Mode is %x", mode);
     746    interpolationKernelBiquadratic(kernel, xFrac, yFrac);
     747
     748    float sumKernel = 0.0;              // Sum of the kernel
     749    double sumBad = 0.0;                // Sum of bad kernel-squared contributions
     750    double sumImage = 0.0;              // Sum of image multiplied by kernel
     751    double sumVariance = 0.0;           // Sum of variance multiplied by kernel-squared
     752    float sumKernel2 = 0.0;             // Sum of kernel^2
     753
     754    // Add contributions in an area outside the image
     755#define INTERPOLATE_KERNEL_ADD_OFFIMAGE() { \
     756        sumBad += PS_SQR(kernel[j][i]); \
     757    }
     758
     759    // Measure kernel contribution from outside image
     760    if (offImage) {
     761        if (!yExact) {
     762            // Bottom rows
     763            for (int j = 0; j < jMin; j++) {
     764                for (int i = 0; i < size; i++) {
     765                    INTERPOLATE_KERNEL_ADD_OFFIMAGE();
     766                }
     767            }
     768        }
     769        if (!xExact) {
     770            // Two sides
     771            for (int j = jMin; j < jMax; j++) {
     772                for (int i = 0; i < iMin; i++) {
     773                    INTERPOLATE_KERNEL_ADD_OFFIMAGE();
     774                }
     775                for (int i = iMax; i < size; i++) {
     776                    INTERPOLATE_KERNEL_ADD_OFFIMAGE();
     777                }
     778            }
     779        }
     780        if (!yExact) {
     781            // Top rows
     782            for (int j = jMax; j < size; j++) {
     783                for (int i = 0; i < size; i++) {
     784                    INTERPOLATE_KERNEL_ADD_OFFIMAGE();
     785                }
     786            }
     787        }
     788    }
     789
     790#define INTERPOLATE_KERNEL_CASE(TYPE) \
     791  case PS_TYPE_##TYPE: { \
     792      if (wantVariance) { \
     793          if (haveMask) { \
     794              /* Variance and mask */ \
     795              for (int j = jMin, yPix = yMin; j < jMax; j++, yPix++) { \
     796                  for (int i = iMin, xPix = xMin; i < iMax; i++, xPix++) { \
     797                      float kernelValue = kernel[j][i]; /* Value of kernel */ \
     798                      float kernelValue2 = PS_SQR(kernelValue); /* Square of kernel */ \
     799                      if (mask->data.PS_TYPE_MASK_DATA[yPix][xPix] & maskVal) { \
     800                          sumBad += kernelValue2; \
     801                      } else { \
     802                          sumImage += kernelValue * image->data.TYPE[yPix][xPix]; \
     803                          sumVariance += kernelValue2 * variance->data.TYPE[yPix][xPix]; \
     804                          sumKernel += kernelValue; \
     805                          sumKernel2 += kernelValue2; \
     806                      } \
     807                  } \
     808              } \
     809              \
     810          } else { \
     811              /* Variance, no mask */ \
     812              for (int j = jMin, yPix = yMin; j < jMax; j++, yPix++) { \
     813                  for (int i = iMin, xPix = xMin; i < iMax; i++, xPix++) { \
     814                      float kernelValue = kernel[j][i]; /* Value of kernel */ \
     815                      float kernelValue2 = PS_SQR(kernelValue); /* Square of kernel */ \
     816                      sumImage += kernelValue * image->data.TYPE[yPix][xPix]; \
     817                      sumVariance += kernelValue2 * variance->data.TYPE[yPix][xPix]; \
     818                      sumKernel += kernelValue; \
     819                      sumKernel2 += kernelValue2; \
     820                  } \
     821              } \
     822          } \
     823      } else if (haveMask) { \
     824          /* Mask, no variance */ \
     825          for (int j = jMin, yPix = yMin; j < jMax; j++, yPix++) { \
     826              for (int i = iMin, xPix = xMin; i < iMax; i++, xPix++) { \
     827                  float kernelValue = kernel[j][i]; /* Value of kernel */ \
     828                  float kernelValue2 = PS_SQR(kernelValue); /* Square of kernel */ \
     829                  if (mask->data.PS_TYPE_MASK_DATA[yPix][xPix] & maskVal) { \
     830                      sumBad += kernelValue2; \
     831                  } else { \
     832                      sumImage += kernelValue * image->data.TYPE[yPix][xPix]; \
     833                      sumKernel += kernelValue; \
     834                      sumKernel2 += kernelValue2; \
     835                  } \
     836              } \
     837          } \
     838      } else { \
     839          /* Neither variance nor mask */ \
     840          for (int j = jMin, yPix = yMin; j < jMax; j++, yPix++) { \
     841              for (int i = iMin, xPix = xMin; i < iMax; i++, xPix++) { \
     842                  float kernelValue = kernel[j][i]; /* Value of kernel */ \
     843                  sumImage += kernelValue * image->data.TYPE[yPix][xPix]; \
     844                  sumKernel += kernelValue; \
     845              } \
     846          } \
     847      } \
     848    } \
     849    break;
     850
     851    switch (image->type.type) {
    632852        INTERPOLATE_KERNEL_CASE(U8);
    633853        INTERPOLATE_KERNEL_CASE(U16);
     
    646866    }
    647867
    648     psImageInterpolateStatus status = PS_INTERPOLATE_STATUS_ERROR; // Status of interpolation
    649     *imageValue = sumImage / sumKernel;
    650     if (wantVariance) {
    651         *varianceValue = sumVariance / sumKernel2;
    652     }
    653     if (sumBad == 0) {
    654         // Completely good pixel
    655         status = PS_INTERPOLATE_STATUS_GOOD;
    656     } else if (sumBad < PS_SQR(interp->poorFrac) * (sumBad + sumKernel2) ) {
    657         // Some pixels masked: poor pixel
    658         if (haveMask && maskValue) {
    659             *maskValue |= interp->poorMask;
    660         }
    661         status = PS_INTERPOLATE_STATUS_POOR;
    662     } else {
    663         // Many pixels (or a few important pixels) masked: bad pixel
    664         if (haveMask && maskValue) {
    665             *maskValue |= interp->badMask;
    666         }
    667         status = PS_INTERPOLATE_STATUS_BAD;
    668     }
    669 
    670     psFree(xKernelNew);
    671     psFree(yKernelNew);
    672     psFree(xKernel2New);
    673     psFree(yKernel2New);
     868    INTERPOLATE_RESULT();
    674869
    675870    return status;
     
    706901      case PS_INTERPOLATE_FLAT:
    707902        return interpolateFlat(imageValue, varianceValue, maskValue, x, y, interp);
     903      case PS_INTERPOLATE_BIQUADRATIC:
     904        return interpolateKernel(imageValue, varianceValue, maskValue, x, y, interp);
    708905      case PS_INTERPOLATE_BILINEAR:
    709       case PS_INTERPOLATE_BIQUADRATIC:
    710906      case PS_INTERPOLATE_GAUSS:
    711907      case PS_INTERPOLATE_LANCZOS2:
    712908      case PS_INTERPOLATE_LANCZOS3:
    713909      case PS_INTERPOLATE_LANCZOS4:
    714         return interpolateKernel(imageValue, varianceValue, maskValue, x, y, interp);
     910        return interpolateSeparable(imageValue, varianceValue, maskValue, x, y, interp);
    715911      default:
    716912        psError(PS_ERR_BAD_PARAMETER_VALUE, true,
     
    730926
    731927    // Kernel basics
    732     INTERPOLATE_KERNEL_SETUP(x, y);
     928    INTERPOLATE_SETUP(x, y);
     929    xExact = yExact = false;
    733930
    734931    psF32 xKernel[size], yKernel[size]; // Interpolation kernels
Note: See TracChangeset for help on using the changeset viewer.