process.c 8.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250
  1. /*****************************************************************************/
  2. /*** ***/
  3. /*** process() ***/
  4. /*** Reads in a raster map row by row for processing. ***/
  5. /*** Jo Wood, Project ASSIST, V1.0 7th February 1993 ***/
  6. /*** ***/
  7. /*****************************************************************************/
  8. #include <stdlib.h>
  9. #include <grass/gis.h>
  10. #include <grass/glocale.h>
  11. #include <grass/gmath.h>
  12. #include "param.h"
  13. #include "nrutil.h"
  14. void process(void)
  15. {
  16. /*--------------------------------------------------------------------------*/
  17. /* INITIALISE */
  18. /*--------------------------------------------------------------------------*/
  19. DCELL *row_in, /* Buffer large enough to hold `wsize' */
  20. *row_out = NULL, /* raster rows. When GRASS reads in a */
  21. /* raster row, each element is of type */
  22. /* DCELL */
  23. *window_ptr, /* Stores local terrain window. */
  24. centre; /* Elevation of central cell in window. */
  25. CELL *featrow_out = NULL; /* store features in CELL */
  26. struct Cell_head region; /* Structure to hold region information */
  27. int nrows, /* Will store the current number of */
  28. ncols, /* rows and columns in the raster. */
  29. row, col, /* Counts through each row and column */
  30. /* of the input raster. */
  31. wind_row, /* Counts through each row and column */
  32. wind_col, /* of the local neighbourhood window. */
  33. *index_ptr; /* Row permutation vector for LU decomp. */
  34. double **normal_ptr, /* Cross-products matrix. */
  35. *obs_ptr, /* Observed vector. */
  36. temp; /* Unused */
  37. double *weight_ptr; /* Weighting matrix for observed values. */
  38. /*--------------------------------------------------------------------------*/
  39. /* GET RASTER AND WINDOW DETAILS */
  40. /*--------------------------------------------------------------------------*/
  41. G_get_window(&region); /* Fill out the region structure (the */
  42. /* geographical limits etc.) */
  43. nrows = G_window_rows(); /* Find out the number of rows and */
  44. ncols = G_window_cols(); /* columns of the raster. */
  45. if ((region.ew_res / region.ns_res >= 1.01) || /* If EW and NS resolns are */
  46. (region.ns_res / region.ew_res >= 1.01)) { /* >1% different, warn user. */
  47. G_warning(_("E-W and N-S grid resolutions are different. Taking average."));
  48. resoln = (region.ns_res + region.ew_res) / 2;
  49. }
  50. else
  51. resoln = region.ns_res;
  52. /*--------------------------------------------------------------------------*/
  53. /* RESERVE MEMORY TO HOLD Z VALUES AND MATRICES */
  54. /*--------------------------------------------------------------------------*/
  55. row_in = (DCELL *) G_malloc(ncols * sizeof(DCELL) * wsize);
  56. /* Reserve `wsize' rows of memory. */
  57. if (mparam != FEATURE)
  58. row_out = G_allocate_raster_buf(DCELL_TYPE); /* Initialise output row buffer. */
  59. else
  60. featrow_out = G_allocate_raster_buf(CELL_TYPE); /* Initialise output row buffer. */
  61. window_ptr = (DCELL *) G_malloc(SQR(wsize) * sizeof(DCELL));
  62. /* Reserve enough memory for local wind. */
  63. weight_ptr = (double *)G_malloc(SQR(wsize) * sizeof(double));
  64. /* Reserve enough memory weights matrix. */
  65. normal_ptr = dmatrix(0, 5, 0, 5); /* Allocate memory for 6*6 matrix */
  66. index_ptr = ivector(0, 5); /* and for 1D vector holding indices */
  67. obs_ptr = dvector(0, 5); /* and for 1D vector holding observed z */
  68. /* ---------------------------------------------------------------- */
  69. /* - CALCULATE LEAST SQUARES COEFFICIENTS - */
  70. /* ---------------------------------------------------------------- */
  71. /*--- Calculate weighting matrix. ---*/
  72. find_weight(weight_ptr);
  73. /* Initial coefficients need only be found once since they are
  74. constant for any given window size. The only element that
  75. changes is the observed vector (RHS of normal equations). */
  76. /*--- Find normal equations in matrix form. ---*/
  77. find_normal(normal_ptr, weight_ptr);
  78. /*--- Apply LU decomposition to normal equations. ---*/
  79. if (constrained) {
  80. G_ludcmp(normal_ptr, 5, index_ptr, &temp);
  81. /* To constrain the quadtratic
  82. through the central cell, ignore
  83. the calculations involving the
  84. coefficient f. Since these are
  85. all in the last row and column of
  86. the matrix, simply redimension. */
  87. /* disp_matrix(normal_ptr,obs_ptr,obs_ptr,5);
  88. */
  89. }
  90. else {
  91. G_ludcmp(normal_ptr, 6, index_ptr, &temp);
  92. /* disp_matrix(normal_ptr,obs_ptr,obs_ptr,6);
  93. */
  94. }
  95. /*--------------------------------------------------------------------------*/
  96. /* PROCESS INPUT RASTER AND WRITE OUT RASTER LINE BY LINE */
  97. /*--------------------------------------------------------------------------*/
  98. if (mparam != FEATURE)
  99. for (wind_row = 0; wind_row < EDGE; wind_row++)
  100. G_put_raster_row(fd_out, row_out, DCELL_TYPE); /* Write out the edge cells as NULL. */
  101. else
  102. for (wind_row = 0; wind_row < EDGE; wind_row++)
  103. G_put_raster_row(fd_out, featrow_out, CELL_TYPE); /* Write out the edge cells as NULL. */
  104. for (wind_row = 0; wind_row < wsize - 1; wind_row++)
  105. G_get_raster_row(fd_in, row_in + (wind_row * ncols), wind_row,
  106. DCELL_TYPE);
  107. /* Read in enough of the first rows to */
  108. /* allow window to be examined. */
  109. for (row = EDGE; row < (nrows - EDGE); row++) {
  110. G_percent(row + 1, nrows - EDGE, 2);
  111. G_get_raster_row(fd_in, row_in + ((wsize - 1) * ncols), row + EDGE,
  112. DCELL_TYPE);
  113. for (col = EDGE; col < (ncols - EDGE); col++) {
  114. /* Find central z value */
  115. centre = *(row_in + EDGE * ncols + col);
  116. for (wind_row = 0; wind_row < wsize; wind_row++)
  117. for (wind_col = 0; wind_col < wsize; wind_col++)
  118. /* Express all window values relative */
  119. /* to the central elevation. */
  120. *(window_ptr + (wind_row * wsize) + wind_col) =
  121. *(row_in + (wind_row * ncols) + col + wind_col -
  122. EDGE) - centre;
  123. /*--- Use LU back substitution to solve normal equations. ---*/
  124. find_obs(window_ptr, obs_ptr, weight_ptr);
  125. /* disp_wind(window_ptr);
  126. disp_matrix(normal_ptr,obs_ptr,obs_ptr,6);
  127. */
  128. if (constrained) {
  129. G_lubksb(normal_ptr, 5, index_ptr, obs_ptr);
  130. /*
  131. disp_matrix(normal_ptr,obs_ptr,obs_ptr,5);
  132. */
  133. }
  134. else {
  135. G_lubksb(normal_ptr, 6, index_ptr, obs_ptr);
  136. /*
  137. disp_matrix(normal_ptr,obs_ptr,obs_ptr,6);
  138. */
  139. }
  140. /*--- Calculate terrain parameter based on quad. coefficients. ---*/
  141. if (mparam == FEATURE)
  142. *(featrow_out + col) = (CELL) feature(obs_ptr);
  143. else
  144. *(row_out + col) = param(mparam, obs_ptr);
  145. if (mparam == ELEV)
  146. *(row_out + col) += centre; /* Add central elevation back */
  147. }
  148. if (mparam != FEATURE)
  149. G_put_raster_row(fd_out, row_out, DCELL_TYPE); /* Write the row buffer to the output */
  150. /* raster. */
  151. else /* write FEATURE to CELL */
  152. G_put_raster_row(fd_out, featrow_out, CELL_TYPE); /* Write the row buffer to the output */
  153. /* raster. */
  154. /* 'Shuffle' rows down one, and read in */
  155. /* one new row. */
  156. for (wind_row = 0; wind_row < wsize - 1; wind_row++)
  157. for (col = 0; col < ncols; col++)
  158. *(row_in + (wind_row * ncols) + col) =
  159. *(row_in + ((wind_row + 1) * ncols) + col);
  160. }
  161. for (wind_row = 0; wind_row < EDGE; wind_row++) {
  162. if (mparam != FEATURE)
  163. G_put_raster_row(fd_out, row_out, DCELL_TYPE); /* Write out the edge cells as NULL. */
  164. else
  165. G_put_raster_row(fd_out, featrow_out, CELL_TYPE); /* Write out the edge cells as NULL. */
  166. }
  167. /*--------------------------------------------------------------------------*/
  168. /* FREE MEMORY USED TO STORE RASTER ROWS, LOCAL WINDOW AND MATRICES */
  169. /*--------------------------------------------------------------------------*/
  170. G_free(row_in);
  171. if (mparam != FEATURE)
  172. G_free(row_out);
  173. else
  174. G_free(featrow_out);
  175. G_free(window_ptr);
  176. free_dmatrix(normal_ptr, 0, 5, 0, 5);
  177. free_dvector(obs_ptr, 0, 5);
  178. free_ivector(index_ptr, 0, 5);
  179. }