source: CPL/oasis3-mct/branches/OASIS3-MCT_5.0_branch/pyoasis/examples/7-multiple-puts/fortran/receiver_one.F90 @ 6331

Last change on this file since 6331 was 6331, checked in by aclsce, 15 months ago

Moved oasis-mct_5.0 in oasis3-mct/branches directory.

File size: 4.1 KB
Line 
1program receiver_one
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=12) :: comp_name = "receiver_one"
10   character(len=10), dimension(2) :: var_name = ["FRECVATM_1","FRECVATM_2"]
11   real(kind=8) :: field(n_points), error
12   real(kind=8) :: epsilon = 1.e-8
13   integer :: ncpl
14   integer, dimension(:), allocatable :: cpl_freqs_1
15   integer, dimension(:), allocatable :: cpl_freqs_2
16   integer :: intra_one, intra_rank, intra_size
17
18   print '(2A)', "Component name: ", comp_name
19
20   call oasis_init_comp(comp_id, comp_name, kinfo)
21   if(kinfo<0) &
22      & call oasis_abort(comp_id, comp_name, &
23      & "Error in oasis_init_comp: ", rcode=kinfo)
24
25   part_params(OASIS_Strategy) = OASIS_Serial
26   part_params(OASIS_Length)   = n_points
27   call oasis_def_partition(part_id, part_params, kinfo)
28   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
29      & "Error in oasis_def_partition: ", rcode=kinfo)
30
31   var_nodims=[1, 1]
32
33   do  i = 1,2
34      print '(2A)', "Recv_one: defining var: ", var_name(i)
35      call oasis_def_var(var_id(i), var_name(i), part_id, var_nodims, OASIS_IN, &
36         &               [1], OASIS_REAL, kinfo)
37      if(kinfo<0 .or. var_id(i)<0) &
38         & call oasis_abort(comp_id, comp_name, &
39         & "Error in oasis_def_var: "//trim(var_name(i)), rcode=kinfo)
40   end do
41
42   call oasis_enddef(kinfo)
43   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
44      & "Error in oasis_enddef: ", rcode=kinfo)
45
46   call oasis_get_intracomm(intra_one, "sender-serial", kinfo)
47   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
48      &"Error in oasis_get_intracomm: ", rcode=kinfo)
49
50   call mpi_comm_size(intra_one, intra_size, kinfo)
51   call mpi_comm_rank(intra_one, intra_rank, kinfo)
52   print '(A,I0,A,I0)', "Recv_one intra_one: rank = ",intra_rank, " of ",intra_size
53
54   call oasis_get_ncpl(var_id(1), ncpl, kinfo)
55   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
56      & "Error in oasis_get_ncpl: ", rcode=kinfo)
57
58   allocate(cpl_freqs_1(ncpl))
59   call oasis_get_freqs(var_id(1), OASIS_IN, ncpl, cpl_freqs_1, kinfo)
60   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
61      & "Error in oasis_get_freqs: ", rcode=kinfo)
62
63   call oasis_get_ncpl(var_id(2), ncpl, kinfo)
64   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
65      & "Error in oasis_get_ncpl: ", rcode=kinfo)
66
67   allocate(cpl_freqs_2(ncpl))
68   call oasis_get_freqs(var_id(2), OASIS_IN, ncpl, cpl_freqs_2, kinfo)
69   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
70      & "Error in oasis_get_freqs: ", rcode=kinfo)
71
72   do date = 0, 43199
73
74      if (any(mod(date,cpl_freqs_1) == 0)) then
75         field(:) = 0
76         call oasis_get(var_id(1), date, field, kinfo)
77         if(kinfo<0) &
78            & call oasis_abort(comp_id, comp_name, &
79            & "Error in oasis_get: "//trim(var_name(1)), rcode=kinfo)
80
81         error = sum(abs(field(:)-date))
82         if(error<epsilon) then
83            print '(A,I0)', &
84               & "Recv_one: field 1 received successfully at time ",date
85         else
86            print '(A,I0,A,F6.0,A,F6.0)', &
87               & "Warning: Recv_one at time ",date," got ",field(1)," instead of ",1.*date
88         endif
89      end if
90
91      if (any(mod(date,cpl_freqs_2) == 0)) then
92         field(:) = 0
93         call oasis_get(var_id(2), date, field, kinfo)
94         if(kinfo<0) &
95            & call oasis_abort(comp_id, comp_name, &
96            & "Error in oasis_get: "//trim(var_name(1)), rcode=kinfo)
97
98         error = sum(abs(field(:)+date))
99         if(error<epsilon) then
100            print '(A,I0)', &
101               & "Recv_one: field 2 received successfully at time ",date
102         else
103            print '(A,I0,A,F6.0,A,F6.0)', &
104               & "Warning: Recv_one at time ",date," got ",field(1)," instead of ",-1.*date
105         endif
106      end if
107
108   end do
109
110   deallocate(cpl_freqs_1, cpl_freqs_2)
111
112   call oasis_terminate(kinfo)
113   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
114      & "Error in oasis_terminate: ", rcode=kinfo)
115
116end program receiver_one
Note: See TracBrowser for help on using the repository browser.