YetAnotherCoupler 3.2.0_a
Loading...
Searching...
No Matches
couple_config.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#include "config.h"
7#endif
8
9#include <stdlib.h>
10#include <stdio.h>
11#include <string.h>
12#include <math.h>
13
14#include "couple_config.h"
15#include "yac_mpi_common.h"
16#include "yac.h"
17#include "mtime_calendar.h"
18#include "utils_common.h"
19
21 char const * name;
22 char* metadata;
23};
24
26 char const * name;
27 size_t grid_idx;
28 char const * timestep;
29 char * metadata;
32};
33
35
37 size_t num_fields;
38
39 char const * name;
40 char * metadata;
41};
42
67
75
90
92
93 struct yac_couple_config * couple_config =
94 xmalloc(1 * sizeof(*couple_config));
95
96 couple_config->start_datetime = NULL;
97 couple_config->end_datetime = NULL;
98
99 couple_config->grids = NULL;
100 couple_config->num_grids = 0;
101
102 couple_config->components = NULL;
103 couple_config->num_components = 0;
104
105 couple_config->couples = NULL;
106 couple_config->num_couples = 0;
107
108 return couple_config;
109}
110
111static char const * string_dup(char const * string) {
112 return (string != NULL)?strdup(string):NULL;
113}
114
115static void yac_couple_config_field_free(void * field_){
116 struct yac_couple_config_field * field = field_;
117 free((void*)field->name);
118 free((void*)field->timestep);
119 free(field->metadata);
120}
121
122
123static void yac_couple_config_component_free(void * component_) {
124
125 struct yac_couple_config_component * component = component_;
126
127 for (size_t i = 0; i < component->num_fields; ++i)
128 yac_couple_config_field_free(component->fields + i);
129 free(component->fields);
130 free((void*)(component->name));
131 free(component->metadata);
132}
133
134static void yac_couple_config_field_couple_free(void * field_couple_) {
135
136 struct yac_couple_config_field_couple * field_couple = field_couple_;
138 free((void*)field_couple->coupling_period);
139 free((void*)field_couple->weight_file_name);
140 for (size_t i = 0; i < field_couple->num_src_mask_names; ++i)
141 free((void*)field_couple->src_mask_names[i]);
142 free(field_couple->src_mask_names);
143 free(field_couple->tgt_mask_name);
144}
145
146static void yac_couple_config_couple_free(void * couple_) {
147
148 struct yac_couple_config_couple * couple = couple_;
149 for (size_t i = 0; i < couple->num_field_couples; ++i)
151 free(couple->field_couples);
152 couple->field_couples = NULL;
153 couple->num_field_couples = 0;
154}
155
156static void yac_couple_config_grid_free(void * grid_) {
157
158 struct yac_couple_config_grid * grid = grid_;
159
160 free((void*)grid->name);
161 free(grid->metadata);
162}
163
164static int yac_couple_config_grid_compare(void const * a, void const * b) {
165 YAC_ASSERT(((struct yac_couple_config_grid const *)a)->name != NULL &&
166 ((struct yac_couple_config_grid const *)b)->name != NULL,
167 "ERROR(yac_couple_config_grid_compare): "
168 "invalid name (NULL is not allowed)");
169 return strcmp(((struct yac_couple_config_grid const *)a)->name,
170 ((struct yac_couple_config_grid const *)b)->name);
171}
172
173static int yac_couple_config_field_compare(void const * a_, void const * b_) {
174 struct yac_couple_config_field const * a = a_;
175 struct yac_couple_config_field const * b = b_;
176 int ret = (a->grid_idx > b->grid_idx) - (a->grid_idx < b->grid_idx);
177 if (ret) return ret;
179 a->name != NULL && b->name != NULL,
180 "ERROR(yac_couple_config_field_compare): "
181 "invalid name (NULL is not allowed)");
182 return strcmp(a->name, b->name);
183}
184
186 void const * a_, void const * b_) {
187 struct yac_couple_config_component const * a = a_;
188 struct yac_couple_config_component const * b = b_;
189 YAC_ASSERT(a->name != NULL && b->name != NULL,
190 "ERROR(yac_couple_config_component_compare): "
191 "invalid name (NULL is not allowed)");
192 return strcmp(a->name, b->name);
193}
194
196 void const * a_, void const * b_) {
197 struct yac_couple_config_field_couple const * a = a_;
198 struct yac_couple_config_field_couple const * b = b_;
199 int ret;
200 ret = (a->source.component_idx > b->source.component_idx) -
202 if (ret) return ret;
203 ret = (a->target.component_idx > b->target.component_idx) -
205 if (ret) return ret;
206 ret = (a->source.field_idx > b->source.field_idx) -
208 if (ret) return ret;
209 return (a->target.field_idx > b->target.field_idx) -
211}
212
214 void const * a_, void const * b_) {
215 struct yac_couple_config_couple const * a = a_;
216 struct yac_couple_config_couple const * b = b_;
217 int ret;
218 ret = (a->component_indices[0] > b->component_indices[0]) -
219 (a->component_indices[0] < b->component_indices[0]);
220 if (ret) return ret;
221 return (a->component_indices[1] > b->component_indices[1]) -
222 (a->component_indices[1] < b->component_indices[1]);
223}
224
226 char const * string_name, char ** string, MPI_Comm comm) {
227
228 int rank;
229 MPI_Comm_rank(comm, &rank);
230 struct {
231 int len;
232 int rank;
233 } data_pair;
234 size_t len = (*string != NULL)?(strlen(*string) + 1):0;
236 len <= INT_MAX,
237 "ERROR(couple_config_sync_string): \"%s\" too long", string_name);
238 data_pair.len = (int)len;
239 data_pair.rank = rank;
240
241 // determine the broadcaster
243 MPI_Allreduce(
244 MPI_IN_PLACE, &data_pair, 1, MPI_2INT, MPI_MAXLOC, comm), comm);
245 if (data_pair.len == 0) return;
246
247 // broadcast string
248 char const * string_bak = NULL;
249 if (data_pair.rank != rank) {
250 string_bak = *string;
251 *string = xmalloc((size_t)data_pair.len * sizeof(**string));
252 }
254 MPI_Bcast(*string, data_pair.len, MPI_CHAR, data_pair.rank, comm), comm);
255
256 // check for consistency
257 if (data_pair.rank != rank) {
259 (string_bak == NULL) ||
260 !strcmp(string_bak, *string),
261 "ERROR(couple_config_sync_string): inconsistent \"%s\" definition "
262 "(\"%s\" != \"%s\")", string_name, string_bak, *string);
263 free((void*)string_bak);
264 }
265}
266
267static void couple_config_sync_calendar(MPI_Comm comm) {
268
269 int calendar = (int)getCalendarType();
270 if (calendar == CALENDAR_NOT_SET) calendar = INT_MAX;
271
272 // broadcast calendar
274 MPI_Allreduce(
275 MPI_IN_PLACE, &calendar, 1, MPI_INT, MPI_MIN, comm), comm);
276
277 // if no process has defined a calendar
278 if (calendar == INT_MAX) return;
279
280 // set calendar (if local process has already defined a calendar,
281 // the definition is checked for consistency)
282 yac_cdef_calendar(calendar);
283}
284
286 void * a_, void * b_, MPI_Comm comm) {
287
288 struct yac_couple_config_component * a = a_;
289 struct yac_couple_config_component * b = b_;
290
291 if (b) {
292 a->num_fields = b->num_fields;
293 a->fields = b->fields;
294 b->num_fields = 0;
295 b->fields = NULL;
296 a->metadata = b->metadata;
297 b->metadata = NULL;
298 }
299
301 "component metadata", (char **)&(a->metadata), comm);
302}
303
305 void * a_, void * b_, MPI_Comm comm) {
306
307 struct yac_couple_config_grid * a = a_;
308 struct yac_couple_config_grid * b = b_;
309
310 if (b!= NULL) {
311 a->metadata = b->metadata;
312 b->metadata = NULL;
313 }
314
315 couple_config_sync_string("grid metadata", (char **)&(a->metadata), comm);
316}
317
319 void * a_, void * b_, MPI_Comm comm) {
320
321 struct yac_couple_config_field * a = a_;
322 struct yac_couple_config_field * b = b_;
323
324 enum {
325 TIMESTEP_IDX,
326 METADATA_IDX,
327 FRAC_MASK_IDX,
328 COLLECTION_SIZE_IDX,
329 DATA_PAIR_COUNT
330 };
331
332 int rank;
333 MPI_Comm_rank(comm, &rank);
334 struct {
335 int len;
336 int rank;
337 } data_pairs[DATA_PAIR_COUNT];
338 size_t timestep_len =
339 ((b != NULL) && (b->timestep != NULL))?(strlen(b->timestep) + 1):0;
340 size_t metadata_len =
341 ((b != NULL) && (b->metadata != NULL))?(strlen(b->metadata) + 1):0;
343 timestep_len <= INT_MAX,
344 "ERROR(yac_couple_config_field_merge): timestep string too long");
346 metadata_len <= INT_MAX,
347 "ERROR(yac_couple_config_field_merge): metadata string too long");
349 (b == NULL) || (b->collection_size <= INT_MAX) ||
350 (b->collection_size == SIZE_MAX),
351 "ERROR(yac_couple_config_field_merge): invalid collection size \"%zu\"",
352 b->collection_size);
353 data_pairs[TIMESTEP_IDX].len = (int)timestep_len;
354 data_pairs[TIMESTEP_IDX].rank = rank;
355 data_pairs[METADATA_IDX].len = (int)metadata_len;
356 data_pairs[METADATA_IDX].rank = rank;
357 data_pairs[FRAC_MASK_IDX].len =
359 data_pairs[FRAC_MASK_IDX].rank = rank;
360 data_pairs[COLLECTION_SIZE_IDX].len =
361 ((b == NULL) || (b->collection_size == SIZE_MAX))?
362 -1:(int)b->collection_size;
363 data_pairs[COLLECTION_SIZE_IDX].rank = rank;
364
365 // determine the broadcaster
367 MPI_Allreduce(
368 MPI_IN_PLACE, data_pairs, DATA_PAIR_COUNT, MPI_2INT,
369 MPI_MAXLOC, comm), comm);
370
371 // if at least one process has a valid timestep
372 if (data_pairs[TIMESTEP_IDX].len > 0) {
373
374 // broadcast timestep
375 char * timestep_buffer = xmalloc((size_t)data_pairs[TIMESTEP_IDX].len);
376 if (data_pairs[TIMESTEP_IDX].rank == rank)
377 strcpy(timestep_buffer, b->timestep);
379 MPI_Bcast(
380 timestep_buffer, data_pairs[TIMESTEP_IDX].len, MPI_CHAR,
381 data_pairs[TIMESTEP_IDX].rank, comm), comm);
382
383 // check for consistency
385 (b == NULL) || (b->timestep == NULL) ||
386 !strcmp(b->timestep, timestep_buffer),
387 "ERROR(yac_couple_config_field_merge): "
388 "inconsistent timestep definition (\"%s\" != \"%s\")",
389 b->timestep, timestep_buffer);
390
391 // update timestep
392 free((void*)(a->timestep));
393 a->timestep = timestep_buffer;
394 }
395
396 // if at least one process has a valid metadata
397 if (data_pairs[METADATA_IDX].len > 0) {
398
399 // broadcast metadata
400 char * metadata_buffer = xmalloc((size_t)data_pairs[METADATA_IDX].len);
401 if (data_pairs[METADATA_IDX].rank == rank)
402 strcpy(metadata_buffer, b->metadata);
404 MPI_Bcast(
405 metadata_buffer, data_pairs[METADATA_IDX].len, MPI_CHAR,
406 data_pairs[METADATA_IDX].rank, comm), comm);
407
408 // check for consistency
410 (b == NULL) || (b->metadata == NULL) ||
411 !strcmp(b->metadata, metadata_buffer),
412 "ERROR(yac_couple_config_field_merge): "
413 "inconsistent metadata definition (\"%s\" != \"%s\")",
414 b->metadata, metadata_buffer);
415
416 // update metadata
417 free(a->metadata);
418 a->metadata = metadata_buffer;
419 }
420
421 // if at least one process has a valid fractional mask fallback value
422 if (data_pairs[FRAC_MASK_IDX].len != 0) {
423
424 // broadcast fractional mask fallback value
425 double frac_mask_fallback_value;
426 if (data_pairs[FRAC_MASK_IDX].rank == rank)
427 frac_mask_fallback_value = b->frac_mask_fallback_value;
429 MPI_Bcast(
430 &frac_mask_fallback_value, 1, MPI_DOUBLE,
431 data_pairs[FRAC_MASK_IDX].rank, comm), comm);
432
433 // check for consistency
434 // (use memcmp for comparison, because it can be nan)
436 (b == NULL) || (b->frac_mask_fallback_value == YAC_FRAC_MASK_NO_VALUE) ||
437 !memcmp(
438 &b->frac_mask_fallback_value, &frac_mask_fallback_value,
439 sizeof(double)),
440 "ERROR(yac_couple_config_field_merge): "
441 "inconsistent fractional mask fallback value definition "
442 "(%lf != %lf)",
443 b->frac_mask_fallback_value, frac_mask_fallback_value);
444
445 // update fractional mask fallback value
446 a->frac_mask_fallback_value = frac_mask_fallback_value;
447 }
448
449
450 // if at least one process has a valid collection size
451 if (data_pairs[COLLECTION_SIZE_IDX].len > -1) {
452
453 // check for consistency
455 (b == NULL) || (b->collection_size == SIZE_MAX) ||
456 (b->collection_size == (size_t)(data_pairs[COLLECTION_SIZE_IDX].len)),
457 "ERROR(yac_couple_config_field_merge): "
458 "inconsistent collection size definition (\"%zu\" != \"%d\")",
459 b->collection_size, data_pairs[COLLECTION_SIZE_IDX].len);
460
461 // update collection size
462 a->collection_size = (size_t)(data_pairs[COLLECTION_SIZE_IDX].len);
463 }
464}
465
467 void * a_, void * b_, MPI_Comm comm) {
468
469 if (b_ == NULL) return;
470
471 struct yac_couple_config_field_couple * a = a_;
472 struct yac_couple_config_field_couple * b = b_;
473
475 (a->source.lag == b->source.lag),
476 "ERROR(yac_couple_config_field_couple_merge): "
477 "inconsistent source lag (%d != %d)", a->source.lag, b->source.lag)
479 (a->target.lag == b->target.lag),
480 "ERROR(yac_couple_config_field_couple_merge): "
481 "inconsistent target lag (%d != %d)", a->target.lag, b->target.lag)
484 "ERROR(yac_couple_config_field_couple_merge): "
485 "inconsistent mapping side (%d != %d)",
489 "ERROR(yac_couple_config_field_couple_merge): "
490 "inconsistent interpolation stack")
493 "ERROR(yac_couple_config_field_couple_merge): "
494 "inconsistent coupling period operation (%d != %d)",
497 !strcmp(a->coupling_period, b->coupling_period),
498 "ERROR(yac_couple_config_field_couple_merge): "
499 "inconsistent coupling period (%s != %s)",
503 "ERROR(yac_couple_config_field_couple_merge): "
504 "inconsistent defintion of enforce_write_weight_file (%d != %d)",
508 !strcmp(a->weight_file_name, b->weight_file_name),
509 "ERROR(yac_couple_config_field_couple_merge): "
510 "inconsistent weight_file_name (%s != %s)",
513 a->scale_factor == b->scale_factor,
514 "ERROR(yac_couple_config_field_couple_merge): "
515 "inconsistent scale factor (%lf != %lf)",
519 "ERROR(yac_couple_config_field_couple_merge): "
520 "inconsistent scale summand (%lf != %lf)",
524 "ERROR(yac_couple_config_field_couple_merge): "
525 "inconsistent number of source mask names (%zu != %zu)",
528 (a->src_mask_names != NULL) == (b->src_mask_names != NULL),
529 "ERROR(yac_couple_config_field_couple_merge): "
530 "inconsistent availability of source mask names (%s != %s)",
531 (a->src_mask_names != NULL)?"TRUE":"FALSE",
532 (b->src_mask_names != NULL)?"TRUE":"FALSE")
533 for (size_t i = 0; i < a->num_src_mask_names; ++i)
535 !strcmp(a->src_mask_names[i], b->src_mask_names[i]),
536 "ERROR(yac_couple_config_field_couple_merge): "
537 "inconsistent source mask names at index %zu (\"%s\" != \"%s\")",
538 i, a->src_mask_names[i], b->src_mask_names[i])
540 (a->tgt_mask_name != NULL) == (b->tgt_mask_name != NULL),
541 "ERROR(yac_couple_config_field_couple_merge): "
542 "inconsistent availability of target mask name (%s != %s)",
543 (a->tgt_mask_name != NULL)?"TRUE":"FALSE",
544 (b->tgt_mask_name != NULL)?"TRUE":"FALSE")
546 (a->tgt_mask_name == NULL) ||
547 !strcmp(a->tgt_mask_name, b->tgt_mask_name),
548 "ERROR(yac_couple_config_field_couple_merge): "
549 "inconsistent target mask name (\"%s\" != \"%s\")",
551}
552
553static void merge_field_couples(
554 size_t * num_field_couples,
555 struct yac_couple_config_field_couple ** field_couples, MPI_Comm comm);
556
558 void * a_, void * b_, MPI_Comm comm) {
559
560 struct yac_couple_config_couple * a = a_;
561 struct yac_couple_config_couple * b = b_;
562
563 if (b) {
566 b->num_field_couples = 0;
567 b->field_couples = NULL;
568 }
569
570 // distribute and merge field couples
572}
573
574void yac_couple_config_delete(struct yac_couple_config * couple_config) {
575
576 free((void*)couple_config->start_datetime);
577 free((void*)couple_config->end_datetime);
578
579 for (size_t i = 0; i < couple_config->num_grids; ++i){
580 yac_couple_config_grid_free(couple_config->grids + i);
581 }
582 free(couple_config->grids);
583
584 for (size_t i = 0;
585 i < couple_config->num_components; ++i)
587 free(couple_config->components);
588
589 for (size_t couple_idx = 0; couple_idx < couple_config->num_couples;
590 ++couple_idx)
591 yac_couple_config_couple_free(couple_config->couples + couple_idx);
592 free(couple_config->couples);
593 free(couple_config);
594}
595
597 struct yac_couple_config * couple_config, char const * name) {
598
600 name != NULL,
601 "ERROR(yac_couple_config_add_grid_): "
602 "invalid name (NULL is not allowed)")
603
604 for (size_t i = 0; i < couple_config->num_grids; ++i)
605 if (!strcmp(couple_config->grids[i].name, name)) return i;
606
607 size_t grid_idx = couple_config->num_grids;
608 couple_config->num_grids++;
609 couple_config->grids =
610 xrealloc(
611 couple_config->grids,
612 couple_config->num_grids * sizeof(*(couple_config->grids)));
613 couple_config->grids[grid_idx].name = string_dup(name);
614 couple_config->grids[grid_idx].metadata = NULL;
615 return grid_idx;
616}
617
619 struct yac_couple_config * couple_config, char const * name) {
620
621 yac_couple_config_add_grid_(couple_config, name);
622}
623
625 struct yac_couple_config * couple_config, char const * name) {
626
628 name != NULL,
629 "ERROR(yac_couple_config_add_component_): "
630 "invalid name (NULL is not allowed)")
631
632 for (size_t i = 0; i < couple_config->num_components; ++i)
633 if (!strcmp(couple_config->components[i].name, name)) return i;
634
635 size_t component_idx = couple_config->num_components;
636 couple_config->num_components++;
637 couple_config->components =
638 xrealloc(
639 couple_config->components,
640 couple_config->num_components * sizeof(*(couple_config->components)));
641 couple_config->components[component_idx].fields = NULL;
642 couple_config->components[component_idx].num_fields = 0;
643 couple_config->components[component_idx].name = strdup(name);
644 couple_config->components[component_idx].metadata = NULL;
645 return component_idx;
646}
647
649 struct yac_couple_config * couple_config, char const * name) {
650
651 yac_couple_config_add_component_(couple_config, name);
652}
653
655 struct yac_couple_config * couple_config,
656 char const * comp_name, const char* metadata) {
657 size_t comp_idx = yac_couple_config_get_component_idx(couple_config, comp_name);
658 free(couple_config->components[comp_idx].metadata);
659 couple_config->components[comp_idx].metadata
660 = metadata==NULL?NULL:strdup(metadata);
661}
662
664 struct yac_couple_config * couple_config,
665 char const * grid_name, const char* metadata) {
666 if(!yac_couple_config_contains_grid_name(couple_config, grid_name))
667 yac_couple_config_add_grid(couple_config, grid_name);
668 size_t grid_idx = yac_couple_config_get_grid_idx(couple_config, grid_name);
669 free(couple_config->grids[grid_idx].metadata);
670 couple_config->grids[grid_idx].metadata
671 = metadata==NULL?NULL:strdup(metadata);
672}
673
675 struct yac_couple_config * couple_config,
676 const char* comp_name, const char * grid_name, const char* field_name,
677 const char* metadata) {
678 size_t comp_idx = yac_couple_config_get_component_idx(couple_config, comp_name);
679 size_t grid_idx = yac_couple_config_get_grid_idx(couple_config, grid_name);
680 size_t field_idx = yac_couple_config_get_field_idx(couple_config,
681 comp_idx, grid_idx, field_name);
682 free(couple_config->components[comp_idx].fields[field_idx].metadata);
683 couple_config->components[comp_idx].fields[field_idx].metadata
684 = metadata==NULL?NULL:strdup(metadata);
685}
686
688 struct yac_couple_config * couple_config,
689 const char * comp_name) {
690 size_t comp_idx = yac_couple_config_get_component_idx(couple_config, comp_name);
691 return couple_config->components[comp_idx].metadata;
692}
693
695 struct yac_couple_config * couple_config,
696 const char * grid_name) {
697 size_t grid_idx = yac_couple_config_get_grid_idx(couple_config, grid_name);
698 return couple_config->grids[grid_idx].metadata;
699}
700
702 struct yac_couple_config * couple_config,
703 const char* comp_name, const char * grid_name, const char* field_name) {
704 size_t comp_idx = yac_couple_config_get_component_idx(couple_config, comp_name);
705 size_t grid_idx = yac_couple_config_get_grid_idx(couple_config, grid_name);
706 size_t field_idx = yac_couple_config_get_field_idx(couple_config,
707 comp_idx, grid_idx, field_name);
708 return couple_config->components[comp_idx].fields[field_idx].metadata;
709}
710
712 struct yac_couple_config * couple_config, size_t component_idx,
713 char const * routine_name, int line) {
714
716 component_idx < couple_config->num_components,
717 "ERROR(%s:%d:%s): invalid component_idx", __FILE__, line, routine_name)
718}
719
720static void check_grid_idx(
721 struct yac_couple_config * couple_config, size_t grid_idx,
722 char const * routine_name, int line) {
723
725 grid_idx < couple_config->num_grids,
726 "ERROR(%s:%d:%s): invalid grid_idx", __FILE__, line, routine_name)
727}
728
730 struct yac_couple_config * couple_config,
731 size_t comp_idx, size_t grid_idx, char const * name,
732 char const * timestep, size_t collection_size) {
733
735 couple_config, comp_idx,
736 "yac_couple_config_component_add_field_", __LINE__);
738 couple_config, grid_idx,
739 "yac_couple_config_component_add_field_", __LINE__);
740
741 // check whether the field already exists
742 struct yac_couple_config_component * component =
743 couple_config->components + comp_idx;
744 for (size_t i = 0; i < component->num_fields; i++) {
745 if((strcmp(component->fields[i].name, name) == 0) &&
746 (component->fields[i].grid_idx == grid_idx)) {
747 // if no timestep is defined for the field
748 if (!component->fields[i].timestep && timestep)
749 component->fields[i].timestep = strdup(timestep);
750 // if no collection size is defined for the field
751 if (component->fields[i].collection_size == SIZE_MAX)
752 component->fields[i].collection_size = collection_size;
754 !timestep ||
755 !strcmp(timestep, component->fields[i].timestep),
756 "ERROR(yac_couple_config_component_add_field): "
757 "inconsistent timestep definition (\"%s\" != \"%s\")",
758 timestep, component->fields[i].timestep);
760 collection_size == SIZE_MAX ||
761 collection_size == component->fields[i].collection_size,
762 "ERROR(yac_couple_config_component_add_field): "
763 "inconsistent collection_size definition (%zu != %zu)",
764 collection_size, component->fields[i].collection_size);
765 return i;
766 }
767 }
768
769 size_t field_idx = component->num_fields;
770 component->num_fields++;
771
772 component->fields =
773 xrealloc(
774 component->fields,
775 component->num_fields *
776 sizeof(*(component->fields)));
777 struct yac_couple_config_field * field =
778 component->fields + field_idx;
779
780 field->name = strdup(name);
781 field->grid_idx = grid_idx;
782 field->timestep = timestep?strdup(timestep):NULL;
783 field->metadata = NULL;
786 return field_idx;
787}
788
790 struct yac_couple_config * couple_config, const char* component_name,
791 const char* grid_name, const char* name, char const * timestep,
792 size_t collection_size) {
793
795 couple_config,
796 yac_couple_config_get_component_idx(couple_config, component_name),
797 yac_couple_config_get_grid_idx(couple_config, grid_name),
799}
800
802 struct yac_couple_config * couple_config) {
803
804 return couple_config->num_couples;
805}
806
808 struct yac_couple_config * couple_config, size_t couple_idx,
809 char const * routine_name, int line) {
810
812 couple_idx < couple_config->num_couples,
813 "ERROR(%s:%d:%s): invalid couple_idx", __FILE__, line, routine_name);
814}
815
817 struct yac_couple_config * couple_config, size_t couple_idx) {
818
820 couple_config, couple_idx, "yac_couple_config_get_num_couple_fields",
821 __LINE__);
822
823 return couple_config->couples[couple_idx].num_field_couples;
824}
825
827 struct yac_couple_config * couple_config, size_t couple_idx,
828 char const * couple_component_names[2]) {
829
831 couple_config, couple_idx, "yac_couple_config_get_couple_component_names",
832 __LINE__);
833
834 for (int i = 0; i < 2; ++i)
835 couple_component_names[i] =
836 couple_config->components[
837 couple_config->couples[couple_idx].component_indices[i]].name;
838}
839
841 struct yac_couple_config * couple_config, char const * component_name) {
842
844 component_name,
845 "ERROR(yac_couple_config_component_name_is_valid): component name is NULL")
847 strlen(component_name) <= YAC_MAX_CHARLEN,
848 "ERROR(yac_couple_config_component_name_is_valid): "
849 "component name is too long (maximum is YAC_MAX_CHARLEN)")
850
851 for (size_t component_idx = 0; component_idx < couple_config->num_components;
852 ++component_idx)
853 if (!strcmp(component_name, couple_config->components[component_idx].name))
854 return 1;
855 return 0;
856}
857
859 struct yac_couple_config * couple_config) {
860
861 return couple_config->num_components;
862}
863
864
866 struct yac_couple_config * couple_config) {
867
868 return couple_config->num_grids;
869}
870
872 struct yac_couple_config * couple_config, size_t component_idx) {
873 check_component_idx(couple_config, component_idx,
874 "yac_couple_config_get_num_fields", __LINE__);
875 return couple_config->components[component_idx].num_fields;
876}
877
879 struct yac_couple_config * couple_config, char const * component_name) {
880
881 size_t component_idx = SIZE_MAX;
882 for (size_t i = 0; (i < couple_config->num_components) &&
883 (component_idx == SIZE_MAX); ++i)
884 if (!strcmp(couple_config->components[i].name, component_name))
885 component_idx = i;
886
888 component_idx != SIZE_MAX,
889 "ERROR(yac_couple_config_get_component_idx): "
890 "Component \"%s\" not found in coupling configuration",
891 component_name);
892
893 return component_idx;
894}
895
897 struct yac_couple_config * couple_config, char const * grid_name) {
898
899 size_t grid_idx = SIZE_MAX;
900 for (size_t i = 0;
901 (i < couple_config->num_grids) && (grid_idx == SIZE_MAX); ++i)
902 if (!strcmp(couple_config->grids[i].name, grid_name))
903 grid_idx = i;
904
906 grid_idx != SIZE_MAX,
907 "ERROR(yac_couple_config_get_grid_idx): "
908 "grid name \"%s\" not in list of grids", grid_name)
909
910 return grid_idx;
911}
912
914 struct yac_couple_config * couple_config, size_t component_idx,
915 size_t grid_idx, char const * field_name) {
917 couple_config, component_idx,
918 "yac_couple_config_get_component_name", __LINE__);
919
920 size_t field_idx = SIZE_MAX;
921 struct yac_couple_config_component * component =
922 couple_config->components + component_idx;
923 size_t nbr_fields = component->num_fields;
924 for(size_t i=0;
925 (i<nbr_fields) && (field_idx == SIZE_MAX); ++i)
926 if((component->fields[i].grid_idx == grid_idx) &&
927 !strcmp(
928 field_name,
929 component->fields[i].name))
930 field_idx = i;
931
933 field_idx != INT_MAX,
934 "ERROR(yac_couple_config_get_field_idx): "
935 "field not found "
936 "(component_idx %zu grid_idx %zu field_name \"%s\"",
937 component_idx, grid_idx, field_name);
938
939 return field_idx;
940}
941
943 struct yac_couple_config * couple_config, size_t component_idx) {
944
946 couple_config, component_idx,
947 "yac_couple_config_get_component_name", __LINE__);
948
949 return couple_config->components[component_idx].name;
950}
951
952static void check_field_idx(
953 struct yac_couple_config * couple_config, size_t component_idx,
954 size_t field_idx, char const * routine_name, int line) {
955
957 couple_config, component_idx, routine_name, __LINE__);
958
960 field_idx <
961 couple_config->components[component_idx].num_fields,
962 "ERROR(%s:%d:%s): invalid field_idx", __FILE__,
963 line, routine_name)
964}
965
967 struct yac_couple_config * couple_config, size_t component_idx,
968 size_t field_idx) {
969
971 couple_config, component_idx, field_idx,
972 "yac_couple_config_get_field_grid_name", __LINE__);
973
974 return
975 couple_config->grids[
976 couple_config->
977 components[component_idx].
978 fields[field_idx].grid_idx].name;
979}
980
982 struct yac_couple_config * couple_config, size_t component_idx,
983 size_t field_idx) {
984
986 couple_config, component_idx, field_idx,
987 "yac_couple_config_get_field_name", __LINE__);
988
989 return
990 couple_config->
991 components[component_idx].
992 fields[field_idx].name;
993}
994
996 struct yac_couple_config * couple_config,
997 char const * component_name, char const * grid_name,
998 char const * field_name) {
999
1000 size_t component_idx =
1001 yac_couple_config_get_component_idx(couple_config, component_name);
1002 size_t grid_idx =
1003 yac_couple_config_get_grid_idx(couple_config, grid_name);
1004 size_t field_idx =
1006 couple_config, component_idx, grid_idx, field_name);
1007
1008 struct yac_couple_config_field * field =
1009 couple_config->components[component_idx].fields +
1010 field_idx;
1011
1013 field->timestep,
1014 "ERROR(yac_couple_config_get_field_timestep): "
1015 "no valid timestep defined (component: \"%s\" field \"%s\")",
1016 couple_config->components[component_idx].name, field->name);
1017
1018 return field->timestep;
1019}
1020
1022 struct yac_couple_config * couple_config,
1023 char const * component_name, char const * grid_name,
1024 char const * field_name) {
1025
1026 size_t component_idx =
1027 yac_couple_config_get_component_idx(couple_config, component_name);
1028 size_t grid_idx =
1029 yac_couple_config_get_grid_idx(couple_config, grid_name);
1030 size_t field_idx =
1032 couple_config, component_idx, grid_idx, field_name);
1033
1034 size_t nbr_couples = couple_config->num_couples;
1035 for(size_t couple_idx = 0; couple_idx<nbr_couples; ++couple_idx){
1036 struct yac_couple_config_couple * couple = couple_config->couples + couple_idx;
1037 if(couple->component_indices[0] != component_idx &&
1038 couple->component_indices[1] != component_idx)
1039 continue;
1040 size_t nbr_trans_couples = couple->num_field_couples;
1041 for(size_t trans_couple_idx = 0; trans_couple_idx < nbr_trans_couples;
1042 ++trans_couple_idx){
1043 struct yac_couple_config_field_couple * transcouple =
1044 couple->field_couples + trans_couple_idx;
1045 if(transcouple->source.component_idx == component_idx &&
1046 transcouple->source.field_idx == field_idx)
1048 if(transcouple->target.component_idx == component_idx &&
1049 transcouple->target.field_idx == field_idx)
1051 }
1052 }
1054}
1055
1057 struct yac_couple_config * couple_config,
1058 size_t component_idx, size_t field_idx) {
1059
1061 couple_config, component_idx, field_idx,
1062 "yac_couple_config_field_is_valid", __LINE__);
1063
1064 struct yac_couple_config_field * field =
1065 couple_config->components[component_idx].fields +
1066 field_idx;
1067
1068 return (field->collection_size != SIZE_MAX) &&
1069 (field->timestep != NULL);
1070}
1071
1073 struct yac_couple_config * couple_config, size_t couple_idx,
1074 size_t field_couple_idx, char const * routine_name, int line) {
1075
1076 check_couple_idx(couple_config, couple_idx, routine_name, __LINE__);
1077
1079 field_couple_idx <
1080 couple_config->couples[couple_idx].num_field_couples,
1081 "ERROR(%s:%d:%s): invalid field_couple_idx",
1082 __FILE__, line, routine_name)
1083}
1084
1086 struct yac_couple_config * couple_config,
1087 size_t couple_idx, size_t field_couple_idx) {
1088
1090 couple_config, couple_idx, field_couple_idx,
1091 "yac_couple_config_get_interp_stack", __LINE__);
1092
1093 return
1094 couple_config->
1095 couples[couple_idx].
1096 field_couples[field_couple_idx].interp_stack;
1097}
1098
1100 struct yac_couple_config * couple_config,
1101 size_t couple_idx, size_t field_couple_idx,
1102 char const ** src_grid_name, char const ** tgt_grid_name) {
1103
1105 couple_config, couple_idx, field_couple_idx,
1106 "yac_couple_config_get_field_grid_names", __LINE__);
1107
1108 size_t src_component_idx =
1109 couple_config->couples[couple_idx].
1110 field_couples[field_couple_idx].
1111 source.component_idx;
1112 size_t src_field_idx =
1113 couple_config->couples[couple_idx].
1114 field_couples[field_couple_idx].
1115 source.field_idx;
1116
1117 size_t tgt_component_idx =
1118 couple_config->couples[couple_idx].
1119 field_couples[field_couple_idx].
1120 target.component_idx;
1121 size_t tgt_field_idx =
1122 couple_config->couples[couple_idx].
1123 field_couples[field_couple_idx].
1124 target.field_idx;
1125
1126 *src_grid_name =
1127 couple_config->grids[
1128 couple_config->components[src_component_idx].
1129 fields[src_field_idx].grid_idx].name;
1130 *tgt_grid_name =
1131 couple_config->grids[
1132 couple_config->components[tgt_component_idx].
1133 fields[tgt_field_idx].grid_idx].name;
1134}
1135
1137 struct yac_couple_config * couple_config,
1138 char const * comp_name, char const * grid_name, char const * field_name,
1139 double frac_mask_fallback_value) {
1140
1142 (frac_mask_fallback_value != YAC_FRAC_MASK_UNDEF) &&
1143 (frac_mask_fallback_value != YAC_FRAC_MASK_NO_VALUE),
1144 "ERROR(yac_couple_config_field_enable_frac_mask): "
1145 "\"%lf\" is not a valid fractional mask fallback value "
1146 "(component: \"%s\" grid: \"%s\" field \"%s\")",
1147 frac_mask_fallback_value, comp_name, grid_name, field_name);
1148
1149 size_t comp_idx = yac_couple_config_get_component_idx(couple_config, comp_name);
1150 size_t grid_idx = yac_couple_config_get_grid_idx(couple_config, grid_name);
1151 size_t field_idx = yac_couple_config_get_field_idx(couple_config,
1152 comp_idx, grid_idx, field_name);
1153
1154 double old_frac_mask_fallback_value =
1155 couple_config->components[comp_idx].fields[field_idx].
1156 frac_mask_fallback_value;
1157
1159 (old_frac_mask_fallback_value == YAC_FRAC_MASK_UNDEF) ||
1160 (old_frac_mask_fallback_value == YAC_FRAC_MASK_NO_VALUE) ||
1161 (old_frac_mask_fallback_value == frac_mask_fallback_value),
1162 "ERROR(yac_couple_config_field_enable_frac_mask): "
1163 "fractional mask fallback value was already set:\n"
1164 "\told value \"%lf\" new value \"%lf\" "
1165 "(component: \"%s\" grid: \"%s\" field \"%s\")",
1166 old_frac_mask_fallback_value, frac_mask_fallback_value,
1167 comp_name, grid_name, field_name);
1168
1169 couple_config->components[comp_idx].fields[field_idx].
1170 frac_mask_fallback_value = frac_mask_fallback_value;
1171}
1172
1174 struct yac_couple_config * couple_config,
1175 char const * component_name, char const * grid_name,
1176 char const * field_name) {
1177
1178 size_t component_idx =
1179 yac_couple_config_get_component_idx(couple_config, component_name);
1180 size_t grid_idx =
1181 yac_couple_config_get_grid_idx(couple_config, grid_name);
1182 size_t field_idx =
1184 couple_config, component_idx, grid_idx, field_name);
1185
1186 struct yac_couple_config_field * field =
1187 couple_config->components[component_idx].fields +
1188 field_idx;
1189
1192 "ERROR(yac_couple_config_get_frac_mask_fallback_value): "
1193 "no valid fractional mask fallback value defined "
1194 "(component: \"%s\" grid: \"%s\" field \"%s\")",
1195 couple_config->components[component_idx].name, grid_name, field->name);
1196
1197 return field->frac_mask_fallback_value;
1198}
1199
1201 struct yac_couple_config * couple_config,
1202 char const * component_name, char const * grid_name,
1203 char const * field_name) {
1204
1205 size_t component_idx =
1206 yac_couple_config_get_component_idx(couple_config, component_name);
1207 size_t grid_idx =
1208 yac_couple_config_get_grid_idx(couple_config, grid_name);
1209 size_t field_idx =
1211 couple_config, component_idx, grid_idx, field_name);
1212
1213 struct yac_couple_config_field * field =
1214 couple_config->components[component_idx].fields +
1215 field_idx;
1216
1218 field->collection_size != SIZE_MAX,
1219 "ERROR(yac_couple_config_get_collection_size): "
1220 "no valid collection size defined "
1221 "(component: \"%s\" grid: \"%s\" field \"%s\")",
1222 couple_config->components[component_idx].name, grid_name, field->name);
1223
1224 return field->collection_size;
1225}
1226
1228 struct yac_couple_config * couple_config,
1229 size_t couple_idx, size_t field_couple_idx,
1230 char const ** src_component_name, char const ** tgt_component_name) {
1231
1233 couple_config, couple_idx, field_couple_idx,
1234 "yac_couple_config_get_field_couple_component_names", __LINE__);
1235
1236 *src_component_name =
1237 couple_config->components[
1238 couple_config->couples[couple_idx].
1239 field_couples[field_couple_idx].source.component_idx].name;
1240 *tgt_component_name =
1241 couple_config->components[
1242 couple_config->couples[couple_idx].
1243 field_couples[field_couple_idx].target.component_idx].name;
1244}
1245
1247 struct yac_couple_config * couple_config,
1248 size_t couple_idx, size_t field_couple_idx,
1249 char const ** src_field_name, const char ** tgt_field_name) {
1250
1252 couple_config, couple_idx, field_couple_idx,
1253 "yac_couple_config_get_field_names", __LINE__);
1254
1255 size_t src_component_idx =
1256 couple_config->couples[couple_idx].
1257 field_couples[field_couple_idx].source.component_idx;
1258 size_t tgt_component_idx =
1259 couple_config->couples[couple_idx].
1260 field_couples[field_couple_idx].target.component_idx;
1261 size_t src_field_idx =
1262 couple_config->couples[couple_idx].
1263 field_couples[field_couple_idx].source.field_idx;
1264 size_t tgt_field_idx =
1265 couple_config->couples[couple_idx].
1266 field_couples[field_couple_idx].target.field_idx;
1267
1268 *src_field_name =
1269 couple_config->components[src_component_idx].
1270 fields[src_field_idx].name;
1271 *tgt_field_name =
1272 couple_config->components[tgt_component_idx].
1273 fields[tgt_field_idx].name;
1274}
1275
1277 struct yac_couple_config * couple_config,
1278 size_t couple_idx, size_t field_couple_idx) {
1279
1281 couple_config, couple_idx, field_couple_idx,
1282 "yac_couple_config_mapping_on_source", __LINE__);
1283
1284 return
1285 couple_config->
1286 couples[couple_idx].
1287 field_couples[field_couple_idx].
1288 mapping_on_source;
1289}
1290
1292 struct yac_couple_config * couple_config,
1293 size_t couple_idx, size_t field_couple_idx) {
1294
1296 couple_config, couple_idx, field_couple_idx,
1297 "yac_couple_config_get_source_lag", __LINE__);
1298
1299 return
1300 couple_config->
1301 couples[couple_idx].
1302 field_couples[field_couple_idx].
1303 source.lag;
1304}
1305
1307 struct yac_couple_config * couple_config,
1308 size_t couple_idx, size_t field_couple_idx) {
1309
1311 couple_config, couple_idx, field_couple_idx,
1312 "yac_couple_config_get_target_lag", __LINE__);
1313
1314 return
1315 couple_config->
1316 couples[couple_idx].
1317 field_couples[field_couple_idx].
1318 target.lag;
1319}
1320
1322 struct yac_couple_config * couple_config,
1323 size_t couple_idx, size_t field_couple_idx) {
1324
1326 couple_config, couple_idx, field_couple_idx,
1327 "yac_couple_config_get_coupling_period", __LINE__);
1328
1329 return
1330 couple_config->
1331 couples[couple_idx].
1332 field_couples[field_couple_idx].
1333 coupling_period;
1334}
1335
1337 struct yac_couple_config * couple_config,
1338 size_t couple_idx, size_t field_couple_idx) {
1339
1341 couple_config, couple_idx, field_couple_idx,
1342 "yac_couple_config_get_source_timestep", __LINE__);
1343
1344 struct yac_couple_config_field_couple * tcouple =
1345 couple_config->couples[couple_idx].field_couples + field_couple_idx;
1346 char const * timestep =
1347 couple_config->components[tcouple->source.component_idx].
1348 fields[tcouple->source.field_idx].timestep;
1349
1351 timestep,
1352 "ERROR(yac_couple_config_get_source_timestep): "
1353 "no valid timestep defined (component: \"%s\" field \"%s\")",
1354 couple_config->components[tcouple->source.component_idx].name,
1355 couple_config->components[tcouple->source.component_idx].
1356 fields[tcouple->source.field_idx].name);
1357
1358 return timestep;
1359}
1360
1362 struct yac_couple_config * couple_config,
1363 size_t couple_idx, size_t field_couple_idx) {
1364
1366 couple_config, couple_idx, field_couple_idx,
1367 "yac_couple_config_get_target_timestep", __LINE__);
1368
1369 struct yac_couple_config_field_couple * tcouple =
1370 couple_config->couples[couple_idx].field_couples + field_couple_idx;
1371 char const * timestep =
1372 couple_config->components[tcouple->target.component_idx].
1373 fields[tcouple->target.field_idx].timestep;
1374
1376 timestep,
1377 "ERROR(yac_couple_config_get_target_timestep): "
1378 "no valid timestep defined (component: \"%s\" field \"%s\")",
1379 couple_config->components[tcouple->target.component_idx].name,
1380 couple_config->components[tcouple->target.component_idx].
1381 fields[tcouple->target.field_idx].name);
1382
1383 return timestep;
1384}
1385
1387 struct yac_couple_config * couple_config,
1388 size_t couple_idx, size_t field_couple_idx) {
1389
1391 couple_config, couple_idx, field_couple_idx,
1392 "yac_couple_config_get_coupling_period_operation", __LINE__);
1393
1394 return
1395 couple_config->
1396 couples[couple_idx].
1397 field_couples[field_couple_idx].
1399}
1400
1402 struct yac_couple_config * couple_config,
1403 char const * start, char const * end) {
1404
1405 if ((start != NULL) && (strlen(start) > 0)) {
1406
1407 // in case start points to couple_config->start_datetime,
1408 // we have to first strdup start, before freeing it
1409 char const * old_start = couple_config->start_datetime;
1410 couple_config->start_datetime = strdup(start);
1411 free((void*)old_start);
1412 }
1413
1414 if ((end != NULL) && (strlen(end) > 0)) {
1415
1416 // in case end points to couple_config->end_datetime,
1417 // we have to first strdup end, before freeing it
1418 char const * old_end = couple_config->end_datetime;
1419 couple_config->end_datetime = strdup(end);
1420 free((void*)old_end);
1421 }
1422}
1423
1425 struct yac_couple_config * couple_config) {
1426
1427 YAC_ASSERT(couple_config->start_datetime,
1428 "ERROR(yac_couple_config_get_start_datetime): "
1429 "start_datetime not yet defined");
1430
1431 return couple_config->start_datetime;
1432}
1433
1435 struct yac_couple_config * couple_config) {
1436
1437 YAC_ASSERT(couple_config->end_datetime,
1438 "ERROR(yac_couple_config_get_start_datetime): "
1439 "start_datetime not yet defined");
1440
1441 return couple_config->end_datetime;
1442}
1443
1445 struct yac_couple_config * couple_config, size_t grid_idx) {
1446
1448 grid_idx < couple_config->num_grids,
1449 "ERROR(yac_couple_config_get_grid_name): "
1450 "Invalid grid idx %zu", grid_idx);
1451
1452 return couple_config->grids[grid_idx].name;
1453}
1454
1456 struct yac_couple_config * couple_config,
1457 size_t couple_idx, size_t field_couple_idx) {
1458
1460 couple_config, couple_idx, field_couple_idx,
1461 "yac_couple_config_enforce_write_weight_file", __LINE__);
1462
1463 return
1464 couple_config->
1465 couples[couple_idx].
1466 field_couples[field_couple_idx].
1468}
1469
1471 struct yac_couple_config * couple_config,
1472 size_t couple_idx, size_t field_couple_idx) {
1473
1475 couple_config, couple_idx, field_couple_idx,
1476 "yac_couple_config_get_weight_file_name", __LINE__);
1477
1478 static char dummy[] = "\0";
1479 char const * weight_file_name =
1480 couple_config->
1481 couples[couple_idx].
1482 field_couples[field_couple_idx].
1484
1485 return (weight_file_name != NULL)?weight_file_name:dummy;
1486}
1487
1489 struct yac_couple_config * couple_config,
1490 size_t couple_idx, size_t field_couple_idx) {
1491
1493 couple_config, couple_idx, field_couple_idx,
1494 "yac_couple_config_get_scale_factor", __LINE__);
1495
1496 return
1497 couple_config->
1498 couples[couple_idx].
1499 field_couples[field_couple_idx].
1501}
1502
1504 struct yac_couple_config * couple_config,
1505 size_t couple_idx, size_t field_couple_idx) {
1506
1508 couple_config, couple_idx, field_couple_idx,
1509 "yac_couple_config_get_scale_summand", __LINE__);
1510
1511 return
1512 couple_config->
1513 couples[couple_idx].
1514 field_couples[field_couple_idx].
1516}
1517
1519 struct yac_couple_config * couple_config,
1520 size_t couple_idx, size_t field_couple_idx,
1521 char const * const ** mask_names, size_t * num_mask_names) {
1522
1524 couple_config, couple_idx, field_couple_idx,
1525 "yac_couple_config_get_src_mask_names", __LINE__);
1526
1527 *mask_names =
1528 (char const * const *)(
1529 couple_config->
1530 couples[couple_idx].
1531 field_couples[field_couple_idx].
1533 *num_mask_names =
1534 couple_config->
1535 couples[couple_idx].
1536 field_couples[field_couple_idx].
1538}
1539
1541 struct yac_couple_config * couple_config,
1542 size_t couple_idx, size_t field_couple_idx) {
1543
1545 couple_config, couple_idx, field_couple_idx,
1546 "yac_couple_config_get_tgt_mask_name", __LINE__);
1547
1548 return
1549 couple_config->
1550 couples[couple_idx].
1551 field_couples[field_couple_idx].
1553}
1554
1556 struct yac_couple_config * couple_config, char const * grid_name) {
1557
1558 for (size_t grid_idx = 0; grid_idx < couple_config->num_grids;
1559 ++grid_idx)
1560 if (!strcmp(couple_config->grids[grid_idx].name, grid_name))
1561 return 1;
1562 return 0;
1563}
1564
1566 char const * string, MPI_Comm comm) {
1567
1568 int strlen_pack_size, string_pack_size;
1569 yac_mpi_call(MPI_Pack_size(1, MPI_INT, comm, &strlen_pack_size), comm);
1570
1571 if (string != NULL) {
1573 MPI_Pack_size(
1574 (int)(strlen(string)), MPI_CHAR, comm, &string_pack_size), comm);
1575 } else {
1576 string_pack_size = 0;
1577 }
1578
1579 return (size_t)strlen_pack_size + (size_t)string_pack_size;
1580}
1581
1583 void * grid_, MPI_Comm comm) {
1584
1585 struct yac_couple_config_grid * grid = grid_;
1586
1587 return yac_couple_config_get_string_pack_size(grid->name, comm);
1588}
1589
1591 size_t num_grids, void * grids_, MPI_Comm comm) {
1592
1593 struct yac_couple_config_grid * grids = grids_;
1594
1595 int num_grids_pack_size;
1596 yac_mpi_call(MPI_Pack_size(1, MPI_INT, comm, &num_grids_pack_size), comm);
1597
1598 size_t grids_pack_size = 0;
1599 for (size_t i = 0; i < num_grids; ++i)
1600 grids_pack_size +=
1602
1603 return (size_t)num_grids_pack_size + grids_pack_size;
1604}
1605
1607 struct yac_couple_config_field * field,
1608 MPI_Comm comm) {
1609
1610 YAC_ASSERT(
1611 field->grid_idx <= INT_MAX,
1612 "ERROR(yac_couple_config_get_field_pack_size):"
1613 "grid_idx is too big")
1614
1615 int int_pack_size;
1616 yac_mpi_call(MPI_Pack_size(1, MPI_INT, comm, &int_pack_size), comm);
1617 int double_pack_size;
1618 yac_mpi_call(MPI_Pack_size(1, MPI_DOUBLE, comm, &double_pack_size), comm);
1619 size_t name_pack_size = yac_couple_config_get_string_pack_size(
1620 field->name, comm);
1621 size_t timestep_pack_size = yac_couple_config_get_string_pack_size(
1622 field->timestep, comm);
1623 size_t metadata_pack_size = yac_couple_config_get_string_pack_size(
1624 field->metadata, comm);
1625
1626 return int_pack_size + // grid_idx
1627 int_pack_size + // collection_size
1628 double_pack_size + // frac_mask_fallback_value
1629 name_pack_size +
1630 timestep_pack_size +
1631 metadata_pack_size;
1632}
1633
1635 size_t num_fields, void * fields_, MPI_Comm comm) {
1636
1637 struct yac_couple_config_field * fields = fields_;
1638
1639 int num_fields_pack_size;
1641 MPI_Pack_size(1, MPI_INT, comm, &num_fields_pack_size), comm);
1642
1643 size_t fields_pack_size = 0;
1644 for (size_t i = 0; i < num_fields; ++i)
1645 fields_pack_size +=
1647 fields + i, comm);
1648
1649 return (size_t)num_fields_pack_size +
1650 fields_pack_size;
1651}
1652
1654 struct yac_couple_config_component * component, MPI_Comm comm) {
1655
1656 return
1658}
1659
1661 size_t num_components, void * components_, MPI_Comm comm) {
1662
1663 struct yac_couple_config_component * components = components_;
1664
1665 int num_components_pack_size;
1667 MPI_Pack_size(1, MPI_INT, comm, &num_components_pack_size), comm);
1668
1669 size_t components_pack_size = 0;
1670 for (size_t i = 0; i < num_components; ++i)
1671 components_pack_size +=
1673
1674 return (size_t)num_components_pack_size + components_pack_size;
1675}
1676
1678 struct yac_couple_config_field_couple * field_couple, MPI_Comm comm) {
1679
1680 int ints_pack_size;
1681 yac_mpi_call(MPI_Pack_size(10, MPI_INT, comm, &ints_pack_size), comm);
1682 int doubles_pack_size;
1683 yac_mpi_call(MPI_Pack_size(2, MPI_DOUBLE, comm, &doubles_pack_size), comm);
1684 int src_mask_names_pack_size = 0;
1685 if (field_couple->num_src_mask_names > 0)
1686 for (size_t i = 0; i < field_couple->num_src_mask_names; ++i)
1687 src_mask_names_pack_size +=
1689 field_couple->src_mask_names[i], comm);
1690
1691 return
1692 (size_t)ints_pack_size + // source.component_idx
1693 // source.field_idx
1694 // source.lag
1695 // target.component_idx
1696 // target.field_idx
1697 // target.lag
1698 // mapping_on_source
1699 // coupling_period_operation
1700 // enforce_write_weight_file
1701 // num_src_mask_names
1703 field_couple->interp_stack, comm) +
1705 field_couple->coupling_period, comm) +
1707 field_couple->weight_file_name, comm) +
1708 doubles_pack_size + // scale_factor
1709 // scale_summand
1710 src_mask_names_pack_size +
1712 field_couple->tgt_mask_name, comm);
1713}
1714
1716 size_t num_field_couples, void * field_couples_, MPI_Comm comm) {
1717
1718 struct yac_couple_config_field_couple * field_couples = field_couples_;
1719
1720 int num_field_couples_pack_size;
1722 MPI_Pack_size(1, MPI_INT, comm, &num_field_couples_pack_size), comm);
1723
1724 size_t field_couples_pack_size = 0;
1725 for (size_t i = 0; i < num_field_couples; ++i)
1726 field_couples_pack_size +=
1727 yac_couple_config_get_field_couple_pack_size(field_couples + i, comm);
1728
1729 return (size_t)num_field_couples_pack_size +
1730 field_couples_pack_size;
1731}
1732
1734 struct yac_couple_config_couple * couple, MPI_Comm comm) {
1735
1736 int component_indices_pack_size;
1738 MPI_Pack_size(2, MPI_INT, comm, &component_indices_pack_size), comm);
1739
1740 return (size_t)component_indices_pack_size;
1741}
1742
1744 size_t num_couples, void * couples_, MPI_Comm comm) {
1745
1746 struct yac_couple_config_couple * couples = couples_;
1747
1748 int num_couples_pack_size;
1750 MPI_Pack_size(1, MPI_INT, comm, &num_couples_pack_size), comm);
1751
1752 size_t couples_pack_size = 0;
1753 for (size_t i = 0; i < num_couples; ++i)
1754 couples_pack_size +=
1756
1757 return (size_t)num_couples_pack_size + couples_pack_size;
1758}
1759
1761 char const * string, void * buffer, int buffer_size, int * position,
1762 MPI_Comm comm) {
1763
1764 size_t len = (string == NULL)?0:strlen(string);
1765
1766 YAC_ASSERT(
1767 len <= INT_MAX, "ERROR(yac_couple_config_pack_string): string too long")
1768
1769 int len_int = (int)len;
1770
1772 MPI_Pack(
1773 &len_int, 1, MPI_INT, buffer, buffer_size, position, comm), comm);
1774
1775 if (len > 0)
1777 MPI_Pack(
1778 string, len_int, MPI_CHAR, buffer, buffer_size, position, comm),
1779 comm);
1780}
1781
1783 struct yac_couple_config_grid * grid,
1784 void * buffer, int buffer_size, int * position, MPI_Comm comm) {
1785
1787 grid->name, buffer, buffer_size, position, comm);
1788}
1789
1791 size_t num_grids, void * grids_,
1792 void * buffer, int buffer_size, int * position, MPI_Comm comm) {
1793
1794 struct yac_couple_config_grid * grids = grids_;
1795
1796 YAC_ASSERT(
1797 num_grids <= INT_MAX,
1798 "ERROR(yac_couple_config_pack_grids): num_grids bigger than INT_MAX")
1799
1800 int num_grids_int = (int)num_grids;
1802 MPI_Pack(
1803 &num_grids_int, 1, MPI_INT, buffer, buffer_size, position, comm), comm);
1804
1805 for (size_t i = 0; i < num_grids; ++i)
1807 grids + i, buffer, buffer_size, position, comm);
1808}
1809
1811 struct yac_couple_config_field * field,
1812 void * buffer, int buffer_size, int * position, MPI_Comm comm) {
1813
1814 YAC_ASSERT(
1815 field->grid_idx <= INT_MAX,
1816 "ERROR(yac_couple_config_pack_field): grid_idx is too big")
1817 YAC_ASSERT(
1818 (field->collection_size < INT_MAX) ||
1819 (field->collection_size == SIZE_MAX),
1820 "ERROR(yac_couple_config_pack_field): invalid collection size")
1821
1822 int grid_idx = (int)field->grid_idx;
1823 double frac_mask_fallback_value = field->frac_mask_fallback_value;
1824 int collection_size =
1825 (field->collection_size == SIZE_MAX)?
1826 INT_MAX:(int)(field->collection_size);
1827
1829 MPI_Pack(
1830 &grid_idx, 1, MPI_INT, buffer, buffer_size, position, comm), comm);
1832 MPI_Pack(
1833 &frac_mask_fallback_value, 1, MPI_DOUBLE, buffer, buffer_size,
1834 position, comm), comm);
1836 MPI_Pack(
1837 &collection_size, 1, MPI_INT, buffer, buffer_size, position, comm), comm);
1838 yac_couple_config_pack_string(field->name, buffer,
1839 buffer_size, position, comm);
1841 buffer_size, position, comm);
1843 buffer_size, position, comm);
1844}
1845
1847 size_t num_fields, void * fields_, void * buffer, int buffer_size,
1848 int * position, MPI_Comm comm) {
1849
1850 struct yac_couple_config_field * fields = fields_;
1851
1852 YAC_ASSERT(
1853 num_fields <= INT_MAX,
1854 "ERROR(yac_couple_config_pack_fields): "
1855 "num_fields bigger than INT_MAX")
1856
1857 int num_fields_int = (int)num_fields;
1859 MPI_Pack(
1860 &num_fields_int, 1, MPI_INT,
1861 buffer, buffer_size, position, comm), comm);
1862
1863 for (size_t i = 0; i < num_fields; ++i)
1865 fields + i, buffer, buffer_size, position, comm);
1866}
1867
1869 struct yac_couple_config_component * component,
1870 void * buffer, int buffer_size, int * position, MPI_Comm comm) {
1871
1873 component->name, buffer, buffer_size, position, comm);
1874}
1875
1877 size_t num_components, void * components_,
1878 void * buffer, int buffer_size, int * position, MPI_Comm comm) {
1879
1880 struct yac_couple_config_component * components = components_;
1881
1882 YAC_ASSERT(
1883 num_components <= INT_MAX,
1884 "ERROR(yac_couple_config_pack_components): "
1885 "num_components bigger than INT_MAX")
1886
1887 int num_components_int = (int)num_components;
1889 MPI_Pack(
1890 &num_components_int, 1, MPI_INT, buffer, buffer_size, position, comm), comm);
1891
1892 for (size_t i = 0; i < num_components; ++i)
1894 components + i, buffer, buffer_size, position, comm);
1895}
1896
1898 struct yac_couple_config_field_couple * field_couple,
1899 void * buffer, int buffer_size, int * position, MPI_Comm comm) {
1900
1901 YAC_ASSERT(
1902 field_couple->source.component_idx <= INT_MAX,
1903 "ERROR(yac_couple_config_pack_field_couple): "
1904 "source.component_idx bigger than INT_MAX")
1905 YAC_ASSERT(
1906 field_couple->source.field_idx <= INT_MAX,
1907 "ERROR(yac_couple_config_pack_field_couple): "
1908 "source.field_idx bigger than INT_MAX")
1909 YAC_ASSERT(
1910 field_couple->target.component_idx <= INT_MAX,
1911 "ERROR(yac_couple_config_pack_field_couple): "
1912 "target.component_idx bigger than INT_MAX")
1913 YAC_ASSERT(
1914 field_couple->target.field_idx <= INT_MAX,
1915 "ERROR(yac_couple_config_pack_field_couple): "
1916 "target.field_idx bigger than INT_MAX")
1917 YAC_ASSERT(
1918 field_couple->mapping_on_source <= INT_MAX,
1919 "ERROR(yac_couple_config_pack_field_couple): "
1920 "mapping_on_source bigger than INT_MAX")
1921 YAC_ASSERT(
1922 field_couple->coupling_period_operation <= INT_MAX,
1923 "ERROR(yac_couple_config_pack_field_couple): "
1924 "coupling_period_operation bigger than INT_MAX")
1925 YAC_ASSERT(
1926 field_couple->enforce_write_weight_file <= INT_MAX,
1927 "ERROR(yac_couple_config_pack_field_couple): "
1928 "enforce_write_weight_file bigger than INT_MAX")
1929 YAC_ASSERT(
1930 field_couple->num_src_mask_names <= INT_MAX,
1931 "ERROR(yac_couple_config_pack_field_couple): "
1932 "num_src_mask_names bigger than INT_MAX")
1933
1934 int ints[10] = {
1935 field_couple->source.component_idx,
1936 field_couple->source.field_idx,
1937 field_couple->source.lag,
1938 field_couple->target.component_idx,
1939 field_couple->target.field_idx,
1940 field_couple->target.lag,
1941 field_couple->mapping_on_source,
1942 field_couple->coupling_period_operation,
1943 field_couple->enforce_write_weight_file,
1944 field_couple->num_src_mask_names};
1945
1947 MPI_Pack(ints, 10, MPI_INT, buffer, buffer_size, position, comm), comm);
1948
1950 field_couple->coupling_period, buffer, buffer_size, position, comm);
1952 field_couple->weight_file_name, buffer, buffer_size, position, comm);
1953
1954 double doubles[2] = {
1955 field_couple->scale_factor,
1956 field_couple->scale_summand};
1957
1959 MPI_Pack(doubles, 2, MPI_DOUBLE, buffer, buffer_size, position, comm),
1960 comm);
1961
1963 field_couple->interp_stack, buffer, buffer_size, position, comm);
1964
1965 if (field_couple->num_src_mask_names > 0)
1966 for (size_t i = 0; i < field_couple->num_src_mask_names; ++i)
1968 field_couple->src_mask_names[i], buffer, buffer_size, position, comm);
1969
1971 field_couple->tgt_mask_name, buffer, buffer_size, position, comm);
1972}
1973
1975 size_t num_field_couples, void * field_couples_,
1976 void * buffer, int buffer_size, int * position, MPI_Comm comm) {
1977
1978 struct yac_couple_config_field_couple * field_couples = field_couples_;
1979
1980 YAC_ASSERT(
1981 num_field_couples <= INT_MAX,
1982 "ERROR(yac_couple_config_pack_field_couples): "
1983 "num_field_couples bigger than INT_MAX")
1984
1985 int num_field_couples_int = (int)num_field_couples;
1987 MPI_Pack(
1988 &num_field_couples_int, 1, MPI_INT,
1989 buffer, buffer_size, position, comm), comm);
1990
1991 for (size_t i = 0; i < num_field_couples; ++i)
1993 field_couples + i, buffer, buffer_size, position, comm);
1994}
1995
1997 struct yac_couple_config_couple * couple,
1998 void * buffer, int buffer_size, int * position, MPI_Comm comm) {
1999
2000 YAC_ASSERT(
2001 (couple->component_indices[0] <= INT_MAX) &&
2002 (couple->component_indices[1] <= INT_MAX),
2003 "ERROR(yac_couple_config_pack_couple_basic): "
2004 "component_indices bigger than INT_MAX")
2005
2006 int component_indices[2] = {(int)(couple->component_indices[0]),
2007 (int)(couple->component_indices[1])};
2008
2010 MPI_Pack(
2011 component_indices, 2, MPI_INT, buffer, buffer_size, position, comm),
2012 comm);
2013}
2014
2016 size_t num_couples, void * couples_,
2017 void * buffer, int buffer_size, int * position, MPI_Comm comm) {
2018
2019 struct yac_couple_config_couple * couples = couples_;
2020
2021 YAC_ASSERT(
2022 num_couples <= INT_MAX,
2023 "ERROR(yac_couple_config_pack_couples_basic): "
2024 "num_couples bigger than INT_MAX")
2025
2026 int num_couples_int = (int)num_couples;
2028 MPI_Pack(
2029 &num_couples_int, 1, MPI_INT, buffer, buffer_size, position, comm), comm);
2030
2031 for (size_t i = 0; i < num_couples; ++i)
2033 couples + i, buffer, buffer_size, position, comm);
2034}
2035
2037 void * buffer, int buffer_size, int * position, MPI_Comm comm) {
2038
2039 int string_len;
2041 MPI_Unpack(
2042 buffer, buffer_size, position, &string_len, 1, MPI_INT, comm), comm);
2043
2044 if (string_len <= 0) return NULL;
2045
2046 char * string = xmalloc(((size_t)string_len + 1) * sizeof(*string));
2048 MPI_Unpack(
2049 buffer, buffer_size, position, string, string_len, MPI_CHAR, comm), comm);
2050 string[string_len] = '\0';
2051 return string;
2052}
2053
2055 void * buffer, int buffer_size, int * position,
2056 struct yac_couple_config_grid * grid, MPI_Comm comm) {
2057
2058 grid->name =
2059 yac_couple_config_unpack_string(buffer, buffer_size, position, comm);
2060 grid->metadata = NULL;
2061}
2062
2064 void * buffer, int buffer_size, int * position,
2065 size_t * num_grids, void * grids_, MPI_Comm comm) {
2066
2067 struct yac_couple_config_grid ** grids = grids_;
2068
2069 int num_grids_int;
2071 MPI_Unpack(
2072 buffer, buffer_size, position, &num_grids_int, 1, MPI_INT, comm), comm);
2073
2074 YAC_ASSERT(
2075 num_grids_int >= 0,
2076 "ERROR(yac_couple_config_unpack_grids): invalid number of grids")
2077
2078 *num_grids = (size_t)num_grids_int;
2079 *grids = xmalloc(*num_grids * sizeof(**grids));
2080
2081 for (size_t i = 0; i < *num_grids; ++i)
2083 buffer, buffer_size, position, *grids + i, comm);
2084}
2085
2087 void * buffer, int buffer_size, int * position,
2088 struct yac_couple_config_field * field, MPI_Comm comm) {
2089
2090 int grid_idx;
2091 double frac_mask_fallback_value;
2092 int collection_size;
2093
2095 MPI_Unpack(
2096 buffer, buffer_size, position, &grid_idx, 1, MPI_INT, comm), comm);
2098 MPI_Unpack(
2099 buffer, buffer_size, position, &frac_mask_fallback_value, 1,
2100 MPI_DOUBLE, comm), comm);
2102 MPI_Unpack(
2103 buffer, buffer_size, position, &collection_size, 1, MPI_INT, comm),
2104 comm);
2105
2106 YAC_ASSERT(
2107 grid_idx >= 0,
2108 "ERROR(yac_couple_config_unpack_field): invalid number of grid_idx")
2109 YAC_ASSERT(
2110 collection_size >= 0,
2111 "ERROR(yac_couple_config_unpack_field): invalid collection_size")
2112
2113 field->grid_idx = (size_t)grid_idx;
2114 field->frac_mask_fallback_value = frac_mask_fallback_value;
2115 field->collection_size =
2116 (collection_size == INT_MAX)?SIZE_MAX:((int)collection_size);
2117 field->name = yac_couple_config_unpack_string(buffer,
2118 buffer_size, position, comm);
2120 buffer_size, position, comm);
2122 buffer_size, position, comm);
2123}
2124
2126 void * buffer, int buffer_size, int * position,
2127 size_t * num_fields, void * fields_, MPI_Comm comm) {
2128
2129 struct yac_couple_config_field ** fields = fields_;
2130
2131 int num_fields_int;
2133 MPI_Unpack(
2134 buffer, buffer_size, position, &num_fields_int, 1,
2135 MPI_INT, comm), comm);
2136
2137 YAC_ASSERT(
2138 num_fields_int >= 0,
2139 "ERROR(yac_couple_config_unpack_fields): "
2140 "invalid number of fields")
2141
2142 *num_fields = (size_t)num_fields_int;
2143 *fields =
2144 xmalloc(*num_fields * sizeof(**fields));
2145
2146 for (size_t i = 0; i < *num_fields; ++i)
2148 buffer, buffer_size, position, *fields + i, comm);
2149}
2150
2152 void * buffer, int buffer_size, int * position,
2153 struct yac_couple_config_component * component, MPI_Comm comm) {
2154
2155 component->name =
2156 yac_couple_config_unpack_string(buffer, buffer_size, position, comm);
2157 component->metadata = NULL;
2158 component->num_fields = 0;
2159 component->fields = NULL;
2160}
2161
2163 void * buffer, int buffer_size, int * position,
2164 size_t * num_components, void * components_, MPI_Comm comm) {
2165
2166 struct yac_couple_config_component ** components = components_;
2167
2168 int num_components_int;
2170 MPI_Unpack(
2171 buffer, buffer_size, position, &num_components_int, 1, MPI_INT, comm),
2172 comm);
2173
2174 YAC_ASSERT(
2175 num_components_int >= 0,
2176 "ERROR(yac_couple_config_unpack_components): "
2177 "invalid number of components")
2178
2179 *num_components = (size_t)num_components_int;
2181
2182 for (size_t i = 0; i < *num_components; ++i)
2184 buffer, buffer_size, position, *components + i, comm);
2185}
2186
2188 void * buffer, int buffer_size, int * position,
2189 struct yac_couple_config_field_couple * field_couple,
2190 MPI_Comm comm) {
2191
2192 int ints[10];
2194 MPI_Unpack(
2195 buffer, buffer_size, position, ints, 10, MPI_INT, comm), comm);
2196
2197 YAC_ASSERT(
2198 ints[0] >= 0,
2199 "ERROR(yac_couple_config_unpack_field_couple): "
2200 "invalid source.component_idx")
2201 YAC_ASSERT(
2202 ints[1] >= 0,
2203 "ERROR(yac_couple_config_unpack_field_couple): "
2204 "invalid source.field_idx")
2205 YAC_ASSERT(
2206 ints[3] >= 0,
2207 "ERROR(yac_couple_config_unpack_field_couple): "
2208 "target.component_idx bigger than INT_MAX")
2209 YAC_ASSERT(
2210 ints[4] >= 0,
2211 "ERROR(yac_couple_config_unpack_field_couple): "
2212 "invalid target.field_idx")
2213 YAC_ASSERT(
2214 ints[6] >= 0,
2215 "ERROR(yac_couple_config_unpack_field_couple): "
2216 "invalid mapping_on_source")
2217 YAC_ASSERT(
2218 ints[7] >= 0,
2219 "ERROR(yac_couple_config_unpack_field_couple): "
2220 "invalid coupling_period_operation")
2221 YAC_ASSERT(
2222 ints[8] >= 0,
2223 "ERROR(yac_couple_config_unpack_field_couple): "
2224 "invalid enforce_write_weight_file")
2225 YAC_ASSERT(
2226 ints[9] >= 0,
2227 "ERROR(yac_couple_config_unpack_field_couple): "
2228 "invalid num_src_mask_names")
2229
2230 field_couple->source.component_idx = (size_t)(ints[0]);
2231 field_couple->source.field_idx = (size_t)(ints[1]);
2232 field_couple->source.lag = ints[2];
2233 field_couple->target.component_idx = (size_t)(ints[3]);
2234 field_couple->target.field_idx = (size_t)(ints[4]);
2235 field_couple->target.lag = ints[5];
2236 field_couple->mapping_on_source = ints[6];
2237 field_couple->coupling_period_operation =
2238 (enum yac_reduction_type)(ints[7]);
2239 field_couple->enforce_write_weight_file = ints[8];
2240 field_couple->num_src_mask_names = (size_t)(ints[9]);
2241
2242 field_couple->coupling_period =
2243 yac_couple_config_unpack_string(buffer, buffer_size, position, comm);
2244 field_couple->weight_file_name =
2245 yac_couple_config_unpack_string(buffer, buffer_size, position, comm);
2246
2247 double doubles[2];
2249 MPI_Unpack(
2250 buffer, buffer_size, position, doubles, 2, MPI_DOUBLE, comm), comm);
2251
2252 field_couple->scale_factor = doubles[0];
2253 field_couple->scale_summand = doubles[1];
2254
2255 field_couple->interp_stack =
2256 yac_interp_stack_config_unpack(buffer, buffer_size, position, comm);
2257
2258 if (field_couple->num_src_mask_names > 0) {
2259 field_couple->src_mask_names =
2260 xmalloc(
2261 field_couple->num_src_mask_names *
2262 sizeof(*(field_couple->src_mask_names)));
2263 for (size_t i = 0; i < field_couple->num_src_mask_names; ++i)
2264 field_couple->src_mask_names[i] =
2266 buffer, buffer_size, position, comm);
2267 } else {
2268 field_couple->src_mask_names = NULL;
2269 }
2270
2271 field_couple->tgt_mask_name =
2273 buffer, buffer_size, position, comm);
2274}
2275
2277 void * buffer, int buffer_size, int * position,
2278 size_t * num_field_couples, void * field_couples_, MPI_Comm comm) {
2279
2280 struct yac_couple_config_field_couple ** field_couples = field_couples_;
2281
2282 int num_field_couples_int;
2284 MPI_Unpack(
2285 buffer, buffer_size, position, &num_field_couples_int, 1,
2286 MPI_INT, comm), comm);
2287
2288 YAC_ASSERT(
2289 num_field_couples_int >= 0,
2290 "ERROR(yac_couple_config_unpack_field_couples): "
2291 "invalid number of field_couples")
2292
2293 *num_field_couples = (size_t)num_field_couples_int;
2294 *field_couples =
2295 xmalloc(*num_field_couples * sizeof(**field_couples));
2296
2297 for (size_t i = 0; i < *num_field_couples; ++i)
2299 buffer, buffer_size, position, *field_couples + i, comm);
2300}
2301
2303 void * buffer, int buffer_size, int * position,
2304 struct yac_couple_config_couple * couple, MPI_Comm comm) {
2305
2306 int component_indices[2];
2308 MPI_Unpack(
2309 buffer, buffer_size, position, component_indices, 2, MPI_INT, comm),
2310 comm);
2311
2312 YAC_ASSERT(
2313 (component_indices[0] >= 0) && (component_indices[1] >= 0),
2314 "ERROR(yac_couple_config_unpack_couple_basic): invalid component indices")
2315
2316 couple->component_indices[0] = (size_t)(component_indices[0]);
2317 couple->component_indices[1] = (size_t)(component_indices[1]);
2318 couple->num_field_couples = 0;
2319 couple->field_couples = NULL;
2320}
2321
2323 void * buffer, int buffer_size, int * position,
2324 size_t * num_couples, void * couples_, MPI_Comm comm) {
2325
2326 struct yac_couple_config_couple ** couples = couples_;
2327
2328 int num_couples_int;
2330 MPI_Unpack(
2331 buffer, buffer_size, position, &num_couples_int, 1, MPI_INT, comm),
2332 comm);
2333
2334 YAC_ASSERT(
2335 num_couples_int >= 0,
2336 "ERROR(yac_couple_config_unpack_couples_basic): "
2337 "invalid number of couples")
2338
2339 *num_couples = (size_t)num_couples_int;
2340 *couples = xmalloc(*num_couples * sizeof(**couples));
2341
2342 for (size_t i = 0; i < *num_couples; ++i)
2344 buffer, buffer_size, position, *couples + i, comm);
2345}
2346
2348 struct yac_couple_config * couple_config,
2349 char const * src_comp_name, char const * src_grid_name, char const * src_field_name,
2350 char const * tgt_comp_name, char const * tgt_grid_name, char const * tgt_field_name,
2351 char const * coupling_period, int time_reduction,
2352 struct yac_interp_stack_config * interp_stack,
2353 int src_lag, int tgt_lag,
2354 const char* weight_file_name, int mapping_on_source,
2355 double scale_factor, double scale_summand,
2356 size_t num_src_mask_names, char const * const * src_mask_names,
2357 char const * tgt_mask_name) {
2358
2359 YAC_ASSERT(src_comp_name && src_comp_name[0] != '\0',
2360 "ERROR(yac_couple_config_def_couple): invalid parameter: src_comp_name");
2361 YAC_ASSERT(src_grid_name && src_grid_name[0] != '\0',
2362 "ERROR(yac_couple_config_def_couple): invalid parameter: src_grid_name");
2363 YAC_ASSERT(src_field_name && src_field_name[0] != '\0',
2364 "ERROR(yac_couple_config_def_couple): invalid parameter: src_field_name");
2365 YAC_ASSERT(tgt_comp_name && tgt_comp_name[0] != '\0',
2366 "ERROR(yac_couple_config_def_couple): invalid parameter: tgt_comp_name");
2367 YAC_ASSERT(tgt_grid_name && tgt_grid_name[0] != '\0',
2368 "ERROR(yac_couple_config_def_couple): invalid parameter: tgt_grid_name");
2369 YAC_ASSERT(tgt_field_name && tgt_field_name[0] != '\0',
2370 "ERROR(yac_couple_config_def_couple): invalid parameter: tgt_field_name");
2371 YAC_ASSERT(coupling_period && coupling_period[0] != '\0',
2372 "ERROR(yac_couple_config_def_couple): invalid parameter: coupling_period");
2373 YAC_ASSERT(
2374 (time_reduction == TIME_NONE) ||
2375 (time_reduction == TIME_ACCUMULATE) ||
2376 (time_reduction == TIME_AVERAGE) ||
2377 (time_reduction == TIME_MINIMUM) ||
2378 (time_reduction == TIME_MAXIMUM),
2379 "ERROR(yac_couple_config_def_couple): invalid parameter: time_reduction");
2380 YAC_ASSERT_F(isnormal(scale_factor),
2381 "ERROR(yac_couple_config_def_couple): \"%lf\" is not a valid scale factor",
2382 scale_factor);
2383 YAC_ASSERT_F(isnormal(scale_summand) || (scale_summand == 0.0),
2384 "ERROR(yac_couple_config_def_couple): \"%lf\" is not a valid scale summand",
2385 scale_summand);
2386
2387 // get component indices
2388 size_t src_comp_idx =
2389 yac_couple_config_add_component_(couple_config, src_comp_name);
2390 size_t tgt_comp_idx =
2391 yac_couple_config_add_component_(couple_config, tgt_comp_name);
2392 size_t src_grid_idx =
2393 yac_couple_config_add_grid_(couple_config, src_grid_name);
2394 size_t tgt_grid_idx =
2395 yac_couple_config_add_grid_(couple_config, tgt_grid_name);
2396
2397 // check if couple exists
2398 size_t component_indices[2];
2399 if(src_comp_idx < tgt_comp_idx){
2400 component_indices[0] = src_comp_idx;
2401 component_indices[1] = tgt_comp_idx;
2402 }else{
2403 component_indices[0] = tgt_comp_idx;
2404 component_indices[1] = src_comp_idx;
2405 }
2406 struct yac_couple_config_couple * couple = NULL;
2407 for(size_t i = 0; (i < couple_config->num_couples) && !couple; ++i)
2408 if(couple_config->couples[i].component_indices[0] == component_indices[0] &&
2409 couple_config->couples[i].component_indices[1] == component_indices[1])
2410 couple = &couple_config->couples[i];
2411
2412 // create if couple does not exists
2413 if(!couple){
2414 couple_config->couples =
2415 xrealloc(
2416 couple_config->couples,
2417 (couple_config->num_couples + 1) * sizeof(*couple_config->couples));
2418 couple = &couple_config->couples[couple_config->num_couples];
2419 couple_config->num_couples++;
2420 couple->component_indices[0] = component_indices[0];
2421 couple->component_indices[1] = component_indices[1];
2422 couple->num_field_couples = 0;
2423 couple->field_couples = NULL;
2424 }
2425
2426 // get field indices
2427 size_t src_field_idx =
2429 couple_config, src_comp_idx, src_grid_idx,
2430 src_field_name, NULL, SIZE_MAX);
2431 size_t tgt_field_idx =
2433 couple_config, tgt_comp_idx, tgt_grid_idx,
2434 tgt_field_name, NULL, SIZE_MAX);
2435
2436 // check if field_couple exists
2437 for(size_t i = 0; i < couple->num_field_couples; ++i)
2439 !(couple->field_couples[i].source.component_idx == src_comp_idx &&
2440 couple->field_couples[i].target.component_idx == tgt_comp_idx &&
2441 couple->field_couples[i].source.field_idx == src_field_idx &&
2442 couple->field_couples[i].target.field_idx == tgt_field_idx),
2443 "ERROR(yac_couple_config_def_couple): Coupling is already defined \n"
2444 "source (comp_name: \"%s\" grid_name: \"%s\" field_name: \"%s\")\n"
2445 "target (comp_name: \"%s\" grid_name: \"%s\" field_name: \"%s\")\n",
2446 src_comp_name, src_grid_name, src_field_name,
2447 tgt_comp_name, tgt_grid_name, tgt_field_name);
2448
2449 couple->field_couples = xrealloc(couple->field_couples,
2450 (couple->num_field_couples + 1) * sizeof(*couple->field_couples));
2451 struct yac_couple_config_field_couple * field_couple =
2452 couple->field_couples + couple->num_field_couples;
2453 couple->num_field_couples++;
2454 field_couple->source.component_idx = src_comp_idx;
2455 field_couple->source.field_idx = src_field_idx;
2456 field_couple->source.lag = src_lag;
2457 field_couple->target.component_idx = tgt_comp_idx;
2458 field_couple->target.field_idx = tgt_field_idx;
2459 field_couple->target.lag = tgt_lag;
2460 field_couple->coupling_period = strdup(coupling_period);
2461 field_couple->coupling_period_operation =
2462 (enum yac_reduction_type)time_reduction;
2464 field_couple->mapping_on_source = mapping_on_source;
2465 field_couple->weight_file_name =
2467 field_couple->scale_factor = scale_factor;
2468 field_couple->scale_summand = scale_summand;
2469 field_couple->enforce_write_weight_file = weight_file_name != NULL;
2470 field_couple->num_src_mask_names = num_src_mask_names;
2471 if (num_src_mask_names > 0) {
2472 field_couple->src_mask_names =
2474 for (size_t i = 0; i < num_src_mask_names; ++i)
2475 field_couple->src_mask_names[i] = strdup(src_mask_names[i]);
2476 } else {
2477 field_couple->src_mask_names = NULL;
2478 }
2479 field_couple->tgt_mask_name =
2480 (tgt_mask_name != NULL)?strdup(tgt_mask_name):NULL;
2481}
2482
2484 struct yac_couple_config * couple_config, MPI_Comm comm) {
2485
2488 "start_datetime", (char**)&(couple_config->start_datetime), comm);
2490 "end_datetime", (char**)&(couple_config->end_datetime), comm);
2491}
2492
2494 size_t(*get_pack_size)(size_t, void*, MPI_Comm);
2495 void(*pack)(size_t, void*, void *, int, int *, MPI_Comm);
2496 void(*unpack)(void *, int, int *, size_t *, void *, MPI_Comm);
2497 int(*compare)(void const*, void const*);
2498 void(*merge)(void*, void *, MPI_Comm);
2499 void(*free_data)(void*);
2535
2536static void dist_merge(size_t* len, void** arr, size_t element_size, MPI_Comm comm,
2537 struct dist_merge_vtable* vtable, size_t** idx_old_to_new){
2538 int rank;
2539 MPI_Comm_rank(comm, &rank);
2540
2541 // initialise
2542 unsigned char * input = *arr;
2543 size_t input_len = *len;
2544 *idx_old_to_new = xmalloc(input_len * sizeof(size_t));
2545 size_t* idx = xmalloc(input_len * sizeof(size_t));
2546 for(size_t i = 0;i<input_len;++i) idx[i] = i;
2547 unsigned char * arr_new = NULL;
2548 size_t len_new = 0;
2549
2550 // sort input data
2551 yac_qsort_index(input, input_len, element_size, vtable->compare, idx);
2552
2553 void * buffer = NULL;
2554 while(1){
2555
2556 // determine pack size of local data
2557 size_t pack_size =
2558 (input_len > 0)?vtable->get_pack_size(input_len, input, comm):0;
2559 YAC_ASSERT(
2560 pack_size <= LONG_MAX, "ERROR(dist_merge): packing size too big");
2561
2562 // determine rank with most amount of data
2563 struct {
2564 long pack_size;
2565 int rank;
2566 } data_pair = {.pack_size = (long)pack_size, .rank = rank};
2568 MPI_Allreduce(
2569 MPI_IN_PLACE, &data_pair, 1, MPI_LONG_INT, MPI_MAXLOC, comm), comm);
2570
2571 // if there is no more data to exchange
2572 if (data_pair.pack_size == 0) break;
2573
2574 // pack data into buffer buffer
2575 pack_size = (size_t)data_pair.pack_size;
2576 if (!buffer) buffer = xmalloc(pack_size);
2577 int position = 0;
2578 if(data_pair.rank == rank)
2579 vtable->pack(input_len, input, buffer, pack_size, &position, comm);
2580
2581 // broadcast and unpack data
2583 MPI_Bcast(buffer, pack_size, MPI_PACKED, data_pair.rank, comm), comm);
2584 void * recved = NULL;
2585 size_t num_recved;
2586 position = 0;
2587 vtable->unpack(buffer, pack_size, &position, &num_recved, &recved, comm);
2588
2589 // add received elements to list
2590 arr_new = xrealloc(arr_new, (len_new + num_recved)*element_size);
2591 memcpy(
2592 arr_new + len_new*element_size, recved, num_recved*element_size);
2593 free(recved);
2594
2595 // for all received elements
2596 size_t input_idx, input_len_new, i = len_new;
2597 len_new += num_recved;
2598 for(input_idx = 0, input_len_new = 0; i < len_new; ++i) {
2599 void* recved_element = arr_new + i*element_size;
2600
2601 // search for matching element in input list
2602 int cmp = 0;
2603 void * input_element = NULL;
2604 while ((input_idx < input_len) &&
2605 (((cmp =
2606 vtable->compare(
2607 ((input_element = input + input_idx * element_size)),
2608 recved_element))) < 0)) {
2609
2610 if (input_idx != input_len_new) {
2611 memcpy(
2612 input + input_len_new * element_size,
2613 input_element, element_size);
2614 idx[input_len_new] = idx[input_idx];
2615 }
2616 ++input_len_new;
2617 ++input_idx;
2618 }
2619
2620 // if the end of the local input list was reached
2621 if (input_idx == input_len) break;
2622
2623 // if a matching element was found in the input list
2624 if (!cmp) {
2625
2626 // merge input list element with received element
2627 if (vtable->merge)
2628 vtable->merge(recved_element, input_element, comm);
2629 // remove element
2630 vtable->free_data(input_element);
2631 // update index
2632 (*idx_old_to_new)[idx[input_idx]] = i;
2633 // upate input idx
2634 input_idx++;
2635
2636 // if no matching element was found in the input list
2637 } else if (vtable->merge) {
2638 // since the merge operation can potentially be collective, we have
2639 // to call it, even if no matching element was found
2640 vtable->merge(recved_element, NULL, comm);
2641 }
2642 }
2643 // process remaining received elements
2644 if (vtable->merge)
2645 for(; i < len_new; ++i)
2646 vtable->merge(arr_new + i*element_size, NULL, comm);
2647 // pack remaining elements in input list
2648 for (; input_idx < input_len; ++input_idx, ++input_len_new) {
2649 if (input_idx != input_len_new) {
2650 memcpy(
2651 input + input_len_new * element_size,
2652 input + input_idx * element_size, element_size);
2653 idx[input_len_new] = idx[input_idx];
2654 }
2655 }
2656 // update length of input list
2657 input_len = input_len_new;
2658 }
2659 free(buffer);
2660 free(input);
2661 free(idx);
2662 *arr = arr_new;
2663 *len = len_new;
2664}
2665
2666static void merge_grids(
2667 struct yac_couple_config * couple_config, MPI_Comm comm) {
2668
2669 size_t* old_to_new_idx;
2670 void * p_grids = couple_config->grids;
2671 dist_merge(
2672 &couple_config->num_grids, &p_grids,
2673 sizeof(couple_config->grids[0]),
2674 comm, &dist_merge_vtable_grid, &old_to_new_idx);
2675 couple_config->grids = p_grids;
2676
2677 // set new grid_idx in fields
2678 for(size_t comp_idx = 0; comp_idx < couple_config->num_components;
2679 ++comp_idx) {
2680 struct yac_couple_config_component * component =
2681 couple_config->components + comp_idx;
2682 for(size_t field_idx = 0; field_idx < component->num_fields; ++field_idx)
2683 component->fields[field_idx].grid_idx =
2684 old_to_new_idx[component->fields[field_idx].grid_idx];
2685 }
2686
2687 free(old_to_new_idx);
2688}
2689
2690static void merge_fields(
2691 struct yac_couple_config * couple_config, size_t comp_idx, MPI_Comm comm) {
2692
2693 struct yac_couple_config_component * component =
2694 couple_config->components + comp_idx;
2695
2696 size_t* old_to_new_idx;
2697 void * p_fields = component->fields;
2698 dist_merge(
2699 &component->num_fields, &p_fields,
2700 sizeof(component->fields[0]),
2701 comm, &dist_merge_vtable_field, &old_to_new_idx);
2702 component->fields = p_fields;
2703
2704 // set new field_idx in all field_couples
2705 for(size_t couple_idx = 0; couple_idx < couple_config->num_couples;
2706 ++couple_idx) {
2707 struct yac_couple_config_couple * couple =
2708 couple_config->couples + couple_idx;
2709 for(size_t field_couple_idx = 0;
2710 field_couple_idx < couple->num_field_couples; ++field_couple_idx){
2711 struct yac_couple_config_field_couple * field_couple =
2712 couple->field_couples + field_couple_idx;
2713 if(field_couple->source.component_idx == comp_idx)
2714 field_couple->source.field_idx =
2715 old_to_new_idx[field_couple->source.field_idx];
2716 if(field_couple->target.component_idx == comp_idx)
2717 field_couple->target.field_idx =
2718 old_to_new_idx[field_couple->target.field_idx];
2719 }
2720 }
2721
2722 free(old_to_new_idx);
2723}
2724
2726 struct yac_couple_config * couple_config, MPI_Comm comm) {
2727
2728 // distribute and merge basic component information while keeping the
2729 // individual fields
2730 size_t* old_to_new_idx;
2731 void * p_components = couple_config->components;
2732 dist_merge(
2733 &couple_config->num_components, &p_components,
2734 sizeof(couple_config->components[0]),
2735 comm, &dist_merge_vtable_component, &old_to_new_idx);
2736 couple_config->components = p_components;
2737
2738 // set new component_idx in couples
2739 for(size_t couple_idx = 0; couple_idx < couple_config->num_couples;
2740 ++couple_idx) {
2741 struct yac_couple_config_couple * couple =
2742 couple_config->couples + couple_idx;
2743 couple->component_indices[0] = old_to_new_idx[couple->component_indices[0]];
2744 couple->component_indices[1] = old_to_new_idx[couple->component_indices[1]];
2745 for(size_t field_couple_idx = 0;
2746 field_couple_idx < couple->num_field_couples; ++field_couple_idx) {
2747 couple->field_couples[field_couple_idx].source.component_idx =
2748 old_to_new_idx[
2749 couple->field_couples[field_couple_idx].source.component_idx];
2750 couple->field_couples[field_couple_idx].target.component_idx =
2751 old_to_new_idx[
2752 couple->field_couples[field_couple_idx].target.component_idx];
2753 }
2754 }
2755 free(old_to_new_idx);
2756
2757 // merge the fields of each component
2758 for (size_t comp_idx = 0; comp_idx < couple_config->num_components; ++comp_idx)
2759 merge_fields(couple_config, comp_idx, comm);
2760}
2761
2763 size_t * num_field_couples,
2764 struct yac_couple_config_field_couple ** field_couples, MPI_Comm comm) {
2765
2766 size_t* old_to_new_idx;
2767 void * p_field_couples = *field_couples;
2768 dist_merge(
2769 num_field_couples, &p_field_couples, sizeof(**field_couples),
2770 comm, &dist_merge_vtable_field_couple, &old_to_new_idx);
2771 free(old_to_new_idx);
2772 *field_couples = p_field_couples;
2773}
2774
2775static void merge_couples(
2776 struct yac_couple_config * couple_config, MPI_Comm comm) {
2777
2778 // distribute and merge basic couple information while keeping the
2779 // individual field couples
2780 size_t* old_to_new_idx;
2781 void * p_couples = couple_config->couples;
2782 dist_merge(
2783 &couple_config->num_couples, &p_couples,
2784 sizeof(couple_config->couples[0]),
2785 comm, &dist_merge_vtable_couple, &old_to_new_idx);
2786 couple_config->couples = p_couples;
2787 free(old_to_new_idx);
2788}
2789
2791 struct yac_couple_config * couple_config, MPI_Comm comm){
2792
2793 // sync time stuff
2794 couple_config_sync_time(couple_config, comm);
2795
2796 merge_grids(couple_config, comm);
2797 merge_components(couple_config, comm);
2798 merge_couples(couple_config, comm);
2799}
void yac_couple_config_grid_set_metadata(struct yac_couple_config *couple_config, char const *grid_name, const char *metadata)
static size_t yac_couple_config_get_field_pack_size(struct yac_couple_config_field *field, MPI_Comm comm)
static void check_field_couple_idx(struct yac_couple_config *couple_config, size_t couple_idx, size_t field_couple_idx, char const *routine_name, int line)
static void yac_couple_config_unpack_field_couples(void *buffer, int buffer_size, int *position, size_t *num_field_couples, void *field_couples_, MPI_Comm comm)
static void yac_couple_config_grid_merge(void *a_, void *b_, MPI_Comm comm)
int yac_couple_config_field_is_valid(struct yac_couple_config *couple_config, size_t component_idx, size_t field_idx)
static void dist_merge(size_t *len, void **arr, size_t element_size, MPI_Comm comm, struct dist_merge_vtable *vtable, size_t **idx_old_to_new)
static void couple_config_sync_time(struct yac_couple_config *couple_config, MPI_Comm comm)
static void yac_couple_config_field_couple_free(void *field_couple_)
static int yac_couple_config_grid_compare(void const *a, void const *b)
static size_t yac_couple_config_get_grid_pack_size(void *grid_, MPI_Comm comm)
static void yac_couple_config_component_merge(void *a_, void *b_, MPI_Comm comm)
void yac_couple_config_field_enable_frac_mask(struct yac_couple_config *couple_config, char const *comp_name, char const *grid_name, char const *field_name, double frac_mask_fallback_value)
static void yac_couple_config_unpack_component(void *buffer, int buffer_size, int *position, struct yac_couple_config_component *component, MPI_Comm comm)
static void couple_config_sync_string(char const *string_name, char **string, MPI_Comm comm)
static void merge_fields(struct yac_couple_config *couple_config, size_t comp_idx, MPI_Comm comm)
static char const * string_dup(char const *string)
static void yac_couple_config_couple_free(void *couple_)
static void yac_couple_config_pack_components(size_t num_components, void *components_, void *buffer, int buffer_size, int *position, MPI_Comm comm)
size_t yac_couple_config_get_num_fields(struct yac_couple_config *couple_config, size_t component_idx)
void yac_couple_config_get_field_grid_names(struct yac_couple_config *couple_config, size_t couple_idx, size_t field_couple_idx, char const **src_grid_name, char const **tgt_grid_name)
int yac_couple_config_contains_grid_name(struct yac_couple_config *couple_config, char const *grid_name)
static void yac_couple_config_pack_field_couples(size_t num_field_couples, void *field_couples_, void *buffer, int buffer_size, int *position, MPI_Comm comm)
void yac_couple_config_sync(struct yac_couple_config *couple_config, MPI_Comm comm)
void yac_couple_config_set_datetime(struct yac_couple_config *couple_config, char const *start, char const *end)
static void yac_couple_config_unpack_couple_basic(void *buffer, int buffer_size, int *position, struct yac_couple_config_couple *couple, MPI_Comm comm)
void yac_couple_config_component_add_field(struct yac_couple_config *couple_config, const char *component_name, const char *grid_name, const char *name, char const *timestep, size_t collection_size)
char const * yac_couple_config_get_field_timestep(struct yac_couple_config *couple_config, char const *component_name, char const *grid_name, char const *field_name)
static void yac_couple_config_pack_couples_basic(size_t num_couples, void *couples_, void *buffer, int buffer_size, int *position, MPI_Comm comm)
int yac_couple_config_mapping_on_source(struct yac_couple_config *couple_config, size_t couple_idx, size_t field_couple_idx)
int yac_couple_config_get_field_role(struct yac_couple_config *couple_config, char const *component_name, char const *grid_name, char const *field_name)
int yac_couple_config_component_name_is_valid(struct yac_couple_config *couple_config, char const *component_name)
static void yac_couple_config_pack_couple_basic(struct yac_couple_config_couple *couple, void *buffer, int buffer_size, int *position, MPI_Comm comm)
static size_t yac_couple_config_get_field_couples_pack_size(size_t num_field_couples, void *field_couples_, MPI_Comm comm)
static char * yac_couple_config_unpack_string(void *buffer, int buffer_size, int *position, MPI_Comm comm)
char const * yac_couple_config_get_coupling_period(struct yac_couple_config *couple_config, size_t couple_idx, size_t field_couple_idx)
struct dist_merge_vtable dist_merge_vtable_component
void yac_couple_config_field_set_metadata(struct yac_couple_config *couple_config, const char *comp_name, const char *grid_name, const char *field_name, const char *metadata)
struct yac_interp_stack_config * yac_couple_config_get_interp_stack(struct yac_couple_config *couple_config, size_t couple_idx, size_t field_couple_idx)
static void yac_couple_config_unpack_grid(void *buffer, int buffer_size, int *position, struct yac_couple_config_grid *grid, MPI_Comm comm)
static size_t yac_couple_config_get_field_couple_pack_size(struct yac_couple_config_field_couple *field_couple, MPI_Comm comm)
char const * yac_couple_config_get_field_name(struct yac_couple_config *couple_config, size_t component_idx, size_t field_idx)
size_t yac_couple_config_get_field_idx(struct yac_couple_config *couple_config, size_t component_idx, size_t grid_idx, char const *field_name)
static size_t yac_couple_config_get_components_pack_size(size_t num_components, void *components_, MPI_Comm comm)
static void couple_config_sync_calendar(MPI_Comm comm)
static size_t yac_couple_config_add_component_(struct yac_couple_config *couple_config, char const *name)
void yac_couple_config_component_set_metadata(struct yac_couple_config *couple_config, char const *comp_name, const char *metadata)
static size_t yac_couple_config_add_grid_(struct yac_couple_config *couple_config, char const *name)
struct dist_merge_vtable dist_merge_vtable_field
static void merge_components(struct yac_couple_config *couple_config, MPI_Comm comm)
static void check_couple_idx(struct yac_couple_config *couple_config, size_t couple_idx, char const *routine_name, int line)
size_t yac_couple_config_get_num_couple_fields(struct yac_couple_config *couple_config, size_t couple_idx)
void yac_couple_config_get_src_mask_names(struct yac_couple_config *couple_config, size_t couple_idx, size_t field_couple_idx, char const *const **mask_names, size_t *num_mask_names)
static int yac_couple_config_component_compare(void const *a_, void const *b_)
static void yac_couple_config_unpack_field(void *buffer, int buffer_size, int *position, struct yac_couple_config_field *field, MPI_Comm comm)
static void yac_couple_config_unpack_couples_basic(void *buffer, int buffer_size, int *position, size_t *num_couples, void *couples_, MPI_Comm comm)
size_t yac_couple_config_get_num_couples(struct yac_couple_config *couple_config)
static size_t yac_couple_config_get_couple_pack_size_basic(struct yac_couple_config_couple *couple, MPI_Comm comm)
static int yac_couple_config_field_compare(void const *a_, void const *b_)
static void yac_couple_config_couple_merge(void *a_, void *b_, MPI_Comm comm)
int yac_couple_config_get_target_lag(struct yac_couple_config *couple_config, size_t couple_idx, size_t field_couple_idx)
static void yac_couple_config_pack_field_couple(struct yac_couple_config_field_couple *field_couple, void *buffer, int buffer_size, int *position, MPI_Comm comm)
double yac_couple_config_get_frac_mask_fallback_value(struct yac_couple_config *couple_config, char const *component_name, char const *grid_name, char const *field_name)
static size_t yac_couple_config_get_fields_pack_size(size_t num_fields, void *fields_, MPI_Comm comm)
int yac_couple_config_get_source_lag(struct yac_couple_config *couple_config, size_t couple_idx, size_t field_couple_idx)
static int yac_couple_config_couple_compare_basic(void const *a_, void const *b_)
const char * yac_couple_config_field_get_metadata(struct yac_couple_config *couple_config, const char *comp_name, const char *grid_name, const char *field_name)
static void merge_grids(struct yac_couple_config *couple_config, MPI_Comm comm)
void yac_couple_config_get_couple_component_names(struct yac_couple_config *couple_config, size_t couple_idx, char const *couple_component_names[2])
const char * yac_couple_config_grid_get_metadata(struct yac_couple_config *couple_config, const char *grid_name)
static void merge_field_couples(size_t *num_field_couples, struct yac_couple_config_field_couple **field_couples, MPI_Comm comm)
static void yac_couple_config_unpack_grids(void *buffer, int buffer_size, int *position, size_t *num_grids, void *grids_, MPI_Comm comm)
const char * yac_couple_config_component_get_metadata(struct yac_couple_config *couple_config, const char *comp_name)
double yac_couple_config_get_scale_factor(struct yac_couple_config *couple_config, size_t couple_idx, size_t field_couple_idx)
char const * yac_couple_config_get_end_datetime(struct yac_couple_config *couple_config)
static void yac_couple_config_pack_fields(size_t num_fields, void *fields_, void *buffer, int buffer_size, int *position, MPI_Comm comm)
size_t yac_couple_config_get_grid_idx(struct yac_couple_config *couple_config, char const *grid_name)
size_t yac_couple_config_get_component_idx(struct yac_couple_config *couple_config, char const *component_name)
char const * yac_couple_config_get_tgt_mask_name(struct yac_couple_config *couple_config, size_t couple_idx, size_t field_couple_idx)
void yac_couple_config_add_grid(struct yac_couple_config *couple_config, char const *name)
char const * yac_couple_config_get_grid_name(struct yac_couple_config *couple_config, size_t grid_idx)
double yac_couple_config_get_scale_summand(struct yac_couple_config *couple_config, size_t couple_idx, size_t field_couple_idx)
static void yac_couple_config_unpack_fields(void *buffer, int buffer_size, int *position, size_t *num_fields, void *fields_, MPI_Comm comm)
enum yac_reduction_type yac_couple_config_get_coupling_period_operation(struct yac_couple_config *couple_config, size_t couple_idx, size_t field_couple_idx)
static void yac_couple_config_pack_grids(size_t num_grids, void *grids_, void *buffer, int buffer_size, int *position, MPI_Comm comm)
static size_t yac_couple_config_get_string_pack_size(char const *string, MPI_Comm comm)
int yac_couple_config_enforce_write_weight_file(struct yac_couple_config *couple_config, size_t couple_idx, size_t field_couple_idx)
static void yac_couple_config_grid_free(void *grid_)
static size_t yac_couple_config_get_couples_pack_size_basic(size_t num_couples, void *couples_, MPI_Comm comm)
static void yac_couple_config_unpack_field_couple(void *buffer, int buffer_size, int *position, struct yac_couple_config_field_couple *field_couple, MPI_Comm comm)
void yac_couple_config_def_couple(struct yac_couple_config *couple_config, char const *src_comp_name, char const *src_grid_name, char const *src_field_name, char const *tgt_comp_name, char const *tgt_grid_name, char const *tgt_field_name, char const *coupling_period, int time_reduction, struct yac_interp_stack_config *interp_stack, int src_lag, int tgt_lag, const char *weight_file_name, int mapping_on_source, double scale_factor, double scale_summand, size_t num_src_mask_names, char const *const *src_mask_names, char const *tgt_mask_name)
static size_t yac_couple_config_get_grids_pack_size(size_t num_grids, void *grids_, MPI_Comm comm)
struct dist_merge_vtable dist_merge_vtable_grid
void yac_couple_config_get_field_names(struct yac_couple_config *couple_config, size_t couple_idx, size_t field_couple_idx, char const **src_field_name, const char **tgt_field_name)
struct yac_couple_config * yac_couple_config_new()
static size_t yac_couple_config_get_component_pack_size(struct yac_couple_config_component *component, MPI_Comm comm)
char const * yac_couple_config_get_weight_file_name(struct yac_couple_config *couple_config, size_t couple_idx, size_t field_couple_idx)
struct dist_merge_vtable dist_merge_vtable_couple
static void check_grid_idx(struct yac_couple_config *couple_config, size_t grid_idx, char const *routine_name, int line)
static size_t yac_couple_config_component_add_field_(struct yac_couple_config *couple_config, size_t comp_idx, size_t grid_idx, char const *name, char const *timestep, size_t collection_size)
static void yac_couple_config_pack_component(struct yac_couple_config_component *component, void *buffer, int buffer_size, int *position, MPI_Comm comm)
void yac_couple_config_get_field_couple_component_names(struct yac_couple_config *couple_config, size_t couple_idx, size_t field_couple_idx, char const **src_component_name, char const **tgt_component_name)
static void yac_couple_config_field_free(void *field_)
char const * yac_couple_config_get_component_name(struct yac_couple_config *couple_config, size_t component_idx)
static void yac_couple_config_field_couple_merge(void *a_, void *b_, MPI_Comm comm)
static int yac_couple_config_field_couple_compare(void const *a_, void const *b_)
char const * yac_couple_config_get_source_timestep(struct yac_couple_config *couple_config, size_t couple_idx, size_t field_couple_idx)
static void yac_couple_config_component_free(void *component_)
static void yac_couple_config_pack_string(char const *string, void *buffer, int buffer_size, int *position, MPI_Comm comm)
static void yac_couple_config_pack_grid(struct yac_couple_config_grid *grid, void *buffer, int buffer_size, int *position, MPI_Comm comm)
void yac_couple_config_add_component(struct yac_couple_config *couple_config, char const *name)
size_t yac_couple_config_get_num_grids(struct yac_couple_config *couple_config)
struct dist_merge_vtable dist_merge_vtable_field_couple
char const * yac_couple_config_get_field_grid_name(struct yac_couple_config *couple_config, size_t component_idx, size_t field_idx)
char const * yac_couple_config_get_start_datetime(struct yac_couple_config *couple_config)
static void check_field_idx(struct yac_couple_config *couple_config, size_t component_idx, size_t field_idx, char const *routine_name, int line)
size_t yac_couple_config_get_collection_size(struct yac_couple_config *couple_config, char const *component_name, char const *grid_name, char const *field_name)
static void check_component_idx(struct yac_couple_config *couple_config, size_t component_idx, char const *routine_name, int line)
static void merge_couples(struct yac_couple_config *couple_config, MPI_Comm comm)
static void yac_couple_config_unpack_components(void *buffer, int buffer_size, int *position, size_t *num_components, void *components_, MPI_Comm comm)
char const * yac_couple_config_get_target_timestep(struct yac_couple_config *couple_config, size_t couple_idx, size_t field_couple_idx)
static void yac_couple_config_field_merge(void *a_, void *b_, MPI_Comm comm)
void yac_couple_config_delete(struct yac_couple_config *couple_config)
static void yac_couple_config_pack_field(struct yac_couple_config_field *field, void *buffer, int buffer_size, int *position, MPI_Comm comm)
size_t yac_couple_config_get_num_components(struct yac_couple_config *couple_config)
yac_reduction_type
@ TIME_NONE
@ TIME_ACCUMULATE
@ TIME_MAXIMUM
@ TIME_MINIMUM
@ TIME_AVERAGE
int yac_interp_stack_config_compare(void const *a_, void const *b_)
struct yac_interp_stack_config * yac_interp_stack_config_unpack(void *buffer, int buffer_size, int *position, MPI_Comm comm)
void yac_interp_stack_config_delete(struct yac_interp_stack_config *interp_stack_config)
struct yac_interp_stack_config * yac_interp_stack_config_copy(struct yac_interp_stack_config *interp_stack)
size_t yac_interp_stack_config_get_pack_size(struct yac_interp_stack_config *interp_stack, MPI_Comm comm)
void yac_interp_stack_config_pack(struct yac_interp_stack_config *interp_stack, void *buffer, int buffer_size, int *position, MPI_Comm comm)
double const YAC_FRAC_MASK_UNDEF
double const YAC_FRAC_MASK_NO_VALUE
#define xrealloc(ptr, size)
Definition ppm_xfuncs.h:67
#define xmalloc(size)
Definition ppm_xfuncs.h:66
void yac_qsort_index(void *a_, size_t count, size_t size, int(*compare)(void const *, void const *), size_t *idx)
Definition quicksort.c:185
void(* pack)(size_t, void *, void *, int, int *, MPI_Comm)
size_t(* get_pack_size)(size_t, void *, MPI_Comm)
void(* free_data)(void *)
void(* unpack)(void *, int, int *, size_t *, void *, MPI_Comm)
int(* compare)(void const *, void const *)
void(* merge)(void *, void *, MPI_Comm)
struct yac_couple_config_field * fields
struct yac_couple_config_field_couple * field_couples
struct yac_interp_stack_config * interp_stack
enum yac_reduction_type coupling_period_operation
struct yac_couple_config_field_couple::@58 source
struct yac_couple_config_field_couple::@58 target
struct yac_couple_config_couple * couples
char const * end_datetime
struct yac_couple_config_grid * grids
char const * start_datetime
struct yac_couple_config_component * components
static struct user_input_data_component ** components
Definition yac.c:117
size_t num_grids
Definition yac.c:115
static size_t num_components
Definition yac.c:118
int const YAC_EXCHANGE_TYPE_SOURCE
Definition yac.c:32
int const YAC_EXCHANGE_TYPE_NONE
Definition yac.c:31
int const YAC_EXCHANGE_TYPE_TARGET
Definition yac.c:33
void yac_cdef_calendar(int calendar)
Definition yac.c:562
struct yac_basic_grid ** grids
Definition yac.c:114
#define YAC_MAX_CHARLEN
Definition yac.h:74
#define YAC_ASSERT_F(exp, format,...)
Definition yac_assert.h:18
#define YAC_ASSERT(exp, msg)
Definition yac_assert.h:15
#define yac_mpi_call(call, comm)