timesync2¶
Two models, A & B, run simulations over time and exchange data via a timestep synchronization communication pattern. The timesync model statesync is explicit in the YAML and allows for specific control of the synchronization (e.g. relationship between variables, interpolation, aggregation). The models compute different variables, have different timesteps, and output data at each time step to a file. This example demonstrates more advanced use of the timesync communication pattern and YAML parameters.
C Version¶
Model Code:
1#define _USE_MATH_DEFINES // Required to use M_PI with MSVC
2#include <math.h>
3#include <stdio.h>
4#include "YggInterface.h"
5
6
7int timestep_calc(double t, const char* t_units, generic_t state,
8 const char* model) {
9 double x_period = 10.0; // Days
10 double y_period = 5.0; // Days
11 double z_period = 20.0; // Days
12 double o_period = 2.5; // Days
13 int ret = 0;
14 if (strcmp(t_units, "day") == 0) {
15 // No conversion necessary
16 } else if (strcmp(t_units, "hr") == 0) {
17 x_period = x_period * 24.0;
18 y_period = y_period * 24.0;
19 z_period = z_period * 24.0;
20 o_period = o_period * 24.0;
21 } else {
22 printf("timestep_calc: Unsupported unit '%s'\n", t_units);
23 ret = -1;
24 }
25 if (ret >= 0) {
26 if (strcmp(model, "A") == 0) {
27 ret = generic_map_set_double(state, "x",
28 sin(2.0 * M_PI * t / x_period),
29 "");
30 ret = generic_map_set_double(state, "y",
31 cos(2.0 * M_PI * t / y_period),
32 "");
33 ret = generic_map_set_double(state, "z1",
34 -cos(2.0 * M_PI * t / z_period),
35 "");
36 ret = generic_map_set_double(state, "z2",
37 -cos(2.0 * M_PI * t / z_period),
38 "");
39 ret = generic_map_set_double(state, "a",
40 sin(2.0 * M_PI * t / o_period),
41 "");
42 } else {
43 ret = generic_map_set_double(state, "xvar",
44 sin(2.0 * M_PI * t / x_period) / 2.0,
45 "");
46 ret = generic_map_set_double(state, "yvar",
47 cos(2.0 * M_PI * t / y_period),
48 "");
49 ret = generic_map_set_double(state, "z",
50 -2.0 * cos(2.0 * M_PI * t / z_period),
51 "");
52 ret = generic_map_set_double(state, "b",
53 cos(2.0 * M_PI * t / o_period),
54 "");
55 }
56 }
57 return ret;
58}
59
60
61int main(int argc, char *argv[]) {
62
63 double t_step = atof(argv[1]);
64 char* t_units = argv[2];
65 char* model = argv[3];
66 int exit_code = 0;
67 printf("Hello from C timesync: timestep %f %s\n", t_step, t_units);
68 double t_start = 0.0;
69 double t_end = 5.0;
70 size_t nkeys, ikey;
71 char** keys = NULL;
72 if (strcmp(t_units, "hr") == 0) {
73 t_end = 24.0 * t_end;
74 }
75 int ret;
76 generic_t state_send = init_generic_map();
77 generic_t state_recv = init_generic_map();
78 ret = timestep_calc(t_start, t_units, state_send, model);
79 if (ret < 0) {
80 printf("timesync(C): Error in initial timestep calculation.");
81 return -1;
82 }
83
84 // Set up connections matching yaml
85 // Timestep synchronization connection will be 'statesync'
86 comm_t* timesync = yggTimesync("statesync", t_units);
87 dtype_t* out_dtype = create_dtype_json_object(0, NULL, NULL, true);
88 comm_t* out = yggOutputType("output", out_dtype);
89
90 // Initialize state and synchronize with other models
91 double t = t_start;
92 ret = rpcCall(timesync, t, state_send, &state_recv);
93 if (ret < 0) {
94 printf("timesync(C): Initial sync failed.\n");
95 return -1;
96 }
97 printf("timesync(C): t = %5.1f %-3s", t, t_units);
98 nkeys = generic_map_get_keys(state_recv, &keys);
99 for (ikey = 0; ikey < nkeys; ikey++) {
100 printf(", %s = %+ 5.2f", keys[ikey],
101 generic_map_get_double(state_recv, keys[ikey]));
102 }
103 printf("\n");
104
105 // Send initial state to output
106 generic_t msg = copy_generic(state_recv);
107 ret = generic_map_set_double(msg, "time", t, t_units);
108 if (ret < 0) {
109 printf("timesync(C): Failed to set time in initial output map.\n");
110 return -1;
111 }
112 ret = yggSend(out, msg);
113 if (ret < 0) {
114 printf("timesync(C): Failed to send initial output for t=%f.\n", t);
115 return -1;
116 }
117 destroy_generic(&msg);
118
119 // Iterate until end
120 while (t < t_end) {
121
122 // Perform calculations to update the state
123 t = t + t_step;
124 ret = timestep_calc(t, t_units, state_send, model);
125 if (ret < 0) {
126 printf("timesync(C): Error in timestep calculation for t = %f.\n", t);
127 return -1;
128 }
129
130 // Synchronize the state
131 ret = rpcCall(timesync, t, state_send, &state_recv);
132 if (ret < 0) {
133 printf("timesync(C): sync for t=%f failed.\n", t);
134 return -1;
135 }
136 printf("timesync(C): t = %5.1f %-3s", t, t_units);
137 nkeys = generic_map_get_keys(state_recv, &keys);
138 for (ikey = 0; ikey < nkeys; ikey++) {
139 printf(", %s = %+ 5.2f", keys[ikey],
140 generic_map_get_double(state_recv, keys[ikey]));
141 }
142 printf("\n");
143
144 // Send output
145 msg = copy_generic(state_recv);
146 ret = generic_map_set_double(msg, "time", t, t_units);
147 if (ret < 0) {
148 printf("timesync(C): Failed to set time in output map.\n");
149 return -1;
150 }
151 ret = yggSend(out, msg);
152 if (ret < 0) {
153 printf("timesync(C): Failed to send output for t=%f.\n", t);
154 return -1;
155 }
156 destroy_generic(&msg);
157
158 }
159
160 printf("Goodbye from C timesync\n");
161 destroy_generic(&state_send);
162 destroy_generic(&state_recv);
163 return 0;
164
165}
Model YAML:
1---
2
3models:
4 - name: statesync
5 language: timesync
6 synonyms:
7 modelB:
8 x:
9 alt: xvar
10 alt2base: ./src/timesync.py:xvar2x
11 base2alt: ./src/timesync.py:x2xvar
12 y: yvar
13 modelA:
14 z:
15 alt: [z1, z2]
16 alt2base: ./src/timesync.py:merge_z
17 base2alt: ./src/timesync.py:split_z
18 interpolation:
19 modelA: krogh
20 modelB:
21 method: spline
22 order: 5
23 aggregation:
24 x: ./src/timesync.py:xagg
25 y: sum
26 additional_variables:
27 modelA: [b]
28 modelB: [a]
29 - name: modelA
30 language: c
31 args:
32 - ./src/timesync.c
33 - 7 # Pass the timestep in hours
34 - hr
35 - A
36 timesync: statesync
37 outputs:
38 name: output
39 default_file:
40 name: modelA_output.txt
41 in_temp: True
42 filetype: table
43 - name: modelB
44 language: c
45 args:
46 - ./src/timesync.c
47 - 1 # Pass the timestep in days
48 - day
49 - B
50 timesync: statesync
51 outputs:
52 name: output
53 default_file:
54 name: modelB_output.txt
55 in_temp: True
56 filetype: table
C++ Version¶
Model Code:
1#define _USE_MATH_DEFINES // Required to use M_PI with MSVC
2#include <math.h>
3#include <stdio.h>
4#include "YggInterface.hpp"
5
6
7void timestep_calc(rapidjson::units::Quantity<double>& t,
8 rapidjson::Document& state, std::string model) {
9 rapidjson::units::Quantity<double> x_period(10.0, "days");
10 rapidjson::units::Quantity<double> y_period(5.0, "days");
11 rapidjson::units::Quantity<double> z_period(10.0, "days");
12 rapidjson::units::Quantity<double> o_period(2.5, "days");
13#define SET_(key, val) \
14 if (!state.HasMember(key)) \
15 state.AddMember(key,rapidjson::Value(val).Move(), \
16 state.GetAllocator()); \
17 else \
18 state[key].SetDouble(val)
19 if (model == "A") {
20 SET_("x", sin(2.0 * M_PI * (t / x_period).value()));
21 SET_("y", cos(2.0 * M_PI * (t / y_period).value()));
22 SET_("z1", -cos(2.0 * M_PI * (t / z_period).value()));
23 SET_("z2", -cos(2.0 * M_PI * (t / z_period).value()));
24 SET_("a", sin(2.0 * M_PI * (t / o_period).value()));
25 } else {
26 SET_("xvar", sin(2.0 * M_PI * (t / x_period).value()) / 2.0);
27 SET_("yvar", cos(2.0 * M_PI * (t / y_period).value()));
28 SET_("z", -2.0 * cos(2.0 * M_PI * (t / z_period).value()));
29 SET_("b", cos(2.0 * M_PI * (t / o_period).value()));
30 }
31}
32
33
34int main(int argc, char *argv[]) {
35
36 char* t_units = argv[2];
37 rapidjson::units::Quantity<double> t_step(atof(argv[1]), t_units);
38 std::cout << "Hello from C++ timesync: timestep " << t_step << std::endl;
39 rapidjson::units::Quantity<double> t_start(0.0, t_units);
40 rapidjson::units::Quantity<double> t_end(5.0, "days");
41
42 std::string model(argv[3]);
43 int ret;
44 rapidjson::Document state_send(rapidjson::kObjectType);
45 rapidjson::Document state_recv(rapidjson::kObjectType);
46 timestep_calc(t_start, state_send, model);
47
48 // Set up connections matching yaml
49 // Timestep synchronization connection will be 'statesync'
50 YggTimesync timesync("statesync", t_units);
51 dtype_t* out_dtype = create_dtype_json_object(0, NULL, NULL, true);
52 YggOutput out("output", out_dtype);
53
54 // Initialize state and synchronize with other models
55 rapidjson::units::Quantity<double> t = t_start;
56 ret = timesync.call(3, t.value(), &state_send, &state_recv);
57 if (ret < 0) {
58 std::cerr << "timesync(C++): Initial sync failed." << std::endl;
59 return -1;
60 }
61 std::cout << "timesync(C++): t = " << t;
62 for (rapidjson::Value::MemberIterator it = state_recv.MemberBegin();
63 it != state_recv.MemberEnd(); it++)
64 std::cout << ", " << it->name.GetString() <<
65 " = " << it->value.GetDouble();
66 std::cout << std::endl;
67
68 // Send initial state to output
69 rapidjson::Document msg;
70 msg.CopyFrom(state_recv, msg.GetAllocator());
71 msg.AddMember("time",
72 rapidjson::Value(t, msg.GetAllocator()).Move(),
73 msg.GetAllocator());
74 ret = out.send(1, &msg);
75 if (ret < 0) {
76 std::cerr << "timesync(C++): Failed to send initial output for t=" <<
77 t << std::endl;
78 return -1;
79 }
80
81 // Iterate until end
82 while (t < t_end) {
83
84 // Perform calculations to update the state
85 t = t + t_step;
86 timestep_calc(t, state_send, model);
87
88 // Synchronize the state
89 ret = timesync.call(3, t.value(), &state_send, &state_recv);
90 if (ret < 0) {
91 std::cerr << "timesync(C++): sync for t=" << t << " failed" << std::endl;
92 return -1;
93 }
94 std::cout << "timesync(C++): t = " << t;
95 for (rapidjson::Value::MemberIterator it = state_recv.MemberBegin();
96 it != state_recv.MemberEnd(); it++)
97 std::cout << ", " << it->name.GetString() <<
98 " = " << it->value.GetDouble();
99 std::cout << std::endl;
100
101 // Send output
102 msg.CopyFrom(state_recv, msg.GetAllocator());
103 msg.AddMember("time",
104 rapidjson::Value(t, msg.GetAllocator()).Move(),
105 msg.GetAllocator());
106 ret = out.send(1, &msg);
107 if (ret < 0) {
108 std::cerr << "timesync(C++): Failed to send output for t=" << t << std::endl;
109 return -1;
110 }
111 }
112
113 std::cout << "Goodbye from C++ timesync" << std::endl;
114 return 0;
115
116}
Model YAML:
1---
2
3models:
4 - name: statesync
5 language: timesync
6 synonyms:
7 modelB:
8 x:
9 alt: xvar
10 alt2base: ./src/timesync.py:xvar2x
11 base2alt: ./src/timesync.py:x2xvar
12 y: yvar
13 modelA:
14 z:
15 alt: [z1, z2]
16 alt2base: ./src/timesync.py:merge_z
17 base2alt: ./src/timesync.py:split_z
18 interpolation:
19 modelA: krogh
20 modelB:
21 method: spline
22 order: 5
23 aggregation:
24 x: ./src/timesync.py:xagg
25 y: sum
26 additional_variables:
27 modelA: [b]
28 modelB: [a]
29 - name: modelA
30 language: c++
31 args:
32 - ./src/timesync.cpp
33 - 7 # Pass the timestep in hours
34 - hr
35 - A
36 timesync: statesync
37 outputs:
38 name: output
39 default_file:
40 name: modelA_output.txt
41 in_temp: True
42 filetype: table
43 - name: modelB
44 language: c++
45 args:
46 - ./src/timesync.cpp
47 - 1 # Pass the timestep in days
48 - day
49 - B
50 timesync: statesync
51 outputs:
52 name: output
53 default_file:
54 name: modelB_output.txt
55 in_temp: True
56 filetype: table
Fortran Version¶
Model Code:
1program main
2 use fygg
3
4 ! Declare resulting variables and create buffer for received message
5 logical :: timestep_calc
6 logical :: ret = .true.
7 type(yggcomm) :: timesync, out
8 character(len=32) :: arg
9 character(len=0), dimension(0) :: dtype_keys
10 type(yggdtype), dimension(0) :: dtype_vals
11 type(yggdtype) :: out_dtype
12 real(kind=8) :: t_step, t_start, t_end, t
13 character(len=32) :: t_units, model
14 character(len=20), dimension(:), pointer :: keys
15 type(ygggeneric) :: state_send, state_recv, msg
16 real(kind=8), pointer :: x
17 integer :: i
18
19 call get_command_argument(1, arg)
20 read(arg, *) t_step
21 call get_command_argument(2, arg)
22 read(arg, *) t_units
23 call get_command_argument(3, arg)
24 read(arg, *) model
25 write (*, '("Hello from Fortran timesync: timestep ",F10.5," ",A3,&
26 &" (model = ",A,")")') t_step, trim(t_units), trim(model)
27 t_start = 0.0
28 t_end = 5.0
29 if (t_units.eq."hr") then
30 t_end = 24.0 * t_end
31 end if
32 state_send = init_generic_map()
33 state_recv = init_generic_map()
34 ret = timestep_calc(t_start, t_units, state_send, model)
35 if (.not.ret) then
36 write (*, '("timesync(Fortran): Error in initial timestep &
37 &calculation.")')
38 end if
39
40 ! Set up connections matching yaml
41 ! Timestep synchronization connection will be 'statesync'
42 timesync = ygg_timesync("statesync", t_units)
43 out_dtype = create_dtype_json_object(0, dtype_keys, dtype_vals, .true.)
44 out = ygg_output_type("output", out_dtype)
45
46 ! Initialize state and synchronize with other models
47 t = t_start
48 ret = ygg_rpc_call(timesync, [yggarg(t), yggarg(state_send)], &
49 yggarg(state_recv))
50 if (.not.ret) then
51 write (*, '("timesync(Fortran): Initial sync failed.")')
52 stop 1
53 end if
54 write (*, '("timesync(Fortran): t = ",F5.1," ",A3)', advance="no") &
55 t, adjustl(t_units)
56 call generic_map_get_keys(state_recv, keys)
57 do i = 1, size(keys)
58 call generic_map_get(state_recv, trim(keys(i)), x)
59 write (*, '(SP,", ",A," = ",F5.2)', advance="no") &
60 trim(keys(i)), x
61 end do
62 write (*, '("")')
63
64 ! Send initial state to output
65 msg = copy_generic(state_recv)
66 call generic_map_set(msg, "time", t, t_units)
67 ret = ygg_send_var(out, yggarg(msg))
68 if (.not.ret) then
69 write (*, '("timesync(Fortran): Failed to send initial output &
70 &for t=",F10.5,".")') t
71 stop 1
72 end if
73 call free_generic(msg)
74
75 ! Iterate until end
76 do while (t.lt.t_end)
77
78 ! Perform calculations to update the state
79 t = t + t_step
80 ret = timestep_calc(t, t_units, state_send, model)
81 if (.not.ret) then
82 write (*, '("timesync(Fortran): Error in timestep &
83 &calculation for t = ",F10.5,".")') t
84 stop 1
85 end if
86
87 ! Synchronize the state
88 ret = ygg_rpc_call(timesync, [yggarg(t), yggarg(state_send)], &
89 yggarg(state_recv))
90 if (.not.ret) then
91 write (*, '("timesync(Fortran): sync failed for t=",F10.5,&
92 &".")') t
93 stop 1
94 end if
95 write (*, '("timesync(Fortran): t = ",F5.1," ",A3)', advance="no") &
96 t, adjustl(t_units)
97 call generic_map_get_keys(state_recv, keys)
98 do i = 1, size(keys)
99 call generic_map_get(state_recv, keys(i), x)
100 write (*, '(SP,", ",A," = ",F5.2)', advance="no") &
101 trim(keys(i)), x
102 end do
103 write (*, '("")')
104
105 ! Send output
106 msg = copy_generic(state_recv)
107 call generic_map_set(msg, "time", t, t_units)
108 ret = ygg_send_var(out, yggarg(msg))
109 if (.not.ret) then
110 write (*, '("timesync(Fortran): Failed to send output for &
111 &t=",F10.5,".")') t
112 stop 1
113 end if
114 call free_generic(msg)
115
116 end do
117
118 write (*, '("Goodbye from Fortran timesync")')
119 call free_generic(state_send)
120 call free_generic(state_recv)
121
122end program main
123
124
125function timestep_calc(t, t_units, state, model) result (ret)
126 use fygg
127 implicit none
128 real(kind=8) :: t
129 character(len=*) :: t_units
130 type(ygggeneric) :: state
131 character(len=*) :: model
132 logical :: ret
133 real(kind=8) :: x_period = 10.0 ! Days
134 real(kind=8) :: y_period = 5.0 ! Days
135 real(kind=8) :: z_period = 20.0 ! Days
136 real(kind=8) :: o_period = 2.5 ! Days
137 ret = .true.
138 if (t_units.eq."day") then
139 ! No conversion necessary
140 else if (t_units.eq."hr") then
141 x_period = x_period * 24.0
142 y_period = y_period * 24.0
143 z_period = z_period * 24.0
144 o_period = o_period * 24.0
145 else
146 write (*, '("timestep_calc: Unsupported unit ''",A,"''")') t_units
147 ret = .false.
148 end if
149 if (ret) then
150 if (model.eq."A") then
151 call generic_map_set(state, "x", &
152 sin(2.0 * PI_8 * t / x_period))
153 call generic_map_set(state, "y", &
154 cos(2.0 * PI_8 * t / y_period))
155 call generic_map_set(state, "z1", &
156 -cos(2.0 * PI_8 * t / z_period))
157 call generic_map_set(state, "z2", &
158 -cos(2.0 * PI_8 * t / z_period))
159 call generic_map_set(state, "a", &
160 sin(2.0 * PI_8 * t / o_period))
161 else
162 call generic_map_set(state, "xvar", &
163 sin(2.0 * PI_8 * t / x_period))
164 call generic_map_set(state, "yvar", &
165 cos(2.0 * PI_8 * t / y_period))
166 call generic_map_set(state, "z", &
167 -2.0 * cos(2.0 * PI_8 * t / z_period))
168 call generic_map_set(state, "b", &
169 cos(2.0 * PI_8 * t / o_period))
170 end if
171 end if
172end function timestep_calc
Model YAML:
1---
2
3models:
4 - name: statesync
5 language: timesync
6 synonyms:
7 modelB:
8 x:
9 alt: xvar
10 alt2base: ./src/timesync.py:xvar2x
11 base2alt: ./src/timesync.py:x2xvar
12 y: yvar
13 modelA:
14 z:
15 alt: [z1, z2]
16 alt2base: ./src/timesync.py:merge_z
17 base2alt: ./src/timesync.py:split_z
18 interpolation:
19 modelA: krogh
20 modelB:
21 method: spline
22 order: 5
23 aggregation:
24 x: ./src/timesync.py:xagg
25 y: sum
26 additional_variables:
27 modelA: [b]
28 modelB: [a]
29 - name: modelA
30 language: fortran
31 args:
32 - ./src/timesync.f90
33 - 7 # Pass the timestep in hours
34 - hr
35 - A
36 timesync: statesync
37 outputs:
38 name: output
39 default_file:
40 name: modelA_output.txt
41 in_temp: True
42 filetype: table
43 - name: modelB
44 language: fortran
45 args:
46 - ./src/timesync.f90
47 - 1 # Pass the timestep in days
48 - day
49 - B
50 timesync: statesync
51 outputs:
52 name: output
53 default_file:
54 name: modelB_output.txt
55 in_temp: True
56 filetype: table
Julia Version¶
Model Code:
1using Yggdrasil
2using Unitful
3using Printf
4
5function timestep_calc(t::Unitful.Quantity, model::String)
6 state = Dict()
7 if (model == "A")
8 state["x"] = sin(2.0 * t / Unitful.Quantity(10.0, u"d"))
9 state["y"] = cos(2.0 * t / Unitful.Quantity(5.0, u"d"))
10 state["z1"] = -cos(2.0 * t / Unitful.Quantity(20.0, u"d"))
11 state["z2"] = -cos(2.0 * t / Unitful.Quantity(20.0, u"d"))
12 state["a"] = sin(2.0 * t / Unitful.Quantity(2.5, u"d"))
13 else
14 state["xvar"] = sin(2.0 * t / Unitful.Quantity(10.0, u"d"))
15 state["yvar"] = cos(2.0 * t / Unitful.Quantity(5.0, u"d"))
16 state["z"] = -2.0 * cos(2.0 * t / Unitful.Quantity(20.0, u"d"))
17 state["b"] = cos(2.0 * t / Unitful.Quantity(2.5, u"d"))
18 end
19 return state
20end
21
22function main(t_step::Float64, t_units::String, model::String)
23
24 @printf("Hello from Julia timesync: timestep = %f %s\n", t_step, t_units)
25 t_step = Unitful.Quantity(t_step, Unitful.uparse(t_units))
26 t_start = Unitful.Quantity(0.0, Unitful.uparse(t_units))
27 t_end = Unitful.Quantity(1.0, u"d")
28 state = timestep_calc(t_start, model)
29
30 # Set up connections matching yaml
31 # Timestep synchronization connection will default to 'timesync'
32 timesync = Yggdrasil.YggInterface("YggTimesync", "statesync")
33 out = Yggdrasil.YggInterface("YggOutput", "output")
34
35 # Initialize state and synchronize with other models
36 t = t_start
37 ret, state = timesync.call(t, state)
38 if (!ret)
39 error("timesync(Julia): Initial sync failed.")
40 end
41 @printf("timesync(Julia): t = %s", t)
42 for (k, v) in state
43 @printf(", %s = %s", k, v)
44 end
45 @printf("\n")
46
47 # Send initial state to output
48 msg = state
49 msg["time"] = t
50 flag = out.send(msg)
51 if (!flag)
52 error(@sprintf("timesync(Julia): Failed to send initial output for t=%s\n", t))
53 end
54
55 # Iterate until end
56 while (t < t_end)
57
58 # Perform calculations to update the state
59 t = t + t_step
60 state = timestep_calc(t, model)
61
62 # Synchronize the state
63 ret, state = timesync.call(t, state)
64 if (!ret)
65 error(@sprintf("timesync(Julia): sync for t=%s failed.", t))
66 end
67 @printf("timesync(Julia): t = %s", t)
68 for (k, v) in state
69 @printf(", %s = %s", k, v)
70 end
71 @printf("\n")
72
73 # Send output
74 msg = state
75 msg["time"] = t
76 flag = out.send(msg)
77 if (!flag)
78 error(@sprintf("timesync(Julia): Failed to send output for t=%s.", t))
79 end
80
81 end
82
83 println("Goodbye from Julia timesync")
84
85end
86
87main(parse(Float64, ARGS[1]), ARGS[2], ARGS[3])
Model YAML:
1---
2
3models:
4 - name: statesync
5 language: timesync
6 synonyms:
7 modelB:
8 x:
9 alt: xvar
10 alt2base: ./src/timesync.py:xvar2x
11 base2alt: ./src/timesync.py:x2xvar
12 y: yvar
13 modelA:
14 z:
15 alt: [z1, z2]
16 alt2base: ./src/timesync.py:merge_z
17 base2alt: ./src/timesync.py:split_z
18 interpolation:
19 modelA: krogh
20 modelB:
21 method: spline
22 order: 5
23 aggregation:
24 x: ./src/timesync.py:xagg
25 y: sum
26 additional_variables:
27 modelA: [b]
28 modelB: [a]
29 - name: modelA
30 language: julia
31 args:
32 - ./src/timesync.jl
33 - 7 # Pass the timestep in hours
34 - hr
35 - A
36 timesync: statesync
37 outputs:
38 name: output
39 default_file:
40 name: modelA_output.txt
41 in_temp: True
42 filetype: table
43 - name: modelB
44 language: julia
45 args:
46 - ./src/timesync.jl
47 - 1 # Pass the timestep in days
48 - d
49 - B
50 timesync: statesync
51 outputs:
52 name: output
53 default_file:
54 name: modelB_output.txt
55 in_temp: True
56 filetype: table
Matlab Version¶
Model Code:
1function timesync(t_step, t_units, model)
2
3 t_step = str2num(t_step);
4 fprintf('Hello from Matlab timesync: timestep = %f %s\n', t_step, t_units);
5 t_step = t_step * str2symunit(t_units);
6 t_start = 0.0000000000000001 * str2symunit(t_units);
7 t_end = 5.0 * str2symunit('day');
8 state = timestep_calc(t_start, model);
9
10 % Set up connections matching yaml
11 % Timestep synchronization connection will be 'statesync'
12 timesync = YggInterface('YggTimesync', 'statesync');
13 out = YggInterface('YggOutput', 'output');
14
15 % Initialize state and synchronize with other models
16 t = t_start;
17 [ret, state] = timesync.call(t, state);
18 if (~ret);
19 error('timesync(Matlab): Initial sync failed.');
20 end;
21 [t_data, t_unit] = separateUnits(t);
22 fprintf('timesync(Matlab): t = %5.1f %-1s', ...
23 t_data, symunit2str(t_unit));
24 state_keys = keys(state);
25 state_vals = values(state);
26 for i = 1:length(state_keys)
27 fprintf(', %s = %+ 5.2f', state_keys{i}, state_vals{i});
28 end;
29 fprintf('\n');
30
31 % Send initial state to output
32 msg_keys = keys(state);
33 msg_keys{length(msg_keys) + 1} = 'time';
34 msg_vals = values(state);
35 msg_vals{length(msg_vals) + 1} = t;
36 msg = containers.Map(msg_keys, msg_vals, 'UniformValues', false);
37 flag = out.send(msg);
38
39 % Iterate until end
40 while (simplify(t/t_end) < 1)
41
42 % Perform calculations to update the state
43 t = t + t_step;
44 state = timestep_calc(t, model);
45
46 % Synchronize the state
47 [ret, state] = timesync.call(t, state);
48 if (~ret);
49 error(sprintf('timesync(Matlab): sync for t=%f failed.\n', t));
50 end;
51 [t_data, t_unit] = separateUnits(t);
52 fprintf('timesync(Matlab): t = %5.1f %-1s', ...
53 t_data, symunit2str(t_unit));
54 state_keys = keys(state);
55 state_vals = values(state);
56 for i = 1:length(state_keys)
57 fprintf(', %s = %+ 5.2f', state_keys{i}, state_vals{i});
58 end;
59 fprintf('\n');
60
61 % Send output
62 msg_keys = keys(state);
63 msg_keys{length(msg_keys) + 1} = 'time';
64 msg_vals = values(state);
65 msg_vals{length(msg_vals) + 1} = t;
66 msg = containers.Map(msg_keys, msg_vals, 'UniformValues', false);
67 flag = out.send(msg);
68 if (~flag);
69 error(sprintf('timesync(Matlab): Failed to send output for t=%s.\n', t));
70 end;
71 end;
72
73 disp('Goodbye from Matlab timesync');
74
75end
76
77
78function state = timestep_calc(t, model)
79 state = containers.Map('UniformValues', false, 'ValueType', 'any');
80 if (model == 'A')
81 state('x') = sin(2.0 * pi * t / (10.0 * str2symunit('day')));
82 state('y') = cos(2.0 * pi * t / (5.0 * str2symunit('day')));
83 state('z1') = -cos(2.0 * pi * t / (20.0 * str2symunit('day')));
84 state('z2') = -cos(2.0 * pi * t / (20.0 * str2symunit('day')));
85 state('a') = sin(2.0 * pi * t / (2.5 * str2symunit('day')));
86 else
87 state('xvar') = sin(2.0 * pi * t / (10.0 * str2symunit('day'))) / 2.0;
88 state('yvar') = cos(2.0 * pi * t / (5.0 * str2symunit('day')));
89 state('z') = -2.0 * cos(2.0 * pi * t / (20.0 * str2symunit('day')));
90 state('b') = cos(2.0 * pi * t / (2.5 * str2symunit('day')));
91 end;
92end
Model YAML:
1---
2
3models:
4 - name: statesync
5 language: timesync
6 synonyms:
7 modelB:
8 x:
9 alt: xvar
10 alt2base: ./src/timesync.py:xvar2x
11 base2alt: ./src/timesync.py:x2xvar
12 y: yvar
13 modelA:
14 z:
15 alt: [z1, z2]
16 alt2base: ./src/timesync.py:merge_z
17 base2alt: ./src/timesync.py:split_z
18 interpolation:
19 modelA: krogh
20 modelB:
21 method: spline
22 order: 5
23 aggregation:
24 x: ./src/timesync.py:xagg
25 y: sum
26 additional_variables:
27 modelA: [b]
28 modelB: [a]
29 - name: modelA
30 language: matlab
31 args:
32 - ./src/timesync.m
33 - 7 # Pass the timestep in hours
34 - hr
35 - A
36 timesync: statesync
37 outputs:
38 name: output
39 default_file:
40 name: modelA_output.txt
41 in_temp: True
42 filetype: table
43 - name: modelB
44 language: matlab
45 args:
46 - ./src/timesync.m
47 - 1 # Pass the timestep in days
48 - day
49 - B
50 timesync: statesync
51 outputs:
52 name: output
53 default_file:
54 name: modelB_output.txt
55 in_temp: True
56 filetype: table
Python Version¶
Model Code:
1import sys
2import numpy as np
3from yggdrasil import units
4from yggdrasil.interface.YggInterface import (
5 YggTimesync, YggOutput)
6
7
8def xvar2x(xvar):
9 r"""Convert xvar to x."""
10 x = 2 * xvar
11 return x
12
13
14def x2xvar(x):
15 r"""Convert x to xvar."""
16 xvar = x / 2
17 return xvar
18
19
20def xagg(series):
21 r"""Aggregate x variables."""
22 return max(series, key=abs)
23
24
25def merge_z(z1, z2):
26 r"""Merge the z1 and z2 variables."""
27 return z1 + z2
28
29
30def split_z(z):
31 r"""Split z into z1 and z2 variables."""
32 return (z / 2.0, z / 2.0)
33
34
35def timestep_calc(t, model):
36 r"""Updates the state based on the time where x is a sine wave
37 with period of 10 days and y is a cosine wave with a period of 5 days.
38 If model is 'A', the forth state variable will be 'a', a sine
39 with a period of 2.5 days. If model is 'B', the forth state
40 variable will be 'b', a cosine with a period of 2.5 days.
41
42 Args:
43 t (float): Current time.
44 model (str): Identifier for the model (A or B).
45
46 Returns:
47 dict: Map of state parameters.
48
49 """
50 if model == 'A':
51 state = {
52 'x': np.sin(2.0 * np.pi * t / units.add_units(10, 'day')),
53 'y': np.cos(2.0 * np.pi * t / units.add_units(5, 'day')),
54 'z1': -np.cos(2.0 * np.pi * t / units.add_units(20, 'day')),
55 'z2': -np.cos(2.0 * np.pi * t / units.add_units(20, 'day')),
56 'a': np.sin(2.0 * np.pi * t / units.add_units(2.5, 'day'))}
57 else:
58 state = {
59 'xvar': x2xvar(np.sin(2.0 * np.pi * t / units.add_units(10, 'day'))),
60 'yvar': np.cos(2.0 * np.pi * t / units.add_units(5, 'day')),
61 'z': -2.0 * np.cos(2.0 * np.pi * t / units.add_units(20, 'day')),
62 'b': np.cos(2.0 * np.pi * t / units.add_units(2.5, 'day'))}
63 return state
64
65
66def main(t_step, t_units, model):
67 r"""Function to execute integration.
68
69 Args:
70 t_step (float): The time step that should be used.
71 t_units (str): Units of the time step.
72 model (str): Identifier for the model (A or B).
73
74 """
75 print('Hello from Python timesync: timestep = %s %s' % (t_step, t_units))
76 t_step = units.add_units(t_step, t_units)
77 t_start = units.add_units(0.0, t_units)
78 t_end = units.add_units(5.0, 'day')
79 state = timestep_calc(t_start, model)
80
81 # Set up connections matching yaml
82 # Timestep synchronization connection will be 'statesync'
83 timesync = YggTimesync("statesync")
84 out = YggOutput('output')
85
86 # Initialize state and synchronize with other models
87 t = t_start
88 ret, state = timesync.call(t, state)
89 if not ret:
90 raise RuntimeError("timesync(Python): Initial sync failed.")
91 print('timesync(Python): t = % 8s' % t, end='')
92 for k, v in state.items():
93 print(', %s = %+ 5.2f' % (k, v), end='')
94 print('')
95
96 # Send initial state to output
97 flag = out.send(dict(state, time=t))
98 if not flag:
99 raise RuntimeError("timesync(Python): Failed to send "
100 "initial output for t=%s." % t)
101
102 # Iterate until end
103 while t < t_end:
104
105 # Perform calculations to update the state
106 t = t + t_step
107 state = timestep_calc(t, model)
108
109 # Synchronize the state
110 ret, state = timesync.call(t, state)
111 if not ret:
112 raise RuntimeError("timesync(Python): sync for t=%f failed." % t)
113 print('timesync(Python): t = % 8s' % t, end='')
114 for k, v in state.items():
115 print(', %s = %+ 5.2f' % (k, v), end='')
116 print('')
117
118 # Send output
119 flag = out.send(dict(state, time=t))
120 if not flag:
121 raise RuntimeError("timesync(Python): Failed to send output for t=%s." % t)
122
123 print('Goodbye from Python timesync')
124
125
126if __name__ == '__main__':
127 # Take time step from the first argument
128 main(float(sys.argv[1]), sys.argv[2], sys.argv[3])
Model YAML:
1---
2
3models:
4 - name: statesync
5 language: timesync
6 synonyms:
7 modelB:
8 x:
9 alt: xvar
10 alt2base: ./src/timesync.py:xvar2x
11 base2alt: ./src/timesync.py:x2xvar
12 y: yvar
13 modelA:
14 z:
15 alt: [z1, z2]
16 alt2base: ./src/timesync.py:merge_z
17 base2alt: ./src/timesync.py:split_z
18 interpolation:
19 modelA: krogh
20 modelB:
21 method: spline
22 order: 5
23 aggregation:
24 x: ./src/timesync.py:xagg
25 y: sum
26 additional_variables:
27 modelA: [b]
28 modelB: [a]
29 - name: modelA
30 language: python
31 args:
32 - ./src/timesync.py
33 - 7 # Pass the timestep in hours
34 - hr
35 - A
36 timesync: statesync
37 outputs:
38 name: output
39 default_file:
40 name: modelA_output.txt
41 in_temp: True
42 filetype: table
43 - name: modelB
44 language: python
45 args:
46 - ./src/timesync.py
47 - 1 # Pass the timestep in days
48 - day
49 - B
50 timesync: statesync
51 outputs:
52 name: output
53 default_file:
54 name: modelB_output.txt
55 in_temp: True
56 filetype: table
R Version¶
Model Code:
1library(yggdrasil)
2
3
4timestep_calc <- function(t, model) {
5 state = list()
6 if (model == 'A') {
7 state[['x']] = sinpi(2.0 * t / units::set_units(10.0, 'day', mode="standard"))
8 state[['y']] = cospi(2.0 * t / units::set_units(5.0, 'day', mode="standard"))
9 state[['z1']] = -cospi(2.0 * t / units::set_units(20.0, 'day', mode="standard"))
10 state[['z2']] = -cospi(2.0 * t / units::set_units(20.0, 'day', mode="standard"))
11 state[['a']] = sinpi(2.0 * t / units::set_units(2.5, 'day', model="standard"))
12 } else {
13 state[['xvar']] = sinpi(2.0 * t / units::set_units(10.0, 'day', mode="standard")) / 2.0
14 state[['yvar']] = cospi(2.0 * t / units::set_units(5.0, 'day', mode="standard"))
15 state[['z']] = -2.0 * cospi(2.0 * t / units::set_units(20.0, 'day', mode="standard"))
16 state[['b']] = cospi(2.0 * t / units::set_units(2.5, 'day', model="standard"))
17 }
18 return(state)
19}
20
21main <- function(t_step, t_units, model) {
22
23 fprintf('Hello from R timesync: timestep = %f %s', t_step, t_units)
24 t_step <- units::set_units(t_step, t_units, mode="standard")
25 t_start <- units::set_units(0.0, t_units, mode="standard")
26 t_end <- units::set_units(5.0, 'day', mode="standard")
27 state <- timestep_calc(t_start, model)
28
29 # Set up connections matching yaml
30 # Timestep synchronization connection will be 'statesync'
31 timesync <- YggInterface('YggTimesync', 'statesync')
32 out <- YggInterface('YggOutput', 'output')
33
34 # Initialize state and synchronize with other models
35 t <- t_start
36 c(ret, state) %<-% timesync$call(t, state)
37 if (!ret) {
38 stop('timesync(R): Initial sync failed.')
39 }
40 fprintf('timesync(R): t = %5.1f %-1s',
41 units::drop_units(t), units::deparse_unit(t))
42 for (key in names(state)) {
43 fprintf(', %s = %+ 5.2f', key, state[[key]])
44 }
45 fprintf('\n')
46
47 # Send initial state to output
48 msg = state
49 msg[['time']] = t
50 flag <- out$send(msg)
51 if (!flag) {
52 stop(sprintf('timesync(R): Failed to send initial output for t=%s', t))
53 }
54
55 # Iterate until end
56 while (t < t_end) {
57
58 # Perform calculations to update the state
59 t <- t + t_step
60 state <- timestep_calc(t, model)
61
62 # Synchronize the state
63 c(ret, state) %<-% timesync$call(t, state)
64 if (!ret) {
65 stop(sprintf('timesync(R): sync for t=%f failed.', t))
66 }
67 fprintf('timesync(R): t = %5.1f %-1s',
68 units::drop_units(t), units::deparse_unit(t))
69 for (key in names(state)) {
70 fprintf(', %s = %+ 5.2f', key, state[[key]])
71 }
72 fprintf('\n')
73
74 # Send output
75 msg = state
76 msg[['time']] = t
77 flag <- out$send(msg)
78 if (!flag) {
79 stop(sprintf('timesync(R): Failed to send output for t=%s.', t))
80 }
81 }
82
83 print('Goodbye from R timesync')
84
85}
86
87
88args = commandArgs(trailingOnly=TRUE)
89main(as.double(args[[1]]), args[[2]], args[[3]])
Model YAML:
1---
2
3models:
4 - name: statesync
5 language: timesync
6 synonyms:
7 modelB:
8 x:
9 alt: xvar
10 alt2base: ./src/timesync.py:xvar2x
11 base2alt: ./src/timesync.py:x2xvar
12 y: yvar
13 modelA:
14 z:
15 alt: [z1, z2]
16 alt2base: ./src/timesync.py:merge_z
17 base2alt: ./src/timesync.py:split_z
18 interpolation:
19 modelA: krogh
20 modelB:
21 method: spline
22 order: 5
23 aggregation:
24 x: ./src/timesync.py:xagg
25 y: sum
26 additional_variables:
27 modelA: [b]
28 modelB: [a]
29 - name: modelA
30 language: R
31 args:
32 - ./src/timesync.R
33 - 7 # Pass the timestep in hours
34 - hr
35 - A
36 timesync: statesync
37 outputs:
38 name: output
39 default_file:
40 name: modelA_output.txt
41 in_temp: True
42 filetype: table
43 - name: modelB
44 language: R
45 args:
46 - ./src/timesync.R
47 - 1 # Pass the timestep in days
48 - day
49 - B
50 timesync: statesync
51 outputs:
52 name: output
53 default_file:
54 name: modelB_output.txt
55 in_temp: True
56 filetype: table