source: CPL/oasis3-mct_5.0/pyoasis/examples/11-test-interpolation/fortran/receiver.F90 @ 6328

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

First import of oasis3-mct_5.0 (from oasis git server, branch OASIS3-MCT_5.0)

File size: 3.3 KB
Line 
1program receiver
2   use mpi
3   use mod_oasis
4   use netcdf
5   implicit none
6   integer :: i, j, k, kinfo, date
7   integer :: comp_id, part_id, var_id
8   integer :: part_params(OASIS_Serial_Params)
9   integer :: var_nodims(2)
10   character(len=8) :: comp_name = "receiver"
11   character(len=8) :: var_name = "FRECVANA"
12   real(kind=8) :: error, epsilon
13   integer, parameter :: nx_global = 144, ny_global = 143
14   real(kind=8) ::  bundle(nx_global,ny_global,2)
15   real(kind=8) ::  expected(nx_global,ny_global,2)
16   integer :: n_points = nx_global*ny_global
17   integer :: ncid, varid
18   real(kind=8) :: lon(nx_global,ny_global), lat(nx_global,ny_global)
19   integer :: imsk(nx_global,ny_global)
20   real(kind=8) :: dp_conv
21   logical :: success
22
23   call oasis_init_comp(comp_id, comp_name, kinfo)
24   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
25      & "Error in oasis_init_comp: ", rcode=kinfo)
26   print '(A,I0)', "receiver: Component ID: ", comp_id
27
28   kinfo = nf90_open('grids.nc',NF90_NOWRITE,ncid)
29   kinfo = nf90_inq_varid(ncid,'bggd.lon',varid)
30   kinfo = nf90_get_var(ncid,varid,lon)
31   kinfo = nf90_inq_varid(ncid,'bggd.lat',varid)
32   kinfo = nf90_get_var(ncid,varid,lat)
33   kinfo = nf90_close(ncid)
34
35   kinfo = nf90_open('masks.nc',NF90_NOWRITE,ncid)
36   kinfo = nf90_inq_varid(ncid,'bggd.msk',varid)
37   kinfo = nf90_get_var(ncid,varid,imsk)
38   kinfo = nf90_close(ncid)
39
40   part_params(OASIS_Strategy) = OASIS_Serial
41   part_params(OASIS_Length)   = n_points
42   call oasis_def_partition(part_id, part_params, kinfo)
43   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
44      & "Error in oasis_def_partition: ", rcode=kinfo)
45
46   var_nodims=[1, 2]
47   print '(2A)', "Receiver: var_name: ", var_name
48   call oasis_def_var(var_id, var_name, part_id, var_nodims, OASIS_IN, &
49      &              [1], OASIS_DOUBLE, kinfo)
50   if(kinfo<0 .or. var_id<0) call oasis_abort(comp_id, comp_name, &
51      & "Error in oasis_def_var: ", rcode=kinfo)
52
53   call oasis_enddef(kinfo)
54   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
55      & "Error in oasis_enddef: ", rcode=kinfo)
56
57   date=0
58
59   bundle(:,:,:)=0
60
61   call oasis_get(var_id, date, bundle, kinfo)
62   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
63      & "Error in oasis_get: ", rcode=kinfo)
64
65   call oasis_terminate(kinfo)
66   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
67      & "Error in oasis_terminate: ", rcode=kinfo)
68
69   dp_conv = atan(1.)/45.0
70   do j = 1, ny_global
71      do i = 1, nx_global
72         expected(i,j,1) = 2.0 + (sin(2.*lat(i,j)*dp_conv))**4 * &
73            & cos(4.*lon(i,j)*dp_conv)
74         expected(i,j,2) = 2.0 - cos(atan(1.0)*4.* &
75            & (acos(cos(lon(i,j)*dp_conv)*cos(lat(i,j)*dp_conv))/ &
76            & (1.2*atan(1.)*4)))
77      end do
78   end do
79
80   epsilon=1.e-3
81   success = .true.
82   do k = 1, 2
83      error=0.
84      do j = 1, ny_global
85         do i = 1, nx_global
86            if (imsk(i,j) == 0) &
87               & error = error + abs((bundle(i,j,k)-expected(i,j,k))/expected(i,j,k))
88         end do
89      end do
90      success = success .and. (error/dble(n_points) < epsilon)
91      if (success) then
92         print '(A,I0,A)',"Receiver: Data for bundle ",k," is ok"
93      else
94         print '(A,I0,A,E12.5)', "Receiver: Error for bundle ",k," is ",error
95      end if
96   end do
97
98   if(success) print '(A)', "Receiver: Data received successfully"
99
100end program receiver
Note: See TracBrowser for help on using the repository browser.