source: LMDZ6/trunk/tools/netcdf95/Groups/nf95_inq_grps.f90 @ 5410

Last change on this file since 5410 was 5084, checked in by Laurent Fairhead, 5 months ago

Reverting to r4065. Updating fortran standard broke too much stuff. Will do it by smaller chunks
AB, LF

File size: 2.3 KB
Line 
1module nf95_inq_grps_m
2
3  implicit none
4
5contains
6
7    Integer(C_INT) function inq_numgrps(ncid, numgrps)
8    use, intrinsic:: ISO_C_BINDING
9    Integer(C_INT) numgrps, ncid
10    Integer(C_INT) cstatus
11      Interface
12        Integer(C_INT) Function nc_inq_grps(ncid, numgrps, ncids) BIND(C)
13          import c_int,c_ptr, c_null_ptr
14          Integer(C_INT), VALUE :: ncid
15          Integer(C_INT), Intent(OUT):: numgrps
16          Type(C_PTR), VALUE :: ncids
17        End Function nc_inq_grps
18      End Interface
19   
20    cstatus = nc_inq_grps(ncid, numgrps, c_null_ptr );
21    inq_numgrps = cstatus
22  end function inq_numgrps
23 
24  subroutine nf95_inq_grps(ncid, ncids, ncerr)
25
26    use, intrinsic:: ISO_C_BINDING
27
28    use netcdf, only: nf90_noerr
29
30    use nc_constants, only: nc_noerr
31    use nf95_abort_m, only: nf95_abort
32
33    integer, intent(in):: ncid ! can be the file id or a group id
34    integer, allocatable, intent(out):: ncids(:)
35    integer, intent(out), optional:: ncerr
36
37    ! Local:
38
39    Integer(C_INT) numgrps, cstatus, cncid
40    integer(C_INT), allocatable:: cncids(:)
41   
42    Interface
43!       Integer(C_INT) Function inq_numgrps(ncid, numgrps) BIND(C)
44!         import c_int
45!         Integer(C_INT), VALUE:: ncid
46!         Integer(C_INT), Intent(OUT):: numgrps
47!       End Function inq_numgrps
48
49       Integer(C_INT) Function nc_inq_grps(ncid, numgrps, ncids) BIND(C)
50         import c_int
51         Integer(C_INT), VALUE:: ncid
52         Integer(C_INT), Intent(OUT):: numgrps
53         Integer(C_INT), Intent(OUT):: ncids(*)
54       End Function nc_inq_grps
55    End Interface
56
57    !------------------------------------------------------------
58
59    cncid = int(ncid, c_int)
60    cstatus = inq_numgrps(cncid, numgrps)
61    if (cstatus /= nc_noerr) call nf95_abort("nf95_inq_grps -- inq_numgrps", &
62         int(cstatus), ncid)
63
64    if (numgrps >= 1) then
65       allocate(cncids(numgrps))
66       cstatus = nc_inq_grps(cncid, numgrps, cncids)
67
68       if (present(ncerr)) then
69          ncerr = cstatus
70       else
71          if (cstatus /= nc_noerr) call &
72               nf95_abort("nf95_inq_grps -- nc_inq_grps", int(cstatus), ncid)
73       end if
74
75       if (cstatus == nf90_noerr) ncids = cncids
76    else
77       allocate(ncids(0))
78    end if
79
80  end subroutine nf95_inq_grps
81
82end module nf95_inq_grps_m
Note: See TracBrowser for help on using the repository browser.