YetAnotherCoupler 3.5.2
Loading...
Searching...
No Matches
interp_method_spmap.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// Get the definition of the 'restrict' keyword.
7#include "config.h"
8#endif
9
10#include <float.h>
11#include <string.h>
12
14#include "interp_method_spmap.h"
15#include "ensure_array_size.h"
16#include "area.h"
17#include "io_utils.h"
18#include "yac_mpi_internal.h"
19
20static size_t do_search_spmap(struct interp_method * method,
21 struct yac_interp_grid * interp_grid,
22 size_t * tgt_points, size_t count,
23 struct yac_interp_weights * weights);
24static void delete_spmap(struct interp_method * method);
25
26static struct interp_method_vtable
31
40
41static inline int compare_size_t(const void * a, const void * b) {
42
43 size_t const * a_ = a, * b_ = b;
44
45 return (*a_ > *b_) - (*b_ > *a_);
46}
47
48static int lists_overlap(
49 size_t * list_a, size_t count_a, size_t * list_b, size_t count_b) {
50
51 if ((count_a == 0) || (count_b == 0)) return 0;
52
53 size_t i = 0, j = 0;
54 size_t curr_a = SIZE_MAX, curr_b = list_b[0];
55
56 do {
57 while ((i < count_a) && (((curr_a = list_a[i++])) < curr_b));
58 if (curr_a == curr_b) return 1;
59 while ((j < count_b) && (((curr_b = list_b[j++])) < curr_a));
60 if (curr_a == curr_b) return 1;
61 } while ((i < count_a) || (j < count_b));
62
63 return 0;
64}
65
66static void merge_lists(
67 size_t * list, size_t * list_size, size_t * insert, size_t insert_size) {
68
69 if (insert_size == 0) return;
70
71 size_t new_list_size = *list_size;
72 size_t old_list_size = *list_size;
73
74 for (size_t i = 0, j = 0; i < insert_size; ++i) {
75 size_t curr_insert = insert[i];
76 while ((j < old_list_size) && (list[j] < curr_insert)) ++j;
77 if ((j >= old_list_size) || (list[j] != curr_insert))
78 list[new_list_size++] = curr_insert;
79 }
80
81 if (new_list_size != old_list_size) {
82 qsort(list, new_list_size, sizeof(*list), compare_size_t);
83 *list_size = new_list_size;
84 }
85}
86
88 yac_const_coordinate_pointer tgt_field_coords, double max_distance,
89 size_t tgt_start_point, size_t * tgt_points, size_t * count) {
90
91 double const * start_coord = tgt_field_coords[tgt_start_point];
92 size_t new_count = 0;
93 for (size_t i = 0, old_count = *count; i < old_count; ++i) {
95 start_coord, tgt_field_coords[tgt_points[i]]) <= max_distance) {
96 if (new_count != i) tgt_points[new_count] = tgt_points[i];
97 ++new_count;
98 }
99 }
100
101 *count = new_count;
102}
103
105 struct yac_interp_grid * interp_grid, size_t tgt_start_point,
106 size_t * from_tgt_points, size_t * to_tgt_points,
107 size_t * count, int * flag, size_t * temp_cell_edges,
108 size_t ** edges_buffer, size_t * edges_buffer_array_size) {
109
110 struct yac_const_basic_grid_data * tgt_grid_data =
114 tgt_grid_data->cell_to_edge_offsets;
115 const int * num_vertices_per_cell =
116 tgt_grid_data->num_vertices_per_cell;
117
118 size_t old_count = *count;
119 memset(flag, 0, old_count * sizeof(*flag));
120
121 for (size_t i = 0; i < old_count; ++i) {
122 if (from_tgt_points[i] == tgt_start_point) {
123 flag[i] = 1;
124 break;
125 }
126 }
127
128 size_t * edges = *edges_buffer;
129 size_t num_edges = 0;
130 size_t edges_array_size = *edges_buffer_array_size;
131
132 const_size_t_pointer curr_edges =
133 cell_to_edge + cell_to_edge_offsets[tgt_start_point];
134 size_t curr_num_edges = num_vertices_per_cell[tgt_start_point];
135
136 ENSURE_ARRAY_SIZE(edges, edges_array_size, curr_num_edges);
137 num_edges = curr_num_edges;
138 memcpy(edges, curr_edges, num_edges * sizeof(*edges));
139 qsort(edges, num_edges, sizeof(*edges), compare_size_t);
140
141 int change_flag = 0;
142 do {
143
144 change_flag = 0;
145
146 for (size_t i = 0; i < old_count; ++i) {
147
148 if (flag[i]) continue;
149
150 const_size_t_pointer curr_edges =
151 tgt_grid_data->cell_to_edge + cell_to_edge_offsets[from_tgt_points[i]];
152 curr_num_edges = num_vertices_per_cell[from_tgt_points[i]];
153 memcpy(
154 temp_cell_edges, curr_edges, curr_num_edges * sizeof(*temp_cell_edges));
155 qsort(temp_cell_edges, curr_num_edges, sizeof(*edges), compare_size_t);
156 ENSURE_ARRAY_SIZE(edges, edges_array_size, num_edges + curr_num_edges);
157
158 if (lists_overlap(edges, num_edges, temp_cell_edges, curr_num_edges)) {
159 merge_lists(edges, &num_edges, temp_cell_edges, curr_num_edges);
160 flag[i] = 1;
161 change_flag = 1;
162 }
163 }
164 } while (change_flag);
165
166 *edges_buffer = edges;
167 *edges_buffer_array_size = edges_array_size;
168
169 size_t new_count = 0;
170 for (size_t i = 0; i < old_count; ++i)
171 if (flag[i]) to_tgt_points[new_count++] = from_tgt_points[i];
172
173 *count = new_count;
174}
175
177 struct yac_interp_grid * interp_grid, double spread_distance,
178 size_t num_src_points, size_t const * const tgt_result_points,
179 size_t * num_tgt_per_src, size_t * spread_tgt_result_points) {
180
181 size_t max_num_tgt_per_src = 0;
182 for (size_t i = 0; i < num_src_points; ++i)
183 if (num_tgt_per_src[i] > max_num_tgt_per_src)
184 max_num_tgt_per_src = num_tgt_per_src[i];
185
186 int * flag = xmalloc(max_num_tgt_per_src * sizeof(*flag));
187
188 size_t new_offset = 0;
189 size_t * cell_edge_buffer = NULL;
190 size_t cell_edge_buffer_array_size = 0;
191 size_t * edge_buffer = NULL;
192 size_t edge_buffer_array_size = 0;
193 int max_num_vertice_per_tgt = 0;
194 const int * num_vertices_per_tgt =
197
198 yac_const_coordinate_pointer tgt_field_coords =
200
201 // for all source points
202 for (size_t i = 0, old_offset = 0; i < num_src_points; ++i) {
203
204 size_t * old_results = spread_tgt_result_points + old_offset;
205 size_t * new_results = spread_tgt_result_points + new_offset;
206 old_offset += num_tgt_per_src[i];
207
208 // remove all target points, which exceed the spread distance from
209 // the original tgt
211 tgt_field_coords, spread_distance, tgt_result_points[i],
212 old_results, num_tgt_per_src + i);
213
214 // check buffer sizes required by routine remove_disconnected_points
215 for (size_t j = 0, curr_num_tgt_per_src = num_tgt_per_src[i];
216 j < curr_num_tgt_per_src; ++j)
217 if (num_vertices_per_tgt[old_results[j]] > max_num_vertice_per_tgt)
218 max_num_vertice_per_tgt = num_vertices_per_tgt[old_results[j]];
219
221 cell_edge_buffer, cell_edge_buffer_array_size,
222 max_num_vertice_per_tgt);
223
224 // remove all tgts that are not directly connected to the original tgt
226 interp_grid, tgt_result_points[i],
227 old_results, new_results, num_tgt_per_src + i, flag,
228 cell_edge_buffer, &edge_buffer, &edge_buffer_array_size);
229
230 new_offset += num_tgt_per_src[i];
231 }
232
233 free(edge_buffer);
234 free(cell_edge_buffer);
235 free(flag);
236
237 return new_offset;
238}
239
240static double * compute_weights(
241 struct yac_interp_grid * interp_grid,
242 enum yac_interp_spmap_weight_type weight_type, size_t num_src_points,
243 size_t const * const src_points, size_t const * const num_tgt_per_src,
244 size_t total_num_tgt, size_t const * const tgt_result_points) {
245
246 double * weights = xmalloc(total_num_tgt * sizeof(*weights));
250 "ERROR(do_search_spmap): invalid weight_type")
251 switch (weight_type) {
252 case (YAC_INTERP_SPMAP_AVG): {
253 for (size_t i = 0, offset = 0; i < num_src_points; ++i) {
254 size_t curr_num_tgt = num_tgt_per_src[i];
255 if (curr_num_tgt == 0) continue;
256 if (curr_num_tgt > 1) {
257 double curr_weight_data = 1.0 / (double)(curr_num_tgt);
258 for (size_t j = 0; j < curr_num_tgt; ++j, ++offset) {
259 weights[offset] = curr_weight_data;
260 }
261 } else {
262 weights[offset] = 1.0;
263 ++offset;
264 }
265 }
266 break;
267 }
268 default:
269 case (YAC_INTERP_SPMAP_DIST): {
270
271 yac_const_coordinate_pointer src_field_coords =
273 yac_const_coordinate_pointer tgt_field_coords =
275
276 for (size_t i = 0, offset = 0; i < num_src_points; ++i) {
277
278 size_t curr_num_tgt = num_tgt_per_src[i];
279
280 if (curr_num_tgt == 0) continue;
281
282 double * curr_weights = weights + offset;
283
284 if (curr_num_tgt > 1) {
285
286 size_t const * const curr_result_points =
287 tgt_result_points + offset;
288 double const * curr_src_coord = src_field_coords[src_points[offset]];
289 offset += curr_num_tgt;
290
291 int match_flag = 0;
292
293 for (size_t j = 0; j < curr_num_tgt; ++j) {
294
295 double distance =
297 (double*)curr_src_coord,
298 (double*)tgt_field_coords[curr_result_points[j]]);
299
300 if (distance < yac_angle_tol) {
301 for (size_t k = 0; k < curr_num_tgt; ++k) curr_weights[k] = 0.0;
302 curr_weights[j] = 1.0;
303 match_flag = 1;
304 break;
305 }
306 curr_weights[j] = 1.0 / distance;
307 }
308
309 if (!match_flag) {
310
311 // compute scaling factor for the weights
312 double inv_distance_sum = 0.0;
313 for (size_t j = 0; j < curr_num_tgt; ++j)
314 inv_distance_sum += curr_weights[j];
315 double scale = 1.0 / inv_distance_sum;
316
317 for (size_t j = 0; j < curr_num_tgt; ++j) curr_weights[j] *= scale;
318 }
319 } else {
320 *curr_weights = 1.0;
321 ++offset;
322 }
323 }
324 break;
325 }
326 };
327
328 return weights;
329}
330
332 struct yac_const_basic_grid_data * basic_grid_data,
333 int * required_cell_areas, double area_scale, char const * type,
334 double * cell_areas) {
335
336 struct yac_grid_cell grid_cell;
337 yac_init_grid_cell(&grid_cell);
338
339 size_t num_cells = basic_grid_data->count[YAC_LOC_CELL];
340
341 for (size_t i = 0; i < num_cells; ++i) {
342 if (required_cell_areas[i]) {
343 yac_const_basic_grid_data_get_grid_cell(basic_grid_data, i, &grid_cell);
344 double cell_area = yac_huiliers_area(grid_cell) * area_scale;
346 cell_area > YAC_AREA_TOL,
347 "ERROR(get_cell_areas): "
348 "area of %s cell (global id %"XT_INT_FMT") is close to zero (%e)",
349 type, basic_grid_data->ids[YAC_LOC_CELL][i], cell_area);
350 cell_areas[i] = cell_area;
351 } else {
352 cell_areas[i] = 0.0;
353 }
354 }
355
356 yac_free_grid_cell(&grid_cell);
357}
358
360 char const * filename, char const * varname, MPI_Comm comm,
361 int * io_ranks, int io_rank_idx, int num_io_ranks,
362 double ** dist_cell_areas, size_t * dist_cell_areas_global_size) {
363
364#ifndef YAC_NETCDF_ENABLED
365
366 UNUSED(filename);
367 UNUSED(varname);
368 UNUSED(comm);
369 UNUSED(io_ranks);
370 UNUSED(io_rank_idx);
371 UNUSED(num_io_ranks);
372
373 *dist_cell_areas = NULL;
374 *dist_cell_areas_global_size = 0;
375
376 die(
377 "ERROR(interp_method_spmap::dist_read_cell_areas): "
378 "YAC is built without the NetCDF support");
379#else
380
381 if ((io_rank_idx >= 0) && (io_rank_idx < num_io_ranks)) {
382
383 // open file
384 int ncid;
385 yac_nc_open(filename, NC_NOWRITE, &ncid);
386
387 // get variable id
388 int varid;
389 yac_nc_inq_varid(ncid, varname, &varid);
390
391 // get dimension ids
392 int ndims;
393 int dimids[NC_MAX_VAR_DIMS];
395 nc_inq_var(ncid, varid, NULL, NULL, &ndims, dimids, NULL));
396
398 (ndims == 1) || (ndims == 2),
399 "ERROR(dist_read_cell_areas): "
400 "invalid number of dimensions for cell area variable \"%s\" from "
401 "file \"%s\" (is %d, but should be either 1 or 2)",
402 varname, filename, ndims);
403
404 // get size of dimensions
405 size_t dimlen[2];
406 *dist_cell_areas_global_size = 1;
407 for (int i = 0; i < ndims; ++i) {
408 YAC_HANDLE_ERROR(nc_inq_dimlen(ncid, dimids[i], &dimlen[i]));
410 dimlen[i] > 0,
411 "ERROR(dist_read_cell_areas): "
412 "invalid dimension size for cell area variable \"%s\" from "
413 "file \"%s\" (is %zu, should by > 0)",
414 varname, filename, dimlen[i]);
415 *dist_cell_areas_global_size *= dimlen[i];
416 }
417
418 // compute start/count
419 // (in 2D case we have to round a little bit and then adjust the
420 // data afterwards)
421 size_t start[2], count[2], offset, read_size;
422 size_t local_start =
423 (size_t)(
424 ((long long)*dist_cell_areas_global_size * (long long)io_rank_idx) /
425 (long long)num_io_ranks);
426 size_t local_count =
427 (size_t)(
428 ((long long)*dist_cell_areas_global_size *
429 (long long)(io_rank_idx+1)) / (long long)num_io_ranks) - local_start;
430 if (ndims == 1) {
431 start[0] = local_start;
432 count[0] = local_count;
433 offset = 0;
434 read_size = local_count;
435 } else {
436 start[0] = local_start / dimlen[1];
437 count[0] =
438 (local_start + local_count + dimlen[1] - 1) / dimlen[1] - start[0];
439 start[1] = 0;
440 count[1] = dimlen[1];
441 offset = local_start - start[0] * dimlen[1];
442 read_size = count[0] * count[1];
443 }
444
445 // read data
446 *dist_cell_areas = xmalloc(read_size * sizeof(**dist_cell_areas));
448 nc_get_vara_double(ncid, varid, start, count, *dist_cell_areas));
449
450 // adjust data if necessary
451 if (ndims == 2)
452 memmove(
453 *dist_cell_areas, *dist_cell_areas + offset,
454 local_count * sizeof(**dist_cell_areas));
455
456 // close file
457 YAC_HANDLE_ERROR(nc_close(ncid));
458
459 } else {
460 *dist_cell_areas = xmalloc(1*sizeof(**dist_cell_areas));
461 *dist_cell_areas_global_size = 0;
462 }
463
465 MPI_Bcast(
466 dist_cell_areas_global_size, 1, YAC_MPI_SIZE_T, io_ranks[0], comm),
467 comm);
468
469#endif // YAC_NETCDF_ENABLED
470}
471
472static void read_cell_areas(
473 struct yac_const_basic_grid_data * basic_grid_data,
474 struct yac_spmap_cell_area_file_config file_config, MPI_Comm comm,
475 int * required_cell_areas, char const * type,
476 double * cell_areas) {
477
478 size_t num_cells = basic_grid_data->count[YAC_LOC_CELL];
479
480 // get io configuration
481 int local_is_io, * io_ranks, num_io_ranks;
482 yac_get_io_ranks(comm, &local_is_io, &io_ranks, &num_io_ranks);
483
484 int comm_rank, comm_size;
485 yac_mpi_call(MPI_Comm_rank(comm, &comm_rank), comm);
486 yac_mpi_call(MPI_Comm_size(comm, &comm_size), comm);
487
488 int io_rank_idx = 0;
489 while ((io_rank_idx < num_io_ranks) &&
490 (comm_rank != io_ranks[io_rank_idx]))
491 ++io_rank_idx;
493 !local_is_io || (io_rank_idx < num_io_ranks),
494 "ERROR(read_cell_areas): unable to determine io_rank_idx");
495
496 double * dist_cell_areas;
497 size_t dist_cell_areas_global_size;
498
499 // read the data on the io processes
501 file_config.filename, file_config.varname, comm,
502 io_ranks, io_rank_idx, num_io_ranks,
503 &dist_cell_areas, &dist_cell_areas_global_size);
504
505 // count the number of locally required cell areas
506 size_t num_required_cell_areas = 0;
507 for (size_t i = 0; i < num_cells; ++i)
508 if (required_cell_areas[i]) ++num_required_cell_areas;
509
510 size_t * global_idx = xmalloc(num_required_cell_areas * sizeof(*global_idx));
511 size_t * reorder_idx =
512 xmalloc(num_required_cell_areas * sizeof(*reorder_idx));
513
514 size_t * sendcounts, * recvcounts, * sdispls, * rdispls;
516 1, &sendcounts, &recvcounts, &sdispls, &rdispls, comm);
517
518 // determine global indices of locally required cell areas
519 yac_int const * global_cell_ids = basic_grid_data->ids[YAC_LOC_CELL];
520 yac_int min_global_id = file_config.min_global_id;
521 for (size_t i = 0, j = 0; i < num_cells; ++i) {
522 if (required_cell_areas[i]) {
524 global_cell_ids[i] >= min_global_id,
525 "ERROR(read_cell_areas): %s global id %" XT_INT_FMT " is smaller than "
526 "the minimum global id provided by the user (%" XT_INT_FMT ")",
527 type, global_cell_ids[i], min_global_id);
528 size_t idx = (size_t)(global_cell_ids[i] - min_global_id);
530 idx < dist_cell_areas_global_size,
531 "ERROR(read_cell_areas): %s global id %" XT_INT_FMT " exceeds "
532 "available size of array \"%s\" in file \"%s\" "
533 "(min_global_id %" XT_INT_FMT ")", type, global_cell_ids[i],
534 file_config.varname, file_config.filename, file_config.min_global_id);
535 global_idx[j] = idx;
536 reorder_idx[j] = i;
537 int dist_rank_idx =
538 ((long long)idx * (long long)num_io_ranks +
539 (long long)num_io_ranks - 1) / (long long)dist_cell_areas_global_size;
540 sendcounts[io_ranks[dist_rank_idx]]++;
541 ++j;
542 }
543 }
544 free(io_ranks);
545
547 1, sendcounts, recvcounts, sdispls, rdispls, comm);
548
549 // sort required data by their global index
551 global_idx, num_required_cell_areas, reorder_idx);
552
553 // send required points to io processes
554 size_t request_count = recvcounts[comm_size - 1] + rdispls[comm_size - 1];
555 size_t * request_global_idx =
556 xmalloc(request_count * sizeof(*request_global_idx));
558 global_idx, sendcounts, sdispls+1,
559 request_global_idx, recvcounts, rdispls,
560 sizeof(*global_idx), YAC_MPI_SIZE_T, comm);
561 free(global_idx);
562
563 // pack requested cell areas
564 double * requested_cell_areas =
565 xmalloc(request_count * sizeof(*requested_cell_areas));
566 size_t global_idx_offset =
567 ((long long)io_rank_idx * (long long)dist_cell_areas_global_size) /
568 (long long)num_io_ranks;
569 for (size_t i = 0; i < request_count; ++i)
570 requested_cell_areas[i] =
571 dist_cell_areas[request_global_idx[i] - global_idx_offset];
572 free(request_global_idx);
573 free(dist_cell_areas);
574
575 // return cell areas
576 double * temp_cell_areas =
577 xmalloc(num_required_cell_areas * sizeof(*temp_cell_areas));
579 requested_cell_areas, recvcounts, rdispls,
580 temp_cell_areas, sendcounts, sdispls+1,
581 sizeof(*requested_cell_areas), MPI_DOUBLE, comm);
582 free(requested_cell_areas);
583
584 yac_free_comm_buffers(sendcounts, recvcounts, sdispls, rdispls);
585
586 // unpack cell areas
587 for (size_t i = 0; i < num_required_cell_areas; ++i)
588 cell_areas[reorder_idx[i]] = temp_cell_areas[i];
589
590 free(temp_cell_areas);
591 free(reorder_idx);
592}
593
594static double * get_cell_areas(
595 struct yac_const_basic_grid_data * basic_grid_data,
596 char const * type, struct yac_spmap_cell_area_config cell_area_config,
597 MPI_Comm comm, size_t const * required_points, size_t num_required_points) {
598
599 size_t num_cells = basic_grid_data->count[YAC_LOC_CELL];
600
601 // determine which cell areas are actually required
602 int * required_cell_areas = xcalloc(num_cells, sizeof(*required_cell_areas));
603 for (size_t i = 0; i < num_required_points; ++i)
604 required_cell_areas[required_points[i]] = 1;
605
606 // compute and scale cell areas
607 double * cell_areas = xmalloc(num_cells * sizeof(*cell_areas));
608
610 (cell_area_config.cell_area_provider == YAC_INTERP_SPMAP_CELL_AREA_YAC) ||
611 (cell_area_config.cell_area_provider == YAC_INTERP_SPMAP_CELL_AREA_FILE),
612 "ERROR(get_cell_areas): unsupported %s cell area origin", type);
613
614 switch (cell_area_config.cell_area_provider) {
615
616 default:
618
619 double area_scale =
620 cell_area_config.sphere_radius * cell_area_config.sphere_radius;
622 basic_grid_data, required_cell_areas, area_scale, type,
623 cell_areas);
624 break;
625 }
627
629 basic_grid_data, cell_area_config.file_config, comm,
630 required_cell_areas, type, cell_areas);
631 }
632 };
633
634 free(required_cell_areas);
635
636 return cell_areas;
637}
638
639static void scale_weights(
640 struct yac_interp_grid * interp_grid,
641 struct yac_spmap_scale_config scale_config,
642 size_t num_src_points, size_t const * src_points,
643 size_t const * num_tgt_per_src, size_t const * tgt_points,
644 size_t total_num_weights, double * weights) {
645
647 (scale_config.type == YAC_INTERP_SPMAP_NONE) ||
648 (scale_config.type == YAC_INTERP_SPMAP_SRCAREA) ||
649 (scale_config.type == YAC_INTERP_SPMAP_INVTGTAREA) ||
650 (scale_config.type == YAC_INTERP_SPMAP_FRACAREA),
651 "ERROR(scale_weights): invalid scale_type")
652
653 // if there is no scaling
654 if (scale_config.type == YAC_INTERP_SPMAP_NONE) return;
655
656 // get cell areas if they are required
657 double const * src_cell_areas =
658 ((scale_config.type == YAC_INTERP_SPMAP_SRCAREA) ||
659 (scale_config.type == YAC_INTERP_SPMAP_FRACAREA))?
662 "source", scale_config.src, yac_interp_grid_get_MPI_Comm(interp_grid),
663 src_points, total_num_weights):NULL;
664 double const * tgt_cell_areas =
665 ((scale_config.type == YAC_INTERP_SPMAP_INVTGTAREA) ||
666 (scale_config.type == YAC_INTERP_SPMAP_FRACAREA))?
669 "target", scale_config.tgt, yac_interp_grid_get_MPI_Comm(interp_grid),
670 tgt_points, total_num_weights):NULL;
671
672#define SCALE_WEIGHTS( \
673 SRC_CELL_AREA, TGT_CELL_AREA) \
674{ \
675 for (size_t i = 0, offset = 0; i < num_src_points; ++i) { \
676 size_t curr_num_tgt = num_tgt_per_src[i]; \
677 for (size_t j = 0; j < curr_num_tgt; ++j, ++offset) { \
678 weights[offset] *= SRC_CELL_AREA / TGT_CELL_AREA; \
679 } \
680 } \
681}
682
683 switch (scale_config.type) {
685 SCALE_WEIGHTS(src_cell_areas[src_points[offset]], 1.0)
686 break;
688 SCALE_WEIGHTS(1.0, tgt_cell_areas[tgt_points[offset]])
689 break;
690 default:
693 src_cell_areas[src_points[offset]], tgt_cell_areas[tgt_points[offset]])
694 break;
695 }
696
697 free((void*)src_cell_areas);
698 free((void*)tgt_cell_areas);
699}
700
701static void spread_src_data(
702 struct yac_interp_grid * interp_grid, double spread_distance,
704 struct yac_spmap_scale_config scale_config,
705 size_t num_src_points, size_t ** src_points_,
706 size_t ** tgt_result_points_, double ** weights_,
707 size_t * total_num_weights_) {
708
709 // shortcut in case there is only a single "1.0" weight for all source points
710 if ((spread_distance <= 0) && (scale_config.type == YAC_INTERP_SPMAP_NONE)) {
711
712 *weights_ = NULL;
713 *total_num_weights_ = num_src_points;
714 return;
715 }
716
717 size_t * src_points = *src_points_;
718 size_t * tgt_result_points = *tgt_result_points_;
719
720 size_t * num_tgt_per_src =
721 xmalloc(num_src_points * sizeof(*num_tgt_per_src));
722 size_t total_num_weights;
723
724 // search for additional target points if spread distance is bigger than 0.0
725 if (spread_distance > 0.0) {
726
727 struct sin_cos_angle inc_angle =
728 sin_cos_angle_new(sin(spread_distance), cos(spread_distance));
729 yac_const_coordinate_pointer tgt_field_coords =
731
732 struct bounding_circle * search_bnd_circles =
733 xmalloc(num_src_points * sizeof(*search_bnd_circles));
734 for (size_t i = 0; i < num_src_points; ++i) {
735 memcpy(
736 search_bnd_circles[i].base_vector,
737 tgt_field_coords[tgt_result_points[i]], sizeof(*tgt_field_coords));
738 search_bnd_circles[i].inc_angle = inc_angle;
739 search_bnd_circles[i].sq_crd = DBL_MAX;
740 }
741 size_t * spread_tgt_result_points = NULL;
742
743 // do bounding circle search around found tgt points
745 interp_grid, search_bnd_circles, num_src_points,
746 &spread_tgt_result_points, num_tgt_per_src);
747 free(search_bnd_circles);
748
749 // remove target points which exceed the spread distance and only keep
750 // targets that are connected to the original target point or other
751 // target that have already been selected
752 total_num_weights =
754 interp_grid, spread_distance, num_src_points, tgt_result_points,
755 num_tgt_per_src, spread_tgt_result_points);
756 free((void*)tgt_result_points);
757 tgt_result_points =
758 xrealloc(
759 spread_tgt_result_points,
760 total_num_weights * sizeof(*spread_tgt_result_points));
761
762 // adjust src_points (one source per target)
763 size_t * new_src_points =
764 xmalloc(total_num_weights * sizeof(*new_src_points));
765 for (size_t i = 0, offset = 0; i < num_src_points; ++i)
766 for (size_t j = 0, curr_src_point = src_points[i];
767 j < num_tgt_per_src[i]; ++j, ++offset)
768 new_src_points[offset] = curr_src_point;
769 free((void*)src_points);
770 src_points = new_src_points;
771
772 } else {
773
774 for (size_t i = 0; i < num_src_points; ++i) num_tgt_per_src[i] = 1;
775 total_num_weights = num_src_points;
776 }
777
778 // compute weights
779 double * weights =
781 interp_grid, weight_type, num_src_points, src_points,
782 num_tgt_per_src, total_num_weights, tgt_result_points);
783
784 // scale weights
786 interp_grid, scale_config, num_src_points, src_points, num_tgt_per_src,
787 tgt_result_points, total_num_weights, weights);
788
789 free(num_tgt_per_src);
790
791 // set return values
792 *tgt_result_points_ = tgt_result_points;
793 *src_points_ = src_points;
794 *weights_ = weights;
795 *total_num_weights_ = total_num_weights;
796}
797
798static size_t do_search_spmap (struct interp_method * method,
799 struct yac_interp_grid * interp_grid,
800 size_t * tgt_points, size_t count,
801 struct yac_interp_weights * weights) {
802
804 yac_interp_grid_get_num_src_fields(interp_grid) == 1,
805 "ERROR(do_search_spmap): invalid number of source fields")
808 "ERROR(do_search_spmap): "
809 "invalid source field location (has to be YAC_LOC_CELL)")
812 "ERROR(do_search_spmap): "
813 "invalid target field location (has to be YAC_LOC_CELL)")
814
815 // get coordinates of all source points
816 size_t * src_points;
817 size_t num_src_points;
819 interp_grid, 0, &src_points, &num_src_points);
820 yac_coordinate_pointer src_coords = xmalloc(num_src_points * sizeof(*src_coords));
822 interp_grid, src_points, num_src_points, 0, src_coords);
823
824 // search for matching tgt points
825 size_t * tgt_result_points =
826 xmalloc(num_src_points * sizeof(*tgt_result_points));
828 interp_grid, src_coords, num_src_points, 1, tgt_result_points,
829 ((struct interp_method_spmap*)method)->max_search_distance);
830
831 free(src_coords);
832
833 // remove source points for which matching target point was found
834 {
835 size_t new_num_src_points = 0;
836 for (size_t i = 0; i < num_src_points; ++i) {
837 if (tgt_result_points[i] != SIZE_MAX) {
838 if (i != new_num_src_points) {
839 src_points[new_num_src_points] = src_points[i];
840 tgt_result_points[new_num_src_points] =
841 tgt_result_points[i];
842 }
843 ++new_num_src_points;
844 }
845 }
846 num_src_points = new_num_src_points;
847 }
848
849 // spread the data from each source point to multiple target points
850 // (weight_data is set to NULL if no spreading was applied)
851 double * weight_data;
852 size_t total_num_weights;
854 interp_grid, ((struct interp_method_spmap*)method)->spread_distance,
855 ((struct interp_method_spmap*)method)->weight_type,
856 ((struct interp_method_spmap*)method)->scale_config, num_src_points,
857 &src_points, &tgt_result_points, &weight_data, &total_num_weights);
858
859 // relocate source-target-point-pairs to dist owners of the respective
860 // target points
861 size_t result_count = total_num_weights;
862 int to_tgt_owner = 1;
864 interp_grid, to_tgt_owner,
865 0, &src_points, &tgt_result_points, &weight_data, &result_count);
866 total_num_weights = result_count;
867
868 // sort source-target-point-pairs by target points
870 tgt_result_points, result_count, src_points, weight_data);
871
872 // generate num_src_per_tgt and compact tgt_result_points
873 size_t * num_src_per_tgt = xmalloc(result_count * sizeof(*num_src_per_tgt));
874 size_t num_unique_tgt_result_points = 0;
875 for (size_t i = 0; i < result_count;) {
876 size_t prev_i = i;
877 size_t curr_tgt = tgt_result_points[i];
878 while ((i < result_count) && (curr_tgt == tgt_result_points[i])) ++i;
879 num_src_per_tgt[num_unique_tgt_result_points] = i - prev_i;
880 tgt_result_points[num_unique_tgt_result_points] = curr_tgt;
881 ++num_unique_tgt_result_points;
882 }
883 result_count = num_unique_tgt_result_points;
884 num_src_per_tgt =
885 xrealloc(num_src_per_tgt, result_count * sizeof(*num_src_per_tgt));
886 tgt_result_points =
887 xrealloc(tgt_result_points, result_count * sizeof(*tgt_result_points));
888
889 // match tgt_result_points with available target points
890 qsort(tgt_points, count, sizeof(*tgt_points), compare_size_t);
891 int * reorder_flag = xmalloc(count * sizeof(*tgt_points));
892 {
893 size_t j = 0;
894 for (size_t i = 0; i < result_count; ++i) {
895 size_t curr_result_tgt = tgt_result_points[i];
896 while ((j < count) && (tgt_points[j] < curr_result_tgt))
897 reorder_flag[j++] = 1;
899 (j < count) && (curr_result_tgt == tgt_points[j]),
900 "ERROR(do_search_spmap): "
901 "required target points already in use or not available")
902 reorder_flag[j++] = 0;
903 }
904 for (; j < count; ++j) reorder_flag[j] = 1;
905 }
906
907 // sort used target points to the beginning of the array
908 yac_quicksort_index_int_size_t(reorder_flag, count, tgt_points);
909 free(reorder_flag);
910
911 struct remote_point * srcs =
913 interp_grid, 0, src_points, total_num_weights);
914 struct remote_points tgts = {
915 .data =
917 interp_grid, tgt_result_points, result_count),
918 .count = result_count};
919 free(tgt_result_points);
920
921 // store results
922 if (weight_data == NULL)
924 weights, &tgts, num_src_per_tgt, srcs);
925 else
927 weights, &tgts, num_src_per_tgt, srcs, weight_data);
928
929 free(weight_data);
930 free(src_points);
931 free(tgts.data);
932 free(srcs);
933 free(num_src_per_tgt);
934
935 return result_count;
936}
937
939 struct yac_spmap_cell_area_config * cell_area_config, char const * type) {
940
942 (cell_area_config->cell_area_provider == YAC_INTERP_SPMAP_CELL_AREA_YAC) ||
943 (cell_area_config->cell_area_provider == YAC_INTERP_SPMAP_CELL_AREA_FILE),
944 "ERROR(check_cell_area_config): invalid %s cell_area_provider type", type);
945
946 switch (cell_area_config->cell_area_provider) {
947 default:
950 cell_area_config->sphere_radius > 0.0,
951 "ERROR(check_cell_area_config): invalid %s sphere_radius "
952 "(has to be >= 0.0", type);
953 cell_area_config->file_config.filename = NULL;
954 cell_area_config->file_config.varname = NULL;
955 cell_area_config->file_config.min_global_id = XT_INT_MAX;
956 break;
957 }
959#ifndef YAC_NETCDF_ENABLED
960 die(
961 "ERROR(interp_method_spmap::check_cell_area_config): "
962 "YAC is built without the NetCDF support");
963#endif
965 (cell_area_config->file_config.filename != NULL) &&
966 (strlen(cell_area_config->file_config.filename) > 0),
967 "ERROR(check_cell_area_config): invalid filename for %s areas", type);
969 (cell_area_config->file_config.varname!= NULL) &&
970 (strlen(cell_area_config->file_config.varname) > 0),
971 "ERROR(check_cell_area_config): invalid varname for %s areas", type);
972 cell_area_config->sphere_radius = 0.0;
973 cell_area_config->file_config.filename =
974 strdup(cell_area_config->file_config.filename);
975 cell_area_config->file_config.varname =
976 strdup(cell_area_config->file_config.varname);
977 break;
978 }
979 }
980}
981
983 double spread_distance, double max_search_distance,
985 struct yac_spmap_scale_config scale_config) {
986
987 struct interp_method_spmap * method = xmalloc(1 * sizeof(*method));
988
991 method->max_search_distance =
993 method->weight_type = weight_type;
994 method->scale_config = scale_config;
995
997 (spread_distance >= 0.0) && (spread_distance <= M_PI_2),
998 "ERROR(yac_interp_method_spmap_new): invalid spread_distance "
999 "(has to be >= 0 and <= PI/2")
1000
1001 YAC_ASSERT(
1002 (max_search_distance >= 0.0) && (max_search_distance <= M_PI),
1003 "ERROR(yac_interp_method_spmap_new): invalid max_search_distance "
1004 "(has to be >= 0 and <= PI")
1005
1006 check_cell_area_config(&(method->scale_config.src), "source");
1007 check_cell_area_config(&(method->scale_config.tgt), "target");
1008
1009 return (struct interp_method*)method;
1010}
1011
1012static void delete_spmap(struct interp_method * method) {
1013
1014 struct interp_method_spmap * method_spmap =
1015 (struct interp_method_spmap*)(method);
1016
1017 free((void*)method_spmap->scale_config.src.file_config.filename);
1018 free((void*)method_spmap->scale_config.src.file_config.varname);
1019 free((void*)method_spmap->scale_config.tgt.file_config.filename);
1020 free((void*)method_spmap->scale_config.tgt.file_config.varname);
1021 free(method);
1022}
double yac_huiliers_area(struct yac_grid_cell cell)
Area calculation on a unit sphere taken from ESMF based on L'Huilier's Theorem.
Definition area.c:396
Structs and interfaces for area calculations.
#define YAC_AREA_TOL
Definition area.h:23
enum yac_interp_ncc_weight_type weight_type
#define UNUSED(x)
Definition core.h:73
void yac_const_basic_grid_data_get_grid_cell(struct yac_const_basic_grid_data *grid_data, size_t cell_idx, struct yac_grid_cell *buffer_cell)
Definition dist_grid.c:2649
size_t const *const const_size_t_pointer
#define ENSURE_ARRAY_SIZE(arrayp, curr_array_size, req_size)
#define yac_angle_tol
Definition geometry.h:34
static double get_vector_angle(double const a[3], double const b[3])
Definition geometry.h:364
static struct sin_cos_angle sin_cos_angle_new(double sin, double cos)
Definition geometry.h:405
void yac_init_grid_cell(struct yac_grid_cell *cell)
Definition grid_cell.c:14
void yac_free_grid_cell(struct yac_grid_cell *cell)
Definition grid_cell.c:44
enum yac_location yac_interp_grid_get_tgt_field_location(struct yac_interp_grid *interp_grid)
void yac_interp_grid_relocate_src_tgt_pairs(struct yac_interp_grid *interp_grid, int to_tgt_owner, size_t src_field_idx, size_t **src_points, size_t **tgt_points, double **weights, size_t *count)
size_t yac_interp_grid_get_num_src_fields(struct yac_interp_grid *interp_grid)
void yac_interp_grid_do_bnd_circle_search_tgt(struct yac_interp_grid *interp_grid, const_bounding_circle_pointer bnd_circles, size_t count, size_t **tgt_cells, size_t *num_tgt_per_bnd_circle)
MPI_Comm yac_interp_grid_get_MPI_Comm(struct yac_interp_grid *interp_grid)
struct remote_point * yac_interp_grid_get_tgt_remote_points(struct yac_interp_grid *interp_grid, size_t *tgt_points, size_t count)
void yac_interp_grid_get_src_points(struct yac_interp_grid *interp_grid, size_t src_field_idx, size_t **src_indices, size_t *count)
Definition interp_grid.c:87
enum yac_location yac_interp_grid_get_src_field_location(struct yac_interp_grid *interp_grid, size_t src_field_idx)
void yac_interp_grid_get_src_coordinates(struct yac_interp_grid *interp_grid, size_t *src_points, size_t count, size_t src_field_idx, yac_coordinate_pointer src_coordinates)
struct remote_point * yac_interp_grid_get_src_remote_points(struct yac_interp_grid *interp_grid, size_t src_field_idx, size_t *src_points, size_t count)
yac_const_coordinate_pointer yac_interp_grid_get_src_field_coords(struct yac_interp_grid *interp_grid, size_t src_field_idx)
yac_const_coordinate_pointer yac_interp_grid_get_tgt_field_coords(struct yac_interp_grid *interp_grid)
struct yac_const_basic_grid_data * yac_interp_grid_get_basic_grid_data_tgt(struct yac_interp_grid *interp_grid)
void yac_interp_grid_do_nnn_search_tgt(struct yac_interp_grid *interp_grid, yac_coordinate_pointer search_coords, size_t count, size_t n, size_t *tgt_points, double max_search_distance)
struct yac_const_basic_grid_data * yac_interp_grid_get_basic_grid_data_src(struct yac_interp_grid *interp_grid)
enum callback_type type
static int lists_overlap(size_t *list_a, size_t count_a, size_t *list_b, size_t count_b)
static void delete_spmap(struct interp_method *method)
static int compare_size_t(const void *a, const void *b)
static size_t do_search_spmap(struct interp_method *method, struct yac_interp_grid *interp_grid, size_t *tgt_points, size_t count, struct yac_interp_weights *weights)
static struct interp_method_vtable interp_method_spmap_vtable
static void merge_lists(size_t *list, size_t *list_size, size_t *insert, size_t insert_size)
static void remove_disconnected_points(struct yac_interp_grid *interp_grid, size_t tgt_start_point, size_t *from_tgt_points, size_t *to_tgt_points, size_t *count, int *flag, size_t *temp_cell_edges, size_t **edges_buffer, size_t *edges_buffer_array_size)
static void compute_cell_areas(struct yac_const_basic_grid_data *basic_grid_data, int *required_cell_areas, double area_scale, char const *type, double *cell_areas)
static void check_cell_area_config(struct yac_spmap_cell_area_config *cell_area_config, char const *type)
static void spread_src_data(struct yac_interp_grid *interp_grid, double spread_distance, enum yac_interp_spmap_weight_type weight_type, struct yac_spmap_scale_config scale_config, size_t num_src_points, size_t **src_points_, size_t **tgt_result_points_, double **weights_, size_t *total_num_weights_)
#define SCALE_WEIGHTS(SRC_CELL_AREA, TGT_CELL_AREA)
static void read_cell_areas(struct yac_const_basic_grid_data *basic_grid_data, struct yac_spmap_cell_area_file_config file_config, MPI_Comm comm, int *required_cell_areas, char const *type, double *cell_areas)
struct interp_method * yac_interp_method_spmap_new(double spread_distance, double max_search_distance, enum yac_interp_spmap_weight_type weight_type, struct yac_spmap_scale_config scale_config)
static size_t check_tgt_result_points(struct yac_interp_grid *interp_grid, double spread_distance, size_t num_src_points, size_t const *const tgt_result_points, size_t *num_tgt_per_src, size_t *spread_tgt_result_points)
static void check_spread_distance(yac_const_coordinate_pointer tgt_field_coords, double max_distance, size_t tgt_start_point, size_t *tgt_points, size_t *count)
static void dist_read_cell_areas(char const *filename, char const *varname, MPI_Comm comm, int *io_ranks, int io_rank_idx, int num_io_ranks, double **dist_cell_areas, size_t *dist_cell_areas_global_size)
static double * compute_weights(struct yac_interp_grid *interp_grid, enum yac_interp_spmap_weight_type weight_type, size_t num_src_points, size_t const *const src_points, size_t const *const num_tgt_per_src, size_t total_num_tgt, size_t const *const tgt_result_points)
static void scale_weights(struct yac_interp_grid *interp_grid, struct yac_spmap_scale_config scale_config, size_t num_src_points, size_t const *src_points, size_t const *num_tgt_per_src, size_t const *tgt_points, size_t total_num_weights, double *weights)
static double * get_cell_areas(struct yac_const_basic_grid_data *basic_grid_data, char const *type, struct yac_spmap_cell_area_config cell_area_config, MPI_Comm comm, size_t const *required_points, size_t num_required_points)
@ YAC_INTERP_SPMAP_NONE
weights are not scaled
@ YAC_INTERP_SPMAP_INVTGTAREA
@ YAC_INTERP_SPMAP_SRCAREA
@ YAC_INTERP_SPMAP_FRACAREA
@ YAC_INTERP_SPMAP_CELL_AREA_FILE
@ YAC_INTERP_SPMAP_CELL_AREA_YAC
yac_interp_spmap_weight_type
@ YAC_INTERP_SPMAP_AVG
@ YAC_INTERP_SPMAP_DIST
void yac_interp_weights_add_wsum(struct yac_interp_weights *weights, struct remote_points *tgts, size_t *num_src_per_tgt, struct remote_point *srcs, double *w)
void yac_interp_weights_add_sum(struct yac_interp_weights *weights, struct remote_points *tgts, size_t *num_src_per_tgt, struct remote_point *srcs)
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
#define YAC_HANDLE_ERROR(exp)
Definition io_utils.h:30
@ 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
struct sin_cos_angle inc_angle
angle between the middle point and the boundary of the spherical cap
Definition geometry.h:61
double base_vector[3]
Definition geometry.h:59
double sq_crd
Definition geometry.h:64
struct yac_spmap_scale_config scale_config
enum yac_interp_spmap_weight_type weight_type
struct interp_method_vtable * vtable
size_t(* do_search)(struct interp_method *method, struct yac_interp_grid *grid, size_t *tgt_points, size_t count, struct yac_interp_weights *weights)
information (global id and location) about a point that
structure containing the information (global id and location)
struct remote_point * data
double sin
Definition geometry.h:41
double cos
Definition geometry.h:41
const const_int_pointer num_vertices_per_cell
const_size_t_pointer cell_to_edge
const const_yac_int_pointer ids[3]
const_size_t_pointer cell_to_edge_offsets
struct yac_spmap_scale_config::yac_spmap_cell_area_config::yac_spmap_cell_area_file_config file_config
struct yac_spmap_scale_config::yac_spmap_cell_area_config src
struct yac_spmap_scale_config::yac_spmap_cell_area_config tgt
enum yac_interp_spmap_scale_type type
void yac_quicksort_index_size_t_size_t_double(size_t *a, size_t n, size_t *b, double *c)
void yac_quicksort_index_int_size_t(int *a, size_t n, size_t *idx)
void yac_quicksort_index_size_t_size_t(size_t *a, size_t n, size_t *idx)
#define YAC_ASSERT_F(exp, format,...)
Definition yac_assert.h:19
#define die(msg)
Definition yac_assert.h:12
#define YAC_ASSERT(exp, msg)
Definition yac_assert.h:16
void yac_generate_alltoallv_args(int count, size_t const *sendcounts, size_t *recvcounts, size_t *sdispls, size_t *rdispls, MPI_Comm comm)
Definition yac_mpi.c:571
void yac_free_comm_buffers(size_t *sendcounts, size_t *recvcounts, size_t *sdispls, size_t *rdispls)
Definition yac_mpi.c:627
void yac_get_comm_buffers(int count, size_t **sendcounts, size_t **recvcounts, size_t **sdispls, size_t **rdispls, MPI_Comm comm)
Definition yac_mpi.c:596
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)
Definition yac_mpi.c:131
#define yac_mpi_call(call, comm)
#define YAC_MPI_SIZE_T
Xt_int yac_int
Definition yac_types.h:15
double const (* yac_const_coordinate_pointer)[3]
Definition yac_types.h:20
double(* yac_coordinate_pointer)[3]
Definition yac_types.h:19