[4918] | 1 | module nf95_get_missing_m |
---|
| 2 | |
---|
[5084] | 3 | use netcdf, only: nf90_noerr |
---|
[4918] | 4 | use nf95_get_att_m, only: nf95_get_att |
---|
| 5 | |
---|
| 6 | implicit none |
---|
| 7 | |
---|
| 8 | interface nf95_get_missing |
---|
| 9 | module procedure nf95_get_missing_real, nf95_get_missing_dble,& |
---|
| 10 | nf95_get_missing_short_int, nf95_get_missing_int, & |
---|
| 11 | nf95_get_missing_char |
---|
| 12 | end interface nf95_get_missing |
---|
| 13 | private |
---|
| 14 | public nf95_get_missing |
---|
| 15 | |
---|
| 16 | contains |
---|
| 17 | |
---|
| 18 | subroutine nf95_get_missing_real(ncid, varid, missing) |
---|
| 19 | |
---|
[5084] | 20 | use netcdf, only: NF90_FILL_REAL |
---|
[4918] | 21 | use typesizes, only: FourByteReal |
---|
| 22 | |
---|
| 23 | integer, intent(in):: ncid, varid |
---|
| 24 | real(kind=FourByteReal), intent(out):: missing ! missing or fill value |
---|
| 25 | |
---|
| 26 | ! Local: |
---|
| 27 | integer ncerr |
---|
| 28 | |
---|
| 29 | !------------------------------------------------------------------- |
---|
| 30 | |
---|
| 31 | call nf95_get_att(ncid, varid, name = "missing_value", values = missing, & |
---|
| 32 | ncerr = ncerr) |
---|
| 33 | |
---|
| 34 | if (ncerr /= nf90_noerr) then |
---|
| 35 | call nf95_get_att(ncid, varid, name = "_FillValue", values = missing, & |
---|
| 36 | ncerr = ncerr) |
---|
| 37 | if (ncerr /= nf90_noerr) missing = NF90_FILL_REAL |
---|
| 38 | end if |
---|
| 39 | |
---|
| 40 | end subroutine nf95_get_missing_real |
---|
| 41 | |
---|
| 42 | !************************************************************************** |
---|
| 43 | |
---|
| 44 | subroutine nf95_get_missing_dble(ncid, varid, missing) |
---|
| 45 | |
---|
[5084] | 46 | use netcdf, only: NF90_FILL_double |
---|
[4918] | 47 | use typesizes, only: EightByteReal |
---|
| 48 | |
---|
| 49 | integer, intent(in):: ncid, varid |
---|
| 50 | real(kind=EightByteReal), intent(out):: missing ! missing or fill value |
---|
| 51 | |
---|
| 52 | ! Local: |
---|
| 53 | integer ncerr |
---|
| 54 | |
---|
| 55 | !------------------------------------------------------------------- |
---|
| 56 | |
---|
| 57 | call nf95_get_att(ncid, varid, name = "missing_value", values = missing, & |
---|
| 58 | ncerr = ncerr) |
---|
| 59 | |
---|
| 60 | if (ncerr /= nf90_noerr) then |
---|
| 61 | call nf95_get_att(ncid, varid, name = "_FillValue", values = missing, & |
---|
| 62 | ncerr = ncerr) |
---|
| 63 | if (ncerr /= nf90_noerr) missing = NF90_FILL_double |
---|
| 64 | end if |
---|
| 65 | |
---|
| 66 | end subroutine nf95_get_missing_dble |
---|
| 67 | |
---|
| 68 | !************************************************************************** |
---|
| 69 | |
---|
| 70 | subroutine nf95_get_missing_short_int(ncid, varid, missing) |
---|
| 71 | |
---|
[5084] | 72 | use netcdf, only: NF90_FILL_short |
---|
[4918] | 73 | use typesizes, only: TwoByteInt |
---|
| 74 | |
---|
| 75 | integer, intent(in):: ncid, varid |
---|
| 76 | integer(kind = TwoByteInt), intent(out):: missing ! missing or fill value |
---|
| 77 | |
---|
| 78 | ! Local: |
---|
| 79 | integer ncerr |
---|
| 80 | |
---|
| 81 | !------------------------------------------------------------------- |
---|
| 82 | |
---|
| 83 | call nf95_get_att(ncid, varid, name = "missing_value", values = missing, & |
---|
| 84 | ncerr = ncerr) |
---|
| 85 | |
---|
| 86 | if (ncerr /= nf90_noerr) then |
---|
| 87 | call nf95_get_att(ncid, varid, name = "_FillValue", values = missing, & |
---|
| 88 | ncerr = ncerr) |
---|
| 89 | if (ncerr /= nf90_noerr) missing = NF90_FILL_short |
---|
| 90 | end if |
---|
| 91 | |
---|
| 92 | end subroutine nf95_get_missing_short_int |
---|
| 93 | |
---|
| 94 | !************************************************************************** |
---|
| 95 | |
---|
| 96 | subroutine nf95_get_missing_int(ncid, varid, missing) |
---|
| 97 | |
---|
[5084] | 98 | use netcdf, only: NF90_FILL_INT |
---|
[4918] | 99 | |
---|
| 100 | integer, intent(in):: ncid, varid |
---|
| 101 | integer, intent(out):: missing ! missing or fill value |
---|
| 102 | |
---|
| 103 | ! Local: |
---|
| 104 | integer ncerr |
---|
| 105 | |
---|
| 106 | !------------------------------------------------------------------- |
---|
| 107 | |
---|
| 108 | call nf95_get_att(ncid, varid, name = "missing_value", values = missing, & |
---|
| 109 | ncerr = ncerr) |
---|
| 110 | |
---|
| 111 | if (ncerr /= nf90_noerr) then |
---|
| 112 | call nf95_get_att(ncid, varid, name = "_FillValue", values = missing, & |
---|
| 113 | ncerr = ncerr) |
---|
| 114 | if (ncerr /= nf90_noerr) missing = NF90_FILL_INT |
---|
| 115 | end if |
---|
| 116 | |
---|
| 117 | end subroutine nf95_get_missing_int |
---|
| 118 | |
---|
| 119 | !************************************************************************** |
---|
| 120 | |
---|
| 121 | subroutine nf95_get_missing_char(ncid, varid, missing) |
---|
| 122 | |
---|
[5084] | 123 | use netcdf, only: NF90_FILL_char |
---|
[4918] | 124 | |
---|
| 125 | integer, intent(in):: ncid, varid |
---|
| 126 | |
---|
| 127 | character, intent(out):: missing |
---|
| 128 | ! Missing or fill value. Note that the missing_value attribute for |
---|
| 129 | ! a NetCDF character variable should be a single character |
---|
| 130 | |
---|
| 131 | ! Local: |
---|
| 132 | integer ncerr |
---|
| 133 | |
---|
| 134 | !------------------------------------------------------------------- |
---|
| 135 | |
---|
| 136 | call nf95_get_att(ncid, varid, name = "missing_value", values = missing, & |
---|
| 137 | ncerr = ncerr) |
---|
| 138 | if (ncerr /= nf90_noerr) then |
---|
| 139 | call nf95_get_att(ncid, varid, name = "_FillValue", values = missing, & |
---|
| 140 | ncerr = ncerr) |
---|
| 141 | if (ncerr /= nf90_noerr) missing = NF90_FILL_char |
---|
| 142 | end if |
---|
| 143 | |
---|
| 144 | end subroutine nf95_get_missing_char |
---|
| 145 | |
---|
| 146 | end module nf95_get_missing_m |
---|