1 | program receiver |
---|
2 | use mod_oasis |
---|
3 | implicit none |
---|
4 | integer :: i, kinfo, date |
---|
5 | integer :: comp_id, part_id, var_id |
---|
6 | integer, parameter :: n_points = 16 |
---|
7 | integer :: part_params(OASIS_Serial_Params) |
---|
8 | integer :: var_nodims(2) |
---|
9 | character(len=8) :: comp_name = "receiver" |
---|
10 | character(len=8) :: var_name = "FRECVATM" |
---|
11 | real(kind=8) :: field(n_points), error, epsilon |
---|
12 | integer :: intra_comm, intra_rank, intra_size |
---|
13 | integer :: inter_comm, inter_rank, inter_size, inter_rsize |
---|
14 | |
---|
15 | print '(2A)', "Component name: ", comp_name |
---|
16 | |
---|
17 | call oasis_init_comp(comp_id, comp_name, kinfo) |
---|
18 | if(kinfo<0) call oasis_abort(comp_id, comp_name, & |
---|
19 | & "Error in oasis_init_comp: ", rcode=kinfo) |
---|
20 | print '(A,I0)', "Receiver: Component ID: ", comp_id |
---|
21 | |
---|
22 | part_params(OASIS_Strategy) = OASIS_Serial |
---|
23 | part_params(OASIS_Length) = n_points |
---|
24 | call oasis_def_partition(part_id, part_params, kinfo) |
---|
25 | if(kinfo<0) call oasis_abort(comp_id, comp_name, & |
---|
26 | & "Error in oasis_def_partition: ", rcode=kinfo) |
---|
27 | print '(A,I0)', "Receiver: part_id: ", part_id |
---|
28 | |
---|
29 | var_nodims=[1, 1] |
---|
30 | print '(2A)', "Receiver: var_name: ", var_name |
---|
31 | call oasis_def_var(var_id, var_name, part_id, var_nodims, OASIS_IN, & |
---|
32 | & [1], OASIS_REAL, kinfo) |
---|
33 | if(kinfo<0 .or. var_id<0) call oasis_abort(comp_id, comp_name, & |
---|
34 | & "Error in oasis_def_var: ", rcode=kinfo) |
---|
35 | print '(A,I0)', "Receiver: var_id: ", var_id |
---|
36 | |
---|
37 | call oasis_enddef(kinfo) |
---|
38 | if(kinfo<0) call oasis_abort(comp_id, comp_name, & |
---|
39 | & "Error in oasis_enddef: ", rcode=kinfo) |
---|
40 | |
---|
41 | call oasis_get_intracomm(intra_comm,"sender-box",kinfo) |
---|
42 | if(kinfo<0) call oasis_abort(comp_id, comp_name, & |
---|
43 | &"Error in oasis_get_intracomm: ", rcode=kinfo) |
---|
44 | |
---|
45 | call mpi_comm_size(intra_comm, intra_size, kinfo) |
---|
46 | call mpi_comm_rank(intra_comm, intra_rank, kinfo) |
---|
47 | print '(A,I0,A,I0)', "Receiver intra_comm: rank = ",intra_rank, " of ",intra_size |
---|
48 | |
---|
49 | call oasis_get_intercomm(inter_comm,"sender-box",kinfo) |
---|
50 | if(kinfo<0) call oasis_abort(comp_id, comp_name, & |
---|
51 | &"Error in oasis_get_intercomm: ", rcode=kinfo) |
---|
52 | |
---|
53 | call mpi_comm_size(inter_comm, inter_size, kinfo) |
---|
54 | call mpi_comm_rank(inter_comm, inter_rank, kinfo) |
---|
55 | call mpi_comm_remote_size(inter_comm, inter_rsize, kinfo) |
---|
56 | print '(A,I0,A,I0,A,I0)', "Receiver inter_comm: rank = ",inter_rank, " of ",inter_size, & |
---|
57 | & " Remote size = ",inter_rsize |
---|
58 | |
---|
59 | date=0 |
---|
60 | |
---|
61 | field(:)=0 |
---|
62 | |
---|
63 | call oasis_get(var_id, date, field, kinfo) |
---|
64 | if(kinfo<0) call oasis_abort(comp_id, comp_name, & |
---|
65 | & "Error in oasis_get: ", rcode=kinfo) |
---|
66 | |
---|
67 | call oasis_terminate(kinfo) |
---|
68 | if(kinfo<0) call oasis_abort(comp_id, comp_name, & |
---|
69 | & "Error in oasis_terminate: ", rcode=kinfo) |
---|
70 | |
---|
71 | epsilon=1e-8 |
---|
72 | error=0. |
---|
73 | do i = 1, n_points |
---|
74 | error=error+abs(field(i)-i) |
---|
75 | end do |
---|
76 | if(error<epsilon) print '(A)', "Receiver: Data received successfully" |
---|
77 | |
---|
78 | end program receiver |
---|