A test for internal interpolation routines.
#include <stdlib.h>
#include <stdio.h>
#include <math.h>
#include <unistd.h>
#include <string.h>
#include <mpi.h>
#include "tests.h"
#include "test_common.h"
#include "weight_file_common.h"
#include "dist_grid_utils.h"
#define YAC_RAD (0.01745329251994329576923690768489)
static void
generate_ref_weights();
static void submain_src(
static void submain_tgt(
unsigned const num_links_file = 48;
int src_address_file[48] = {
9,10,16,17, 10,11,17,18, 15,16,22,23, 16,17,23,24, 17,18,24,25, 18,19,25,26,
22,23,29,30, 23,24,30,31, 24,25,31,32, 25,26,32,33, 30,31,37,38, 31,32,38,39};
int tgt_address_file[48] = {
8,8,8,8, 9,9,9,9, 13,13,13,13, 14,14,14,14, 15,15,15,15, 16,16,16,16,
19,19,19,19, 20,20,20,20, 21,21,21,21, 22,22,22,22, 26,26,26,26, 27,27,27,27};
double weights_file[48] = {
0.1,0.2,0.3,0.4, 0.1,0.2,0.3,0.4, 0.1,0.2,0.3,0.4, 0.1,0.2,0.3,0.4,
0.1,0.2,0.3,0.4, 0.1,0.2,0.3,0.4, 0.1,0.2,0.3,0.4, 0.1,0.2,0.3,0.4,
0.1,0.2,0.3,0.4, 0.1,0.2,0.3,0.4, 0.1,0.2,0.3,0.4, 0.1,0.2,0.3,0.4};
char const weight_file_in[] =
"test_interpolation_parallel3_weight_file_in.nc";
char const weight_file_out[] =
"test_interpolation_parallel3_weight_file_out.nc";
char const src_grid_name[] = "src_grid";
char const tgt_grid_name[] = "tgt_grid";
unsigned ref_num_links = 4 * 36;
int ref_src_address[4 * 36];
int ref_tgt_address[4 * 36];
double ref_weights[4 * 36];
int const * ref_tgt_address_fixed = NULL;
double const * ref_fixed_values = NULL;
int const * ref_num_tgt_per_fixed_value = NULL;
unsigned ref_num_fixed_values = 0;
size_t num_cells[2] = {6,6};
double coordinates_x[] = {0,1,2,3,4,5,6};
double coordinates_y[] = {0,1,2,3,4,5,6};
double cell_coordinates_x[] = {0.5,1.5,2.5,3.5,4.5,5.5};
double cell_coordinates_y[] = {0.5,1.5,2.5,3.5,4.5,5.5};
double cell_coords[36][3];
int with_halo = 1;
int main(
int argc,
char *argv[]) {
if (argc != 2) {
PUT_ERR("wrong number of arguments\n");
return TEST_EXIT_CODE;
}
PUT_ERR("invalid argument (has to be either \"src\" or \"tgt\")\n");
return TEST_EXIT_CODE;
}
MPI_Init(NULL, NULL);
xt_initialize(MPI_COMM_WORLD);
int comm_rank, comm_size;
MPI_Comm_rank(MPI_COMM_WORLD, &comm_rank);
MPI_Comm_size(MPI_COMM_WORLD, &comm_size);
MPI_Barrier(MPI_COMM_WORLD);
if (comm_size != 8) {
PUT_ERR("ERROR: wrong number of processes");
xt_finalize();
MPI_Finalize();
return TEST_EXIT_CODE;
}
for (
size_t i = 0; i <= num_cells[0]; ++i) coordinates_x[i] *=
YAC_RAD;
for (
size_t i = 0; i <= num_cells[1]; ++i) coordinates_y[i] *=
YAC_RAD;
for (size_t i = 0; i < 6; ++i)
for (size_t j = 0; j < 6; ++j)
cell_coordinates_x[j], cell_coordinates_y[i], cell_coords[6 * i + j]);
generate_ref_weights();
int comp_flag = comm_rank < 4;
MPI_Comm split_comm;
MPI_Comm_split(
MPI_COMM_WORLD, comp_flag, 0, &split_comm);
if (comp_flag) submain_src(split_comm, reorder_type);
else submain_tgt(split_comm, reorder_type);
MPI_Comm_free(&split_comm);
xt_finalize();
MPI_Finalize();
return TEST_EXIT_CODE;
}
static void submain_tgt(
int my_rank;
MPI_Comm_rank(comp_comm, &my_rank);
size_t local_start[4][2] = {{0,0},{3,0},{0,3},{3,3}};
size_t local_count[4][2] = {{3,3},{3,3},{3,3},{3,3}};
yac_generate_basic_grid_data_reg2d(
local_start[my_rank], local_count[my_rank], with_halo);
memcpy(cell_field_coords[i], cell_coords[grid_data.
cell_ids[i]],
3 * sizeof(double));
{
size_t num_src_fields = sizeof(src_fields) / sizeof(src_fields[0]);
NULL};
weights, reorder_type, 1,
double * target_data_field =
xmalloc(num_cells *
sizeof(*target_data_field));
double * target_data[1] = {target_data_field};
for (size_t i = 0; i < num_cells; ++i) target_data_field[i] = -1.0;
weights, weight_file_out, src_grid_name, tgt_grid_name, 0, 0);
double ref_global_target_data[6*6];
for (size_t i = 0; i < 6*6; ++i) {
ref_global_target_data[i] = 0;
for (size_t j = 0; j < 4; ++j)
ref_global_target_data[i] +=
ref_weights[4 * i + j] * (double)ref_src_address[4 * i + j];
}
for (size_t i = 0; i < num_cells; ++i) {
if (fabs(
target_data[0][i] -
ref_global_target_data[grid_data.
cell_ids[i]]) > 1e-10)
PUT_ERR("error in interpolated data on target side\n");
} else {
if (target_data[0][i] != -1.0)
PUT_ERR("error in interpolated data on target side\n");
}
}
free(target_data_field);
}
}
static void submain_src(
int my_rank;
MPI_Comm_rank(comp_comm, &my_rank);
if (my_rank == 0) {
int * tgt_id_fixed = NULL;
unsigned num_fixed_tgt = 0;
double * fixed_values = NULL;
int * num_tgt_per_fixed_value = NULL;
unsigned num_fixed_values = 0;
write_weight_file(weight_file_in, src_address_file, tgt_address_file,
weights_file, num_links_file, src_locations,
1, (int*)&num_links_file, tgt_id_fixed, num_fixed_tgt,
fixed_values, num_tgt_per_fixed_value, num_fixed_values,
tgt_location, src_grid_name, tgt_grid_name);
}
size_t local_start[4][2] = {{0,0},{1,0},{3,0},{5,0}};
size_t local_count[4][2] = {{1,6},{2,6},{2,6},{1,6}};
yac_generate_basic_grid_data_reg2d(
local_start[my_rank], local_count[my_rank], with_halo);
{
size_t num_src_fields = sizeof(src_fields) / sizeof(src_fields[0]);
NULL};
weights, reorder_type, 1,
double * source_data_field =
xmalloc(num_vertices *
sizeof(*source_data_field));
double * source_data_pointset[1] = {source_data_field};
double ** source_data[1] = {source_data_pointset};
for (size_t i = 0; i < num_vertices; ++i)
source_data_field[i] =
weights, weight_file_out, src_grid_name, tgt_grid_name, 0, 0);
if (my_rank == 0) {
check_weight_file(weight_file_out, ref_src_address, ref_tgt_address,
ref_weights, ref_num_links, ref_src_locations,
1, (int*)&ref_num_links, ref_tgt_address_fixed,
ref_fixed_values, ref_num_tgt_per_fixed_value,
ref_num_fixed_values, ref_tgt_location,
src_grid_name, tgt_grid_name);
}
free(source_data_field);
}
if (my_rank == 0) {
unlink(weight_file_in);
unlink(weight_file_out);
}
}
static void
generate_ref_weights() {
for (unsigned i = 0; i < 6; ++i) {
for (unsigned j = 0; j < 6; ++j) {
ref_src_address[4 * (i * 6 + j) + 0] = 0 + i * 7 + j;
ref_src_address[4 * (i * 6 + j) + 1] = 1 + i * 7 + j;
ref_src_address[4 * (i * 6 + j) + 2] = 7 + i * 7 + j;
ref_src_address[4 * (i * 6 + j) + 3] = 8 + i * 7 + j;
}
}
for (unsigned i = 0; i < ref_num_links; ++i) ref_tgt_address[i] = i / 4;
for (unsigned i = 0; i < ref_num_links; ++i) ref_weights[i] = 0.25;
for (unsigned i = 0; i < num_links_file; ++i)
ref_weights[tgt_address_file[i] * 4 + (i & 3)] = weights_file[i];
}
struct yac_basic_grid * yac_basic_grid_new(char const *name, struct yac_basic_grid_data grid_data)
size_t yac_basic_grid_add_coordinates_nocpy(struct yac_basic_grid *grid, enum yac_location location, yac_coordinate_pointer coordinates)
struct yac_basic_grid * yac_basic_grid_empty_new(char const *name)
void yac_basic_grid_delete(struct yac_basic_grid *grid)
int main(int argc, char **argv)
void yac_dist_grid_pair_delete(struct yac_dist_grid_pair *grid_pair)
struct yac_dist_grid_pair * yac_dist_grid_pair_new(struct yac_basic_grid *grid_a, struct yac_basic_grid *grid_b, MPI_Comm comm)
static void LLtoXYZ_deg(double lon, double lat, double p_out[])
void yac_interp_grid_delete(struct yac_interp_grid *interp_grid)
struct yac_interp_grid * yac_interp_grid_new(struct yac_dist_grid_pair *grid_pair, char const *src_grid_name, char const *tgt_grid_name, size_t num_src_fields, struct yac_interp_field const *src_fields, struct yac_interp_field const tgt_field)
void yac_interp_method_delete(struct interp_method **method)
struct yac_interp_weights * yac_interp_method_do_search(struct interp_method **method, struct yac_interp_grid *interp_grid)
struct interp_method * yac_interp_method_avg_new(enum yac_interp_avg_weight_type weight_type, int partial_coverage)
@ YAC_INTERP_AVG_ARITHMETIC
struct interp_method * yac_interp_method_file_new(char const *weight_file_name)
void yac_interp_weights_delete(struct yac_interp_weights *weights)
void yac_interp_weights_write_to_file(struct yac_interp_weights *weights, char const *filename, char const *src_grid_name, char const *tgt_grid_name, size_t src_grid_size, size_t tgt_grid_size)
struct yac_interpolation * yac_interp_weights_get_interpolation(struct yac_interp_weights *weights, enum yac_interp_weights_reorder_type reorder, size_t collection_size, double frac_mask_fallback_value, double scaling_factor, double scaling_summand, char const *yaxt_exchanger_name)
yac_interp_weights_reorder_type
@ YAC_MAPPING_ON_TGT
weights will be applied at target processes
@ YAC_MAPPING_ON_SRC
weights will be appied at source processes
void yac_interpolation_delete(struct yac_interpolation *interp)
void yac_interpolation_execute_get(struct yac_interpolation *interp, double **tgt_field)
void yac_interpolation_execute_put(struct yac_interpolation *interp, double ***src_fields)
double const YAC_FRAC_MASK_NO_VALUE
enum yac_location location
struct yac_interp_field tgt_field
struct yac_dist_grid_pair * grid_pair
struct yac_interp_field src_fields[]
double(* yac_coordinate_pointer)[3]