1 | program sender_apple |
---|
2 | use mod_oasis |
---|
3 | use mpi |
---|
4 | implicit none |
---|
5 | integer :: i, kinfo, date |
---|
6 | integer :: comp_id, part_id, var_id |
---|
7 | integer :: n_points = 16 |
---|
8 | integer :: part_params(OASIS_Apple_Params), offset, local_size |
---|
9 | integer :: local_comm, local_comm_size, local_comm_rank |
---|
10 | integer :: icpl, coupl_comm, comm_size, comm_rank |
---|
11 | integer :: var_nodims(2) |
---|
12 | character(len=13) :: comp_name = "sender-apple" |
---|
13 | character(len=8) :: var_name = "FSENDOCN" |
---|
14 | real, allocatable :: field(:) |
---|
15 | |
---|
16 | print '(2A)', "Component name: ", comp_name |
---|
17 | |
---|
18 | call oasis_init_comp(comp_id, comp_name, kinfo) |
---|
19 | if(kinfo<0) call oasis_abort(comp_id, comp_name, & |
---|
20 | & "Error in oasis_init_comp: ", rcode=kinfo) |
---|
21 | print '(A,I0)', "Sender: Component ID: ", comp_id |
---|
22 | |
---|
23 | call oasis_get_localcomm(local_comm, kinfo) |
---|
24 | if(kinfo<0) call oasis_abort(comp_id, comp_name, & |
---|
25 | & "Error in oasis_get_localcomm: ", rcode=kinfo) |
---|
26 | |
---|
27 | call mpi_comm_size(local_comm, local_comm_size, kinfo) |
---|
28 | call mpi_comm_rank(local_comm, local_comm_rank, kinfo) |
---|
29 | print '(A,I0,A,I0)', "Sender: local_comm_rank = ",local_comm_rank, & |
---|
30 | & " of ",local_comm_size |
---|
31 | |
---|
32 | icpl = 1 |
---|
33 | if (local_comm_size .gt. 3) then |
---|
34 | if (local_comm_rank .ge. local_comm_size - 2) icpl = 0 |
---|
35 | end if |
---|
36 | |
---|
37 | call mpi_comm_split(local_comm, icpl, local_comm_rank, coupl_comm, kinfo) |
---|
38 | if(kinfo<0) call oasis_abort(comp_id, comp_name, & |
---|
39 | & "Error in mpi_comm_split: ", rcode=kinfo) |
---|
40 | |
---|
41 | if (icpl == 0) coupl_comm = MPI_COMM_NULL |
---|
42 | |
---|
43 | call oasis_set_couplcomm(coupl_comm, kinfo) |
---|
44 | if(kinfo<0) call oasis_abort(comp_id, comp_name, & |
---|
45 | & "Error in oasis_set_couplcomm: ", rcode=kinfo) |
---|
46 | |
---|
47 | if (icpl == 1) then |
---|
48 | |
---|
49 | call mpi_comm_size(coupl_comm, comm_size, kinfo) |
---|
50 | call mpi_comm_rank(coupl_comm, comm_rank, kinfo) |
---|
51 | print '(A,I0,A,I0)', "Sender: comm_rank = ",comm_rank, " of ",comm_size |
---|
52 | |
---|
53 | if ( mod(n_points,comm_size) /= 0) & |
---|
54 | & call oasis_abort(comp_id, comp_name, & |
---|
55 | & "Sender: comm_size has to divide n_points exaclty", rcode=kinfo) |
---|
56 | |
---|
57 | local_size=n_points/comm_size |
---|
58 | offset=comm_rank*local_size |
---|
59 | |
---|
60 | part_params(OASIS_Strategy) = OASIS_Apple |
---|
61 | part_params(OASIS_Offset) = offset |
---|
62 | part_params(OASIS_Length) = local_size |
---|
63 | call oasis_def_partition(part_id, part_params, kinfo) |
---|
64 | if(kinfo<0) call oasis_abort(comp_id, comp_name, & |
---|
65 | & "Error in oasis_def_partition: ", rcode=kinfo) |
---|
66 | print '(A,I0,A,I0)', "Sender rank(",comm_rank,"): part_id: ", part_id |
---|
67 | |
---|
68 | var_nodims=[1, 1] |
---|
69 | print '(A,I0,2A)', "Sender rank(",comm_rank,"): var_name: ", var_name |
---|
70 | call oasis_def_var(var_id, var_name, part_id, var_nodims, OASIS_OUT, & |
---|
71 | & [1], OASIS_REAL, kinfo) |
---|
72 | if(kinfo<0 .or. var_id<0) call oasis_abort(comp_id, comp_name, & |
---|
73 | & "Error in oasis_def_var: ", rcode=kinfo) |
---|
74 | print '(A,I0,A,I0)', "Sender rank(",comm_rank,"): var_id: ", var_id |
---|
75 | end if |
---|
76 | |
---|
77 | call oasis_enddef(kinfo) |
---|
78 | if(kinfo<0) call oasis_abort(comp_id, comp_name, & |
---|
79 | & "Error in oasis_enddef: ", rcode=kinfo) |
---|
80 | |
---|
81 | if (icpl == 1) then |
---|
82 | |
---|
83 | allocate(field(local_size)) |
---|
84 | field(:) = [(i+offset, i=1,local_size)] |
---|
85 | |
---|
86 | date=0 |
---|
87 | |
---|
88 | call oasis_put(var_id, date, field, kinfo) |
---|
89 | if(kinfo<0) call oasis_abort(comp_id, comp_name, & |
---|
90 | & "Error in oasis_put: ", rcode=kinfo) |
---|
91 | |
---|
92 | end if |
---|
93 | |
---|
94 | call oasis_terminate(kinfo) |
---|
95 | if(kinfo<0) call oasis_abort(comp_id, comp_name, & |
---|
96 | & "Error in oasis_terminate: ", rcode=kinfo) |
---|
97 | |
---|
98 | end program sender_apple |
---|