23#ifdef YAC_NETCDF_ENABLED
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);
58 size_t nbr_cells,
size_t ** duplicated_cell_idx,
59 size_t ** orig_cell_idx,
size_t * nbr_duplicated_cells);
63 MPI_Comm comm,
size_t ** duplicated_cell_idx,
64 yac_int ** orig_cell_ids,
size_t * num_duplicated_cells_);
67 const void * a,
const void * b) {
68 return *(
int const *)a - *(
int const*)b;
72 const void * a,
const void * b) {
78 const void * a_,
const void * b_) {
87 int ret = (num_vertices_a > num_vertices_b) -
88 (num_vertices_a < num_vertices_b);
90 for (
size_t i = 0; (
i < num_vertices_a) && !ret; ++
i)
98 const void * a_,
const void * b_) {
112 const void * a_,
const void * b_) {
121 return (a->
id > b->
id) - (a->
id < b->
id);
125 const void * a_,
const void * b_) {
136 const void * a_,
const void * b_) {
145 const void * a_,
const void * b_) {
150 return (a->
i > b->
i) - (a->
i < b->
i);
154 const char * filename,
const char * grid_name,
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) {
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];
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);
189 int dim_ids[NC_MAX_VAR_DIMS];
193 nc_inq_var(ncid, cla_var_id, NULL, NULL, &ndims, dim_ids, NULL));
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);
202 nc_inq_var(ncid, lat_var_id, NULL, NULL, &ndims, dim_ids, NULL));
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];
212 int crn_dim_id = dim_ids[0];
213 int x_dim_id = dim_ids[2];
214 int y_dim_id = dim_ids[1];
217 size_t crn_dim_len, x_dim_len, y_dim_len;
225 size_t num_cells = x_dim_len * y_dim_len;
226 size_t num_vertices =
num_cells * crn_dim_len;
230 "ERROR(yac_read_scrip_basic_grid_information): "
231 "cell mask size is inconsistent with number of grid cells")
233 if ((x_cells != NULL) && (y_cells != NULL)) {
243 double * cla =
xmalloc(num_vertices *
sizeof(*cla));
244 double * clo =
xmalloc(num_vertices *
sizeof(*clo));
255 (x_cells != NULL) && (y_cells != NULL),
256 "ERROR(yac_read_scrip_basic_grid_information): "
257 "no pointers for cell coordinates were provided")
259 memcpy(cla, *y_cells,
num_cells *
sizeof(*cla));
260 memcpy(clo, *x_cells,
num_cells *
sizeof(*clo));
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;
281 int * num_vertices_per_cell =
283 size_t total_num_cell_vertices = 0;
286 if (crn_dim_len > 1) {
288 size_t curr_num_vertices = 0;
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;
303 from_vertices += crn_dim_len;
305 num_vertices_per_cell[
i] = (int)curr_num_vertices;
306 total_num_cell_vertices += curr_num_vertices;
311 if (to_vertices != from_vertices) *to_vertices = *from_vertices;
313 ++total_num_cell_vertices;
314 num_vertices_per_cell[
i] = 1;
316 num_vertices_per_cell[
i] = 0;
322 if (total_num_cell_vertices !=
num_cells * crn_dim_len)
329 int * cell_to_vertex_copy =
330 xmalloc(total_num_cell_vertices *
sizeof(*cell_to_vertex_copy));
333 total_num_cell_vertices *
sizeof(*cell_to_vertex_copy));
336 cell_to_vertex_copy + offset, (
size_t)(num_vertices_per_cell[
i]),
338 offset += num_vertices_per_cell[
i];
342 size_t * duplicated_cell_idx;
343 size_t * orig_cell_idx;
344 size_t nbr_duplicated_cells;
347 &duplicated_cell_idx, &orig_cell_idx, &nbr_duplicated_cells);
348 free(cell_to_vertex_copy);
351 for (
size_t i = 0;
i < nbr_duplicated_cells; ++
i)
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;
360 free(duplicated_cell_idx);
366 *num_vertices_ = num_vertices;
367 *num_vertices_per_cell_ = num_vertices_per_cell;
374 free(num_vertices_per_cell);
384 const char * filename,
const char * grid_name,
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) {
395 size_t num_vertices = 0;
396 size_t * reorder_idx = NULL;
397 size_t max_num_vertices_per_cell = 0;
399 double * x_vertices = NULL;
400 double * y_vertices = NULL;
402 double * x_cells = NULL;
403 double * y_cells = NULL;
406 if ((io_rank_idx >= 0) && (io_rank_idx < num_io_ranks)) {
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];
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);
423 int cla_var_id, clo_var_id, lat_var_id, lon_var_id;
432 int dim_ids[NC_MAX_VAR_DIMS];
436 nc_inq_var(ncid, cla_var_id, NULL, NULL, &ndims, dim_ids, NULL));
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);
445 nc_inq_var(ncid, lat_var_id, NULL, NULL, &ndims, dim_ids, NULL));
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];
455 int crn_dim_id = dim_ids[0];
456 int x_dim_id = dim_ids[2];
457 int y_dim_id = dim_ids[1];
460 size_t crn_dim_len, x_dim_len, y_dim_len;
469 size_t x_start_local =
471 (((
unsigned long)x_dim_len * (
unsigned long)io_rank_idx) /
472 (
unsigned long)num_io_ranks);
473 size_t x_count_local =
475 (((
unsigned long)x_dim_len * (io_rank_idx + 1)) /
476 (
unsigned long)num_io_ranks) - x_start_local;
483 max_num_vertices_per_cell = crn_dim_len;
488 max_num_vertices_per_cell = 1;
493 "ERROR(yac_read_part_scrip_basic_grid_information): "
494 "cell mask size is inconsistent with number of grid cells")
496 if ((x_cells_ != NULL) && (y_cells_ != NULL)) {
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));
510 x_vertices =
xmalloc(num_vertices *
sizeof(*x_vertices));
511 y_vertices =
xmalloc(num_vertices *
sizeof(*y_vertices));
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));
527 (x_cells_ != NULL) && (y_cells_ != NULL),
528 "ERROR(yac_read_part_scrip_basic_grid_information): "
529 "no pointers for cell coordinates were provided")
531 memcpy(x_vertices, x_cells,
num_cells *
sizeof(*x_vertices));
532 memcpy(y_vertices, y_cells,
num_cells *
sizeof(*y_vertices));
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)
542 x + y * x_count_local + n * x_count_local * y_dim_len] = l;
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;
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);
559 for (
size_t n = 0, l = 0; n < crn_dim_len; ++n)
561 vertex_ids[l] = cell_ids[
i] * crn_dim_len + n;
567 &x_vertices, &y_vertices, &num_vertices,
576 int * num_vertices_per_cell =
578 size_t total_num_cell_vertices = 0;
581 if (max_num_vertices_per_cell > 1) {
583 size_t curr_num_vertices = 0;
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;
598 from_vertices += max_num_vertices_per_cell;
600 num_vertices_per_cell[
i] = (int)curr_num_vertices;
601 total_num_cell_vertices += curr_num_vertices;
606 if (to_vertices != from_vertices) *to_vertices = *from_vertices;
608 ++total_num_cell_vertices;
609 num_vertices_per_cell[
i] = 1;
611 num_vertices_per_cell[
i] = 0;
617 if (total_num_cell_vertices !=
num_cells * max_num_vertices_per_cell)
623 size_t * duplicated_cell_idx;
625 size_t nbr_duplicated_cells;
628 vertex_ids, comm, &duplicated_cell_idx, &orig_cell_ids,
629 &nbr_duplicated_cells);
632 for (
size_t i = 0;
i < nbr_duplicated_cells; ++
i)
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;
641 free(duplicated_cell_idx);
647 *num_vertices_ = num_vertices;
648 *num_vertices_per_cell_ = num_vertices_per_cell;
649 *x_vertices_ = x_vertices;
650 *y_vertices_ = y_vertices;
652 *vertex_ids_ = vertex_ids;
656 free(num_vertices_per_cell);
663 if ((x_cells_ != NULL) && (y_cells_ != NULL)) {
667 *cell_ids_ = cell_ids;
671 const char * filename,
const char * grid_name,
size_t * num_cells_,
674 size_t grid_name_len = strlen(grid_name) + 1;
675 char msk_var_name[4 + grid_name_len];
677 snprintf(msk_var_name, 4 + grid_name_len,
"%s.msk", grid_name);
688 int dim_ids[NC_MAX_VAR_DIMS];
690 nc_inq_var(ncid, msk_var_id, NULL, NULL, &ndims, dim_ids, NULL));
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];
705 size_t num_cells = x_dim_len * y_dim_len;
719 const char * filename,
const char * grid_name,
size_t * num_cells_,
720 int **
cell_mask,
int io_rank_idx,
int num_io_ranks) {
722 if ((io_rank_idx >= 0) && (io_rank_idx < num_io_ranks)) {
724 size_t grid_name_len = strlen(grid_name) + 1;
725 char msk_var_name[4 + grid_name_len];
727 snprintf(msk_var_name, 4 + grid_name_len,
"%s.msk", grid_name);
738 int dim_ids[NC_MAX_VAR_DIMS];
740 nc_inq_var(ncid, msk_var_id, NULL, NULL, &ndims, dim_ids, NULL));
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];
756 size_t x_start_local =
758 (((
unsigned long)x_dim_len * (
unsigned long)io_rank_idx) /
759 (
unsigned long)num_io_ranks);
760 size_t x_count_local =
762 (((
unsigned long)x_dim_len * (io_rank_idx+1)) /
763 (
unsigned long)num_io_ranks) - x_start_local;
765 size_t num_cells = x_count_local * y_dim_len;
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));
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,
793 size_t ** orig_cell_idx,
size_t * nbr_duplicated_cells,
int with_vertices) {
795 size_t cell_mask_size;
798 mask_filename, grid_name, &cell_mask_size, &
cell_mask);
800 for (
size_t i = 0;
i < cell_mask_size; ++
i)
804 grid_filename, grid_name, cell_mask_size,
cell_mask,
806 x_vertices, y_vertices, x_cells, y_cells, duplicated_cell_idx,
807 orig_cell_idx, nbr_duplicated_cells, with_vertices);
810 for (
size_t i = 0;
i < *num_vertices; ++
i) {
815 if ((x_cells != NULL) && (y_cells != NULL)) {
823 !with_vertices || (*num_vertices <= INT_MAX),
824 "ERROR(yac_read_scrip_grid_information): "
825 "too man vertices in grid")
828 "ERROR(yac_read_scrip_grid_information): "
829 "too man cells in grid")
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,
841 size_t ** orig_cell_idx,
size_t * nbr_duplicated_cells) {
844 grid_filename, mask_filename, grid_name, valid_mask_value,
845 num_vertices,
num_cells, num_vertices_per_cell, x_vertices, y_vertices,
847 duplicated_cell_idx, orig_cell_idx, nbr_duplicated_cells, 1);
851 MPI_Comm comm,
int * io_rank_idx,
int * num_io_ranks) {
853 int local_is_io, * io_ranks;
860 while ((*io_rank_idx < *num_io_ranks) &&
861 (comm_rank != io_ranks[*io_rank_idx]))
864 !local_is_io || (*io_rank_idx < *num_io_ranks),
865 "ERROR(get_io_configuration): unable to determine io_rank_idx");
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,
878 yac_int ** orig_cell_ids,
size_t * nbr_duplicated_cells,
int with_vertices) {
881 if (comm == MPI_COMM_NULL) {
883 size_t * orig_cell_idx = NULL;
886 grid_filename, mask_filename, grid_name, valid_mask_value,
887 num_vertices,
num_cells, num_vertices_per_cell,
890 (orig_cell_ids != NULL)?&orig_cell_idx:NULL,
891 nbr_duplicated_cells, 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;
901 if (orig_cell_ids != NULL) {
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]];
911 int io_rank_idx, num_io_ranks;
915 size_t cell_mask_size;
917 mask_filename, grid_name, &cell_mask_size, &
cell_mask,
918 io_rank_idx, num_io_ranks);
920 for (
size_t i = 0;
i < cell_mask_size; ++
i)
924 grid_filename, grid_name, cell_mask_size,
cell_mask,
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);
932 for (
size_t i = 0;
i < *num_vertices; ++
i) {
937 if ((x_cells != NULL) && (y_cells != NULL)) {
945 !with_vertices || (*num_vertices <= INT_MAX),
946 "ERROR(yac_read_part_scrip_grid_information): "
947 "too man vertices in grid")
950 "ERROR(yac_read_part_scrip_grid_information): "
951 "too man cells in grid")
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)
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) {
981 grid_filename, mask_filename, comm, grid_name, valid_mask_value,
985 orig_cell_global_ids, nbr_duplicated_cells, 1);
990 for (
size_t i = 0; (i <
num_cells) && use_ll; ++i)
993 if (comm != MPI_COMM_NULL) {
996 MPI_IN_PLACE, &use_ll, 1, MPI_INT, MPI_MIN, comm), comm);
999 if (!use_ll && (comm_rank == 0))
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);
1008 (*generate_basic_grid_data_unstruct_ptr)(
1009 size_t,
size_t,
int *,
double *,
double *,
int *) =
1015 generate_basic_grid_data_unstruct_ptr(
1036 char const * grid_filename,
char const * mask_filename,
1037 char const * grid_name,
int valid_mask_value,
int use_ll_edges) {
1041 grid_filename, mask_filename, MPI_COMM_NULL,
1042 grid_name, valid_mask_value, use_ll_edges,
1043 NULL, NULL, NULL, NULL, NULL);
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) {
1054 double * x_cells, * y_cells;
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));
1064 if (cell_coord_idx != NULL) {
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) {
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);
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) {
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);
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) {
1118 size_t file_num_cells;
1120 double * file_x_cells;
1121 double * file_y_cells;
1123 int * file_cell_core_mask;
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);
1134 file_num_cells, file_x_cells, file_y_cells);
1140 grid_data.core_vertex_mask = file_cell_core_mask;
1149 char const * grid_filename,
char const * mask_filename,
1150 char const * grid_name,
int valid_mask_value) {
1154 grid_filename, mask_filename, MPI_COMM_NULL,
1155 grid_name, valid_mask_value, NULL, NULL, NULL);
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) {
1168 grid_filename, mask_filename, comm, grid_name, valid_mask_value,
1169 duplicated_vertex_idx, orig_vertex_global_ids, nbr_duplicated_vertices);
1173 if (vertex_coord_idx != NULL) {
1177 xmalloc(num_vertices *
sizeof(*vertex_coords));
1180 num_vertices *
sizeof(*vertex_coords));
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) {
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);
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) {
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);
1217 char const * grid_filename, MPI_Comm comm,
char const * grid_name) {
1219 int local_is_io, * io_ranks, num_io_ranks;
1221 int root_rank = io_ranks[0];
1226 int contains_corner_information = 0;
1228 if (local_is_io && (comm_rank == root_rank)) {
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];
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);
1241 int cla_var_id, clo_var_id;
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));
1252 &contains_corner_information, 1, MPI_INT, root_rank, comm), comm);
1254 return contains_corner_information;
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) {
1264 int contains_corner_information =
1267 if (contains_corner_information) {
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);
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);
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) {
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,
1302 const void * a,
const void * b) {
1306 if (ret)
return ret;
1308 for (
int i = 0; (
i < count) && !ret; ++
i)
1315 double ** vertex_lon,
double ** vertex_lat,
yac_int ** vertex_ids,
1316 size_t * nbr_vertices,
int * old_to_new_id) {
1319 xmalloc(*nbr_vertices *
sizeof(*sort_array));
1321 double const scale = (double)(2 << 20);
1322 int32_t
const periode = (int32_t)(360.0 * scale);
1324 for (
size_t i = 0;
i < *nbr_vertices; ++
i) {
1326 int32_t curr_lon = (int32_t)round((*vertex_lon)[
i] * scale);
1327 int32_t curr_lat = (int32_t)round((*vertex_lat)[
i] * scale);
1329 if ((curr_lat == (int32_t)(90.0 * scale)) ||
1330 (curr_lat == (int32_t)(-90.0 * scale))) {
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);
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;
1347 sort_array[
i].
i =
i;
1350 yac_mergesort(sort_array, *nbr_vertices,
sizeof(*sort_array),
1354 {.
lon = INT32_MAX, .lat = INT32_MAX, .i = SIZE_MAX};
1356 size_t new_nbr_vertices = 0;
1357 for (
size_t i = 0;
i < *nbr_vertices; ++
i, ++curr) {
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;
1369 old_to_new_id[curr->i] = (int)(new_nbr_vertices - 1);
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;
1382 double ** vertex_lon,
double ** vertex_lat,
1383 size_t * nbr_vertices,
int * old_to_new_id) {
1387 vertex_lon, vertex_lat, NULL, nbr_vertices, old_to_new_id));
1391unsigned long djb2_hash(
unsigned char * values,
size_t count) {
1393 unsigned long hash = 5381;
1395 for (
size_t i = 0;
i < count; ++
i) {
1396 unsigned int value = values[
i];
1397 hash = ((hash << 5) + hash) +
value;
1404 void *
data,
size_t data_size,
int comm_size) {
1407 djb2_hash((
unsigned char*)
data, data_size) % (
unsigned long)comm_size;
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);
1435 double ** vertex_lon,
double ** vertex_lat,
1436 size_t * nbr_vertices_,
int * old_to_new_id,
yac_int ** vertex_ids,
1439 char const * routine =
"remove_duplicated_vertices_parallel";
1444 vertex_lon, vertex_lat, vertex_ids, nbr_vertices_, old_to_new_id);
1446 sort_array =
xrealloc(sort_array, *nbr_vertices_ *
sizeof(*sort_array));
1448 size_t nbr_vertices = *nbr_vertices_;
1450 size_t * sendcounts, * recvcounts, * sdispls, * rdispls;
1452 1, &sendcounts, &recvcounts, &sdispls, &rdispls, comm);
1454 int comm_rank, comm_size;
1459 for (
size_t i = 0;
i < nbr_vertices; ++
i) {
1462 &(sort_array[
i].
lon), 2 *
sizeof(sort_array[
i].
lon), comm_size);
1465 sort_array[
i].
i =
i;
1469 yac_mergesort(sort_array, nbr_vertices,
sizeof(*sort_array),
1473 1, sendcounts, recvcounts, sdispls, rdispls, comm);
1475 size_t recv_count = recvcounts[comm_size - 1] + rdispls[comm_size - 1];
1477 xmalloc(recv_count *
sizeof(*dist_sort_array));
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__);
1488 for (
size_t i = 0;
i < recv_count; ++
i) dist_sort_array[
i].
i =
i;
1491 yac_mergesort(dist_sort_array, recv_count,
sizeof(*dist_sort_array),
1495 {.
lon = INT32_MAX, .lat = INT32_MAX};
1497 for (
size_t i = 0;
i < recv_count; ++
i, ++curr) {
1504 curr->
id = prev->
id;
1505 curr->dlon = prev->
dlon;
1506 curr->dlat = prev->
dlat;
1511 yac_mergesort(dist_sort_array, recv_count,
sizeof(*dist_sort_array),
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);
1521 yac_mpi_call(MPI_Type_free(&point_with_index_dt), comm);
1524 for (
size_t i = 0;
i < nbr_vertices; ++
i) {
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;
1536 size_t nbr_cells,
size_t ** duplicated_cell_idx,
1537 size_t ** orig_cell_idx,
size_t * nbr_duplicated_cells) {
1540 size_t nbr_unmasked_cells = 0;
1541 for (
size_t i = 0;
i < nbr_cells; ++
i)
1545 xmalloc(nbr_unmasked_cells *
sizeof(*sort_array));
1548 for (
size_t i = 0, offset = 0, j = 0;
i < nbr_cells; ++
i) {
1554 sort_array[j].
i =
i;
1558 offset += (size_t)(num_vertices_per_cell[
i]);
1562 yac_mergesort(sort_array, nbr_unmasked_cells,
sizeof(*sort_array),
1566 *nbr_duplicated_cells = 0;
1567 for (
size_t i = 1;
i < nbr_unmasked_cells; ++
i)
1569 ++*nbr_duplicated_cells;
1572 *duplicated_cell_idx =
1573 xmalloc(*nbr_duplicated_cells *
sizeof(**duplicated_cell_idx));
1575 xmalloc(*nbr_duplicated_cells *
sizeof(**orig_cell_idx));
1577 for (
size_t i = 1, j = 0;
i < nbr_unmasked_cells; ++
i, ++curr) {
1581 (*duplicated_cell_idx)[j] = (size_t)(curr->i);
1582 (*orig_cell_idx)[j] = (size_t)(prev->
i);
1591 MPI_Comm comm,
size_t count) {
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;
1603 MPI_Comm comm,
size_t ** duplicated_cell_idx_,
1604 yac_int ** orig_cell_ids_,
size_t * num_duplicated_cells_) {
1606 char const * routine =
"detect_duplicated_cells_parallel";
1613 size_t max_num_vertices_per_cell = 0;
1614 size_t local_num_cells = 0;
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]);
1626 MPI_IN_PLACE, &max_num_vertices_per_cell, 1,
YAC_MPI_SIZE_T, MPI_MAX,
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;
1636 yac_int * curr_local_cell_ids = local_cell_ids;
1637 yac_int * curr_local_cell_to_vertex_ids = local_cell_to_vertex_ids;
1641 size_t curr_num_vertices_per_cell = (size_t)(num_vertices_per_cell[
i]);
1644 *curr_local_cell_ids = cell_ids[
i];
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]];
1653 qsort(curr_local_cell_to_vertex_ids, curr_num_vertices_per_cell,
1657 for (; j < max_num_vertices_per_cell; ++j)
1658 curr_local_cell_to_vertex_ids[j] = YAC_INT_MAX;
1660 curr_local_cell_ids++;
1661 curr_local_cell_to_vertex_ids += max_num_vertices_per_cell;
1663 curr_cell_to_vertex += curr_num_vertices_per_cell;
1666 size_t * sendcounts, * recvcounts, * sdispls, * rdispls;
1668 1, &sendcounts, &recvcounts, &sdispls, &rdispls, comm);
1671 int * dist_cell_ranks =
xmalloc(local_num_cells *
sizeof(*dist_cell_ranks));
1672 for (
size_t i = 0;
i < local_num_cells; ++
i) {
1675 local_cell_to_vertex_ids +
i * max_num_vertices_per_cell,
1676 max_num_vertices_per_cell *
sizeof(*local_cell_to_vertex_ids),
1679 dist_cell_ranks[
i] = rank;
1683 1, sendcounts, recvcounts, sdispls, rdispls, comm);
1685 size_t dist_num_cells = recvcounts[comm_size-1] + rdispls[comm_size-1];
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);
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];
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));
1707 xrealloc(local_cell_ids, local_num_cells *
sizeof(*local_cell_ids));
1710 MPI_Datatype cell_to_vertex_ids_dt =
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);
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;
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];
1738 dist_cell_vertex_ids +
i * max_num_vertices_per_cell;
1740 sort_array[
i].
i =
i;
1745 sort_array, dist_num_cells,
sizeof(*sort_array),
1751 for (
size_t i = 1;
i < dist_num_cells; ++
i, ++curr) {
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;
1772 yac_alltoallv_yac_int_p2p(
1773 dist_cell_ids, recvcounts, rdispls,
1774 recv_local_cell_ids, sendcounts, sdispls, comm, routine, __LINE__);
1779 size_t num_duplicated_cells = 0;
1780 for (
size_t i = 0;
i < local_num_cells; ++
i) {
1782 size_t pos = sdispls[dist_cell_ranks[
i]]++;
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];
1788 local_cell_ids[
i] = YAC_INT_MAX;
1791 free(yac_int_buffer);
1792 free(dist_cell_ranks);
1794 size_t * duplicated_cell_idx =
1795 xmalloc(num_duplicated_cells *
sizeof(*duplicated_cell_idx));
1797 xmalloc(num_duplicated_cells *
sizeof(*orig_cell_ids));
1804 if (local_cell_ids[j] != YAC_INT_MAX) {
1805 duplicated_cell_idx[k] =
i;
1806 orig_cell_ids[k] = local_cell_ids[j];
1812 free(local_cell_ids);
1814 *duplicated_cell_idx_ = duplicated_cell_idx;
1815 *orig_cell_ids_ = orig_cell_ids;
1816 *num_duplicated_cells_ = num_duplicated_cells;
1822 char const * grid_filename,
char const * mask_filename,
1823 char const * grid_name,
int valid_mask_value,
int use_ll_edges) {
1828 UNUSED(valid_mask_value);
1831 "ERROR(yac_read_scrip_basic_grid_data): "
1832 "YAC is built without the NetCDF support");
1836 (
size_t[]){0,0}, (
int[]){0,0}, NULL, NULL);
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) {
1849 UNUSED(valid_mask_value);
1853 UNUSED(duplicated_cell_idx);
1854 UNUSED(orig_cell_global_ids);
1855 UNUSED(nbr_duplicated_cells);
1857 "ERROR(yac_read_scrip_basic_grid): "
1858 "YAC is built without the NetCDF support");
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) {
1874 UNUSED(valid_mask_value);
1878 UNUSED(duplicated_cell_idx);
1879 UNUSED(orig_cell_global_ids);
1880 UNUSED(nbr_duplicated_cells);
1882 "ERROR(yac_read_scrip_basic_grid_parallel): "
1883 "YAC is built without the NetCDF support");
1889 char const * grid_filename,
char const * mask_filename,
1890 char const * grid_name,
int valid_mask_value) {
1895 UNUSED(valid_mask_value);
1897 "ERROR(yac_read_scrip_cloud_basic_grid_data): "
1898 "YAC is built without the NetCDF support");
1902 (
size_t[]){0,0}, (
int[]){0,0}, NULL, NULL);
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) {
1914 UNUSED(valid_mask_value);
1916 UNUSED(vertex_coord_idx);
1917 UNUSED(duplicated_vertex_idx);
1918 UNUSED(orig_vertex_global_ids);
1919 UNUSED(nbr_duplicated_vertices);
1921 "ERROR(yac_read_scrip_cloud_basic_grid): "
1922 "YAC is built without the NetCDF support");
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) {
1938 UNUSED(valid_mask_value);
1940 UNUSED(vertex_coord_idx);
1941 UNUSED(duplicated_vertex_idx);
1942 UNUSED(orig_vertex_global_ids);
1943 UNUSED(nbr_duplicated_vertices);
1945 "ERROR(yac_read_scrip_cloud_basic_grid_parallel): "
1946 "YAC is built without the NetCDF support");
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) {
1962 UNUSED(valid_mask_value);
1966 UNUSED(duplicated_point_idx);
1967 UNUSED(orig_point_global_ids);
1968 UNUSED(nbr_duplicated_points);
1971 "ERROR(yac_read_scrip_generic_basic_grid_parallel): "
1972 "YAC is built without the NetCDF support");
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,
1984 size_t ** orig_cell_idx,
size_t * nbr_duplicated_cells) {
1989 UNUSED(valid_mask_value);
1992 UNUSED(num_vertices_per_cell);
1999 UNUSED(duplicated_cell_idx);
2001 UNUSED(nbr_duplicated_cells);
2003 "ERROR(yac_read_scrip_grid_information): "
2004 "YAC is built without the NetCDF support");
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) {
2018 UNUSED(valid_mask_value);
2022 UNUSED(duplicated_cell_idx);
2023 UNUSED(orig_cell_global_ids);
2024 UNUSED(nbr_duplicated_cells);
2026 "ERROR(yac_read_scrip_basic_grid_parallel_f2c): "
2027 "YAC is built without the NetCDF support");
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) {
2042 UNUSED(valid_mask_value);
2044 UNUSED(vertex_coord_idx);
2045 UNUSED(duplicated_vertex_idx);
2046 UNUSED(orig_vertex_global_ids);
2047 UNUSED(nbr_duplicated_vertices);
2049 "ERROR(yac_read_scrip_cloud_basic_grid_parallel_f2c): "
2050 "YAC is built without the NetCDF support");
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) {
2065 UNUSED(valid_mask_value);
2069 UNUSED(duplicated_point_idx);
2070 UNUSED(orig_point_global_ids);
2071 UNUSED(nbr_duplicated_points);
2074 "ERROR(yac_read_scrip_generic_basic_grid_parallel_f2c): "
2075 "YAC is built without the NetCDF support");
#define YAC_ASSERT(exp, msg)
struct yac_basic_grid * yac_basic_grid_new(char const *name, struct yac_basic_grid_data grid_data)
size_t yac_basic_grid_add_coordinates_nocpy(struct yac_basic_grid *grid, enum yac_location location, yac_coordinate_pointer coordinates)
size_t yac_basic_grid_get_data_size(struct yac_basic_grid *grid, enum yac_location location)
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)
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)
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)
void yac_get_io_ranks(MPI_Comm comm, int *local_is_io_, int **io_ranks_, int *num_io_ranks_)
void yac_nc_inq_varid(int ncid, char const *name, int *varidp)
void yac_nc_open(const char *path, int omode, int *ncidp)
void yac_mergesort(void *base, size_t num, size_t size, int(*compar)(const void *, const void *))
#define xrealloc(ptr, size)
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 * num_vertices_per_cell
double cell_coords[36][3]
#define YAC_HANDLE_ERROR(exp)
static void LLtoXYZ(double lon, double lat, double p_out[])
void yac_quicksort_index_size_t_int(size_t *a, size_t n, int *idx)
#define YAC_ASSERT_F(exp, format,...)
void yac_generate_alltoallv_args(int count, size_t const *sendcounts, size_t *recvcounts, size_t *sdispls, size_t *rdispls, MPI_Comm comm)
void yac_free_comm_buffers(size_t *sendcounts, size_t *recvcounts, size_t *sdispls, size_t *rdispls)
void yac_get_comm_buffers(int count, size_t **sendcounts, size_t **recvcounts, size_t **sdispls, size_t **rdispls, MPI_Comm comm)
MPI_Datatype yac_create_resized(MPI_Datatype dt, size_t new_size, MPI_Comm comm)
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)
#define yac_mpi_call(call, comm)
double(* yac_coordinate_pointer)[3]