source: CPL/oasis3-mct/branches/OASIS3-MCT_5.0_branch/pyoasis/examples/11-test-interpolation/fortran/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: 3.0 KB
Line 
1program sender_apple
2   use mpi
3   use mod_oasis
4   use netcdf
5   implicit none
6   integer :: i, kinfo, date
7   integer :: comp_id, part_id, var_id
8   integer :: part_params(OASIS_Apple_Params), offset, local_size
9   integer :: local_comm, comm_size, comm_rank
10   integer :: var_nodims(2)
11   character(len=13) :: comp_name = "sender-apple"
12   character(len=8) :: var_name = "FSENDANA"
13   real(kind=8), allocatable :: bundle(:,:)
14   integer, parameter :: nx_global = 362, ny_global = 294
15   integer :: n_points = nx_global*ny_global
16   integer :: ncid, varid
17   real(kind=8) :: lon(nx_global,ny_global), lat(nx_global,ny_global)
18   integer :: ll_i, ll_j
19   real(kind=8) :: dp_conv
20
21   call oasis_init_comp(comp_id, comp_name, kinfo)
22   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
23      & "Error in oasis_init_comp: ", rcode=kinfo)
24   print '(A,I0)', "sender-apple: Component ID: ", comp_id
25
26   call oasis_get_localcomm(local_comm, kinfo)
27   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
28      & "Error in oasis_get_localcomm: ", rcode=kinfo)
29
30   call mpi_comm_size(local_comm, comm_size, kinfo)
31   call mpi_comm_rank(local_comm, comm_rank, kinfo)
32
33   kinfo = nf90_open('grids.nc',NF90_NOWRITE,ncid)
34   kinfo = nf90_inq_varid(ncid,'nogt.lon',varid)
35   kinfo = nf90_get_var(ncid,varid,lon)
36   kinfo = nf90_inq_varid(ncid,'nogt.lat',varid)
37   kinfo = nf90_get_var(ncid,varid,lat)
38   kinfo = nf90_close(ncid)
39
40   local_size=n_points/comm_size
41   offset=comm_rank*local_size
42   if (comm_rank == comm_size - 1) &
43      & local_size = n_points - offset
44
45   part_params(OASIS_Strategy) = OASIS_Apple
46   part_params(OASIS_Offset)   = offset
47   part_params(OASIS_Length)   = local_size
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
52   var_nodims=[1, 2]
53   print '(A,I0,2A)', "Sender rank(",comm_rank,"): var_name: ", var_name
54   call oasis_def_var(var_id, var_name, part_id, var_nodims, OASIS_OUT, &
55      &               [1], OASIS_DOUBLE, kinfo)
56   if(kinfo<0 .or. var_id<0) call oasis_abort(comp_id, comp_name, &
57      & "Error in oasis_def_var: ", rcode=kinfo)
58
59   call oasis_enddef(kinfo)
60   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
61      & "Error in oasis_enddef: ", rcode=kinfo)
62
63   dp_conv = atan(1.0)/45.0
64   allocate(bundle(local_size,2))
65   do i = 1, local_size
66      ll_j = int((offset+i-1)/nx_global)+1
67      ll_i = mod(offset+i-1,nx_global)+1
68      bundle(i,1) = 2.0 + (sin(2.*lat(ll_i,ll_j)*dp_conv))**4 * &
69         & cos(4.*lon(ll_i,ll_j)*dp_conv)
70      bundle(i,2) = 2.0 - cos(atan(1.0)*4.* &
71          & (acos(cos(lon(ll_i,ll_j)*dp_conv)*cos(lat(ll_i,ll_j)*dp_conv))/ &
72          & (1.2*atan(1.)*4)))
73   end do
74
75   date=0
76
77   call oasis_put(var_id, date, bundle, kinfo)
78   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
79      & "Error in oasis_put: ", rcode=kinfo)
80
81   call oasis_terminate(kinfo)
82   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
83      & "Error in oasis_terminate: ", rcode=kinfo)
84
85end program sender_apple
Note: See TracBrowser for help on using the repository browser.