YetAnotherCoupler 3.5.2
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#include "config_yaml.h"
20#include "dist_merge.h"
21#include "mtime_datetime.h"
22
28
30 char const * name; //< name of the file
31 enum yac_text_filetype type; //< type of the file
32 char const * ref; //< name of the synchronization location
33 //< at which this file is to be written
34 //< out
35 enum flag_value include_definitions; //< include user definitions (components,
36 //< grids, and fields)
37};
38
40 char const * name;
41 char* metadata;
42 char const * output_filename;
43};
44
46 char const * name;
47 size_t grid_idx;
48 char const * timestep;
49 char * metadata;
52};
53
55
57 size_t num_fields;
58
59 char const * name;
60 char * metadata;
61};
62
89
97
99
100 struct _datetime * start_datetime;
101 struct _datetime * end_datetime;
102
105
107 size_t num_grids;
108
111
114
115 // this flag indicates whether YAC aborts if for a defined couple at least
116 // one associated field was not defined by the user
118};
119
121
122 struct yac_couple_config * couple_config =
123 xmalloc(1 * sizeof(*couple_config));
124
125 couple_config->start_datetime = NULL;
126 couple_config->end_datetime = NULL;
127
128 couple_config->config_outputs = NULL;
129 couple_config->num_config_outputs = 0;
130
131 couple_config->grids = NULL;
132 couple_config->num_grids = 0;
133
134 couple_config->components = NULL;
135 couple_config->num_components = 0;
136
137 couple_config->couples = NULL;
138 couple_config->num_couples = 0;
139
141
142 return couple_config;
143}
144
145static char * string_dup(char const * string) {
146 return (string != NULL)?strdup(string):NULL;
147}
148
149static void yac_couple_config_field_free(void * field_){
150 struct yac_couple_config_field * field = field_;
151 free((void*)field->name);
152 free((void*)field->timestep);
153 free(field->metadata);
154}
155
156
157static void yac_couple_config_component_free(void * component_) {
158
159 struct yac_couple_config_component * component = component_;
160
161 for (size_t i = 0; i < component->num_fields; ++i)
162 yac_couple_config_field_free(component->fields + i);
163 free(component->fields);
164 free((void*)(component->name));
165 free(component->metadata);
166}
167
168static void yac_couple_config_field_couple_free(void * field_couple_) {
169
170 struct yac_couple_config_field_couple * field_couple = field_couple_;
172 free((void*)field_couple->coupling_period);
173 free((void*)field_couple->weight_file_name);
174 for (size_t i = 0; i < field_couple->num_src_mask_names; ++i)
175 free((void*)field_couple->src_mask_names[i]);
176 free(field_couple->src_mask_names);
177 free(field_couple->tgt_mask_name);
178 free(field_couple->yaxt_exchanger_name);
179}
180
181static void yac_couple_config_couple_free(void * couple_) {
182
183 struct yac_couple_config_couple * couple = couple_;
184 for (size_t i = 0; i < couple->num_field_couples; ++i)
186 free(couple->field_couples);
187 couple->field_couples = NULL;
188 couple->num_field_couples = 0;
189}
190
191static void yac_couple_config_grid_free(void * grid_) {
192
193 struct yac_couple_config_grid * grid = grid_;
194
195 free((void*)grid->name);
196 free((void*)grid->output_filename);
197 free(grid->metadata);
198}
199
201 void * config_output_) {
202
203 struct yac_couple_config_config_output * config_output = config_output_;
204
205 free((void*)config_output->name);
206 free((void*)config_output->ref);
207}
208
209static int yac_couple_config_grid_compare(void const * a, void const * b) {
210 YAC_ASSERT(((struct yac_couple_config_grid const *)a)->name != NULL &&
211 ((struct yac_couple_config_grid const *)b)->name != NULL,
212 "ERROR(yac_couple_config_grid_compare): "
213 "invalid name (NULL is not allowed)");
214 return strcmp(((struct yac_couple_config_grid const *)a)->name,
215 ((struct yac_couple_config_grid const *)b)->name);
216}
217
218static int yac_couple_config_field_compare(void const * a_, void const * b_) {
219 struct yac_couple_config_field const * a = a_;
220 struct yac_couple_config_field const * b = b_;
221 int ret = (a->grid_idx > b->grid_idx) - (a->grid_idx < b->grid_idx);
222 if (ret) return ret;
224 a->name != NULL && b->name != NULL,
225 "ERROR(yac_couple_config_field_compare): "
226 "invalid name (NULL is not allowed)");
227 return strcmp(a->name, b->name);
228}
229
231 void const * a_, void const * b_) {
232 struct yac_couple_config_component const * a = a_;
233 struct yac_couple_config_component const * b = b_;
234 YAC_ASSERT(a->name != NULL && b->name != NULL,
235 "ERROR(yac_couple_config_component_compare): "
236 "invalid name (NULL is not allowed)");
237 return strcmp(a->name, b->name);
238}
239
241 void const * a_, void const * b_) {
242 struct yac_couple_config_field_couple const * a = a_;
243 struct yac_couple_config_field_couple const * b = b_;
244 int ret;
245 ret = (a->source.component_idx > b->source.component_idx) -
247 if (ret) return ret;
248 ret = (a->target.component_idx > b->target.component_idx) -
250 if (ret) return ret;
251 ret = (a->source.field_idx > b->source.field_idx) -
253 if (ret) return ret;
254 return (a->target.field_idx > b->target.field_idx) -
256}
257
259 void const * a_, void const * b_) {
260 struct yac_couple_config_couple const * a = a_;
261 struct yac_couple_config_couple const * b = b_;
262 int ret;
263 ret = (a->component_indices[0] > b->component_indices[0]) -
264 (a->component_indices[0] < b->component_indices[0]);
265 if (ret) return ret;
266 return (a->component_indices[1] > b->component_indices[1]) -
267 (a->component_indices[1] < b->component_indices[1]);
268}
269
271 void const * a, void const * b) {
273 ((struct yac_couple_config_config_output const *)a)->ref != NULL &&
274 ((struct yac_couple_config_config_output const *)b)->ref != NULL,
275 "ERROR(yac_couple_config_config_output_compare): "
276 "invalid ref (NULL is not allowed)");
277 return strcmp(((struct yac_couple_config_config_output const *)a)->ref,
278 ((struct yac_couple_config_config_output const *)b)->ref);
279}
280
282 char const * string_name, char ** string, MPI_Comm comm) {
283
284 int rank;
285 MPI_Comm_rank(comm, &rank);
286 struct {
287 int len;
288 int rank;
289 } data_pair;
290 size_t len = (*string != NULL)?(strlen(*string) + 1):0;
292 len <= INT_MAX,
293 "ERROR(couple_config_sync_string): \"%s\" too long", string_name);
294 data_pair.len = (int)len;
295 data_pair.rank = rank;
296
297 // determine the broadcaster
299 MPI_Allreduce(
300 MPI_IN_PLACE, &data_pair, 1, MPI_2INT, MPI_MAXLOC, comm), comm);
301 if (data_pair.len == 0) return;
302
303 // broadcast string
304 char const * string_bak = NULL;
305 if (data_pair.rank != rank) {
306 string_bak = *string;
307 *string = xmalloc((size_t)data_pair.len * sizeof(**string));
308 }
310 MPI_Bcast(*string, data_pair.len, MPI_CHAR, data_pair.rank, comm), comm);
311
312 // check for consistency
313 if (data_pair.rank != rank) {
315 (string_bak == NULL) ||
316 !strcmp(string_bak, *string),
317 "ERROR(couple_config_sync_string): inconsistent \"%s\" definition "
318 "(\"%s\" != \"%s\")", string_name, string_bak, *string);
319 free((void*)string_bak);
320 }
321}
322
324 char const * flag_value_name, enum flag_value * flag, MPI_Comm comm) {
325
326 int int_flag = (*flag == FLAG_UNSET)?INT_MAX:((int)*flag);
327
328 // broadcast flag value
330 MPI_Allreduce(
331 MPI_IN_PLACE, &int_flag, 1, MPI_INT, MPI_MIN, comm), comm);
332
333 // if no process has defined a flag value
334 if (int_flag == INT_MAX) return;
335
336 // check for consistency
338 (*flag == FLAG_UNSET) || ((int)(*flag) == int_flag),
339 "ERROR(couple_config_sync_flag_value): inconsistent \"%s\"-flag definition "
340 "(%s != %s)", flag_value_name, ((int_flag == FLAG_TRUE)?"TRUE":"FALSE"),
341 ((*flag == FLAG_TRUE)?"TRUE":"FALSE"));
342
343 // update flag
344 *flag = (enum flag_value)int_flag;
345}
346
347static void couple_config_sync_calendar(MPI_Comm comm) {
348
349 int calendar = (int)getCalendarType();
350 if (calendar == CALENDAR_NOT_SET) calendar = INT_MAX;
351
352 // broadcast calendar
354 MPI_Allreduce(
355 MPI_IN_PLACE, &calendar, 1, MPI_INT, MPI_MIN, comm), comm);
356
357 // if no process has defined a calendar
358 if (calendar == INT_MAX) return;
359
360 // set calendar (if local process has already defined a calendar,
361 // the definition is checked for consistency)
362 yac_cdef_calendar(calendar);
363}
364
366 void * a_, void * b_, MPI_Comm comm) {
367
368 struct yac_couple_config_component * a = a_;
369 struct yac_couple_config_component * b = b_;
370
371 if (b) {
372 a->num_fields = b->num_fields;
373 a->fields = b->fields;
374 b->num_fields = 0;
375 b->fields = NULL;
376 a->metadata = b->metadata;
377 b->metadata = NULL;
378 }
379
381 "component metadata", (char **)&(a->metadata), comm);
382}
383
385 void * a_, void * b_, MPI_Comm comm) {
386
387 struct yac_couple_config_grid * a = a_;
388 struct yac_couple_config_grid * b = b_;
389
390 if (b!= NULL) {
392 b->output_filename = NULL;
393 a->metadata = b->metadata;
394 b->metadata = NULL;
395 }
396
398 "grid output_filename", (char **)&(a->output_filename), comm);
399 couple_config_sync_string("grid metadata", (char **)&(a->metadata), comm);
400}
401
403 void * a_, void * b_, MPI_Comm comm) {
404
405 struct yac_couple_config_field * a = a_;
406 struct yac_couple_config_field * b = b_;
407
408 enum {
409 TIMESTEP_IDX,
410 METADATA_IDX,
411 FRAC_MASK_IDX,
412 COLLECTION_SIZE_IDX,
413 DATA_PAIR_COUNT
414 };
415
416 int rank;
417 MPI_Comm_rank(comm, &rank);
418 struct {
419 int len;
420 int rank;
421 } data_pairs[DATA_PAIR_COUNT];
422 size_t timestep_len =
423 ((b != NULL) && (b->timestep != NULL))?(strlen(b->timestep) + 1):0;
424 size_t metadata_len =
425 ((b != NULL) && (b->metadata != NULL))?(strlen(b->metadata) + 1):0;
427 timestep_len <= INT_MAX,
428 "ERROR(yac_couple_config_field_merge): timestep string too long");
430 metadata_len <= INT_MAX,
431 "ERROR(yac_couple_config_field_merge): metadata string too long");
433 (b == NULL) || (b->collection_size <= INT_MAX) ||
434 (b->collection_size == SIZE_MAX),
435 "ERROR(yac_couple_config_field_merge): invalid collection size \"%zu\"",
436 b->collection_size);
437 data_pairs[TIMESTEP_IDX].len = (int)timestep_len;
438 data_pairs[TIMESTEP_IDX].rank = rank;
439 data_pairs[METADATA_IDX].len = (int)metadata_len;
440 data_pairs[METADATA_IDX].rank = rank;
441 data_pairs[FRAC_MASK_IDX].len =
443 data_pairs[FRAC_MASK_IDX].rank = rank;
444 data_pairs[COLLECTION_SIZE_IDX].len =
445 ((b == NULL) || (b->collection_size == SIZE_MAX))?
446 -1:(int)b->collection_size;
447 data_pairs[COLLECTION_SIZE_IDX].rank = rank;
448
449 // determine the broadcaster
451 MPI_Allreduce(
452 MPI_IN_PLACE, data_pairs, DATA_PAIR_COUNT, MPI_2INT,
453 MPI_MAXLOC, comm), comm);
454
455 // if at least one process has a valid timestep
456 if (data_pairs[TIMESTEP_IDX].len > 0) {
457
458 // broadcast timestep
459 char * timestep_buffer = xmalloc((size_t)data_pairs[TIMESTEP_IDX].len);
460 if (data_pairs[TIMESTEP_IDX].rank == rank)
461 strcpy(timestep_buffer, b->timestep);
463 MPI_Bcast(
464 timestep_buffer, data_pairs[TIMESTEP_IDX].len, MPI_CHAR,
465 data_pairs[TIMESTEP_IDX].rank, comm), comm);
466
467 // check for consistency
469 (b == NULL) || (b->timestep == NULL) ||
470 !strcmp(b->timestep, timestep_buffer),
471 "ERROR(yac_couple_config_field_merge): "
472 "inconsistent timestep definition (\"%s\" != \"%s\")",
473 b->timestep, timestep_buffer);
474
475 // update timestep
476 free((void*)(a->timestep));
477 a->timestep = timestep_buffer;
478 }
479
480 // if at least one process has a valid metadata
481 if (data_pairs[METADATA_IDX].len > 0) {
482
483 // broadcast metadata
484 char * metadata_buffer = xmalloc((size_t)data_pairs[METADATA_IDX].len);
485 if (data_pairs[METADATA_IDX].rank == rank)
486 strcpy(metadata_buffer, b->metadata);
488 MPI_Bcast(
489 metadata_buffer, data_pairs[METADATA_IDX].len, MPI_CHAR,
490 data_pairs[METADATA_IDX].rank, comm), comm);
491
492 // check for consistency
494 (b == NULL) || (b->metadata == NULL) ||
495 !strcmp(b->metadata, metadata_buffer),
496 "ERROR(yac_couple_config_field_merge): "
497 "inconsistent metadata definition (\"%s\" != \"%s\")",
498 b->metadata, metadata_buffer);
499
500 // update metadata
501 free(a->metadata);
502 a->metadata = metadata_buffer;
503 }
504
505 // if at least one process has a valid fractional mask fallback value
506 if (data_pairs[FRAC_MASK_IDX].len != 0) {
507
508 // broadcast fractional mask fallback value
509 double frac_mask_fallback_value;
510 if (data_pairs[FRAC_MASK_IDX].rank == rank)
511 frac_mask_fallback_value = b->frac_mask_fallback_value;
513 MPI_Bcast(
514 &frac_mask_fallback_value, 1, MPI_DOUBLE,
515 data_pairs[FRAC_MASK_IDX].rank, comm), comm);
516
517 // check for consistency
518 // (use memcmp for comparison, because it can be nan)
520 (b == NULL) ||
522 !memcmp(
523 &b->frac_mask_fallback_value, &frac_mask_fallback_value,
524 sizeof(frac_mask_fallback_value)),
525 "ERROR(yac_couple_config_field_merge): "
526 "inconsistent fractional mask fallback value definition "
527 "(%lf != %lf)",
528 b->frac_mask_fallback_value, frac_mask_fallback_value);
529
530 // update fractional mask fallback value
531 a->frac_mask_fallback_value = frac_mask_fallback_value;
532 }
533
534
535 // if at least one process has a valid collection size
536 if (data_pairs[COLLECTION_SIZE_IDX].len > -1) {
537
538 // check for consistency
540 (b == NULL) || (b->collection_size == SIZE_MAX) ||
541 (b->collection_size == (size_t)(data_pairs[COLLECTION_SIZE_IDX].len)),
542 "ERROR(yac_couple_config_field_merge): "
543 "inconsistent collection size definition (\"%zu\" != \"%d\")",
544 b->collection_size, data_pairs[COLLECTION_SIZE_IDX].len);
545
546 // update collection size
547 a->collection_size = (size_t)(data_pairs[COLLECTION_SIZE_IDX].len);
548 }
549}
550
552 void * a_, void * b_, MPI_Comm comm) {
553
554 UNUSED(comm);
555
556 if (b_ == NULL) return;
557
558 struct yac_couple_config_field_couple * a = a_;
559 struct yac_couple_config_field_couple * b = b_;
560
562 (a->source.lag == b->source.lag),
563 "ERROR(yac_couple_config_field_couple_merge): "
564 "inconsistent source lag (%d != %d)", a->source.lag, b->source.lag)
566 (a->target.lag == b->target.lag),
567 "ERROR(yac_couple_config_field_couple_merge): "
568 "inconsistent target lag (%d != %d)", a->target.lag, b->target.lag)
571 "ERROR(yac_couple_config_field_couple_merge): "
572 "inconsistent mapping side (%d != %d)",
576 "ERROR(yac_couple_config_field_couple_merge): "
577 "inconsistent interpolation stack")
580 "ERROR(yac_couple_config_field_couple_merge): "
581 "inconsistent coupling period operation (%d != %d)",
584 !strcmp(a->coupling_period, b->coupling_period),
585 "ERROR(yac_couple_config_field_couple_merge): "
586 "inconsistent coupling period (%s != %s)",
590 "ERROR(yac_couple_config_field_couple_merge): "
591 "inconsistent defintion of enforce_write_weight_file (%d != %d)",
595 !strcmp(a->weight_file_name, b->weight_file_name),
596 "ERROR(yac_couple_config_field_couple_merge): "
597 "inconsistent weight_file_name (%s != %s)",
600 a->scale_factor == b->scale_factor,
601 "ERROR(yac_couple_config_field_couple_merge): "
602 "inconsistent scale factor (%lf != %lf)",
606 "ERROR(yac_couple_config_field_couple_merge): "
607 "inconsistent scale summand (%lf != %lf)",
611 "ERROR(yac_couple_config_field_couple_merge): "
612 "inconsistent number of source mask names (%zu != %zu)",
615 (a->src_mask_names != NULL) == (b->src_mask_names != NULL),
616 "ERROR(yac_couple_config_field_couple_merge): "
617 "inconsistent availability of source mask names (%s != %s)",
618 (a->src_mask_names != NULL)?"TRUE":"FALSE",
619 (b->src_mask_names != NULL)?"TRUE":"FALSE")
620 for (size_t i = 0; i < a->num_src_mask_names; ++i)
622 !strcmp(a->src_mask_names[i], b->src_mask_names[i]),
623 "ERROR(yac_couple_config_field_couple_merge): "
624 "inconsistent source mask names at index %zu (\"%s\" != \"%s\")",
625 i, a->src_mask_names[i], b->src_mask_names[i])
627 (a->tgt_mask_name != NULL) == (b->tgt_mask_name != NULL),
628 "ERROR(yac_couple_config_field_couple_merge): "
629 "inconsistent availability of target mask name (%s != %s)",
630 (a->tgt_mask_name != NULL)?"TRUE":"FALSE",
631 (b->tgt_mask_name != NULL)?"TRUE":"FALSE")
633 (a->tgt_mask_name == NULL) ||
634 !strcmp(a->tgt_mask_name, b->tgt_mask_name),
635 "ERROR(yac_couple_config_field_couple_merge): "
636 "inconsistent target mask name (\"%s\" != \"%s\")",
639 (a->yaxt_exchanger_name != NULL) == (b->yaxt_exchanger_name != NULL),
640 "ERROR(yac_couple_config_field_couple_merge): "
641 "inconsistent availability of yaxt exchanger name (%s != %s)",
642 (a->yaxt_exchanger_name != NULL)?"TRUE":"FALSE",
643 (b->yaxt_exchanger_name != NULL)?"TRUE":"FALSE")
645 (a->yaxt_exchanger_name == NULL) ||
647 "ERROR(yac_couple_config_field_couple_merge): "
648 "inconsistent yaxt exchanger name (\"%s\" != \"%s\")",
650}
651
652static void merge_field_couples(
653 size_t * num_field_couples,
654 struct yac_couple_config_field_couple ** field_couples, MPI_Comm comm);
655
657 void * a_, void * b_, MPI_Comm comm) {
658
659 struct yac_couple_config_couple * a = a_;
660 struct yac_couple_config_couple * b = b_;
661
662 if (b) {
665 b->num_field_couples = 0;
666 b->field_couples = NULL;
667 }
668
669 // distribute and merge field couples
671}
672
674 void * a_, void * b_, MPI_Comm comm) {
675
676 struct yac_couple_config_config_output * a = a_;
677 struct yac_couple_config_config_output * b = b_;
678
679 if (b) {
680
682 !strcmp(a->name, b->name),
683 "ERROR(yac_couple_config_config_output_merge): "
684 "inconsistent file names for ref \"%s\" (\"%s\" != \"%s\")",
685 a->ref, a->name, b->name)
687 a->type == b->type,
688 "ERROR(yac_couple_config_config_output_merge): "
689 "inconsistent file types for ref \"%s\" (%d != %d)",
690 a->ref, a->type, b->type)
691
693 }
694
696 "include definitions", &(a->include_definitions), comm);
697}
698
699void yac_couple_config_delete(struct yac_couple_config * couple_config) {
700
701 deallocateDateTime((void*)couple_config->start_datetime);
702 deallocateDateTime((void*)couple_config->end_datetime);
703
704 for (size_t i = 0; i < couple_config->num_config_outputs; ++i)
706 couple_config->config_outputs + i);
707 free((void*)couple_config->config_outputs);
708
709 for (size_t i = 0; i < couple_config->num_grids; ++i){
710 yac_couple_config_grid_free(couple_config->grids + i);
711 }
712 free(couple_config->grids);
713
714 for (size_t i = 0;
715 i < couple_config->num_components; ++i)
717 free(couple_config->components);
718
719 for (size_t couple_idx = 0; couple_idx < couple_config->num_couples;
720 ++couple_idx)
721 yac_couple_config_couple_free(couple_config->couples + couple_idx);
722 free(couple_config->couples);
723 free(couple_config);
724}
725
727 struct yac_couple_config * couple_config, char const * name) {
728
730 name != NULL,
731 "ERROR(yac_couple_config_add_grid_): "
732 "invalid name (NULL is not allowed)")
733
734 for (size_t i = 0; i < couple_config->num_grids; ++i)
735 if (!strcmp(couple_config->grids[i].name, name)) return i;
736
737 size_t grid_idx = couple_config->num_grids;
738 couple_config->num_grids++;
739 couple_config->grids =
740 xrealloc(
741 couple_config->grids,
742 couple_config->num_grids * sizeof(*(couple_config->grids)));
743 couple_config->grids[grid_idx].name = string_dup(name);
744 couple_config->grids[grid_idx].output_filename = NULL;
745 couple_config->grids[grid_idx].metadata = NULL;
746 return grid_idx;
747}
748
750 struct yac_couple_config * couple_config, char const * name) {
751
752 yac_couple_config_add_grid_(couple_config, name);
753}
754
756 struct yac_couple_config * couple_config, char const * name) {
757
759 name != NULL,
760 "ERROR(yac_couple_config_add_component_): "
761 "invalid name (NULL is not allowed)")
762
763 for (size_t i = 0; i < couple_config->num_components; ++i)
764 if (!strcmp(couple_config->components[i].name, name)) return i;
765
766 size_t component_idx = couple_config->num_components;
767 couple_config->num_components++;
768 couple_config->components =
769 xrealloc(
770 couple_config->components,
771 couple_config->num_components * sizeof(*(couple_config->components)));
772 couple_config->components[component_idx].fields = NULL;
773 couple_config->components[component_idx].num_fields = 0;
774 couple_config->components[component_idx].name = strdup(name);
775 couple_config->components[component_idx].metadata = NULL;
776 return component_idx;
777}
778
780 struct yac_couple_config * couple_config, char const * name) {
781
783}
784
786 struct yac_couple_config * couple_config,
787 char const * comp_name, const char* metadata) {
788 size_t comp_idx = yac_couple_config_get_component_idx(couple_config, comp_name);
789 free(couple_config->components[comp_idx].metadata);
790 couple_config->components[comp_idx].metadata
791 = metadata==NULL?NULL:strdup(metadata);
792}
793
795 struct yac_couple_config * couple_config,
796 char const * grid_name, const char* output_filename) {
797 if(!yac_couple_config_contains_grid_name(couple_config, grid_name))
798 yac_couple_config_add_grid(couple_config, grid_name);
799 size_t grid_idx = yac_couple_config_get_grid_idx(couple_config, grid_name);
801 output_filename != NULL,
802 "ERROR(yac_couple_config_grid_set_output_filename): invalid output filename "
803 "for grid \"%s\" (NULL is not allowed)",
804 grid_name);
806 (couple_config->grids[grid_idx].output_filename == NULL) ||
807 (!strcmp(couple_config->grids[grid_idx].output_filename, output_filename)),
808 "ERROR(yac_couple_config_grid_set_output_filename): output file name for grid "
809 "\"%s\" has already been set (\"%s\" != \"%s\")",
810 grid_name, couple_config->grids[grid_idx].output_filename, output_filename);
811 if (couple_config->grids[grid_idx].output_filename == NULL)
812 couple_config->grids[grid_idx].output_filename = strdup(output_filename);
813}
814
816 struct yac_couple_config * couple_config,
817 char const * grid_name, const char* metadata) {
818 if(!yac_couple_config_contains_grid_name(couple_config, grid_name))
819 yac_couple_config_add_grid(couple_config, grid_name);
820 size_t grid_idx = yac_couple_config_get_grid_idx(couple_config, grid_name);
821 free(couple_config->grids[grid_idx].metadata);
822 couple_config->grids[grid_idx].metadata
823 = metadata==NULL?NULL:strdup(metadata);
824}
825
827 struct yac_couple_config * couple_config,
828 const char* comp_name, const char * grid_name, const char* field_name,
829 const char* metadata) {
830 size_t comp_idx = yac_couple_config_get_component_idx(couple_config, comp_name);
831 size_t grid_idx = yac_couple_config_get_grid_idx(couple_config, grid_name);
832 size_t field_idx = yac_couple_config_get_field_idx(couple_config,
833 comp_idx, grid_idx, field_name);
834 free(couple_config->components[comp_idx].fields[field_idx].metadata);
835 couple_config->components[comp_idx].fields[field_idx].metadata
836 = metadata==NULL?NULL:strdup(metadata);
837}
838
840 struct yac_couple_config * couple_config,
841 const char * comp_name) {
842 size_t comp_idx = yac_couple_config_get_component_idx(couple_config, comp_name);
843 return couple_config->components[comp_idx].metadata;
844}
845
847 struct yac_couple_config * couple_config,
848 const char * grid_name) {
849 size_t grid_idx = yac_couple_config_get_grid_idx(couple_config, grid_name);
850 return couple_config->grids[grid_idx].output_filename;
851}
852
854 struct yac_couple_config * couple_config,
855 const char * grid_name) {
856 size_t grid_idx = yac_couple_config_get_grid_idx(couple_config, grid_name);
857 return couple_config->grids[grid_idx].metadata;
858}
859
861 struct yac_couple_config * couple_config,
862 const char* comp_name, const char * grid_name, const char* field_name) {
863 size_t comp_idx = yac_couple_config_get_component_idx(couple_config, comp_name);
864 size_t grid_idx = yac_couple_config_get_grid_idx(couple_config, grid_name);
865 size_t field_idx = yac_couple_config_get_field_idx(couple_config,
866 comp_idx, grid_idx, field_name);
867 return couple_config->components[comp_idx].fields[field_idx].metadata;
868}
869
871 struct yac_couple_config * couple_config, size_t component_idx,
872 char const * routine_name, int line) {
873
875 component_idx < couple_config->num_components,
876 "ERROR(%s:%d:%s): invalid component_idx", __FILE__, line, routine_name)
877}
878
879static void check_grid_idx(
880 struct yac_couple_config * couple_config, size_t grid_idx,
881 char const * routine_name, int line) {
882
884 grid_idx < couple_config->num_grids,
885 "ERROR(%s:%d:%s): invalid grid_idx", __FILE__, line, routine_name)
886}
887
889 struct yac_couple_config * couple_config,
890 size_t comp_idx, size_t grid_idx, char const * name,
891 char const * timestep, size_t collection_size) {
892
894 couple_config, comp_idx,
895 "yac_couple_config_component_add_field_", __LINE__);
897 couple_config, grid_idx,
898 "yac_couple_config_component_add_field_", __LINE__);
899
900 // check whether the field already exists
901 struct yac_couple_config_component * component =
902 couple_config->components + comp_idx;
903 for (size_t i = 0; i < component->num_fields; i++) {
904 if((strcmp(component->fields[i].name, name) == 0) &&
905 (component->fields[i].grid_idx == grid_idx)) {
906 // if no timestep is defined for the field
907 if (!component->fields[i].timestep && timestep)
908 component->fields[i].timestep = strdup(timestep);
909 // if no collection size is defined for the field
910 if (component->fields[i].collection_size == SIZE_MAX)
911 component->fields[i].collection_size = collection_size;
913 !timestep ||
914 !strcmp(timestep, component->fields[i].timestep),
915 "ERROR(yac_couple_config_component_add_field): "
916 "inconsistent timestep definition (\"%s\" != \"%s\")",
917 timestep, component->fields[i].timestep);
919 collection_size == SIZE_MAX ||
920 collection_size == component->fields[i].collection_size,
921 "ERROR(yac_couple_config_component_add_field): "
922 "inconsistent collection_size definition (%zu != %zu)",
923 collection_size, component->fields[i].collection_size);
924 return i;
925 }
926 }
927
928 size_t field_idx = component->num_fields;
929 component->num_fields++;
930
931 component->fields =
932 xrealloc(
933 component->fields,
934 component->num_fields *
935 sizeof(*(component->fields)));
936 struct yac_couple_config_field * field =
937 component->fields + field_idx;
938
939 field->name = strdup(name);
940 field->grid_idx = grid_idx;
941 field->timestep = timestep?strdup(timestep):NULL;
942 field->metadata = NULL;
945 return field_idx;
946}
947
949 struct yac_couple_config * couple_config, const char* component_name,
950 const char* grid_name, const char* name, char const * timestep,
951 size_t collection_size) {
952
954 couple_config,
955 yac_couple_config_get_component_idx(couple_config, component_name),
956 yac_couple_config_get_grid_idx(couple_config, grid_name),
958}
959
961 struct yac_couple_config * couple_config) {
962
963 return couple_config->num_couples;
964}
965
967 struct yac_couple_config * couple_config) {
968
969 switch(couple_config->missing_definition_is_fatal) {
970 case(FLAG_FALSE): return 0;
971 case(FLAG_TRUE): return 1;
972 default:
974 }
975}
976
978 struct yac_couple_config * couple_config,
979 int missing_definition_is_fatal) {
980
981 enum flag_value user_value =
982 (missing_definition_is_fatal != 0)?FLAG_TRUE:FLAG_FALSE;
983
985 (couple_config->missing_definition_is_fatal == FLAG_UNSET) ||
986 (couple_config->missing_definition_is_fatal == user_value),
987 "ERROR(yac_couple_config_set_missing_definition_is_fatal): "
988 "inconsistent setting of \"missing_definition_is_fatal\"-flag "
989 "(old: %s new: %s)",
990 (couple_config->missing_definition_is_fatal == FLAG_TRUE)?"TRUE":"FALSE",
991 (user_value == FLAG_TRUE)?"TRUE":"FALSE");
992
993 couple_config->missing_definition_is_fatal = user_value;
994}
995
997 struct yac_couple_config * couple_config, size_t couple_idx,
998 char const * routine_name, int line) {
999
1001 couple_idx < couple_config->num_couples,
1002 "ERROR(%s:%d:%s): invalid couple_idx", __FILE__, line, routine_name);
1003}
1004
1006 struct yac_couple_config * couple_config, size_t couple_idx) {
1007
1009 couple_config, couple_idx, "yac_couple_config_get_num_couple_fields",
1010 __LINE__);
1011
1012 return couple_config->couples[couple_idx].num_field_couples;
1013}
1014
1016 struct yac_couple_config * couple_config, size_t couple_idx,
1017 char const * couple_component_names[2]) {
1018
1020 couple_config, couple_idx, "yac_couple_config_get_couple_component_names",
1021 __LINE__);
1022
1023 for (int i = 0; i < 2; ++i)
1024 couple_component_names[i] =
1025 couple_config->components[
1026 couple_config->couples[couple_idx].component_indices[i]].name;
1027}
1028
1030 struct yac_couple_config * couple_config, char const * component_name) {
1031
1032 YAC_ASSERT(
1033 component_name,
1034 "ERROR(yac_couple_config_component_name_is_valid): component name is NULL")
1035 YAC_ASSERT(
1036 strlen(component_name) <= YAC_MAX_CHARLEN,
1037 "ERROR(yac_couple_config_component_name_is_valid): "
1038 "component name is too long (maximum is YAC_MAX_CHARLEN)")
1039
1040 for (size_t component_idx = 0; component_idx < couple_config->num_components;
1041 ++component_idx)
1042 if (!strcmp(component_name, couple_config->components[component_idx].name))
1043 return 1;
1044 return 0;
1045}
1046
1048 struct yac_couple_config * couple_config) {
1049
1050 return couple_config->num_components;
1051}
1052
1054 struct yac_couple_config * couple_config) {
1055
1056 return couple_config->num_grids;
1057}
1058
1060 struct yac_couple_config * couple_config, size_t component_idx) {
1061 check_component_idx(couple_config, component_idx,
1062 "yac_couple_config_get_num_fields", __LINE__);
1063 return couple_config->components[component_idx].num_fields;
1064}
1065
1067 struct yac_couple_config * couple_config, char const * component_name) {
1068
1069 size_t component_idx = SIZE_MAX;
1070 for (size_t i = 0; (i < couple_config->num_components) &&
1071 (component_idx == SIZE_MAX); ++i)
1072 if (!strcmp(couple_config->components[i].name, component_name))
1073 component_idx = i;
1074
1076 component_idx != SIZE_MAX,
1077 "ERROR(yac_couple_config_get_component_idx): "
1078 "Component \"%s\" not found in coupling configuration",
1079 component_name);
1080
1081 return component_idx;
1082}
1083
1085 struct yac_couple_config * couple_config, char const * grid_name) {
1086
1087 size_t grid_idx = SIZE_MAX;
1088 for (size_t i = 0;
1089 (i < couple_config->num_grids) && (grid_idx == SIZE_MAX); ++i)
1090 if (!strcmp(couple_config->grids[i].name, grid_name))
1091 grid_idx = i;
1092
1094 grid_idx != SIZE_MAX,
1095 "ERROR(yac_couple_config_get_grid_idx): "
1096 "grid name \"%s\" not in list of grids", grid_name)
1097
1098 return grid_idx;
1099}
1100
1102 struct yac_couple_config * couple_config, size_t component_idx,
1103 size_t grid_idx, char const * field_name) {
1105 couple_config, component_idx,
1106 "yac_couple_config_get_component_name", __LINE__);
1107
1108 size_t field_idx = SIZE_MAX;
1109 struct yac_couple_config_component * component =
1110 couple_config->components + component_idx;
1111 size_t nbr_fields = component->num_fields;
1112 for(size_t i=0;
1113 (i<nbr_fields) && (field_idx == SIZE_MAX); ++i)
1114 if((component->fields[i].grid_idx == grid_idx) &&
1115 !strcmp(
1116 field_name,
1117 component->fields[i].name))
1118 field_idx = i;
1119
1121 field_idx != SIZE_MAX,
1122 "ERROR(yac_couple_config_get_field_idx): "
1123 "field not found "
1124 "(component_idx %zu grid_idx %zu field_name \"%s\"",
1125 component_idx, grid_idx, field_name);
1126
1127 return field_idx;
1128}
1129
1131 struct yac_couple_config * couple_config, size_t component_idx) {
1132
1134 couple_config, component_idx,
1135 "yac_couple_config_get_component_name", __LINE__);
1136
1137 return couple_config->components[component_idx].name;
1138}
1139
1141 struct yac_couple_config * couple_config, size_t component_idx,
1142 size_t field_idx, char const * routine_name, int line) {
1143
1145 couple_config, component_idx, routine_name, __LINE__);
1146
1148 field_idx <
1149 couple_config->components[component_idx].num_fields,
1150 "ERROR(%s:%d:%s): invalid field_idx", __FILE__,
1151 line, routine_name)
1152}
1153
1155 struct yac_couple_config * couple_config, size_t component_idx,
1156 size_t field_idx) {
1157
1159 couple_config, component_idx, field_idx,
1160 "yac_couple_config_get_field_grid_name", __LINE__);
1161
1162 return
1163 couple_config->grids[
1164 couple_config->
1165 components[component_idx].
1166 fields[field_idx].grid_idx].name;
1167}
1168
1170 struct yac_couple_config * couple_config, size_t component_idx,
1171 size_t field_idx) {
1172
1174 couple_config, component_idx, field_idx,
1175 "yac_couple_config_get_field_name", __LINE__);
1176
1177 return
1178 couple_config->
1179 components[component_idx].
1180 fields[field_idx].name;
1181}
1182
1184 struct yac_couple_config * couple_config,
1185 char const * component_name, char const * grid_name,
1186 char const * field_name) {
1187
1188 size_t component_idx =
1189 yac_couple_config_get_component_idx(couple_config, component_name);
1190 size_t grid_idx =
1191 yac_couple_config_get_grid_idx(couple_config, grid_name);
1192 size_t field_idx =
1194 couple_config, component_idx, grid_idx, field_name);
1195
1196 struct yac_couple_config_field * field =
1197 couple_config->components[component_idx].fields +
1198 field_idx;
1199
1200 return field->timestep;
1201}
1202
1204 struct yac_couple_config * couple_config,
1205 char const * component_name, char const * grid_name,
1206 char const * field_name) {
1207
1208 size_t component_idx =
1209 yac_couple_config_get_component_idx(couple_config, component_name);
1210 size_t grid_idx =
1211 yac_couple_config_get_grid_idx(couple_config, grid_name);
1212 size_t field_idx =
1214 couple_config, component_idx, grid_idx, field_name);
1215
1216 size_t nbr_couples = couple_config->num_couples;
1217 for(size_t couple_idx = 0; couple_idx<nbr_couples; ++couple_idx){
1218 struct yac_couple_config_couple * couple = couple_config->couples + couple_idx;
1219 if(couple->component_indices[0] != component_idx &&
1220 couple->component_indices[1] != component_idx)
1221 continue;
1222 size_t nbr_trans_couples = couple->num_field_couples;
1223 for(size_t trans_couple_idx = 0; trans_couple_idx < nbr_trans_couples;
1224 ++trans_couple_idx){
1225 struct yac_couple_config_field_couple * transcouple =
1226 couple->field_couples + trans_couple_idx;
1227 if(transcouple->source.component_idx == component_idx &&
1228 transcouple->source.field_idx == field_idx)
1230 if(transcouple->target.component_idx == component_idx &&
1231 transcouple->target.field_idx == field_idx)
1233 }
1234 }
1236}
1237
1239 struct yac_couple_config * couple_config,
1240 size_t component_idx, size_t field_idx) {
1241
1243 couple_config, component_idx,
1244 "yac_couple_config_field_is_valid", __LINE__);
1246 couple_config, component_idx, field_idx,
1247 "yac_couple_config_field_is_valid", __LINE__);
1248
1249 struct yac_couple_config_field * field =
1250 couple_config->components[component_idx].fields +
1251 field_idx;
1252
1253 return (field->collection_size != SIZE_MAX) &&
1254 (field->timestep != NULL);
1255}
1256
1258 struct yac_couple_config * couple_config, size_t couple_idx,
1259 size_t field_couple_idx, char const * routine_name, int line) {
1260
1261 check_couple_idx(couple_config, couple_idx, routine_name, __LINE__);
1262
1264 field_couple_idx <
1265 couple_config->couples[couple_idx].num_field_couples,
1266 "ERROR(%s:%d:%s): invalid field_couple_idx",
1267 __FILE__, line, routine_name)
1268}
1269
1271 struct yac_couple_config * couple_config,
1272 size_t couple_idx, size_t field_couple_idx) {
1273
1275 couple_config, couple_idx, field_couple_idx,
1276 "yac_couple_config_get_interp_stack", __LINE__);
1277
1278 return
1279 couple_config->
1280 couples[couple_idx].
1281 field_couples[field_couple_idx].interp_stack;
1282}
1283
1285 struct yac_couple_config * couple_config,
1286 size_t couple_idx, size_t field_couple_idx,
1287 char const ** src_grid_name, char const ** tgt_grid_name) {
1288
1290 couple_config, couple_idx, field_couple_idx,
1291 "yac_couple_config_get_field_grid_names", __LINE__);
1292
1293 size_t src_component_idx =
1294 couple_config->couples[couple_idx].
1295 field_couples[field_couple_idx].
1296 source.component_idx;
1297 size_t src_field_idx =
1298 couple_config->couples[couple_idx].
1299 field_couples[field_couple_idx].
1300 source.field_idx;
1301
1302 size_t tgt_component_idx =
1303 couple_config->couples[couple_idx].
1304 field_couples[field_couple_idx].
1305 target.component_idx;
1306 size_t tgt_field_idx =
1307 couple_config->couples[couple_idx].
1308 field_couples[field_couple_idx].
1309 target.field_idx;
1310
1311 *src_grid_name =
1312 couple_config->grids[
1313 couple_config->components[src_component_idx].
1314 fields[src_field_idx].grid_idx].name;
1315 *tgt_grid_name =
1316 couple_config->grids[
1317 couple_config->components[tgt_component_idx].
1318 fields[tgt_field_idx].grid_idx].name;
1319}
1320
1322 struct yac_couple_config * couple_config,
1323 char const * comp_name, char const * grid_name, char const * field_name,
1324 double frac_mask_fallback_value) {
1325
1327 YAC_FRAC_MASK_VALUE_IS_VALID(frac_mask_fallback_value),
1328 "ERROR(yac_couple_config_field_enable_frac_mask): "
1329 "\"%lf\" is not a valid fractional mask fallback value "
1330 "(component: \"%s\" grid: \"%s\" field \"%s\")",
1331 frac_mask_fallback_value, comp_name, grid_name, field_name);
1332
1333 size_t comp_idx = yac_couple_config_get_component_idx(couple_config, comp_name);
1334 size_t grid_idx = yac_couple_config_get_grid_idx(couple_config, grid_name);
1335 size_t field_idx = yac_couple_config_get_field_idx(couple_config,
1336 comp_idx, grid_idx, field_name);
1337
1338 double old_frac_mask_fallback_value =
1339 couple_config->components[comp_idx].fields[field_idx].
1340 frac_mask_fallback_value;
1341
1342 // (use memcmp for comparison, because it can be nan)
1344 !YAC_FRAC_MASK_VALUE_IS_VALID(old_frac_mask_fallback_value) ||
1345 !memcmp(
1346 &old_frac_mask_fallback_value, &frac_mask_fallback_value,
1347 sizeof(frac_mask_fallback_value)),
1348 "ERROR(yac_couple_config_field_enable_frac_mask): "
1349 "fractional mask fallback value was already set:\n"
1350 "\told value \"%lf\" new value \"%lf\" "
1351 "(component: \"%s\" grid: \"%s\" field \"%s\")",
1352 old_frac_mask_fallback_value, frac_mask_fallback_value,
1353 comp_name, grid_name, field_name);
1354
1355 couple_config->components[comp_idx].fields[field_idx].
1356 frac_mask_fallback_value = frac_mask_fallback_value;
1357}
1358
1360 struct yac_couple_config * couple_config,
1361 char const * component_name, char const * grid_name,
1362 char const * field_name) {
1363
1364 size_t component_idx =
1365 yac_couple_config_get_component_idx(couple_config, component_name);
1366 size_t grid_idx =
1367 yac_couple_config_get_grid_idx(couple_config, grid_name);
1368 size_t field_idx =
1370 couple_config, component_idx, grid_idx, field_name);
1371
1372 struct yac_couple_config_field * field =
1373 couple_config->components[component_idx].fields +
1374 field_idx;
1375
1378 "ERROR(yac_couple_config_get_frac_mask_fallback_value): "
1379 "no valid fractional mask fallback value defined "
1380 "(component: \"%s\" grid: \"%s\" field \"%s\")",
1381 couple_config->components[component_idx].name, grid_name, field->name);
1382
1383 return field->frac_mask_fallback_value;
1384}
1385
1387 struct yac_couple_config * couple_config,
1388 char const * component_name, char const * grid_name,
1389 char const * field_name) {
1390
1391 size_t component_idx =
1392 yac_couple_config_get_component_idx(couple_config, component_name);
1393 size_t grid_idx =
1394 yac_couple_config_get_grid_idx(couple_config, grid_name);
1395 size_t field_idx =
1397 couple_config, component_idx, grid_idx, field_name);
1398
1399 struct yac_couple_config_field * field =
1400 couple_config->components[component_idx].fields +
1401 field_idx;
1402
1403 return field->collection_size;
1404}
1405
1407 struct yac_couple_config * couple_config,
1408 size_t couple_idx, size_t field_couple_idx,
1409 char const ** src_component_name, char const ** tgt_component_name) {
1410
1412 couple_config, couple_idx, field_couple_idx,
1413 "yac_couple_config_get_field_couple_component_names", __LINE__);
1414
1415 *src_component_name =
1416 couple_config->components[
1417 couple_config->couples[couple_idx].
1418 field_couples[field_couple_idx].source.component_idx].name;
1419 *tgt_component_name =
1420 couple_config->components[
1421 couple_config->couples[couple_idx].
1422 field_couples[field_couple_idx].target.component_idx].name;
1423}
1424
1426 struct yac_couple_config * couple_config,
1427 size_t couple_idx, size_t field_couple_idx,
1428 char const ** src_field_name, const char ** tgt_field_name) {
1429
1431 couple_config, couple_idx, field_couple_idx,
1432 "yac_couple_config_get_field_names", __LINE__);
1433
1434 size_t src_component_idx =
1435 couple_config->couples[couple_idx].
1436 field_couples[field_couple_idx].source.component_idx;
1437 size_t tgt_component_idx =
1438 couple_config->couples[couple_idx].
1439 field_couples[field_couple_idx].target.component_idx;
1440 size_t src_field_idx =
1441 couple_config->couples[couple_idx].
1442 field_couples[field_couple_idx].source.field_idx;
1443 size_t tgt_field_idx =
1444 couple_config->couples[couple_idx].
1445 field_couples[field_couple_idx].target.field_idx;
1446
1447 *src_field_name =
1448 couple_config->components[src_component_idx].
1449 fields[src_field_idx].name;
1450 *tgt_field_name =
1451 couple_config->components[tgt_component_idx].
1452 fields[tgt_field_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_mapping_on_source", __LINE__);
1462
1463 return
1464 couple_config->
1465 couples[couple_idx].
1466 field_couples[field_couple_idx].
1467 mapping_on_source;
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_source_lag", __LINE__);
1477
1478 return
1479 couple_config->
1480 couples[couple_idx].
1481 field_couples[field_couple_idx].
1482 source.lag;
1483}
1484
1486 struct yac_couple_config * couple_config,
1487 size_t couple_idx, size_t field_couple_idx) {
1488
1490 couple_config, couple_idx, field_couple_idx,
1491 "yac_couple_config_get_target_lag", __LINE__);
1492
1493 return
1494 couple_config->
1495 couples[couple_idx].
1496 field_couples[field_couple_idx].
1497 target.lag;
1498}
1499
1501 struct yac_couple_config * couple_config,
1502 size_t couple_idx, size_t field_couple_idx) {
1503
1505 couple_config, couple_idx, field_couple_idx,
1506 "yac_couple_config_get_coupling_period", __LINE__);
1507
1508 return
1509 couple_config->
1510 couples[couple_idx].
1511 field_couples[field_couple_idx].
1512 coupling_period;
1513}
1514
1516 struct yac_couple_config * couple_config,
1517 size_t couple_idx, size_t field_couple_idx) {
1518
1520 couple_config, couple_idx, field_couple_idx,
1521 "yac_couple_config_get_source_timestep", __LINE__);
1522
1523 struct yac_couple_config_field_couple * tcouple =
1524 couple_config->couples[couple_idx].field_couples + field_couple_idx;
1525 char const * timestep =
1526 couple_config->components[tcouple->source.component_idx].
1527 fields[tcouple->source.field_idx].timestep;
1528
1530 timestep,
1531 "ERROR(yac_couple_config_get_source_timestep): "
1532 "no valid timestep defined (component: \"%s\" field \"%s\")",
1533 couple_config->components[tcouple->source.component_idx].name,
1534 couple_config->components[tcouple->source.component_idx].
1535 fields[tcouple->source.field_idx].name);
1536
1537 return timestep;
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_target_timestep", __LINE__);
1547
1548 struct yac_couple_config_field_couple * tcouple =
1549 couple_config->couples[couple_idx].field_couples + field_couple_idx;
1550 char const * timestep =
1551 couple_config->components[tcouple->target.component_idx].
1552 fields[tcouple->target.field_idx].timestep;
1553
1554 return timestep;
1555}
1556
1558 struct yac_couple_config * couple_config,
1559 size_t couple_idx, size_t field_couple_idx) {
1560
1562 couple_config, couple_idx, field_couple_idx,
1563 "yac_couple_config_get_coupling_period_operation", __LINE__);
1564
1565 return
1566 couple_config->
1567 couples[couple_idx].
1568 field_couples[field_couple_idx].
1570}
1571
1572static void set_datetime(
1573 char const * type_datetime, struct _datetime ** old, char const * str_new) {
1574
1575 if ((str_new == NULL) || (strlen(str_new) == 0)) return;
1576
1577 YAC_ASSERT(
1579 "ERROR(set_datetime): calendar has not yet been set");
1580
1581 struct _datetime * new = newDateTime(str_new);
1582
1584 new != NULL,
1585 "ERROR(set_datetime): failed to parse datetime \"%s\"", str_new);
1586
1587 char old_datetime_buffer[MAX_DATETIME_STR_LEN],
1588 new_datetime_buffer[MAX_DATETIME_STR_LEN];
1589
1591 (*old == NULL) || (equal_to == compareDatetime(*old, new)),
1592 "ERROR(set_datetime): inconsistent %s datetime "
1593 "(old: \"%s\" new: \"%s\")", type_datetime,
1594 datetimeToString(*old, old_datetime_buffer),
1595 datetimeToString(new, new_datetime_buffer));
1596
1597 deallocateDateTime(*old);
1598
1599 *old = new;
1600}
1601
1603 struct yac_couple_config * couple_config,
1604 char const * start, char const * end) {
1605
1606 set_datetime("start", &(couple_config->start_datetime), start);
1607 set_datetime("end", &(couple_config->end_datetime), end);
1608}
1609
1611 struct yac_couple_config * couple_config) {
1612
1613 YAC_ASSERT(couple_config->start_datetime,
1614 "ERROR(yac_couple_config_get_start_datetime): "
1615 "start_datetime not yet defined");
1616
1617 char datetime_buffer[MAX_DATETIME_STR_LEN];
1618
1619 return
1620 strdup(datetimeToString(couple_config->start_datetime, datetime_buffer));
1621}
1622
1624 struct yac_couple_config * couple_config) {
1625
1626 YAC_ASSERT(couple_config->end_datetime,
1627 "ERROR(yac_couple_config_get_start_datetime): "
1628 "start_datetime not yet defined");
1629
1630 char datetime_buffer[MAX_DATETIME_STR_LEN];
1631
1632 return
1633 strdup(datetimeToString(couple_config->end_datetime, datetime_buffer));
1634}
1635
1637 struct yac_couple_config * couple_config, size_t grid_idx) {
1638
1640 grid_idx < couple_config->num_grids,
1641 "ERROR(yac_couple_config_get_grid_name): "
1642 "Invalid grid idx %zu", grid_idx);
1643
1644 return couple_config->grids[grid_idx].name;
1645}
1646
1648 struct yac_couple_config * couple_config,
1649 size_t couple_idx, size_t field_couple_idx) {
1650
1652 couple_config, couple_idx, field_couple_idx,
1653 "yac_couple_config_enforce_write_weight_file", __LINE__);
1654
1655 return
1656 couple_config->
1657 couples[couple_idx].
1658 field_couples[field_couple_idx].
1659 enforce_write_weight_file;
1660}
1661
1663 struct yac_couple_config * couple_config,
1664 size_t couple_idx, size_t field_couple_idx) {
1665
1667 couple_config, couple_idx, field_couple_idx,
1668 "yac_couple_config_get_weight_file_name", __LINE__);
1669
1670 static char dummy[] = "\0";
1671 char const * weight_file_name =
1672 couple_config->
1673 couples[couple_idx].
1674 field_couples[field_couple_idx].
1675 weight_file_name;
1676
1677 return (weight_file_name != NULL)?weight_file_name:dummy;
1678}
1679
1681 struct yac_couple_config * couple_config,
1682 size_t couple_idx, size_t field_couple_idx) {
1683
1685 couple_config, couple_idx, field_couple_idx,
1686 "yac_couple_config_get_scale_factor", __LINE__);
1687
1688 return
1689 couple_config->
1690 couples[couple_idx].
1691 field_couples[field_couple_idx].
1692 scale_factor;
1693}
1694
1696 struct yac_couple_config * couple_config,
1697 size_t couple_idx, size_t field_couple_idx) {
1698
1700 couple_config, couple_idx, field_couple_idx,
1701 "yac_couple_config_get_scale_summand", __LINE__);
1702
1703 return
1704 couple_config->
1705 couples[couple_idx].
1706 field_couples[field_couple_idx].
1707 scale_summand;
1708}
1709
1711 struct yac_couple_config * couple_config,
1712 size_t couple_idx, size_t field_couple_idx,
1713 char const * const ** mask_names, size_t * num_mask_names) {
1714
1716 couple_config, couple_idx, field_couple_idx,
1717 "yac_couple_config_get_src_mask_names", __LINE__);
1718
1719 *mask_names =
1720 (char const * const *)(
1721 couple_config->
1722 couples[couple_idx].
1723 field_couples[field_couple_idx].
1724 src_mask_names);
1725 *num_mask_names =
1726 couple_config->
1727 couples[couple_idx].
1728 field_couples[field_couple_idx].
1729 num_src_mask_names;
1730}
1731
1733 struct yac_couple_config * couple_config,
1734 size_t couple_idx, size_t field_couple_idx) {
1735
1737 couple_config, couple_idx, field_couple_idx,
1738 "yac_couple_config_get_tgt_mask_name", __LINE__);
1739
1740 return
1741 couple_config->
1742 couples[couple_idx].
1743 field_couples[field_couple_idx].
1744 tgt_mask_name;
1745}
1746
1748 struct yac_couple_config * couple_config,
1749 size_t couple_idx, size_t field_couple_idx) {
1750
1752 couple_config, couple_idx, field_couple_idx,
1753 "yac_couple_config_get_yaxt_exchanger_name", __LINE__);
1754
1755 return
1756 couple_config->
1757 couples[couple_idx].
1758 field_couples[field_couple_idx].
1759 yaxt_exchanger_name;
1760}
1761
1763 struct yac_couple_config * couple_config, char const * grid_name) {
1764
1765 for (size_t grid_idx = 0; grid_idx < couple_config->num_grids;
1766 ++grid_idx)
1767 if (!strcmp(couple_config->grids[grid_idx].name, grid_name))
1768 return 1;
1769 return 0;
1770}
1771
1773 char const * string, MPI_Comm comm) {
1774
1775 int strlen_pack_size, string_pack_size;
1776 yac_mpi_call(MPI_Pack_size(1, MPI_INT, comm, &strlen_pack_size), comm);
1777
1778 if (string != NULL) {
1780 MPI_Pack_size(
1781 (int)(strlen(string)), MPI_CHAR, comm, &string_pack_size), comm);
1782 } else {
1783 string_pack_size = 0;
1784 }
1785
1786 return (size_t)strlen_pack_size + (size_t)string_pack_size;
1787}
1788
1790 void * grid_, MPI_Comm comm) {
1791
1792 struct yac_couple_config_grid * grid = grid_;
1793
1794 return yac_couple_config_get_string_pack_size(grid->name, comm);
1795}
1796
1798 void * field_, MPI_Comm comm) {
1799
1800 struct yac_couple_config_field * field =
1801 (struct yac_couple_config_field *)field_;
1802
1803 YAC_ASSERT(
1804 field->grid_idx <= INT_MAX,
1805 "ERROR(yac_couple_config_get_field_pack_size):"
1806 "grid_idx is too big")
1807
1808 int int_pack_size;
1809 yac_mpi_call(MPI_Pack_size(1, MPI_INT, comm, &int_pack_size), comm);
1810 size_t name_pack_size = yac_couple_config_get_string_pack_size(
1811 field->name, comm);
1812
1813 return int_pack_size + // grid_idx
1814 name_pack_size;
1815}
1816
1818 void * component_, MPI_Comm comm) {
1819
1820 struct yac_couple_config_component * component =
1821 (struct yac_couple_config_component *)component_;
1822
1823 return
1825}
1826
1828 void * field_couple_, MPI_Comm comm) {
1829
1830 struct yac_couple_config_field_couple * field_couple =
1831 (struct yac_couple_config_field_couple *)field_couple_;
1832
1833 int ints_pack_size;
1834 yac_mpi_call(MPI_Pack_size(10, MPI_INT, comm, &ints_pack_size), comm);
1835 int doubles_pack_size;
1836 yac_mpi_call(MPI_Pack_size(2, MPI_DOUBLE, comm, &doubles_pack_size), comm);
1837 int src_mask_names_pack_size = 0;
1838 if (field_couple->num_src_mask_names > 0)
1839 for (size_t i = 0; i < field_couple->num_src_mask_names; ++i)
1840 src_mask_names_pack_size +=
1842 field_couple->src_mask_names[i], comm);
1843
1844 return
1845 (size_t)ints_pack_size + // source.component_idx
1846 // source.field_idx
1847 // source.lag
1848 // target.component_idx
1849 // target.field_idx
1850 // target.lag
1851 // mapping_on_source
1852 // coupling_period_operation
1853 // enforce_write_weight_file
1854 // num_src_mask_names
1856 field_couple->interp_stack, comm) +
1858 field_couple->coupling_period, comm) +
1860 field_couple->weight_file_name, comm) +
1861 doubles_pack_size + // scale_factor
1862 // scale_summand
1863 src_mask_names_pack_size +
1865 field_couple->tgt_mask_name, comm) +
1867 field_couple->yaxt_exchanger_name, comm);
1868}
1869
1871 void * couple_, MPI_Comm comm) {
1872
1873 UNUSED(couple_);
1874
1875 int component_indices_pack_size;
1877 MPI_Pack_size(2, MPI_INT, comm, &component_indices_pack_size), comm);
1878
1879 return (size_t)component_indices_pack_size;
1880}
1881
1883 void * config_output_, MPI_Comm comm) {
1884
1885 struct yac_couple_config_config_output * config_output = config_output_;
1886
1887 int filetype_pack_size;
1889 MPI_Pack_size(1, MPI_INT, comm, &filetype_pack_size), comm);
1890
1891 return
1892 filetype_pack_size +
1893 yac_couple_config_get_string_pack_size(config_output->name, comm) +
1894 yac_couple_config_get_string_pack_size(config_output->ref, comm);
1895}
1896
1898 char const * string, void * buffer, int buffer_size, int * position,
1899 MPI_Comm comm) {
1900
1901 size_t len = (string == NULL)?0:strlen(string);
1902
1903 YAC_ASSERT(
1904 len <= INT_MAX, "ERROR(yac_couple_config_pack_string): string too long")
1905
1906 int len_int = (int)len;
1907
1909 MPI_Pack(
1910 &len_int, 1, MPI_INT, buffer, buffer_size, position, comm), comm);
1911
1912 if (len > 0)
1914 MPI_Pack(
1915 string, len_int, MPI_CHAR, buffer, buffer_size, position, comm),
1916 comm);
1917}
1918
1920 void * grid_, void * buffer, int buffer_size, int * position, MPI_Comm comm) {
1921
1922 struct yac_couple_config_grid * grid =
1923 (struct yac_couple_config_grid *)grid_;
1925 grid->name, buffer, buffer_size, position, comm);
1926}
1927
1929 void * field_, void * buffer, int buffer_size, int * position, MPI_Comm comm) {
1930
1931 struct yac_couple_config_field * field =
1932 (struct yac_couple_config_field *)field_;
1933
1934 YAC_ASSERT(
1935 field->grid_idx <= INT_MAX,
1936 "ERROR(yac_couple_config_pack_field): grid_idx is too big")
1937
1938 int grid_idx = (int)field->grid_idx;
1939
1941 MPI_Pack(
1942 &grid_idx, 1, MPI_INT, buffer, buffer_size, position, comm), comm);
1943 yac_couple_config_pack_string(field->name, buffer,
1944 buffer_size, position, comm);
1945}
1946
1948 void * component_, void * buffer, int buffer_size, int * position, MPI_Comm comm) {
1949
1950 struct yac_couple_config_component * component =
1951 (struct yac_couple_config_component *)component_;
1952
1954 component->name, buffer, buffer_size, position, comm);
1955}
1956
1958 void * field_couple_, void * buffer, int buffer_size, int * position,
1959 MPI_Comm comm) {
1960
1961 struct yac_couple_config_field_couple * field_couple =
1962 (struct yac_couple_config_field_couple *)field_couple_;
1963
1964 YAC_ASSERT(
1965 field_couple->source.component_idx <= INT_MAX,
1966 "ERROR(yac_couple_config_pack_field_couple): "
1967 "source.component_idx bigger than INT_MAX")
1968 YAC_ASSERT(
1969 field_couple->source.field_idx <= INT_MAX,
1970 "ERROR(yac_couple_config_pack_field_couple): "
1971 "source.field_idx bigger than INT_MAX")
1972 YAC_ASSERT(
1973 field_couple->target.component_idx <= INT_MAX,
1974 "ERROR(yac_couple_config_pack_field_couple): "
1975 "target.component_idx bigger than INT_MAX")
1976 YAC_ASSERT(
1977 field_couple->target.field_idx <= INT_MAX,
1978 "ERROR(yac_couple_config_pack_field_couple): "
1979 "target.field_idx bigger than INT_MAX")
1980 YAC_ASSERT(
1981 field_couple->mapping_on_source <= INT_MAX,
1982 "ERROR(yac_couple_config_pack_field_couple): "
1983 "mapping_on_source bigger than INT_MAX")
1984 YAC_ASSERT(
1985 field_couple->coupling_period_operation <= INT_MAX,
1986 "ERROR(yac_couple_config_pack_field_couple): "
1987 "coupling_period_operation bigger than INT_MAX")
1988 YAC_ASSERT(
1989 field_couple->enforce_write_weight_file <= INT_MAX,
1990 "ERROR(yac_couple_config_pack_field_couple): "
1991 "enforce_write_weight_file bigger than INT_MAX")
1992 YAC_ASSERT(
1993 field_couple->num_src_mask_names <= INT_MAX,
1994 "ERROR(yac_couple_config_pack_field_couple): "
1995 "num_src_mask_names bigger than INT_MAX")
1996
1997 int ints[10] = {
1998 field_couple->source.component_idx,
1999 field_couple->source.field_idx,
2000 field_couple->source.lag,
2001 field_couple->target.component_idx,
2002 field_couple->target.field_idx,
2003 field_couple->target.lag,
2004 field_couple->mapping_on_source,
2005 field_couple->coupling_period_operation,
2006 field_couple->enforce_write_weight_file,
2007 field_couple->num_src_mask_names};
2008
2010 MPI_Pack(ints, 10, MPI_INT, buffer, buffer_size, position, comm), comm);
2011
2013 field_couple->coupling_period, buffer, buffer_size, position, comm);
2015 field_couple->weight_file_name, buffer, buffer_size, position, comm);
2016
2017 double doubles[2] = {
2018 field_couple->scale_factor,
2019 field_couple->scale_summand};
2020
2022 MPI_Pack(doubles, 2, MPI_DOUBLE, buffer, buffer_size, position, comm),
2023 comm);
2024
2026 field_couple->interp_stack, buffer, buffer_size, position, comm);
2027
2028 if (field_couple->num_src_mask_names > 0)
2029 for (size_t i = 0; i < field_couple->num_src_mask_names; ++i)
2031 field_couple->src_mask_names[i], buffer, buffer_size, position, comm);
2032
2034 field_couple->tgt_mask_name, buffer, buffer_size, position, comm);
2035
2037 field_couple->yaxt_exchanger_name, buffer, buffer_size, position, comm);
2038}
2039
2041 void * couple_, void * buffer, int buffer_size, int * position, MPI_Comm comm) {
2042
2043 struct yac_couple_config_couple * couple =
2044 (struct yac_couple_config_couple *)couple_;
2045
2046 YAC_ASSERT(
2047 (couple->component_indices[0] <= INT_MAX) &&
2048 (couple->component_indices[1] <= INT_MAX),
2049 "ERROR(yac_couple_config_pack_couple): "
2050 "component_indices bigger than INT_MAX")
2051
2052 int component_indices[2] = {(int)(couple->component_indices[0]),
2053 (int)(couple->component_indices[1])};
2054
2056 MPI_Pack(
2057 component_indices, 2, MPI_INT, buffer, buffer_size, position, comm),
2058 comm);
2059}
2060
2062 void * config_output_, void * buffer, int buffer_size, int * position,
2063 MPI_Comm comm) {
2064
2065 struct yac_couple_config_config_output * config_output =
2066 (struct yac_couple_config_config_output *)config_output_;
2067
2069 config_output->name, buffer, buffer_size, position, comm);
2070 int filetype_int = (int)config_output->type;
2072 MPI_Pack(
2073 &filetype_int, 1, MPI_INT, buffer, buffer_size, position, comm), comm);
2075 config_output->ref, buffer, buffer_size, position, comm);
2076}
2077
2079 void * buffer, int buffer_size, int * position, MPI_Comm comm) {
2080
2081 int string_len;
2083 MPI_Unpack(
2084 buffer, buffer_size, position, &string_len, 1, MPI_INT, comm), comm);
2085
2086 if (string_len <= 0) return NULL;
2087
2088 char * string = xmalloc(((size_t)string_len + 1) * sizeof(*string));
2090 MPI_Unpack(
2091 buffer, buffer_size, position, string, string_len, MPI_CHAR, comm), comm);
2092 string[string_len] = '\0';
2093 return string;
2094}
2095
2097 void * buffer, int buffer_size, int * position, void * grid_,
2098 MPI_Comm comm) {
2099
2100 struct yac_couple_config_grid * grid =
2101 (struct yac_couple_config_grid *)grid_;
2102
2103 grid->name =
2104 yac_couple_config_unpack_string(buffer, buffer_size, position, comm);
2105 grid->output_filename = NULL;
2106 grid->metadata = NULL;
2107}
2108
2110 void * buffer, int buffer_size, int * position, void * field_,
2111 MPI_Comm comm) {
2112
2113 struct yac_couple_config_field * field =
2114 (struct yac_couple_config_field *)field_;
2115
2116 int grid_idx_int;
2118 MPI_Unpack(
2119 buffer, buffer_size, position, &grid_idx_int, 1, MPI_INT, comm), comm);
2120
2121 YAC_ASSERT(
2122 grid_idx_int >= 0,
2123 "ERROR(yac_couple_config_unpack_field): invalid number of grid_idx_int")
2124
2125 field->grid_idx = (size_t)grid_idx_int;
2127 field->collection_size = SIZE_MAX;
2128 field->name = yac_couple_config_unpack_string(buffer,
2129 buffer_size, position, comm);
2130 field->timestep = NULL;
2131 field->metadata = NULL;
2132}
2133
2135 void * buffer, int buffer_size, int * position, void * component_,
2136 MPI_Comm comm) {
2137
2138 struct yac_couple_config_component * component =
2139 (struct yac_couple_config_component *)component_;
2140
2141 component->name =
2142 yac_couple_config_unpack_string(buffer, buffer_size, position, comm);
2143 component->metadata = NULL;
2144 component->num_fields = 0;
2145 component->fields = NULL;
2146}
2147
2149 void * buffer, int buffer_size, int * position, void * field_couple_,
2150 MPI_Comm comm) {
2151
2152 struct yac_couple_config_field_couple * field_couple =
2153 (struct yac_couple_config_field_couple *)field_couple_;
2154
2155 int ints[10];
2157 MPI_Unpack(
2158 buffer, buffer_size, position, ints, 10, MPI_INT, comm), comm);
2159
2160 YAC_ASSERT(
2161 ints[0] >= 0,
2162 "ERROR(yac_couple_config_unpack_field_couple): "
2163 "invalid source.component_idx")
2164 YAC_ASSERT(
2165 ints[1] >= 0,
2166 "ERROR(yac_couple_config_unpack_field_couple): "
2167 "invalid source.field_idx")
2168 YAC_ASSERT(
2169 ints[3] >= 0,
2170 "ERROR(yac_couple_config_unpack_field_couple): "
2171 "target.component_idx bigger than INT_MAX")
2172 YAC_ASSERT(
2173 ints[4] >= 0,
2174 "ERROR(yac_couple_config_unpack_field_couple): "
2175 "invalid target.field_idx")
2176 YAC_ASSERT(
2177 ints[6] >= 0,
2178 "ERROR(yac_couple_config_unpack_field_couple): "
2179 "invalid mapping_on_source")
2180 YAC_ASSERT(
2181 ints[7] >= 0,
2182 "ERROR(yac_couple_config_unpack_field_couple): "
2183 "invalid coupling_period_operation")
2184 YAC_ASSERT(
2185 ints[8] >= 0,
2186 "ERROR(yac_couple_config_unpack_field_couple): "
2187 "invalid enforce_write_weight_file")
2188 YAC_ASSERT(
2189 ints[9] >= 0,
2190 "ERROR(yac_couple_config_unpack_field_couple): "
2191 "invalid num_src_mask_names")
2192
2193 field_couple->source.component_idx = (size_t)(ints[0]);
2194 field_couple->source.field_idx = (size_t)(ints[1]);
2195 field_couple->source.lag = ints[2];
2196 field_couple->target.component_idx = (size_t)(ints[3]);
2197 field_couple->target.field_idx = (size_t)(ints[4]);
2198 field_couple->target.lag = ints[5];
2199 field_couple->mapping_on_source = ints[6];
2200 field_couple->coupling_period_operation =
2201 (enum yac_reduction_type)(ints[7]);
2202 field_couple->enforce_write_weight_file = ints[8];
2203 field_couple->num_src_mask_names = (size_t)(ints[9]);
2204
2205 field_couple->coupling_period =
2206 yac_couple_config_unpack_string(buffer, buffer_size, position, comm);
2207 field_couple->weight_file_name =
2208 yac_couple_config_unpack_string(buffer, buffer_size, position, comm);
2209
2210 double doubles[2];
2212 MPI_Unpack(
2213 buffer, buffer_size, position, doubles, 2, MPI_DOUBLE, comm), comm);
2214
2215 field_couple->scale_factor = doubles[0];
2216 field_couple->scale_summand = doubles[1];
2217
2218 field_couple->interp_stack =
2219 yac_interp_stack_config_unpack(buffer, buffer_size, position, comm);
2220
2221 if (field_couple->num_src_mask_names > 0) {
2222 field_couple->src_mask_names =
2223 xmalloc(
2224 field_couple->num_src_mask_names *
2225 sizeof(*(field_couple->src_mask_names)));
2226 for (size_t i = 0; i < field_couple->num_src_mask_names; ++i)
2227 field_couple->src_mask_names[i] =
2229 buffer, buffer_size, position, comm);
2230 } else {
2231 field_couple->src_mask_names = NULL;
2232 }
2233
2234 field_couple->tgt_mask_name =
2236 buffer, buffer_size, position, comm);
2237
2238 field_couple->yaxt_exchanger_name =
2240 buffer, buffer_size, position, comm);
2241}
2242
2244 void * buffer, int buffer_size, int * position, void * couple_,
2245 MPI_Comm comm) {
2246
2247 struct yac_couple_config_couple * couple =
2248 (struct yac_couple_config_couple *)couple_;
2249
2250 int component_indices[2];
2252 MPI_Unpack(
2253 buffer, buffer_size, position, component_indices, 2, MPI_INT, comm),
2254 comm);
2255
2256 YAC_ASSERT(
2257 (component_indices[0] >= 0) && (component_indices[1] >= 0),
2258 "ERROR(yac_couple_config_unpack_couple): invalid component indices")
2259
2260 couple->component_indices[0] = (size_t)(component_indices[0]);
2261 couple->component_indices[1] = (size_t)(component_indices[1]);
2262 couple->num_field_couples = 0;
2263 couple->field_couples = NULL;
2264}
2265
2267 void * buffer, int buffer_size, int * position, void * config_output_,
2268 MPI_Comm comm) {
2269
2270 struct yac_couple_config_config_output * config_output =
2271 (struct yac_couple_config_config_output *)config_output_;
2272
2273 config_output->name =
2274 yac_couple_config_unpack_string(buffer, buffer_size, position, comm);
2275 int filetype_int;
2277 MPI_Unpack(
2278 buffer, buffer_size, position, &filetype_int, 1, MPI_INT, comm),
2279 comm);
2280 config_output->type = (enum yac_text_filetype)filetype_int;
2281 config_output->ref =
2282 yac_couple_config_unpack_string(buffer, buffer_size, position, comm);
2283 config_output->include_definitions = FLAG_UNSET;
2284}
2285
2287 struct yac_couple_config * couple_config,
2288 char const * src_comp_name, char const * src_grid_name, char const * src_field_name,
2289 char const * tgt_comp_name, char const * tgt_grid_name, char const * tgt_field_name,
2290 char const * coupling_period, int time_reduction,
2291 struct yac_interp_stack_config * interp_stack,
2292 int src_lag, int tgt_lag,
2293 const char* weight_file_name, int mapping_on_source,
2294 double scale_factor, double scale_summand,
2295 size_t num_src_mask_names, char const * const * src_mask_names,
2296 char const * tgt_mask_name, char const * yaxt_exchanger_name) {
2297
2298 YAC_ASSERT(src_comp_name && src_comp_name[0] != '\0',
2299 "ERROR(yac_couple_config_def_couple): invalid parameter: src_comp_name");
2300 YAC_ASSERT(src_grid_name && src_grid_name[0] != '\0',
2301 "ERROR(yac_couple_config_def_couple): invalid parameter: src_grid_name");
2302 YAC_ASSERT(src_field_name && src_field_name[0] != '\0',
2303 "ERROR(yac_couple_config_def_couple): invalid parameter: src_field_name");
2304 YAC_ASSERT(tgt_comp_name && tgt_comp_name[0] != '\0',
2305 "ERROR(yac_couple_config_def_couple): invalid parameter: tgt_comp_name");
2306 YAC_ASSERT(tgt_grid_name && tgt_grid_name[0] != '\0',
2307 "ERROR(yac_couple_config_def_couple): invalid parameter: tgt_grid_name");
2308 YAC_ASSERT(tgt_field_name && tgt_field_name[0] != '\0',
2309 "ERROR(yac_couple_config_def_couple): invalid parameter: tgt_field_name");
2310 YAC_ASSERT(coupling_period && coupling_period[0] != '\0',
2311 "ERROR(yac_couple_config_def_couple): invalid parameter: coupling_period");
2312 YAC_ASSERT(
2313 (time_reduction == TIME_NONE) ||
2314 (time_reduction == TIME_ACCUMULATE) ||
2315 (time_reduction == TIME_AVERAGE) ||
2316 (time_reduction == TIME_MINIMUM) ||
2317 (time_reduction == TIME_MAXIMUM),
2318 "ERROR(yac_couple_config_def_couple): invalid parameter: time_reduction");
2319 YAC_ASSERT_F(isnormal(scale_factor),
2320 "ERROR(yac_couple_config_def_couple): \"%lf\" is not a valid scale factor",
2321 scale_factor);
2322 YAC_ASSERT_F(isnormal(scale_summand) || (scale_summand == 0.0),
2323 "ERROR(yac_couple_config_def_couple): \"%lf\" is not a valid scale summand",
2324 scale_summand);
2325
2326 // get component indices
2327 size_t src_comp_idx =
2328 yac_couple_config_add_component_(couple_config, src_comp_name);
2329 size_t tgt_comp_idx =
2330 yac_couple_config_add_component_(couple_config, tgt_comp_name);
2331 size_t src_grid_idx =
2332 yac_couple_config_add_grid_(couple_config, src_grid_name);
2333 size_t tgt_grid_idx =
2334 yac_couple_config_add_grid_(couple_config, tgt_grid_name);
2335
2336 // check if couple exists
2337 size_t component_indices[2];
2338 if(src_comp_idx < tgt_comp_idx){
2339 component_indices[0] = src_comp_idx;
2340 component_indices[1] = tgt_comp_idx;
2341 }else{
2342 component_indices[0] = tgt_comp_idx;
2343 component_indices[1] = src_comp_idx;
2344 }
2345 struct yac_couple_config_couple * couple = NULL;
2346 for(size_t i = 0; (i < couple_config->num_couples) && !couple; ++i)
2347 if(couple_config->couples[i].component_indices[0] == component_indices[0] &&
2348 couple_config->couples[i].component_indices[1] == component_indices[1])
2349 couple = &couple_config->couples[i];
2350
2351 // create if couple does not exists
2352 if(!couple){
2353 couple_config->couples =
2354 xrealloc(
2355 couple_config->couples,
2356 (couple_config->num_couples + 1) * sizeof(*couple_config->couples));
2357 couple = &couple_config->couples[couple_config->num_couples];
2358 couple_config->num_couples++;
2359 couple->component_indices[0] = component_indices[0];
2360 couple->component_indices[1] = component_indices[1];
2361 couple->num_field_couples = 0;
2362 couple->field_couples = NULL;
2363 }
2364
2365 // get field indices
2366 size_t src_field_idx =
2368 couple_config, src_comp_idx, src_grid_idx,
2369 src_field_name, NULL, SIZE_MAX);
2370 size_t tgt_field_idx =
2372 couple_config, tgt_comp_idx, tgt_grid_idx,
2373 tgt_field_name, NULL, SIZE_MAX);
2374
2375 struct yac_couple_config_field_couple field_couple;
2376 field_couple.source.component_idx = src_comp_idx;
2377 field_couple.source.field_idx = src_field_idx;
2378 field_couple.source.lag = src_lag;
2379 field_couple.target.component_idx = tgt_comp_idx;
2380 field_couple.target.field_idx = tgt_field_idx;
2381 field_couple.target.lag = tgt_lag;
2382 field_couple.coupling_period = strdup(coupling_period);
2383 field_couple.coupling_period_operation =
2384 (enum yac_reduction_type)time_reduction;
2386 field_couple.mapping_on_source = mapping_on_source;
2387 field_couple.weight_file_name =
2389 field_couple.scale_factor = scale_factor;
2390 field_couple.scale_summand = scale_summand;
2391 field_couple.enforce_write_weight_file = weight_file_name != NULL;
2393 if (num_src_mask_names > 0) {
2394 field_couple.src_mask_names =
2396 for (size_t i = 0; i < num_src_mask_names; ++i)
2397 field_couple.src_mask_names[i] = strdup(src_mask_names[i]);
2398 } else {
2399 field_couple.src_mask_names = NULL;
2400 }
2401 field_couple.tgt_mask_name = string_dup(tgt_mask_name);
2403
2404 // check if field_couple exists
2405 int found_flag = 0;
2406 for(size_t i = 0; (i < couple->num_field_couples) && !found_flag; ++i) {
2408 couple->field_couples + i, &field_couple)) {
2410 couple->field_couples + i, &field_couple, MPI_COMM_SELF);
2412 found_flag = 1;
2413 }
2414 }
2415
2416 if (!found_flag) {
2417 couple->field_couples = xrealloc(couple->field_couples,
2418 (couple->num_field_couples + 1) * sizeof(*couple->field_couples));
2419 couple->field_couples[couple->num_field_couples] = field_couple;
2420 couple->num_field_couples++;
2421 }
2422}
2423
2425 char const * type_datetime, struct _datetime ** datetime, MPI_Comm comm) {
2426
2427 // convert datetime to string
2428 char datetime_buffer[MAX_DATETIME_STR_LEN];
2429 char * str_datetime = datetimeToString(*datetime, datetime_buffer);
2430
2431 // copy from static to dynamic memory
2432 if (str_datetime != NULL) str_datetime = strdup(str_datetime);
2433
2434 // delete old datetime
2435 deallocateDateTime(*datetime);
2436
2437 // synchronize datetime string across processes
2438 couple_config_sync_string(type_datetime, &str_datetime, comm);
2439
2440 // convert datetime string to _datetime struct
2441 *datetime = newDateTime(str_datetime);
2442
2443 free(str_datetime);
2444}
2445
2447 struct yac_couple_config * couple_config, MPI_Comm comm) {
2448
2451 "start_datetime", &(couple_config->start_datetime), comm);
2453 "end_datetime", &(couple_config->end_datetime), comm);
2454}
2455
2499
2500static void merge_grids(
2501 struct yac_couple_config * couple_config, MPI_Comm comm) {
2502
2503 size_t* old_to_new_idx;
2504 void * p_grids = couple_config->grids;
2506 &couple_config->num_grids, &p_grids,
2507 sizeof(couple_config->grids[0]),
2508 comm, &dist_merge_vtable_grid, &old_to_new_idx);
2509 couple_config->grids = p_grids;
2510
2511 // set new grid_idx in fields
2512 for(size_t comp_idx = 0; comp_idx < couple_config->num_components;
2513 ++comp_idx) {
2514 struct yac_couple_config_component * component =
2515 couple_config->components + comp_idx;
2516 for(size_t field_idx = 0; field_idx < component->num_fields; ++field_idx)
2517 component->fields[field_idx].grid_idx =
2518 old_to_new_idx[component->fields[field_idx].grid_idx];
2519 }
2520
2521 free(old_to_new_idx);
2522}
2523
2524static void merge_fields(
2525 struct yac_couple_config * couple_config, size_t comp_idx, MPI_Comm comm) {
2526
2527 struct yac_couple_config_component * component =
2528 couple_config->components + comp_idx;
2529
2530 size_t* old_to_new_idx;
2531 void * p_fields = component->fields;
2533 &component->num_fields, &p_fields,
2534 sizeof(component->fields[0]),
2535 comm, &dist_merge_vtable_field, &old_to_new_idx);
2536 component->fields = p_fields;
2537
2538 // set new field_idx in all field_couples
2539 for(size_t couple_idx = 0; couple_idx < couple_config->num_couples;
2540 ++couple_idx) {
2541 struct yac_couple_config_couple * couple =
2542 couple_config->couples + couple_idx;
2543 for(size_t field_couple_idx = 0;
2544 field_couple_idx < couple->num_field_couples; ++field_couple_idx){
2545 struct yac_couple_config_field_couple * field_couple =
2546 couple->field_couples + field_couple_idx;
2547 if(field_couple->source.component_idx == comp_idx)
2548 field_couple->source.field_idx =
2549 old_to_new_idx[field_couple->source.field_idx];
2550 if(field_couple->target.component_idx == comp_idx)
2551 field_couple->target.field_idx =
2552 old_to_new_idx[field_couple->target.field_idx];
2553 }
2554 }
2555
2556 free(old_to_new_idx);
2557}
2558
2560 struct yac_couple_config * couple_config, MPI_Comm comm) {
2561
2562 // distribute and merge basic component information while keeping the
2563 // individual fields
2564 size_t* old_to_new_idx;
2565 void * p_components = couple_config->components;
2567 &couple_config->num_components, &p_components,
2568 sizeof(couple_config->components[0]),
2569 comm, &dist_merge_vtable_component, &old_to_new_idx);
2570 couple_config->components = p_components;
2571
2572 // set new component_idx in couples
2573 for(size_t couple_idx = 0; couple_idx < couple_config->num_couples;
2574 ++couple_idx) {
2575 struct yac_couple_config_couple * couple =
2576 couple_config->couples + couple_idx;
2577 couple->component_indices[0] = old_to_new_idx[couple->component_indices[0]];
2578 couple->component_indices[1] = old_to_new_idx[couple->component_indices[1]];
2579 if (couple->component_indices[1] < couple->component_indices[0]) {
2580 size_t temp_component_index = couple->component_indices[0];
2581 couple->component_indices[0] = couple->component_indices[1];
2582 couple->component_indices[1] = temp_component_index;
2583 }
2584 for(size_t field_couple_idx = 0;
2585 field_couple_idx < couple->num_field_couples; ++field_couple_idx) {
2586 couple->field_couples[field_couple_idx].source.component_idx =
2587 old_to_new_idx[
2588 couple->field_couples[field_couple_idx].source.component_idx];
2589 couple->field_couples[field_couple_idx].target.component_idx =
2590 old_to_new_idx[
2591 couple->field_couples[field_couple_idx].target.component_idx];
2592 }
2593 }
2594 free(old_to_new_idx);
2595
2596 // merge the fields of each component
2597 for (size_t comp_idx = 0; comp_idx < couple_config->num_components; ++comp_idx)
2598 merge_fields(couple_config, comp_idx, comm);
2599}
2600
2602 size_t * num_field_couples,
2603 struct yac_couple_config_field_couple ** field_couples, MPI_Comm comm) {
2604
2605 void * p_field_couples = *field_couples;
2607 num_field_couples, &p_field_couples, sizeof(**field_couples),
2608 comm, &dist_merge_vtable_field_couple, NULL);
2609 *field_couples = p_field_couples;
2610}
2611
2612static void merge_couples(
2613 struct yac_couple_config * couple_config, MPI_Comm comm) {
2614
2615 // distribute and merge basic couple information while keeping the
2616 // individual field couples
2617 void * p_couples = couple_config->couples;
2619 &couple_config->num_couples, &p_couples,
2620 sizeof(couple_config->couples[0]),
2621 comm, &dist_merge_vtable_couple, NULL);
2622 couple_config->couples = p_couples;
2623}
2624
2626 struct yac_couple_config * couple_config, MPI_Comm comm) {
2627
2628 // distribute and merge basic couple information while keeping the
2629 // individual field couples
2630 void * p_config_outputs = couple_config->config_outputs;
2632 &couple_config->num_config_outputs, &p_config_outputs,
2633 sizeof(couple_config->config_outputs[0]),
2634 comm, &dist_merge_vtable_config_output, NULL);
2635 couple_config->config_outputs = p_config_outputs;
2636}
2637
2639 struct yac_couple_config * couple_config, MPI_Comm comm,
2640 char const * output_ref){
2641
2642 // sync time stuff
2643 couple_config_sync_time(couple_config, comm);
2644
2645 merge_config_output(couple_config, comm);
2646 merge_grids(couple_config, comm);
2647 merge_components(couple_config, comm);
2648 merge_couples(couple_config, comm);
2649
2651 "missing_definition_is_fatal",
2652 &(couple_config->missing_definition_is_fatal), comm);
2653
2654 struct yac_couple_config_config_output * config_output = NULL;
2655 if (output_ref != NULL)
2656 for (size_t i = 0; (i < couple_config->num_config_outputs) &&
2657 (config_output == NULL); ++i)
2658 if (!strcmp(output_ref, couple_config->config_outputs[i].ref))
2659 config_output = couple_config->config_outputs + i;
2660
2661 if (config_output != NULL) {
2662
2663 int rank;
2664 MPI_Comm_rank(comm, &rank);
2665
2666 // rank 0 writes the coupling configuration to file
2667 if (rank == 0) {
2668
2669 FILE * config_file = fopen(config_output->name, "w");
2670
2672 config_file != NULL,
2673 "ERROR(yac_couple_config_sync): "
2674 "failed to create coupling configuration file \"%s\"",
2675 config_output->name);
2676
2678 (config_output->type == YAC_TEXT_FILETYPE_YAML) ||
2679 (config_output->type == YAC_TEXT_FILETYPE_JSON),
2680 "ERROR(yac_couple_config_sync): "
2681 "invalid coupling configuration filetype (type = %d)",
2682 config_output->type);
2683
2684 int emit_flags;
2685 switch(config_output->type) {
2686 default:
2688 emit_flags = YAC_YAML_EMITTER_DEFAULT;
2689 break;
2691 emit_flags = YAC_YAML_EMITTER_JSON;
2692 break;
2693 }
2694
2695 enum flag_value include_definitions = FLAG_FALSE; // default value
2696 if (config_output->include_definitions != FLAG_UNSET)
2698 char * str_couple_config =
2700 couple_config, emit_flags, (int)include_definitions);
2701
2702 fputs(str_couple_config, config_file);
2703 free(str_couple_config);
2704 fclose(config_file);
2705 }
2706 yac_mpi_call(MPI_Barrier(comm), comm);
2707 }
2708}
2709
2711 struct yac_couple_config * couple_config,
2712 char const * filename, enum yac_text_filetype filetype, char const * ref,
2713 int include_definitions) {
2714
2715 YAC_ASSERT(
2716 filename != NULL,
2717 "ERROR(yac_couple_config_set_config_output_filename): filename is NULL");
2718
2720 (filetype == YAC_TEXT_FILETYPE_YAML) ||
2721 (filetype == YAC_TEXT_FILETYPE_JSON),
2722 "ERROR(yac_couple_config_set_config_output_filename): "
2723 "invalid output configuration filetype (type = %d)",
2724 (int)filetype);
2725
2726 YAC_ASSERT(
2727 ref != NULL,
2728 "ERROR(yac_couple_config_set_config_output_filename): ref is NULL");
2729
2730 // check whether there already exists a config output with the same reference
2731 for (size_t i = 0; i < couple_config->num_config_outputs; ++i) {
2732 if (!strcmp(couple_config->config_outputs[i].ref, ref)) {
2734 !strcmp(couple_config->config_outputs[i].name, filename) &&
2735 couple_config->config_outputs[i].type == filetype,
2736 "ERROR(yac_couple_config_set_config_output_filename): "
2737 "an filename has already been set for reference "
2738 "(ref \"%s\" curr filename \"%s\" filetype %d; "
2739 "new filename \"%s\" filetype %d",
2740 ref, couple_config->config_outputs[i].name,
2741 couple_config->config_outputs[i].type, filename, filetype);
2742 return;
2743 }
2744 }
2745
2746 couple_config->config_outputs =
2747 xrealloc(
2748 couple_config->config_outputs,
2749 (couple_config->num_config_outputs + 1) *
2750 sizeof(*(couple_config->config_outputs)));
2751
2752 couple_config->config_outputs[
2753 couple_config->num_config_outputs].name = strdup(filename);
2754 couple_config->config_outputs[
2755 couple_config->num_config_outputs].type = filetype;
2756 couple_config->config_outputs[
2757 couple_config->num_config_outputs].ref = strdup(ref);
2758 couple_config->config_outputs[
2759 couple_config->num_config_outputs].include_definitions =
2760 (enum flag_value)(include_definitions != 0);
2761 couple_config->num_config_outputs++;
2762}
2763
2765 struct yac_couple_config * couple_config,
2766 char const * tgt_component_name, char const * tgt_grid_name,
2767 char const * tgt_field_name, char const ** src_component_name,
2768 char const ** src_grid_name, char const ** src_field_name) {
2769
2770 size_t tgt_comp_idx =
2771 yac_couple_config_get_component_idx(couple_config, tgt_component_name);
2772 size_t tgt_grid_idx =
2773 yac_couple_config_get_grid_idx(couple_config, tgt_grid_name);
2774 size_t tgt_field_idx =
2776 couple_config, tgt_comp_idx, tgt_grid_idx, tgt_field_name);
2777
2778 struct yac_couple_config_field_couple * field_couple = NULL;
2779 for (size_t couple_idx = 0; couple_idx < couple_config->num_couples;
2780 ++couple_idx) {
2781
2782 struct yac_couple_config_couple * couple =
2783 couple_config->couples + couple_idx;
2784
2785 if ((couple->component_indices[0] != tgt_comp_idx) &&
2786 (couple->component_indices[1] != tgt_comp_idx))
2787 continue;
2788
2789 for (size_t temp_field_couple_idx = 0;
2790 temp_field_couple_idx < couple->num_field_couples;
2791 ++temp_field_couple_idx) {
2792
2793 struct yac_couple_config_field_couple * temp_field_couple =
2794 couple->field_couples + temp_field_couple_idx;
2795
2796 if ((temp_field_couple->target.component_idx == tgt_comp_idx) &&
2797 (temp_field_couple->target.field_idx == tgt_field_idx)) {
2798
2800 field_couple == NULL,
2801 "ERROR(yac_couple_config_get_field_source): "
2802 "multiple couples with the same target field is not supported by "
2803 "this routine (target component: \"%s\" target grid: \"%s\""
2804 "target field: \"%s\")", tgt_component_name, tgt_grid_name,
2805 tgt_field_name);
2806
2807 field_couple = temp_field_couple;
2808 }
2809 }
2810 }
2811
2813 field_couple != NULL,
2814 "ERROR(yac_couple_config_get_field_source): "
2815 "provided field is not defined as a target in any coupling "
2816 "(target component: \"%s\" target grid: \"%s\" target field: \"%s\")",
2817 tgt_component_name, tgt_grid_name, tgt_field_name);
2818
2819 struct yac_couple_config_component * src_component =
2820 couple_config->components + field_couple->source.component_idx;
2821 struct yac_couple_config_field * src_field =
2822 src_component->fields + field_couple->source.field_idx;
2823 struct yac_couple_config_grid * src_grid =
2824 couple_config->grids + src_field->grid_idx;
2825
2826 *src_component_name = src_component->name;
2827 *src_grid_name = src_grid->name;
2828 *src_field_name = src_field->name;
2829}
char * yac_yaml_emit_coupling(struct yac_couple_config *couple_config, int emit_flags, int include_definitions)
int const YAC_YAML_EMITTER_DEFAULT
emit to YAML format
Definition config_yaml.c:63
int const YAC_YAML_EMITTER_JSON
emit to JSON format
Definition config_yaml.c:64
#define UNUSED(x)
Definition core.h:73
void yac_couple_config_grid_set_metadata(struct yac_couple_config *couple_config, char const *grid_name, const char *metadata)
struct yac_dist_merge_vtable dist_merge_vtable_couple
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_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)
void yac_couple_config_get_field_source(struct yac_couple_config *couple_config, char const *tgt_component_name, char const *tgt_grid_name, char const *tgt_field_name, char const **src_component_name, char const **src_grid_name, char const **src_field_name)
static void set_datetime(char const *type_datetime, struct _datetime **old, char const *str_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_)
char const * yac_couple_config_get_yaxt_exchanger_name(struct yac_couple_config *couple_config, size_t couple_idx, size_t field_couple_idx)
struct yac_dist_merge_vtable dist_merge_vtable_config_output
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)
void yac_couple_config_sync(struct yac_couple_config *couple_config, MPI_Comm comm, char const *output_ref)
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 couple_config_sync_string(char const *string_name, char **string, MPI_Comm comm)
static void yac_couple_config_config_output_free(void *config_output_)
static void couple_config_sync_datetime(char const *type_datetime, struct _datetime **datetime, MPI_Comm comm)
static size_t yac_couple_config_get_field_couple_pack_size(void *field_couple_, MPI_Comm comm)
static void merge_fields(struct yac_couple_config *couple_config, size_t comp_idx, MPI_Comm comm)
static void yac_couple_config_couple_free(void *couple_)
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_unpack_field(void *buffer, int buffer_size, int *position, void *field_, MPI_Comm comm)
const char * yac_couple_config_grid_get_output_filename(struct yac_couple_config *couple_config, const char *grid_name)
char * yac_couple_config_get_start_datetime(struct yac_couple_config *couple_config)
void yac_couple_config_set_datetime(struct yac_couple_config *couple_config, char const *start, char const *end)
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_unpack_config_output(void *buffer, int buffer_size, int *position, void *config_output_, 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)
struct yac_dist_merge_vtable dist_merge_vtable_grid
static void yac_couple_config_pack_field(void *field_, void *buffer, int buffer_size, int *position, 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)
flag_value
@ FLAG_UNSET
@ FLAG_FALSE
@ FLAG_TRUE
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)
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 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)
size_t yac_couple_config_get_field_collection_size(struct yac_couple_config *couple_config, char const *component_name, char const *grid_name, char const *field_name)
static size_t yac_couple_config_add_grid_(struct yac_couple_config *couple_config, char const *name)
static size_t yac_couple_config_get_config_output_pack_size(void *config_output_, MPI_Comm comm)
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)
static void couple_config_sync_flag_value(char const *flag_value_name, enum flag_value *flag, MPI_Comm comm)
static void yac_couple_config_unpack_grid(void *buffer, int buffer_size, int *position, void *grid_, MPI_Comm comm)
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 size_t yac_couple_config_get_couple_pack_size(void *couple_, MPI_Comm comm)
static char * string_dup(char const *string)
static void yac_couple_config_pack_config_output(void *config_output_, void *buffer, int buffer_size, int *position, MPI_Comm comm)
static void yac_couple_config_pack_field_couple(void *field_couple_, void *buffer, int buffer_size, int *position, MPI_Comm comm)
struct yac_dist_merge_vtable dist_merge_vtable_field
static size_t yac_couple_config_get_field_pack_size(void *field_, MPI_Comm comm)
static void yac_couple_config_pack_grid(void *grid_, void *buffer, int buffer_size, int *position, MPI_Comm comm)
size_t yac_couple_config_get_num_couples(struct yac_couple_config *couple_config)
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)
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, char const *yaxt_exchanger_name)
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)
int yac_couple_config_get_source_lag(struct yac_couple_config *couple_config, size_t couple_idx, size_t field_couple_idx)
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)
static int yac_couple_config_couple_compare(void const *a_, void const *b_)
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)
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)
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)
struct yac_dist_merge_vtable dist_merge_vtable_component
void yac_couple_config_add_grid(struct yac_couple_config *couple_config, char const *name)
static void yac_couple_config_unpack_couple(void *buffer, int buffer_size, int *position, void *couple_, MPI_Comm comm)
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)
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_unpack_field_couple(void *buffer, int buffer_size, int *position, void *field_couple_, MPI_Comm comm)
static size_t yac_couple_config_get_string_pack_size(char const *string, MPI_Comm comm)
static void yac_couple_config_config_output_merge(void *a_, void *b_, MPI_Comm comm)
static size_t yac_couple_config_get_component_pack_size(void *component_, 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_pack_couple(void *couple_, void *buffer, int buffer_size, int *position, MPI_Comm comm)
static void yac_couple_config_grid_free(void *grid_)
void yac_couple_config_grid_set_output_filename(struct yac_couple_config *couple_config, char const *grid_name, const char *output_filename)
static void merge_config_output(struct yac_couple_config *couple_config, MPI_Comm comm)
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)
static void yac_couple_config_unpack_component(void *buffer, int buffer_size, int *position, void *component_, MPI_Comm comm)
struct yac_couple_config * yac_couple_config_new()
char const * yac_couple_config_get_weight_file_name(struct yac_couple_config *couple_config, size_t couple_idx, size_t field_couple_idx)
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)
int yac_couple_config_get_missing_definition_is_fatal(struct yac_couple_config *couple_config)
struct yac_dist_merge_vtable dist_merge_vtable_field_couple
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)
void yac_couple_config_set_missing_definition_is_fatal(struct yac_couple_config *couple_config, int missing_definition_is_fatal)
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)
void yac_couple_config_add_component(struct yac_couple_config *couple_config, char const *name)
static int yac_couple_config_config_output_compare(void const *a, void const *b)
size_t yac_couple_config_get_num_grids(struct yac_couple_config *couple_config)
char const * yac_couple_config_get_field_grid_name(struct yac_couple_config *couple_config, size_t component_idx, size_t field_idx)
void yac_couple_config_set_config_output_filename(struct yac_couple_config *couple_config, char const *filename, enum yac_text_filetype filetype, char const *ref, int include_definitions)
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)
char * yac_couple_config_get_end_datetime(struct yac_couple_config *couple_config)
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_pack_component(void *component_, void *buffer, int buffer_size, int *position, 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)
size_t yac_couple_config_get_num_components(struct yac_couple_config *couple_config)
#define MISSING_DEFINITION_IS_FATAL_DEFAULT_VALUE
yac_text_filetype
@ YAC_TEXT_FILETYPE_YAML
YAML format.
@ YAC_TEXT_FILETYPE_JSON
JSON format.
yac_reduction_type
@ TIME_NONE
@ TIME_ACCUMULATE
@ TIME_MAXIMUM
@ TIME_MINIMUM
@ TIME_AVERAGE
void yac_dist_merge(size_t *count, void **array, size_t element_size, MPI_Comm comm, struct yac_dist_merge_vtable *vtable, size_t **idx_old_to_new)
Definition dist_merge.c:64
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_NO_VALUE
#define YAC_FRAC_MASK_VALUE_IS_SET(value)
#define YAC_FRAC_MASK_VALUE_IS_VALID(value)
#define xrealloc(ptr, size)
Definition ppm_xfuncs.h:67
#define xmalloc(size)
Definition ppm_xfuncs.h:66
struct yac_couple_config_field * fields
enum yac_text_filetype type
struct yac_couple_config_field_couple * field_couples
struct yac_couple_config_field_couple::@60 source
struct yac_interp_stack_config * interp_stack
enum yac_reduction_type coupling_period_operation
struct yac_couple_config_field_couple::@60 target
char const * output_filename
struct yac_couple_config_couple * couples
struct yac_couple_config_config_output * config_outputs
struct _datetime * end_datetime
struct yac_couple_config_grid * grids
enum flag_value missing_definition_is_fatal
struct _datetime * start_datetime
struct yac_couple_config_component * components
size_t(* get_pack_size)(void *element, MPI_Comm comm)
Definition dist_merge.h:20
static struct user_input_data_component ** components
Definition yac.c:133
size_t num_grids
Definition yac.c:131
static size_t num_components
Definition yac.c:134
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:648
int yac_cget_calendar()
Definition yac.c:665
int const YAC_CALENDAR_NOT_SET
Definition yac.c:60
#define YAC_MAX_CHARLEN
Definition yac.h:81
#define YAC_ASSERT_F(exp, format,...)
Definition yac_assert.h:19
#define YAC_ASSERT(exp, msg)
Definition yac_assert.h:16
#define yac_mpi_call(call, comm)