source: LMDZ6/trunk/tools/netcdf95/Attributes/nf95_get_att.f90 @ 5130

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

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

File size: 7.2 KB
RevLine 
[4918]1module nf95_get_att_m
2
3  use nf95_abort_m, only: nf95_abort
[5084]4  use netcdf, only: nf90_get_att, nf90_noerr
[4918]5  use nf95_inquire_attribute_m, only: nf95_inquire_attribute
6  use nf95_constants, only: nf95_noerr
7
8  implicit none
9
10  interface nf95_get_att
11     module procedure nf95_get_att_text, nf95_get_att_one_twoByteInt, &
12          nf95_get_att_one_FourByteInt, nf95_get_att_one_FourByteReal, &
13          nf95_get_att_one_eightByteReal
14
15     ! The difference between the specific procedures is the type of
16     ! argument "values".
17  end interface nf95_get_att
18
19  private
20  public nf95_get_att
21
22contains
23
24  subroutine nf95_get_att_text(ncid, varid, name, values, ncerr)
25
26    integer, intent(in):: ncid, varid
27    character(len = *), intent(in):: name
28    character(len = *), intent(out):: values
29    integer, intent(out), optional:: ncerr
30
31    ! Variables local to the procedure:
32    integer ncerr_not_opt
33    integer att_len
34
35    !-------------------
36
37    ! Check that the length of "values" is large enough:
38    call nf95_inquire_attribute(ncid, varid, name, nclen=att_len, &
39         ncerr=ncerr_not_opt)
40    if (ncerr_not_opt == nf90_noerr) then
41       if (len(values) < att_len) then
42          print *, "nf95_get_att_text"
43          print *, "varid = ", varid
44          print *, "attribute name: ", name
45          print *, 'length of "values" is not large enough'
46          print *, "len(values) = ", len(values)
47          print *, "number of characters in attribute: ", att_len
48          stop 1
49       end if
50    end if
51
52    values = "" ! useless in NetCDF version 3.6.2 or better
53    ncerr_not_opt = nf90_get_att(ncid, varid, name, values)
54    if (present(ncerr)) then
55       ncerr = ncerr_not_opt
56    else
57       if (ncerr_not_opt /= nf95_noerr) call nf95_abort("nf95_get_att_text " &
58            // trim(name), ncerr_not_opt, ncid, varid)
59    end if
60
61    if (att_len >= 1 .and. ncerr_not_opt == nf90_noerr) then
62       ! Remove null terminator, if any:
63       if (iachar(values(att_len:att_len)) == 0) values(att_len:att_len) = " "
64    end if
65
66  end subroutine nf95_get_att_text
67
68  !***********************
69
70  subroutine nf95_get_att_one_TwoByteInt(ncid, varid, name, values, ncerr)
71
72    use typesizes, only: TwoByteInt
73
74    integer, intent(in):: ncid, varid
75    character(len = *), intent(in):: name
76    integer (kind = TwoByteInt), intent(out):: values
77    integer, intent(out), optional:: ncerr
78
79    ! Variables local to the procedure:
80    integer ncerr_not_opt
81    integer att_len
82
83    !-------------------
84
85    ! Check that the attribute contains a single value:
86    call nf95_inquire_attribute(ncid, varid, name, nclen=att_len, &
87         ncerr=ncerr_not_opt)
88    if (ncerr_not_opt == nf90_noerr) then
89       if (att_len /= 1) then
90          print *, "nf95_get_att_one_TwoByteInt"
91          print *, "varid = ", varid
92          print *, "attribute name: ", name
93          print *, 'the attribute does not contain a single value'
94          print *, "number of values in attribute: ", att_len
95          stop 1
96       end if
97    end if
98
99    ncerr_not_opt = nf90_get_att(ncid, varid, name, values)
100    if (present(ncerr)) then
101       ncerr = ncerr_not_opt
102    else
103       if (ncerr_not_opt /= nf95_noerr) call &
104            nf95_abort("nf95_get_att_one_TwoByteInt " // trim(name), &
105            ncerr_not_opt, ncid, varid)
106    end if
107
108  end subroutine nf95_get_att_one_TwoByteInt
109
110  !***********************
111
112  subroutine nf95_get_att_one_FourByteInt(ncid, varid, name, values, ncerr)
113
114    use typesizes, only: FourByteInt
115
116    integer, intent(in):: ncid, varid
117    character(len = *), intent(in):: name
118    integer (kind = FourByteInt), intent(out):: values
119    integer, intent(out), optional:: ncerr
120
121    ! Variables local to the procedure:
122    integer ncerr_not_opt
123    integer att_len
124
125    !-------------------
126
127    ! Check that the attribute contains a single value:
128    call nf95_inquire_attribute(ncid, varid, name, nclen=att_len, &
129         ncerr=ncerr_not_opt)
130    if (ncerr_not_opt == nf90_noerr) then
131       if (att_len /= 1) then
132          print *, "nf95_get_att_one_FourByteInt"
133          print *, "varid = ", varid
134          print *, "attribute name: ", name
135          print *, 'the attribute does not contain a single value'
136          print *, "number of values in attribute: ", att_len
137          stop 1
138       end if
139    end if
140
141    ncerr_not_opt = nf90_get_att(ncid, varid, name, values)
142    if (present(ncerr)) then
143       ncerr = ncerr_not_opt
144    else
145       if (ncerr_not_opt /= nf95_noerr) call &
146            nf95_abort("nf95_get_att_one_FourByteInt " // trim(name), &
147            ncerr_not_opt, ncid, varid)
148    end if
149
150  end subroutine nf95_get_att_one_FourByteInt
151
152  !***********************
153
154  subroutine nf95_get_att_one_FourByteReal(ncid, varid, name, values, ncerr)
155
156    use typesizes, only: FourByteReal
157
158    integer, intent(in):: ncid, varid
159    character(len = *), intent(in):: name
160    real (kind = FourByteReal), intent(out):: values
161    integer, intent(out), optional:: ncerr
162
163    ! Variables local to the procedure:
164    integer ncerr_not_opt
165    integer att_len
166
167    !-------------------
168
169    ! Check that the attribute contains a single value:
170    call nf95_inquire_attribute(ncid, varid, name, nclen=att_len, &
171         ncerr=ncerr_not_opt)
172    if (ncerr_not_opt == nf90_noerr) then
173       if (att_len /= 1) then
174          print *, "nf95_get_att_one_Fourbytereal"
175          print *, "varid = ", varid
176          print *, "attribute name: ", name
177          print *, 'the attribute does not contain a single value'
178          print *, "number of values in attribute: ", att_len
179          stop 1
180       end if
181    end if
182
183    ncerr_not_opt = nf90_get_att(ncid, varid, name, values)
184    if (present(ncerr)) then
185       ncerr = ncerr_not_opt
186    else
187       if (ncerr_not_opt /= nf95_noerr) call &
188            nf95_abort("nf95_get_att_one_Fourbytereal " // trim(name), &
189            ncerr_not_opt, ncid, varid)
190    end if
191
192  end subroutine nf95_get_att_one_Fourbytereal
193
194  !***********************
195
196  subroutine nf95_get_att_one_EightByteReal(ncid, varid, name, values, ncerr)
197
198    use typesizes, only: EightByteReal
199
200    integer, intent(in):: ncid, varid
201    character(len = *), intent(in):: name
202    real (kind = EightByteReal), intent(out):: values
203    integer, intent(out), optional:: ncerr
204
205    ! Variables local to the procedure:
206    integer ncerr_not_opt
207    integer att_len
208
209    !-------------------
210
211    ! Check that the attribute contains a single value:
212    call nf95_inquire_attribute(ncid, varid, name, nclen=att_len, &
213         ncerr=ncerr_not_opt)
214    if (ncerr_not_opt == nf90_noerr) then
215       if (att_len /= 1) then
216          print *, "nf95_get_att_one_Eightbytereal"
217          print *, "varid = ", varid
218          print *, "attribute name: ", name
219          print *, 'the attribute does not contain a single value'
220          print *, "number of values in attribute: ", att_len
221          stop 1
222       end if
223    end if
224
225    ncerr_not_opt = nf90_get_att(ncid, varid, name, values)
226    if (present(ncerr)) then
227       ncerr = ncerr_not_opt
228    else
229       if (ncerr_not_opt /= nf95_noerr) call &
230            nf95_abort("nf95_get_att_one_Eightbytereal " // trim(name), &
231            ncerr_not_opt, ncid, varid)
232    end if
233
234  end subroutine nf95_get_att_one_Eightbytereal
235
236end module nf95_get_att_m
Note: See TracBrowser for help on using the repository browser.