source: CPL/oasis3-mct/branches/OASIS3-MCT_5.0_branch/pyoasis/examples/8-interoperability/fortran_and_C/sender-apple.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: 2.7 KB
Line 
1program sender_apple
2   use mod_oasis
3   implicit none
4   integer :: i, lat, kinfo, date
5   integer :: comp_id, part_id, var_id
6   integer :: n_points = 16
7   integer :: part_params(3), offset, local_size
8   integer :: local_comm, comm_size, comm_rank
9   integer :: var_nodims(2)
10   character(len=13) :: comp_name = "sender-apple"
11   character(len=8) :: var_name = "FSENDOCN"
12   real(kind=8), allocatable, dimension(:,:,:) :: bundle
13
14   call oasis_init_comp(comp_id, comp_name, kinfo)
15   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
16      & "Error in oasis_init_comp: ", rcode=kinfo)
17
18   call oasis_get_localcomm(local_comm, kinfo)
19   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
20      & "Error in oasis_get_localcomm: ", rcode=kinfo)
21
22   call mpi_comm_size(local_comm, comm_size, kinfo)
23   call mpi_comm_rank(local_comm, comm_rank, kinfo)
24
25   if ( comm_rank == 0 ) &
26      & print '(3A,I0)', "Component name: ", trim(comp_name), &
27      & " = Component ID: ", comp_id
28
29   local_size=n_points/comm_size
30   offset=comm_rank*local_size
31   if (comm_rank == comm_size - 1) &
32      & local_size = n_points - offset
33
34   part_params=[1, offset, local_size]
35   call oasis_def_partition(part_id, part_params, kinfo)
36   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
37      & "Error in oasis_def_partition: ", rcode=kinfo)
38   if ( comm_rank == 0 ) &
39      & print '(2A,I0)', trim(comp_name),": part_id: ", part_id
40
41   var_nodims=[2, 2]
42
43   call oasis_def_var(var_id, var_name, part_id, var_nodims, OASIS_OUT, &
44      &               [1], OASIS_REAL, kinfo)
45   if(kinfo<0 .or. var_id<0) call oasis_abort(comp_id, comp_name, &
46      & "Error in oasis_def_var: ", rcode=kinfo)
47   if ( comm_rank == 0 ) &
48      & print '(4A,I0)', trim(comp_name),": var_name: ", trim(var_name), &
49      & " = var_id: ", var_id
50
51   call oasis_enddef(kinfo)
52   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
53      & "Error in oasis_enddef: ", rcode=kinfo)
54
55   allocate(bundle(local_size,1,2))
56   lat = int(comm_rank/2)+1
57   do i = 1, 2
58      bundle(:,:,i)=i
59   end do
60   do i=1, local_size
61      bundle(i,:,:)=bundle(i,:,:)+(i+offset-(lat-1)*2*local_size)*100
62   end do
63   bundle(:,:,:)=bundle(:,:,:)+lat*10
64
65   print '(A,I2,A)', "On sender side rank", comm_rank, ": bundle(1)"
66   print '(4F6.0)', bundle(:,:,1)
67   print '(A,I2,A)', "On sender side rank", comm_rank, ": bundle(2)"
68   print '(4F6.0)', bundle(:,:,2)
69
70   date=0
71
72   call oasis_put(var_id, date, bundle, kinfo)
73   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
74      & "Error in oasis_put: ", rcode=kinfo)
75
76   deallocate(bundle)
77
78   call oasis_terminate(kinfo)
79   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
80      & "Error in oasis_terminate: ", rcode=kinfo)
81
82end program sender_apple
83
Note: See TracBrowser for help on using the repository browser.