YAC 3.13.2
Yet Another Coupler
Loading...
Searching...
No Matches
read_icon_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 <string.h>
12#include <assert.h>
13#include <math.h>
14#include <limits.h>
15
16#include <mpi.h>
17
18#include "read_icon_grid.h"
19#include "utils_common.h"
20#include "io_utils.h"
21#include "geometry.h"
22#include "read_grid.h"
23
24#ifdef YAC_NETCDF_ENABLED
25#include <netcdf.h>
26
27static int * get_icon_cell_mask( int ncid, size_t nbr_cells );
28static void get_icon_connect(
29 int ncid, size_t nbr_cells, int ** vertex_of_cell, size_t * nv);
30
32 const char * filename, int * nbr_vertices, int * nbr_cells,
33 int ** num_vertices_per_cell, int ** cell_to_vertex,
34 double ** x_vertices, double ** y_vertices,
35 double ** x_cells, double ** y_cells, int ** global_cell_id,
36 int ** cell_mask, int ** cell_core_mask,
37 int ** global_corner_id, int ** corner_core_mask, int rank, int size) {
38
39 // read the global grid
41 filename, nbr_vertices, nbr_cells, num_vertices_per_cell, cell_to_vertex,
42 x_vertices, y_vertices, x_cells, y_cells, cell_mask);
43
44 // determine local range
45 int local_start =
46 ((unsigned long)(*nbr_cells) * (unsigned long)rank) /
47 (unsigned long)size;
48 int num_local_cells =
49 ((unsigned long)(*nbr_cells) * ((unsigned long)rank+1)) /
50 (unsigned long)size - (unsigned long)local_start;
51
52 // mask for required vertices and cells
53 int * required_vertices = xcalloc(*nbr_vertices, sizeof(*required_vertices));
54 int * required_cells = xcalloc(*nbr_cells, sizeof(*required_cells));
55
56 int offset = 0;
57 for (int i = 0; i < local_start; ++i)
58 offset += (*num_vertices_per_cell)[i];
59
60 // mark all local cells and their vertices as required
61 for (int i = local_start; i < local_start + num_local_cells; ++i) {
62
63 required_cells[i] = 2;
64 for (int j = 0; j < (*num_vertices_per_cell)[i]; ++j)
65 required_vertices[(*cell_to_vertex)[offset+j]] = 2;
66
67 offset += (*num_vertices_per_cell)[i];
68 }
69
70 // mark all halo cells as required
71 offset = 0;
72 for (int i = 0; i < *nbr_cells; ++i) {
73
74 if (!required_cells[i]) {
75
76 for (int j = 0; j < (*num_vertices_per_cell)[i]; ++j) {
77
78 if (required_vertices[(*cell_to_vertex)[offset+j]]) {
79
80 required_cells[i] = 1;
81 break;
82 }
83 }
84
85 }
86
87 offset += (*num_vertices_per_cell)[i];
88 }
89
90 // mark all halo vertices as required
91 offset = 0;
92 for (int i = 0; i < *nbr_cells; ++i) {
93
94 if (required_cells[i] == 1) {
95
96 for (int j = 0; j < (*num_vertices_per_cell)[i]; ++j) {
97
98 if (!required_vertices[(*cell_to_vertex)[offset+j]]) {
99
100 required_vertices[(*cell_to_vertex)[offset+j]] = 1;
101 }
102 }
103 }
104
105 offset += (*num_vertices_per_cell)[i];
106 }
107
108 // count the number of cells and vertices
109 int part_num_vertices = 0;
110 int part_num_cells = 0;
111 for (int i = 0; i < *nbr_vertices; ++i)
112 if (required_vertices[i])
113 part_num_vertices++;
114 for (int i = 0; i < *nbr_cells; ++i)
115 if(required_cells[i])
116 part_num_cells++;
117
118 *global_cell_id = xmalloc(part_num_cells * sizeof(**global_cell_id));
119 *cell_core_mask = xmalloc(part_num_cells * sizeof(**cell_core_mask));
120 *global_corner_id = xmalloc(part_num_vertices * sizeof(**global_corner_id));
121 *corner_core_mask = xmalloc(part_num_vertices * sizeof(**corner_core_mask));
122
123 // generate final vertex data
124 part_num_vertices = 0;
125 int * global_to_local_vertex =
126 xmalloc(*nbr_vertices * sizeof(*global_to_local_vertex));
127 for (int i = 0; i < *nbr_vertices; ++i) {
128
129 if (required_vertices[i]) {
130
131 (*global_corner_id)[part_num_vertices] = i;
132 (*corner_core_mask)[part_num_vertices] = required_vertices[i] == 2;
133 (*x_vertices)[part_num_vertices] = (*x_vertices)[i];
134 (*y_vertices)[part_num_vertices] = (*y_vertices)[i];
135 global_to_local_vertex[i] = part_num_vertices;
136 part_num_vertices++;
137 }
138 }
139
140 *x_vertices = xrealloc(*x_vertices, part_num_vertices * sizeof(**x_vertices));
141 *y_vertices = xrealloc(*y_vertices, part_num_vertices * sizeof(**y_vertices));
142 *nbr_vertices = part_num_vertices;
143 free(required_vertices);
144
145 // generate final cell data
146 int num_cell_vertex_dependencies = 0;
147 part_num_cells = 0;
148 offset = 0;
149 for (int i = 0; i < *nbr_cells; ++i) {
150
151 if (required_cells[i]) {
152
153 (*global_cell_id)[part_num_cells] = i;
154 (*cell_core_mask)[part_num_cells] = required_cells[i] == 2;
155 (*x_cells)[part_num_cells] = (*x_cells)[i];
156 (*y_cells)[part_num_cells] = (*y_cells)[i];
157 (*cell_mask)[part_num_cells] = (*cell_mask)[i];
158
159 for (int j = 0; j < (*num_vertices_per_cell)[i]; ++j)
160 (*cell_to_vertex)[num_cell_vertex_dependencies++] =
161 global_to_local_vertex[(*cell_to_vertex)[offset+j]];
162
163 (*num_vertices_per_cell)[part_num_cells] = (*num_vertices_per_cell)[i];
164
165 part_num_cells++;
166 }
167
168 offset += (*num_vertices_per_cell)[i];
169 }
170
171 *x_cells = xrealloc(*x_cells, part_num_cells * sizeof(**x_cells));
172 *y_cells = xrealloc(*y_cells, part_num_cells * sizeof(**y_cells));
173 *cell_mask = xrealloc(*cell_mask, part_num_cells * sizeof(**cell_mask));
174
175 *num_vertices_per_cell = xrealloc(*num_vertices_per_cell, part_num_cells *
176 sizeof(**num_vertices_per_cell));
177 *cell_to_vertex = xrealloc(*cell_to_vertex, num_cell_vertex_dependencies *
178 sizeof(**cell_to_vertex));
179 *nbr_cells = part_num_cells;
180 free(required_cells);
181 free(global_to_local_vertex);
182}
183
184void yac_read_icon_grid_information(const char * filename, int * nbr_vertices,
185 int * nbr_cells, int ** num_vertices_per_cell,
186 int ** cell_to_vertex, double ** x_vertices,
187 double ** y_vertices, double ** x_cells,
188 double ** y_cells, int ** cell_mask) {
189
190 /* Open file */
191 int ncid;
192 yac_nc_open(filename, NC_NOWRITE, &ncid);
193
194 /* Get vertex longitudes and latitudes of cells */
195 size_t nbr_vertices_;
197 ncid, "vlon", "vlat", x_vertices, y_vertices, &nbr_vertices_);
198 *nbr_vertices = (int)nbr_vertices_;
199
200 /* Get cell center longitudes and latitudes of cells */
201 size_t nbr_cells_;
202 yac_read_coords(ncid, "clon", "clat", x_cells, y_cells, &nbr_cells_);
203 *nbr_cells = (int)nbr_cells_;
204
205 /* Get mask of cells */
206 *cell_mask = get_icon_cell_mask ( ncid, *nbr_cells );
207
208 /* Get relations between vertices and cells */
209 int * vertex_of_cell;
210 size_t nv;
211 get_icon_connect(ncid, *nbr_cells, &vertex_of_cell, &nv);
212
213 /* Close file */
214 YAC_HANDLE_ERROR(nc_close(ncid));
215
216 //-------------------------------------------------------------------------//
217
218 *num_vertices_per_cell = xmalloc(*nbr_cells * sizeof(**num_vertices_per_cell));
219 for (int i = 0; i < *nbr_cells; (*num_vertices_per_cell)[i++] = (int)nv);
220
221 *cell_to_vertex = xmalloc(*nbr_cells * nv * sizeof(**cell_to_vertex));
222
223 // Unfortunately the data is only available in Fortran order
224 for (int i = 0; i < *nbr_cells; ++i) {
225
226 for (size_t j = 0; j < nv; ++j)
227 (*cell_to_vertex)[nv * i + j] =
228 vertex_of_cell[i + j * (*nbr_cells)] - 1;
229 }
230
231 free(vertex_of_cell);
232}
233
234// taken from scales-ppm library
235// https://www.dkrz.de/redmine/projects/scales-ppm
236static inline int
237partition_idx_from_element_idx(unsigned element_idx, unsigned num_elements,
238 int num_partitions) {
239
240 return (int)((((unsigned long)element_idx) * ((unsigned long)num_partitions) +
241 (unsigned long)num_partitions - 1) /
242 ((unsigned long)num_elements));
243}
244
245static void convert_to_rad(int ncid, int varid, double * array, size_t count) {
246
247 if (yac_check_coord_units(ncid, varid))
248 for (size_t i = 0; i < count; ++i) array[i] *= YAC_RAD;
249}
250
252 const char * filename, MPI_Comm comm, int * nbr_vertices, int * nbr_cells,
253 int ** num_vertices_per_cell, int ** cell_to_vertex, int ** global_cell_ids,
254 int ** cell_owner, int ** global_vertex_ids, int ** vertex_owner,
255 double ** x_vertices, double ** y_vertices,
256 double ** x_cells, double ** y_cells, int ** cell_msk) {
257
259 ((x_cells == NULL) && (y_cells == NULL)) ||
260 ((x_cells != NULL) && (y_cells != NULL)),
261 "ERROR(yac_read_icon_grid_information_parallel): "
262 "arguments x_cells and y_cells have to be either both NULL or "
263 "both different from NULL");
264
265 const int with_cell_coords = x_cells != NULL;
266 const int with_cell_mask = cell_msk != NULL;
267
268 int comm_rank, comm_size;
269
270 MPI_Comm_rank(comm, &comm_rank);
271 MPI_Comm_size(comm, &comm_size);
272
273 int local_is_io, * io_ranks, num_io_ranks;
274 yac_get_io_ranks(comm, &local_is_io, &io_ranks, &num_io_ranks);
275
276 size_t ncells, nvertices, nv, ne;
277
278 size_t read_local_start_cell = 0;
279 size_t read_num_local_cells = 0;
280 size_t read_local_start_vertex = 0;
281 size_t read_num_local_vertices = 0;
282 double * read_cell_lon = NULL;
283 double * read_cell_lat = NULL;
284 int * read_cell_mask = NULL;
285 double * read_vertex_lon = NULL;
286 double * read_vertex_lat = NULL;
287 int * read_dist_vertex_of_cell = NULL;
288 int * read_dist_cells_of_vertex = NULL;
289
290 if (local_is_io) {
291
292 unsigned long io_proc_idx = ULONG_MAX;
293 for (int i = 0; (i < num_io_ranks) && (io_proc_idx == ULONG_MAX); ++i)
294 if (io_ranks[i] == comm_rank)
295 io_proc_idx = (unsigned long)i;
296
297 // open file
298 int ncid;
299 yac_nc_open(filename, NC_NOWRITE, &ncid);
300
301 // get number of cells and vertices
302 int dim_id;
303 yac_nc_inq_dimid(ncid, "cell", &dim_id);
304 YAC_HANDLE_ERROR(nc_inq_dimlen(ncid, dim_id, &ncells));
305 yac_nc_inq_dimid(ncid, "vertex", &dim_id);
306 YAC_HANDLE_ERROR(nc_inq_dimlen(ncid, dim_id, &nvertices));
307 yac_nc_inq_dimid(ncid, "nv", &dim_id);
308 YAC_HANDLE_ERROR(nc_inq_dimlen(ncid, dim_id, &nv));
309 yac_nc_inq_dimid(ncid, "ne", &dim_id);
310 YAC_HANDLE_ERROR(nc_inq_dimlen(ncid, dim_id, &ne));
311
312 // determine local range for cell and vertex data
313 read_local_start_cell =
314 ((unsigned long)ncells * io_proc_idx) / (unsigned long)num_io_ranks;
315 read_num_local_cells =
316 ((unsigned long)ncells * (io_proc_idx+1)) / (unsigned long)num_io_ranks -
317 (unsigned long)read_local_start_cell;
318 read_local_start_vertex =
319 ((unsigned long)nvertices * io_proc_idx) / (unsigned long)num_io_ranks;
320 read_num_local_vertices =
321 ((unsigned long)nvertices * (io_proc_idx+1)) / (unsigned long)num_io_ranks -
322 (unsigned long)read_local_start_vertex;
323
324 // read basic grid data (each process its individual part)
325 read_cell_lon =
326 with_cell_coords?
327 xmalloc(read_num_local_cells * sizeof(*read_cell_lon)):NULL;
328 read_cell_lat =
329 with_cell_coords?
330 xmalloc(read_num_local_cells * sizeof(*read_cell_lat)):NULL;
331 read_cell_mask =
332 with_cell_mask?
333 xmalloc(read_num_local_cells * sizeof(*read_cell_mask)):NULL;
334 read_vertex_lon = xmalloc(read_num_local_vertices * sizeof(*read_vertex_lon));
335 read_vertex_lat = xmalloc(read_num_local_vertices * sizeof(*read_vertex_lat));
336 read_dist_vertex_of_cell =
337 xmalloc(read_num_local_cells * nv * sizeof(*read_dist_vertex_of_cell));
338 read_dist_cells_of_vertex =
339 xmalloc(read_num_local_vertices * ne * sizeof(*read_dist_cells_of_vertex));
340 int varid;
341 if (with_cell_coords) {
342 yac_nc_inq_varid(ncid, "clon", &varid);
343 YAC_HANDLE_ERROR(nc_get_vara_double(ncid, varid, &read_local_start_cell,
344 &read_num_local_cells, read_cell_lon));
345 convert_to_rad(ncid, varid, read_cell_lon, read_num_local_cells);
346 yac_nc_inq_varid(ncid, "clat", &varid);
347 YAC_HANDLE_ERROR(nc_get_vara_double(ncid, varid, &read_local_start_cell,
348 &read_num_local_cells, read_cell_lat));
349 convert_to_rad(ncid, varid, read_cell_lat, read_num_local_cells);
350 }
351 if (with_cell_mask) {
352 yac_nc_inq_varid(ncid, "cell_sea_land_mask", &varid);
353 YAC_HANDLE_ERROR(nc_get_vara_int(ncid, varid, &read_local_start_cell,
354 &read_num_local_cells, read_cell_mask));
355 }
356 yac_nc_inq_varid(ncid, "vlon", &varid);
357 YAC_HANDLE_ERROR(nc_get_vara_double(ncid, varid, &read_local_start_vertex,
358 &read_num_local_vertices, read_vertex_lon));
359 convert_to_rad(ncid, varid, read_vertex_lon, read_num_local_vertices);
360 yac_nc_inq_varid(ncid, "vlat", &varid);
361 YAC_HANDLE_ERROR(nc_get_vara_double(ncid, varid, &read_local_start_vertex,
362 &read_num_local_vertices, read_vertex_lat));
363 convert_to_rad(ncid, varid, read_vertex_lat, read_num_local_vertices);
364 {
365 int * buffer = xmalloc(MAX(nv*read_num_local_cells,
366 ne*read_num_local_vertices) * sizeof(*buffer));
367 {
368 size_t tmp_start[2] = {0, read_local_start_cell};
369 size_t tmp_count[2] = {nv, read_num_local_cells};
370 yac_nc_inq_varid(ncid, "vertex_of_cell", &varid);
371 YAC_HANDLE_ERROR(nc_get_vara_int(ncid, varid, tmp_start, tmp_count, buffer));
372 for (size_t i = 0; i < read_num_local_cells; ++i)
373 for (size_t j = 0; j < nv; ++j)
374 read_dist_vertex_of_cell[i * nv + j] =
375 buffer[i + j * read_num_local_cells];
376 }
377 {
378 size_t tmp_start[2] = {0, read_local_start_vertex};
379 size_t tmp_count[2] = {ne, read_num_local_vertices};
380 yac_nc_inq_varid(ncid, "cells_of_vertex", &varid);
381 YAC_HANDLE_ERROR(nc_get_vara_int(ncid, varid, tmp_start, tmp_count, buffer));
382 for (size_t i = 0; i < read_num_local_vertices; ++i)
383 for (size_t j = 0; j < ne; ++j)
384 read_dist_cells_of_vertex[i * ne + j] =
385 buffer[i + j * read_num_local_vertices];
386 }
387 free(buffer);
388
389 // adjust for c indices
390 for (size_t i = 0; i < read_num_local_cells * nv; ++i)
391 if (read_dist_vertex_of_cell[i] > 0) read_dist_vertex_of_cell[i]--;
392 for (size_t i = 0; i < read_num_local_vertices * ne; ++i)
393 read_dist_cells_of_vertex[i]--;
394 }
395
396 YAC_HANDLE_ERROR(nc_close(ncid));
397 } else {
398 read_cell_lon =
399 with_cell_coords?xmalloc(1 * sizeof(*read_cell_lon)):NULL;
400 read_cell_lat =
401 with_cell_coords?xmalloc(1 * sizeof(*read_cell_lat)):NULL;
402 read_cell_mask =
403 with_cell_mask?xmalloc(1 * sizeof(*read_cell_mask)):NULL;
404 read_vertex_lon = xmalloc(1 * sizeof(*read_vertex_lon));
405 read_vertex_lat = xmalloc(1 * sizeof(*read_vertex_lat));
406 read_dist_vertex_of_cell = xmalloc(1 * sizeof(*read_dist_vertex_of_cell));
407 read_dist_cells_of_vertex = xmalloc(1 * sizeof(*read_dist_cells_of_vertex));
408 }
409
410 free(io_ranks);
411
412 {
413 int tmp;
414 if (comm_rank == 0) tmp = (int)ncells;
415 MPI_Bcast(&tmp, 1, MPI_INT, 0, comm);
416 ncells = (size_t)tmp;
417 if (comm_rank == 0) tmp = (int)nvertices;
418 MPI_Bcast(&tmp, 1, MPI_INT, 0, comm);
419 nvertices = (size_t)tmp;
420 if (comm_rank == 0) tmp = (int)nv;
421 MPI_Bcast(&tmp, 1, MPI_INT, 0, comm);
422 nv = (size_t)tmp;
423 if (comm_rank == 0) tmp = (int)ne;
424 MPI_Bcast(&tmp, 1, MPI_INT, 0, comm);
425 ne = (size_t)tmp;
426 }
427
428 // determine local range for cell and vertex data
429 size_t local_start_cell =
430 ((unsigned long)ncells * (unsigned long)comm_rank) /
431 (unsigned long)comm_size;
432 size_t num_local_cells =
433 ((unsigned long)ncells * ((unsigned long)comm_rank+1)) /
434 (unsigned long)comm_size - (unsigned long)local_start_cell;
435 size_t local_start_vertex =
436 ((unsigned long)nvertices * (unsigned long)comm_rank) /
437 (unsigned long)comm_size;
438 size_t num_local_vertices =
439 ((unsigned long)nvertices * ((unsigned long)comm_rank+1)) /
440 (unsigned long)comm_size - (unsigned long)local_start_vertex;
441
442 // redistribute basic cell data (from io decomposition)
443 double * cell_lon =
444 with_cell_coords?xmalloc(num_local_cells * sizeof(*cell_lon)):NULL;
445 double * cell_lat =
446 with_cell_coords?xmalloc(num_local_cells * sizeof(*cell_lat)):NULL;
447 int * cell_mask =
448 with_cell_mask?xmalloc(num_local_cells * sizeof(*cell_mask)):NULL;
449 int * dist_vertex_of_cell =
450 xmalloc(num_local_cells * nv * sizeof(*dist_vertex_of_cell));
451 {
452 int * send_count = xcalloc(comm_size, sizeof(*send_count));
453 int * recv_count = xcalloc(comm_size, sizeof(*recv_count));
454
455 for (size_t i = 0; i < read_num_local_cells; ++i)
456 send_count[
458 read_local_start_cell + i, ncells, comm_size)]++;
459
460 MPI_Alltoall(send_count, 1, MPI_INT, recv_count, 1, MPI_INT, comm);
461
462 int * send_displ = xmalloc(comm_size * sizeof(*send_displ));
463 int * recv_displ = xmalloc(comm_size * sizeof(*recv_displ));
464 int send_accum = 0, recv_accum = 0;
465 for (int i = 0; i < comm_size; ++i) {
466 send_displ[i] = send_accum;
467 recv_displ[i] = recv_accum;
468 send_accum += send_count[i];
469 recv_accum += recv_count[i];
470 }
471
472 if (with_cell_coords) {
473 MPI_Alltoallv(read_cell_lon, send_count, send_displ, MPI_DOUBLE,
474 cell_lon, recv_count, recv_displ, MPI_DOUBLE, comm);
475 MPI_Alltoallv(read_cell_lat, send_count, send_displ, MPI_DOUBLE,
476 cell_lat, recv_count, recv_displ, MPI_DOUBLE, comm);
477 }
478 if (with_cell_mask) {
479 MPI_Alltoallv(read_cell_mask, send_count, send_displ, MPI_INT,
480 cell_mask, recv_count, recv_displ, MPI_INT, comm);
481 }
482
483 for (int i = 0; i < comm_size; ++i) {
484 send_count[i] *= nv;
485 send_displ[i] *= nv;
486 recv_count[i] *= nv;
487 recv_displ[i] *= nv;
488 }
489
490 MPI_Alltoallv(read_dist_vertex_of_cell, send_count, send_displ, MPI_INT,
491 dist_vertex_of_cell, recv_count, recv_displ, MPI_INT, comm);
492
493 free(recv_displ);
494 free(send_displ);
495 free(recv_count);
496 free(send_count);
497 free(read_cell_lon);
498 free(read_cell_lat);
499 free(read_cell_mask);
500 free(read_dist_vertex_of_cell);
501 }
502
503 // redistribute basic vertex data (from io decomposition)
504 double * vertex_lon = xmalloc(num_local_vertices * sizeof(*vertex_lon));
505 double * vertex_lat = xmalloc(num_local_vertices * sizeof(*vertex_lat));
506 int * dist_cells_of_vertex =
507 xmalloc(num_local_vertices * ne * sizeof(*dist_cells_of_vertex));
508 {
509 int * send_count = xcalloc(comm_size, sizeof(*send_count));
510 int * recv_count = xcalloc(comm_size, sizeof(*recv_count));
511
512 for (size_t i = 0; i < read_num_local_vertices; ++i)
513 send_count[
515 read_local_start_vertex + i, nvertices, comm_size)]++;
516
517 MPI_Alltoall(send_count, 1, MPI_INT, recv_count, 1, MPI_INT, comm);
518
519 int * send_displ = xmalloc(comm_size * sizeof(*send_displ));
520 int * recv_displ = xmalloc(comm_size * sizeof(*recv_displ));
521 int send_accum = 0, recv_accum = 0;
522 for (int i = 0; i < comm_size; ++i) {
523 send_displ[i] = send_accum;
524 recv_displ[i] = recv_accum;
525 send_accum += send_count[i];
526 recv_accum += recv_count[i];
527 }
528
529 MPI_Alltoallv(read_vertex_lon, send_count, send_displ, MPI_DOUBLE,
530 vertex_lon, recv_count, recv_displ, MPI_DOUBLE, comm);
531 MPI_Alltoallv(read_vertex_lat, send_count, send_displ, MPI_DOUBLE,
532 vertex_lat, recv_count, recv_displ, MPI_DOUBLE, comm);
533
534 for (int i = 0; i < comm_size; ++i) {
535 send_count[i] *= ne;
536 send_displ[i] *= ne;
537 recv_count[i] *= ne;
538 recv_displ[i] *= ne;
539 }
540
541 MPI_Alltoallv(read_dist_cells_of_vertex, send_count, send_displ, MPI_INT,
542 dist_cells_of_vertex, recv_count, recv_displ, MPI_INT, comm);
543
544 free(recv_displ);
545 free(send_displ);
546 free(recv_count);
547 free(send_count);
548 free(read_vertex_lon);
549 free(read_vertex_lat);
550 free(read_dist_cells_of_vertex);
551 }
552
553 // determine required vertices for core cells
554 size_t num_core_vertices = num_local_cells * nv;
555 int * core_vertices = xmalloc(num_core_vertices * sizeof(*core_vertices));
556 {
557 memcpy(core_vertices, dist_vertex_of_cell,
558 num_core_vertices * sizeof(*core_vertices));
559 yac_quicksort_index(core_vertices, num_core_vertices, NULL);
560 yac_remove_duplicates_int(core_vertices, &num_core_vertices);
561 core_vertices =
562 xrealloc(core_vertices, num_core_vertices * sizeof(*core_vertices));
563 }
564
565 // get cells_of_vertex for core vertices and compute dist_vertex_owner
566 int * cells_of_vertex_core =
567 xmalloc(num_core_vertices * ne * sizeof(*cells_of_vertex_core));
568 int * dist_vertex_owner =
569 xmalloc(num_local_vertices * sizeof(*dist_vertex_owner));
570 {
571 int * send_count = xcalloc(comm_size, sizeof(*send_count));
572 int * recv_count = xcalloc(comm_size, sizeof(*recv_count));
573
574 for (size_t i = 0; i < num_core_vertices; ++i)
575 send_count[((unsigned long)(core_vertices[i]) *
576 (unsigned long)comm_size + (unsigned long)comm_size - 1) /
577 (unsigned long)nvertices]++;
578
579 MPI_Alltoall(send_count, 1, MPI_INT, recv_count, 1, MPI_INT, comm);
580
581 int * send_displ = xmalloc(comm_size * sizeof(*send_displ));
582 int * recv_displ = xmalloc(comm_size * sizeof(*recv_displ));
583 int send_accum = 0, recv_accum = 0;
584 for (int i = 0; i < comm_size; ++i) {
585 send_displ[i] = send_accum;
586 recv_displ[i] = recv_accum;
587 send_accum += send_count[i];
588 recv_accum += recv_count[i];
589 }
590
591 int num_core_vertices_remote = 0;
592 for (int i = 0; i < comm_size; ++i)
593 num_core_vertices_remote += recv_count[i];
594
595 int * remote_vertex_buffer =
596 xmalloc(num_core_vertices_remote * sizeof(*remote_vertex_buffer));
597
598 MPI_Alltoallv(core_vertices, send_count, send_displ, MPI_INT,
599 remote_vertex_buffer, recv_count, recv_displ, MPI_INT, comm);
600
601 for (size_t i = 0; i < num_local_vertices; ++i) dist_vertex_owner[i] = -1;
602
603 for (int i = 0, j = 0; i < comm_size; ++i)
604 for (int k = 0; k < recv_count[i]; ++k, ++j)
605 dist_vertex_owner[remote_vertex_buffer[j] - local_start_vertex] = i;
606
607 int * send_cell_of_vertex =
608 xmalloc(num_core_vertices_remote * ne * sizeof(*send_cell_of_vertex));
609
610 for (int i = 0, l = 0, m = 0; i < comm_size; ++i) {
611 for (int j = 0; j < recv_count[i]; ++j, ++l) {
612 int idx = remote_vertex_buffer[l] - local_start_vertex;
613 for (size_t k = 0; k < ne; ++k, ++m)
614 send_cell_of_vertex[m] = dist_cells_of_vertex[idx * ne + k];
615 }
616 }
617 free(dist_cells_of_vertex);
618
619 for (int i = 0; i < comm_size; ++i) {
620 send_count[i] *= ne;
621 send_displ[i] *= ne;
622 recv_count[i] *= ne;
623 recv_displ[i] *= ne;
624 }
625
626 MPI_Alltoallv(send_cell_of_vertex, recv_count, recv_displ, MPI_INT,
627 cells_of_vertex_core, send_count, send_displ, MPI_INT, comm);
628
629 free(send_cell_of_vertex);
630 free(remote_vertex_buffer);
631 free(recv_displ);
632 free(send_displ);
633 free(recv_count);
634 free(send_count);
635 }
636
637 // determine halo cells
638 int * halo_cells;
639 size_t num_halo_cells = 0;
640 {
641 int * tmp_cells = cells_of_vertex_core;
642 size_t num_tmp_cells = num_core_vertices * ne;
643
644 yac_quicksort_index(tmp_cells, num_tmp_cells, NULL);
645 yac_remove_duplicates_int(tmp_cells, &num_tmp_cells);
646
647 if ((num_tmp_cells > 0) && (tmp_cells[0] == -1)) {
648 num_tmp_cells--;
649 tmp_cells++;
650 }
651
652 halo_cells = xmalloc((num_tmp_cells - num_local_cells) * sizeof(*halo_cells));
653
654 size_t i = 0;
655 for (; i < num_tmp_cells && tmp_cells[i] < (int)local_start_cell; ++i)
656 halo_cells[num_halo_cells++] = tmp_cells[i];
657 i += num_local_cells;
658 for (; i < num_tmp_cells; ++i) halo_cells[num_halo_cells++] = tmp_cells[i];
659
660 assert(num_halo_cells == num_tmp_cells - num_local_cells);
661
662 free(cells_of_vertex_core);
663 }
664
665 // determine all vertices and get coordinates of halo cells
666 size_t num_all_local_vertices = num_halo_cells * nv + num_core_vertices;
667 int * all_cell_to_vertex;
668 int * all_local_vertices = xrealloc(core_vertices, num_all_local_vertices *
669 sizeof(*all_local_vertices));
670 if (with_cell_coords) {
671 cell_lon = xrealloc(cell_lon, (num_local_cells + num_halo_cells) *
672 sizeof(*cell_lon));
673 cell_lat = xrealloc(cell_lat, (num_local_cells + num_halo_cells) *
674 sizeof(*cell_lat));
675 }
676 if (with_cell_mask) {
677 cell_mask = xrealloc(cell_mask, (num_local_cells + num_halo_cells) *
678 sizeof(*cell_mask));
679 }
680 {
681 int * send_count = xcalloc(comm_size, sizeof(*send_count));
682 int * recv_count = xcalloc(comm_size, sizeof(*recv_count));
683
684 for (size_t i = 0; i < num_halo_cells; ++i)
685 send_count[((unsigned long)(halo_cells[i]) *
686 (unsigned long)comm_size + (unsigned long)comm_size - 1) /
687 (unsigned long)ncells]++;
688
689 MPI_Alltoall(send_count, 1, MPI_INT, recv_count, 1, MPI_INT, comm);
690
691 int * send_displ = xmalloc(comm_size * sizeof(*send_displ));
692 int * recv_displ = xmalloc(comm_size * sizeof(*recv_displ));
693 int send_accum = 0, recv_accum = 0;
694 for (int i = 0; i < comm_size; ++i) {
695 send_displ[i] = send_accum;
696 recv_displ[i] = recv_accum;
697 send_accum += send_count[i];
698 recv_accum += recv_count[i];
699 }
700
701 int num_halo_cells_remote = 0;
702 for (int i = 0; i < comm_size; ++i)
703 num_halo_cells_remote += recv_count[i];
704
705 int * remote_halo_cell_buffer =
706 xmalloc(num_halo_cells_remote * sizeof(*remote_halo_cell_buffer));
707
708 MPI_Alltoallv(halo_cells, send_count, send_displ, MPI_INT,
709 remote_halo_cell_buffer, recv_count, recv_displ, MPI_INT,
710 comm);
711
712 int * send_halo_cell_vertices =
713 xmalloc(num_halo_cells_remote * nv * sizeof(*send_halo_cell_vertices));
714 double * send_cell_lon =
715 with_cell_coords?
716 xmalloc(num_halo_cells_remote * sizeof(*send_cell_lon)):NULL;
717 double * send_cell_lat =
718 with_cell_coords?
719 xmalloc(num_halo_cells_remote * sizeof(*send_cell_lat)):NULL;
720 int * send_cell_mask =
721 with_cell_mask?
722 xmalloc(num_halo_cells_remote * sizeof(*send_cell_mask)):NULL;
723
724 for (int i = 0, l = 0, m = 0; i < comm_size; ++i) {
725 for (int j = 0; j < recv_count[i]; ++j, ++l) {
726 int idx = remote_halo_cell_buffer[l] - local_start_cell;
727 for (size_t k = 0; k < nv; ++k, ++m)
728 send_halo_cell_vertices[m] = dist_vertex_of_cell[idx * nv + k];
729 }
730 }
731 if (with_cell_coords) {
732 for (int i = 0, l = 0; i < comm_size; ++i) {
733 for (int j = 0; j < recv_count[i]; ++j, ++l) {
734 int idx = remote_halo_cell_buffer[l] - local_start_cell;
735 send_cell_lon[l] = cell_lon[idx];
736 send_cell_lat[l] = cell_lat[idx];
737 }
738 }
739 }
740 if (with_cell_mask) {
741 for (int i = 0, l = 0; i < comm_size; ++i) {
742 for (int j = 0; j < recv_count[i]; ++j, ++l) {
743 int idx = remote_halo_cell_buffer[l] - local_start_cell;
744 send_cell_mask[l] = cell_mask[idx];
745 }
746 }
747 }
748
749 if (with_cell_coords) {
750 MPI_Alltoallv(send_cell_lon, recv_count, recv_displ, MPI_DOUBLE,
751 cell_lon + num_local_cells, send_count, send_displ,
752 MPI_DOUBLE, comm);
753 MPI_Alltoallv(send_cell_lat, recv_count, recv_displ, MPI_DOUBLE,
754 cell_lat + num_local_cells, send_count, send_displ,
755 MPI_DOUBLE, comm);
756 }
757 if (with_cell_mask) {
758 MPI_Alltoallv(send_cell_mask, recv_count, recv_displ, MPI_INT,
759 cell_mask + num_local_cells, send_count, send_displ,
760 MPI_INT, comm);
761 }
762
763 for (int i = 0; i < comm_size; ++i) {
764 send_count[i] *= nv;
765 send_displ[i] *= nv;
766 recv_count[i] *= nv;
767 recv_displ[i] *= nv;
768 }
769
770 all_cell_to_vertex =
771 xrealloc(dist_vertex_of_cell, (num_local_cells + num_halo_cells) * nv *
772 sizeof(*all_cell_to_vertex));
773
774 MPI_Alltoallv(send_halo_cell_vertices, recv_count, recv_displ, MPI_INT,
775 all_cell_to_vertex + num_local_cells * nv, send_count,
776 send_displ, MPI_INT, comm);
777
778 memcpy(all_local_vertices + num_core_vertices,
779 all_cell_to_vertex + num_local_cells * nv,
780 num_halo_cells * nv * sizeof(*all_local_vertices));
781
782 free(send_cell_mask);
783 free(send_cell_lat);
784 free(send_cell_lon);
785 free(send_halo_cell_vertices);
786 free(remote_halo_cell_buffer);
787 free(recv_displ);
788 free(send_displ);
789 free(recv_count);
790 free(send_count);
791
793 all_local_vertices, num_all_local_vertices, NULL);
794 yac_remove_duplicates_int(all_local_vertices, &num_all_local_vertices);
795 }
796
797 // determine owner and coordinates for all vertices
798 double * all_vertex_lon =
799 xmalloc(num_all_local_vertices * sizeof(*all_vertex_lon));
800 double * all_vertex_lat =
801 xmalloc(num_all_local_vertices * sizeof(*all_vertex_lat));
802 int * all_local_vertices_owner =
803 xmalloc(num_all_local_vertices * sizeof(*all_local_vertices_owner));
804 {
805 int * send_count = xcalloc(comm_size, sizeof(*send_count));
806 int * recv_count = xcalloc(comm_size, sizeof(*recv_count));
807
808 for (size_t i = 0; i < num_all_local_vertices; ++i)
809 send_count[((unsigned long)(all_local_vertices[i]) *
810 (unsigned long)comm_size + (unsigned long)comm_size - 1) /
811 (unsigned long)nvertices]++;
812
813 MPI_Alltoall(send_count, 1, MPI_INT, recv_count, 1, MPI_INT, comm);
814
815 int * send_displ = xmalloc(comm_size * sizeof(*send_displ));
816 int * recv_displ = xmalloc(comm_size * sizeof(*recv_displ));
817 int send_accum = 0, recv_accum = 0;
818 for (int i = 0; i < comm_size; ++i) {
819 send_displ[i] = send_accum;
820 recv_displ[i] = recv_accum;
821 send_accum += send_count[i];
822 recv_accum += recv_count[i];
823 }
824
825 int num_all_local_vertices_remote = 0;
826 for (int i = 0; i < comm_size; ++i)
827 num_all_local_vertices_remote += recv_count[i];
828
829 int * remote_vertex_buffer =
830 xmalloc(num_all_local_vertices_remote * sizeof(*remote_vertex_buffer));
831
832 MPI_Alltoallv(all_local_vertices, send_count, send_displ, MPI_INT,
833 remote_vertex_buffer, recv_count, recv_displ, MPI_INT, comm);
834
835 int * send_vertex_owner = remote_vertex_buffer;
836 double * send_vertex_lon =
837 xmalloc(num_all_local_vertices_remote * sizeof(*send_vertex_lon));
838 double * send_vertex_lat =
839 xmalloc(num_all_local_vertices_remote * sizeof(*send_vertex_lat));
840
841 for (int i = 0, l = 0; i < comm_size; ++i) {
842 for (int j = 0; j < recv_count[i]; ++j, ++l) {
843 int idx = remote_vertex_buffer[l] - local_start_vertex;
844 send_vertex_owner[l] = dist_vertex_owner[idx];
845 send_vertex_lon[l] = vertex_lon[idx];
846 send_vertex_lat[l] = vertex_lat[idx];
847 }
848 }
849
850 free(dist_vertex_owner);
851 free(vertex_lon);
852 free(vertex_lat);
853
854 MPI_Alltoallv(send_vertex_owner, recv_count, recv_displ, MPI_INT,
855 all_local_vertices_owner, send_count, send_displ, MPI_INT,
856 comm);
857 MPI_Alltoallv(send_vertex_lon, recv_count, recv_displ, MPI_DOUBLE,
858 all_vertex_lon, send_count, send_displ, MPI_DOUBLE, comm);
859 MPI_Alltoallv(send_vertex_lat, recv_count, recv_displ, MPI_DOUBLE,
860 all_vertex_lat, send_count, send_displ, MPI_DOUBLE, comm);
861
862 free(send_vertex_lat);
863 free(send_vertex_lon);
864 free(send_vertex_owner);
865 free(recv_displ);
866 free(send_displ);
867 free(recv_count);
868 free(send_count);
869 }
870
871 // convert global ids within all_cell_to_vertex into local ids
872 {
873 size_t vertex_count = nv * (num_local_cells + num_halo_cells);
874 int * permutation = xmalloc(vertex_count * sizeof(*permutation));
875 for (size_t i = 0; i < vertex_count; ++i) permutation[i] = (int)i;
876
877 yac_quicksort_index(all_cell_to_vertex, vertex_count, permutation);
878
879 for (size_t i = 0, j = 0; i < vertex_count; ++i) {
880 while (all_local_vertices[j] != all_cell_to_vertex[i]) ++j;
881 all_cell_to_vertex[i] = (int)j;
882 }
883
884 yac_quicksort_index(permutation, vertex_count, all_cell_to_vertex);
885 free(permutation);
886 }
887
888 *nbr_vertices = num_all_local_vertices;
889 *nbr_cells = num_local_cells + num_halo_cells;
890
891 *num_vertices_per_cell =
892 xmalloc((num_local_cells + num_halo_cells) * sizeof(**num_vertices_per_cell));
893 for (size_t i = 0; i < num_local_cells + num_halo_cells; ++i)
894 (*num_vertices_per_cell)[i] = nv;
895
896 *cell_to_vertex = all_cell_to_vertex;
897
898 *global_cell_ids =
899 xmalloc((num_local_cells + num_halo_cells) * sizeof(**global_cell_ids));
900 for (size_t i = 0; i < num_local_cells; ++i)
901 (*global_cell_ids)[i] = local_start_cell + i;
902 memcpy(*global_cell_ids + num_local_cells, halo_cells,
903 num_halo_cells * sizeof(*halo_cells));
904 *global_vertex_ids = xrealloc(all_local_vertices, num_all_local_vertices *
905 sizeof(*all_local_vertices));
906
907 *cell_owner =
908 xmalloc((num_local_cells + num_halo_cells) * sizeof(**cell_owner));
909 for (size_t i = 0; i < num_local_cells; ++i)
910 (*cell_owner)[i] = -1;
911 for (size_t i = 0; i < num_halo_cells; ++i)
912 (*cell_owner)[num_local_cells + i] =
913 ((unsigned long)(halo_cells[i]) * (unsigned long)comm_size +
914 (unsigned long)comm_size - 1) / (unsigned long)ncells;
915 free(halo_cells);
916
917 for (size_t i = 0; i < num_all_local_vertices; ++i)
918 if (all_local_vertices_owner[i] == comm_rank)
919 all_local_vertices_owner[i] = -1;
920 *vertex_owner = all_local_vertices_owner;
921
922 *x_vertices = all_vertex_lon;
923 *y_vertices = all_vertex_lat;
924 if (with_cell_coords) {
925 *x_cells = cell_lon;
926 *y_cells = cell_lat;
927 }
928 if (with_cell_mask) *cell_msk = cell_mask;
929}
930
932 const char * filename, MPI_Fint comm, int * nbr_vertices, int * nbr_cells,
933 int ** num_vertices_per_cell, int ** cell_to_vertex, int ** global_cell_ids,
934 int ** cell_owner, int ** global_vertex_ids, int ** vertex_owner,
935 double ** x_vertices, double ** y_vertices,
936 double ** x_cells, double ** y_cells, int ** cell_msk) {
937
939 filename, MPI_Comm_f2c(comm), nbr_vertices, nbr_cells,
940 num_vertices_per_cell, cell_to_vertex, global_cell_ids,
941 cell_owner, global_vertex_ids, vertex_owner,
942 x_vertices, y_vertices, x_cells, y_cells, cell_msk);
943}
944
945static int * generate_simple_core_mask(size_t N) {
946 int * mask = xmalloc(N * sizeof(*mask));
947 for (size_t i = 0; i < N; ++i) mask[i] = 1;
948 return mask;
949}
950
951static size_t * generate_offsets(size_t N, int * counts) {
952
953 size_t * offsets = xmalloc(N * sizeof(*offsets));
954 for (size_t i = 0, accu = 0; i < N; ++i) {
955 offsets[i] = accu;
956 accu += (size_t)(counts[i]);
957 }
958 return offsets;
959}
960
962 const char * filename, MPI_Comm comm,
963 double ** x_vertices, double ** y_vertices,
964 yac_int ** cell_ids, yac_int ** vertex_ids, yac_int ** edge_ids,
965 size_t * num_cells, size_t * num_vertices, size_t * num_edges,
966 int ** num_vertices_per_cell, int ** num_cells_per_vertex,
967 size_t ** cell_to_vertex, size_t ** cell_to_edge, size_t ** vertex_to_cell,
968 size_t ** edge_to_vertex, enum yac_edge_type ** edge_type,
969 double ** x_cells, double ** y_cells, int ** cell_msk) {
970
972 ((x_cells == NULL) && (y_cells == NULL)) ||
973 ((x_cells != NULL) && (y_cells != NULL)),
974 "ERROR(yac_read_icon_grid_information_parallel_2): "
975 "arguments x_cells and y_cells have to be either both NULL or "
976 "both different from NULL");
977
978 const int with_cell_coords = x_cells != NULL;
979 const int with_cell_mask = cell_msk != NULL;
980
981 int comm_rank, comm_size;
982
983 MPI_Comm_rank(comm, &comm_rank);
984 MPI_Comm_size(comm, &comm_size);
985
986 int local_is_io, * io_ranks, num_io_ranks;
987 yac_get_io_ranks(comm, &local_is_io, &io_ranks, &num_io_ranks);
988
989 size_t ncells, nvertices, nedges, nv, ne;
990
991 size_t read_local_start_cell = 0;
992 size_t read_num_local_cells = 0;
993 size_t read_local_start_vertex = 0;
994 size_t read_num_local_vertices = 0;
995 size_t read_local_start_edge = 0;
996 size_t read_num_local_edges = 0;
997 double * read_cell_lon = NULL;
998 double * read_cell_lat = NULL;
999 int * read_cell_mask = NULL;
1000 double * read_vertex_lon = NULL;
1001 double * read_vertex_lat = NULL;
1002 int * read_dist_vertex_of_cell = NULL;
1003 int * read_dist_edge_of_cell = NULL;
1004 int * read_dist_edge_vertices = NULL;
1005
1006 if (local_is_io) {
1007
1008 unsigned long io_proc_idx = ULONG_MAX;
1009 for (int i = 0; (i < num_io_ranks) && (io_proc_idx == ULONG_MAX); ++i)
1010 if (io_ranks[i] == comm_rank)
1011 io_proc_idx = (unsigned long)i;
1012
1013 // open file
1014 int ncid;
1015 yac_nc_open(filename, NC_NOWRITE, &ncid);
1016
1017 // get number of cells and vertices
1018 int dim_id;
1019 yac_nc_inq_dimid(ncid, "cell", &dim_id);
1020 YAC_HANDLE_ERROR(nc_inq_dimlen(ncid, dim_id, &ncells));
1021 yac_nc_inq_dimid(ncid, "vertex", &dim_id);
1022 YAC_HANDLE_ERROR(nc_inq_dimlen(ncid, dim_id, &nvertices));
1023 yac_nc_inq_dimid(ncid, "edge", &dim_id);
1024 YAC_HANDLE_ERROR(nc_inq_dimlen(ncid, dim_id, &nedges));
1025 yac_nc_inq_dimid(ncid, "nv", &dim_id);
1026 YAC_HANDLE_ERROR(nc_inq_dimlen(ncid, dim_id, &nv));
1027 yac_nc_inq_dimid(ncid, "ne", &dim_id);
1028 YAC_HANDLE_ERROR(nc_inq_dimlen(ncid, dim_id, &ne));
1029
1030 // determine local range for cell and vertex data
1031 read_local_start_cell =
1032 ((unsigned long)ncells * io_proc_idx) / (unsigned long)num_io_ranks;
1033 read_num_local_cells =
1034 ((unsigned long)ncells * (io_proc_idx+1)) / (unsigned long)num_io_ranks -
1035 (unsigned long)read_local_start_cell;
1036 read_local_start_vertex =
1037 ((unsigned long)nvertices * io_proc_idx) / (unsigned long)num_io_ranks;
1038 read_num_local_vertices =
1039 ((unsigned long)nvertices * (io_proc_idx+1)) / (unsigned long)num_io_ranks -
1040 (unsigned long)read_local_start_vertex;
1041 read_local_start_edge =
1042 ((unsigned long)nedges * io_proc_idx) / (unsigned long)num_io_ranks;
1043 read_num_local_edges =
1044 ((unsigned long)nedges * (io_proc_idx+1)) / (unsigned long)num_io_ranks -
1045 (unsigned long)read_local_start_edge;
1046
1047 // read basic grid data (each process its individual part)
1048 read_vertex_lon = xmalloc(read_num_local_vertices * sizeof(*read_vertex_lon));
1049 read_vertex_lat = xmalloc(read_num_local_vertices * sizeof(*read_vertex_lat));
1050 read_cell_lon =
1051 with_cell_coords?
1052 xmalloc(read_num_local_cells * sizeof(*read_cell_lon)):NULL;
1053 read_cell_lat =
1054 with_cell_coords?
1055 xmalloc(read_num_local_cells * sizeof(*read_cell_lat)):NULL;
1056 read_cell_mask =
1057 with_cell_mask?
1058 xmalloc(read_num_local_cells * sizeof(*read_cell_mask)):NULL;
1059 read_dist_vertex_of_cell =
1060 xmalloc(read_num_local_cells * nv * sizeof(*read_dist_vertex_of_cell));
1061 read_dist_edge_of_cell =
1062 xmalloc(read_num_local_cells * nv * sizeof(*read_dist_edge_of_cell));
1063 read_dist_edge_vertices =
1064 xmalloc(read_num_local_edges * 2 * sizeof(*read_dist_edge_vertices));
1065 int varid;
1066 yac_nc_inq_varid(ncid, "vlon", &varid);
1067 YAC_HANDLE_ERROR(nc_get_vara_double(ncid, varid, &read_local_start_vertex,
1068 &read_num_local_vertices, read_vertex_lon));
1069 convert_to_rad(ncid, varid, read_vertex_lon, read_num_local_vertices);
1070 yac_nc_inq_varid(ncid, "vlat", &varid);
1071 YAC_HANDLE_ERROR(nc_get_vara_double(ncid, varid, &read_local_start_vertex,
1072 &read_num_local_vertices, read_vertex_lat));
1073 convert_to_rad(ncid, varid, read_vertex_lat, read_num_local_vertices);
1074 if (with_cell_coords) {
1075 yac_nc_inq_varid(ncid, "clon", &varid);
1077 nc_get_vara_double(
1078 ncid, varid, &read_local_start_cell,
1079 &read_num_local_cells, read_cell_lon));
1080 convert_to_rad(ncid, varid, read_cell_lon, read_num_local_cells);
1081 yac_nc_inq_varid(ncid, "clat", &varid);
1083 nc_get_vara_double(
1084 ncid, varid, &read_local_start_cell,
1085 &read_num_local_cells, read_cell_lat));
1086 convert_to_rad(ncid, varid, read_cell_lat, read_num_local_cells);
1087 }
1088 if (read_cell_mask) {
1089 yac_nc_inq_varid(ncid, "cell_sea_land_mask", &varid);
1091 nc_get_vara_int(
1092 ncid, varid, &read_local_start_cell,
1093 &read_num_local_cells, read_cell_mask));
1094 }
1095 {
1096 int * buffer =
1097 xmalloc(
1098 MAX(nv*read_num_local_cells, 2*read_num_local_edges) *
1099 sizeof(*buffer));
1100 {
1101 size_t tmp_start[2] = {0, read_local_start_cell};
1102 size_t tmp_count[2] = {nv, read_num_local_cells};
1103 yac_nc_inq_varid(ncid, "vertex_of_cell", &varid);
1104 YAC_HANDLE_ERROR(nc_get_vara_int(ncid, varid, tmp_start, tmp_count, buffer));
1105 for (size_t i = 0; i < read_num_local_cells; ++i)
1106 for (size_t j = 0; j < nv; ++j)
1107 read_dist_vertex_of_cell[i * nv + j] =
1108 buffer[i + j * read_num_local_cells];
1109 yac_nc_inq_varid(ncid, "edge_of_cell", &varid);
1110 YAC_HANDLE_ERROR(nc_get_vara_int(ncid, varid, tmp_start, tmp_count, buffer));
1111 for (size_t i = 0; i < read_num_local_cells; ++i)
1112 for (size_t j = 0; j < nv; ++j)
1113 read_dist_edge_of_cell[i * nv + j] =
1114 buffer[i + j * read_num_local_cells];
1115 }
1116 {
1117 size_t tmp_start[2] = {0, read_local_start_edge};
1118 size_t tmp_count[2] = {2, read_num_local_edges};
1119 yac_nc_inq_varid(ncid, "edge_vertices", &varid);
1120 YAC_HANDLE_ERROR(nc_get_vara_int(ncid, varid, tmp_start, tmp_count, buffer));
1121 for (size_t i = 0; i < read_num_local_edges; ++i)
1122 for (int j = 0; j < 2; ++j)
1123 read_dist_edge_vertices[i * 2 + j] =
1124 buffer[i + j * read_num_local_edges];
1125 }
1126 free(buffer);
1127
1128 // adjust for c indices
1129 for (size_t i = 0; i < read_num_local_cells * nv; ++i)
1130 if (read_dist_vertex_of_cell[i] > 0) read_dist_vertex_of_cell[i]--;
1131 for (size_t i = 0; i < read_num_local_cells * nv; ++i)
1132 if (read_dist_edge_of_cell[i] > 0) read_dist_edge_of_cell[i]--;
1133 for (size_t i = 0; i < read_num_local_edges * 2; ++i)
1134 if (read_dist_edge_vertices[i] > 0) read_dist_edge_vertices[i]--;
1135 }
1136
1137 YAC_HANDLE_ERROR(nc_close(ncid));
1138 } else {
1139 read_cell_lon = with_cell_coords?xmalloc(1 * sizeof(*read_cell_lon)):NULL;
1140 read_cell_lat = with_cell_coords?xmalloc(1 * sizeof(*read_cell_lat)):NULL;
1141 read_cell_mask = with_cell_mask?xmalloc(1 * sizeof(*read_cell_mask)):NULL;
1142 read_vertex_lon = xmalloc(1 * sizeof(*read_vertex_lon));
1143 read_vertex_lat = xmalloc(1 * sizeof(*read_vertex_lat));
1144 read_dist_vertex_of_cell = xmalloc(1 * sizeof(*read_dist_vertex_of_cell));
1145 read_dist_edge_of_cell = xmalloc(1 * sizeof(*read_dist_edge_of_cell));
1146 read_dist_edge_vertices = xmalloc(1 * sizeof(*read_dist_edge_vertices));
1147 }
1148
1149 free(io_ranks);
1150
1151 {
1152 int tmp;
1153 if (comm_rank == 0) tmp = (int)ncells;
1154 MPI_Bcast(&tmp, 1, MPI_INT, 0, comm);
1155 ncells = (size_t)tmp;
1156 if (comm_rank == 0) tmp = (int)nvertices;
1157 MPI_Bcast(&tmp, 1, MPI_INT, 0, comm);
1158 nvertices = (size_t)tmp;
1159 if (comm_rank == 0) tmp = (int)nedges;
1160 MPI_Bcast(&tmp, 1, MPI_INT, 0, comm);
1161 nedges = (size_t)tmp;
1162 if (comm_rank == 0) tmp = (int)nv;
1163 MPI_Bcast(&tmp, 1, MPI_INT, 0, comm);
1164 nv = (size_t)tmp;
1165 if (comm_rank == 0) tmp = (int)ne;
1166 MPI_Bcast(&tmp, 1, MPI_INT, 0, comm);
1167 ne = (size_t)tmp;
1168 }
1169
1170 // determine local range for cell and vertex data
1171 size_t local_start_cell =
1172 ((unsigned long)ncells * (unsigned long)comm_rank) /
1173 (unsigned long)comm_size;
1174 size_t num_local_cells =
1175 ((unsigned long)ncells * ((unsigned long)comm_rank+1)) /
1176 (unsigned long)comm_size - (unsigned long)local_start_cell;
1177 size_t local_start_vertex =
1178 ((unsigned long)nvertices * (unsigned long)comm_rank) /
1179 (unsigned long)comm_size;
1180 size_t num_local_vertices =
1181 ((unsigned long)nvertices * ((unsigned long)comm_rank+1)) /
1182 (unsigned long)comm_size - (unsigned long)local_start_vertex;
1183 size_t local_start_edge =
1184 ((unsigned long)nedges * (unsigned long)comm_rank) /
1185 (unsigned long)comm_size;
1186 size_t num_local_edges =
1187 ((unsigned long)nedges * ((unsigned long)comm_rank+1)) /
1188 (unsigned long)comm_size - (unsigned long)local_start_edge;
1189
1190 // redistribute basic cell data (from io decomposition)
1191 if (with_cell_coords) {
1192 *x_cells = xmalloc(num_local_cells * sizeof(**x_cells));
1193 *y_cells = xmalloc(num_local_cells * sizeof(**y_cells));
1194 }
1195 if (with_cell_mask) {
1196 *cell_msk = xmalloc(num_local_cells * sizeof(**cell_msk));
1197 }
1198 int * dist_vertex_of_cell =
1199 xmalloc(num_local_cells * nv * sizeof(*dist_vertex_of_cell));
1200 int * dist_edge_of_cell =
1201 xmalloc(num_local_cells * nv * sizeof(*dist_edge_of_cell));
1202 {
1203 int * send_count = xcalloc(comm_size, sizeof(*send_count));
1204 int * recv_count = xcalloc(comm_size, sizeof(*recv_count));
1205
1206 for (unsigned i = 0; i < read_num_local_cells; ++i)
1207 send_count[
1209 read_local_start_cell + i, ncells, comm_size)]++;
1210
1211 MPI_Alltoall(send_count, 1, MPI_INT, recv_count, 1, MPI_INT, comm);
1212
1213 int * send_displ = xmalloc(comm_size * sizeof(*send_displ));
1214 int * recv_displ = xmalloc(comm_size * sizeof(*recv_displ));
1215 int send_accum = 0, recv_accum = 0;
1216 for (int i = 0; i < comm_size; ++i) {
1217 send_displ[i] = send_accum;
1218 recv_displ[i] = recv_accum;
1219 send_accum += send_count[i];
1220 recv_accum += recv_count[i];
1221 }
1222
1223 if (with_cell_coords) {
1224 MPI_Alltoallv(read_cell_lon, send_count, send_displ, MPI_DOUBLE,
1225 *x_cells, recv_count, recv_displ, MPI_DOUBLE, comm);
1226 MPI_Alltoallv(read_cell_lat, send_count, send_displ, MPI_DOUBLE,
1227 *y_cells, recv_count, recv_displ, MPI_DOUBLE, comm);
1228 }
1229 if (with_cell_mask) {
1230 MPI_Alltoallv(read_cell_mask, send_count, send_displ, MPI_INT,
1231 *cell_msk, recv_count, recv_displ, MPI_INT, comm);
1232 }
1233
1234 for (int i = 0; i < comm_size; ++i) {
1235 send_count[i] *= nv;
1236 send_displ[i] *= nv;
1237 recv_count[i] *= nv;
1238 recv_displ[i] *= nv;
1239 }
1240
1241 MPI_Alltoallv(read_dist_vertex_of_cell, send_count, send_displ, MPI_INT,
1242 dist_vertex_of_cell, recv_count, recv_displ, MPI_INT, comm);
1243 MPI_Alltoallv(read_dist_edge_of_cell, send_count, send_displ, MPI_INT,
1244 dist_edge_of_cell, recv_count, recv_displ, MPI_INT, comm);
1245
1246 free(recv_displ);
1247 free(send_displ);
1248 free(recv_count);
1249 free(send_count);
1250 free(read_dist_vertex_of_cell);
1251 free(read_dist_edge_of_cell);
1252 free(read_cell_mask);
1253 free(read_cell_lat);
1254 free(read_cell_lon);
1255 }
1256
1257 // redistribute basic vertex data (from io decomposition)
1258 double * vertex_lon = xmalloc(num_local_vertices * sizeof(*vertex_lon));
1259 double * vertex_lat = xmalloc(num_local_vertices * sizeof(*vertex_lat));
1260 {
1261 int * send_count = xcalloc(comm_size, sizeof(*send_count));
1262 int * recv_count = xcalloc(comm_size, sizeof(*recv_count));
1263
1264 for (unsigned i = 0; i < read_num_local_vertices; ++i)
1265 send_count[
1267 read_local_start_vertex + i, nvertices, comm_size)]++;
1268
1269 MPI_Alltoall(send_count, 1, MPI_INT, recv_count, 1, MPI_INT, comm);
1270
1271 int * send_displ = xmalloc(comm_size * sizeof(*send_displ));
1272 int * recv_displ = xmalloc(comm_size * sizeof(*recv_displ));
1273 int send_accum = 0, recv_accum = 0;
1274 for (int i = 0; i < comm_size; ++i) {
1275 send_displ[i] = send_accum;
1276 recv_displ[i] = recv_accum;
1277 send_accum += send_count[i];
1278 recv_accum += recv_count[i];
1279 }
1280
1281 MPI_Alltoallv(read_vertex_lon, send_count, send_displ, MPI_DOUBLE,
1282 vertex_lon, recv_count, recv_displ, MPI_DOUBLE, comm);
1283 MPI_Alltoallv(read_vertex_lat, send_count, send_displ, MPI_DOUBLE,
1284 vertex_lat, recv_count, recv_displ, MPI_DOUBLE, comm);
1285
1286 free(recv_displ);
1287 free(send_displ);
1288 free(recv_count);
1289 free(send_count);
1290 free(read_vertex_lon);
1291 free(read_vertex_lat);
1292 }
1293
1294 // redistribute basic edge data (from io decomposition)
1295 int * dist_edge_vertices =
1296 xmalloc(num_local_edges * 2 * sizeof(*dist_edge_vertices));
1297 {
1298 int * send_count = xcalloc(comm_size, sizeof(*send_count));
1299 int * recv_count = xcalloc(comm_size, sizeof(*recv_count));
1300
1301 for (unsigned i = 0; i < read_num_local_edges; ++i)
1302 send_count[
1304 read_local_start_edge + i, nedges, comm_size)] += 2;
1305
1306 MPI_Alltoall(send_count, 1, MPI_INT, recv_count, 1, MPI_INT, comm);
1307
1308 int * send_displ = xmalloc(comm_size * sizeof(*send_displ));
1309 int * recv_displ = xmalloc(comm_size * sizeof(*recv_displ));
1310 int send_accum = 0, recv_accum = 0;
1311 for (int i = 0; i < comm_size; ++i) {
1312 send_displ[i] = send_accum;
1313 recv_displ[i] = recv_accum;
1314 send_accum += send_count[i];
1315 recv_accum += recv_count[i];
1316 }
1317
1318 MPI_Alltoallv(read_dist_edge_vertices, send_count, send_displ, MPI_INT,
1319 dist_edge_vertices, recv_count, recv_displ, MPI_INT, comm);
1320
1321 free(recv_displ);
1322 free(send_displ);
1323 free(recv_count);
1324 free(send_count);
1325 free(read_dist_edge_vertices);
1326 }
1327
1328 // determine required vertices for core cells
1329 // in additional compute num_cells_per_vertex, vertex_to_cell,
1330 // and cell_to_vertex
1331 size_t num_core_vertices;
1332 {
1333 size_t N = num_local_cells * nv;
1334 *vertex_ids = xmalloc(N * sizeof(**vertex_ids));
1335 *num_cells_per_vertex = xmalloc(N * sizeof(**num_cells_per_vertex));
1336 *vertex_to_cell = xmalloc(N * sizeof(**vertex_to_cell));
1337 *cell_to_vertex = xmalloc(N * sizeof(**cell_to_vertex));
1338 for (size_t i = 0; i < N; ++i)
1339 (*vertex_ids)[i] = (yac_int)dist_vertex_of_cell[i];
1340 size_t * permutation = *vertex_to_cell;
1341 for (size_t i = 0; i < N; ++i) permutation[i] = i;
1342 yac_quicksort_index_yac_int_size_t(*vertex_ids, N, permutation);
1343 // remove duplicated core vertices and count number of cells per vertex
1344 yac_int prev_vertex_id = YAC_INT_MAX;
1345 num_core_vertices = 0;
1346 for (size_t i = 0; i < N; ++i) {
1347 yac_int curr_vertex_id = (*vertex_ids)[i];
1348 if (prev_vertex_id == curr_vertex_id) {
1349 (*num_cells_per_vertex)[num_core_vertices-1]++;
1350 } else {
1351 (*num_cells_per_vertex)[num_core_vertices] = 1;
1352 (*vertex_ids)[num_core_vertices] = (prev_vertex_id = curr_vertex_id);
1353 ++num_core_vertices;
1354 }
1355 (*cell_to_vertex)[permutation[i]] = num_core_vertices-1;
1356 permutation[i] /= nv;
1357 }
1358 *vertex_ids =
1359 xrealloc(*vertex_ids, num_core_vertices * sizeof(**vertex_ids));
1360 *num_cells_per_vertex =
1361 xrealloc(*num_cells_per_vertex,
1362 num_core_vertices * sizeof(**num_cells_per_vertex));
1363 free(dist_vertex_of_cell);
1364 }
1365
1366 // determine required edges for core cells
1367 // in additional compute edge_to_cell and cell_to_edge
1368 size_t num_core_edges;
1369 {
1370 size_t N = num_local_cells * nv;
1371 *edge_ids = xmalloc(N * sizeof(**edge_ids));
1372 size_t * permutation = xmalloc(N * sizeof(*permutation));
1373 *cell_to_edge = xmalloc(N * sizeof(**cell_to_edge));
1374 for (size_t i = 0; i < N; ++i)
1375 (*edge_ids)[i] = (yac_int)dist_edge_of_cell[i];
1376 for (size_t i = 0; i < N; ++i) permutation[i] = i;
1377 yac_quicksort_index_yac_int_size_t(*edge_ids, N, permutation);
1378 // remove duplicated core edges and count number of cells per edge
1379 yac_int prev_edge_id = YAC_INT_MAX;
1380 num_core_edges = 0;
1381 for (size_t i = 0; i < N; ++i) {
1382 yac_int curr_edge_id = (*edge_ids)[i];
1383 if (prev_edge_id != curr_edge_id)
1384 (*edge_ids)[num_core_edges++] = (prev_edge_id = curr_edge_id);
1385 (*cell_to_edge)[permutation[i]] = num_core_edges-1;
1386 permutation[i] /= nv;
1387 }
1388 free(permutation);
1389 *edge_ids =
1390 xrealloc(*edge_ids, num_core_edges * sizeof(**edge_ids));
1391 free(dist_edge_of_cell);
1392 }
1393
1394 // generate vertex coordinate data
1395 {
1396 *x_vertices = xmalloc(num_core_vertices * sizeof(**x_vertices));
1397 *y_vertices = xmalloc(num_core_vertices * sizeof(**y_vertices));
1398 int * send_count = xcalloc(comm_size, sizeof(*send_count));
1399 int * recv_count = xcalloc(comm_size, sizeof(*recv_count));
1400
1401 for (size_t i = 0; i < num_core_vertices; ++i)
1402 send_count[((unsigned long)((*vertex_ids)[i]) *
1403 (unsigned long)comm_size + (unsigned long)comm_size - 1) /
1404 (unsigned long)nvertices]++;
1405
1406 MPI_Alltoall(send_count, 1, MPI_INT, recv_count, 1, MPI_INT, comm);
1407
1408 int * send_displ = xmalloc(comm_size * sizeof(*send_displ));
1409 int * recv_displ = xmalloc(comm_size * sizeof(*recv_displ));
1410 int send_accum = 0, recv_accum = 0;
1411 for (int i = 0; i < comm_size; ++i) {
1412 send_displ[i] = send_accum;
1413 recv_displ[i] = recv_accum;
1414 send_accum += send_count[i];
1415 recv_accum += recv_count[i];
1416 }
1417
1418 int num_all_local_vertices_remote = 0;
1419 for (int i = 0; i < comm_size; ++i)
1420 num_all_local_vertices_remote += recv_count[i];
1421
1422 yac_int * remote_vertex_buffer =
1423 xmalloc(num_all_local_vertices_remote * sizeof(*remote_vertex_buffer));
1424
1425 MPI_Alltoallv(
1426 *vertex_ids, send_count, send_displ, yac_int_dt,
1427 remote_vertex_buffer, recv_count, recv_displ, yac_int_dt, comm);
1428
1429 double * send_vertex_lon =
1430 xmalloc(num_all_local_vertices_remote * sizeof(*send_vertex_lon));
1431 double * send_vertex_lat =
1432 xmalloc(num_all_local_vertices_remote * sizeof(*send_vertex_lat));
1433
1434 for (int i = 0, l = 0; i < comm_size; ++i) {
1435 for (int j = 0; j < recv_count[i]; ++j, ++l) {
1436 size_t idx = (size_t)(remote_vertex_buffer[l]) - local_start_vertex;
1437 send_vertex_lon[l] = vertex_lon[idx];
1438 send_vertex_lat[l] = vertex_lat[idx];
1439 }
1440 }
1441
1442 free(remote_vertex_buffer);
1443 free(vertex_lon);
1444 free(vertex_lat);
1445
1446 MPI_Alltoallv(send_vertex_lon, recv_count, recv_displ, MPI_DOUBLE,
1447 *x_vertices, send_count, send_displ, MPI_DOUBLE, comm);
1448 MPI_Alltoallv(send_vertex_lat, recv_count, recv_displ, MPI_DOUBLE,
1449 *y_vertices, send_count, send_displ, MPI_DOUBLE, comm);
1450
1451 free(send_vertex_lat);
1452 free(send_vertex_lon);
1453 free(recv_displ);
1454 free(send_displ);
1455 free(recv_count);
1456 free(send_count);
1457 }
1458
1459 // generate edge vertex data
1460 *edge_to_vertex =
1461 xmalloc(2 * num_core_edges * sizeof(**edge_to_vertex));
1462 {
1463 int * local_edge_to_vertex =
1464 xmalloc(2 * num_core_edges * sizeof(*local_edge_to_vertex));
1465 int * send_count = xcalloc(comm_size, sizeof(*send_count));
1466 int * recv_count = xcalloc(comm_size, sizeof(*recv_count));
1467
1468 for (size_t i = 0; i < num_core_edges; ++i)
1469 send_count[((unsigned long)((*edge_ids)[i]) *
1470 (unsigned long)comm_size + (unsigned long)comm_size - 1) /
1471 (unsigned long)nedges]++;
1472
1473 MPI_Alltoall(send_count, 1, MPI_INT, recv_count, 1, MPI_INT, comm);
1474
1475 int * send_displ = xmalloc(comm_size * sizeof(*send_displ));
1476 int * recv_displ = xmalloc(comm_size * sizeof(*recv_displ));
1477 int send_accum = 0, recv_accum = 0;
1478 for (int i = 0; i < comm_size; ++i) {
1479 send_displ[i] = send_accum;
1480 recv_displ[i] = recv_accum;
1481 send_accum += send_count[i];
1482 recv_accum += recv_count[i];
1483 }
1484
1485 int num_all_local_edges_remote = 0;
1486 for (int i = 0; i < comm_size; ++i)
1487 num_all_local_edges_remote += recv_count[i];
1488
1489 yac_int * remote_edge_buffer =
1490 xmalloc(num_all_local_edges_remote * sizeof(*remote_edge_buffer));
1491
1492 MPI_Alltoallv(
1493 *edge_ids, send_count, send_displ, yac_int_dt,
1494 remote_edge_buffer, recv_count, recv_displ, yac_int_dt, comm);
1495
1496 int * send_edge_vertices =
1497 xmalloc(2 * num_all_local_edges_remote * sizeof(*send_edge_vertices));
1498
1499 for (int i = 0, l = 0; i < comm_size; ++i) {
1500 for (int j = 0; j < recv_count[i]; ++j, ++l) {
1501 size_t idx = (size_t)(remote_edge_buffer[l]) - local_start_edge;
1502 send_edge_vertices[2*l+0] = dist_edge_vertices[2*idx+0];
1503 send_edge_vertices[2*l+1] = dist_edge_vertices[2*idx+1];
1504
1505 }
1506 send_count[i] *= 2;
1507 recv_count[i] *= 2;
1508 send_displ[i] *= 2;
1509 recv_displ[i] *= 2;
1510 }
1511
1512 free(remote_edge_buffer);
1513 free(dist_edge_vertices);
1514
1515 MPI_Alltoallv(send_edge_vertices, recv_count, recv_displ, MPI_INT,
1516 local_edge_to_vertex, send_count, send_displ, MPI_INT, comm);
1517
1518 size_t * permutation = xmalloc(2 * num_core_edges * sizeof(*permutation));
1519 for (size_t i = 0; i < 2 * num_core_edges; ++i) permutation[i] = i;
1520
1522 local_edge_to_vertex, 2 * num_core_edges, permutation);
1523
1524 for (size_t i = 0, j = 0; i < 2 * num_core_edges; ++i) {
1525 yac_int curr_vertex = (yac_int)(local_edge_to_vertex[i]);
1526 while ((j < nvertices) && ((*vertex_ids)[j] < curr_vertex)) ++j;
1527 YAC_ASSERT(
1528 (j < nvertices) && ((*vertex_ids)[j] == curr_vertex),
1529 "ERROR(yac_read_icon_basic_grid_data_parallel): vertex id missmatch")
1530 (*edge_to_vertex)[permutation[i]] = j;
1531 }
1532
1533 free(permutation);
1534 free(send_edge_vertices);
1535 free(recv_displ);
1536 free(send_displ);
1537 free(recv_count);
1538 free(send_count);
1539 free(local_edge_to_vertex);
1540 }
1541
1542 // generate cell ids for local partition
1543 *cell_ids = xmalloc(num_local_cells * sizeof(**cell_ids));
1544 for (size_t i = 0; i < num_local_cells; ++i)
1545 (*cell_ids)[i] = (yac_int)(local_start_cell + i);
1546
1547 // generate num_vertices_per_cell
1548 *num_vertices_per_cell =
1549 xmalloc(num_local_cells * sizeof(**num_vertices_per_cell));
1550 for (size_t i = 0; i < num_local_cells; ++i) (*num_vertices_per_cell)[i] = nv;
1551
1552 // generate edge_type (for icon grids this is always YAC_GREAT_CIRCLE_EDGE)
1553 *edge_type = xmalloc(num_core_edges * sizeof(**edge_type));
1554 for (size_t i = 0; i < num_core_edges; ++i)
1555 (*edge_type)[i] = YAC_GREAT_CIRCLE_EDGE;
1556
1557 *num_cells = num_local_cells;
1558 *num_vertices = num_core_vertices;
1559 *num_edges = num_core_edges;
1560}
1561
1563 const char * filename, MPI_Comm comm) {
1564
1565 double * x_vertices;
1566 double * y_vertices;
1567 yac_int * cell_ids;
1569 yac_int * edge_ids;
1570 size_t num_cells;
1571 size_t num_vertices;
1572 size_t num_edges;
1575 size_t * cell_to_vertex;
1576 size_t * cell_to_edge;
1577 size_t * vertex_to_cell;
1578 size_t * edge_to_vertex;
1579 enum yac_edge_type * edge_type;
1580
1582 filename, comm, &x_vertices, &y_vertices, &cell_ids, &vertex_ids, &edge_ids,
1585 &edge_to_vertex, &edge_type, NULL, NULL, NULL);
1586
1589 for (size_t i = 0; i < num_vertices; ++i)
1590 LLtoXYZ(x_vertices[i], y_vertices[i], vertex_coordinates[i]);
1591 free(x_vertices);
1592 free(y_vertices);
1593
1595 grid_data.vertex_coordinates = vertex_coordinates;
1596 grid_data.cell_ids = cell_ids;
1597 grid_data.vertex_ids = vertex_ids;
1598 grid_data.edge_ids = edge_ids;
1599 grid_data.num_cells = num_cells;
1600 grid_data.num_vertices = num_vertices;
1601 grid_data.num_edges = num_edges;
1606 grid_data.num_cells_per_vertex = num_cells_per_vertex;
1608 grid_data.cell_to_vertex_offsets = generate_offsets(num_cells, num_vertices_per_cell);
1609 grid_data.cell_to_edge = cell_to_edge;
1610 grid_data.cell_to_edge_offsets = grid_data.cell_to_vertex_offsets;
1611 grid_data.vertex_to_cell = vertex_to_cell;
1613 grid_data.edge_to_vertex = (yac_size_t_2_pointer)&(edge_to_vertex[0]);
1614 grid_data.edge_type = edge_type;
1615 grid_data.num_total_cells = num_cells;
1616 grid_data.num_total_vertices = num_vertices;
1617 grid_data.num_total_edges = num_edges;
1618
1619 return grid_data;
1620}
1621
1623 char const * filename, char const * gridname, MPI_Comm comm,
1624 struct yac_basic_grid ** basic_grid, size_t * cell_coordinate_idx,
1625 int ** cell_mask) {
1626
1627 double * x_vertices;
1628 double * y_vertices;
1629 yac_int * cell_ids;
1631 yac_int * edge_ids;
1632 size_t num_cells;
1633 size_t num_vertices;
1634 size_t num_edges;
1637 size_t * cell_to_vertex;
1638 size_t * cell_to_edge;
1639 size_t * vertex_to_cell;
1640 size_t * edge_to_vertex;
1641 enum yac_edge_type * edge_type;
1642 double * x_cells;
1643 double * y_cells;
1644
1646 filename, comm, &x_vertices, &y_vertices, &cell_ids, &vertex_ids, &edge_ids,
1649 &edge_to_vertex, &edge_type, &x_cells, &y_cells, cell_mask);
1650
1653 for (size_t i = 0; i < num_vertices; ++i)
1654 LLtoXYZ(x_vertices[i], y_vertices[i], vertex_coordinates[i]);
1655 free(x_vertices);
1656 free(y_vertices);
1657
1658 struct yac_basic_grid_data basic_grid_data;
1659 basic_grid_data.vertex_coordinates = vertex_coordinates;
1660 basic_grid_data.cell_ids = cell_ids;
1661 basic_grid_data.vertex_ids = vertex_ids;
1662 basic_grid_data.edge_ids = edge_ids;
1663 basic_grid_data.num_cells = num_cells;
1664 basic_grid_data.num_vertices = num_vertices;
1665 basic_grid_data.num_edges = num_edges;
1671 basic_grid_data.cell_to_vertex = cell_to_vertex;
1673 basic_grid_data.cell_to_edge = cell_to_edge;
1674 basic_grid_data.cell_to_edge_offsets = basic_grid_data.cell_to_vertex_offsets;
1675 basic_grid_data.vertex_to_cell = vertex_to_cell;
1677 basic_grid_data.edge_to_vertex = (yac_size_t_2_pointer)&(edge_to_vertex[0]);
1678 basic_grid_data.edge_type = edge_type;
1679 basic_grid_data.num_total_cells = num_cells;
1680 basic_grid_data.num_total_vertices = num_vertices;
1681 basic_grid_data.num_total_edges = num_edges;
1682
1683 *basic_grid = yac_basic_grid_new(gridname, basic_grid_data);
1684
1685 yac_coordinate_pointer cell_coordinates =
1686 xmalloc(num_cells * sizeof(*cell_coordinates));
1687 for (size_t i = 0; i < num_cells; ++i)
1688 LLtoXYZ(x_cells[i], y_cells[i], cell_coordinates[i]);
1689 *cell_coordinate_idx =
1691 *basic_grid, YAC_LOC_CELL, cell_coordinates);
1692 free(x_cells);
1693 free(y_cells);
1694}
1695
1697 char const * filename, char const * gridname, MPI_Comm comm) {
1698
1699 return
1701 gridname,
1703}
1704
1706 char const * filename, char const * gridname, MPI_Fint comm) {
1707
1708 return
1710 filename, gridname, MPI_Comm_f2c(comm));
1711}
1712
1714 char const * filename) {
1715
1716 int nbr_vertices;
1717 int nbr_cells;
1718 int * num_vertices_per_cell = NULL;
1719 int * cell_to_vertex = NULL;
1720 int * cell_mask = NULL;
1721
1722 double * x_vertices = NULL;
1723 double * y_vertices = NULL;
1724 double * x_cells = NULL;
1725 double * y_cells = NULL;
1726
1727 yac_read_icon_grid_information(filename, &nbr_vertices, &nbr_cells,
1729 &x_vertices, &y_vertices,
1730 &x_cells, &y_cells,
1731 &cell_mask);
1732
1733 free(x_cells);
1734 free(y_cells);
1735 free(cell_mask);
1736
1737 struct yac_basic_grid_data grid =
1739 (size_t)nbr_vertices, (int)nbr_cells, num_vertices_per_cell,
1740 x_vertices, y_vertices, cell_to_vertex);
1742 free(x_vertices);
1743 free(y_vertices);
1744 free(cell_to_vertex);
1745
1746 return grid;
1747}
1748
1750 char const * filename, char const * gridname) {
1751
1752 return
1754 gridname, yac_read_icon_basic_grid_data(filename));
1755}
1756
1757/* ---------------------------------------------------------------- */
1758
1759static int * get_icon_cell_mask ( int ncid, size_t nbr_cells ) {
1760
1761 // get variable id
1762 int mask_id;
1763 yac_nc_inq_varid (ncid, "cell_sea_land_mask", &mask_id);
1764
1765 // check number of dimension (has to be 1)
1766 int ndims;
1767 YAC_HANDLE_ERROR(nc_inq_varndims(ncid, mask_id, &ndims));
1768 YAC_ASSERT(
1769 ndims == 1,
1770 "ERROR(get_icon_cell_mask): mask array has more than one dimension")
1771
1772 // get id of dimension
1773 int dimid;
1774 YAC_HANDLE_ERROR(nc_inq_vardimid(ncid, mask_id, &dimid));
1775
1776 // check size of mask (has to be equal to nbr_cells)
1777 size_t dimlen;
1778 YAC_HANDLE_ERROR(nc_inq_dimlen(ncid, dimid, &dimlen));
1779 YAC_ASSERT(
1780 dimlen == nbr_cells,
1781 "ERROR(get_icon_cell_mask): invalid size of mask array")
1782
1783 // check mask type (has to be NC_INT)
1784 nc_type mask_type;
1785 YAC_HANDLE_ERROR(nc_inq_vartype(ncid, mask_id, &mask_type));
1786 YAC_ASSERT(
1787 mask_type == NC_INT, "ERROR(get_icon_cell_mask): invalid mask type")
1788
1789 // get and return mask
1790 int * cell_mask = xmalloc(nbr_cells * sizeof(*cell_mask));
1791 YAC_HANDLE_ERROR(nc_get_var_int (ncid, mask_id, cell_mask));
1792 return cell_mask;
1793}
1794
1795/* ---------------------------------------------------------------- */
1796
1798 int ncid, size_t nbr_cells, int ** vertex_of_cell, size_t * nv) {
1799
1800 // get variable id
1801 int conn_id;
1802 yac_nc_inq_varid(ncid, "vertex_of_cell", &conn_id);
1803
1804 // check number of dimension (has to be 1)
1805 int ndims;
1806 YAC_HANDLE_ERROR(nc_inq_varndims(ncid, conn_id, &ndims));
1807 YAC_ASSERT(
1808 ndims == 2,
1809 "ERROR(get_icon_connect): "
1810 "connectivity array has invalid number of dimensions")
1811
1812 // get ids of dimensions
1813 int dimids[2];
1814 YAC_HANDLE_ERROR(nc_inq_vardimid(ncid, conn_id, dimids));
1815
1816 // check size of dimensions
1817 size_t dimlen;
1818 YAC_HANDLE_ERROR(nc_inq_dimlen(ncid, dimids[0], nv));
1819 YAC_HANDLE_ERROR(nc_inq_dimlen(ncid, dimids[1], &dimlen));
1820 YAC_ASSERT(
1821 dimlen == nbr_cells,
1822 "ERROR(get_icon_connect): invalid size of second dimension of "
1823 "connectivity array (has to be nbr_cells)")
1824
1825 // get and return connectivity array
1826 *vertex_of_cell = xmalloc(*nv * nbr_cells * sizeof(**vertex_of_cell));
1827 YAC_HANDLE_ERROR(nc_get_var_int (ncid, conn_id, *vertex_of_cell));
1828}
1829
1831 int ** global_cell_id,
1832 int ** global_cell_id_rank,
1833 int ** num_vertices_per_cell,
1834 int ** global_corner_id,
1835 int ** global_corner_id_rank,
1836 int ** cell_to_vertex,
1837 double ** x_cells,
1838 double ** y_cells,
1839 double ** x_vertices,
1840 double ** y_vertices) {
1841
1842 free (*cell_mask);
1843 free (*global_cell_id);
1844 free (*global_cell_id_rank);
1845 free (*num_vertices_per_cell);
1846 free (*global_corner_id);
1847 free (*global_corner_id_rank);
1848 free (*cell_to_vertex);
1849 free (*x_cells);
1850 free (*y_cells);
1851 free (*x_vertices);
1852 free (*y_vertices);
1853}
1854
1855#else
1856
1858 char const * filename) {
1859
1860 UNUSED(filename);
1861 die(
1862 "ERROR(yac_basic_grid_data): "
1863 "YAC is built without the NetCDF support");
1864
1865 return
1867 (size_t[]){0,0}, (int[]){0,0}, NULL, NULL);
1868}
1869
1871 char const * filename, char const * gridname) {
1872
1873 UNUSED(filename);
1874 UNUSED(gridname);
1875 die(
1876 "ERROR(yac_read_icon_basic_grid): "
1877 "YAC is built without the NetCDF support");
1878
1879 return NULL;
1880}
1881
1882void yac_read_icon_grid_information(const char * filename, int * nbr_vertices,
1883 int * nbr_cells, int ** num_vertices_per_cell,
1884 int ** cell_to_vertex, double ** x_vertices,
1885 double ** y_vertices, double ** x_cells,
1886 double ** y_cells, int ** cell_mask) {
1887
1888 UNUSED(filename);
1889 UNUSED(nbr_vertices);
1890 UNUSED(nbr_cells);
1891 UNUSED(num_vertices_per_cell);
1893 UNUSED(x_vertices);
1894 UNUSED(y_vertices);
1895 UNUSED(x_cells);
1896 UNUSED(y_cells);
1898 die(
1899 "ERROR(yac_read_icon_grid_information): "
1900 "YAC is built without the NetCDF support");
1901}
1902
1904 const char * filename, int * num_vertices, int * num_cells,
1905 int ** num_vertices_per_cell, int ** cell_to_vertex,
1906 double ** x_vertices, double ** y_vertices,
1907 double ** x_cells, double ** y_cells, int ** global_cell_id,
1908 int ** cell_mask, int ** cell_core_mask,
1909 int ** global_corner_id, int ** corner_core_mask, int rank, int size) {
1910
1911 UNUSED(filename);
1912 UNUSED(num_vertices);
1914 UNUSED(num_vertices_per_cell);
1916 UNUSED(x_vertices);
1917 UNUSED(y_vertices);
1918 UNUSED(x_cells);
1919 UNUSED(y_cells);
1920 UNUSED(global_cell_id);
1923 UNUSED(global_corner_id);
1924 UNUSED(corner_core_mask);
1925 UNUSED(rank);
1926 UNUSED(size);
1927 die(
1928 "ERROR(yac_read_part_icon_grid_information): "
1929 "YAC is built without the NetCDF support");
1930}
1931
1933 const char * filename, MPI_Comm comm, int * nbr_vertices, int * nbr_cells,
1934 int ** num_vertices_per_cell, int ** cell_to_vertex, int ** global_cell_ids,
1935 int ** cell_owner, int ** global_vertex_ids, int ** vertex_owner,
1936 double ** x_vertices, double ** y_vertices,
1937 double ** x_cells, double ** y_cells, int ** cell_msk) {
1938
1939 UNUSED(filename);
1940 UNUSED(comm);
1941 UNUSED(nbr_vertices);
1942 UNUSED(nbr_cells);
1943 UNUSED(num_vertices_per_cell);
1945 UNUSED(global_cell_ids);
1946 UNUSED(cell_owner);
1947 UNUSED(global_vertex_ids);
1948 UNUSED(vertex_owner);
1949 UNUSED(x_vertices);
1950 UNUSED(y_vertices);
1951 UNUSED(x_cells);
1952 UNUSED(y_cells);
1953 UNUSED(cell_msk);
1954 die(
1955 "ERROR(yac_read_icon_grid_information_parallel): "
1956 "YAC is built without the NetCDF support");
1957}
1958
1960 const char * filename, MPI_Comm comm) {
1961
1962 UNUSED(filename);
1963 UNUSED(comm);
1964 die(
1965 "ERROR(yac_read_icon_basic_grid_data_parallel): "
1966 "YAC is built without the NetCDF support");
1967
1968 return
1970 (size_t[]){0,0}, (int[]){0,0}, NULL, NULL);
1971}
1972
1974 char const * filename, char const * gridname, MPI_Comm comm) {
1975
1976 UNUSED(filename);
1977 UNUSED(gridname);
1978 UNUSED(comm);
1979 die(
1980 "ERROR(yac_read_icon_basic_grid_parallel): "
1981 "YAC is built without the NetCDF support");
1982
1983 return NULL;
1984}
1985
1987 int ** global_cell_id,
1988 int ** cell_core_mask,
1989 int ** num_vertices_per_cell,
1990 int ** global_corner_id,
1991 int ** corner_core_mask,
1992 int ** cell_to_vertex,
1993 double ** x_cells,
1994 double ** y_cells,
1995 double ** x_vertices,
1996 double ** y_vertices) {
1997
1999 UNUSED(global_cell_id);
2001 UNUSED(num_vertices_per_cell);
2002 UNUSED(global_corner_id);
2003 UNUSED(corner_core_mask);
2005 UNUSED(x_vertices);
2006 UNUSED(y_vertices);
2007 UNUSED(x_cells);
2008 UNUSED(y_cells);
2009 die(
2010 "ERROR(yac_delete_icon_grid_data): "
2011 "YAC is built without the NetCDF support");
2012}
2013
2014#endif // YAC_NETCDF_ENABLED
#define YAC_ASSERT(exp, msg)
struct yac_basic_grid * yac_basic_grid_new(char const *name, struct yac_basic_grid_data grid_data)
Definition basic_grid.c:53
size_t yac_basic_grid_add_coordinates_nocpy(struct yac_basic_grid *grid, enum yac_location location, yac_coordinate_pointer coordinates)
Definition basic_grid.c:211
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(size_t nbr_vertices, size_t nbr_cells, int *num_vertices_per_cell, double *x_vertices, double *y_vertices, int *cell_to_vertex)
#define UNUSED(x)
Definition core.h:72
#define YAC_RAD
yac_edge_type
Definition grid_cell.h:12
@ YAC_GREAT_CIRCLE_EDGE
great circle
Definition grid_cell.h:13
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
void yac_nc_inq_dimid(int ncid, char const *name, int *dimidp)
Definition io_utils.c:385
@ YAC_LOC_CELL
Definition location.h:14
#define xrealloc(ptr, size)
Definition ppm_xfuncs.h:67
#define xcalloc(nmemb, size)
Definition ppm_xfuncs.h:64
#define xmalloc(size)
Definition ppm_xfuncs.h:66
int yac_check_coord_units(int ncid, int varid)
Definition read_grid.c:45
void yac_read_coords(int ncid, char const *lon_name, char const *lat_name, double **lon, double **lat, size_t *len)
Definition read_grid.c:94
struct yac_basic_grid * yac_read_icon_basic_grid_parallel_f2c(char const *filename, char const *gridname, MPI_Fint comm)
static void convert_to_rad(int ncid, int varid, double *array, size_t count)
void yac_delete_icon_grid_data(int **cell_mask, int **global_cell_id, int **global_cell_id_rank, int **num_vertices_per_cell, int **global_corner_id, int **global_corner_id_rank, int **cell_to_vertex, double **x_cells, double **y_cells, double **x_vertices, double **y_vertices)
struct yac_basic_grid * yac_read_icon_basic_grid(char const *filename, char const *gridname)
void yac_read_part_icon_grid_information(const char *filename, int *nbr_vertices, int *nbr_cells, int **num_vertices_per_cell, int **cell_to_vertex, double **x_vertices, double **y_vertices, double **x_cells, double **y_cells, int **global_cell_id, int **cell_mask, int **cell_core_mask, int **global_corner_id, int **corner_core_mask, int rank, int size)
struct yac_basic_grid_data yac_read_icon_basic_grid_data_parallel(const char *filename, MPI_Comm comm)
static int * get_icon_cell_mask(int ncid, size_t nbr_cells)
void yac_read_icon_grid_information_parallel(const char *filename, MPI_Comm comm, int *nbr_vertices, int *nbr_cells, int **num_vertices_per_cell, int **cell_to_vertex, int **global_cell_ids, int **cell_owner, int **global_vertex_ids, int **vertex_owner, double **x_vertices, double **y_vertices, double **x_cells, double **y_cells, int **cell_msk)
void yac_read_icon_basic_grid_parallel_2(char const *filename, char const *gridname, MPI_Comm comm, struct yac_basic_grid **basic_grid, size_t *cell_coordinate_idx, int **cell_mask)
void yac_read_icon_grid_information_parallel_f2c(const char *filename, MPI_Fint comm, int *nbr_vertices, int *nbr_cells, int **num_vertices_per_cell, int **cell_to_vertex, int **global_cell_ids, int **cell_owner, int **global_vertex_ids, int **vertex_owner, double **x_vertices, double **y_vertices, double **x_cells, double **y_cells, int **cell_msk)
static int partition_idx_from_element_idx(unsigned element_idx, unsigned num_elements, int num_partitions)
void yac_read_icon_grid_information_parallel_2(const char *filename, MPI_Comm comm, double **x_vertices, double **y_vertices, yac_int **cell_ids, yac_int **vertex_ids, yac_int **edge_ids, size_t *num_cells, size_t *num_vertices, size_t *num_edges, int **num_vertices_per_cell, int **num_cells_per_vertex, size_t **cell_to_vertex, size_t **cell_to_edge, size_t **vertex_to_cell, size_t **edge_to_vertex, enum yac_edge_type **edge_type, double **x_cells, double **y_cells, int **cell_msk)
static size_t * generate_offsets(size_t N, int *counts)
static int * generate_simple_core_mask(size_t N)
void yac_read_icon_grid_information(const char *filename, int *nbr_vertices, int *nbr_cells, int **num_vertices_per_cell, int **cell_to_vertex, double **x_vertices, double **y_vertices, double **x_cells, double **y_cells, int **cell_mask)
static void get_icon_connect(int ncid, size_t nbr_cells, int **vertex_of_cell, size_t *nv)
struct yac_basic_grid * yac_read_icon_basic_grid_parallel(char const *filename, char const *gridname, MPI_Comm comm)
struct yac_basic_grid_data yac_read_icon_basic_grid_data(char const *filename)
yac_coordinate_pointer vertex_coordinates
size_t * vertex_to_cell_offsets
yac_size_t_2_pointer edge_to_vertex
enum yac_edge_type * edge_type
size_t * cell_to_vertex_offsets
int * cell_to_vertex
int * cell_mask
double * cell_lat
double * cell_lon
#define N
size_t num_cells[2]
int vertex_of_cell[3][16]
static int mask[16]
int * cell_core_mask
double * buffer
#define YAC_HANDLE_ERROR(exp)
Definition toy_output.c:13
static void LLtoXYZ(double lon, double lat, double p_out[])
Definition toy_scrip.c:587
static void yac_remove_duplicates_int(int *array, size_t *n)
#define MAX(a, b)
void yac_quicksort_index_yac_int_size_t(yac_int *a, size_t n, size_t *idx)
void yac_quicksort_index_int_size_t(int *a, size_t n, size_t *idx)
void yac_quicksort_index(int *a, size_t n, int *idx)
#define die(msg)
Definition yac_assert.h:12
YAC_INT yac_int
Definition yac_types.h:15
size_t(* yac_size_t_2_pointer)[2]
Definition yac_types.h:25
#define yac_int_dt
Definition yac_types.h:18
double(* yac_coordinate_pointer)[3]
Definition yac_types.h:21