source: LMDZ6/branches/contrails/tools/netcdf95/Datasets/nf95_find_coord.f90 @ 5418

Last change on this file since 5418 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: 3.5 KB
Line 
1module nf95_find_coord_m
2
3  implicit none
4
5contains
6
7  subroutine nf95_find_coord(ncid, name, dimid, varid, std_name)
8
9    ! This procedure returns the name, dimension id or variable id of
10    ! the NetCDF coordinate with standard name "std_name", if such a
11    ! coordinate exists. The standard name is only used to know what
12    ! to search, it is not used for the search itself. The search
13    ! itself is done via a string match on the attribute "units". So
14    ! the NetCDF variable one looks for does not need to have the
15    ! attribute "std_name".
16
17    use netcdf, only: NF90_MAX_NAME, NF90_NOERR
18    use nf95_get_att_m, only: nf95_get_att
19    use nf95_inq_varid_m, only: nf95_inq_varid
20    use nf95_inquire_dimension_m, only: nf95_inquire_dimension
21    use nf95_inquire_m, only: nf95_inquire
22    use nf95_inquire_variable_m, only: nf95_inquire_variable
23
24    integer, intent(in):: ncid
25
26    character(len=*), intent(out), optional:: name ! blanks if not found
27    ! The actual character argument should normally have the length
28    ! "NF95_MAX_NAME".
29
30    integer, intent(out), optional:: dimid ! 0 if not found
31    integer, intent(out), optional:: varid ! 0 if not found
32
33    character(len=*), intent(in):: std_name
34    ! standard name : "plev", "latitude", "longitude" or "time"
35
36    ! Variables local to the procedure:
37
38    character(len=13), allocatable:: units(:)
39    logical exact ! "units" must be matched exactly
40
41    integer ncerr, nDimensions, dimid_local, varid_local
42    character(len=NF90_MAX_NAME) name_local
43    integer, allocatable:: dimids(:)
44    character(len=80) values
45    logical found
46
47    !----------------------------------------------
48
49    select case (std_name)
50    case("longitude")
51       allocate(units(1))
52       units(1)="degrees_east"
53       exact=.true.
54    case("latitude")
55       allocate(units(1))
56       units(1)="degrees_north"
57       exact=.true.
58    case("time")
59       allocate(units(1))
60       units(1)=" since"
61       exact=.false.
62    case("plev")
63       allocate(units(4))
64       units = ["Pa           ", "hPa          ", "millibar     ", &
65            "mbar         "]
66       exact = .true.
67    case default
68       print *, "nf95_find_coord: bad value of std_name"
69       print *, "std_name = ", std_name
70       stop 1
71    end select
72
73    call nf95_inquire(ncid, nDimensions)
74    dimid_local = 0
75    found = .false.
76
77    ! Loop on dimensions:
78    do while (.not. found .and. dimid_local < nDimensions)
79       dimid_local = dimid_local + 1
80       call nf95_inquire_dimension(ncid, dimid_local, name_local)
81       call nf95_inq_varid(ncid, name_local, varid_local, ncerr)
82       if (ncerr == NF90_NOERR) then
83          call nf95_inquire_variable(ncid, varid_local, dimids=dimids)
84          if (size(dimids) == 1) then
85             if (dimids(1) == dimid_local) then
86                ! We have found a coordinate
87                call nf95_get_att(ncid, varid_local, "units", values, ncerr)
88                if (ncerr == NF90_NOERR)then
89                   if (exact) then
90                      found = any(values == units)
91                   else
92                      found = index(values, trim(units(1))) /= 0
93                   end if
94                end if
95             end if
96          end if
97       end if
98    end do
99
100    if (found) then
101       if (present(name)) name = name_local
102       if (present(dimid)) dimid = dimid_local
103       if (present(varid)) varid = varid_local
104    else
105       if (present(name)) name = ""
106       if (present(dimid)) dimid = 0
107       if (present(varid)) varid = 0
108    end if
109
110  end subroutine nf95_find_coord
111
112end module nf95_find_coord_m
Note: See TracBrowser for help on using the repository browser.