source: CPL/oasis3-mct/branches/OASIS3-MCT_5.0_branch/pyoasis/examples/6-apple_and_orange/fortran/receiver-orange.F90 @ 6331

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

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

File size: 3.0 KB
Line 
1program receiver
2   use mod_oasis
3   implicit none
4   integer :: i, kinfo, date
5   integer :: comp_id, part_id, var_id
6   integer :: local_comm, comm_size, comm_rank
7   integer, parameter :: n_points = 16
8   integer :: offset, local_size, n_segments
9   integer, allocatable, dimension(:) :: part_params
10   integer :: var_nodims(2)
11   character(len=8) :: comp_name = "receiver"
12   character(len=8) :: var_name = "FRECVATM"
13   real(kind=8), allocatable :: field(:)
14   real(kind=8) :: error, epsilon
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)', "Receiver: 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, comm_size, kinfo)
28   call mpi_comm_rank(local_comm, comm_rank, kinfo)
29   print '(A,I0,A,I0)', "Receiver: local_comm_rank = ",comm_rank, &
30      &                 " of ",comm_size
31
32   if ( mod(n_points,comm_size) /= 0) &
33      &  call oasis_abort(comp_id, comp_name, &
34      & "Receiver: comm_size has to divide n_points exaclty", rcode=kinfo)
35
36   local_size=n_points/comm_size
37   offset=comm_rank*local_size
38
39   n_segments = 1
40   allocate(part_params(OASIS_Orange_Params(n_segments)))
41   part_params(OASIS_Strategy) = OASIS_Orange
42   part_params(OASIS_Segments) = n_segments
43   do i = 0, n_segments-1
44      part_params(OASIS_Segments + 2*i + 1) = offset
45      part_params(OASIS_Segments + 2*i + 2) = local_size
46   end do
47
48   call oasis_def_partition(part_id, part_params, kinfo)
49   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
50      & "Error in oasis_def_partition: ", rcode=kinfo)
51   print '(A,I0,A,I0)', "Receiver rank(",comm_rank,"): part_id: ", part_id
52   deallocate(part_params)
53
54   var_nodims=[1, 1]
55   print '(A,I0,2A)', "Receiver rank(",comm_rank,"): var_name: ", var_name
56   call oasis_def_var(var_id, var_name, part_id, var_nodims, OASIS_IN, &
57      &               [1], OASIS_REAL, kinfo)
58   if(kinfo<0 .or. var_id<0) call oasis_abort(comp_id, comp_name, &
59      & "Error in oasis_def_var: ", rcode=kinfo)
60   print '(A,I0,A,I0)', "Receiver rank(",comm_rank,"): var_id: ", var_id
61
62   call oasis_enddef(kinfo)
63   if(kinfo<0) then
64      print *, "Error in oasis_enddef: ", kinfo
65      stop
66   endif
67
68   date=0
69
70   allocate(field(local_size))
71   field(:)=0
72
73   call oasis_get(var_id, date, field, kinfo)
74   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
75      & "Error in oasis_get: ", rcode=kinfo)
76
77   call oasis_terminate(kinfo)
78   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
79      & "Error in oasis_terminate: ", rcode=kinfo)
80
81   epsilon=1e-8
82   error=0
83   do i = 1, local_size
84      error=error+abs(field(i)-i-offset)
85   end do
86   if(error<epsilon)  print '(A,I0,A)', "Receiver rank(",comm_rank,"): Data received successfully"
87
88end program receiver
Note: See TracBrowser for help on using the repository browser.