YAC 3.13.0
Yet Another Coupler
Loading...
Searching...
No Matches
interp_method_ncc.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 <string.h>
11
14
15static size_t do_search_ncc(struct interp_method * method,
16 struct yac_interp_grid * interp_grid,
17 size_t * tgt_points, size_t count,
18 struct yac_interp_weights * weights,
19 int * interpolation_compelte);
20static void delete_ncc(struct interp_method * method);
21
22typedef int (*func_compute_weights)(
23 double[3], size_t, yac_const_coordinate_pointer, int *, double *);
24
25static struct interp_method_vtable
29
35
36// determines for a given coordinate which corner of a specific
37// cell is the closest one
39 struct yac_const_basic_grid_data * grid_data, size_t cell,
40 double const coord[3]) {
41
42 size_t const * cell_to_vertex =
43 grid_data->cell_to_vertex + grid_data->cell_to_vertex_offsets[cell];
44 yac_const_coordinate_pointer corner_coords =
45 grid_data->vertex_coordinates;
46 size_t num_vertices = grid_data->num_vertices_per_cell[cell];
47 yac_int const * vertex_ids = grid_data->ids[YAC_LOC_CORNER];
48
49 double min_distance = DBL_MAX;
50 size_t min_vertex = SIZE_MAX;
51 yac_int min_vertex_id = YAC_INT_MAX;
52
53 for (size_t i = 0; i < num_vertices; ++i) {
54 size_t curr_vertex = cell_to_vertex[i];
55 double curr_distance =
56 get_vector_angle(coord, corner_coords[curr_vertex]);
57 if (curr_distance < min_distance) {
58 min_distance = curr_distance;
59 min_vertex = curr_vertex;
60 min_vertex_id = vertex_ids[curr_vertex];
61 } else if ((curr_distance == min_distance) &&
62 (vertex_ids[curr_vertex] < min_vertex_id)) {
63
64 min_vertex = curr_vertex;
65 min_vertex_id = vertex_ids[curr_vertex];
66 }
67 }
68 return min_vertex;
69}
70
71static size_t do_search_ncc (struct interp_method * method,
72 struct yac_interp_grid * interp_grid,
73 size_t * tgt_points, size_t count,
74 struct yac_interp_weights * weights,
75 int * interpolation_complete) {
76
77 if (*interpolation_complete) return 0;
78
79 struct interp_method_ncc * method_ncc = (struct interp_method_ncc *)method;
80
82 yac_interp_grid_get_num_src_fields(interp_grid) == 1,
83 "ERROR(do_search_ncc): invalid number of source fields")
84
87 "ERROR(do_search_ncc): unsupported source field location type "
88 "(has to be CELL)")
89
90 // get coordinates of target points
91 yac_coordinate_pointer tgt_coords = xmalloc(count * sizeof(*tgt_coords));
93 interp_grid, tgt_points, count, tgt_coords);
94
95 size_t * size_t_buffer = xmalloc(3 * count * sizeof(*size_t_buffer));
96 size_t * src_cells = size_t_buffer;
97 size_t * src_corners = src_cells;
98 size_t * reorder_idx = size_t_buffer + count;
99 size_t * interp_flag = size_t_buffer + 2 * count;
100
101 // get matching source cells for all target points
103 interp_grid, tgt_coords, count, src_cells);
104
105 // sort target points, for which we found a source cell, to the beginning of
106 // the array
107 for (size_t i = 0; i < count; ++i) reorder_idx[i] = i;
108 yac_quicksort_index_size_t_size_t(src_cells, count, reorder_idx);
109
110 size_t result_count = 0;
111 for (result_count = 0; result_count < count; ++result_count)
112 if (src_cells[result_count] == SIZE_MAX) break;
113 for (size_t i = result_count; i < count; ++i)
114 interp_flag[reorder_idx[i]] = SIZE_MAX;
115
116 // for each target point find the closest corner of its source cell
117 {
118 struct yac_const_basic_grid_data * src_grid_data =
120 for (size_t i = 0; i < result_count; ++i)
121 src_corners[i] =
123 src_grid_data, src_cells[i], tgt_coords[reorder_idx[i]]);
124 }
125
126 // sort source corners
128 src_corners, result_count, reorder_idx);
129
130 // count number of unique source corners
131 size_t num_unique_src_corners = 0;
132 for (size_t i = 0, prev_src_corner = SIZE_MAX; i < result_count; ++i) {
133 size_t curr_src_corner = src_corners[i];
134 if (prev_src_corner != curr_src_corner) {
135 prev_src_corner = curr_src_corner;
136 ++num_unique_src_corners;
137 }
138 }
139
140 // count number of target points per unique source corner and
141 // remove duplicated corners
142 size_t * num_tgt_per_corner =
143 xcalloc(num_unique_src_corners, sizeof(*num_tgt_per_corner));
144 num_unique_src_corners = 0;
145 for (size_t i = 0, prev_src_corner = SIZE_MAX; i < result_count; ++i) {
146 size_t curr_src_corner = src_corners[i];
147 if (prev_src_corner != curr_src_corner) {
148 prev_src_corner = curr_src_corner;
149 src_corners[num_unique_src_corners] = curr_src_corner;
150 ++num_unique_src_corners;
151 }
152 num_tgt_per_corner[num_unique_src_corners-1]++;
153 }
154
155 // get all source cells associated with each unique source corner;
156 // results cells for each corner are sorted by their global id
157 size_t * src_corner_cells = NULL;
158 size_t * num_cells_per_corner =
159 xmalloc(num_unique_src_corners * sizeof(*num_cells_per_corner));
161 interp_grid, src_corners, num_unique_src_corners, &src_corner_cells,
162 num_cells_per_corner);
163
164 // get maximum number of source cells per target
165 size_t max_num_cell_per_corner = 0;
166 size_t total_num_weights = 0;
167 for (size_t i = 0; i < num_unique_src_corners; ++i) {
168 if (max_num_cell_per_corner < num_cells_per_corner[i])
169 max_num_cell_per_corner = num_cells_per_corner[i];
170 total_num_weights += num_cells_per_corner[i] * num_tgt_per_corner[i];
171 }
172
173 yac_const_coordinate_pointer src_field_coordinates =
175 const_int_pointer src_field_mask =
177
178 yac_coordinate_pointer src_coord_buffer =
179 xmalloc(max_num_cell_per_corner * sizeof(*src_coord_buffer));
180 int * mask_buffer =
181 (src_field_mask != NULL)?
182 xmalloc(max_num_cell_per_corner * sizeof(*mask_buffer)):NULL;
183 double * w = xmalloc(total_num_weights * sizeof(*w));
184 size_t * src_points = xmalloc(total_num_weights * sizeof(*src_points));
185 size_t * num_weights_per_tgt =
186 xmalloc(result_count * sizeof(*num_weights_per_tgt));
187
188 // for each target point, extract relevant source data and compute the weights
189 // based on that
190 func_compute_weights compute_weights = method_ncc->compute_weights;
191 total_num_weights = 0;
192 result_count = 0;
193 for (size_t i = 0, l = 0, src_corner_cell_offset = 0;
194 i < num_unique_src_corners; ++i) {
195
196 size_t curr_num_cells = num_cells_per_corner[i];
197
198 const_size_t_pointer curr_src_corner_cells =
199 src_corner_cells + src_corner_cell_offset;
200 src_corner_cell_offset += curr_num_cells;
201
202 // get field coordinates and mask values for current source cells
203 for (size_t j = 0; j < curr_num_cells; ++j) {
204 size_t const curr_src_cell = curr_src_corner_cells[j];
205 for (size_t k = 0; k < 3; ++k)
206 src_coord_buffer[j][k] = src_field_coordinates[curr_src_cell][k];
207 if (src_field_mask != NULL)
208 mask_buffer[j] = src_field_mask[curr_src_cell];
209 }
210
211 for (size_t j = 0; j < num_tgt_per_corner[i]; ++j, ++l) {
212
213 // compute weights
214 if (compute_weights(
215 tgt_coords[reorder_idx[l]], curr_num_cells,
216 (yac_const_coordinate_pointer)src_coord_buffer,
217 mask_buffer, w + total_num_weights)) {
218
219 interp_flag[reorder_idx[l]] = result_count;
220 memcpy(
221 src_points + total_num_weights, curr_src_corner_cells,
222 curr_num_cells * sizeof(*src_points));
223 num_weights_per_tgt[result_count] = curr_num_cells;
224 total_num_weights += curr_num_cells;
225 result_count++;
226 } else {
227 interp_flag[reorder_idx[l]] = SIZE_MAX;
228 }
229 }
230 }
231
232 free(src_corner_cells);
233 free(num_tgt_per_corner);
234 free(num_cells_per_corner);
235
236 // move target points for which the interpolation failed to the end of
237 // the array while keeping the origin order for the remaining target points
238 yac_quicksort_index_size_t_size_t(interp_flag, count, tgt_points);
239
240 free(tgt_coords);
241 free(mask_buffer);
242 free(src_coord_buffer);
243 free(size_t_buffer);
244
245 struct remote_points tgts = {
246 .data =
248 interp_grid, tgt_points, result_count),
249 .count = result_count};
250 struct remote_point * srcs =
252 interp_grid, 0, src_points, total_num_weights);
253
254 // store weights
256 weights, &tgts, num_weights_per_tgt, srcs, w);
257
258 free(tgts.data);
259 free(srcs);
260 free(w);
261 free(num_weights_per_tgt);
262 free(src_points);
263
264 return result_count;
265}
266
268 double tgt_coords[3], size_t num_src,
269 yac_const_coordinate_pointer src_coords, int * src_mask, double * weights) {
270
271 UNUSED(tgt_coords);
272 UNUSED(src_coords);
273
274 size_t num_unmasked_points = 0;
275
276 // if we have a source mask, check whether any points is masked
277 if (src_mask != NULL) {
278 for (size_t i = 0; i < num_src; ++i)
279 if (src_mask[i]) num_unmasked_points++;
280 if (num_unmasked_points == 0) return 0;
281
282 double weight = 1.0 / (double)num_unmasked_points;
283 for (size_t i = 0; i < num_src; ++i)
284 weights[i] = (src_mask[i])?weight:0.0;
285
286 } else {
287 double weight = 1.0 / (double)num_src;
288 for (size_t i = 0; i < num_src; ++i) weights[i] = weight;
289 }
290 return 1;
291}
292
294 double tgt_coords[3], size_t num_src,
295 yac_const_coordinate_pointer src_coords, int * src_mask, double * weights) {
296
297 UNUSED(tgt_coords);
298 UNUSED(src_coords);
299
300 // if we have a source mask, check whether any points is masked
301 if (src_mask != NULL)
302 for (size_t i = 0; i < num_src; ++i)
303 if (!src_mask[i]) return 0;
304
305 double weight = 1.0 / (double)num_src;
306 for (size_t i = 0; i < num_src; ++i) weights[i] = weight;
307
308 return 1;
309}
310
312 double tgt_coords[3], size_t num_src,
313 yac_const_coordinate_pointer src_coords, int * src_mask, double * weights) {
314
315 // if there is a source mask, check if there are unmasked vertices
316 if (src_mask != NULL) {
317 int has_unmasked = 0;
318 for (size_t i = 0; i < num_src; ++i) has_unmasked |= src_mask[i];
319 if (!has_unmasked) return 0;
320 }
321
322 if (src_mask != NULL) {
323 for (size_t i = 0; i < num_src; ++i) {
324
325 if (src_mask[i]) {
326
327 double distance =
328 get_vector_angle(tgt_coords, (double*)(src_coords[i]));
329
330 // if the target and source point are nearly identical
331 if (distance < yac_angle_tol) {
332 for (size_t j = 0; j < num_src; ++j) weights[j] = 0.0;
333 weights[i] = 1.0;
334 return 1;
335 }
336
337 weights[i] = 1.0 / distance;
338 } else {
339 weights[i] = 0.0;
340 }
341 }
342 } else {
343 for (size_t i = 0; i < num_src; ++i) {
344
345 double distance =
346 get_vector_angle(tgt_coords, (double*)(src_coords[i]));
347
348 // if the target and source point are nearly identical
349 if (distance < yac_angle_tol) {
350 for (size_t j = 0; j < num_src; ++j) weights[j] = 0.0;
351 weights[i] = 1.0;
352 return 1;
353 }
354
355 weights[i] = 1.0 / distance;
356 }
357 }
358
359 // compute scaling factor for the weights
360 double inv_distance_sum = 0.0;
361 for (size_t i = 0; i < num_src; ++i)
362 inv_distance_sum += weights[i];
363 double scale = 1.0 / inv_distance_sum;
364
365 for (size_t i = 0; i < num_src; ++i) weights[i] *= scale;
366
367 return 1;
368}
369
371 double tgt_coords[3], size_t num_src,
372 yac_const_coordinate_pointer src_coords, int * src_mask, double * weights) {
373
374 for (size_t i = 0; i < num_src; ++i) {
375
376 double distance =
377 get_vector_angle(tgt_coords, (double*)(src_coords[i]));
378
379 // if the target and source point are nearly identical
380 if (distance < yac_angle_tol) {
381 if ((src_mask != NULL) && !src_mask[i]) return 0;
382 for (size_t j = 0; j < num_src; ++j) weights[j] = 0.0;
383 weights[i] = 1.0;
384 return 1;
385 }
386
387 weights[i] = 1.0 / distance;
388 }
389
390 // if there is a source mask, check if there are masked vertices
391 // (we do this here because there may be a matching source point
392 // that is not masked)
393 if (src_mask != NULL)
394 for (size_t i = 0; i < num_src; ++i) if(!src_mask[i]) return 0;
395
396 // compute scaling factor for the weights
397 double inv_distance_sum = 0.0;
398 for (size_t i = 0; i < num_src; ++i)
399 inv_distance_sum += weights[i];
400 double scale = 1.0 / inv_distance_sum;
401
402 for (size_t i = 0; i < num_src; ++i) weights[i] *= scale;
403
404 return 1;
405}
406
408 enum yac_interp_ncc_weight_type weight_type, int partial_coverage) {
409
413 "ERROR(select_compute_weight_routine): invalid weight type")
414
415 switch(weight_type) {
416 case(YAC_INTERP_NCC_AVG):
417 return
420 default:
421 return
423 };
424}
425
428 int partial_coverage) {
429
430 struct interp_method_ncc * method = xmalloc(1 * sizeof(*method));
431
433 method->compute_weights =
435
436 return (struct interp_method*)method;
437}
438
439static void delete_ncc(struct interp_method * method) {
440 free(method);
441}
#define YAC_ASSERT(exp, msg)
#define UNUSED(x)
Definition core.h:72
int const * const_int_pointer
size_t const *const const_size_t_pointer
#define yac_angle_tol
Definition geometry.h:26
static double get_vector_angle(double const a[3], double const b[3])
Definition geometry.h:340
void yac_interp_grid_do_points_search(struct yac_interp_grid *interp_grid, yac_coordinate_pointer search_coords, size_t count, size_t *src_cells)
size_t yac_interp_grid_get_num_src_fields(struct yac_interp_grid *interp_grid)
const_int_pointer yac_interp_grid_get_src_field_mask(struct yac_interp_grid *interp_grid, size_t src_field_idx)
struct remote_point * yac_interp_grid_get_tgt_remote_points(struct yac_interp_grid *interp_grid, size_t *tgt_points, size_t count)
enum yac_location yac_interp_grid_get_src_field_location(struct yac_interp_grid *interp_grid, size_t src_field_idx)
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)
void yac_interp_grid_get_src_corner_cells(struct yac_interp_grid *interp_grid, size_t *src_corners, size_t count, size_t **src_cells, size_t *num_cells_per_corner)
void yac_interp_grid_get_tgt_coordinates(struct yac_interp_grid *interp_grid, size_t *tgt_points, size_t count, yac_coordinate_pointer tgt_coordinates)
struct yac_const_basic_grid_data * yac_interp_grid_get_basic_grid_data_src(struct yac_interp_grid *interp_grid)
static void compute_weights(struct tgt_point_search_data *tgt_point_data, size_t num_tgt_points, struct edge_interp_data *edge_data, size_t num_edges, struct triangle_interp_data *triangle_data, size_t num_triangles, struct weight_vector_data **weights, size_t **num_weights_per_tgt, size_t *total_num_weights)
static int compute_weights_dist_no(double tgt_coords[3], size_t num_src, yac_const_coordinate_pointer src_coords, int *src_mask, double *weights)
static void delete_ncc(struct interp_method *method)
static int compute_weights_avg_yes(double tgt_coords[3], size_t num_src, yac_const_coordinate_pointer src_coords, int *src_mask, double *weights)
struct interp_method * yac_interp_method_ncc_new(enum yac_interp_ncc_weight_type weight_type, int partial_coverage)
static size_t do_search_ncc(struct interp_method *method, struct yac_interp_grid *interp_grid, size_t *tgt_points, size_t count, struct yac_interp_weights *weights, int *interpolation_compelte)
static size_t get_closest_src_corner(struct yac_const_basic_grid_data *grid_data, size_t cell, double const coord[3])
static func_compute_weights select_compute_weight_routine(enum yac_interp_ncc_weight_type weight_type, int partial_coverage)
static struct interp_method_vtable interp_method_ncc_vtable
int(* func_compute_weights)(double[3], size_t, yac_const_coordinate_pointer, int *, double *)
static int compute_weights_avg_no(double tgt_coords[3], size_t num_src, yac_const_coordinate_pointer src_coords, int *src_mask, double *weights)
static int compute_weights_dist_yes(double tgt_coords[3], size_t num_src, yac_const_coordinate_pointer src_coords, int *src_mask, double *weights)
yac_interp_ncc_weight_type
@ YAC_INTERP_NCC_DIST
distance weighted average of n source points
@ YAC_INTERP_NCC_AVG
average of n source points
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)
@ YAC_LOC_CORNER
Definition location.h:15
@ YAC_LOC_CELL
Definition location.h:14
#define xcalloc(nmemb, size)
Definition ppm_xfuncs.h:64
#define xmalloc(size)
Definition ppm_xfuncs.h:66
struct interp_method_vtable * vtable
func_compute_weights compute_weights
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, int *interpolation_complete)
information (global id and location) about a point that
structure containing the information (global id and location)
struct remote_point * data
int * cell_to_vertex
void yac_quicksort_index_size_t_size_t(size_t *a, size_t n, size_t *idx)
double const (* yac_const_coordinate_pointer)[3]
Definition yac_types.h:22
YAC_INT yac_int
Definition yac_types.h:15
double(* yac_coordinate_pointer)[3]
Definition yac_types.h:21