YAC 3.13.2
Yet Another Coupler
Loading...
Searching...
No Matches
read_scrip_grid.c
Go to the documentation of this file.
1// Copyright (c) 2024 The YAC Authors
2//
3// SPDX-License-Identifier: BSD-3-Clause
4
5#ifdef HAVE_CONFIG_H
6#include "config.h"
7#endif
8
9#include <stdlib.h>
10#include <stdio.h>
11#include <stdint.h>
12#include <string.h>
13#include <limits.h>
14
15#include <mpi.h>
16
17#include "read_scrip_grid.h"
18#include "geometry.h"
19#include "utils_common.h"
20#include "io_utils.h"
21#include "yac_mpi_internal.h"
22
23#ifdef YAC_NETCDF_ENABLED
24
25#include <netcdf.h>
26
27struct point_with_index {
28
29 int32_t lon, lat;
30 int rank;
32 size_t i;
33 double dlon, dlat;
34};
35
41
48
50 double ** vertex_lon, double ** vertex_lat,
51 size_t * nbr_vertices, int * old_to_new_id);
53 double ** vertex_lon, double ** vertex_lat,
54 size_t * nbr_vertices, int * old_to_new_id,
55 yac_int ** vertex_ids, MPI_Comm comm);
56static void detect_duplicated_cells(
57 int * cell_to_vertex, int * num_vertices_per_cell, int * cell_mask,
58 size_t nbr_cells, size_t ** duplicated_cell_idx,
59 size_t ** orig_cell_idx, size_t * nbr_duplicated_cells);
61 int * cell_to_vertex, int * num_vertices_per_cell, int * cell_mask,
62 yac_int * cell_ids, size_t num_cells, yac_int * vertex_ids,
63 MPI_Comm comm, size_t ** duplicated_cell_idx,
64 yac_int ** orig_cell_ids, size_t * num_duplicated_cells_);
65
66static inline int compare_int(
67 const void * a,const void * b) {
68 return *(int const *)a - *(int const*)b;
69}
70
71static inline int compare_yac_int(
72 const void * a, const void * b) {
73 return (*(yac_int const *)a > *(yac_int const*)b) -
74 (*(yac_int const *)a < *(yac_int const*)b);
75}
76
78 const void * a_, const void * b_) {
79
80 struct cell_vertex_ids_with_index const * a =
81 (struct cell_vertex_ids_with_index const *)a_;
82 struct cell_vertex_ids_with_index const * b =
83 (struct cell_vertex_ids_with_index const *)b_;
84
85 size_t num_vertices_a = a->num_vertices;
86 size_t num_vertices_b = b->num_vertices;
87 int ret = (num_vertices_a > num_vertices_b) -
88 (num_vertices_a < num_vertices_b);
89
90 for (size_t i = 0; (i < num_vertices_a) && !ret; ++i)
91 ret = (a->vertex_ids[i] > b->vertex_ids[i]) -
92 (a->vertex_ids[i] < b->vertex_ids[i]);
93
94 return ret;
95}
96
98 const void * a_, const void * b_) {
99
101 if (ret) return ret;
102
103 struct cell_vertex_ids_with_index const * a =
104 (struct cell_vertex_ids_with_index const *)a_;
105 struct cell_vertex_ids_with_index const * b =
106 (struct cell_vertex_ids_with_index const *)b_;
107
108 return (a->cell_id > b->cell_id) - (a->cell_id < b->cell_id);
109}
110
112 const void * a_, const void * b_) {
113
114 struct point_with_index * a = (struct point_with_index *)a_;
115 struct point_with_index * b = (struct point_with_index *)b_;
116
117 int ret = (a->lon > b->lon) - (a->lon < b->lon);
118 if (ret) return ret;
119 ret = (a->lat > b->lat) - (a->lat < b->lat);
120 if (ret) return ret;
121 return (a->id > b->id) - (a->id < b->id);
122}
123
125 const void * a_, const void * b_) {
126
127 struct point_with_index * a = (struct point_with_index *)a_;
128 struct point_with_index * b = (struct point_with_index *)b_;
129
130 int ret = (a->lon > b->lon) - (a->lon < b->lon);
131 if (ret) return ret;
132 return (a->lat > b->lat) - (a->lat < b->lat);
133}
134
136 const void * a_, const void * b_) {
137
138 struct point_with_index * a = (struct point_with_index *)a_;
139 struct point_with_index * b = (struct point_with_index *)b_;
140
141 return (a->rank > b->rank) - (a->rank < b->rank);
142}
143
145 const void * a_, const void * b_) {
146
147 struct point_with_index * a = (struct point_with_index *)a_;
148 struct point_with_index * b = (struct point_with_index *)b_;
149
150 return (a->i > b->i) - (a->i < b->i);
151}
152
154 const char * filename, const char * grid_name,
155 size_t cell_mask_size, int * cell_mask,
156 size_t * num_vertices_, size_t * num_cells_, int ** num_vertices_per_cell_,
157 int ** cell_to_vertex_, double ** x_vertices, double ** y_vertices,
158 double ** x_cells, double ** y_cells,
159 size_t ** duplicated_cell_idx_, size_t ** orig_cell_idx_,
160 size_t * nbr_duplicated_cells_, int with_vertices) {
161
162 size_t grid_name_len = strlen(grid_name) + 1;
163 char cla_var_name[4 + grid_name_len];
164 char clo_var_name[4 + grid_name_len];
165 char lat_var_name[4 + grid_name_len];
166 char lon_var_name[4 + grid_name_len];
167
168 snprintf(cla_var_name, 4 + grid_name_len, "%s.cla", grid_name);
169 snprintf(clo_var_name, 4 + grid_name_len, "%s.clo", grid_name);
170 snprintf(lat_var_name, 4 + grid_name_len, "%s.lat", grid_name);
171 snprintf(lon_var_name, 4 + grid_name_len, "%s.lon", grid_name);
172
173 int ncid;
174 yac_nc_open(filename, NC_NOWRITE, &ncid);
175
176 // get variable ids
177 int cla_var_id;
178 int clo_var_id;
179 int lat_var_id;
180 int lon_var_id;
181 if (with_vertices) {
182 yac_nc_inq_varid(ncid, cla_var_name, &cla_var_id);
183 yac_nc_inq_varid(ncid, clo_var_name, &clo_var_id);
184 }
185 yac_nc_inq_varid(ncid, lat_var_name, &lat_var_id);
186 yac_nc_inq_varid(ncid, lon_var_name, &lon_var_id);
187
188 // get dimension ids
189 int dim_ids[NC_MAX_VAR_DIMS];
190 if (with_vertices) {
191 int ndims;
193 nc_inq_var(ncid, cla_var_id, NULL, NULL, &ndims, dim_ids, NULL));
195 ndims == 3,
196 "ERROR(yac_read_scrip_basic_grid_information): "
197 "invalid number of dimensions for variable \"%s\" (has to be 3 not %d)",
198 cla_var_name, ndims);
199 } else {
200 int ndims;
202 nc_inq_var(ncid, lat_var_id, NULL, NULL, &ndims, dim_ids, NULL));
204 ndims == 2,
205 "ERROR(yac_read_scrip_basic_grid_information): "
206 "invalid number of dimensions for variable \"%s\" (has to be 2 not %d)",
207 lat_var_name, ndims);
208 dim_ids[2] = dim_ids[1];
209 dim_ids[1] = dim_ids[0];
210 dim_ids[0] = -1;
211 }
212 int crn_dim_id = dim_ids[0];
213 int x_dim_id = dim_ids[2];
214 int y_dim_id = dim_ids[1];
215
216 // get dimension length
217 size_t crn_dim_len, x_dim_len, y_dim_len;
218 if (with_vertices)
219 YAC_HANDLE_ERROR(nc_inq_dimlen(ncid, crn_dim_id, &crn_dim_len));
220 else
221 crn_dim_len = 1;
222 YAC_HANDLE_ERROR(nc_inq_dimlen(ncid, x_dim_id, &x_dim_len));
223 YAC_HANDLE_ERROR(nc_inq_dimlen(ncid, y_dim_id, &y_dim_len));
224
225 size_t num_cells = x_dim_len * y_dim_len;
226 size_t num_vertices = num_cells * crn_dim_len;
227
229 num_cells == cell_mask_size,
230 "ERROR(yac_read_scrip_basic_grid_information): "
231 "cell mask size is inconsistent with number of grid cells")
232
233 if ((x_cells != NULL) && (y_cells != NULL)) {
234 double * lat = xmalloc(num_cells * sizeof(*lat));
235 double * lon = xmalloc(num_cells * sizeof(*lon));
236 YAC_HANDLE_ERROR(nc_get_var_double(ncid, lat_var_id, lat));
237 YAC_HANDLE_ERROR(nc_get_var_double(ncid, lon_var_id, lon));
238 *x_cells = lon;
239 *y_cells = lat;
240 }
241
242 // allocate variables
243 double * cla = xmalloc(num_vertices * sizeof(*cla));
244 double * clo = xmalloc(num_vertices * sizeof(*clo));
245
246 //read variables
247 if (with_vertices) {
248
249 YAC_HANDLE_ERROR(nc_get_var_double(ncid, cla_var_id, cla));
250 YAC_HANDLE_ERROR(nc_get_var_double(ncid, clo_var_id, clo));
251
252 } else {
253
255 (x_cells != NULL) && (y_cells != NULL),
256 "ERROR(yac_read_scrip_basic_grid_information): "
257 "no pointers for cell coordinates were provided")
258
259 memcpy(cla, *y_cells, num_cells * sizeof(*cla));
260 memcpy(clo, *x_cells, num_cells * sizeof(*clo));
261 }
262
263 YAC_HANDLE_ERROR(nc_close(ncid));
264
265 size_t * reorder_idx = xmalloc(num_vertices * sizeof(*reorder_idx));
266 for (size_t y = 0, l = 0; y < y_dim_len; ++y)
267 for (size_t x = 0; x < x_dim_len; ++x)
268 for (size_t n = 0; n < crn_dim_len; ++n, ++l)
269 reorder_idx[x + y * x_dim_len + n * x_dim_len * y_dim_len] = l;
270
271 // remove duplicated vertices
272 int * cell_to_vertex = xmalloc(num_vertices * sizeof(*cell_to_vertex));
273 remove_duplicated_vertices(&clo, &cla, &num_vertices, cell_to_vertex);
274
275 // we have to reorder cell_to_vertex
277 reorder_idx, num_cells * crn_dim_len, cell_to_vertex);
278 free(reorder_idx);
279
280 // determine number of vertices per cell and compact cell_to_vertex
281 int * num_vertices_per_cell =
282 xmalloc(num_cells * sizeof(*num_vertices_per_cell));
283 size_t total_num_cell_vertices = 0;
284 int * to_vertices = cell_to_vertex;
285 int * from_vertices = cell_to_vertex;
286 if (crn_dim_len > 1) {
287 for (size_t i = 0; i < num_cells; ++i) {
288 size_t curr_num_vertices = 0;
289 if (cell_mask[i]) {
290 if (crn_dim_len > 1) {
291 int prev_vertex = from_vertices[crn_dim_len-1];
292 for (size_t j = 0; j < crn_dim_len; ++j, ++from_vertices) {
293 int curr_vertex = *from_vertices;
294 if (prev_vertex != curr_vertex) {
295 prev_vertex = curr_vertex;
296 if (to_vertices != from_vertices) *to_vertices = curr_vertex;
297 ++curr_num_vertices;
298 ++to_vertices;
299 }
300 }
301 }
302 } else {
303 from_vertices += crn_dim_len;
304 }
305 num_vertices_per_cell[i] = (int)curr_num_vertices;
306 total_num_cell_vertices += curr_num_vertices;
307 }
308 } else {
309 for (size_t i = 0; i < num_cells; ++i) {
310 if (cell_mask[i]) {
311 if (to_vertices != from_vertices) *to_vertices = *from_vertices;
312 ++to_vertices;
313 ++total_num_cell_vertices;
314 num_vertices_per_cell[i] = 1;
315 } else {
316 num_vertices_per_cell[i] = 0;
317 }
318 ++from_vertices;
319 }
320 }
321
322 if (total_num_cell_vertices != num_cells * crn_dim_len)
324 xrealloc(
325 cell_to_vertex, total_num_cell_vertices * sizeof(*cell_to_vertex));
326
327 // make a copy of cell_to_vertex and sort the vertex indices of each cell,
328 // this is required by the detection of duplicated cells
329 int * cell_to_vertex_copy =
330 xmalloc(total_num_cell_vertices * sizeof(*cell_to_vertex_copy));
331 memcpy(
332 cell_to_vertex_copy, cell_to_vertex,
333 total_num_cell_vertices * sizeof(*cell_to_vertex_copy));
334 for (size_t i = 0, offset = 0; i < num_cells; ++i) {
335 qsort(
336 cell_to_vertex_copy + offset, (size_t)(num_vertices_per_cell[i]),
337 sizeof(*cell_to_vertex_copy), compare_int);
338 offset += num_vertices_per_cell[i];
339 }
340
341 // detect duplicated cells
342 size_t * duplicated_cell_idx;
343 size_t * orig_cell_idx;
344 size_t nbr_duplicated_cells;
346 cell_to_vertex_copy, num_vertices_per_cell, cell_mask, num_cells,
347 &duplicated_cell_idx, &orig_cell_idx, &nbr_duplicated_cells);
348 free(cell_to_vertex_copy);
349
350 // mask out duplicated cells
351 for (size_t i = 0; i < nbr_duplicated_cells; ++i)
352 cell_mask[duplicated_cell_idx[i]] = 0;
353
354 if ((duplicated_cell_idx_ != NULL) && (orig_cell_idx_ != NULL) &&
355 (nbr_duplicated_cells_ != NULL)) {
356 *duplicated_cell_idx_ = duplicated_cell_idx;
357 *orig_cell_idx_ = orig_cell_idx;
358 *nbr_duplicated_cells_ = nbr_duplicated_cells;
359 } else {
360 free(duplicated_cell_idx);
361 free(orig_cell_idx);
362 }
363
364 if (with_vertices) {
365
366 *num_vertices_ = num_vertices;
367 *num_vertices_per_cell_ = num_vertices_per_cell;
368 *cell_to_vertex_ = cell_to_vertex;
369 *x_vertices = clo;
370 *y_vertices = cla;
371
372 } else {
373
374 free(num_vertices_per_cell);
375 free(cell_to_vertex);
376 free(clo);
377 free(cla);
378
379 }
380 *num_cells_ = num_cells;
381}
382
384 const char * filename, const char * grid_name,
385 size_t cell_mask_size, int * cell_mask,
386 size_t * num_vertices_, size_t * num_cells_, int ** num_vertices_per_cell_,
387 int ** cell_to_vertex_,
388 double ** x_vertices_, double ** y_vertices_, yac_int ** vertex_ids_,
389 double ** x_cells_, double ** y_cells_, yac_int ** cell_ids_,
390 size_t ** duplicated_cell_idx_, yac_int ** orig_cell_ids_,
391 size_t * nbr_duplicated_cells_,
392 int io_rank_idx, int num_io_ranks, MPI_Comm comm, int with_vertices) {
393
394 size_t num_cells = 0;
395 size_t num_vertices = 0;
396 size_t * reorder_idx = NULL;
397 size_t max_num_vertices_per_cell = 0;
398
399 double * x_vertices = NULL;
400 double * y_vertices = NULL;
401 yac_int * vertex_ids = NULL;
402 double * x_cells = NULL;
403 double * y_cells = NULL;
404 yac_int * cell_ids = NULL;
405
406 if ((io_rank_idx >= 0) && (io_rank_idx < num_io_ranks)) {
407
408 size_t grid_name_len = strlen(grid_name) + 1;
409 char cla_var_name[4 + grid_name_len];
410 char clo_var_name[4 + grid_name_len];
411 char lat_var_name[4 + grid_name_len];
412 char lon_var_name[4 + grid_name_len];
413
414 snprintf(cla_var_name, 4 + grid_name_len, "%s.cla", grid_name);
415 snprintf(clo_var_name, 4 + grid_name_len, "%s.clo", grid_name);
416 snprintf(lat_var_name, 4 + grid_name_len, "%s.lat", grid_name);
417 snprintf(lon_var_name, 4 + grid_name_len, "%s.lon", grid_name);
418
419 int ncid;
420 yac_nc_open(filename, NC_NOWRITE, &ncid);
421
422 // get variable ids
423 int cla_var_id, clo_var_id, lat_var_id, lon_var_id;
424 if (with_vertices) {
425 yac_nc_inq_varid(ncid, cla_var_name, &cla_var_id);
426 yac_nc_inq_varid(ncid, clo_var_name, &clo_var_id);
427 }
428 yac_nc_inq_varid(ncid, lat_var_name, &lat_var_id);
429 yac_nc_inq_varid(ncid, lon_var_name, &lon_var_id);
430
431 // get dimension ids
432 int dim_ids[NC_MAX_VAR_DIMS];
433 if (with_vertices) {
434 int ndims;
436 nc_inq_var(ncid, cla_var_id, NULL, NULL, &ndims, dim_ids, NULL));
438 ndims == 3,
439 "ERROR(yac_read_part_scrip_basic_grid_information): "
440 "invalid number of dimensions for variable \"%s\" (has to be 3 not %d)",
441 cla_var_name, ndims);
442 } else {
443 int ndims;
445 nc_inq_var(ncid, lat_var_id, NULL, NULL, &ndims, dim_ids, NULL));
447 ndims == 2,
448 "ERROR(yac_read_part_scrip_basic_grid_information): "
449 "invalid number of dimensions for variable \"%s\" (has to be 2 not %d)",
450 lat_var_name, ndims);
451 dim_ids[2] = dim_ids[1];
452 dim_ids[1] = dim_ids[0];
453 dim_ids[0] = -1;
454 }
455 int crn_dim_id = dim_ids[0];
456 int x_dim_id = dim_ids[2];
457 int y_dim_id = dim_ids[1];
458
459 // get dimension length
460 size_t crn_dim_len, x_dim_len, y_dim_len;
461 if (with_vertices)
462 YAC_HANDLE_ERROR(nc_inq_dimlen(ncid, crn_dim_id, &crn_dim_len));
463 else
464 crn_dim_len = 1;
465 YAC_HANDLE_ERROR(nc_inq_dimlen(ncid, x_dim_id, &x_dim_len));
466 YAC_HANDLE_ERROR(nc_inq_dimlen(ncid, y_dim_id, &y_dim_len));
467
468 // decompose x dimension among io ranks
469 size_t x_start_local =
470 (size_t)
471 (((unsigned long)x_dim_len * (unsigned long)io_rank_idx) /
472 (unsigned long)num_io_ranks);
473 size_t x_count_local =
474 (size_t)
475 (((unsigned long)x_dim_len * (io_rank_idx + 1)) /
476 (unsigned long)num_io_ranks) - x_start_local;
477
478 num_cells = x_count_local * y_dim_len;
479
480 if (with_vertices) {
481
482 num_vertices = num_cells * crn_dim_len;
483 max_num_vertices_per_cell = crn_dim_len;
484
485 } else {
486
487 num_vertices = num_cells;
488 max_num_vertices_per_cell = 1;
489 }
490
492 num_cells == cell_mask_size,
493 "ERROR(yac_read_part_scrip_basic_grid_information): "
494 "cell mask size is inconsistent with number of grid cells")
495
496 if ((x_cells_ != NULL) && (y_cells_ != NULL)) {
497 x_cells = xmalloc(num_cells * sizeof(*x_cells));
498 y_cells = xmalloc(num_cells * sizeof(*y_cells));
499 size_t start[2] = {0, x_start_local};
500 size_t count[2] = {y_dim_len, x_count_local};
502 nc_get_vara_double(ncid, lon_var_id, start, count, x_cells));
504 nc_get_vara_double(ncid, lat_var_id, start, count, y_cells));
505 }
506
507 if (with_vertices) {
508
509 // allocate variables
510 x_vertices = xmalloc(num_vertices * sizeof(*x_vertices));
511 y_vertices = xmalloc(num_vertices * sizeof(*y_vertices));
512
513 //read variables
514 size_t start[3] = {0, 0, x_start_local};
515 size_t count[3] = {crn_dim_len, y_dim_len, x_count_local};
517 nc_get_vara_double(ncid, clo_var_id, start, count, x_vertices));
519 nc_get_vara_double(ncid, cla_var_id, start, count, y_vertices));
520
521 } else {
522
523 x_vertices = xmalloc(num_cells * sizeof(*x_vertices));
524 y_vertices = xmalloc(num_cells * sizeof(*y_vertices));
525
527 (x_cells_ != NULL) && (y_cells_ != NULL),
528 "ERROR(yac_read_part_scrip_basic_grid_information): "
529 "no pointers for cell coordinates were provided")
530
531 memcpy(x_vertices, x_cells, num_cells * sizeof(*x_vertices));
532 memcpy(y_vertices, y_cells, num_cells * sizeof(*y_vertices));
533 }
534
535 YAC_HANDLE_ERROR(nc_close(ncid));
536
537 reorder_idx = xmalloc(num_vertices * sizeof(*reorder_idx));
538 for (size_t y = 0, l = 0; y < y_dim_len; ++y)
539 for (size_t x = 0; x < x_count_local; ++x)
540 for (size_t n = 0; n < crn_dim_len; ++n, ++l)
541 reorder_idx[
542 x + y * x_count_local + n * x_count_local * y_dim_len] = l;
543
544 // generate global cell ids
545 cell_ids = xmalloc(num_cells * sizeof(*cell_ids));
546 for (size_t y = 0, k = 0; y < y_dim_len; ++y)
547 for (size_t x = 0; x < x_count_local; ++x, ++k)
548 cell_ids[k] = y * x_dim_len + x + x_start_local;
549
550
552 (y_dim_len * x_dim_len * crn_dim_len) <= YAC_INT_MAX,
553 "ERROR(yac_read_part_scrip_basic_grid_information): "
554 "number of vertices exceed maximum global id (%zu)",
555 (size_t)YAC_INT_MAX);
556
557 // generate global vertex ids
558 vertex_ids = xmalloc(num_cells * crn_dim_len * sizeof(*vertex_ids));
559 for (size_t n = 0, l = 0; n < crn_dim_len; ++n)
560 for (size_t i = 0; i < num_cells; ++i, ++l)
561 vertex_ids[l] = cell_ids[i] * crn_dim_len + n;
562 }
563
564 // remove duplicated vertices
565 int * cell_to_vertex = xmalloc(num_vertices * sizeof(*cell_to_vertex));
567 &x_vertices, &y_vertices, &num_vertices,
568 cell_to_vertex, &vertex_ids, comm);
569
570 // we have to reorder cell_to_vertex
572 reorder_idx, num_cells * max_num_vertices_per_cell, cell_to_vertex);
573 free(reorder_idx);
574
575 // determine number of vertices per cell and compact cell_to_vertex
576 int * num_vertices_per_cell =
577 xmalloc(num_cells * sizeof(*num_vertices_per_cell));
578 size_t total_num_cell_vertices = 0;
579 int * to_vertices = cell_to_vertex;
580 int * from_vertices = cell_to_vertex;
581 if (max_num_vertices_per_cell > 1) {
582 for (size_t i = 0; i < num_cells; ++i) {
583 size_t curr_num_vertices = 0;
584 if (cell_mask[i]) {
585 if (max_num_vertices_per_cell > 1) {
586 int prev_vertex = from_vertices[max_num_vertices_per_cell-1];
587 for (size_t j = 0; j < max_num_vertices_per_cell; ++j, ++from_vertices) {
588 int curr_vertex = *from_vertices;
589 if (prev_vertex != curr_vertex) {
590 prev_vertex = curr_vertex;
591 if (to_vertices != from_vertices) *to_vertices = curr_vertex;
592 ++curr_num_vertices;
593 ++to_vertices;
594 }
595 }
596 }
597 } else {
598 from_vertices += max_num_vertices_per_cell;
599 }
600 num_vertices_per_cell[i] = (int)curr_num_vertices;
601 total_num_cell_vertices += curr_num_vertices;
602 }
603 } else {
604 for (size_t i = 0; i < num_cells; ++i) {
605 if (cell_mask[i]) {
606 if (to_vertices != from_vertices) *to_vertices = *from_vertices;
607 ++to_vertices;
608 ++total_num_cell_vertices;
609 num_vertices_per_cell[i] = 1;
610 } else {
611 num_vertices_per_cell[i] = 0;
612 }
613 ++from_vertices;
614 }
615 }
616
617 if (total_num_cell_vertices != num_cells * max_num_vertices_per_cell)
619 xrealloc(
620 cell_to_vertex, total_num_cell_vertices * sizeof(*cell_to_vertex));
621
622 // detect duplicated cells
623 size_t * duplicated_cell_idx;
624 yac_int * orig_cell_ids;
625 size_t nbr_duplicated_cells;
627 cell_to_vertex, num_vertices_per_cell, cell_mask, cell_ids, num_cells,
628 vertex_ids, comm, &duplicated_cell_idx, &orig_cell_ids,
629 &nbr_duplicated_cells);
630
631 // mask out duplicated cells
632 for (size_t i = 0; i < nbr_duplicated_cells; ++i)
633 cell_mask[duplicated_cell_idx[i]] = 0;
634
635 if ((duplicated_cell_idx_ != NULL) && (orig_cell_ids_ != NULL) &&
636 (nbr_duplicated_cells_ != NULL)) {
637 *duplicated_cell_idx_ = duplicated_cell_idx;
638 *orig_cell_ids_ = orig_cell_ids;
639 *nbr_duplicated_cells_ = nbr_duplicated_cells;
640 } else {
641 free(duplicated_cell_idx);
642 free(orig_cell_ids);
643 }
644
645 if (with_vertices) {
646
647 *num_vertices_ = num_vertices;
648 *num_vertices_per_cell_ = num_vertices_per_cell;
649 *x_vertices_ = x_vertices;
650 *y_vertices_ = y_vertices;
651 *cell_to_vertex_ = cell_to_vertex;
652 *vertex_ids_ = vertex_ids;
653
654 } else {
655
656 free(num_vertices_per_cell);
657 free(x_vertices);
658 free(y_vertices);
659 free(cell_to_vertex);
660 free(vertex_ids);
661 }
662 *num_cells_ = num_cells;
663 if ((x_cells_ != NULL) && (y_cells_ != NULL)) {
664 *x_cells_ = x_cells;
665 *y_cells_ = y_cells;
666 }
667 *cell_ids_ = cell_ids;
668}
669
671 const char * filename, const char * grid_name, size_t * num_cells_,
672 int ** cell_mask) {
673
674 size_t grid_name_len = strlen(grid_name) + 1;
675 char msk_var_name[4 + grid_name_len];
676
677 snprintf(msk_var_name, 4 + grid_name_len, "%s.msk", grid_name);
678
679 int ncid;
680 yac_nc_open(filename, NC_NOWRITE, &ncid);
681
682 // get variable id
683 int msk_var_id;
684 yac_nc_inq_varid(ncid, msk_var_name, &msk_var_id);
685
686 // get dimension ids
687 int ndims;
688 int dim_ids[NC_MAX_VAR_DIMS];
690 nc_inq_var(ncid, msk_var_id, NULL, NULL, &ndims, dim_ids, NULL));
692 ndims == 2,
693 "ERROR(yac_read_scrip_mask_information): "
694 "invalid number of dimensions for variable \"%s\" (has to be 2 not %d)",
695 msk_var_name, ndims);
696 int x_dim_id = dim_ids[1];
697 int y_dim_id = dim_ids[0];
698
699 // get dimension length
700 size_t x_dim_len;
701 size_t y_dim_len;
702 YAC_HANDLE_ERROR(nc_inq_dimlen(ncid, x_dim_id, &x_dim_len));
703 YAC_HANDLE_ERROR(nc_inq_dimlen(ncid, y_dim_id, &y_dim_len));
704
705 size_t num_cells = x_dim_len * y_dim_len;
706
707 // allocate variable
708 *cell_mask = xmalloc(num_cells * sizeof(**cell_mask));
709
710 //read variable
711 YAC_HANDLE_ERROR(nc_get_var_int(ncid, msk_var_id, *cell_mask));
712
713 YAC_HANDLE_ERROR(nc_close(ncid));
714
715 *num_cells_ = num_cells;
716}
717
719 const char * filename, const char * grid_name, size_t * num_cells_,
720 int ** cell_mask, int io_rank_idx, int num_io_ranks) {
721
722 if ((io_rank_idx >= 0) && (io_rank_idx < num_io_ranks)) {
723
724 size_t grid_name_len = strlen(grid_name) + 1;
725 char msk_var_name[4 + grid_name_len];
726
727 snprintf(msk_var_name, 4 + grid_name_len, "%s.msk", grid_name);
728
729 int ncid;
730 yac_nc_open(filename, NC_NOWRITE, &ncid);
731
732 // get variable id
733 int msk_var_id;
734 yac_nc_inq_varid(ncid, msk_var_name, &msk_var_id);
735
736 // get dimension ids
737 int ndims;
738 int dim_ids[NC_MAX_VAR_DIMS];
740 nc_inq_var(ncid, msk_var_id, NULL, NULL, &ndims, dim_ids, NULL));
742 ndims == 2,
743 "ERROR(yac_read_part_scrip_mask_information): "
744 "invalid number of dimensions for variable \"%s\" (has to be 2 not %d)",
745 msk_var_name, ndims);
746 int x_dim_id = dim_ids[1];
747 int y_dim_id = dim_ids[0];
748
749 // get dimension length
750 size_t x_dim_len;
751 size_t y_dim_len;
752 YAC_HANDLE_ERROR(nc_inq_dimlen(ncid, x_dim_id, &x_dim_len));
753 YAC_HANDLE_ERROR(nc_inq_dimlen(ncid, y_dim_id, &y_dim_len));
754
755 // decompose x dimension among io ranks
756 size_t x_start_local =
757 (size_t)
758 (((unsigned long)x_dim_len * (unsigned long)io_rank_idx) /
759 (unsigned long)num_io_ranks);
760 size_t x_count_local =
761 (size_t)
762 (((unsigned long)x_dim_len * (io_rank_idx+1)) /
763 (unsigned long)num_io_ranks) - x_start_local;
764
765 size_t num_cells = x_count_local * y_dim_len;
766
767 // allocate variable
768 *cell_mask = xmalloc(num_cells * sizeof(**cell_mask));
769
770 //read variable
771 size_t start[2] = {0, x_start_local};
772 size_t count[2] = {y_dim_len, x_count_local};
774 nc_get_vara_int(ncid, msk_var_id, start, count, *cell_mask));
775
776 YAC_HANDLE_ERROR(nc_close(ncid));
777
778 *num_cells_ = num_cells;
779
780 } else {
781 *cell_mask = NULL;
782 *num_cells_ = 0;
783 }
784}
785
787 char const * grid_filename, char const * mask_filename,
788 char const * grid_name, int valid_mask_value,
789 size_t * num_vertices, size_t * num_cells, int ** num_vertices_per_cell,
790 double ** x_vertices, double ** y_vertices,
791 double ** x_cells, double ** y_cells,
792 int ** cell_to_vertex, int ** cell_core_mask, size_t ** duplicated_cell_idx,
793 size_t ** orig_cell_idx, size_t * nbr_duplicated_cells, int with_vertices) {
794
795 size_t cell_mask_size;
796 int * cell_mask;
798 mask_filename, grid_name, &cell_mask_size, &cell_mask);
799
800 for (size_t i = 0; i < cell_mask_size; ++i)
801 cell_mask[i] = cell_mask[i] == valid_mask_value;
802
804 grid_filename, grid_name, cell_mask_size, cell_mask,
805 num_vertices, num_cells, num_vertices_per_cell, cell_to_vertex,
806 x_vertices, y_vertices, x_cells, y_cells, duplicated_cell_idx,
807 orig_cell_idx, nbr_duplicated_cells, with_vertices);
808
809 if (with_vertices) {
810 for (size_t i = 0; i < *num_vertices; ++i) {
811 (*x_vertices)[i] *= YAC_RAD;
812 (*y_vertices)[i] *= YAC_RAD;
813 }
814 }
815 if ((x_cells != NULL) && (y_cells != NULL)) {
816 for (size_t i = 0; i < *num_cells; ++i) {
817 (*x_cells)[i] *= YAC_RAD;
818 (*y_cells)[i] *= YAC_RAD;
819 }
820 }
821
823 !with_vertices || (*num_vertices <= INT_MAX),
824 "ERROR(yac_read_scrip_grid_information): "
825 "too man vertices in grid")
827 *num_cells <= INT_MAX,
828 "ERROR(yac_read_scrip_grid_information): "
829 "too man cells in grid")
830
832}
833
835 char const * grid_filename, char const * mask_filename,
836 char const * grid_name, int valid_mask_value,
837 size_t * num_vertices, size_t * num_cells, int ** num_vertices_per_cell,
838 double ** x_vertices, double ** y_vertices,
839 double ** x_cells, double ** y_cells,
840 int ** cell_to_vertex, int ** cell_core_mask, size_t ** duplicated_cell_idx,
841 size_t ** orig_cell_idx, size_t * nbr_duplicated_cells) {
842
844 grid_filename, mask_filename, grid_name, valid_mask_value,
845 num_vertices, num_cells, num_vertices_per_cell, x_vertices, y_vertices,
846 x_cells, y_cells, cell_to_vertex, cell_core_mask,
847 duplicated_cell_idx, orig_cell_idx, nbr_duplicated_cells, 1);
848}
849
851 MPI_Comm comm, int * io_rank_idx, int * num_io_ranks) {
852
853 int local_is_io, * io_ranks;
854 yac_get_io_ranks(comm, &local_is_io, &io_ranks, num_io_ranks);
855
856 int comm_rank;
857 yac_mpi_call(MPI_Comm_rank(comm, &comm_rank), comm);
858
859 *io_rank_idx = 0;
860 while ((*io_rank_idx < *num_io_ranks) &&
861 (comm_rank != io_ranks[*io_rank_idx]))
862 ++*io_rank_idx;
864 !local_is_io || (*io_rank_idx < *num_io_ranks),
865 "ERROR(get_io_configuration): unable to determine io_rank_idx");
866 free(io_ranks);
867
868 return local_is_io;
869}
870
872 char const * grid_filename, char const * mask_filename,
873 MPI_Comm comm, char const * grid_name, int valid_mask_value,
874 size_t * num_vertices, size_t * num_cells, int ** num_vertices_per_cell,
875 double ** x_vertices, double ** y_vertices, yac_int ** vertex_ids,
876 double ** x_cells, double ** y_cells, yac_int ** cell_ids,
877 int ** cell_to_vertex, int ** cell_core_mask, size_t ** duplicated_cell_idx,
878 yac_int ** orig_cell_ids, size_t * nbr_duplicated_cells, int with_vertices) {
879
880 // if no communicator was provided
881 if (comm == MPI_COMM_NULL) {
882
883 size_t * orig_cell_idx = NULL;
884
886 grid_filename, mask_filename, grid_name, valid_mask_value,
887 num_vertices, num_cells, num_vertices_per_cell,
888 x_vertices, y_vertices, x_cells, y_cells, cell_to_vertex,
889 cell_core_mask, duplicated_cell_idx,
890 (orig_cell_ids != NULL)?&orig_cell_idx:NULL,
891 nbr_duplicated_cells, with_vertices);
892
893 if (with_vertices) {
894 *vertex_ids = xmalloc(*num_vertices * sizeof(**vertex_ids));
895 for (size_t i = 0; i < *num_vertices; ++i) (*vertex_ids)[i] = i;
896 }
897
898 *cell_ids = xmalloc(*num_cells * sizeof(**cell_ids));
899 for (size_t i = 0; i < *num_cells; ++i) (*cell_ids)[i] = i;
900
901 if (orig_cell_ids != NULL) {
902 *orig_cell_ids =
903 xmalloc(*nbr_duplicated_cells * sizeof(**orig_cell_ids));
904 for (size_t i = 0; i < *nbr_duplicated_cells; ++i)
905 (*orig_cell_ids)[i] = (*cell_ids)[orig_cell_idx[i]];
906 free(orig_cell_idx);
907 }
908
909 } else {
910
911 int io_rank_idx, num_io_ranks;
912 get_io_configuration(comm, &io_rank_idx, &num_io_ranks);
913
914 int * cell_mask;
915 size_t cell_mask_size;
917 mask_filename, grid_name, &cell_mask_size, &cell_mask,
918 io_rank_idx, num_io_ranks);
919
920 for (size_t i = 0; i < cell_mask_size; ++i)
921 cell_mask[i] = cell_mask[i] == valid_mask_value;
922
924 grid_filename, grid_name, cell_mask_size, cell_mask,
925 num_vertices, num_cells, num_vertices_per_cell, cell_to_vertex,
926 x_vertices, y_vertices, vertex_ids,
927 x_cells, y_cells, cell_ids, duplicated_cell_idx,
928 orig_cell_ids, nbr_duplicated_cells,
929 io_rank_idx, num_io_ranks, comm, with_vertices);
930
931 if (with_vertices) {
932 for (size_t i = 0; i < *num_vertices; ++i) {
933 (*x_vertices)[i] *= YAC_RAD;
934 (*y_vertices)[i] *= YAC_RAD;
935 }
936 }
937 if ((x_cells != NULL) && (y_cells != NULL)) {
938 for (size_t i = 0; i < *num_cells; ++i) {
939 (*x_cells)[i] *= YAC_RAD;
940 (*y_cells)[i] *= YAC_RAD;
941 }
942 }
943
945 !with_vertices || (*num_vertices <= INT_MAX),
946 "ERROR(yac_read_part_scrip_grid_information): "
947 "too man vertices in grid")
949 *num_cells <= INT_MAX,
950 "ERROR(yac_read_part_scrip_grid_information): "
951 "too man cells in grid")
953 cell_mask_size == *num_cells,
954 "ERROR(yac_read_part_scrip_grid_information): "
955 "mask and grid size do not match "
956 "(mask_file: \"%s\" grid_file: \"%s\" grid_name: \"%s\"",
957 mask_filename, grid_filename, grid_name)
958
960 }
961}
962
964 char const * grid_filename, char const * mask_filename,
965 MPI_Comm comm, char const * grid_name, int valid_mask_value, int use_ll,
966 double ** x_cells, double ** y_cells, size_t ** duplicated_cell_idx,
967 yac_int ** orig_cell_global_ids, size_t * nbr_duplicated_cells) {
968
969 size_t num_vertices;
970 size_t num_cells;
971
973 double * x_vertices;
974 double * y_vertices;
977 int * cell_to_vertex;
978 int * cell_core_mask;
979
981 grid_filename, mask_filename, comm, grid_name, valid_mask_value,
983 &x_vertices, &y_vertices, &vertex_ids, x_cells, y_cells, &cell_ids,
984 &cell_to_vertex, &cell_core_mask, duplicated_cell_idx,
985 orig_cell_global_ids, nbr_duplicated_cells, 1);
986
987 // check whether any cell has more then four vertices
988 // --> write out a warning and switch to gc edges
989 if (use_ll) {
990 for (size_t i = 0; (i < num_cells) && use_ll; ++i)
991 use_ll = !cell_core_mask[i] || (num_vertices_per_cell[i] <= 4);
992 int comm_rank = 0;
993 if (comm != MPI_COMM_NULL) {
995 MPI_Allreduce(
996 MPI_IN_PLACE, &use_ll, 1, MPI_INT, MPI_MIN, comm), comm);
997 yac_mpi_call(MPI_Comm_rank(comm, &comm_rank), comm);
998 }
999 if (!use_ll && (comm_rank == 0))
1000 fprintf(
1001 stderr, "WARNING(yac_read_scrip_basic_grid_data_): "
1002 "grid \"%s\" from \"%s\": required LL but stored as GC "
1003 "(>4 vertices per cell)\n",
1004 grid_name, grid_filename);
1005 }
1006
1007 struct yac_basic_grid_data
1008 (*generate_basic_grid_data_unstruct_ptr)(
1009 size_t, size_t, int *, double *, double *, int *) =
1010 (use_ll)?
1013
1015 generate_basic_grid_data_unstruct_ptr(
1017 x_vertices, y_vertices, cell_to_vertex);
1018
1020 free(x_vertices);
1021 free(y_vertices);
1022 free(cell_to_vertex);
1023
1024 free(grid_data.core_cell_mask);
1025 grid_data.core_cell_mask = cell_core_mask;
1026
1027 free(grid_data.vertex_ids);
1028 grid_data.vertex_ids = vertex_ids;
1029 free(grid_data.cell_ids);
1030 grid_data.cell_ids = cell_ids;
1031
1032 return grid_data;
1033}
1034
1036 char const * grid_filename, char const * mask_filename,
1037 char const * grid_name, int valid_mask_value, int use_ll_edges) {
1038
1039 return
1041 grid_filename, mask_filename, MPI_COMM_NULL,
1042 grid_name, valid_mask_value, use_ll_edges,
1043 NULL, NULL, NULL, NULL, NULL);
1044}
1045
1047 char const * grid_filename, char const * mask_filename,
1048 MPI_Comm comm, char const * grid_name, int valid_mask_value,
1049 char const * name, int use_ll_edges, size_t * cell_coord_idx,
1050 size_t ** duplicated_cell_idx, yac_int ** orig_cell_global_ids,
1051 size_t * nbr_duplicated_cells) {
1052
1053 // generate basic grid
1054 double * x_cells, * y_cells;
1055 struct yac_basic_grid * basic_grid =
1057 name,
1059 grid_filename, mask_filename, comm, grid_name, valid_mask_value,
1060 use_ll_edges, (cell_coord_idx != NULL)?&x_cells:NULL,
1061 (cell_coord_idx != NULL)?&y_cells:NULL, duplicated_cell_idx,
1062 orig_cell_global_ids, nbr_duplicated_cells));
1063
1064 if (cell_coord_idx != NULL) {
1065
1066 // generate 3d cell center coordinates
1067 size_t num_cells =
1070 xmalloc(num_cells * sizeof(*cell_coords));
1071 for (size_t i = 0; i < num_cells; ++i)
1072 LLtoXYZ(x_cells[i], y_cells[i], cell_coords[i]);
1073 free(x_cells);
1074 free(y_cells);
1075
1076 // register cell center coordinates in basic grid
1077 *cell_coord_idx =
1079 basic_grid, YAC_LOC_CELL, cell_coords);
1080 }
1081
1082 return basic_grid;
1083}
1084
1086 char const * grid_filename, char const * mask_filename,
1087 MPI_Fint comm, char const * grid_name, int valid_mask_value,
1088 char const * name, int use_ll_edges, size_t * cell_coord_idx,
1089 size_t ** duplicated_cell_idx, yac_int ** orig_cell_global_ids,
1090 size_t * nbr_duplicated_cells) {
1091
1092 return
1094 grid_filename, mask_filename, MPI_Comm_f2c(comm), grid_name,
1095 valid_mask_value, name, use_ll_edges, cell_coord_idx,
1096 duplicated_cell_idx, orig_cell_global_ids, nbr_duplicated_cells);
1097}
1098
1100 char const * grid_filename, char const * mask_filename,
1101 char const * grid_name, int valid_mask_value, char const * name,
1102 int use_ll_edges, size_t * cell_coord_idx,
1103 size_t ** duplicated_cell_idx, yac_int ** orig_cell_global_ids,
1104 size_t * nbr_duplicated_cells) {
1105
1106 return
1108 grid_filename, mask_filename, MPI_COMM_NULL, grid_name, valid_mask_value,
1109 name, use_ll_edges, cell_coord_idx, duplicated_cell_idx,
1110 orig_cell_global_ids, nbr_duplicated_cells);
1111}
1112
1114 char const * grid_filename, char const * mask_filename, MPI_Comm comm,
1115 char const * grid_name, int valid_mask_value, size_t ** duplicated_vertex_idx,
1116 yac_int ** orig_vertex_global_ids, size_t * nbr_duplicated_vertices) {
1117
1118 size_t file_num_cells;
1119
1120 double * file_x_cells;
1121 double * file_y_cells;
1122 yac_int * file_cell_ids;
1123 int * file_cell_core_mask;
1124
1126 grid_filename, mask_filename, comm, grid_name, valid_mask_value,
1127 NULL, &file_num_cells, NULL, NULL, NULL, NULL,
1128 &file_x_cells, &file_y_cells, &file_cell_ids,
1129 NULL, &file_cell_core_mask, duplicated_vertex_idx,
1130 orig_vertex_global_ids, nbr_duplicated_vertices, 0);
1131
1134 file_num_cells, file_x_cells, file_y_cells);
1135
1136 free(file_x_cells);
1137 free(file_y_cells);
1138
1139 free(grid_data.core_vertex_mask);
1140 grid_data.core_vertex_mask = file_cell_core_mask;
1141
1142 free(grid_data.vertex_ids);
1143 grid_data.vertex_ids = file_cell_ids;
1144
1145 return grid_data;
1146}
1147
1149 char const * grid_filename, char const * mask_filename,
1150 char const * grid_name, int valid_mask_value) {
1151
1152 return
1154 grid_filename, mask_filename, MPI_COMM_NULL,
1155 grid_name, valid_mask_value, NULL, NULL, NULL);
1156}
1157
1159 char const * grid_filename, char const * mask_filename,
1160 MPI_Comm comm, char const * grid_name, int valid_mask_value,
1161 char const * name, size_t * vertex_coord_idx,
1162 size_t ** duplicated_vertex_idx, yac_int ** orig_vertex_global_ids,
1163 size_t * nbr_duplicated_vertices) {
1164
1165 // generate basic grid
1166 struct yac_basic_grid_data basic_grid_data =
1168 grid_filename, mask_filename, comm, grid_name, valid_mask_value,
1169 duplicated_vertex_idx, orig_vertex_global_ids, nbr_duplicated_vertices);
1170 struct yac_basic_grid * basic_grid =
1171 yac_basic_grid_new(name, basic_grid_data);
1172
1173 if (vertex_coord_idx != NULL) {
1174
1175 size_t num_vertices = basic_grid_data.num_vertices;
1176 yac_coordinate_pointer vertex_coords =
1177 xmalloc(num_vertices * sizeof(*vertex_coords));
1178 memcpy(
1179 vertex_coords, basic_grid_data.vertex_coordinates,
1180 num_vertices * sizeof(*vertex_coords));
1181 *vertex_coord_idx =
1183 basic_grid, YAC_LOC_CORNER, vertex_coords);
1184 }
1185
1186 return basic_grid;
1187}
1188
1190 char const * grid_filename, char const * mask_filename,
1191 MPI_Fint comm, char const * grid_name, int valid_mask_value,
1192 char const * name, size_t * vertex_coord_idx,
1193 size_t ** duplicated_vertex_idx, yac_int ** orig_vertex_global_ids,
1194 size_t * nbr_duplicated_vertices) {
1195
1196 return
1198 grid_filename, mask_filename, MPI_Comm_f2c(comm), grid_name,
1199 valid_mask_value, name, vertex_coord_idx, duplicated_vertex_idx,
1200 orig_vertex_global_ids, nbr_duplicated_vertices);
1201}
1202
1204 char const * grid_filename, char const * mask_filename,
1205 char const * grid_name, int valid_mask_value, char const * name,
1206 size_t * vertex_coord_idx, size_t ** duplicated_vertex_idx,
1207 yac_int ** orig_vertex_global_ids, size_t * nbr_duplicated_vertices) {
1208
1209 return
1211 grid_filename, mask_filename, MPI_COMM_NULL, grid_name, valid_mask_value,
1212 name, vertex_coord_idx, duplicated_vertex_idx, orig_vertex_global_ids,
1213 nbr_duplicated_vertices);
1214}
1215
1217 char const * grid_filename, MPI_Comm comm, char const * grid_name) {
1218
1219 int local_is_io, * io_ranks, num_io_ranks;
1220 yac_get_io_ranks(comm, &local_is_io, &io_ranks, &num_io_ranks);
1221 int root_rank = io_ranks[0];
1222
1223 int comm_rank;
1224 yac_mpi_call(MPI_Comm_rank(comm, &comm_rank), comm);
1225
1226 int contains_corner_information = 0;
1227
1228 if (local_is_io && (comm_rank == root_rank)) {
1229
1230 size_t grid_name_len = strlen(grid_name) + 1;
1231 char cla_var_name[4 + grid_name_len];
1232 char clo_var_name[4 + grid_name_len];
1233
1234 snprintf(cla_var_name, 4 + grid_name_len, "%s.cla", grid_name);
1235 snprintf(clo_var_name, 4 + grid_name_len, "%s.clo", grid_name);
1236
1237 int ncid;
1238 yac_nc_open(grid_filename, NC_NOWRITE, &ncid);
1239
1240 // get variable ids
1241 int cla_var_id, clo_var_id;
1242
1243 contains_corner_information =
1244 (NC_NOERR == nc_inq_varid(ncid, cla_var_name, &cla_var_id)) &&
1245 (NC_NOERR == nc_inq_varid(ncid, clo_var_name, &clo_var_id));
1246
1247 YAC_HANDLE_ERROR(nc_close(ncid));
1248 }
1249
1251 MPI_Bcast(
1252 &contains_corner_information, 1, MPI_INT, root_rank, comm), comm);
1253
1254 return contains_corner_information;
1255}
1256
1258 char const * grid_filename, char const * mask_filename,
1259 MPI_Comm comm, char const * grid_name, int valid_mask_value,
1260 char const * name, int use_ll_edges, size_t * point_coord_idx,
1261 size_t ** duplicated_point_idx, yac_int ** orig_point_global_ids,
1262 size_t * nbr_duplicated_points, int * point_location) {
1263
1264 int contains_corner_information =
1265 check_corner_information(grid_filename, comm, grid_name);
1266
1267 if (contains_corner_information) {
1268 *point_location = YAC_LOC_CELL;
1269 return
1271 grid_filename, mask_filename, comm, grid_name,
1272 valid_mask_value, name, use_ll_edges, point_coord_idx,
1273 duplicated_point_idx, orig_point_global_ids, nbr_duplicated_points);
1274 } else {
1275 *point_location = YAC_LOC_CORNER;
1276 return
1278 grid_filename, mask_filename, comm, grid_name,
1279 valid_mask_value, name, point_coord_idx, duplicated_point_idx,
1280 orig_point_global_ids, nbr_duplicated_points);
1281 }
1282}
1283
1285 char const * grid_filename, char const * mask_filename,
1286 MPI_Fint comm, char const * grid_name, int valid_mask_value,
1287 char const * name, int use_ll_edges, size_t * point_coord_idx,
1288 size_t ** duplicated_point_idx, yac_int ** orig_point_global_ids,
1289 size_t * nbr_duplicated_points, int * point_location) {
1290
1291 return
1293 grid_filename, mask_filename, MPI_Comm_f2c(comm), grid_name,
1294 valid_mask_value, name, use_ll_edges, point_coord_idx,
1295 duplicated_point_idx, orig_point_global_ids, nbr_duplicated_points,
1296 point_location);
1297}
1298
1299/* ---------------------------------------------------------------- */
1300
1302 const void * a,const void * b) {
1303
1304 int ret = ((struct cell_vertices_with_index const *)a)->num_vertices -
1305 ((struct cell_vertices_with_index const *)b)->num_vertices;
1306 if (ret) return ret;
1307 int count = ((struct cell_vertices_with_index const *)a)->num_vertices;
1308 for (int i = 0; (i < count) && !ret; ++i)
1309 ret = ((struct cell_vertices_with_index const *)a)->cell_to_vertex[i] -
1310 ((struct cell_vertices_with_index const *)b)->cell_to_vertex[i];
1311 return ret;
1312}
1313
1315 double ** vertex_lon, double ** vertex_lat, yac_int ** vertex_ids,
1316 size_t * nbr_vertices, int * old_to_new_id) {
1317
1318 struct point_with_index * sort_array =
1319 xmalloc(*nbr_vertices * sizeof(*sort_array));
1320
1321 double const scale = (double)(2 << 20);
1322 int32_t const periode = (int32_t)(360.0 * scale);
1323
1324 for (size_t i = 0; i < *nbr_vertices; ++i) {
1325
1326 int32_t curr_lon = (int32_t)round((*vertex_lon)[i] * scale);
1327 int32_t curr_lat = (int32_t)round((*vertex_lat)[i] * scale);
1328
1329 if ((curr_lat == (int32_t)(90.0 * scale)) ||
1330 (curr_lat == (int32_t)(-90.0 * scale))) {
1331
1332 curr_lon = 0;
1333
1334 } else {
1335 if (curr_lon <= 0)
1336 curr_lon = curr_lon - (((curr_lon - periode) / periode) * periode);
1337 else if (curr_lon > periode)
1338 curr_lon = curr_lon - ((curr_lon / periode) * periode);
1339 }
1340
1341 sort_array[i].lon = curr_lon;
1342 sort_array[i].lat = curr_lat;
1343 sort_array[i].dlon = (*vertex_lon)[i];
1344 sort_array[i].dlat = (*vertex_lat)[i];
1345 sort_array[i].id = (vertex_ids != NULL)?(*vertex_ids)[i]:YAC_INT_MAX;
1346
1347 sort_array[i].i = i;
1348 }
1349
1350 yac_mergesort(sort_array, *nbr_vertices, sizeof(*sort_array),
1352
1353 struct point_with_index dummy =
1354 {.lon = INT32_MAX, .lat = INT32_MAX, .i = SIZE_MAX};
1355 struct point_with_index * prev = &dummy, * curr = sort_array;
1356 size_t new_nbr_vertices = 0;
1357 for (size_t i = 0; i < *nbr_vertices; ++i, ++curr) {
1358
1359 if (compare_point_with_index_coord(prev, curr)) {
1360
1361 (*vertex_lon)[new_nbr_vertices] = curr->dlon;
1362 (*vertex_lat)[new_nbr_vertices] = curr->dlat;
1363 if (vertex_ids != NULL)
1364 (*vertex_ids)[new_nbr_vertices] = curr->id;
1365 sort_array[new_nbr_vertices] = *curr;
1366 new_nbr_vertices++;
1367 prev = curr;
1368 }
1369 old_to_new_id[curr->i] = (int)(new_nbr_vertices - 1);
1370 }
1371
1372 (*vertex_lon) = xrealloc(*vertex_lon, new_nbr_vertices * sizeof(**vertex_lon));
1373 (*vertex_lat) = xrealloc(*vertex_lat, new_nbr_vertices * sizeof(**vertex_lat));
1374 if (vertex_ids != NULL)
1375 (*vertex_ids) = xrealloc(*vertex_ids, new_nbr_vertices * sizeof(**vertex_ids));
1376 *nbr_vertices = new_nbr_vertices;
1377
1378 return sort_array;
1379}
1380
1382 double ** vertex_lon, double ** vertex_lat,
1383 size_t * nbr_vertices, int * old_to_new_id) {
1384
1385 free(
1387 vertex_lon, vertex_lat, NULL, nbr_vertices, old_to_new_id));
1388}
1389
1390// djb2 hash function
1391unsigned long djb2_hash(unsigned char * values, size_t count) {
1392
1393 unsigned long hash = 5381;
1394
1395 for (size_t i = 0; i < count; ++i) {
1396 unsigned int value = values[i];
1397 hash = ((hash << 5) + hash) + value; /* hash * 33 + value */
1398 }
1399
1400 return hash;
1401}
1402
1404 void * data, size_t data_size, int comm_size) {
1405
1406 return
1407 djb2_hash((unsigned char*)data, data_size) % (unsigned long)comm_size;
1408}
1409
1410static MPI_Datatype get_point_with_index_mpi_datatype(MPI_Comm comm) {
1411
1412 struct point_with_index dummy;
1413 MPI_Datatype point_with_index_dt;
1414 int array_of_blocklengths[] = {1, 1, 1, 1, 1};
1415 const MPI_Aint array_of_displacements[] =
1416 {(MPI_Aint)(intptr_t)(const void *)&(dummy.lon) -
1417 (MPI_Aint)(intptr_t)(const void *)&dummy,
1418 (MPI_Aint)(intptr_t)(const void *)&(dummy.lat) -
1419 (MPI_Aint)(intptr_t)(const void *)&dummy,
1420 (MPI_Aint)(intptr_t)(const void *)&(dummy.id) -
1421 (MPI_Aint)(intptr_t)(const void *)&dummy,
1422 (MPI_Aint)(intptr_t)(const void *)&(dummy.dlon) -
1423 (MPI_Aint)(intptr_t)(const void *)&dummy,
1424 (MPI_Aint)(intptr_t)(const void *)&(dummy.dlat) -
1425 (MPI_Aint)(intptr_t)(const void *)&dummy};
1426 const MPI_Datatype array_of_types[] =
1427 {MPI_INT32_T, MPI_INT32_T, yac_int_dt, MPI_DOUBLE, MPI_DOUBLE};
1429 MPI_Type_create_struct(5, array_of_blocklengths, array_of_displacements,
1430 array_of_types, &point_with_index_dt), comm);
1431 return yac_create_resized(point_with_index_dt, sizeof(dummy), comm);
1432}
1433
1435 double ** vertex_lon, double ** vertex_lat,
1436 size_t * nbr_vertices_, int * old_to_new_id, yac_int ** vertex_ids,
1437 MPI_Comm comm) {
1438
1439 char const * routine = "remove_duplicated_vertices_parallel";
1440
1441 // remove the locally duplicated vertices
1442 struct point_with_index * sort_array =
1444 vertex_lon, vertex_lat, vertex_ids, nbr_vertices_, old_to_new_id);
1445
1446 sort_array = xrealloc(sort_array, *nbr_vertices_ * sizeof(*sort_array));
1447
1448 size_t nbr_vertices = *nbr_vertices_;
1449
1450 size_t * sendcounts, * recvcounts, * sdispls, * rdispls;
1452 1, &sendcounts, &recvcounts, &sdispls, &rdispls, comm);
1453
1454 int comm_rank, comm_size;
1455 yac_mpi_call(MPI_Comm_rank(comm, &comm_rank), comm);
1456 yac_mpi_call(MPI_Comm_size(comm, &comm_size), comm);
1457
1458 // set up counts and displs for redistribution of points
1459 for (size_t i = 0; i < nbr_vertices; ++i) {
1460 int rank =
1462 &(sort_array[i].lon), 2 * sizeof(sort_array[i].lon), comm_size);
1463 sendcounts[rank]++;
1464 sort_array[i].rank = rank;
1465 sort_array[i].i = i;
1466 }
1467
1468 // sort points by bucket ranks
1469 yac_mergesort(sort_array, nbr_vertices, sizeof(*sort_array),
1471
1473 1, sendcounts, recvcounts, sdispls, rdispls, comm);
1474
1475 size_t recv_count = recvcounts[comm_size - 1] + rdispls[comm_size - 1];
1476 struct point_with_index * dist_sort_array =
1477 xmalloc(recv_count * sizeof(*dist_sort_array));
1478
1479 // redistribute point information
1480 MPI_Datatype point_with_index_dt =
1482 yac_mpi_call(MPI_Type_commit(&point_with_index_dt), comm);
1484 sort_array, sendcounts, sdispls+1,
1485 dist_sort_array, recvcounts, rdispls,
1486 sizeof(*sort_array), point_with_index_dt, comm, routine, __LINE__);
1487
1488 for (size_t i = 0; i < recv_count; ++i) dist_sort_array[i].i = i;
1489
1490 // sort received point data based on coordinates
1491 yac_mergesort(dist_sort_array, recv_count, sizeof(*dist_sort_array),
1493
1494 struct point_with_index dummy =
1495 {.lon = INT32_MAX, .lat = INT32_MAX};
1496 struct point_with_index * prev = &dummy, * curr = dist_sort_array;
1497 for (size_t i = 0; i < recv_count; ++i, ++curr) {
1498 if (compare_point_with_index_coord(prev, curr)) {
1499
1500 prev = curr;
1501
1502 } else {
1503
1504 curr->id = prev->id;
1505 curr->dlon = prev->dlon;
1506 curr->dlat = prev->dlat;
1507 }
1508 }
1509
1510 // bring array into original order
1511 yac_mergesort(dist_sort_array, recv_count, sizeof(*dist_sort_array),
1513
1514 // redistribute point information
1516 dist_sort_array, recvcounts, rdispls,
1517 sort_array, sendcounts, sdispls+1,
1518 sizeof(*sort_array), point_with_index_dt, comm, routine, __LINE__);
1519 free(dist_sort_array);
1520 yac_free_comm_buffers(sendcounts, recvcounts, sdispls, rdispls);
1521 yac_mpi_call(MPI_Type_free(&point_with_index_dt), comm);
1522
1523 // bring vertex indices into original order
1524 for (size_t i = 0; i < nbr_vertices; ++i) {
1525
1526 (*vertex_lon)[sort_array[i].i] = sort_array[i].dlon;
1527 (*vertex_lat)[sort_array[i].i] = sort_array[i].dlat;
1528 (*vertex_ids)[sort_array[i].i] = sort_array[i].id;
1529 }
1530
1531 free(sort_array);
1532}
1533
1535 int * cell_to_vertex, int * num_vertices_per_cell, int * cell_mask,
1536 size_t nbr_cells, size_t ** duplicated_cell_idx,
1537 size_t ** orig_cell_idx, size_t * nbr_duplicated_cells) {
1538
1539 // count number of unmasked cells
1540 size_t nbr_unmasked_cells = 0;
1541 for (size_t i = 0; i < nbr_cells; ++i)
1542 if (cell_mask[i]) ++nbr_unmasked_cells;
1543
1544 struct cell_vertices_with_index * sort_array =
1545 xmalloc(nbr_unmasked_cells * sizeof(*sort_array));
1546
1547 // get all unmasked cells
1548 for (size_t i = 0, offset = 0, j = 0; i < nbr_cells; ++i) {
1549
1550 if (cell_mask[i]) {
1551
1552 sort_array[j].num_vertices = num_vertices_per_cell[i];
1553 sort_array[j].cell_to_vertex = cell_to_vertex + offset;
1554 sort_array[j].i = i;
1555 ++j;
1556 }
1557
1558 offset += (size_t)(num_vertices_per_cell[i]);
1559 }
1560
1561 // sort cells
1562 yac_mergesort(sort_array, nbr_unmasked_cells, sizeof(*sort_array),
1564
1565 // count number of duplicated cells
1566 *nbr_duplicated_cells = 0;
1567 for (size_t i = 1; i < nbr_unmasked_cells; ++i)
1568 if (!compare_cell_vertices_with_index(sort_array + (i-1), sort_array + i))
1569 ++*nbr_duplicated_cells;
1570
1571 // get duplicated cells
1572 *duplicated_cell_idx =
1573 xmalloc(*nbr_duplicated_cells * sizeof(**duplicated_cell_idx));
1574 *orig_cell_idx =
1575 xmalloc(*nbr_duplicated_cells * sizeof(**orig_cell_idx));
1576 struct cell_vertices_with_index *prev = sort_array, *curr = sort_array + 1;
1577 for (size_t i = 1, j = 0; i < nbr_unmasked_cells; ++i, ++curr) {
1578 if (compare_cell_vertices_with_index(prev, curr)) {
1579 prev = curr;
1580 } else {
1581 (*duplicated_cell_idx)[j] = (size_t)(curr->i);
1582 (*orig_cell_idx)[j] = (size_t)(prev->i);
1583 ++j;
1584 }
1585 }
1586
1587 free(sort_array);
1588}
1589
1591 MPI_Comm comm, size_t count) {
1592
1593 MPI_Datatype cell_to_vertex_ids_dt;
1595 MPI_Type_contiguous(
1596 (int)count, yac_int_dt, &cell_to_vertex_ids_dt), comm);
1597 return cell_to_vertex_ids_dt;
1598}
1599
1601 int * cell_to_vertex, int * num_vertices_per_cell, int * cell_mask,
1602 yac_int * cell_ids, size_t num_cells, yac_int * vertex_ids,
1603 MPI_Comm comm, size_t ** duplicated_cell_idx_,
1604 yac_int ** orig_cell_ids_, size_t * num_duplicated_cells_) {
1605
1606 char const * routine = "detect_duplicated_cells_parallel";
1607
1608 int comm_size;
1609 yac_mpi_call(MPI_Comm_size(comm, &comm_size), comm);
1610
1611 // determine the local maximum number of vertices per cell and count
1612 // number of unmasked cells
1613 size_t max_num_vertices_per_cell = 0;
1614 size_t local_num_cells = 0;
1615 for (size_t i = 0; i < num_cells; ++i) {
1616 if (cell_mask[i]) {
1617 if ((size_t)(num_vertices_per_cell[i]) > max_num_vertices_per_cell)
1618 max_num_vertices_per_cell = (size_t)(num_vertices_per_cell[i]);
1619 local_num_cells++;
1620 }
1621 }
1622
1623 // determine the global maximum number of vertices per cell
1625 MPI_Allreduce(
1626 MPI_IN_PLACE, &max_num_vertices_per_cell, 1, YAC_MPI_SIZE_T, MPI_MAX,
1627 comm), comm);
1628
1629 // generate local cell to vertex ids array
1630 yac_int * local_cell_ids =
1631 xmalloc(local_num_cells * (1 + max_num_vertices_per_cell) *
1632 sizeof(*local_cell_ids));
1633 yac_int * local_cell_to_vertex_ids = local_cell_ids + local_num_cells;
1634
1635 // fill local cell to vertex ids array
1636 yac_int * curr_local_cell_ids = local_cell_ids;
1637 yac_int * curr_local_cell_to_vertex_ids = local_cell_to_vertex_ids;
1638 int * curr_cell_to_vertex = cell_to_vertex;
1639 for (size_t i = 0; i < num_cells; ++i) {
1640
1641 size_t curr_num_vertices_per_cell = (size_t)(num_vertices_per_cell[i]);
1642 if (cell_mask[i]) {
1643
1644 *curr_local_cell_ids = cell_ids[i];
1645
1646 // get global vertex ids for current cell
1647 size_t j = 0;
1648 for (; j < curr_num_vertices_per_cell; ++j)
1649 curr_local_cell_to_vertex_ids[j] =
1650 vertex_ids[curr_cell_to_vertex[j]];
1651
1652 // sort global vertex ids
1653 qsort(curr_local_cell_to_vertex_ids, curr_num_vertices_per_cell,
1654 sizeof(*curr_local_cell_to_vertex_ids), compare_yac_int);
1655
1656 // fill remaining entries
1657 for (; j < max_num_vertices_per_cell; ++j)
1658 curr_local_cell_to_vertex_ids[j] = YAC_INT_MAX;
1659
1660 curr_local_cell_ids++;
1661 curr_local_cell_to_vertex_ids += max_num_vertices_per_cell;
1662 }
1663 curr_cell_to_vertex += curr_num_vertices_per_cell;
1664 }
1665
1666 size_t * sendcounts, * recvcounts, * sdispls, * rdispls;
1668 1, &sendcounts, &recvcounts, &sdispls, &rdispls, comm);
1669
1670 // determine distributed owners for all cells
1671 int * dist_cell_ranks = xmalloc(local_num_cells * sizeof(*dist_cell_ranks));
1672 for (size_t i = 0; i < local_num_cells; ++i) {
1673 int rank =
1675 local_cell_to_vertex_ids + i * max_num_vertices_per_cell,
1676 max_num_vertices_per_cell * sizeof(*local_cell_to_vertex_ids),
1677 comm_size);
1678 sendcounts[rank]++;
1679 dist_cell_ranks[i] = rank;
1680 }
1681
1683 1, sendcounts, recvcounts, sdispls, rdispls, comm);
1684
1685 size_t dist_num_cells = recvcounts[comm_size-1] + rdispls[comm_size-1];
1686 yac_int * yac_int_buffer =
1687 xmalloc((local_num_cells + dist_num_cells) *
1688 (1 + max_num_vertices_per_cell) * sizeof(*yac_int_buffer));
1689 yac_int * dist_cell_ids = yac_int_buffer;
1690 yac_int * dist_cell_vertex_ids = yac_int_buffer + dist_num_cells;
1691 yac_int * send_local_cell_ids =
1692 yac_int_buffer + dist_num_cells * (1 + max_num_vertices_per_cell);
1693 yac_int * send_local_cell_to_vertex_ids =
1694 yac_int_buffer + local_num_cells +
1695 dist_num_cells * (1 + max_num_vertices_per_cell);
1696
1697 // pack send data
1698 for (size_t i = 0; i < local_num_cells; ++i) {
1699 size_t pos = sdispls[dist_cell_ranks[i] + 1]++;
1700 send_local_cell_ids[pos] = local_cell_ids[i];
1701 memcpy(
1702 send_local_cell_to_vertex_ids + pos * max_num_vertices_per_cell,
1703 local_cell_to_vertex_ids + i * max_num_vertices_per_cell,
1704 max_num_vertices_per_cell * sizeof(*send_local_cell_to_vertex_ids));
1705 }
1706 local_cell_ids =
1707 xrealloc(local_cell_ids, local_num_cells * sizeof(*local_cell_ids));
1708
1709 // redistribute cell ids and cell vertex ids
1710 MPI_Datatype cell_to_vertex_ids_dt =
1711 get_cell_to_vertex_ids_mpi_datatype(comm, max_num_vertices_per_cell);
1712 yac_mpi_call(MPI_Type_commit(&cell_to_vertex_ids_dt), comm);
1713 yac_alltoallv_yac_int_p2p(
1714 send_local_cell_ids, sendcounts, sdispls,
1715 dist_cell_ids, recvcounts, rdispls, comm, routine, __LINE__);
1717 send_local_cell_to_vertex_ids, sendcounts, sdispls,
1718 dist_cell_vertex_ids, recvcounts, rdispls,
1719 max_num_vertices_per_cell * sizeof(*dist_cell_vertex_ids),
1720 cell_to_vertex_ids_dt, comm, routine, __LINE__);
1721 yac_mpi_call(MPI_Type_free(&cell_to_vertex_ids_dt), comm);
1722
1723 // free some memory
1724 yac_int_buffer =
1725 xrealloc(
1726 yac_int_buffer,
1727 dist_num_cells * (1 + max_num_vertices_per_cell) *
1728 sizeof(*yac_int_buffer));
1729 dist_cell_ids = yac_int_buffer;
1730 dist_cell_vertex_ids = yac_int_buffer + dist_num_cells;
1731
1732 // generate sortable data structure
1733 struct cell_vertex_ids_with_index * sort_array =
1734 xmalloc(dist_num_cells * sizeof(*sort_array));
1735 for (size_t i = 0; i < dist_num_cells; ++i) {
1736 sort_array[i].cell_id = dist_cell_ids[i];
1737 sort_array[i].vertex_ids =
1738 dist_cell_vertex_ids + i * max_num_vertices_per_cell;
1739 sort_array[i].num_vertices = max_num_vertices_per_cell;
1740 sort_array[i].i = i;
1741 }
1742
1743 // sort cells by vertex ids
1745 sort_array, dist_num_cells, sizeof(*sort_array),
1747
1748 // determine duplicated cells
1749 struct cell_vertex_ids_with_index * prev = sort_array;
1750 struct cell_vertex_ids_with_index * curr = sort_array + 1;
1751 for (size_t i = 1; i < dist_num_cells; ++i, ++curr) {
1752
1754 prev = curr;
1755 } else {
1756 curr->cell_id = prev->cell_id;
1757 }
1758 }
1759
1760 // pack cell ids
1761 yac_int_buffer =
1762 xrealloc(
1763 yac_int_buffer,
1764 (local_num_cells + dist_num_cells) * sizeof(*yac_int_buffer));
1765 yac_int * recv_local_cell_ids = yac_int_buffer;
1766 dist_cell_ids = yac_int_buffer + local_num_cells;
1767 for (size_t i = 0; i < dist_num_cells; ++i)
1768 dist_cell_ids[sort_array[i].i] = sort_array[i].cell_id;
1769 free(sort_array);
1770
1771 // redistribute cell ids
1772 yac_alltoallv_yac_int_p2p(
1773 dist_cell_ids, recvcounts, rdispls,
1774 recv_local_cell_ids, sendcounts, sdispls, comm, routine, __LINE__);
1775 yac_free_comm_buffers(sendcounts, recvcounts, sdispls, rdispls);
1776
1777 // bring received local cell ids into original order and
1778 // count number of duplicated cells
1779 size_t num_duplicated_cells = 0;
1780 for (size_t i = 0; i < local_num_cells; ++i) {
1781
1782 size_t pos = sdispls[dist_cell_ranks[i]]++;
1783
1784 if (local_cell_ids[i] != recv_local_cell_ids[pos]) {
1785 num_duplicated_cells++;
1786 local_cell_ids[i] = recv_local_cell_ids[pos];
1787 } else {
1788 local_cell_ids[i] = YAC_INT_MAX;
1789 }
1790 }
1791 free(yac_int_buffer);
1792 free(dist_cell_ranks);
1793
1794 size_t * duplicated_cell_idx =
1795 xmalloc(num_duplicated_cells * sizeof(*duplicated_cell_idx));
1796 yac_int * orig_cell_ids =
1797 xmalloc(num_duplicated_cells * sizeof(*orig_cell_ids));
1798
1799 // get duplicated cells
1800 for (size_t i = 0, j = 0, k = 0; i < num_cells; ++i) {
1801
1802 if (cell_mask[i]) {
1803
1804 if (local_cell_ids[j] != YAC_INT_MAX) {
1805 duplicated_cell_idx[k] = i;
1806 orig_cell_ids[k] = local_cell_ids[j];
1807 ++k;
1808 }
1809 ++j;
1810 }
1811 }
1812 free(local_cell_ids);
1813
1814 *duplicated_cell_idx_ = duplicated_cell_idx;
1815 *orig_cell_ids_ = orig_cell_ids;
1816 *num_duplicated_cells_ = num_duplicated_cells;
1817}
1818
1819#else
1820
1822 char const * grid_filename, char const * mask_filename,
1823 char const * grid_name, int valid_mask_value, int use_ll_edges) {
1824
1825 UNUSED(grid_filename);
1826 UNUSED(mask_filename);
1827 UNUSED(grid_name);
1828 UNUSED(valid_mask_value);
1829 UNUSED(use_ll_edges);
1830 die(
1831 "ERROR(yac_read_scrip_basic_grid_data): "
1832 "YAC is built without the NetCDF support");
1833
1834 return
1836 (size_t[]){0,0}, (int[]){0,0}, NULL, NULL);
1837}
1838
1840 char const * grid_filename, char const * mask_filename,
1841 char const * grid_name, int valid_mask_value, char const * name,
1842 int use_ll_edges, size_t * cell_coord_idx,
1843 size_t ** duplicated_cell_idx, yac_int ** orig_cell_global_ids,
1844 size_t * nbr_duplicated_cells) {
1845
1846 UNUSED(grid_filename);
1847 UNUSED(mask_filename);
1848 UNUSED(grid_name);
1849 UNUSED(valid_mask_value);
1850 UNUSED(name);
1851 UNUSED(use_ll_edges);
1852 UNUSED(cell_coord_idx);
1853 UNUSED(duplicated_cell_idx);
1854 UNUSED(orig_cell_global_ids);
1855 UNUSED(nbr_duplicated_cells);
1856 die(
1857 "ERROR(yac_read_scrip_basic_grid): "
1858 "YAC is built without the NetCDF support");
1859
1860 return NULL;
1861}
1862
1864 char const * grid_filename, char const * mask_filename,
1865 MPI_Comm comm, char const * grid_name, int valid_mask_value,
1866 char const * name, int use_ll_edges, size_t * cell_coord_idx,
1867 size_t ** duplicated_cell_idx, yac_int ** orig_cell_global_ids,
1868 size_t * nbr_duplicated_cells) {
1869
1870 UNUSED(grid_filename);
1871 UNUSED(mask_filename);
1872 UNUSED(comm);
1873 UNUSED(grid_name);
1874 UNUSED(valid_mask_value);
1875 UNUSED(name);
1876 UNUSED(use_ll_edges);
1877 UNUSED(cell_coord_idx);
1878 UNUSED(duplicated_cell_idx);
1879 UNUSED(orig_cell_global_ids);
1880 UNUSED(nbr_duplicated_cells);
1881 die(
1882 "ERROR(yac_read_scrip_basic_grid_parallel): "
1883 "YAC is built without the NetCDF support");
1884
1885 return NULL;
1886}
1887
1889 char const * grid_filename, char const * mask_filename,
1890 char const * grid_name, int valid_mask_value) {
1891
1892 UNUSED(grid_filename);
1893 UNUSED(mask_filename);
1894 UNUSED(grid_name);
1895 UNUSED(valid_mask_value);
1896 die(
1897 "ERROR(yac_read_scrip_cloud_basic_grid_data): "
1898 "YAC is built without the NetCDF support");
1899
1900 return
1902 (size_t[]){0,0}, (int[]){0,0}, NULL, NULL);
1903}
1904
1906 char const * grid_filename, char const * mask_filename,
1907 char const * grid_name, int valid_mask_value, char const * name,
1908 size_t * vertex_coord_idx, size_t ** duplicated_vertex_idx,
1909 yac_int ** orig_vertex_global_ids, size_t * nbr_duplicated_vertices) {
1910
1911 UNUSED(grid_filename);
1912 UNUSED(mask_filename);
1913 UNUSED(grid_name);
1914 UNUSED(valid_mask_value);
1915 UNUSED(name);
1916 UNUSED(vertex_coord_idx);
1917 UNUSED(duplicated_vertex_idx);
1918 UNUSED(orig_vertex_global_ids);
1919 UNUSED(nbr_duplicated_vertices);
1920 die(
1921 "ERROR(yac_read_scrip_cloud_basic_grid): "
1922 "YAC is built without the NetCDF support");
1923
1924 return NULL;
1925}
1926
1928 char const * grid_filename, char const * mask_filename,
1929 MPI_Comm comm, char const * grid_name, int valid_mask_value,
1930 char const * name, size_t * vertex_coord_idx,
1931 size_t ** duplicated_vertex_idx, yac_int ** orig_vertex_global_ids,
1932 size_t * nbr_duplicated_vertices) {
1933
1934 UNUSED(grid_filename);
1935 UNUSED(mask_filename);
1936 UNUSED(comm);
1937 UNUSED(grid_name);
1938 UNUSED(valid_mask_value);
1939 UNUSED(name);
1940 UNUSED(vertex_coord_idx);
1941 UNUSED(duplicated_vertex_idx);
1942 UNUSED(orig_vertex_global_ids);
1943 UNUSED(nbr_duplicated_vertices);
1944 die(
1945 "ERROR(yac_read_scrip_cloud_basic_grid_parallel): "
1946 "YAC is built without the NetCDF support");
1947
1948 return NULL;
1949}
1950
1952 char const * grid_filename, char const * mask_filename,
1953 MPI_Comm comm, char const * grid_name, int valid_mask_value,
1954 char const * name, int use_ll_edges, size_t * point_coord_idx,
1955 size_t ** duplicated_point_idx, yac_int ** orig_point_global_ids,
1956 size_t * nbr_duplicated_points, int * point_location) {
1957
1958 UNUSED(grid_filename);
1959 UNUSED(mask_filename);
1960 UNUSED(comm);
1961 UNUSED(grid_name);
1962 UNUSED(valid_mask_value);
1963 UNUSED(name);
1964 UNUSED(use_ll_edges);
1965 UNUSED(point_coord_idx);
1966 UNUSED(duplicated_point_idx);
1967 UNUSED(orig_point_global_ids);
1968 UNUSED(nbr_duplicated_points);
1969 UNUSED(point_location);
1970 die(
1971 "ERROR(yac_read_scrip_generic_basic_grid_parallel): "
1972 "YAC is built without the NetCDF support");
1973
1974 return NULL;
1975}
1976
1978 char const * grid_filename, char const * mask_filename,
1979 char const * grid_name, int valid_mask_value,
1980 size_t * num_vertices, size_t * num_cells, int ** num_vertices_per_cell,
1981 double ** x_vertices, double ** y_vertices,
1982 double ** x_cells, double ** y_cells,
1983 int ** cell_to_vertex, int ** cell_core_mask, size_t ** duplicated_cell_idx,
1984 size_t ** orig_cell_idx, size_t * nbr_duplicated_cells) {
1985
1986 UNUSED(grid_filename);
1987 UNUSED(mask_filename);
1988 UNUSED(grid_name);
1989 UNUSED(valid_mask_value);
1990 UNUSED(num_vertices);
1992 UNUSED(num_vertices_per_cell);
1993 UNUSED(x_vertices);
1994 UNUSED(y_vertices);
1995 UNUSED(x_cells);
1996 UNUSED(y_cells);
1999 UNUSED(duplicated_cell_idx);
2000 UNUSED(orig_cell_idx);
2001 UNUSED(nbr_duplicated_cells);
2002 die(
2003 "ERROR(yac_read_scrip_grid_information): "
2004 "YAC is built without the NetCDF support");
2005}
2006
2008 char const * grid_filename, char const * mask_filename,
2009 MPI_Fint comm, char const * grid_name, int valid_mask_value,
2010 char const * name, int use_ll_edges, size_t * cell_coord_idx,
2011 size_t ** duplicated_cell_idx, yac_int ** orig_cell_global_ids,
2012 size_t * nbr_duplicated_cells) {
2013
2014 UNUSED(grid_filename);
2015 UNUSED(mask_filename);
2016 UNUSED(comm);
2017 UNUSED(grid_name);
2018 UNUSED(valid_mask_value);
2019 UNUSED(name);
2020 UNUSED(use_ll_edges);
2021 UNUSED(cell_coord_idx);
2022 UNUSED(duplicated_cell_idx);
2023 UNUSED(orig_cell_global_ids);
2024 UNUSED(nbr_duplicated_cells);
2025 die(
2026 "ERROR(yac_read_scrip_basic_grid_parallel_f2c): "
2027 "YAC is built without the NetCDF support");
2028 return NULL;
2029}
2030
2032 char const * grid_filename, char const * mask_filename,
2033 MPI_Fint comm, char const * grid_name, int valid_mask_value,
2034 char const * name, size_t * vertex_coord_idx,
2035 size_t ** duplicated_vertex_idx, yac_int ** orig_vertex_global_ids,
2036 size_t * nbr_duplicated_vertices) {
2037
2038 UNUSED(grid_filename);
2039 UNUSED(mask_filename);
2040 UNUSED(comm);
2041 UNUSED(grid_name);
2042 UNUSED(valid_mask_value);
2043 UNUSED(name);
2044 UNUSED(vertex_coord_idx);
2045 UNUSED(duplicated_vertex_idx);
2046 UNUSED(orig_vertex_global_ids);
2047 UNUSED(nbr_duplicated_vertices);
2048 die(
2049 "ERROR(yac_read_scrip_cloud_basic_grid_parallel_f2c): "
2050 "YAC is built without the NetCDF support");
2051 return NULL;
2052}
2053
2055 char const * grid_filename, char const * mask_filename,
2056 MPI_Fint comm, char const * grid_name, int valid_mask_value,
2057 char const * name, int use_ll_edges, size_t * point_coord_idx,
2058 size_t ** duplicated_point_idx, yac_int ** orig_point_global_ids,
2059 size_t * nbr_duplicated_points, int * point_location) {
2060
2061 UNUSED(grid_filename);
2062 UNUSED(mask_filename);
2063 UNUSED(comm);
2064 UNUSED(grid_name);
2065 UNUSED(valid_mask_value);
2066 UNUSED(name);
2067 UNUSED(use_ll_edges);
2068 UNUSED(point_coord_idx);
2069 UNUSED(duplicated_point_idx);
2070 UNUSED(orig_point_global_ids);
2071 UNUSED(nbr_duplicated_points);
2072 UNUSED(point_location);
2073 die(
2074 "ERROR(yac_read_scrip_generic_basic_grid_parallel_f2c): "
2075 "YAC is built without the NetCDF support");
2076 return NULL;
2077}
2078
2079#endif // YAC_NETCDF_ENABLED
#define YAC_ASSERT(exp, msg)
struct yac_basic_grid * yac_basic_grid_new(char const *name, struct yac_basic_grid_data grid_data)
Definition basic_grid.c:53
size_t yac_basic_grid_add_coordinates_nocpy(struct yac_basic_grid *grid, enum yac_location location, yac_coordinate_pointer coordinates)
Definition basic_grid.c:211
size_t yac_basic_grid_get_data_size(struct yac_basic_grid *grid, enum yac_location location)
Definition basic_grid.c:150
struct yac_basic_grid_data yac_generate_basic_grid_data_reg_2d(size_t nbr_vertices[2], int cyclic[2], double *lon_vertices, double *lat_vertices)
Definition grid_reg2d.c:65
struct yac_basic_grid_data yac_generate_basic_grid_data_unstruct_ll(size_t nbr_vertices, size_t nbr_cells, int *num_vertices_per_cell, double *x_vertices, double *y_vertices, int *cell_to_vertex)
struct yac_basic_grid_data yac_generate_basic_grid_data_cloud(size_t nbr_points, double *x_points, double *y_points)
Definition grid_cloud.c:50
struct yac_basic_grid_data yac_generate_basic_grid_data_unstruct(size_t nbr_vertices, size_t nbr_cells, int *num_vertices_per_cell, double *x_vertices, double *y_vertices, int *cell_to_vertex)
#define UNUSED(x)
Definition core.h:72
#define YAC_RAD
struct @34::@35 value
void yac_get_io_ranks(MPI_Comm comm, int *local_is_io_, int **io_ranks_, int *num_io_ranks_)
Definition io_utils.c:309
void yac_nc_inq_varid(int ncid, char const *name, int *varidp)
Definition io_utils.c:411
void yac_nc_open(const char *path, int omode, int *ncidp)
Definition io_utils.c:350
@ YAC_LOC_CORNER
Definition location.h:15
@ YAC_LOC_CELL
Definition location.h:14
void yac_mergesort(void *base, size_t num, size_t size, int(*compar)(const void *, const void *))
Definition mergesort.c:134
#define xrealloc(ptr, size)
Definition ppm_xfuncs.h:67
#define xmalloc(size)
Definition ppm_xfuncs.h:66
struct yac_basic_grid * yac_read_scrip_cloud_basic_grid_parallel_f2c(char const *grid_filename, char const *mask_filename, MPI_Fint comm, char const *grid_name, int valid_mask_value, char const *name, size_t *vertex_coord_idx, size_t **duplicated_vertex_idx, yac_int **orig_vertex_global_ids, size_t *nbr_duplicated_vertices)
struct yac_basic_grid * yac_read_scrip_basic_grid_parallel(char const *grid_filename, char const *mask_filename, MPI_Comm comm, char const *grid_name, int valid_mask_value, char const *name, int use_ll_edges, size_t *cell_coord_idx, size_t **duplicated_cell_idx, yac_int **orig_cell_global_ids, size_t *nbr_duplicated_cells)
unsigned long djb2_hash(unsigned char *values, size_t count)
static void yac_read_part_scrip_grid_information(char const *grid_filename, char const *mask_filename, MPI_Comm comm, char const *grid_name, int valid_mask_value, size_t *num_vertices, size_t *num_cells, int **num_vertices_per_cell, double **x_vertices, double **y_vertices, yac_int **vertex_ids, double **x_cells, double **y_cells, yac_int **cell_ids, int **cell_to_vertex, int **cell_core_mask, size_t **duplicated_cell_idx, yac_int **orig_cell_ids, size_t *nbr_duplicated_cells, int with_vertices)
static void yac_read_part_scrip_basic_grid_information(const char *filename, const char *grid_name, size_t cell_mask_size, int *cell_mask, size_t *num_vertices_, size_t *num_cells_, int **num_vertices_per_cell_, int **cell_to_vertex_, double **x_vertices_, double **y_vertices_, yac_int **vertex_ids_, double **x_cells_, double **y_cells_, yac_int **cell_ids_, size_t **duplicated_cell_idx_, yac_int **orig_cell_ids_, size_t *nbr_duplicated_cells_, int io_rank_idx, int num_io_ranks, MPI_Comm comm, int with_vertices)
struct yac_basic_grid_data yac_read_scrip_basic_grid_data(char const *grid_filename, char const *mask_filename, char const *grid_name, int valid_mask_value, int use_ll_edges)
static void yac_read_scrip_basic_grid_information(const char *filename, const char *grid_name, size_t cell_mask_size, int *cell_mask, size_t *num_vertices_, size_t *num_cells_, int **num_vertices_per_cell_, int **cell_to_vertex_, double **x_vertices, double **y_vertices, double **x_cells, double **y_cells, size_t **duplicated_cell_idx_, size_t **orig_cell_idx_, size_t *nbr_duplicated_cells_, int with_vertices)
static int compare_point_with_index_coord_id(const void *a_, const void *b_)
static void yac_read_part_scrip_mask_information(const char *filename, const char *grid_name, size_t *num_cells_, int **cell_mask, int io_rank_idx, int num_io_ranks)
void yac_read_scrip_grid_information(char const *grid_filename, char const *mask_filename, char const *grid_name, int valid_mask_value, size_t *num_vertices, size_t *num_cells, int **num_vertices_per_cell, double **x_vertices, double **y_vertices, double **x_cells, double **y_cells, int **cell_to_vertex, int **cell_core_mask, size_t **duplicated_cell_idx, size_t **orig_cell_idx, size_t *nbr_duplicated_cells)
static void yac_read_scrip_mask_information(const char *filename, const char *grid_name, size_t *num_cells_, int **cell_mask)
static int compare_cell_vertices_with_index(const void *a, const void *b)
static void detect_duplicated_cells_parallel(int *cell_to_vertex, int *num_vertices_per_cell, int *cell_mask, yac_int *cell_ids, size_t num_cells, yac_int *vertex_ids, MPI_Comm comm, size_t **duplicated_cell_idx, yac_int **orig_cell_ids, size_t *num_duplicated_cells_)
static MPI_Datatype get_point_with_index_mpi_datatype(MPI_Comm comm)
struct yac_basic_grid * yac_read_scrip_generic_basic_grid_parallel_f2c(char const *grid_filename, char const *mask_filename, MPI_Fint comm, char const *grid_name, int valid_mask_value, char const *name, int use_ll_edges, size_t *point_coord_idx, size_t **duplicated_point_idx, yac_int **orig_point_global_ids, size_t *nbr_duplicated_points, int *point_location)
static void remove_duplicated_vertices_parallel(double **vertex_lon, double **vertex_lat, size_t *nbr_vertices, int *old_to_new_id, yac_int **vertex_ids, MPI_Comm comm)
struct yac_basic_grid_data yac_read_scrip_cloud_basic_grid_data(char const *grid_filename, char const *mask_filename, char const *grid_name, int valid_mask_value)
static struct point_with_index * remove_duplicated_vertices_(double **vertex_lon, double **vertex_lat, yac_int **vertex_ids, size_t *nbr_vertices, int *old_to_new_id)
struct yac_basic_grid * yac_read_scrip_cloud_basic_grid_parallel(char const *grid_filename, char const *mask_filename, MPI_Comm comm, char const *grid_name, int valid_mask_value, char const *name, size_t *vertex_coord_idx, size_t **duplicated_vertex_idx, yac_int **orig_vertex_global_ids, size_t *nbr_duplicated_vertices)
struct yac_basic_grid * yac_read_scrip_cloud_basic_grid(char const *grid_filename, char const *mask_filename, char const *grid_name, int valid_mask_value, char const *name, size_t *vertex_coord_idx, size_t **duplicated_vertex_idx, yac_int **orig_vertex_global_ids, size_t *nbr_duplicated_vertices)
struct yac_basic_grid * yac_read_scrip_generic_basic_grid_parallel(char const *grid_filename, char const *mask_filename, MPI_Comm comm, char const *grid_name, int valid_mask_value, char const *name, int use_ll_edges, size_t *point_coord_idx, size_t **duplicated_point_idx, yac_int **orig_point_global_ids, size_t *nbr_duplicated_points, int *point_location)
static int check_corner_information(char const *grid_filename, MPI_Comm comm, char const *grid_name)
static int compare_cell_vertex_ids_with_index_ids_id(const void *a_, const void *b_)
static int compute_bucket(void *data, size_t data_size, int comm_size)
struct yac_basic_grid * yac_read_scrip_basic_grid_parallel_f2c(char const *grid_filename, char const *mask_filename, MPI_Fint comm, char const *grid_name, int valid_mask_value, char const *name, int use_ll_edges, size_t *cell_coord_idx, size_t **duplicated_cell_idx, yac_int **orig_cell_global_ids, size_t *nbr_duplicated_cells)
static struct yac_basic_grid_data yac_read_scrip_basic_grid_data_(char const *grid_filename, char const *mask_filename, MPI_Comm comm, char const *grid_name, int valid_mask_value, int use_ll, double **x_cells, double **y_cells, size_t **duplicated_cell_idx, yac_int **orig_cell_global_ids, size_t *nbr_duplicated_cells)
static int compare_yac_int(const void *a, const void *b)
static int compare_point_with_index_rank(const void *a_, const void *b_)
static void detect_duplicated_cells(int *cell_to_vertex, int *num_vertices_per_cell, int *cell_mask, size_t nbr_cells, size_t **duplicated_cell_idx, size_t **orig_cell_idx, size_t *nbr_duplicated_cells)
static MPI_Datatype get_cell_to_vertex_ids_mpi_datatype(MPI_Comm comm, size_t count)
static int compare_cell_vertex_ids_with_index_ids(const void *a_, const void *b_)
static int get_io_configuration(MPI_Comm comm, int *io_rank_idx, int *num_io_ranks)
static int compare_point_with_index_i(const void *a_, const void *b_)
static struct yac_basic_grid_data yac_read_scrip_cloud_basic_grid_data_(char const *grid_filename, char const *mask_filename, MPI_Comm comm, char const *grid_name, int valid_mask_value, size_t **duplicated_vertex_idx, yac_int **orig_vertex_global_ids, size_t *nbr_duplicated_vertices)
static int compare_int(const void *a, const void *b)
static void yac_read_scrip_grid_information_(char const *grid_filename, char const *mask_filename, char const *grid_name, int valid_mask_value, size_t *num_vertices, size_t *num_cells, int **num_vertices_per_cell, double **x_vertices, double **y_vertices, double **x_cells, double **y_cells, int **cell_to_vertex, int **cell_core_mask, size_t **duplicated_cell_idx, size_t **orig_cell_idx, size_t *nbr_duplicated_cells, int with_vertices)
static void remove_duplicated_vertices(double **vertex_lon, double **vertex_lat, size_t *nbr_vertices, int *old_to_new_id)
struct yac_basic_grid * yac_read_scrip_basic_grid(char const *grid_filename, char const *mask_filename, char const *grid_name, int valid_mask_value, char const *name, int use_ll_edges, size_t *cell_coord_idx, size_t **duplicated_cell_idx, yac_int **orig_cell_global_ids, size_t *nbr_duplicated_cells)
static int compare_point_with_index_coord(const void *a_, const void *b_)
yac_coordinate_pointer vertex_coordinates
int * cell_to_vertex
int * cell_mask
double * data
size_t num_cells[2]
double cell_coords[36][3]
int * cell_core_mask
#define YAC_HANDLE_ERROR(exp)
Definition toy_output.c:13
char const * name
Definition toy_scrip.c:114
static void LLtoXYZ(double lon, double lat, double p_out[])
Definition toy_scrip.c:587
void yac_quicksort_index_size_t_int(size_t *a, size_t n, int *idx)
#define YAC_ASSERT_F(exp, format,...)
Definition yac_assert.h:30
#define die(msg)
Definition yac_assert.h:12
void yac_generate_alltoallv_args(int count, size_t const *sendcounts, size_t *recvcounts, size_t *sdispls, size_t *rdispls, MPI_Comm comm)
Definition yac_mpi.c:578
void yac_free_comm_buffers(size_t *sendcounts, size_t *recvcounts, size_t *sdispls, size_t *rdispls)
Definition yac_mpi.c:634
void yac_get_comm_buffers(int count, size_t **sendcounts, size_t **recvcounts, size_t **sdispls, size_t **rdispls, MPI_Comm comm)
Definition yac_mpi.c:603
MPI_Datatype yac_create_resized(MPI_Datatype dt, size_t new_size, MPI_Comm comm)
Definition yac_mpi.c:557
void yac_alltoallv_p2p(void const *send_buffer, size_t const *sendcounts, size_t const *sdispls, void *recv_buffer, size_t const *recvcounts, size_t const *rdispls, size_t dt_size, MPI_Datatype dt, MPI_Comm comm, char const *caller, int line)
Definition yac_mpi.c:132
#define yac_mpi_call(call, comm)
#define YAC_MPI_SIZE_T
YAC_INT yac_int
Definition yac_types.h:15
#define yac_int_dt
Definition yac_types.h:18
double(* yac_coordinate_pointer)[3]
Definition yac_types.h:21