1 /*
2   This file is part of CDO. CDO is a collection of Operators to manipulate and analyse Climate model Data.
3 
4   Author: Uwe Schulzweida
5 
6 */
7 /*
8   This is a C library of the Fortran SCRIP version 1.4
9 
10   ===>>> Please send bug reports to <https://mpimet.mpg.de/cdo> <<<===
11 
12   Spherical Coordinate Remapping and Interpolation Package (SCRIP)
13   ================================================================
14 
15   SCRIP is a software package which computes addresses and weights for
16   remapping and interpolating fields between grids in spherical coordinates.
17   It was written originally for remapping fields to other grids in a coupled
18   climate model, but is sufficiently general that it can be used in other
19   applications as well. The package should work for any grid on the surface
20   of a sphere. SCRIP currently supports four remapping options:
21 
22   Conservative remapping
23   ----------------------
24   First- and second-order conservative remapping as described in
25   Jones (1999, Monthly Weather Review, 127, 2204-2210).
26 
27   Bilinear interpolation
28   ----------------------
29   Slightly generalized to use a local bilinear approximation
30   (only logically-rectangular grids).
31 
32   Bicubic interpolation
33   ----------------------
34   Similarly generalized (only logically-rectangular grids).
35 
36   Distance-weighted averaging
37   ---------------------------
38   Distance-weighted average of a user-specified number of nearest neighbor
39   values.
40 
41   Documentation
42   =============
43 
44   http://climate.lanl.gov/Software/SCRIP/SCRIPusers.pdf
45 
46 */
47 /*
48   2013-11-08 Uwe Schulzweida: split remapgrid class to src_grid and tgt_grid
49   2012-01-16 Uwe Schulzweida: alloc grid2_bound_box only for conservative remapping
50   2011-01-07 Uwe Schulzweida: Changed remap weights from 2D to 1D array
51   2009-05-25 Uwe Schulzweida: Changed restrict data type from double to int
52   2009-01-11 Uwe Schulzweida: OpenMP parallelization
53  */
54 
55 #include <cdi.h>
56 
57 #include "cdo_options.h"
58 #include "cimdOmp.h"
59 #include "process_int.h"
60 #include "cdo_wtime.h"
61 #include <mpim_grid.h>
62 #include "gridreference.h"
63 #include "remap.h"
64 
65 #define IS_REG2D_GRID(gridID) (gridInqType(gridID) == GRID_LONLAT || gridInqType(gridID) == GRID_GAUSSIAN)
66 
67 static bool remap_gen_weights = true;
68 static bool remap_write_remap = false;
69 static int remap_num_srch_bins = 180;
70 #define DEFAULT_MAX_ITER 100
71 size_t remap_max_iter = DEFAULT_MAX_ITER;  // Max iteration count for i, j iteration
72 
73 void
remap_set_int(int remapvar,int value)74 remap_set_int(int remapvar, int value)
75 {
76   // clang-format off
77   if      (remapvar == REMAP_WRITE_REMAP)   remap_write_remap = (value > 0);
78   else if (remapvar == REMAP_MAX_ITER)      remap_max_iter = value;
79   else if (remapvar == REMAP_NUM_SRCH_BINS) remap_num_srch_bins = value;
80   else if (remapvar == REMAP_GENWEIGHTS)    remap_gen_weights = (value > 0);
81   else    cdo_abort("Unsupported remap variable (%d)!", remapvar);
82   // clang-format on
83 }
84 
85 /*****************************************************************************/
86 static void
boundboxFromCorners(size_t size,size_t nc,const Varray<double> & corner_lon,const Varray<double> & corner_lat,float * bound_box)87 boundboxFromCorners(size_t size, size_t nc, const Varray<double> &corner_lon, const Varray<double> &corner_lat, float *bound_box)
88 {
89 #ifdef _OPENMP
90 #pragma omp parallel for default(none) shared(bound_box, corner_lat, corner_lon, nc, size)
91 #endif
92   for (size_t i = 0; i < size; ++i)
93     {
94       size_t i4 = i << 2;  // *4
95       size_t inc = i * nc;
96       float clat = corner_lat[inc];
97       float clon = corner_lon[inc];
98       bound_box[i4 + 0] = clat;
99       bound_box[i4 + 1] = clat;
100       bound_box[i4 + 2] = clon;
101       bound_box[i4 + 3] = clon;
102       for (size_t j = 1; j < nc; ++j)
103         {
104           clat = corner_lat[inc + j];
105           clon = corner_lon[inc + j];
106           if (clat < bound_box[i4 + 0]) bound_box[i4 + 0] = clat;
107           if (clat > bound_box[i4 + 1]) bound_box[i4 + 1] = clat;
108           if (clon < bound_box[i4 + 2]) bound_box[i4 + 2] = clon;
109           if (clon > bound_box[i4 + 3]) bound_box[i4 + 3] = clon;
110         }
111     }
112 }
113 
114 static void
boundboxFromCenter(bool lonIsCyclic,size_t size,size_t nx,size_t ny,const Varray<double> & center_lon,const Varray<double> & center_lat,float * bound_box)115 boundboxFromCenter(bool lonIsCyclic, size_t size, size_t nx, size_t ny, const Varray<double> &center_lon,
116                    const Varray<double> &center_lat, float *bound_box)
117 {
118 #ifdef _OPENMP
119 #pragma omp parallel for default(none) shared(lonIsCyclic, size, nx, ny, center_lon, center_lat, bound_box)
120 #endif
121   for (size_t n = 0; n < size; n++)
122     {
123       size_t idx[4];
124       float tmp_lats[4], tmp_lons[4];
125 
126       size_t n4 = n << 2;
127 
128       // Find N,S and NE points to this grid point
129 
130       size_t j = n / nx;
131       size_t i = n - j * nx;
132 
133       size_t ip1 = (i < (nx - 1)) ? i + 1 : lonIsCyclic ? 0 : i;
134       size_t jp1 = (j < (ny - 1)) ? j + 1 : j;
135 
136       idx[0] = n;
137       idx[1] = j * nx + ip1;    // east
138       idx[2] = jp1 * nx + ip1;  // north-east
139       idx[3] = jp1 * nx + i;    // north
140 
141       // Find N,S and NE lat/lon coords and check bounding box
142 
143       for (unsigned k = 0; k < 4; ++k) tmp_lons[k] = center_lon[idx[k]];
144       for (unsigned k = 0; k < 4; ++k) tmp_lats[k] = center_lat[idx[k]];
145 
146       bound_box[n4 + 0] = tmp_lats[0];
147       bound_box[n4 + 1] = tmp_lats[0];
148       bound_box[n4 + 2] = tmp_lons[0];
149       bound_box[n4 + 3] = tmp_lons[0];
150 
151       for (unsigned k = 1; k < 4; k++)
152         {
153           if (tmp_lats[k] < bound_box[n4 + 0]) bound_box[n4 + 0] = tmp_lats[k];
154           if (tmp_lats[k] > bound_box[n4 + 1]) bound_box[n4 + 1] = tmp_lats[k];
155           if (tmp_lons[k] < bound_box[n4 + 2]) bound_box[n4 + 2] = tmp_lons[k];
156           if (tmp_lons[k] > bound_box[n4 + 3]) bound_box[n4 + 3] = tmp_lons[k];
157         }
158     }
159 }
160 
161 LonLatPoint
remapgrid_get_lonlat(RemapGrid * grid,size_t cell_add)162 remapgrid_get_lonlat(RemapGrid *grid, size_t cell_add)
163 {
164   double lon, lat;
165 
166   if (grid->type == RemapGridType::Reg2D)
167     {
168       const auto nx = grid->dims[0];
169       const auto iy = cell_add / nx;
170       const auto ix = cell_add - iy * nx;
171       lat = grid->reg2d_center_lat[iy];
172       lon = grid->reg2d_center_lon[ix];
173       if (lon < 0) lon += PI2;
174     }
175   else
176     {
177       lat = grid->cell_center_lat[cell_add];
178       lon = grid->cell_center_lon[cell_add];
179     }
180 
181   return LonLatPoint(lon, lat);
182 }
183 
184 void
check_lon_range(const char * txt,size_t nlons,Varray<double> & lons)185 check_lon_range(const char *txt, size_t nlons, Varray<double> &lons)
186 {
187   assert(!lons.empty());
188 
189   if (txt)
190     {
191       double minval = 1.e36;
192       double maxval = -1.e36;
193 #ifdef _OPENMP
194 #pragma omp parallel for default(none) shared(nlons, lons) reduction(min : minval) reduction(max : maxval)
195 #endif
196       for (size_t i = 0; i < nlons; ++i)
197         {
198           minval = std::min(minval, lons[i]);
199           maxval = std::max(maxval, lons[i]);
200         }
201       if (minval < -PI2 || maxval > 2 * PI2)
202         cdo_warning("%s grid cell center longitudes out of range (min=%.3g/max=%.3g)!", txt, RAD2DEG * minval, RAD2DEG * maxval);
203     }
204 
205 #ifdef _OPENMP
206 #pragma omp parallel for default(none) shared(nlons, lons)
207 #endif
208   for (size_t i = 0; i < nlons; ++i)
209     {
210       // remove missing values
211       if (lons[i] < -PI2) lons[i] = 0;
212       if (lons[i] > 2 * PI2) lons[i] = PI2;
213 
214       if (lons[i] > PI2) lons[i] -= PI2;
215       if (lons[i] < 0.0) lons[i] += PI2;
216     }
217 }
218 
219 void
check_lat_range(const char * txt,size_t nlats,Varray<double> & lats)220 check_lat_range(const char *txt, size_t nlats, Varray<double> &lats)
221 {
222   assert(!lats.empty());
223 
224   if (txt)
225     {
226       double minval = 1.e36;
227       double maxval = -1.e36;
228 #ifdef _OPENMP
229 #pragma omp parallel for default(none) shared(nlats, lats) reduction(min : minval) reduction(max : maxval)
230 #endif
231       for (size_t i = 0; i < nlats; ++i)
232         {
233           minval = std::min(minval, lats[i]);
234           maxval = std::max(maxval, lats[i]);
235         }
236       if (minval < -(PIH + 0.0001) || maxval > (PIH + 0.0001))
237         cdo_warning("%s grid cell center latitudes out of range (min=%.3g/max=%.3g)!", txt, RAD2DEG * minval, RAD2DEG * maxval);
238     }
239 
240 #ifdef _OPENMP
241 #pragma omp parallel for default(none) shared(nlats, lats)
242 #endif
243   for (size_t i = 0; i < nlats; ++i)
244     {
245       if (lats[i] > PIH) lats[i] = PIH;
246       if (lats[i] < -PIH) lats[i] = -PIH;
247     }
248 }
249 
250 static void
check_lon_boundbox_range(size_t nlons,float * bound_box)251 check_lon_boundbox_range(size_t nlons, float *bound_box)
252 {
253   assert(bound_box != nullptr);
254 
255 #ifdef _OPENMP
256 #pragma omp parallel for default(none) shared(nlons, bound_box)
257 #endif
258   for (size_t n = 0; n < nlons; ++n)
259     {
260       const auto n4 = n << 2;
261       if (fabsf(bound_box[n4 + 3] - bound_box[n4 + 2]) > PI_f)
262         {
263           bound_box[n4 + 2] = 0.0f;
264           bound_box[n4 + 3] = PI2_f;
265         }
266     }
267 }
268 
269 static void
check_lat_boundbox_range(size_t nlats,float * bound_box,Varray<double> & lats)270 check_lat_boundbox_range(size_t nlats, float *bound_box, Varray<double> &lats)
271 {
272   assert(bound_box != nullptr);
273 
274 #ifdef _OPENMP
275 #pragma omp parallel for default(none) shared(nlats, bound_box, lats)
276 #endif
277   for (size_t n = 0; n < nlats; ++n)
278     {
279       const auto n4 = n << 2;
280       if ((float) lats[n] < bound_box[n4]) bound_box[n4] = -PIH_f;
281       if ((float) lats[n] > bound_box[n4 + 1]) bound_box[n4 + 1] = PIH_f;
282     }
283 }
284 
285 /*****************************************************************************/
286 
287 static void
grid_check_lat_borders_rad(size_t n,Varray<double> & ybounds)288 grid_check_lat_borders_rad(size_t n, Varray<double> &ybounds)
289 {
290   constexpr double YMAX = PIH;
291   constexpr double YLIM = 88 * DEG2RAD;
292   const auto lrev = (ybounds[0] > ybounds[n - 1]);
293   if (lrev)
294     {
295       if (ybounds[0] > ybounds[1])
296         {
297           if (ybounds[0] > YLIM) ybounds[0] = YMAX;
298           if (ybounds[n - 1] < -YLIM) ybounds[n - 1] = -YMAX;
299         }
300       else
301         {
302           if (ybounds[1] > YLIM) ybounds[1] = YMAX;
303           if (ybounds[n - 2] < -YLIM) ybounds[n - 2] = -YMAX;
304         }
305     }
306   else
307     {
308       if (ybounds[0] < ybounds[1])
309         {
310           if (ybounds[0] < -YLIM) ybounds[0] = -YMAX;
311           if (ybounds[n - 1] > YLIM) ybounds[n - 1] = YMAX;
312         }
313       else
314         {
315           if (ybounds[1] < -YLIM) ybounds[1] = -YMAX;
316           if (ybounds[n - 2] > YLIM) ybounds[n - 2] = YMAX;
317         }
318     }
319 }
320 
321 static void
convert_bounds_reg2d(size_t n,const Varray<double> & boundsIn,Varray<double> & boundsOut)322 convert_bounds_reg2d(size_t n, const Varray<double> &boundsIn, Varray<double> &boundsOut)
323 {
324   auto lrev = (boundsIn[0] > boundsIn[2 * n - 1]);
325   if (boundsIn[0] > boundsIn[1]) lrev = !lrev;
326   if (lrev)
327     {
328       boundsOut[0] = boundsIn[1];
329       for (size_t i = 0; i < n; ++i) boundsOut[i + 1] = boundsIn[2 * i];
330     }
331   else
332     {
333       boundsOut[0] = boundsIn[0];
334       for (size_t i = 0; i < n; ++i) boundsOut[i + 1] = boundsIn[2 * i + 1];
335     }
336 }
337 
338 static void
remap_define_reg2d(int gridID,RemapGrid & grid)339 remap_define_reg2d(int gridID, RemapGrid &grid)
340 {
341   const auto nx = grid.dims[0];
342   const auto ny = grid.dims[1];
343   const auto nxp1 = nx + 1;
344   const auto nyp1 = ny + 1;
345 
346   auto nxm = nx;
347   if (grid.is_cyclic) nxm++;
348 
349   if (grid.size != nx * ny) cdo_abort("Internal error, wrong dimensions!");
350 
351   grid.reg2d_center_lon.resize(nxm);
352   grid.reg2d_center_lat.resize(ny);
353 
354   grid.reg2d_center_lon[0] = 0.0;
355   grid.reg2d_center_lat[0] = 0.0;
356   gridInqXvals(gridID, grid.reg2d_center_lon.data());
357   gridInqYvals(gridID, grid.reg2d_center_lat.data());
358 
359   // Convert lat/lon units if required
360   cdo_grid_to_radian(gridID, CDI_XAXIS, nx, grid.reg2d_center_lon.data(), "grid reg2d center lon");
361   cdo_grid_to_radian(gridID, CDI_YAXIS, ny, grid.reg2d_center_lat.data(), "grid reg2d center lat");
362 
363   if (grid.reg2d_center_lon[nx - 1] < grid.reg2d_center_lon[0])
364     for (size_t i = 1; i < nx; ++i)
365       if (grid.reg2d_center_lon[i] < grid.reg2d_center_lon[i - 1]) grid.reg2d_center_lon[i] += PI2;
366 
367   if (grid.is_cyclic) grid.reg2d_center_lon[nx] = grid.reg2d_center_lon[0] + PI2;
368 
369   grid.reg2d_corner_lon.resize(nxp1);
370   grid.reg2d_corner_lat.resize(nyp1);
371 
372   if (gridInqXbounds(gridID, nullptr))
373     {
374       Varray<double> xbounds(2 * nx);
375       gridInqXbounds(gridID, xbounds.data());
376       convert_bounds_reg2d(nx, xbounds, grid.reg2d_corner_lon);
377       cdo_grid_to_radian(gridID, CDI_XAXIS, nx + 1, grid.reg2d_corner_lon.data(), "grid reg2d corner lon");
378     }
379   else
380     {
381       grid_gen_corners(nx, grid.reg2d_center_lon.data(), grid.reg2d_corner_lon.data());
382     }
383 
384   if (gridInqYbounds(gridID, nullptr))
385     {
386       Varray<double> ybounds(2 * ny);
387       gridInqYbounds(gridID, ybounds.data());
388       convert_bounds_reg2d(ny, ybounds, grid.reg2d_corner_lat);
389       cdo_grid_to_radian(gridID, CDI_YAXIS, ny + 1, grid.reg2d_corner_lat.data(), "grid reg2d corner lat");
390     }
391   else
392     {
393       grid_gen_corners(ny, grid.reg2d_center_lat.data(), grid.reg2d_corner_lat.data());
394       grid_check_lat_borders_rad(ny + 1, grid.reg2d_corner_lat);
395     }
396 }
397 
398 static void
remapDefineGrid(RemapMethod mapType,int gridID,RemapGrid & grid,const char * txt)399 remapDefineGrid(RemapMethod mapType, int gridID, RemapGrid &grid, const char *txt)
400 {
401   bool lgrid_destroy = false;
402   int gridID_gme = -1;
403 
404   const auto gridtype = gridInqType(grid.gridID);
405   if (gridtype != GRID_UNSTRUCTURED && gridtype != GRID_CURVILINEAR)
406     {
407       if (gridtype == GRID_GME)
408         {
409           gridID_gme = gridToUnstructured(grid.gridID, 1);
410           grid.nvgp = gridInqSize(gridID_gme);
411           gridID = gridDuplicate(gridID_gme);
412           gridCompress(gridID);
413           grid.luse_cell_corners = true;
414         }
415       else if (gridtype == GRID_GAUSSIAN_REDUCED)
416         {
417           lgrid_destroy = true;
418           gridID = gridToUnstructured(grid.gridID, 1);
419         }
420       else if (remap_write_remap || grid.type != RemapGridType::Reg2D)
421         {
422           lgrid_destroy = true;
423           gridID = gridToCurvilinear(grid.gridID, 1);
424         }
425     }
426 
427   auto gridsize = grid.size = gridInqSize(gridID);
428 
429   grid.dims[0] = gridInqXsize(gridID);
430   grid.dims[1] = gridInqYsize(gridID);
431   if (gridInqType(grid.gridID) != GRID_UNSTRUCTURED)
432     {
433       if (grid.dims[0] == 0) cdo_abort("%s grid without longitude coordinates!", gridNamePtr(gridInqType(grid.gridID)));
434       if (grid.dims[1] == 0) cdo_abort("%s grid without latitude coordinates!", gridNamePtr(gridInqType(grid.gridID)));
435     }
436 
437   grid.is_cyclic = (gridIsCircular(gridID) > 0);
438 
439   grid.rank = (gridInqType(gridID) == GRID_UNSTRUCTURED) ? 1 : 2;
440 
441   grid.num_cell_corners = (gridInqType(gridID) == GRID_UNSTRUCTURED) ? gridInqNvertex(gridID) : 4;
442 
443   remap_grid_alloc(mapType, grid);
444 
445   // Initialize logical mask
446 
447 #ifdef _OPENMP
448 #pragma omp parallel for default(none) shared(gridsize, grid)
449 #endif
450   for (size_t i = 0; i < gridsize; ++i) grid.mask[i] = 1;
451 
452   if (gridInqMask(gridID, nullptr))
453     {
454       std::vector<int> mask(gridsize);
455       gridInqMask(gridID, &mask[0]);
456       for (size_t i = 0; i < gridsize; ++i)
457         if (mask[i] == 0) grid.mask[i] = 0;
458     }
459 
460   if (!remap_write_remap && grid.type == RemapGridType::Reg2D) return;
461 
462   if (!gridHasCoordinates(gridID)) cdo_abort("%s grid cell center coordinates missing!", txt);
463 
464   gridInqXvals(gridID, grid.cell_center_lon.data());
465   gridInqYvals(gridID, grid.cell_center_lat.data());
466 
467   if (grid.lneed_cell_corners)
468     {
469       if (!gridHasBounds(gridID)) cdo_abort("%s grid cell corner coordinates missing!", txt);
470 
471       gridInqXbounds(gridID, grid.cell_corner_lon.data());
472       gridInqYbounds(gridID, grid.cell_corner_lat.data());
473     }
474 
475   if (gridInqType(grid.gridID) == GRID_GME) gridInqMaskGME(gridID_gme, &grid.vgpm[0]);
476 
477   // Convert lat/lon units if required
478   cdo_grid_to_radian(gridID, CDI_XAXIS, grid.size, grid.cell_center_lon.data(), "grid center lon");
479   cdo_grid_to_radian(gridID, CDI_YAXIS, grid.size, grid.cell_center_lat.data(), "grid center lat");
480   if (grid.num_cell_corners && grid.lneed_cell_corners)
481     {
482       cdo_grid_to_radian(gridID, CDI_XAXIS, grid.num_cell_corners * grid.size, grid.cell_corner_lon.data(), "grid corner lon");
483       cdo_grid_to_radian(gridID, CDI_YAXIS, grid.num_cell_corners * grid.size, grid.cell_corner_lat.data(), "grid corner lat");
484     }
485 
486   if (lgrid_destroy) gridDestroy(gridID);
487 
488   // Convert longitudes to 0,2pi interval
489   check_lon_range(txt, grid.size, grid.cell_center_lon);
490   if (grid.num_cell_corners && grid.lneed_cell_corners)
491     check_lon_range(nullptr, grid.num_cell_corners * grid.size, grid.cell_corner_lon);
492 
493   // Make sure input latitude range is within the machine values for +/- pi/2
494   check_lat_range(txt, grid.size, grid.cell_center_lat);
495   if (grid.num_cell_corners && grid.lneed_cell_corners)
496     check_lat_range(nullptr, grid.num_cell_corners * grid.size, grid.cell_corner_lat);
497 }
498 
499 // Compute bounding boxes for restricting future grid searches
500 static void
cell_bounding_boxes(RemapGrid & grid,float * cell_bound_box,int remap_grid_basis)501 cell_bounding_boxes(RemapGrid &grid, float *cell_bound_box, int remap_grid_basis)
502 {
503   if (grid.luse_cell_corners)
504     {
505       if (grid.lneed_cell_corners)
506         {
507           if (Options::cdoVerbose) cdo_print("Grid: boundboxFromCorners");
508           boundboxFromCorners(grid.size, grid.num_cell_corners, grid.cell_corner_lon, grid.cell_corner_lat, cell_bound_box);
509         }
510       else  // full grid search
511         {
512           if (Options::cdoVerbose) cdo_print("Grid: bounds missing -> full grid search!");
513 
514           const auto gridsize = grid.size;
515           for (size_t i = 0; i < gridsize; ++i)
516             {
517               cell_bound_box[i * 4] = -PIH_f;
518               cell_bound_box[i * 4 + 1] = PIH_f;
519               cell_bound_box[i * 4 + 2] = 0.0f;
520               cell_bound_box[i * 4 + 3] = PI2_f;
521             }
522         }
523     }
524   else if (remap_grid_basis == REMAP_GRID_BASIS_SRC)
525     {
526       if (Options::cdoVerbose) cdo_print("Grid: boundboxFromCenter");
527       if (grid.rank != 2) cdo_abort("Internal problem, grid rank = %d!", grid.rank);
528 
529       const auto nx = grid.dims[0];
530       const auto ny = grid.dims[1];
531       boundboxFromCenter(grid.is_cyclic, grid.size, nx, ny, grid.cell_center_lon, grid.cell_center_lat, cell_bound_box);
532     }
533 
534   if (remap_grid_basis == REMAP_GRID_BASIS_SRC || grid.lneed_cell_corners) check_lon_boundbox_range(grid.size, cell_bound_box);
535 
536   // Try to check for cells that overlap poles
537   if (remap_grid_basis == REMAP_GRID_BASIS_SRC || grid.lneed_cell_corners)
538     check_lat_boundbox_range(grid.size, cell_bound_box, grid.cell_center_lat);
539 }
540 
541 void
remap_grid_alloc(RemapMethod mapType,RemapGrid & grid)542 remap_grid_alloc(RemapMethod mapType, RemapGrid &grid)
543 {
544   if (grid.nvgp) grid.vgpm.resize(grid.nvgp);
545 
546   grid.mask.resize(grid.size);
547 
548   if (remap_write_remap || grid.type != RemapGridType::Reg2D)
549     {
550       grid.cell_center_lon.resize(grid.size);
551       grid.cell_center_lat.resize(grid.size);
552     }
553 
554   const auto needCellarea = (mapType == RemapMethod::CONSERV_SCRIP || mapType == RemapMethod::CONSERV);
555   if (needCellarea) grid.cell_area.resize(grid.size, 0.0);
556 
557   grid.cell_frac.resize(grid.size, 0.0);
558 
559   if (grid.lneed_cell_corners && grid.num_cell_corners > 0)
560     {
561       const auto nalloc = grid.num_cell_corners * grid.size;
562       grid.cell_corner_lon.resize(nalloc, 0);
563       grid.cell_corner_lat.resize(nalloc, 0);
564     }
565 }
566 
567 void
remap_grid_init(RemapGrid & grid)568 remap_grid_init(RemapGrid &grid)
569 {
570   grid.tmpgridID = -1;
571   grid.type = RemapGridType::Undefined;
572 
573   grid.num_cell_corners = 0;
574   grid.luse_cell_corners = false;
575   grid.lneed_cell_corners = false;
576 
577   grid.nvgp = 0;
578 }
579 
580 void
remap_grid_free(RemapGrid & grid)581 remap_grid_free(RemapGrid &grid)
582 {
583   varray_free(grid.vgpm);
584   varray_free(grid.mask);
585 
586   varray_free(grid.reg2d_center_lat);
587   varray_free(grid.reg2d_center_lon);
588   varray_free(grid.reg2d_corner_lat);
589   varray_free(grid.reg2d_corner_lon);
590 
591   varray_free(grid.cell_center_lat);
592   varray_free(grid.cell_center_lon);
593   varray_free(grid.cell_corner_lat);
594   varray_free(grid.cell_corner_lon);
595 
596   varray_free(grid.cell_area);
597   varray_free(grid.cell_frac);
598 
599   if (grid.tmpgridID != -1) gridDestroy(grid.tmpgridID);
600 }
601 
602 void
remap_search_init(RemapMethod mapType,RemapSearch & search,RemapGrid & src_grid,RemapGrid & tgt_grid)603 remap_search_init(RemapMethod mapType, RemapSearch &search, RemapGrid &src_grid, RemapGrid &tgt_grid)
604 {
605   extern PointSearchMethod pointSearchMethod;
606   extern CellSearchMethod cellSearchMethod;
607 
608   search.srcGrid = &src_grid;
609   search.tgtGrid = &tgt_grid;
610 
611   search.srcBins.ncells = src_grid.size;
612   search.tgtBins.ncells = tgt_grid.size;
613 
614   search.srcBins.nbins = remap_num_srch_bins;
615   search.tgtBins.nbins = remap_num_srch_bins;
616 
617   auto usePointsearch = (mapType == RemapMethod::DISTWGT);
618   if (src_grid.type != RemapGridType::Reg2D && pointSearchMethod != PointSearchMethod::latbins)
619     {
620       usePointsearch |= (mapType == RemapMethod::BILINEAR);
621       usePointsearch |= (mapType == RemapMethod::BICUBIC);
622     }
623 
624   auto useCellsearch = (mapType == RemapMethod::CONSERV)
625                        && (cellSearchMethod == CellSearchMethod::spherepart || src_grid.type == RemapGridType::Reg2D);
626 
627   const char *searchMethodStr = NULL;
628   const auto start = Options::cdoVerbose ? cdo_get_wtime() : 0.0;
629 
630   if (usePointsearch)
631     {
632       searchMethodStr = "Point search";
633       bool xIsCyclic = src_grid.is_cyclic;
634       if (src_grid.type == RemapGridType::Reg2D)
635         grid_point_search_create_reg_2d(search.gps, xIsCyclic, src_grid.dims, src_grid.reg2d_center_lon, src_grid.reg2d_center_lat);
636       else
637         grid_point_search_create(search.gps, xIsCyclic, src_grid.dims, src_grid.size, src_grid.cell_center_lon,
638                                  src_grid.cell_center_lat);
639 
640       if (src_grid.lextrapolate) grid_point_search_extrapolate(search.gps);
641     }
642   else if (useCellsearch)
643     {
644       searchMethodStr = "Cell search";
645       if (src_grid.type == RemapGridType::Reg2D)
646         grid_cell_search_create_reg_2d(search.gcs, src_grid.dims, src_grid.reg2d_corner_lon, src_grid.reg2d_corner_lat);
647       else
648         grid_cell_search_create(search.gcs, src_grid.size, src_grid.num_cell_corners, src_grid.cell_corner_lon,
649                                 src_grid.cell_corner_lat);
650     }
651   else if (!(src_grid.type == RemapGridType::Reg2D || tgt_grid.type == RemapGridType::Reg2D))
652     {
653       searchMethodStr = "Latitude bins";
654       search.srcBins.cell_bound_box.resize(4 * src_grid.size);
655       if (tgt_grid.luse_cell_corners) search.tgtBins.cell_bound_box.resize(4 * tgt_grid.size);
656 
657       cell_bounding_boxes(src_grid, &search.srcBins.cell_bound_box[0], REMAP_GRID_BASIS_SRC);
658       cell_bounding_boxes(tgt_grid, &search.tgtBins.cell_bound_box[0], REMAP_GRID_BASIS_TGT);
659       // Set up and assign address ranges to search bins in order to further restrict later searches
660       calc_lat_bins(search.srcBins);
661       if (mapType == RemapMethod::CONSERV_SCRIP || mapType == RemapMethod::CONSERV)
662         {
663           calc_lat_bins(search.tgtBins);
664           varray_free(search.tgtBins.bin_lats);
665           varray_free(search.srcBins.bin_lats);
666           if (mapType == RemapMethod::CONSERV) varray_free(search.tgtBins.cell_bound_box);
667         }
668     }
669 
670   if (Options::cdoVerbose && searchMethodStr) cdo_print("%s created: %.2f seconds", searchMethodStr, cdo_get_wtime() - start);
671 }
672 
673 void
remap_search_free(RemapSearch & search)674 remap_search_free(RemapSearch &search)
675 {
676   varray_free(search.srcBins.bin_addr);
677   varray_free(search.srcBins.bin_lats);
678   varray_free(search.srcBins.cell_bound_box);
679 
680   varray_free(search.tgtBins.bin_addr);
681   varray_free(search.tgtBins.bin_lats);
682   varray_free(search.tgtBins.cell_bound_box);
683 
684   grid_point_search_delete(search.gps);
685   grid_cell_search_delete(search.gcs);
686 }
687 
688 void
remap_init_grids(RemapMethod mapType,bool lextrapolate,int gridID1,RemapGrid & src_grid,int gridID2,RemapGrid & tgt_grid)689 remap_init_grids(RemapMethod mapType, bool lextrapolate, int gridID1, RemapGrid &src_grid, int gridID2, RemapGrid &tgt_grid)
690 {
691   auto reg2d_src_gridID = gridID1;
692   auto reg2d_tgt_gridID = gridID2;
693 
694   remap_grid_init(src_grid);
695   remap_grid_init(tgt_grid);
696 
697   if (mapType == RemapMethod::BILINEAR || mapType == RemapMethod::BICUBIC || mapType == RemapMethod::DISTWGT
698       || mapType == RemapMethod::CONSERV)
699     {
700       if (IS_REG2D_GRID(gridID1)) src_grid.type = RemapGridType::Reg2D;
701       // src_grid.type = 0;
702     }
703 
704   //#ifdef YAC_CELL_SEARCH
705   // if (IS_REG2D_GRID(gridID2) && mapType == RemapMethod::CONSERV) tgt_grid.type = RemapGridType::Reg2D;
706   //#else
707   if (src_grid.type == RemapGridType::Reg2D)
708     {
709       if (IS_REG2D_GRID(gridID2) && mapType == RemapMethod::CONSERV) tgt_grid.type = RemapGridType::Reg2D;
710       // else src_grid.type = -1;
711     }
712   //#endif
713 
714   if (!remap_gen_weights && IS_REG2D_GRID(gridID2) && tgt_grid.type != RemapGridType::Reg2D)
715     {
716       if (mapType == RemapMethod::DISTWGT) tgt_grid.type = RemapGridType::Reg2D;
717       if (mapType == RemapMethod::BILINEAR && src_grid.type == RemapGridType::Reg2D) tgt_grid.type = RemapGridType::Reg2D;
718     }
719 
720   src_grid.lextrapolate = lextrapolate;
721 
722   if (mapType == RemapMethod::CONSERV_SCRIP || mapType == RemapMethod::CONSERV)
723     {
724       if (src_grid.type != RemapGridType::Reg2D)
725         {
726           src_grid.luse_cell_corners = true;
727           src_grid.lneed_cell_corners = true;
728         }
729 
730       if (tgt_grid.type != RemapGridType::Reg2D)
731         {
732           tgt_grid.luse_cell_corners = true;
733           tgt_grid.lneed_cell_corners = true;
734         }
735     }
736 
737   src_grid.gridID = gridID1;
738   tgt_grid.gridID = gridID2;
739 
740   if (gridInqType(gridID1) == GRID_UNSTRUCTURED && !gridHasCoordinates(gridID1))
741     {
742       const auto reference = dereferenceGrid(gridID1);
743       if (reference.isValid) src_grid.gridID = gridID1 = reference.gridID;
744       if (reference.notFound) cdo_abort("Reference to source grid not found!");
745     }
746 
747   if (gridInqType(gridID2) == GRID_UNSTRUCTURED && !gridHasCoordinates(gridID2))
748     {
749       const auto reference = dereferenceGrid(gridID2);
750       if (reference.isValid) tgt_grid.gridID = gridID2 = reference.gridID;
751       if (reference.notFound) cdo_abort("Reference to target grid not found!");
752     }
753 
754   const auto sgridID = src_grid.gridID;
755   if (gridInqSize(sgridID) > 1 && gridProjIsSupported(sgridID))
756     {
757       bool lbounds = true;
758       src_grid.gridID = gridID1 = gridToCurvilinear(src_grid.gridID, lbounds);
759       src_grid.tmpgridID = src_grid.gridID;
760     }
761 
762   // if (src_grid.type != RemapGridType::Reg2D)
763   remapDefineGrid(mapType, gridID1, src_grid, "Source");
764   remapDefineGrid(mapType, gridID2, tgt_grid, "Target");
765 
766   if (src_grid.type == RemapGridType::Reg2D || tgt_grid.type == RemapGridType::Reg2D)
767     {
768       if (src_grid.type == RemapGridType::Reg2D) remap_define_reg2d(reg2d_src_gridID, src_grid);
769       if (tgt_grid.type == RemapGridType::Reg2D) remap_define_reg2d(reg2d_tgt_gridID, tgt_grid);
770     }
771 }
772 
773 /*****************************************************************************/
774 
775 template <typename T>
776 static void
remap_stat(int remapOrder,RemapGrid & src_grid,RemapGrid & tgt_grid,RemapVars & rv,const Varray<T> & array1,const Varray<T> & array2,T missval)777 remap_stat(int remapOrder, RemapGrid &src_grid, RemapGrid &tgt_grid, RemapVars &rv, const Varray<T> &array1,
778            const Varray<T> &array2, T missval)
779 {
780   if (remapOrder == 2)
781     cdo_print("Second order mapping from grid1 to grid2:");
782   else
783     cdo_print("First order mapping from grid1 to grid2:");
784   cdo_print("----------------------------------------------");
785 
786   auto mmm = varray_min_max_mean_mv(src_grid.size, array1, missval);
787   cdo_print("  Grid1 min,mean,max: %g %g %g", mmm.min, mmm.mean, mmm.max);
788 
789   mmm = varray_min_max_mean_mv(tgt_grid.size, array2, missval);
790   cdo_print("  Grid2 min,mean,max: %g %g %g", mmm.min, mmm.mean, mmm.max);
791 
792   // Conservation Test
793 
794   if (src_grid.cell_area.size())
795     {
796       cdo_print("  Conservation:");
797       double sum = 0;
798       for (size_t n = 0; n < src_grid.size; ++n)
799         if (!DBL_IS_EQUAL(array1[n], missval)) sum += array1[n] * src_grid.cell_area[n] * src_grid.cell_frac[n];
800       cdo_print("  Grid1 Integral = %g", sum);
801 
802       sum = 0;
803       for (size_t n = 0; n < tgt_grid.size; ++n)
804         if (!DBL_IS_EQUAL(array2[n], missval)) sum += array2[n] * tgt_grid.cell_area[n] * tgt_grid.cell_frac[n];
805       cdo_print("  Grid2 Integral = %g", sum);
806       /*
807       for ( n = 0; n < src_grid.size; n++ )
808        fprintf(stderr, "1 %d %g %g %g\n", n, array1[n], src_grid.cell_area[n],
809       src_grid.cell_frac[n]); for ( n = 0; n < tgt_grid.size; n++ )
810         fprintf(stderr, "2 %d %g %g %g\n", n, array2[n], tgt_grid.cell_area[n],
811       tgt_grid.cell_frac[n]);
812       */
813     }
814 
815   cdo_print("  Number of weights %zu", rv.num_wts);
816   cdo_print("  Number of sparse matrix entries %zu", rv.num_links);
817   cdo_print("  Total number of dest cells %zu", tgt_grid.size);
818 
819   std::vector<size_t> tgt_count(tgt_grid.size, 0);
820 
821 #ifdef HAVE_OPENMP4
822 #pragma omp simd
823 #endif
824   for (size_t n = 0; n < rv.num_links; ++n) tgt_count[rv.tgt_cell_add[n]]++;
825 
826   size_t imin = SIZE_MAX;
827   size_t imax = 0;
828   for (size_t n = 0; n < tgt_grid.size; ++n)
829     {
830       if (tgt_count[n] > 0)
831         {
832           if (tgt_count[n] < imin) imin = tgt_count[n];
833           if (tgt_count[n] > imax) imax = tgt_count[n];
834         }
835     }
836 
837   size_t idiff = (imax - imin) / 10 + 1;
838   size_t icount = 0;
839   for (size_t i = 0; i < tgt_grid.size; ++i)
840     if (tgt_count[i] > 0) icount++;
841 
842   cdo_print("  Number of cells participating in remap %zu", icount);
843 
844   if (icount)
845     {
846       cdo_print("  Min no of entries/row = %zu", imin);
847       cdo_print("  Max no of entries/row = %zu", imax);
848 
849       imax = imin + idiff;
850       for (size_t n = 0; n < 10; ++n)
851         {
852           icount = 0;
853           for (size_t i = 0; i < tgt_grid.size; ++i)
854             if (tgt_count[i] >= imin && tgt_count[i] < imax) icount++;
855 
856           if (icount) cdo_print("  Num of rows with entries between %zu - %zu  %zu", imin, imax - 1, icount);
857 
858           imin = imin + idiff;
859           imax = imax + idiff;
860         }
861     }
862 
863   if (rv.sort_add) cdo_print("  Sparse matrix entries are explicitly sorted.");
864 }
865 
866 void
remap_stat(int remapOrder,RemapGrid & src_grid,RemapGrid & tgt_grid,RemapVars & rv,const Field & field1,const Field & field2)867 remap_stat(int remapOrder, RemapGrid &src_grid, RemapGrid &tgt_grid, RemapVars &rv, const Field &field1, const Field &field2)
868 {
869   if (field1.memType == MemType::Float)
870     remap_stat(remapOrder, src_grid, tgt_grid, rv, field1.vec_f, field2.vec_f, (float) field1.missval);
871   else
872     remap_stat(remapOrder, src_grid, tgt_grid, rv, field1.vec_d, field2.vec_d, field1.missval);
873 }
874 
875 /*****************************************************************************/
876 
877 template <typename T>
878 void
remap_gradients(RemapGrid & grid,const Varray<short> & mask,const Varray<T> & array,RemapGradients & gradients)879 remap_gradients(RemapGrid &grid, const Varray<short> &mask, const Varray<T> &array, RemapGradients &gradients)
880 {
881   if (grid.rank != 2) cdo_abort("Internal problem (remap_gradients), grid rank = %d!", grid.rank);
882 
883   auto grid_size = grid.size;
884   auto nx = grid.dims[0];
885   auto ny = grid.dims[1];
886 
887 #ifdef _OPENMP
888 #pragma omp parallel for default(none) shared(grid_size, gradients, grid, nx, ny, array, mask)
889 #endif
890   for (size_t n = 0; n < grid_size; ++n)
891     {
892       if (mask[n])
893         {
894           // clang-format off
895           auto delew = 0.5;
896           auto delns = 0.5;
897 
898           const auto j = n / nx + 1;
899           const auto i = n - (j - 1) * nx + 1;
900 
901           auto ip1 = i + 1;
902           auto im1 = i - 1;
903           auto jp1 = j + 1;
904           auto jm1 = j - 1;
905 
906           if (ip1 > nx) ip1 = ip1 - nx;
907           if (im1 < 1)  im1 = nx;
908           if (jp1 > ny) { jp1 = j; delns = 1.0; }
909           if (jm1 < 1)  { jm1 = j; delns = 1.0; }
910 
911           auto in = (jp1 - 1) * nx + i - 1;
912           auto is = (jm1 - 1) * nx + i - 1;
913           auto ie = (j - 1) * nx + ip1 - 1;
914           auto iw = (j - 1) * nx + im1 - 1;
915 
916           auto ine = (jp1 - 1) * nx + ip1 - 1;
917           auto inw = (jp1 - 1) * nx + im1 - 1;
918           auto ise = (jm1 - 1) * nx + ip1 - 1;
919           auto isw = (jm1 - 1) * nx + im1 - 1;
920 
921           // Compute i-gradient
922           if (!mask[ie]) { ie = n; delew = 1.0; }
923           if (!mask[iw]) { iw = n; delew = 1.0; }
924 
925           gradients.grad_lat[n] = delew * (array[ie] - array[iw]);
926 
927           // Compute j-gradient
928           if (!mask[in]) { in = n; delns = 1.0; }
929           if (!mask[is]) { is = n; delns = 1.0; }
930 
931           gradients.grad_lon[n] = delns * (array[in] - array[is]);
932           // clang-format on
933 
934           // Compute ij-gradient
935           delew = 0.5;
936           delns = (jp1 == j || jm1 == j) ? 1.0 : 0.5;
937 
938           if (!mask[ine])
939             {
940               if (in != n)
941                 {
942                   ine = in;
943                   delew = 1.0;
944                 }
945               else if (ie != n)
946                 {
947                   ine = ie;
948                   inw = iw;
949                   if (inw == n) delew = 1.0;
950                   delns = 1.0;
951                 }
952               else
953                 {
954                   ine = n;
955                   inw = iw;
956                   delew = 1.0;
957                   delns = 1.0;
958                 }
959             }
960 
961           if (!mask[inw])
962             {
963               if (in != n)
964                 {
965                   inw = in;
966                   delew = 1.0;
967                 }
968               else if (iw != n)
969                 {
970                   inw = iw;
971                   ine = ie;
972                   if (ie == n) delew = 1.0;
973                   delns = 1.0;
974                 }
975               else
976                 {
977                   inw = n;
978                   ine = ie;
979                   delew = 1.0;
980                   delns = 1.0;
981                 }
982             }
983 
984           const auto grad_lat_zero = delew * (array[ine] - array[inw]);
985 
986           if (!mask[ise])
987             {
988               if (is != n)
989                 {
990                   ise = is;
991                   delew = 1.0;
992                 }
993               else if (ie != n)
994                 {
995                   ise = ie;
996                   isw = iw;
997                   if (isw == n) delew = 1.0;
998                   delns = 1.0;
999                 }
1000               else
1001                 {
1002                   ise = n;
1003                   isw = iw;
1004                   delew = 1.0;
1005                   delns = 1.0;
1006                 }
1007             }
1008 
1009           if (!mask[isw])
1010             {
1011               if (is != n)
1012                 {
1013                   isw = is;
1014                   delew = 1.0;
1015                 }
1016               else if (iw != n)
1017                 {
1018                   isw = iw;
1019                   ise = ie;
1020                   if (ie == n) delew = 1.0;
1021                   delns = 1.0;
1022                 }
1023               else
1024                 {
1025                   isw = n;
1026                   ise = ie;
1027                   delew = 1.0;
1028                   delns = 1.0;
1029                 }
1030             }
1031 
1032           const auto grad_lon_zero = delew * (array[ise] - array[isw]);
1033           gradients.grad_latlon[n] = delns * (grad_lat_zero - grad_lon_zero);
1034         }
1035       else
1036         {
1037           gradients.grad_lat[n] = 0.0;
1038           gradients.grad_lon[n] = 0.0;
1039           gradients.grad_latlon[n] = 0.0;
1040         }
1041     }
1042 } /* remap_gradients */
1043 
1044 // Explicit instantiation
1045 template void remap_gradients(RemapGrid &grid, const Varray<short> &mask, const Varray<float> &array, RemapGradients &gradients);
1046 template void remap_gradients(RemapGrid &grid, const Varray<short> &mask, const Varray<double> &array, RemapGradients &gradients);
1047 
1048 void
remap_gradients(RemapGrid & grid,const Field & field,RemapGradients & gradients)1049 remap_gradients(RemapGrid &grid, const Field &field, RemapGradients &gradients)
1050 {
1051   auto grid_size = grid.size;
1052   Varray<short> mask(grid_size);
1053 #ifdef _OPENMP
1054 #pragma omp parallel for default(none) schedule(static) shared(grid_size, grid, mask)
1055 #endif
1056   for (size_t i = 0; i < grid_size; ++i) mask[i] = (grid.mask[i] > 0);
1057 
1058   if (field.memType == MemType::Float)
1059     remap_gradients(grid, mask, field.vec_f, gradients);
1060   else
1061     remap_gradients(grid, mask, field.vec_d, gradients);
1062 }
1063 
1064 /*****************************************************************************/
1065 
1066 void
remap_check_area(size_t grid_size,double * cell_area,const char * name)1067 remap_check_area(size_t grid_size, double *cell_area, const char *name)
1068 {
1069   for (size_t n = 0; n < grid_size; ++n)
1070     {
1071       if (cell_area[n] < -.01) cdo_print("%s grid area error: %zu %g", name, n, cell_area[n]);
1072     }
1073 }
1074