44 int * nbr_cells,
int ** num_vertices_per_cell,
45 int ** cell_to_vertex,
double ** x_vertices,
46 double ** y_vertices,
double ** x_cells,
47 double ** y_cells,
int ** cell_mask) {
59 double * temp_vertex_lon;
60 double * temp_vertex_lat;
61 size_t temp_nbr_vertices;
64 ncid, &temp_vertex_lon, &temp_vertex_lat, &temp_nbr_vertices);
66 *nbr_cells = (int)(temp_nbr_vertices / 4);
68 int * old_to_new_id =
xmalloc(temp_nbr_vertices *
sizeof(*old_to_new_id));
71 temp_vertex_lon, temp_vertex_lat, NULL,
72 &temp_nbr_vertices, old_to_new_id);
74 *nbr_vertices = temp_nbr_vertices;
77 xrealloc(temp_vertex_lon, temp_nbr_vertices *
sizeof(*temp_vertex_lon));
79 xrealloc(temp_vertex_lat, temp_nbr_vertices *
sizeof(*temp_vertex_lat));
81 for (
size_t i = 0; i < temp_nbr_vertices; ++i) {
86 *num_vertices_per_cell =
xmalloc(*nbr_cells *
sizeof(**num_vertices_per_cell));
87 for (
int i = 0; i < *nbr_cells; ++i) (*num_vertices_per_cell)[i] = 4;
89 *cell_to_vertex =
xmalloc(*nbr_cells * 4 *
sizeof(**cell_to_vertex));
92 for (
int i = 0; i < *nbr_cells; ++i)
93 for (
int j = 0; j < 4; ++j)
94 (*cell_to_vertex)[4*i+j] = old_to_new_id[4*i+j];
101 double * temp_cell_lon;
102 double * temp_cell_lat;
103 int * temp_cell_mask;
105 size_t temp_nbr_cells;
106 size_t temp_nbr_masks;
109 ncid, &temp_cell_lon, &temp_cell_lat, &temp_nbr_cells);
114 temp_nbr_cells == temp_nbr_masks,
115 "ERROR(yac_read_mpiom_grid_information): "
116 "missmatch between cell count and mask value count")
118 int * old_to_new_id =
xmalloc(temp_nbr_cells *
sizeof(*old_to_new_id));
121 temp_cell_lon, temp_cell_lat, temp_cell_mask,
122 &temp_nbr_cells, old_to_new_id);
124 *nbr_cells = temp_nbr_cells;
126 *x_cells =
xrealloc(temp_cell_lon, temp_nbr_cells *
sizeof(**x_cells));
127 *y_cells =
xrealloc(temp_cell_lat, temp_nbr_cells *
sizeof(**y_cells));
129 for (
size_t i = 0; i < temp_nbr_cells; ++i) {
134 *cell_mask =
xrealloc(temp_cell_mask,temp_nbr_cells *
sizeof(
int));
135 for (
size_t i = 0; i < temp_nbr_cells; ++i)
136 (*cell_mask)[i] = !(*cell_mask)[i];
144 const char * filename,
int * num_vertices,
int * num_cells,
145 int ** num_vertices_per_cell,
int ** cell_to_vertex,
double ** x_vertices,
146 double ** y_vertices,
double ** x_cells,
147 double ** y_cells,
int ** global_cell_id,
148 int ** cell_mask,
int ** cell_core_mask,
149 int ** global_corner_id,
int ** corner_core_mask,
int rank,
int size) {
153 num_vertices_per_cell, cell_to_vertex,
154 x_vertices, y_vertices, x_cells, y_cells, cell_mask);
157 int num_cells_per_process = (*num_cells + size - 1)/size;
158 int local_start = rank * num_cells_per_process;
159 int num_local_cells =
MIN(num_cells_per_process,
160 *num_cells - local_start);
163 int * required_vertices =
xcalloc(*num_vertices,
sizeof(*required_vertices));
164 int * required_cells =
xcalloc(*num_cells,
sizeof(*required_cells));
167 for (
int i = 0; i < local_start; ++i)
168 offset += (*num_vertices_per_cell)[i];
171 for (
int i = local_start; i < local_start + num_local_cells; ++i) {
173 required_cells[i] = 2;
174 for (
int j = 0; j < (*num_vertices_per_cell)[i]; ++j)
175 required_vertices[(*cell_to_vertex)[offset+j]] = 2;
177 offset += (*num_vertices_per_cell)[i];
182 for (
int i = 0; i < *num_cells; ++i) {
184 if (!required_cells[i]) {
186 for (
int j = 0; j < (*num_vertices_per_cell)[i]; ++j) {
188 if (required_vertices[(*cell_to_vertex)[offset+j]]) {
190 required_cells[i] = 1;
195 offset += (*num_vertices_per_cell)[i];
200 for (
int i = 0; i < *num_cells; ++i) {
202 if (required_cells[i] == 1) {
204 for (
int j = 0; j < (*num_vertices_per_cell)[i]; ++j) {
206 if (!required_vertices[(*cell_to_vertex)[offset+j]]) {
208 required_vertices[(*cell_to_vertex)[offset+j]] = 1;
213 offset += (*num_vertices_per_cell)[i];
217 int part_num_vertices = 0;
218 int part_num_cells = 0;
219 for (
int i = 0; i < *num_vertices; ++i)
220 if (required_vertices[i])
222 for (
int i = 0; i < *num_cells; ++i)
223 if(required_cells[i])
226 *global_cell_id =
xmalloc(part_num_cells *
sizeof(**global_cell_id));
227 *cell_core_mask =
xmalloc(part_num_cells *
sizeof(**cell_core_mask));
228 *global_corner_id =
xmalloc(part_num_vertices *
sizeof(**global_corner_id));
229 *corner_core_mask =
xmalloc(part_num_vertices *
sizeof(**corner_core_mask));
232 part_num_vertices = 0;
233 int * global_to_local_vertex =
xmalloc(*num_vertices *
sizeof(*global_to_local_vertex));
234 for (
int i = 0; i < *num_vertices; ++i) {
236 if (required_vertices[i]) {
238 (*global_corner_id)[part_num_vertices] = i;
239 (*corner_core_mask)[part_num_vertices] = required_vertices[i] == 2;
240 (*x_vertices)[part_num_vertices] = (*x_vertices)[i];
241 (*y_vertices)[part_num_vertices] = (*y_vertices)[i];
242 global_to_local_vertex[i] = part_num_vertices;
247 *x_vertices =
xrealloc(*x_vertices, part_num_vertices *
sizeof(**x_vertices));
248 *y_vertices =
xrealloc(*y_vertices, part_num_vertices *
sizeof(**y_vertices));
249 *num_vertices = part_num_vertices;
250 free(required_vertices);
253 int num_cell_vertex_dependencies = 0;
256 for (
int i = 0; i < *num_cells; ++i) {
258 if (required_cells[i]) {
260 (*global_cell_id)[part_num_cells] = i;
261 (*cell_core_mask)[part_num_cells] = required_cells[i] == 2;
262 (*x_cells)[part_num_cells] = (*x_cells)[i];
263 (*y_cells)[part_num_cells] = (*y_cells)[i];
264 (*cell_mask)[part_num_cells] = (*cell_mask)[i];
266 for (
int j = 0; j < (*num_vertices_per_cell)[i]; ++j)
267 (*cell_to_vertex)[num_cell_vertex_dependencies++] =
268 global_to_local_vertex[(*cell_to_vertex)[offset+j]];
270 (*num_vertices_per_cell)[part_num_cells] = (*num_vertices_per_cell)[i];
275 offset += (*num_vertices_per_cell)[i];
278 *x_cells =
xrealloc(*x_cells, part_num_cells *
sizeof(**x_cells));
279 *y_cells =
xrealloc(*y_cells, part_num_cells *
sizeof(**y_cells));
280 *cell_mask =
xrealloc(*cell_mask, part_num_cells *
sizeof(**cell_mask));
282 *num_vertices_per_cell =
xrealloc(*num_vertices_per_cell, part_num_cells *
283 sizeof(**num_vertices_per_cell));
284 *cell_to_vertex =
xrealloc(*cell_to_vertex, num_cell_vertex_dependencies *
285 sizeof(**cell_to_vertex));
286 *num_cells = part_num_cells;
287 free(required_cells);
288 free(global_to_local_vertex);
338 int ncid,
double ** vertex_lon,
double ** vertex_lat,
339 size_t * nbr_vertices) {
345 if (nc_inq_varid(ncid,
"lon_bounds", &glon_id) != NC_NOERR)
348 if (nc_inq_varid(ncid,
"lat_bounds", &glat_id) != NC_NOERR)
358 (glon_ndims == 3) && (glat_ndims == 3),
359 "ERROR(get_mpiom_vertices): "
360 "number of dimensions for lon and lat does not match")
369 size_t dimlen_lat[3];
370 size_t dimlen_lon[3];
373 for (
int i = 0; i < 3; ++i) {
377 dimlen_lon[i] == dimlen_lat[i],
378 "ERROR(get_mpiom_vertices): dimlen_lon[i] != dimlen_lat[i]")
381 dimlen_lon[2] == 4,
"ERROR(get_mpiom_vertices): dimlen_lon[2] has to be 4")
383 *nbr_vertices = dimlen_lon[0] * dimlen_lon[1] * dimlen_lon[2];
385 *vertex_lon =
xmalloc(*nbr_vertices *
sizeof(**vertex_lon));
386 *vertex_lat =
xmalloc(*nbr_vertices *
sizeof(**vertex_lat));
397 size_t *nbr_cells ) {
413 (glon_ndims == 2) && (glat_ndims == 2),
414 "ERROR(get_mpiom_cell_center): (glon_ndims != 2) || (glat_ndims != 2)")
423 size_t dimlen_lon[2];
424 size_t dimlen_lat[2];
427 for (
int i = 0; i < 2; ++i) {
431 dimlen_lon[i] == dimlen_lat[i],
432 "ERROR(get_mpiom_cell_center): dimlen_lon[i] != dimlen_lat[i]")
435 *nbr_cells = dimlen_lon[0] * dimlen_lon[1];
437 *cell_lon =
xmalloc(*nbr_cells *
sizeof(**cell_lon));
438 *cell_lat =
xmalloc(*nbr_cells *
sizeof(**cell_lat));
447 int ncid,
int **cell_mask,
size_t *nbr_cells ) {
453 status = nc_inq_varid (ncid,
"var1", &mask_id);
460 (mask_ndims == 3) || (mask_ndims == 4),
461 "ERROR(get_mpiom_cell_mask): invalid number of dimensions")
468 size_t dimlen_mask[4] = {0, 0, 0, 0};
469 for (
int i = 0; i < mask_ndims; ++i)
472 *nbr_cells = dimlen_mask[mask_ndims - 2] * dimlen_mask[mask_ndims - 1];
474 int * temp_cell_mask =
xmalloc(*nbr_cells *
sizeof(*temp_cell_mask));
479 for (
size_t i = 0; i < *nbr_cells; ++i)
480 temp_cell_mask[i] = !(temp_cell_mask[i] == 1);
482 *cell_mask = temp_cell_mask;
499 double * temp_coords_lon,
double * temp_coords_lat,
int * temp_coords_mask,
500 size_t * temp_nbr_coords,
int * old_to_new_id) {
503 xmalloc(*temp_nbr_coords *
sizeof(*sort_array));
505 double const scale = (double)(2 << 20);
506 int32_t
const periode = (int32_t)(360.0 * scale);
508 for (
size_t i = 0;
i < *temp_nbr_coords; ++
i) {
510 int32_t curr_lon, curr_lat;
512 curr_lon = (int32_t)(temp_coords_lon[
i] * scale);
513 curr_lat = (int32_t)(temp_coords_lat[
i] * scale);
516 curr_lon = curr_lon - (((curr_lon - periode) / periode) * periode);
517 else if (curr_lon > periode)
518 curr_lon = curr_lon - ((curr_lon / periode) * periode);
520 sort_array[
i].
lon = curr_lon;
521 sort_array[
i].
lat = curr_lat;
526 qsort(sort_array, *temp_nbr_coords,
sizeof(*sort_array),
529 old_to_new_id[sort_array[0].
i] = 1;
531 int last_unique_idx = (int)(sort_array[0].
i);
533 for (
size_t i = 1;
i < *temp_nbr_coords; ++
i) {
537 old_to_new_id[sort_array[
i].
i] = 1;
538 last_unique_idx = (int)(sort_array[
i].
i);
542 old_to_new_id[sort_array[
i].
i] = -last_unique_idx;
548 size_t new_nbr_coords = 0;
550 for (
size_t i = 0;
i < *temp_nbr_coords; ++
i) {
552 if (old_to_new_id[
i] == 1) {
554 temp_coords_lon[new_nbr_coords] = temp_coords_lon[
i];
555 temp_coords_lat[new_nbr_coords] = temp_coords_lat[
i];
557 if (temp_coords_mask != NULL)
558 temp_coords_mask[new_nbr_coords] = temp_coords_mask[
i];
560 old_to_new_id[
i] = new_nbr_coords;
566 for (
size_t i = 0;
i < *temp_nbr_coords; ++
i)
567 if (old_to_new_id[
i] <= 0)
568 old_to_new_id[
i] = old_to_new_id[-old_to_new_id[
i]];
570 *temp_nbr_coords = new_nbr_coords;
void yac_read_part_mpiom_grid_information(const char *filename, int *num_vertices, int *num_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)