YetAnotherCoupler 3.5.2
Loading...
Searching...
No Matches
test_interp_method_creep_parallel.c

A test for the parallel creep interpolation method.

// Copyright (c) 2024 The YAC Authors
//
// SPDX-License-Identifier: BSD-3-Clause
#include <stdlib.h>
#include "tests.h"
#include "test_common.h"
#include "geometry.h"
#include "interp_method.h"
#include "dist_grid_utils.h"
#include "yac_mpi.h"
#include <mpi.h>
#include <yaxt.h>
#include <netcdf.h>
static char const * grid_names[2] = {"src_grid", "tgt_grid"};
int main(void) {
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 != 3) {
PUT_ERR("ERROR: wrong number of processes");
xt_finalize();
MPI_Finalize();
return TEST_EXIT_CODE;
}
MPI_Comm split_comm;
MPI_Comm_split(MPI_COMM_WORLD, comm_rank < 2, 0, &split_comm);
int split_comm_rank, split_comm_size;
MPI_Comm_rank(split_comm, &split_comm_rank);
MPI_Comm_size(split_comm, &split_comm_size);
{// stack that starts with linear point interpolation, followed by creep
// fill and fixed as backup with two source processes and a single target
// process
// the global source grid is a 3x2 grid:
// 08--14--09--15--10--16--11
// | | | |
// 08 03 10 04 12 05 13
// | | | |
// 04--07--05--09--06--11--07
// | | | |
// 01 00 03 01 05 02 06
// | | | |
// 00--00--01--02--02--04--03
//
//---------------
// setup
//---------------
int is_tgt = split_comm_size == 1;
double coordinates_x[] = {0.0,1.0,2.0,3.0};
double coordinates_y[] = {0.0,1.0,2.0,3.0,4.0,5.0,6.0,7.0,8.0};
double tgt_cell_coordinates_x[] = {0.5,1.5,2.5};
double tgt_cell_coordinates_y[] = {0.5,1.5,2.5,3.5,4.5,5.5,6.5,7.5};
double tgt_cell_coordinates[24][3];
int tgt_cell_mask[24] = {1, 1, 1,
1, 1, 1,
1, 1, 1,
1, 0, 1,
1, 0, 0,
1, 1, 1,
0, 0, 0,
1, 1, 1};
size_t const num_cells[2][2] = {{3,2}, {3,8}};
size_t local_start[2][2][2] = {{{0,0},{1,0}}, {{0,0}}};
size_t local_count[2][2][2] = {{{2,2},{2,2}}, {{3,8}}};
int with_halo = 0;
for (size_t i = 0; i <= num_cells[is_tgt][0]; ++i)
coordinates_x[i] *= YAC_RAD;
for (size_t i = 0; i <= num_cells[is_tgt][1]; ++i)
coordinates_y[i] *= YAC_RAD;
struct yac_basic_grid_data grid_data =
yac_generate_basic_grid_data_reg2d(
coordinates_x, coordinates_y, num_cells[is_tgt],
local_start[is_tgt][split_comm_rank],
local_count[is_tgt][split_comm_rank], with_halo);
struct yac_basic_grid * grids[2] =
{yac_basic_grid_new(grid_names[is_tgt], grid_data),
yac_basic_grid_empty_new(grid_names[is_tgt^1])};
if (is_tgt) {
for (size_t i = 0, k = 0; i < num_cells[is_tgt][1]; ++i)
for (size_t j = 0; j < num_cells[is_tgt][0]; ++j, ++k)
tgt_cell_coordinates_x[j], tgt_cell_coordinates_y[i],
tgt_cell_coordinates[k]);
grids[0], YAC_LOC_CELL, tgt_cell_coordinates, grid_data.num_cells);
grids[0], YAC_LOC_CELL, tgt_cell_mask, grid_data.num_cells, NULL);
}
struct yac_dist_grid_pair * grid_pair =
yac_dist_grid_pair_new(grids[is_tgt], grids[is_tgt^1], MPI_COMM_WORLD);
struct yac_interp_field src_fields[] =
{{.location = YAC_LOC_CORNER, .coordinates_idx = SIZE_MAX, .masks_idx = SIZE_MAX}};
size_t num_src_fields = sizeof(src_fields) / sizeof(src_fields[0]);
struct yac_interp_field tgt_field =
{.location = YAC_LOC_CELL, .coordinates_idx = 0, .masks_idx = 0};
struct yac_interp_grid * interp_grid =
yac_interp_grid_new(grid_pair, grid_names[0], grid_names[1],
for (int creep_distance = -1; creep_distance <= 7; ++creep_distance) {
struct interp_method * method_stack[] =
yac_interp_method_creep_new(creep_distance),
struct yac_interp_weights * weights =
yac_interp_method_do_search(method_stack, interp_grid);
enum yac_interp_weights_reorder_type reorder_type[2] =
for (size_t i = 0; i < 2; ++i) {
struct yac_interpolation * interpolation =
weights, reorder_type[i], 1,
YAC_FRAC_MASK_NO_VALUE, 1.0, 0.0, NULL);
// check generated interpolation
{
double * src_field = NULL;
double ** src_fields = &src_field;
double * tgt_field = NULL;
double * ref_tgt_field = NULL;
double ref_global_tgt_field[9][24] =
{{ 2.5, 3.5, 4.5, // creep_distance = -1
6.5, 7.5, 8.5,
6.5, 7.5, 8.5,
6.5, -2.0, 8.5,
6.5, -2.0, -2.0,
6.5, 6.5, 6.5,
-2.0, -2.0, -2.0,
-1.0, -1.0, -1.0},
{ 2.5, 3.5, 4.5, // creep_distance = 0
6.5, 7.5, 8.5,
-1.0, -1.0, -1.0,
-1.0, -2.0, -1.0,
-1.0, -2.0, -2.0,
-1.0, -1.0, -1.0,
-2.0, -2.0, -2.0,
-1.0, -1.0, -1.0},
{ 2.5, 3.5, 4.5, // creep_distance = 1
6.5, 7.5, 8.5,
6.5, 7.5, 8.5,
-1.0, -2.0, -1.0,
-1.0, -2.0, -2.0,
-1.0, -1.0, -1.0,
-2.0, -2.0, -2.0,
-1.0, -1.0, -1.0},
{ 2.5, 3.5, 4.5, // creep_distance = 2
6.5, 7.5, 8.5,
6.5, 7.5, 8.5,
6.5, -2.0, 8.5,
-1.0, -2.0, -2.0,
-1.0, -1.0, -1.0,
-2.0, -2.0, -2.0,
-1.0, -1.0, -1.0},
{ 2.5, 3.5, 4.5, // creep_distance = 3
6.5, 7.5, 8.5,
6.5, 7.5, 8.5,
6.5, -2.0, 8.5,
6.5, -2.0, -2.0,
-1.0, -1.0, -1.0,
-2.0, -2.0, -2.0,
-1.0, -1.0, -1.0},
{ 2.5, 3.5, 4.5, // creep_distance = 4
6.5, 7.5, 8.5,
6.5, 7.5, 8.5,
6.5, -2.0, 8.5,
6.5, -2.0, -2.0,
6.5, -1.0, -1.0,
-2.0, -2.0, -2.0,
-1.0, -1.0, -1.0},
{ 2.5, 3.5, 4.5, // creep_distance = 5
6.5, 7.5, 8.5,
6.5, 7.5, 8.5,
6.5, -2.0, 8.5,
6.5, -2.0, -2.0,
6.5, 6.5, -1.0,
-2.0, -2.0, -2.0,
-1.0, -1.0, -1.0},
{ 2.5, 3.5, 4.5, // creep_distance = 6
6.5, 7.5, 8.5,
6.5, 7.5, 8.5,
6.5, -2.0, 8.5,
6.5, -2.0, -2.0,
6.5, 6.5, 6.5,
-2.0, -2.0, -2.0,
-1.0, -1.0, -1.0},
{ 2.5, 3.5, 4.5, // creep_distance = 7
6.5, 7.5, 8.5,
6.5, 7.5, 8.5,
6.5, -2.0, 8.5,
6.5, -2.0, -2.0,
6.5, 6.5, 6.5,
-2.0, -2.0, -2.0,
-1.0, -1.0, -1.0}};
if (is_tgt) {
tgt_field = xmalloc(grid_data.num_cells * sizeof(*tgt_field));
ref_tgt_field = xmalloc(grid_data.num_cells * sizeof(*ref_tgt_field));
for (size_t i = 0; i < grid_data.num_cells; ++i) {
tgt_field[i] = -2.0;
ref_tgt_field[i] =
ref_global_tgt_field[creep_distance+1][grid_data.cell_ids[i]];
}
} else {
src_field = xmalloc(grid_data.num_vertices * sizeof(*src_field));
for (size_t i = 0; i < grid_data.num_vertices; ++i)
src_field[i] = (double)(grid_data.vertex_ids[i]);
}
yac_interpolation_execute(interpolation, &src_fields, &tgt_field);
if (is_tgt)
for (size_t i = 0; i < grid_data.num_cells; ++i)
if (fabs(ref_tgt_field[i] - tgt_field[i]) > 1e-9)
PUT_ERR("wrong interpolation result");
free(src_field);
free(tgt_field);
free(ref_tgt_field);
}
yac_interpolation_delete(interpolation);
}
yac_interp_method_delete(method_stack);
}
yac_interp_grid_delete(interp_grid);
}
{// stack that starts with linear point interpolation, followed by creep
// fill and fixed as backup with two source processes and a single target
// process
// the global source grid is a 7x7 grid
// the global target grid is a 6x6 grid
//
//---------------
// setup
//---------------
int is_tgt = split_comm_size == 1;
double coordinates_x[2][8] = {{0.0,1.0,2.0,3.0,4.0,5.0,6.0,7.0},
{0.5,1.5,2.5,3.5,4.5,5.5,6.5}};
double coordinates_y[2][8] = {{0.0,1.0,2.0,3.0,4.0,5.0,6.0,7.0},
{0.5,1.5,2.5,3.5,4.5,5.5,6.5}};
int vertex_mask[2][8*8] =
{{1,1,1,1,1,1,1,1,
1,1,1,1,1,1,1,1,
1,1,0,0,0,0,1,1,
1,1,0,0,0,0,1,1,
1,1,0,0,0,0,1,1,
1,1,0,0,0,0,1,1,
1,1,1,1,1,1,1,1,
1,1,1,1,1,1,1,1},
{1,1,1,1,1,1,1,
1,1,1,1,1,1,1,
1,1,1,0,1,1,1,
1,1,0,1,0,1,1,
1,1,1,0,1,1,1,
1,1,1,1,1,1,1,
1,1,1,1,1,1,1}};
size_t const num_cells[2][2] = {{7,7}, {6,6}};
size_t local_start[2][2][2] = {{{0,0},{3,0}}, {{0,0}}};
size_t local_count[2][2][2] = {{{3,7},{4,7}}, {{6,6}}};
int with_halo = 0;
for (size_t i = 0; i <= num_cells[is_tgt][0]; ++i)
coordinates_x[is_tgt][i] *= YAC_RAD;
for (size_t i = 0; i <= num_cells[is_tgt][1]; ++i)
coordinates_y[is_tgt][i] *= YAC_RAD;
struct yac_basic_grid_data grid_data =
yac_generate_basic_grid_data_reg2d(
coordinates_x[is_tgt], coordinates_y[is_tgt], num_cells[is_tgt],
local_start[is_tgt][split_comm_rank],
local_count[is_tgt][split_comm_rank], with_halo);
struct yac_basic_grid * grids[2] =
{yac_basic_grid_new(grid_names[is_tgt], grid_data),
yac_basic_grid_empty_new(grid_names[is_tgt^1])};
int temp_vertex_mask[8*8];
for (size_t i = 0; i < grid_data.num_vertices; ++i)
temp_vertex_mask[i] = vertex_mask[is_tgt][grid_data.vertex_ids[i]];
grids[0], YAC_LOC_CORNER, temp_vertex_mask, grid_data.num_vertices, NULL);
struct yac_dist_grid_pair * grid_pair =
yac_dist_grid_pair_new(grids[is_tgt], grids[is_tgt^1], MPI_COMM_WORLD);
struct yac_interp_field src_fields[] =
{{.location = YAC_LOC_CORNER, .coordinates_idx = SIZE_MAX, .masks_idx = 0}};
size_t num_src_fields = sizeof(src_fields) / sizeof(src_fields[0]);
struct yac_interp_field tgt_field =
{.location = YAC_LOC_CORNER, .coordinates_idx = SIZE_MAX, .masks_idx = 0};
struct yac_interp_grid * interp_grid =
yac_interp_grid_new(grid_pair, grid_names[0], grid_names[1],
for (int creep_distance = -1; creep_distance <= 3; ++creep_distance) {
struct interp_method * method_stack[] =
yac_interp_method_creep_new(creep_distance),
struct yac_interp_weights * weights =
yac_interp_method_do_search(method_stack, interp_grid);
enum yac_interp_weights_reorder_type reorder_type[2] =
for (size_t i = 0; i < 2; ++i) {
struct yac_interpolation * interpolation =
weights, reorder_type[i], 1,
YAC_FRAC_MASK_NO_VALUE, 1.0, 0.0, NULL);
// check generated interpolation
{
double * src_field = NULL;
double ** src_fields = &src_field;
double * tgt_field = NULL;
double * ref_tgt_field = NULL;
double ref_global_tgt_field[5][7*7] =
{{ 4.5, 5.5, 6.5, 7.5, 8.5, 9.5, 10.5, // creep_distance = -1
12.5, 9.0, 6.5, 7.5, 8.5, 14.0, 18.5,
20.5, 20.5, 13.5, -2.0, 17.5, 26.5, 26.5,
28.5, 28.5, -2.0, -1.0, -2.0, 34.5, 34.5,
36.5, 36.5, 45.5, -2.0, 49.5, 42.5, 42.5,
44.5, 49.0, 54.5, 55.5, 56.5, 54.0, 50.5,
52.5, 53.5, 54.5, 55.5, 56.5, 57.5, 58.5},
{ 4.5, 5.5, 6.5, 7.5, 8.5, 9.5, 10.5, // creep_distance = 0
12.5, -1.0, -1.0, -1.0, -1.0, -1.0, 18.5,
20.5, -1.0, -1.0, -2.0, -1.0, -1.0, 26.5,
28.5, -1.0, -2.0, -1.0, -2.0, -1.0, 34.5,
36.5, -1.0, -1.0, -2.0, -1.0, -1.0, 42.5,
44.5, -1.0, -1.0, -1.0, -1.0, -1.0, 50.5,
52.5, 53.5, 54.5, 55.5, 56.5, 57.5, 58.5},
{ 4.5, 5.5, 6.5, 7.5, 8.5, 9.5, 10.5, // creep_distance = 1
12.5, 9.0, 6.5, 7.5, 8.5, 14.0, 18.5,
20.5, 20.5, -1.0, -2.0, -1.0, 26.5, 26.5,
28.5, 28.5, -2.0, -1.0, -2.0, 34.5, 34.5,
36.5, 36.5, -1.0, -2.0, -1.0, 42.5, 42.5,
44.5, 49.0, 54.5, 55.5, 56.5, 54.0, 50.5,
52.5, 53.5, 54.5, 55.5, 56.5, 57.5, 58.5},
{ 4.5, 5.5, 6.5, 7.5, 8.5, 9.5, 10.5, // creep_distance = 2
12.5, 9.0, 6.5, 7.5, 8.5, 14.0, 18.5,
20.5, 20.5, 13.5, -2.0, 17.5, 26.5, 26.5,
28.5, 28.5, -2.0, -1.0, -2.0, 34.5, 34.5,
36.5, 36.5, 45.5, -2.0, 49.5, 42.5, 42.5,
44.5, 49.0, 54.5, 55.5, 56.5, 54.0, 50.5,
52.5, 53.5, 54.5, 55.5, 56.5, 57.5, 58.5},
{ 4.5, 5.5, 6.5, 7.5, 8.5, 9.5, 10.5, // creep_distance = 3
12.5, 9.0, 6.5, 7.5, 8.5, 14.0, 18.5,
20.5, 20.5, 13.5, -2.0, 17.5, 26.5, 26.5,
28.5, 28.5, -2.0, -1.0, -2.0, 34.5, 34.5,
36.5, 36.5, 45.5, -2.0, 49.5, 42.5, 42.5,
44.5, 49.0, 54.5, 55.5, 56.5, 54.0, 50.5,
52.5, 53.5, 54.5, 55.5, 56.5, 57.5, 58.5}};
if (is_tgt) {
tgt_field = xmalloc(grid_data.num_vertices * sizeof(*tgt_field));
ref_tgt_field = xmalloc(grid_data.num_vertices * sizeof(*ref_tgt_field));
for (size_t i = 0; i < grid_data.num_vertices; ++i) {
tgt_field[i] = -2.0;
ref_tgt_field[i] =
ref_global_tgt_field[creep_distance+1][grid_data.vertex_ids[i]];
}
} else {
src_field = xmalloc(grid_data.num_vertices * sizeof(*src_field));
for (size_t i = 0; i < grid_data.num_vertices; ++i)
src_field[i] = (double)(grid_data.vertex_ids[i]);
}
yac_interpolation_execute(interpolation, &src_fields, &tgt_field);
if (is_tgt)
for (size_t i = 0; i < grid_data.num_vertices; ++i)
if (fabs(ref_tgt_field[i] - tgt_field[i]) > 1e-9)
PUT_ERR("wrong interpolation result");
free(src_field);
free(tgt_field);
free(ref_tgt_field);
}
yac_interpolation_delete(interpolation);
}
yac_interp_method_delete(method_stack);
}
yac_interp_grid_delete(interp_grid);
}
{// stack that starts with linear point interpolation, followed by creep
// fill and fixed as backup with two source processes and a single target
// process
// the global source grid is a 5x1 grid
// the global target grid is a 5x3 grid
//
//---------------
// setup
//---------------
int is_tgt = split_comm_size == 1;
double coordinates_x[6] = {0.0,1.0,2.0,3.0,4.0,5.0};
double coordinates_y[4] = {0.0,1.0,2.0,3.0};
int cell_mask[2][5*3] =
{{1,0,1,0,1,},
{1,1,1,1,1,
1,1,0,1,1,
1,1,1,1,1}};
size_t const num_cells[2][2] = {{5,1}, {5,3}};
size_t local_start[2][2][2] = {{{0,0},{3,0}}, {{0,0}}};
size_t local_count[2][2][2] = {{{3,1},{2,1}}, {{5,3}}};
int with_halo = 0;
for (size_t i = 0; i <= num_cells[is_tgt][0]; ++i)
coordinates_x[i] *= YAC_RAD;
for (size_t i = 0; i <= num_cells[is_tgt][1]; ++i)
coordinates_y[i] *= YAC_RAD;
struct yac_basic_grid_data grid_data =
yac_generate_basic_grid_data_reg2d(
coordinates_x, coordinates_y, num_cells[is_tgt],
local_start[is_tgt][split_comm_rank],
local_count[is_tgt][split_comm_rank], with_halo);
struct yac_basic_grid * grids[2] =
{yac_basic_grid_new(grid_names[is_tgt], grid_data),
yac_basic_grid_empty_new(grid_names[is_tgt^1])};
int temp_cell_mask[5*3];
for (size_t i = 0; i < grid_data.num_cells; ++i)
temp_cell_mask[i] = cell_mask[is_tgt][grid_data.cell_ids[i]];
grids[0], YAC_LOC_CELL, temp_cell_mask, grid_data.num_cells, NULL);
struct yac_dist_grid_pair * grid_pair =
yac_dist_grid_pair_new(grids[is_tgt], grids[is_tgt^1], MPI_COMM_WORLD);
struct yac_interp_field src_fields[] =
{{.location = YAC_LOC_CELL, .coordinates_idx = SIZE_MAX, .masks_idx = 0}};
size_t num_src_fields = sizeof(src_fields) / sizeof(src_fields[0]);
struct yac_interp_field tgt_field =
{.location = YAC_LOC_CELL, .coordinates_idx = SIZE_MAX, .masks_idx = 0};
struct yac_interp_grid * interp_grid =
yac_interp_grid_new(grid_pair, grid_names[0], grid_names[1],
for (int creep_distance = -1; creep_distance <= 5; ++creep_distance) {
struct interp_method * method_stack[] =
yac_interp_method_creep_new(creep_distance),
struct yac_interp_weights * weights =
yac_interp_method_do_search(method_stack, interp_grid);
enum yac_interp_weights_reorder_type reorder_type[2] =
for (size_t i = 0; i < 2; ++i) {
struct yac_interpolation * interpolation =
weights, reorder_type[i], 1,
YAC_FRAC_MASK_NO_VALUE, 1.0, 0.0, NULL);
// check generated interpolation
{
double * src_field = NULL;
double ** src_fields = &src_field;
double * tgt_field = NULL;
double * ref_tgt_field = NULL;
double ref_global_tgt_field[7][7*7] =
{{ 0.0, 1.0 , 2.0 , 3.0 , 4.0, // creep_distance = -1
0.0, 0.5 , -2.0 , 3.5 , 4.0,
0.0, 0.25, 2.0 , 3.75, 4.0},
{ 0.0, -1.0 , 2.0 , -1.0 , 4.0, // creep_distance = 0
-1.0, -1.0 , -2.0 , -1.0 , -1.0,
-1.0, -1.0 , -1.0 , -1.0 , -1.0},
{ 0.0, 1.0 , 2.0 , 3.0 , 4.0, // creep_distance = 1
0.0, -1.0 , -2.0 , -1.0 , 4.0,
-1.0, -1.0 , -1.0 , -1.0 , -1.0},
{ 0.0, 1.0 , 2.0 , 3.0 , 4.0, // creep_distance = 2
0.0, 0.5 , -2.0 , 3.5 , 4.0,
0.0, -1.0 , -1.0 , -1.0 , 4.0},
{ 0.0, 1.0 , 2.0 , 3.0 , 4.0, // creep_distance = 3
0.0, 0.5 , -2.0 , 3.5 , 4.0,
0.0, 0.25, -1.0 , 3.75, 4.0},
{ 0.0, 1.0 , 2.0 , 3.0 , 4.0, // creep_distance = 4
0.0, 0.5 , -2.0 , 3.5 , 4.0,
0.0, 0.25, 2.0 , 3.75, 4.0},
{ 0.0, 1.0 , 2.0 , 3.0 , 4.0, // creep_distance = 5
0.0, 0.5 , -2.0 , 3.5 , 4.0,
0.0, 0.25, 2.0 , 3.75, 4.0}};
if (is_tgt) {
tgt_field = xmalloc(grid_data.num_cells* sizeof(*tgt_field));
ref_tgt_field = xmalloc(grid_data.num_cells* sizeof(*ref_tgt_field));
for (size_t i = 0; i < grid_data.num_cells; ++i) {
tgt_field[i] = -2.0;
ref_tgt_field[i] =
ref_global_tgt_field[creep_distance+1][grid_data.cell_ids[i]];
}
} else {
src_field = xmalloc(grid_data.num_cells* sizeof(*src_field));
for (size_t i = 0; i < grid_data.num_cells; ++i)
src_field[i] = (double)(grid_data.cell_ids[i]);
}
yac_interpolation_execute(interpolation, &src_fields, &tgt_field);
if (is_tgt)
for (size_t i = 0; i < grid_data.num_cells; ++i)
if (fabs(ref_tgt_field[i] - tgt_field[i]) > 1e-9)
PUT_ERR("wrong interpolation result");
free(src_field);
free(tgt_field);
free(ref_tgt_field);
}
yac_interpolation_delete(interpolation);
}
yac_interp_method_delete(method_stack);
}
yac_interp_grid_delete(interp_grid);
}
MPI_Comm_free(&split_comm);
xt_finalize();
MPI_Finalize();
return TEST_EXIT_CODE;
}
struct yac_basic_grid * yac_basic_grid_new(char const *name, struct yac_basic_grid_data grid_data)
Definition basic_grid.c:50
size_t yac_basic_grid_add_mask(struct yac_basic_grid *grid, enum yac_location location, int const *mask, size_t count, char const *mask_name)
Definition basic_grid.c:284
size_t yac_basic_grid_add_coordinates(struct yac_basic_grid *grid, enum yac_location location, yac_coordinate_pointer coordinates, size_t count)
Definition basic_grid.c:232
struct yac_basic_grid * yac_basic_grid_empty_new(char const *name)
Definition basic_grid.c:63
void yac_basic_grid_delete(struct yac_basic_grid *grid)
Definition basic_grid.c:70
int main(int argc, char **argv)
void yac_dist_grid_pair_delete(struct yac_dist_grid_pair *grid_pair)
Definition dist_grid.c:2608
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)
Definition dist_grid.c:2359
#define YAC_RAD
Definition geometry.h:30
static void LLtoXYZ_deg(double lon, double lat, double p_out[])
Definition geometry.h:304
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)
Definition interp_grid.c:31
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_conserv_new(int order, int enforced_conserv, int partial_coverage, enum yac_interp_method_conserv_normalisation normalisation)
@ YAC_INTERP_CONSERV_DESTAREA
struct interp_method * yac_interp_method_creep_new(int creep_distance)
struct interp_method * yac_interp_method_fixed_new(double value)
void yac_interp_weights_delete(struct yac_interp_weights *weights)
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_execute(struct yac_interpolation *interp, double ***src_fields, double **tgt_field)
void yac_interpolation_delete(struct yac_interpolation *interp)
double const YAC_FRAC_MASK_NO_VALUE
@ YAC_LOC_CORNER
Definition location.h:15
@ YAC_LOC_CELL
Definition location.h:14
#define xmalloc(size)
Definition ppm_xfuncs.h:66
enum yac_location location
Definition basic_grid.h:18
struct yac_interp_field tgt_field
Definition interp_grid.c:26
size_t num_src_fields
Definition interp_grid.c:27
struct yac_dist_grid_pair * grid_pair
Definition interp_grid.c:25
struct yac_interp_field src_fields[]
Definition interp_grid.c:28
struct yac_basic_grid ** grids
Definition yac.c:130