source: CPL/oasis3-mct/branches/OASIS3-MCT_5.0_branch/pyoasis/examples/8-interoperability/fortran_and_python/receiver-serial.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: 1.8 KB
Line 
1program receiver_serial
2   use mod_oasis
3   implicit none
4   integer :: kinfo, date
5   integer :: comp_id, part_id, var_id
6   character(len=15) :: comp_name = "receiver-serial"
7   integer :: n_points = 16
8   integer :: part_params(3)
9   integer ::  var_nodims(2)
10   character(len=8) :: var_name = "FRECVICE"
11   real(kind=4), allocatable, dimension(:) :: field
12
13   call oasis_init_comp(comp_id, comp_name, kinfo)
14   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
15      & "Error in oasis_init_comp: ", rcode=kinfo)
16
17   print '(3A,I0)', "Component name: ", trim(comp_name), &
18      & " = Component ID: ", comp_id
19
20   part_params=[0, 0, n_points]
21   call oasis_def_partition(part_id, part_params, kinfo)
22   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
23      & "Error in oasis_def_partition: ", rcode=kinfo)
24   print '(2A,I0)', trim(comp_name),": part_id: ", part_id
25
26   var_nodims=[1, 1]
27
28   call oasis_def_var(var_id, var_name, part_id, var_nodims, OASIS_IN, &
29      &               [1], OASIS_REAL, kinfo)
30   if(kinfo<0 .or. var_id<0) call oasis_abort(comp_id, comp_name, &
31      & "Error in oasis_def_var: ", rcode=kinfo)
32   print '(4A,I0)', trim(comp_name),": var_name: ", trim(var_name), &
33      & " = var_id: ", var_id
34
35   call oasis_enddef(kinfo)
36   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
37      & "Error in oasis_enddef: ", rcode=kinfo)
38
39   allocate(field(n_points))
40   field(:) = 999.
41
42   date = 0
43
44   call oasis_get(var_id, date, field, kinfo)
45   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
46      & "Error in oasis_get: ", rcode=kinfo)
47
48   print *, 'FRECVICE = '
49   print '(8F6.0)', field
50
51   deallocate(field)
52
53   call oasis_terminate(kinfo)
54   if(kinfo<0) call oasis_abort(comp_id, comp_name, &
55      & "Error in oasis_terminate: ", rcode=kinfo)
56
57end program receiver_serial
Note: See TracBrowser for help on using the repository browser.