IPP Software Navigation Tools IPP Links Communication Pan-STARRS Links

Ignore:
Timestamp:
Jul 13, 2010, 2:34:04 PM (16 years ago)
Author:
Paul Price
Message:

It's bothered me for a while that the noise model that we have been using is only kind of spatially variable. The variance map is intended to track the change in the noise properties as a function of position, and while it does keep track of different noise levels due to, e.g., different detectors, vignetting, etc., it has not tracked the different noise that comes from convolution kernels changing over an image. This had been an intentional oversight --- the prescription we were using for calculating the variance and covariance removed the spatial variation (by prescribing that the variance level under a convolution remain unchanged) and we were sweeping it under the carpet so we could calculate the covariance pseudo-matrix and use that for setting the absolute level of the noise.

I've now figured out how to have our cake and eat it too. By changing where the scalings are applied, the variance map can track the changing noise introduced by a varying convolution kernel, while the covariance pseudo-matrix tracks an absolute change. This fix is purely algorithmic --- it does not affect any of the variances/covariances stored on disk, nor is it affected by them. I hope that this means that any Magic recertification required is either minimal or completely unnecessary.

Again, all that has changed is where the scalings are applied.

File:
1 edited

Legend:

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

    r28405 r28667  
    3434static float imageCovarianceCalculate(const psKernel *covar, // Original covariance matrix
    3535                                      const psKernel *kernel, // Convolution kernel
    36                                       int x, int y            // Coordinates in output covariance matrix
     36                                      int x, int y,           // Coordinates in output covariance matrix
     37                                      float scale             // Scale to apply
    3738                                      )
    3839{
     
    8384    }
    8485
    85     return sum;
     86    return scale * sum;
    8687}
    8788
     
    9192    PS_ASSERT_THREAD_JOB_NON_NULL(job, false);
    9293    psAssert(job->args, "No job arguments");
    93     psAssert(job->args->n == 5, "Wrong number of job arguments: %ld", job->args->n);
     94    psAssert(job->args->n == 6, "Wrong number of job arguments: %ld", job->args->n);
    9495
    9596    psKernel *out = job->args->data[0]; // Output covariance matrix
     
    9899    int x = PS_SCALAR_VALUE(job->args->data[3], S32); // x coordinate in output covariance matrix
    99100    int y = PS_SCALAR_VALUE(job->args->data[4], S32); // y coordinate in output covariance matrix
    100 
    101     out->kernel[y][x] = imageCovarianceCalculate(covar, kernel, x, y);
     101    float scale = PS_SCALAR_VALUE(job->args->data[5], F32); // Scaling to apply
     102
     103    out->kernel[y][x] = imageCovarianceCalculate(covar, kernel, x, y, scale);
    102104
    103105    return true;
     
    127129
    128130    // Check for non-finite elements
     131    double sumKernel = 0.0, sumKernel2 = 0.0; // Sum of the kernel
    129132    for (int y = kernel->yMin; y <= kernel->yMax; y++) {
    130133        for (int x = kernel->xMin; x <= kernel->xMax; x++) {
     
    135138                return NULL;
    136139            }
     140            sumKernel += kernel->kernel[y][x];
     141            sumKernel2 += PS_SQR(kernel->kernel[y][x]);
    137142        }
    138143    }
     
    155160    int yMin = kernel->yMin - kernel->yMax + covar->yMin, yMax = kernel->yMax - kernel->yMin + covar->yMax;
    156161    psKernel *out = psKernelAlloc(xMin, xMax, yMin, yMax); // Covariance matrix for output
     162    float scale = 1.0 / sumKernel2;          // Scaling to apply
    157163
    158164    for (int y = yMin; y <= yMax; y++) {
     
    165171                PS_ARRAY_ADD_SCALAR(job->args, x, PS_TYPE_S32);
    166172                PS_ARRAY_ADD_SCALAR(job->args, y, PS_TYPE_S32);
     173                PS_ARRAY_ADD_SCALAR(job->args, scale, PS_TYPE_F32);
    167174                if (!psThreadJobAddPending(job)) {
    168175                    psFree(covar);
     
    170177                }
    171178            } else {
    172                 out->kernel[y][x] = imageCovarianceCalculate(covar, kernel, x, y);
    173             }
    174         }
    175     }
    176     psFree(covar);
     179                out->kernel[y][x] = imageCovarianceCalculate(covar, kernel, x, y, scale);
     180            }
     181        }
     182    }
    177183
    178184    if (threaded && !psThreadPoolWait(true)) {
     
    180186        return false;
    181187    }
     188
     189    psFree(covar);
    182190
    183191    return out;
     
    194202
    195203    // Check for non-finite elements
     204    double sumKernel2 = 0.0; // Sum of the squared kernel
    196205    for (int y = kernel->yMin; y <= kernel->yMax; y++) {
    197206        for (int x = kernel->xMin; x <= kernel->xMax; x++) {
     
    202211                return NAN;
    203212            }
     213            sumKernel2 += PS_SQR(kernel->kernel[y][x]);
    204214        }
    205215    }
     
    215225    }
    216226
    217     float factor = imageCovarianceCalculate(covar, kernel, 0, 0); // Covariance factor
     227    float scale = 1.0 / sumKernel2;     // Scale to apply
     228    float factor = imageCovarianceCalculate(covar, kernel, 0, 0, scale); // Covariance factor
    218229    psFree(covar);
    219230    return factor;
     
    338349        }
    339350    }
    340     psFree(covar);
    341 
    342351    if (threaded && !psThreadPoolWait(true)) {
    343352        psError(PS_ERR_UNKNOWN, false, "Error waiting for threads.");
    344353        return false;
    345354    }
     355    psFree(covar);
    346356
    347357    return out;
     
    616626    if (set && !threaded) {
    617627        {
    618             psThreadTask *task = psThreadTaskAlloc("PSLIB_IMAGE_COVARIANCE_CALCULATE", 5);
     628            psThreadTask *task = psThreadTaskAlloc("PSLIB_IMAGE_COVARIANCE_CALCULATE", 6);
    619629            task->function = &imageCovarianceCalculateThread;
    620630            psThreadTaskAdd(task);
Note: See TracChangeset for help on using the changeset viewer.