19#define YAC_RAD (0.01745329251994329576923690768489)
46static void dummy_compA (
char const * configFilename,
char const * gridPath);
47static void dummy_compB (
char const * configFilename);
49#define STR_USAGE "Usage: %s -c configFilename -g gridPath\n"
50#define YAC_ASSERT_ARGS(exp, msg) \
53 fprintf(stderr, "ERROR: %s\n" STR_USAGE, msg, argv[0]); \
59 int argc,
char ** argv,
char const ** configFilenamem,
char const ** gridPath);
60static inline void LLtoXYZ(
double lon,
double lat,
double p_out[]);
64int main (
int argc,
char *argv[]) {
69 MPI_Comm_rank ( MPI_COMM_WORLD, &rank );
70 MPI_Comm_size ( MPI_COMM_WORLD, &size );
72 YAC_ASSERT(size == 2,
"Too many processes have been launched")
74 char const * configFilename = "toy_coupling.yaml";
75 char const * gridPath = ".";
78 if (rank == 0)
dummy_compA(configFilename, gridPath);
101 int lon_diff = fabs(a_->
p.
lon - b_->
p.
lon) > 1e-9;
102 int lat_diff = fabs(a_->
p.
lat - b_->
p.
lat) > 1e-9;
106 if (a_->
p.
lon > b_->
p.
lon)
return -1;
109 }
else if (lat_diff) {
111 if (a_->
p.
lat > b_->
p.
lat)
return -1;
118 double * temp_vertex_lat,
119 int * temp_nbr_vertices,
120 int * old_to_new_id) {
123 malloc(*temp_nbr_vertices *
sizeof(*sort_array));
125 for (
int i = 0;
i < *temp_nbr_vertices; ++
i) {
127 double curr_lon, curr_lat;
129 curr_lon = temp_vertex_lon[
i];
130 curr_lat = temp_vertex_lat[
i];
132 while (curr_lon < 0.0) curr_lon += 2.0 * M_PI;
133 while (curr_lon >= 2.0 * M_PI) curr_lon -= 2.0 * M_PI;
135 sort_array[
i].
p.
lon = curr_lon;
136 sort_array[
i].
p.
lat = curr_lat;
141 qsort(sort_array, *temp_nbr_vertices,
sizeof(*sort_array),
144 old_to_new_id[sort_array[0].
i] = 1;
146 int last_unique_idx = sort_array[0].
i;
148 for (
int i = 1;
i < *temp_nbr_vertices; ++
i) {
152 old_to_new_id[sort_array[
i].
i] = 1;
153 last_unique_idx = sort_array[
i].i;
157 old_to_new_id[sort_array[
i].
i] = -last_unique_idx;
163 size_t new_nbr_vertices = 0;
165 for (
int i = 0;
i < *temp_nbr_vertices; ++
i) {
167 if (old_to_new_id[
i] == 1) {
169 temp_vertex_lon[new_nbr_vertices] = temp_vertex_lon[
i];
170 temp_vertex_lat[new_nbr_vertices] = temp_vertex_lat[
i];
172 old_to_new_id[
i] = new_nbr_vertices;
178 for (
int i = 0;
i < *temp_nbr_vertices; ++
i)
179 if (old_to_new_id[
i] <= 0)
180 old_to_new_id[
i] = old_to_new_id[-old_to_new_id[
i]];
182 *temp_nbr_vertices = new_nbr_vertices;
186 int * nbr_vertices,
int * nbr_cells,
188 double ** buffer_lonv,
double ** buffer_latv,
200 char * gridFilename =
202 "mesh_compA_exact.txt";
204 "mesh_compA_convex.txt";
207 char * full_gridFilename = malloc(strlen(gridPath) + strlen(gridFilename) + 2);
208 sprintf(full_gridFilename,
"%s/%s", gridPath, gridFilename);
210 char* inFormat0 =
"%s %*s %*s";
211 char* inFormat1 =
"%*10c %lf";
213 char* inFormat2 =
"%*10c %lf %lf %lf %lf %lf %lf %lf %lf %lf %lf %lf %lf";
214 enum {num_entries_in_inFormat2 = 12};
216 char* inFormat2 =
"%*10c %lf %lf %lf %lf %lf %lf %lf";
217 enum {num_entries_in_inFormat2 = 7};
220 FILE *inData = fopen(full_gridFilename,
"r");
222 inData,
"Cannot open file \"%s\"", full_gridFilename);
223 free(full_gridFilename);
225 char dummy_string[32];
229 for (
int i = 0;
i < 3; ++
i )
230 YAC_ASSERT(1 == fscanf (inData, inFormat0, dummy_string),
231 "fscanf failed reading data");
237 *buffer_lonv = malloc ( (*nbr_vertices) *
sizeof(*buffer_lonv) );
238 *buffer_latv = malloc ( (*nbr_vertices) *
sizeof(*buffer_latv) );
240 for (
int i = 0;
i < *nbr_cells; ++
i ) {
242 "fscanf failed reading data")
243 printf (
"point buffer_lat %i : %lf \n",
i, (*
buffer_lon)[
i]);
246 for (
int i = 0;
i < *nbr_cells; ++
i ) {
248 num_entries_in_inFormat2 ==
249 fscanf ( inData, inFormat2,
262 "fscanf failed reading data")
265 for (
int i = 0;
i < *nbr_cells; ++
i ) {
267 "fscanf failed reading data")
268 printf (
"point buffer_lat %i : %lf \n",
i, (*
buffer_lat)[
i]);
271 for (
int i = 0;
i < *nbr_cells; ++
i ) {
273 num_entries_in_inFormat2 ==
274 fscanf ( inData, inFormat2,
287 "fscanf failed reading data")
290 for (
int i = 0;
i < *nbr_vertices; ++
i )
291 printf (
"%i vertex buffer_lonv : %lf buffer_latv : %lf \n",
i,
292 (*buffer_lonv)[
i], (*buffer_latv)[
i]);
294 for (
int i = 0;
i < *nbr_cells; ++
i ) {
299 for (
int i = 0;
i < *nbr_vertices; ++
i ) {
308 (*num_vertices_per_cell) = malloc ( *nbr_cells *
sizeof(**num_vertices_per_cell) );
310 int * old_to_new_id = malloc((*nbr_vertices) *
sizeof(*old_to_new_id));
316 for (
int i = 0;
i < *nbr_cells; ++
i) {
319 (*num_vertices_per_cell)[
i] = 1;
327 (*num_vertices_per_cell)[
i]++;
337static void dummy_compA (
char const * configFilename,
char const * gridPath) {
343 char * comp_name =
"dummy_compA";
344 char * grid_name =
"dummy_compA_grid";
357 MPI_Comm_rank ( local_comm, &
rank );
358 MPI_Comm_size ( local_comm, &size );
360 int nbr_vertices, nbr_cells, * num_vertices_per_cell, *
cell_to_vertex;
371 grid_name, nbr_vertices, nbr_cells, num_vertices_per_cell,
379 for (
int i = 0;
i < nbr_cells; ++
i ) {
423 double *point_set_data[1];
424 double **collection_data[1];
427 collection_data[0] = &(point_set_data[0]);
456 char vtk_filename[32];
458 sprintf(vtk_filename,
"compA_out_%d.vtk",
rank);
460 double point_data[nbr_vertices][3];
461 for (
int i = 0;
i < nbr_vertices; ++
i)
462 LLtoXYZ(buffer_lonv[
i], buffer_latv[
i], point_data[
i]);
469 (
unsigned*)num_vertices_per_cell, nbr_cells);
482 free (num_vertices_per_cell);
489 free ( buffer_lonv );
490 free ( buffer_latv );
512 char * comp_name =
"dummy_compB";
513 char * grid_name =
"dummy_compB_grid";
526 MPI_Comm_rank ( local_comm, &
rank );
527 MPI_Comm_size ( local_comm, &size );
532 int nbr_cells[2] = {200, 320};
537 nbr_vertices[0] = nbr_cells[0] + 1;
538 nbr_vertices[1] = nbr_cells[1] + 1;
540 double * buffer_lonv = malloc ( nbr_vertices[0] *
sizeof(*buffer_lonv) );
541 double * buffer_latv = malloc ( nbr_vertices[1] *
sizeof(*buffer_latv) );
545 buffer_lonv[0] = 265.6 *
YAC_RAD - 0.5 * dx;
546 buffer_latv[0] = 58.5 *
YAC_RAD - 0.5 * dx;
548 for (
int i = 1;
i < nbr_vertices[0]; ++
i )
549 buffer_lonv[
i] = buffer_lonv[
i-1] + dx;
551 for (
int i = 1;
i < nbr_vertices[1]; ++
i )
552 buffer_latv[
i] = buffer_latv[
i-1] + dx;
563 int total_nbr_cells = nbr_cells[0] * nbr_cells[1];
568 for (
int j = 0; j < nbr_cells[1]; ++j ) {
569 for (
int i = 0;
i < nbr_cells[0]; ++
i ) {
586 for (
int i = 1;
i < nbr_cells[0]; ++
i )
589 for (
int i = 1;
i < nbr_cells[1]; ++
i )
603 for (
int i = 0;
i < 2; ++
i )
619 for (
int j = 0; j < nbr_cells[1]; ++j ) {
620 for (
int i = 0;
i < nbr_cells[0]; ++
i ) {
649 double *collection_data[1];
663 double * point_data = malloc(nbr_vertices[0] * nbr_vertices[1] * 3 *
sizeof(*point_data));
664 for (
int i = 0;
i < nbr_vertices[1]; ++
i)
665 for (
int j = 0; j < nbr_vertices[0]; ++j)
666 LLtoXYZ(buffer_lonv[j], buffer_latv[
i],
667 &point_data[3*(
i * nbr_vertices[0] + j)]);
669 unsigned * cell_corners = malloc(total_nbr_cells * 4 *
sizeof(*cell_corners));
670 unsigned * num_points_per_cell = malloc(total_nbr_cells *
sizeof(*num_points_per_cell));
672 for (
int i = 0;
i < total_nbr_cells; ++
i) {
674 unsigned x_index, y_index;
676 y_index =
i / nbr_cells[0];
677 x_index =
i - y_index * nbr_cells[0];
679 cell_corners[
i*4+0] = y_index * (nbr_cells[0] + 1) + x_index;
680 cell_corners[
i*4+1] = y_index * (nbr_cells[0] + 1) + x_index + 1;
681 cell_corners[
i*4+2] = (y_index + 1) * (nbr_cells[0] + 1) + x_index + 1;
682 cell_corners[
i*4+3] = (y_index + 1) * (nbr_cells[0] + 1) + x_index;
684 num_points_per_cell[
i] = 4;
687 char vtk_filename[32];
689 sprintf(vtk_filename,
"compB_out_%d.vtk",
rank);
694 vtk_file, point_data, nbr_vertices[0]*nbr_vertices[1]);
697 vtk_file, cell_corners, num_points_per_cell, total_nbr_cells);
700 vtk_file,
global_index, total_nbr_cells,
"global_cell_id");
703 vtk_file,
recv_buffer, total_nbr_cells,
"cell_in");
712 free(num_points_per_cell);
721 free ( buffer_lonv );
722 free ( buffer_latv );
737 int argc,
char ** argv,
738 char const ** configFilename,
char const ** gridPath) {
741 while ((opt = getopt(argc, argv,
"c:g:")) != -1) {
742 YAC_ASSERT_ARGS((opt ==
'c') || (opt ==
'g'),
"invalid command argument")
746 *configFilename = optarg;
757 while (
lon < -M_PI)
lon += 2.0 * M_PI;
758 while (
lon >= M_PI)
lon -= 2.0 * M_PI;
760 double cos_lat = cos(
lat);
761 p_out[0] = cos_lat * cos(
lon);
762 p_out[1] = cos_lat * sin(
lon);
#define YAC_ASSERT(exp, msg)
struct point_with_index::@5 p
int nbr_vertices_per_cell[NBR_CELLS]
static int compare_point_with_index(const void *a, const void *b)
static void dummy_compA(char const *configFilename, char const *gridPath)
#define YAC_ASSERT_ARGS(exp, msg)
static void remove_duplicated_vertices(double *temp_vertex_lon, double *temp_vertex_lat, int *temp_nbr_vertices, int *old_to_new_id)
static void read_compA_grid_data(char const *gridPath, int *nbr_vertices, int *nbr_cells, int **num_vertices_per_cell, int **cell_to_vertex, double **buffer_lonv, double **buffer_latv, double **buffer_lon, double **buffer_lat)
static void parse_arguments(int argc, char **argv, char const **configFilenamem, char const **gridPath)
static void LLtoXYZ(double lon, double lat, double p_out[])
static void dummy_compB(char const *configFilename)
double(* p)(double lon, double lat)
void yac_vtk_write_cell_scalars_double(YAC_VTK_FILE *vtk_file, double *scalars, unsigned num_cells, char *name)
void yac_vtk_write_cell_scalars_int(YAC_VTK_FILE *vtk_file, int *scalars, unsigned num_cells, char *name)
void yac_vtk_write_point_data(YAC_VTK_FILE *vtk_file, double *point_data, unsigned num_points)
void yac_vtk_close(YAC_VTK_FILE *vtk_file)
void yac_vtk_write_cell_data(YAC_VTK_FILE *vtk_file, unsigned *cell_corners, unsigned *num_points_per_cell, unsigned num_cells)
YAC_VTK_FILE * yac_vtk_open(const char *filename, const char *title)
void yac_cset_global_index(int const *global_index, int location, int grid_id)
int const YAC_LOCATION_CELL
void yac_cfinalize()
Finalises YAC.
void yac_cget_comp_comm(int comp_id, MPI_Comm *comp_comm)
void yac_cget(int const field_id, int collection_size, double **recv_field, int *info, int *ierr)
void yac_cdef_grid_reg2d(const char *grid_name, int nbr_vertices[2], int cyclic[2], double *x_vertices, double *y_vertices, int *grid_id)
void yac_cdef_points_reg2d(int const grid_id, int const *nbr_points, int const located, double const *x_points, double const *y_points, int *point_id)
int const YAC_TIME_UNIT_SECOND
void yac_cdef_grid_unstruct(const char *grid_name, int nbr_vertices, int nbr_cells, int *num_vertices_per_cell, double *x_vertices, double *y_vertices, int *cell_to_vertex, int *grid_id)
void yac_cdef_points_unstruct(int const grid_id, int const nbr_points, int const located, double const *x_points, double const *y_points, int *point_id)
void yac_cput(int const field_id, int const collection_size, double ***const send_field, int *info, int *ierr)
void yac_cread_config_yaml(const char *yaml_filename)
void yac_cdef_calendar(int calendar)
void yac_cset_core_mask(int const *is_core, int location, int grid_id)
int const YAC_PROLEPTIC_GREGORIAN
void yac_cdef_comp(char const *comp_name, int *comp_id)
void yac_cdef_field(char const *name, int const comp_id, int const *point_ids, int const num_pointsets, int collection_size, const char *timestep, int time_unit, int *field_id)
#define YAC_ASSERT_F(exp, format,...)