source: CONFIG_DEVT/IPSLCM6.5_work_ENSEMBLES/oasis3-mct/pyoasis/examples/7-multiple-puts/fortran/sender-serial.F90

Last change on this file 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: 3.8 KB
Line 
1program sender_serial
2   use mod_oasis
3   implicit none
4   integer :: i, kinfo, date
5   integer :: comp_id, part_id, var_id(2)
6   integer, parameter :: n_points = 1
7   integer :: part_params(OASIS_Serial_Params)
8   integer :: var_nodims(2)
9   character(len=13) :: comp_name = "sender-serial"
10   character(len=10), dimension(2) :: var_name = ["FSENDOCN_1","FSENDOCN_2"]
11   real(kind=8) :: field(n_points)
12   integer :: ncpl
13   integer, dimension(:), allocatable :: cpl_freqs
14   logical :: ll_rst
15   integer :: intra_one, intra_rank, intra_size
16   integer :: inter_two, inter_rank, inter_size
17
18   print '(2A)', "Component name: ", comp_name
19
20   call oasis_init_comp(comp_id, comp_name, kinfo)
21   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
22      & "Error in oasis_init_comp: ", rcode=kinfo)
23
24   part_params(OASIS_Strategy) = OASIS_Serial
25   part_params(OASIS_Length)   = n_points
26   call oasis_def_partition(part_id, part_params, kinfo)
27   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
28      & "Error in oasis_def_partition: ", rcode=kinfo)
29
30   var_nodims=[1, 1]
31   do i = 1, 2
32      print '(2A)', "Sender: var_name: ", var_name(i)
33      call oasis_def_var(var_id(i), var_name(i), part_id, var_nodims, OASIS_OUT, &
34         &               [1], OASIS_REAL, kinfo)
35      if(kinfo<0 .or. var_id(i)<0) &
36         & call oasis_abort(comp_id, comp_name, &
37         & "Error in oasis_def_var: "//trim(var_name(i)), rcode=kinfo)
38   end do
39
40   call oasis_enddef(kinfo)
41   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
42      &"Error in oasis_enddef: ", rcode=kinfo)
43
44   call oasis_get_intracomm(intra_one,"receiver_one",kinfo)
45   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
46      &"Error in oasis_get_intracomm: ", rcode=kinfo)
47
48   call mpi_comm_size(intra_one, intra_size, kinfo)
49   call mpi_comm_rank(intra_one, intra_rank, kinfo)
50   print '(A,I0,A,I0)', "Sender intra_one: rank = ",intra_rank, " of ",intra_size
51
52   call oasis_get_intercomm(inter_two,"receiver_two",kinfo)
53   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
54      &"Error in oasis_get_intercomm: ", rcode=kinfo)
55
56   call mpi_comm_size(inter_two, inter_size, kinfo)
57   call mpi_comm_rank(inter_two, inter_rank, kinfo)
58   print '(A,I0,A,I0)', "Sender inter_two: rank = ",inter_rank, " of ",inter_size
59
60   call oasis_get_ncpl(var_id(2), ncpl, kinfo)
61   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
62      & "Error in oasis_get_ncpl: ", rcode=kinfo)
63
64   allocate(cpl_freqs(ncpl))
65
66   call oasis_get_freqs(var_id(2), OASIS_OUT, ncpl, cpl_freqs, kinfo)
67   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
68      & "Error in oasis_get_freqs: ", rcode=kinfo)
69
70   print '(3A,2I6)', &
71      & "Sender: coupling frequencies for ",trim(var_name(2)), " are ", cpl_freqs(:)
72
73   do date= 0, 43199
74
75      call oasis_put_inquire(var_id(1), date, kinfo)
76      if(kinfo<0) call oasis_abort(comp_id, comp_name, &
77         & "Error in oasis_put_inquire: ", rcode=kinfo)
78
79      if ( kinfo == OASIS_Sent) then
80
81         field(:) = date
82
83         call oasis_put(var_id(1), date, field, kinfo)
84         if(kinfo<0) &
85            & call oasis_abort(comp_id, comp_name, &
86            & "Error in oasis_put: "//trim(var_name(1)), rcode=kinfo)
87
88      endif
89
90      if (any(mod(date,cpl_freqs) == 0)) then
91
92         call oasis_set_debug(2)
93
94         field(:) = (-1.) * date
95         ll_rst = mod(date, 10800) == 0
96         call oasis_put(var_id(2), date, field, kinfo, write_restart=ll_rst)
97         if(kinfo<0) &
98            & call oasis_abort(comp_id, comp_name, &
99            & "Error in oasis_put: "//trim(var_name(2)), rcode=kinfo)
100
101         call oasis_set_debug(0)
102
103      end if
104
105   end do
106
107   deallocate(cpl_freqs)
108
109   call oasis_terminate(kinfo)
110   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
111      & "Error in oasis_terminate: ", rcode=kinfo)
112
113end program sender_serial
Note: See TracBrowser for help on using the repository browser.