22 size_t * tgt_points,
size_t count,
43 size_t const * a_ = a, * b_ = b;
45 return (*a_ > *b_) - (*b_ > *a_);
49 size_t * list_a,
size_t count_a,
size_t * list_b,
size_t count_b) {
51 if ((count_a == 0) || (count_b == 0))
return 0;
54 size_t curr_a = SIZE_MAX, curr_b = list_b[0];
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));
67 size_t * list,
size_t * list_size,
size_t * insert,
size_t insert_size) {
69 if (insert_size == 0)
return;
71 size_t new_list_size = *list_size;
72 size_t old_list_size = *list_size;
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;
81 if (new_list_size != old_list_size) {
83 *list_size = new_list_size;
89 size_t tgt_start_point,
size_t * tgt_points,
size_t * count) {
91 double const * start_coord = tgt_field_coords[tgt_start_point];
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];
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) {
118 size_t old_count = *
count;
119 memset(flag, 0, old_count *
sizeof(*flag));
121 for (
size_t i = 0; i < old_count; ++i) {
122 if (from_tgt_points[i] == tgt_start_point) {
128 size_t * edges = *edges_buffer;
129 size_t num_edges = 0;
130 size_t edges_array_size = *edges_buffer_array_size;
137 num_edges = curr_num_edges;
138 memcpy(edges, curr_edges, num_edges *
sizeof(*edges));
146 for (
size_t i = 0; i < old_count; ++i) {
148 if (flag[i])
continue;
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);
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);
164 }
while (change_flag);
166 *edges_buffer = edges;
167 *edges_buffer_array_size = edges_array_size;
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];
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) {
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];
186 int * flag =
xmalloc(max_num_tgt_per_src *
sizeof(*flag));
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 =
202 for (
size_t i = 0, old_offset = 0; i < num_src_points; ++i) {
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];
211 tgt_field_coords, spread_distance, tgt_result_points[i],
212 old_results, num_tgt_per_src + i);
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]];
221 cell_edge_buffer, cell_edge_buffer_array_size,
222 max_num_vertice_per_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);
230 new_offset += num_tgt_per_src[i];
234 free(cell_edge_buffer);
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) {
246 double * weights =
xmalloc(total_num_tgt *
sizeof(*weights));
250 "ERROR(do_search_spmap): invalid weight_type")
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;
262 weights[offset] = 1.0;
276 for (
size_t i = 0, offset = 0; i < num_src_points; ++i) {
278 size_t curr_num_tgt = num_tgt_per_src[i];
280 if (curr_num_tgt == 0)
continue;
282 double * curr_weights = weights + offset;
284 if (curr_num_tgt > 1) {
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;
293 for (
size_t j = 0; j < curr_num_tgt; ++j) {
297 (
double*)curr_src_coord,
298 (
double*)tgt_field_coords[curr_result_points[j]]);
301 for (
size_t k = 0; k < curr_num_tgt; ++k) curr_weights[k] = 0.0;
302 curr_weights[j] = 1.0;
306 curr_weights[j] = 1.0 / distance;
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;
317 for (
size_t j = 0; j < curr_num_tgt; ++j) curr_weights[j] *= scale;
333 int * required_cell_areas,
double area_scale,
char const *
type,
334 double * cell_areas) {
341 for (
size_t i = 0; i < num_cells; ++i) {
342 if (required_cell_areas[i]) {
347 "ERROR(get_cell_areas): "
348 "area of %s cell (global id %"XT_INT_FMT
") is close to zero (%e)",
350 cell_areas[i] = cell_area;
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) {
364#ifndef YAC_NETCDF_ENABLED
373 *dist_cell_areas = NULL;
374 *dist_cell_areas_global_size = 0;
377 "ERROR(interp_method_spmap::dist_read_cell_areas): "
378 "YAC is built without the NetCDF support");
381 if ((io_rank_idx >= 0) && (io_rank_idx < num_io_ranks)) {
393 int dimids[NC_MAX_VAR_DIMS];
395 nc_inq_var(ncid, varid, NULL, NULL, &ndims, dimids, NULL));
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);
406 *dist_cell_areas_global_size = 1;
407 for (
int i = 0; i < ndims; ++i) {
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];
421 size_t start[2], count[2], offset, read_size;
424 ((
long long)*dist_cell_areas_global_size * (
long long)io_rank_idx) /
425 (
long long)num_io_ranks);
428 ((
long long)*dist_cell_areas_global_size *
429 (
long long)(io_rank_idx+1)) / (
long long)num_io_ranks) - local_start;
431 start[0] = local_start;
432 count[0] = local_count;
434 read_size = local_count;
436 start[0] = local_start / dimlen[1];
438 (local_start + local_count + dimlen[1] - 1) / dimlen[1] - start[0];
440 count[1] = dimlen[1];
441 offset = local_start - start[0] * dimlen[1];
442 read_size = count[0] * count[1];
446 *dist_cell_areas =
xmalloc(read_size *
sizeof(**dist_cell_areas));
448 nc_get_vara_double(ncid, varid, start, count, *dist_cell_areas));
453 *dist_cell_areas, *dist_cell_areas + offset,
454 local_count *
sizeof(**dist_cell_areas));
460 *dist_cell_areas =
xmalloc(1*
sizeof(**dist_cell_areas));
461 *dist_cell_areas_global_size = 0;
466 dist_cell_areas_global_size, 1,
YAC_MPI_SIZE_T, io_ranks[0], comm),
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) {
481 int local_is_io, * io_ranks, num_io_ranks;
484 int comm_rank, comm_size;
489 while ((io_rank_idx < num_io_ranks) &&
490 (comm_rank != io_ranks[io_rank_idx]))
493 !local_is_io || (io_rank_idx < num_io_ranks),
494 "ERROR(read_cell_areas): unable to determine io_rank_idx");
496 double * dist_cell_areas;
497 size_t dist_cell_areas_global_size;
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);
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;
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));
514 size_t * sendcounts, * recvcounts, * sdispls, * rdispls;
516 1, &sendcounts, &recvcounts, &sdispls, &rdispls, comm);
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);
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]]++;
547 1, sendcounts, recvcounts, sdispls, rdispls, comm);
551 global_idx, num_required_cell_areas, reorder_idx);
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,
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);
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);
587 for (
size_t i = 0; i < num_required_cell_areas; ++i)
588 cell_areas[reorder_idx[i]] = temp_cell_areas[i];
590 free(temp_cell_areas);
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) {
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;
607 double * cell_areas =
xmalloc(num_cells *
sizeof(*cell_areas));
612 "ERROR(get_cell_areas): unsupported %s cell area origin",
type);
614 switch (cell_area_config.cell_area_provider) {
620 cell_area_config.sphere_radius * cell_area_config.sphere_radius;
622 basic_grid_data, required_cell_areas, area_scale,
type,
629 basic_grid_data, cell_area_config.file_config, comm,
630 required_cell_areas,
type, cell_areas);
634 free(required_cell_areas);
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) {
651 "ERROR(scale_weights): invalid scale_type")
657 double const * src_cell_areas =
663 src_points, total_num_weights):NULL;
664 double const * tgt_cell_areas =
670 tgt_points, total_num_weights):NULL;
672#define SCALE_WEIGHTS( \
673 SRC_CELL_AREA, TGT_CELL_AREA) \
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; \
683 switch (scale_config.
type) {
693 src_cell_areas[src_points[offset]], tgt_cell_areas[tgt_points[offset]])
697 free((
void*)src_cell_areas);
698 free((
void*)tgt_cell_areas);
705 size_t num_src_points,
size_t ** src_points_,
706 size_t ** tgt_result_points_,
double ** weights_,
707 size_t * total_num_weights_) {
713 *total_num_weights_ = num_src_points;
717 size_t * src_points = *src_points_;
718 size_t * tgt_result_points = *tgt_result_points_;
720 size_t * num_tgt_per_src =
721 xmalloc(num_src_points *
sizeof(*num_tgt_per_src));
722 size_t total_num_weights;
725 if (spread_distance > 0.0) {
733 xmalloc(num_src_points *
sizeof(*search_bnd_circles));
734 for (
size_t i = 0; i < num_src_points; ++i) {
737 tgt_field_coords[tgt_result_points[i]],
sizeof(*tgt_field_coords));
739 search_bnd_circles[i].
sq_crd = DBL_MAX;
741 size_t * spread_tgt_result_points = NULL;
745 interp_grid, search_bnd_circles, num_src_points,
746 &spread_tgt_result_points, num_tgt_per_src);
747 free(search_bnd_circles);
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);
759 spread_tgt_result_points,
760 total_num_weights *
sizeof(*spread_tgt_result_points));
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;
774 for (
size_t i = 0; i < num_src_points; ++i) num_tgt_per_src[i] = 1;
775 total_num_weights = num_src_points;
781 interp_grid,
weight_type, num_src_points, src_points,
782 num_tgt_per_src, total_num_weights, tgt_result_points);
786 interp_grid, scale_config, num_src_points, src_points, num_tgt_per_src,
787 tgt_result_points, total_num_weights, weights);
789 free(num_tgt_per_src);
792 *tgt_result_points_ = tgt_result_points;
793 *src_points_ = src_points;
795 *total_num_weights_ = total_num_weights;
800 size_t * tgt_points,
size_t count,
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)")
817 size_t num_src_points;
819 interp_grid, 0, &src_points, &num_src_points);
822 interp_grid, src_points, num_src_points, 0, src_coords);
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,
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];
843 ++new_num_src_points;
846 num_src_points = new_num_src_points;
851 double * weight_data;
852 size_t total_num_weights;
857 &src_points, &tgt_result_points, &weight_data, &total_num_weights);
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;
870 tgt_result_points, result_count, src_points, weight_data);
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;) {
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;
883 result_count = num_unique_tgt_result_points;
885 xrealloc(num_src_per_tgt, result_count *
sizeof(*num_src_per_tgt));
887 xrealloc(tgt_result_points, result_count *
sizeof(*tgt_result_points));
891 int * reorder_flag =
xmalloc(count *
sizeof(*tgt_points));
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;
904 for (; j < count; ++j) reorder_flag[j] = 1;
913 interp_grid, 0, src_points, total_num_weights);
917 interp_grid, tgt_result_points, result_count),
918 .count = result_count};
919 free(tgt_result_points);
922 if (weight_data == NULL)
924 weights, &tgts, num_src_per_tgt, srcs);
927 weights, &tgts, num_src_per_tgt, srcs, weight_data);
933 free(num_src_per_tgt);
939 struct yac_spmap_cell_area_config * cell_area_config,
char const *
type) {
944 "ERROR(check_cell_area_config): invalid %s cell_area_provider type",
type);
946 switch (cell_area_config->cell_area_provider) {
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;
959#ifndef YAC_NETCDF_ENABLED
961 "ERROR(interp_method_spmap::check_cell_area_config): "
962 "YAC is built without the NetCDF support");
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);
983 double spread_distance,
double max_search_distance,
998 "ERROR(yac_interp_method_spmap_new): invalid spread_distance "
999 "(has to be >= 0 and <= PI/2")
1003 "ERROR(yac_interp_method_spmap_new): invalid max_search_distance "
1004 "(has to be >= 0 and <= PI")
double yac_huiliers_area(struct yac_grid_cell cell)
Area calculation on a unit sphere taken from ESMF based on L'Huilier's Theorem.
Structs and interfaces for area calculations.
enum yac_interp_ncc_weight_type weight_type
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)
size_t const *const const_size_t_pointer
#define ENSURE_ARRAY_SIZE(arrayp, curr_array_size, req_size)
static double get_vector_angle(double const a[3], double const b[3])
static struct sin_cos_angle sin_cos_angle_new(double sin, double cos)
void yac_init_grid_cell(struct yac_grid_cell *cell)
void yac_free_grid_cell(struct yac_grid_cell *cell)
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)
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)
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
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_)
void yac_nc_inq_varid(int ncid, char const *name, int *varidp)
void yac_nc_open(const char *path, int omode, int *ncidp)
#define YAC_HANDLE_ERROR(exp)
#define xrealloc(ptr, size)
#define xcalloc(nmemb, size)
struct sin_cos_angle inc_angle
angle between the middle point and the boundary of the spherical cap
struct yac_spmap_scale_config scale_config
enum yac_interp_spmap_weight_type weight_type
double max_search_distance
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
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,...)
#define YAC_ASSERT(exp, msg)
void yac_generate_alltoallv_args(int count, size_t const *sendcounts, size_t *recvcounts, size_t *sdispls, size_t *rdispls, MPI_Comm comm)
void yac_free_comm_buffers(size_t *sendcounts, size_t *recvcounts, size_t *sdispls, size_t *rdispls)
void yac_get_comm_buffers(int count, size_t **sendcounts, size_t **recvcounts, size_t **sdispls, size_t **rdispls, MPI_Comm comm)
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)
#define yac_mpi_call(call, comm)
double const (* yac_const_coordinate_pointer)[3]
double(* yac_coordinate_pointer)[3]