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