source: CONFIG_DEVT/IPSLCM6.5_work_ENSEMBLES/oasis3-mct/pyoasis/examples/3-box/fortran/sender-box.F90 @ 5725

Last change on this file since 5725 was 5725, checked in by aclsce, 3 years ago

Added new oasis3-MCT version to be used to handle ensembles simulations with XIOS.

File size: 4.1 KB
Line 
1program sender_box
2   use mod_oasis
3   implicit none
4   integer :: i, kinfo, date
5   INTEGER :: comp_id, part_id, var_id, var_id2
6   integer :: n_points = 16
7   integer :: part_params(OASIS_Box_Params), offsets(4)
8   integer :: local_comm, local_size, comm_size, comm_rank
9   integer :: var_nodims(2)
10   character(len=13) :: comp_name = "sender-box"
11   character(len=8) :: var_name = "FSENDOCN"
12   character(len=8) :: var_name2 = "NOTANAME"
13   real, allocatable :: field(:)
14   integer :: intra_comm, intra_rank, intra_size
15   integer :: inter_comm, inter_rank, inter_size, inter_rsize
16
17   print '(2A)', "Component name: ", comp_name
18
19   call oasis_init_comp(comp_id, comp_name, kinfo)
20   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
21      & "Error in oasis_init_comp: ", rcode=kinfo)
22   print '(A,I0)', "Sender: Component ID: ", comp_id
23
24   call oasis_get_localcomm(local_comm, kinfo)
25   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
26      & "Error in oasis_get_localcomm: ", rcode=kinfo)
27
28   call mpi_comm_size(local_comm, comm_size, kinfo)
29   call mpi_comm_rank(local_comm, comm_rank, kinfo)
30   print '(A,I0,A,I0)', "Sender: rank = ",comm_rank, " of ",comm_size
31
32   if ( comm_size /= 4) call oasis_abort(comp_id, comp_name, &
33      & "Sender: comm_size has to be 4 for this example", rcode=kinfo)
34
35   local_size = n_points/comm_size
36
37   offsets=[0, 2, 8, 10]
38
39   part_params(OASIS_Strategy) = OASIS_Box
40   part_params(OASIS_Offset)   = offsets(comm_rank+1)
41   part_params(OASIS_SizeX)    = 2
42   part_params(OASIS_SizeY)    = 2
43   part_params(OASIS_LdX)      = 4
44   call oasis_def_partition(part_id, part_params, kinfo)
45   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
46      & "Error in oasis_def_partition: ", rcode=kinfo)
47   print '(A,I0,A,I0)', "Sender rank(",comm_rank,"): part_id: ", part_id
48
49   var_nodims=[1, 1]
50   print '(A,I0,2A)', "Sender rank(",comm_rank,"): var_name: ", var_name
51   call oasis_def_var(var_id, var_name, part_id, var_nodims, OASIS_OUT, &
52      &               [1], OASIS_REAL, kinfo)
53   if(kinfo<0 .or. var_id<0) call oasis_abort(comp_id, comp_name, &
54      & "Error in oasis_def_var: ", rcode=kinfo)
55   print '(A,I0,A,I0)', "Sender rank(",comm_rank,"): var_id: ", var_id
56
57   print '(A,I0,2A)', "Sender rank(",comm_rank,"): var_name: ", var_name2
58   call oasis_def_var(var_id2, var_name2, part_id, var_nodims, OASIS_OUT, &
59      &               [1], OASIS_REAL, kinfo)
60   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
61      & "Error in oasis_def_var: ", rcode=kinfo)
62   if(var_id2<0) then
63      print '(A,I0,A,I0,A)', "Sender rank(",comm_rank,"): var_id: ", var_id2, " is not active"
64   else
65      print '(A,I0,A,I0)', "Sender rank(",comm_rank,"): var_id: ", var_id2     
66   end if
67   
68   call oasis_enddef(kinfo)
69   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
70      & "Error in oasis_enddef: ", rcode=kinfo)
71
72   call oasis_get_intracomm(intra_comm,"receiver",kinfo)
73   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
74      &"Error in oasis_get_intracomm: ", rcode=kinfo)
75
76   call mpi_comm_size(intra_comm, intra_size, kinfo)
77   call mpi_comm_rank(intra_comm, intra_rank, kinfo)
78   print '(A,I0,A,I0)', "Sender intra_comm: rank = ",intra_rank, " of ",intra_size
79
80   call oasis_get_intercomm(inter_comm,"receiver",kinfo)
81   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
82      &"Error in oasis_get_intercomm: ", rcode=kinfo)
83
84   call mpi_comm_size(inter_comm, inter_size, kinfo)
85   call mpi_comm_rank(inter_comm, inter_rank, kinfo)
86   call mpi_comm_remote_size(inter_comm, inter_rsize, kinfo)
87   print '(A,I0,A,I0,A,I0)', "Sender inter_comm: rank = ",inter_rank, " of ",inter_size, &
88      & " Remote size = ",inter_rsize
89
90   allocate(field(local_size))
91   field=[offsets(comm_rank+1)+1, offsets(comm_rank+1)+2,&
92      &   offsets(comm_rank+1)+5, offsets(comm_rank+1)+6]
93
94   date=0
95
96   call oasis_put(var_id, date, field, kinfo)
97   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
98      & "Error in oasis_put: ", rcode=kinfo)
99
100   call oasis_terminate(kinfo)
101   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
102      & "Error in oasis_terminate: ", rcode=kinfo)
103
104end program sender_box
Note: See TracBrowser for help on using the repository browser.