NAG Library Manual, Mark 29.3
Interfaces:  FL   CL   CPP   AD 

NAG CL Interface Introduction
Example description
/* nag_wav_dim3_coeff_ins (c09fzc) Example Program.
 *
 * Copyright 2023 Numerical Algorithms Group.
 *
 * Mark 29.3, 2023.
 */
#include <math.h>
#include <nag.h>
#include <stdio.h>

#define A(I, J, K) a[(K - 1) * lda * sda + (J - 1) * lda + I - 1]
#define AN(I, J, K) an[(K - 1) * lda * sda + (J - 1) * lda + I - 1]
#define B(I, J, K) b[(K - 1) * ldb * sdb + (J - 1) * ldb + I - 1]
#define D(I, J, K) d[(K - 1) * ldd * sdd + (J - 1) * ldd + I - 1]

int main(void) {
  /* Scalars */
  Integer exit_status = 0;
  Integer lstate = 1, lseed = 1;
  Integer i, j, k, lda, ldb, ldd, lenc, m, n, fr, mnfr, nf;
  Integer nwcn, nwct, nwcfr;
  Integer nwl, subid, denoised, cindex, ilev, sda, sdb, sdd, kk;
  Nag_BaseRNG genid;
  double mse, thresh, var, xmu;
  /* Arrays */
  char mode[25], wavnam[25];
  double *a = 0, *an = 0, *b = 0, *c = 0, *d = 0, *x = 0;
  Integer *dwtlvm = 0, *dwtlvn = 0, *dwtlvfr = 0, *state = 0;
  Integer icomm[260], seed[1];
  /* Nag Types */
  Nag_Wavelet wavnamenum;
  Nag_WaveletMode modenum;
  Nag_MatrixType matrix = Nag_GeneralMatrix;
  Nag_OrderType order = Nag_ColMajor;
  Nag_DiagType diag = Nag_NonUnitDiag;
  NagError fail;

  INIT_FAIL(fail);

  printf("nag_wav_dim3_coeff_ins (c09fzc) Example Program Results\n\n");
  /* Skip heading in data file and read problem parameters. */
  scanf("%*[^\n] %" NAG_IFMT "%" NAG_IFMT "%" NAG_IFMT "%*[^\n] ", &m, &n, &fr);
  scanf("%24s%24s%*[^\n] ", wavnam, mode);

  printf("MLDWT :: Wavelet  : %s\n", wavnam);
  printf("         End mode : %s\n", mode);
  printf("         m  : %4" NAG_IFMT "\n", m);
  printf("         n  : %4" NAG_IFMT "\n", n);
  printf("         fr : %4" NAG_IFMT "\n\n", fr);

  /* Allocate arrays to hold the original data, A, original data plus noise,
   * AN, reconstruction using denoised coefficients, B, and randomly generated
   * noise, X.
   */
  lda = m;
  ldb = m;
  sda = n;
  sdb = n;
  if (!(a = NAG_ALLOC((sda) * (lda) * (fr), double)) ||
      !(an = NAG_ALLOC((sda) * (lda) * (fr), double)) ||
      !(b = NAG_ALLOC((sdb) * (ldb) * (fr), double)) ||
      !(x = NAG_ALLOC((m * n * fr), double))) {
    printf("Allocation failure\n");
    exit_status = 1;
    goto END;
  }

  /* nag_enum_name_to_value (x04nac).
   * Converts NAG enum member name to value.
   */
  wavnamenum = (Nag_Wavelet)nag_enum_name_to_value(wavnam);
  modenum = (Nag_WaveletMode)nag_enum_name_to_value(mode);

  /* Read in the original data. */
  for (k = 1; k <= fr; k++) {
    for (i = 1; i <= m; i++) {
      for (j = 1; j <= n; j++)
        scanf("%lf", &A(i, j, k));
      scanf("%*[^\n]");
    }
    scanf("%*[^\n]");
  }

  /* Output the original data. */
  printf("Input data :\n");
  fflush(stdout);
  for (k = 1; k <= fr; k++) {
    nag_file_print_matrix_real_gen_comp(order, matrix, diag, m, n, &A(1, 1, k),
                                        lda, "%11.4e", " ", Nag_NoLabels, 0,
                                        Nag_NoLabels, 0, 80, 0, 0, &fail);
    if (fail.code != NE_NOERROR) {
      printf("Error from nag_file_print_matrix_real_gen_comp (x04cbc).\n%s\n",
             fail.message);
      exit_status = 2;
      goto END;
    }
    printf("\n");
    fflush(stdout);
  }

  /* Set up call to nag_rand_dist_normal (g05skc) in order to create some random
   * noise from a normal distribution to add to the original data.
   * Initial call to RNG initializer to get size of STATE array.
   */
  seed[0] = 642521;
  genid = Nag_MersenneTwister;
  subid = 0;
  if (!(state = NAG_ALLOC((lstate), Integer))) {
    printf("Allocation failure\n");
    exit_status = -1;
    goto END;
  }

  /* nag_rand_init_repeat (g05kfc).
   * Query the size of state.
   */
  lstate = 0;
  nag_rand_init_repeat(genid, subid, seed, lseed, state, &lstate, &fail);
  if (fail.code != NE_NOERROR) {
    printf("Error from nag_rand_init_repeat (g05kfc).\n%s\n", fail.message);
    exit_status = 3;
    goto END;
  }

  /* Reallocate STATE. */
  NAG_FREE(state);
  if (!(state = NAG_ALLOC((lstate), Integer))) {
    printf("Allocation failure\n");
    exit_status = 4;
    goto END;
  }

  /* nag_rand_init_repeat (g05kfc).
   * Initialize the generator to a repeatable sequence.
   */
  nag_rand_init_repeat(genid, subid, seed, lseed, state, &lstate, &fail);
  if (fail.code != NE_NOERROR) {
    printf("Error from nag_rand_init_repeat (g05kfc).\n%s\n", fail.message);
    exit_status = 5;
    goto END;
  }

  /* Set the distribution parameters for the random noise. */
  xmu = 0.0;
  var = 0.1E-3;

  /* Generate the noise variates. */

  /* nag_rand_dist_normal (g05skc).
   * Generates a vector of pseudorandom numbers from a Normal distribution.
   */
  mnfr = n * m * fr;
  nag_rand_dist_normal(mnfr, xmu, var, state, x, &fail);
  if (fail.code != NE_NOERROR) {
    printf("Error from nag_rand_dist_normal (g05skc).\n%s\n", fail.message);
    exit_status = 6;
    goto END;
  }

  /* Add the noise to the original input and save in AN */
  kk = 0;
  for (k = 1; k <= n; k++) {
    for (j = 1; j <= n; j++) {
      for (i = 1; i <= m; i++) {
        AN(i, j, k) = A(i, j, k) + x[kk];
        kk = kk + 1;
      }
    }
  }

  /* Output the noisy data */
  printf("Input data plus noise :\n");
  fflush(stdout);
  for (k = 1; k <= fr; k++) {
    nag_file_print_matrix_real_gen_comp(order, matrix, diag, m, n, &AN(1, 1, k),
                                        lda, "%11.4e", " ", Nag_NoLabels, 0,
                                        Nag_NoLabels, 0, 80, 0, 0, &fail);
    if (fail.code != NE_NOERROR) {
      printf("Error from nag_file_print_matrix_real_gen_comp (x04cbc).\n%s\n",
             fail.message);
      exit_status = 7;
      goto END;
    }
    printf("\n");
    fflush(stdout);
  }

  /* nag_wav_dim3_init (c09acc).
   * Three-dimensional wavelet filter initialization.
   */
  nag_wav_dim3_init(wavnamenum, Nag_MultiLevel, modenum, m, n, fr, &nwl, &nf,
                    &nwct, &nwcn, &nwcfr, icomm, &fail);
  if (fail.code != NE_NOERROR) {
    printf("Error from nag_wav_dim3_init (c09acc).\n%s\n", fail.message);
    exit_status = 8;
    goto END;
  }

  /* Allocate arrays to hold the coefficients, c, and the dimensions
   * of the coefficients at each level, dwtlvm, dwtlvn.
   */
  lenc = nwct;
  if (!(c = NAG_ALLOC((lenc), double)) ||
      !(dwtlvm = NAG_ALLOC((nwl), Integer)) ||
      !(dwtlvn = NAG_ALLOC((nwl), Integer)) ||
      !(dwtlvfr = NAG_ALLOC((nwl), Integer))) {
    printf("Allocation failure\n");
    exit_status = 9;
    goto END;
  }

  /* Perform a forwards multi-level transform on the noisy data. */

  /* nag_wav_dim3_multi_fwd (c09fcc).
   * Two-dimensional multi-level discrete wavelet transform.
   */
  nag_wav_dim3_multi_fwd(m, n, fr, an, lda, sda, lenc, c, nwl, dwtlvm, dwtlvn,
                         dwtlvfr, icomm, &fail);
  if (fail.code != NE_NOERROR) {
    printf("Error from nag_wav_dim3_multi_fwd (c09fcc).\n%s\n", fail.message);
    exit_status = 10;
    goto END;
  }

  /* Reconstruct without thresholding of detail coefficients. */

  /* nag_wav_dim3_mxolap_multi_inv (c09fdc).
   * Two-dimensional inverse multi-level discrete wavelet transform.
   */
  nag_wav_dim3_mxolap_multi_inv(nwl, lenc, c, m, n, fr, b, ldb, sdb, icomm,
                                &fail);
  if (fail.code != NE_NOERROR) {
    printf("Error from nag_wav_dim3_mxolap_multi_inv (c09fdc).\n%s\n",
           fail.message);
    exit_status = 11;
    goto END;
  }

  /* Calculate the Mean Square Error of the noisy reconstruction. */
  mse = 0.0;
  for (k = 1; k <= fr; k++)
    for (j = 1; j <= n; j++)
      for (i = 1; i <= m; i++)
        mse = mse + pow((A(i, j, k) - B(i, j, k)), 2);
  mse = mse / (double)(m * n * fr);
  printf("Without denoising Mean Square Error is %11.4e\n\n", mse);

  /* Now perform the denoising by extracting each of the detail
   * coefficients at each level and applying hard thresholding
   * Allocate a 2D array to hold the detail coefficients
   */
  ldd = dwtlvm[nwl - 1];
  sdd = dwtlvn[nwl - 1];
  if (!(d = NAG_ALLOC((ldd) * (sdd) * (dwtlvfr[nwl - 1]), double))) {
    printf("Allocation failure\n");
    exit_status = 12;
    goto END;
  }

  /* Calculate the threshold based on VisuShrink denoising. */
  thresh = sqrt(var) * sqrt(2. * log((double)(m * n * fr)));
  denoised = 0;
  /* For each level */
  for (ilev = nwl; ilev >= 1; ilev -= 1) {
    /* Select detail coefficients */
    for (cindex = 1; cindex <= 7; cindex++) {
      /* Extract coefficients into the 2D array d */

      /* nag_wav_dim3_coeff_ext (c09fyc).
       * Three-dimensional discrete wavelet transform coefficient extraction.
       */
      nag_wav_dim3_coeff_ext(ilev, cindex, lenc, c, d, ldd, sdd, icomm, &fail);
      if (fail.code != NE_NOERROR) {
        printf("Error from nag_wav_dim3_coeff_ext (c09fyc).\n%s\n",
               fail.message);
        exit_status = 13;
        goto END;
      }

      /* Perform the hard thresholding operation */
      for (k = 1; k <= dwtlvfr[nwl - ilev]; k++)
        for (j = 1; j <= dwtlvn[nwl - ilev]; j++)
          for (i = 1; i <= dwtlvm[nwl - ilev]; i++)
            if (fabs(D(i, j, k)) < thresh) {
              D(i, j, k) = 0.0;
              denoised = denoised + 1;
            }

      /* Insert the denoised coefficients back into c. */

      /* nag_wav_dim3_coeff_ins (c09fzc).
       * Three-dimensional discrete wavelet transform coefficient insertion.
       */
      nag_wav_dim3_coeff_ins(ilev, cindex, lenc, c, d, ldd, sdd, icomm, &fail);
      if (fail.code != NE_NOERROR) {
        printf("Error from nag_wav_dim3_coeff_ins (c09fzc).\n%s\n",
               fail.message);
        exit_status = 14;
        goto END;
      }
    }
  }

  /* Output the number of coefficients that were set to zero */
  printf("Number of coefficients denoised is %4" NAG_IFMT " out of %4" NAG_IFMT
         "\n\n",
         denoised, nwct - dwtlvm[0] * dwtlvn[0] * dwtlvfr[0]);
  fflush(stdout);

  /* Reconstruct original data following thresholding of detail coefficients */

  /* nag_wav_dim3_mxolap_multi_inv (c09fdc).
   * Three-dimensional inverse multi-level discrete wavelet transform.
   */
  nag_wav_dim3_mxolap_multi_inv(nwl, lenc, c, m, n, fr, b, ldb, sdb, icomm,
                                &fail);
  if (fail.code != NE_NOERROR) {
    printf("Error from nag_wav_dim3_mxolap_multi_inv (c09fdc).\n%s\n",
           fail.message);
    exit_status = 15;
    goto END;
  }

  /* Calculate the Mean Square Error of the denoised reconstruction. */
  mse = 0.0;
  for (k = 1; k <= n; k++)
    for (j = 1; j <= n; j++)
      for (i = 1; i <= m; i++)
        mse = mse + pow((A(i, j, k) - B(i, j, k)), 2);
  mse = mse / (double)(m * n * fr);
  printf("With denoising Mean Square Error is %11.4e \n\n", mse);

  /* Output the denoised reconstruction. */
  printf("Reconstruction of denoised input :\n");
  fflush(stdout);
  for (k = 1; k <= fr; k++) {
    nag_file_print_matrix_real_gen_comp(order, matrix, diag, m, n, &B(1, 1, k),
                                        ldb, "%11.4e", " ", Nag_NoLabels, 0,
                                        Nag_NoLabels, 0, 80, 0, 0, &fail);
    if (fail.code != NE_NOERROR) {
      printf("Error from nag_file_print_matrix_real_gen_comp (x04cbc).\n%s\n",
             fail.message);
      exit_status = 16;
      goto END;
    }
    if (k < fr)
      printf("\n");
    fflush(stdout);
  }

END:
  NAG_FREE(a);
  NAG_FREE(an);
  NAG_FREE(b);
  NAG_FREE(c);
  NAG_FREE(d);
  NAG_FREE(x);
  NAG_FREE(dwtlvm);
  NAG_FREE(dwtlvn);
  NAG_FREE(dwtlvfr);
  NAG_FREE(state);
  return exit_status;
}