FV3 Bundle
sim_nc_nlm_mod.F90
Go to the documentation of this file.
1 !***********************************************************************
2 !* GNU General Public License *
3 !* This file is a part of fvGFS. *
4 !* *
5 !* fvGFS is free software; you can redistribute it and/or modify it *
6 !* and are expected to follow the terms of the GNU General Public *
7 !* License as published by the Free Software Foundation; either *
8 !* version 2 of the License, or (at your option) any later version. *
9 !* *
10 !* fvGFS is distributed in the hope that it will be useful, but *
11 !* WITHOUT ANY WARRANTY; without even the implied warranty of *
12 !* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU *
13 !* General Public License for more details. *
14 !* *
15 !* For the full text of the GNU General Public License, *
16 !* write to: Free Software Foundation, Inc., *
17 !* 675 Mass Ave, Cambridge, MA 02139, USA. *
18 !* or see: http://www.gnu.org/licenses/gpl.html *
19 !***********************************************************************
21 
22 ! This is S-J Lin's private netcdf file reader
23 ! This code is needed because FMS utilitty (read_data) led to too much
24 ! memory usage and too many files openned. Perhaps lower-level FMS IO
25 ! calls should be used instead.
26 
27 #if defined(OLD_PT_TO_T) || defined(OLD_COS_SG)
28 #error
29 #error Compile time options -DOLD_PT_TO_T and -DOLD_COS_SG are no longer supported. Please remove them from your XML.
30 #error
31 #endif
32 
33  use mpp_mod, only: mpp_error, fatal
34 
35  implicit none
36 #include <netcdf.inc>
37 
38  private
42 
43  contains
44 
45  subroutine open_ncfile( iflnm, ncid )
46  character(len=*), intent(in):: iflnm
47  integer, intent(out):: ncid
48  integer:: status
49 
50  status = nf_open(iflnm, nf_nowrite, ncid)
51  if (status .ne. nf_noerr) call handle_err(status)
52 
53 
54  end subroutine open_ncfile
55 
56 
57  subroutine close_ncfile( ncid )
58  integer, intent(in):: ncid
59  integer:: status
60 
61  status = nf_close(ncid)
62  if (status .ne. nf_noerr) call handle_err(status)
63 
64 
65  end subroutine close_ncfile
66 
67 
68  subroutine get_ncdim1( ncid, var1_name, im )
69  integer, intent(in):: ncid
70  character(len=*), intent(in):: var1_name
71  integer, intent(out):: im
72  integer:: status, var1id
73 
74  status = nf_inq_dimid(ncid, var1_name, var1id)
75  if (status .ne. nf_noerr) call handle_err(status)
76 
77  status = nf_inq_dimlen(ncid, var1id, im)
78  if (status .ne. nf_noerr) call handle_err(status)
79 
80  end subroutine get_ncdim1
81 
82 
83 
84 
85  subroutine get_var1_double( ncid, var1_name, im, var1, var_exist )
86  integer, intent(in):: ncid
87  character(len=*), intent(in):: var1_name
88  integer, intent(in):: im
89  logical, intent(out), optional:: var_exist
90  real(kind=8), intent(out):: var1(im)
91  integer:: status, var1id
92 
93  status = nf_inq_varid(ncid, var1_name, var1id)
94  if (status .ne. nf_noerr) then
95 ! call handle_err(status)
96  if(present(var_exist) ) var_exist = .false.
97  else
98  status = nf_get_var_double(ncid, var1id, var1)
99  if (status .ne. nf_noerr) call handle_err(status)
100  if(present(var_exist) ) var_exist = .true.
101  endif
102 
103 
104  end subroutine get_var1_double
105 
106 
107 ! 4-byte data:
108  subroutine get_var1_real( ncid, var1_name, im, var1, var_exist )
109  integer, intent(in):: ncid
110  character(len=*), intent(in):: var1_name
111  integer, intent(in):: im
112  logical, intent(out), optional:: var_exist
113  real(kind=4), intent(out):: var1(im)
114  integer:: status, var1id
115 
116  status = nf_inq_varid(ncid, var1_name, var1id)
117  if (status .ne. nf_noerr) then
118 ! call handle_err(status)
119  if(present(var_exist) ) var_exist = .false.
120  else
121  status = nf_get_var_real(ncid, var1id, var1)
122  if (status .ne. nf_noerr) call handle_err(status)
123  if(present(var_exist) ) var_exist = .true.
124  endif
125 
126 
127  end subroutine get_var1_real
128 
129  subroutine get_var2_real( ncid, var_name, im, jm, var2 )
130  integer, intent(in):: ncid
131  character(len=*), intent(in):: var_name
132  integer, intent(in):: im, jm
133  real(kind=4), intent(out):: var2(im)
134 
135  integer:: status, var1id
136 
137  status = nf_inq_varid(ncid, var_name, var1id)
138  if (status .ne. nf_noerr) call handle_err(status)
139 
140  status = nf_get_var_real(ncid, var1id, var2)
141  if (status .ne. nf_noerr) call handle_err(status)
142 
143  end subroutine get_var2_real
144 
145  subroutine get_var2_r4( ncid, var2_name, is,ie, js,je, var2, time_slice )
146  integer, intent(in):: ncid
147  character(len=*), intent(in):: var2_name
148  integer, intent(in):: is, ie, js, je
149  real(kind=4), intent(out):: var2(is:ie,js:je)
150  integer, intent(in), optional :: time_slice
151 !
152  real(kind=4), dimension(1) :: time
153  integer, dimension(3):: start, nreco
154  integer:: status, var2id
155 
156  status = nf_inq_varid(ncid, var2_name, var2id)
157  if (status .ne. nf_noerr) call handle_err(status)
158 
159  start(1) = is; start(2) = js; start(3) = 1
160  if ( present(time_slice) ) then
161  start(3) = time_slice
162  end if
163 
164  nreco(1) = ie - is + 1
165  nreco(2) = je - js + 1
166  nreco(3) = 1
167 
168  status = nf_get_vara_real(ncid, var2id, start, nreco, var2)
169  if (status .ne. nf_noerr) call handle_err(status)
170 
171  end subroutine get_var2_r4
172 
173  subroutine get_var2_double( ncid, var2_name, im, jm, var2 )
174  integer, intent(in):: ncid
175  character(len=*), intent(in):: var2_name
176  integer, intent(in):: im, jm
177  real(kind=8), intent(out):: var2(im,jm)
178 
179  integer:: status, var2id
180 
181  status = nf_inq_varid(ncid, var2_name, var2id)
182  if (status .ne. nf_noerr) call handle_err(status)
183 
184  status = nf_get_var_double(ncid, var2id, var2)
185  if (status .ne. nf_noerr) call handle_err(status)
186 
187 
188  end subroutine get_var2_double
189 
190 
191  subroutine get_var3_double( ncid, var3_name, im, jm, km, var3 )
192  integer, intent(in):: ncid
193  character(len=*), intent(in):: var3_name
194  integer, intent(in):: im, jm, km
195  real(kind=8), intent(out):: var3(im,jm,km)
196 
197  integer:: status, var3id
198 
199  status = nf_inq_varid(ncid, var3_name, var3id)
200 
201  if (status .ne. nf_noerr) call handle_err(status)
202 
203  status = nf_get_var_double(ncid, var3id, var3)
204  if (status .ne. nf_noerr) call handle_err(status)
205 
206  end subroutine get_var3_double
207 
208  subroutine get_var3_real( ncid, var3_name, im, jm, km, var3 )
209  integer, intent(in):: ncid
210  character(len=*), intent(in):: var3_name
211  integer, intent(in):: im, jm, km
212  real(kind=4), intent(out):: var3(im,jm,km)
213 
214  integer:: status, var3id
215 
216  status = nf_inq_varid(ncid, var3_name, var3id)
217 
218  if (status .ne. nf_noerr) call handle_err(status)
219  status = nf_get_var_real(ncid, var3id, var3)
220 
221  if (status .ne. nf_noerr) call handle_err(status)
222 
223  end subroutine get_var3_real
224 
225 
226  subroutine get_var3_r4( ncid, var3_name, is,ie, js,je, ks,ke, var3, time_slice )
227  integer, intent(in):: ncid
228  character(len=*), intent(in):: var3_name
229  integer, intent(in):: is, ie, js, je, ks,ke
230  real(kind=4), intent(out):: var3(is:ie,js:je,ks:ke)
231  integer, intent(in), optional :: time_slice
232 !
233  real(kind=4), dimension(1) :: time
234  integer, dimension(4):: start, nreco
235  integer:: status, var3id
236 
237  status = nf_inq_varid(ncid, var3_name, var3id)
238  if (status .ne. nf_noerr) call handle_err(status)
239 
240  start(1) = is; start(2) = js; start(3) = ks; start(4) = 1
241  if ( present(time_slice) ) then
242  start(4) = time_slice
243  end if
244 
245  nreco(1) = ie - is + 1
246  nreco(2) = je - js + 1
247  nreco(3) = ke - ks + 1
248  nreco(4) = 1
249 
250  status = nf_get_vara_real(ncid, var3id, start, nreco, var3)
251  if (status .ne. nf_noerr) call handle_err(status)
252 
253  end subroutine get_var3_r4
254  subroutine get_var4_real( ncid, var4_name, im, jm, km, nt, var4 )
255  implicit none
256 #include <netcdf.inc>
257  integer, intent(in):: ncid
258  character*(*), intent(in):: var4_name
259  integer, intent(in):: im, jm, km, nt
260  real*4:: wk4(im,jm,km,4)
261  real*4, intent(out):: var4(im,jm)
262  integer:: status, var4id
263  integer:: start(4), icount(4)
264  integer:: i,j
265 
266  start(1) = 1
267  start(2) = 1
268  start(3) = 1
269  start(4) = nt
270 
271  icount(1) = im ! all range
272  icount(2) = jm ! all range
273  icount(3) = km ! all range
274  icount(4) = 1 ! one time level at a time
275 
276 ! write(*,*) nt, 'Within get_var4_double: ', var4_name
277 
278  status = nf_inq_varid(ncid, var4_name, var4id)
279 ! write(*,*) '#1', status, ncid, var4id
280 
281  status = nf_get_vara_real(ncid, var4id, start, icount, var4)
282 ! status = nf_get_vara_real(ncid, var4id, start, icount, wk4)
283 ! write(*,*) '#2', status, ncid, var4id
284 
285  do j=1,jm
286  do i=1,im
287 ! var4(i,j) = wk4(i,j,1,nt)
288  enddo
289  enddo
290 
291  if (status .ne. nf_noerr) call handle_err(status)
292 
293  end subroutine get_var4_real
294 
295 
296  subroutine get_var4_double( ncid, var4_name, im, jm, km, nt, var4 )
297  integer, intent(in):: ncid
298  character(len=*), intent(in):: var4_name
299  integer, intent(in):: im, jm, km, nt
300  real(kind=8), intent(out):: var4(im,jm,km,1)
301  integer:: status, var4id
302 !
303  integer:: start(4), icount(4)
304 
305  start(1) = 1
306  start(2) = 1
307  start(3) = 1
308  start(4) = nt
309 
310  icount(1) = im ! all range
311  icount(2) = jm ! all range
312  icount(3) = km ! all range
313  icount(4) = 1 ! one time level at a time
314 
315  status = nf_inq_varid(ncid, var4_name, var4id)
316  status = nf_get_vara_double(ncid, var4id, start, icount, var4)
317 
318  if (status .ne. nf_noerr) call handle_err(status)
319 
320  end subroutine get_var4_double
321 !------------------------------------------------------------------------
322 
323  subroutine get_real3( ncid, var4_name, im, jm, nt, var4 )
324 ! This is for multi-time-level 2D var
325  integer, intent(in):: ncid
326  character(len=*), intent(in):: var4_name
327  integer, intent(in):: im, jm, nt
328  real(kind=4), intent(out):: var4(im,jm)
329  integer:: status, var4id
330  integer:: start(3), icount(3)
331  integer:: i,j
332 
333  start(1) = 1
334  start(2) = 1
335  start(3) = nt
336 
337  icount(1) = im
338  icount(2) = jm
339  icount(3) = 1
340 
341  status = nf_inq_varid(ncid, var4_name, var4id)
342  status = nf_get_vara_real(ncid, var4id, start, icount, var4)
343 
344  if (status .ne. nf_noerr) call handle_err(status)
345 
346  end subroutine get_real3
347 !------------------------------------------------------------------------
348 
349  logical function check_var( ncid, var3_name)
350  integer, intent(in):: ncid
351  character(len=*), intent(in):: var3_name
352 
353  integer:: status, var3id
354 
355  status = nf_inq_varid(ncid, var3_name, var3id)
356  check_var = (status == nf_noerr)
357 
358  end function check_var
359 
360  subroutine get_var_att_str(ncid, var_name, att_name, att)
361  implicit none
362 #include <netcdf.inc>
363  integer, intent(in):: ncid
364  character*(*), intent(in):: var_name, att_name
365  character*(*), intent(out):: att
366 
367  integer:: status, varid
368 
369  status = nf_inq_varid(ncid, var_name, varid)
370  status = nf_get_att_text(ncid, varid, att_name, att)
371 
372  if (status .ne. nf_noerr) call handle_err(status)
373 
374  end subroutine get_var_att_str
375 
376  subroutine get_var_att_double(ncid, var_name, att_name, value)
377  implicit none
378 #include <netcdf.inc>
379  integer, intent(in):: ncid
380  character*(*), intent(in):: var_name, att_name
381  real(kind=8), intent(out):: value
382 
383  integer:: status, varid
384 
385  status = nf_inq_varid(ncid, var_name, varid)
386  status = nf_get_att(ncid, varid, att_name, value)
387 
388  if (status .ne. nf_noerr) call handle_err(status)
389 
390  end subroutine get_var_att_double
391 
392 
393  subroutine handle_err(status)
394  integer status
395  character(len=120) :: errstr
396 
397  if (status .ne. nf_noerr) then
398  write(errstr,*) 'Error in handle_err: ', nf_strerror(status)
399  call mpp_error(fatal,errstr)
400  endif
401 
402  end subroutine handle_err
403 
404  subroutine calendar(year, month, day, hour)
405  integer, intent(inout) :: year ! year
406  integer, intent(inout) :: month ! month
407  integer, intent(inout) :: day ! day
408  integer, intent(inout) :: hour
409 !
410 ! Local variables
411 !
412  integer irem4,irem100
413  integer mdays(12) ! number day of month
414  data mdays /31,28,31,30,31,30,31,31,30,31,30,31/
415 !
416 !***********************************************************************
417 !****** compute current GMT ******
418 !***********************************************************************
419 !
420 !**** consider leap year
421 !
422  irem4 = mod( year, 4 )
423  irem100 = mod( year, 100 )
424  if( irem4 == 0 .and. irem100 /= 0) mdays(2) = 29
425 !
426  if( hour >= 24 ) then
427  day = day + 1
428  hour = hour - 24
429  end if
430 
431  if( day > mdays(month) ) then
432  day = day - mdays(month)
433  month = month + 1
434  end if
435  if( month > 12 ) then
436  year = year + 1
437  month = 1
438  end if
439 
440  end subroutine calendar
441 
442 end module sim_nc_nlm_mod
subroutine, public open_ncfile(iflnm, ncid)
subroutine, public get_var2_real(ncid, var_name, im, jm, var2)
subroutine, public get_var1_real(ncid, var1_name, im, var1, var_exist)
subroutine, public get_var2_r4(ncid, var2_name, is, ie, js, je, var2, time_slice)
subroutine get_var4_double(ncid, var4_name, im, jm, km, nt, var4)
subroutine get_var4_real(ncid, var4_name, im, jm, km, nt, var4)
subroutine, public get_var3_r4(ncid, var3_name, is, ie, js, je, ks, ke, var3, time_slice)
Definition: mpp.F90:39
subroutine get_var_att_str(ncid, var_name, att_name, att)
subroutine, public get_var_att_double(ncid, var_name, att_name, value)
subroutine, public get_var3_double(ncid, var3_name, im, jm, km, var3)
subroutine calendar(year, month, day, hour)
subroutine, public close_ncfile(ncid)
subroutine get_real3(ncid, var4_name, im, jm, nt, var4)
subroutine, public get_var1_double(ncid, var1_name, im, var1, var_exist)
subroutine, public get_ncdim1(ncid, var1_name, im)
subroutine, public handle_err(status)
logical function, public check_var(ncid, var3_name)
subroutine, public get_var3_real(ncid, var3_name, im, jm, km, var3)
subroutine, public get_var2_double(ncid, var2_name, im, jm, var2)