Subversion Repositories tropofold.echam

Rev

Details | Last modification | View Log | RSS feed

Rev Author Line No. Line
3 michaesp 1
module netcdf_tools
2
 
3
  USE netcdf
4
 
5
  implicit none
6
 
7
  INTRINSIC :: TRIM, ADJUSTL, NINT, SIZE
8
 
9
  interface read_file
10
    module procedure  read_file_1d
11
    module procedure  read_file_2d
12
    module procedure  read_file_3d
13
    module procedure  read_file_4d
14
  end interface 
15
 
16
  contains
17
 
18
  SUBROUTINE inquire_file(fname,             &
19
             x_name, y_name, z_name,t_name,  &
20
             nx, ny, nz, nt)
21
 
22
 
23
    IMPLICIT NONE
24
 
25
 
26
    ! I/O
27
    CHARACTER(LEN=*), INTENT(IN)  :: fname   ! filename
28
    CHARACTER(LEN=*), INTENT(IN)  :: x_name  ! name of dimension in file dimension 
29
    CHARACTER(LEN=*), INTENT(IN)  :: y_name  ! name of dimension in file dimension 
30
    CHARACTER(LEN=*), INTENT(IN)  :: z_name  ! name of dimension in file dimension 
31
    CHARACTER(LEN=*), INTENT(IN)  :: t_name ! name of dimension in file dimension 
32
    INTEGER,          INTENT(OUT) :: nx  ! file dimension lenght 
33
    INTEGER,          INTENT(OUT) :: ny  ! file dimension length
34
    INTEGER,          INTENT(OUT) :: nz  ! file dimension length
35
    INTEGER,          INTENT(OUT) :: nt  ! file dimension length
36
 
37
    ! LOCAL
38
    INTEGER                     :: status
39
    CHARACTER(LEN=*), PARAMETER :: substr = 'inquire_file'
40
    INTEGER,SAVE                :: ncid   ! netCDF-ID
41
    INTEGER                     :: dimid, varid
42
    LOGICAL                     :: file_exists  ! checking existence of file
43
 
44
 
45
    INTEGER, DIMENSION(:), ALLOCATABLE  :: date_file
46
 
47
    CHARACTER(LEN=30) :: name_dim   ! line
48
 
49
 
50
    INQUIRE(FILE=TRIM(fname), EXIST=file_exists)
51
    IF (.not.file_exists) THEN
52
      WRITE(*,*) 'File ',TRIM(fname),' does NOT exist! skipping....'
53
      STOP
54
    ENDIF
55
 
56
    ! OPEN FILE
57
    CALL NFERR( status, &
58
         nf90_open(TRIM(fname), NF90_NOWRITE, ncid) &
59
         ,1)
60
    ! latitude dimension check
61
    CALL  NFERR( status, &
62
         nf90_inq_dimid(ncid, x_name, dimid ) &
63
         ,2)
64
    CALL  NFERR( status, &
65
         nf90_Inquire_Dimension(ncid, dimid, name_dim, nx ) &
66
         ,3)
67
    ! longitude dimension check
68
    CALL  NFERR( status, &
69
         nf90_inq_dimid(ncid, y_name, dimid ) &
70
         ,4)
71
    CALL  NFERR( status, &
72
         nf90_Inquire_Dimension(ncid, dimid, name_dim, ny ) &
73
         ,5)
74
    ! vertical dimension check
75
    CALL  NFERR( status, &
76
         nf90_inq_dimid(ncid, z_name, dimid ) &
77
         ,6)
78
    CALL  NFERR( status, &
79
         nf90_Inquire_Dimension(ncid, dimid, name_dim, nz ) &
80
         ,7)
81
    ! time dimension check
82
    CALL  NFERR( status, &
83
         nf90_inq_dimid(ncid, t_name, dimid ) &
84
         ,8)
85
    CALL  NFERR( status, &
86
         nf90_Inquire_Dimension(ncid, dimid, name_dim, nt ) &
87
         ,9)
88
 
89
 
90
    !CLOSE FILE
91
    CALL NFERR( status, &
92
         nf90_close(ncid) &
93
         ,14)
94
 
95
    ! RETURN
96
    status = 0
97
 
98
  END SUBROUTINE inquire_file
99
 
100
  ! ------------------------------------------------------------------
101
 
102
  ! ------------------------------------------------------------------------
103
  SUBROUTINE read_file_1D(fname, varname,  data_file)
104
 
105
 
106
    IMPLICIT NONE
107
 
108
    ! I/O
109
    CHARACTER(LEN=*),          INTENT(IN)  :: fname       ! filename
110
    CHARACTER(LEN=*),          INTENT(IN)  :: varname         ! variable name
111
    REAL, DIMENSION(:), INTENT(OUT):: data_file ! INTENT(OUT)
112
 
113
    ! LOCAL
114
    INTEGER                     :: status
115
    CHARACTER(LEN=*), PARAMETER :: substr = 'read_file_1d'
116
    INTEGER,SAVE                :: ncid   ! netCDF-ID
117
    INTEGER                     :: dimid, varid
118
 
119
    ! OPEN FILE
120
    CALL NFERR( status, &
121
         nf90_open(TRIM(fname), NF90_NOWRITE, ncid) &
122
         ,21)
123
 
124
    CALL  NFERR( status, &
125
         nf90_inq_varid(ncid, TRIM(varname), varid ) &
126
         ,22)
127
    IF (status.ne.0) then
128
        write(*,*) "variable not found in NetCDF file, skipping"
129
        STOP
130
    ENDIF
131
 
132
    CALL  NFERR( status, &
133
         nf90_get_var(ncid, varid, data_file ) &
134
         ,23)
135
 
136
    !CLOSE FILE
137
    CALL NFERR( status, &
138
         nf90_close(ncid) &
139
         ,24)
140
 
141
    ! RETURN
142
    status = 0
143
 
144
  END SUBROUTINE read_file_1D
145
  ! ------------------------------------------------------------------------
146
  ! ------------------------------------------------------------------------
147
  SUBROUTINE read_file_2D(fname, varname,  data_file)
148
 
149
 
150
    IMPLICIT NONE
151
 
152
    ! I/O
153
    CHARACTER(LEN=*),          INTENT(IN)  :: fname       ! filename
154
    CHARACTER(LEN=*),          INTENT(IN)  :: varname         ! variable name
155
    REAL, DIMENSION(:,:), INTENT(OUT):: data_file ! INTENT(OUT)
156
 
157
    ! LOCAL
158
    INTEGER                     :: status
159
    CHARACTER(LEN=*), PARAMETER :: substr = 'read_file_2d'
160
    INTEGER,SAVE                :: ncid   ! netCDF-ID
161
    INTEGER                     :: dimid, varid
162
 
163
    ! OPEN FILE
164
    CALL NFERR( status, &
165
         nf90_open(TRIM(fname), NF90_NOWRITE, ncid) &
166
         ,21)
167
 
168
    CALL  NFERR( status, &
169
         nf90_inq_varid(ncid, TRIM(varname), varid ) &
170
         ,22)
171
    IF (status.ne.0) then
172
        write(*,*) "variable not found in NetCDF file, skipping"
173
        STOP
174
    ENDIF
175
 
176
    CALL  NFERR( status, &
177
         nf90_get_var(ncid, varid, data_file ) &
178
         ,23)
179
 
180
    !CLOSE FILE
181
    CALL NFERR( status, &
182
         nf90_close(ncid) &
183
         ,24)
184
 
185
    ! RETURN
186
    status = 0
187
 
188
  END SUBROUTINE read_file_2D
189
  ! ------------------------------------------------------------------------
190
  ! ------------------------------------------------------------------------
191
  SUBROUTINE read_file_3D(fname, varname,  data_file)
192
 
193
 
194
    IMPLICIT NONE
195
 
196
    ! I/O
197
    CHARACTER(LEN=*),          INTENT(IN)  :: fname       ! filename
198
    CHARACTER(LEN=*),          INTENT(IN)  :: varname         ! variable name
199
    REAL, DIMENSION(:,:,:), INTENT(OUT):: data_file ! INTENT(OUT)
200
 
201
    ! LOCAL
202
    INTEGER                     :: status
203
    CHARACTER(LEN=*), PARAMETER :: substr = 'read_file_3d'
204
    INTEGER,SAVE                :: ncid   ! netCDF-ID
205
    INTEGER                     :: dimid, varid
206
 
207
    ! OPEN FILE
208
    CALL NFERR( status, &
209
         nf90_open(TRIM(fname), NF90_NOWRITE, ncid) &
210
         ,21)
211
 
212
    CALL  NFERR( status, &
213
         nf90_inq_varid(ncid, TRIM(varname), varid ) &
214
         ,22)
215
    IF (status.ne.0) then
216
        write(*,*) "variable not found in NetCDF file, skipping"
217
        STOP
218
    ENDIF
219
 
220
    CALL  NFERR( status, &
221
         nf90_get_var(ncid, varid, data_file ) &
222
         ,23)
223
 
224
    !CLOSE FILE
225
    CALL NFERR( status, &
226
         nf90_close(ncid) &
227
         ,24)
228
 
229
    ! RETURN
230
    status = 0
231
 
232
  END SUBROUTINE read_file_3D
233
  ! ------------------------------------------------------------------------
234
  ! ------------------------------------------------------------------------
235
  SUBROUTINE read_file_4D(fname, varname,  data_file)
236
 
237
 
238
    IMPLICIT NONE
239
 
240
    ! I/O
241
    CHARACTER(LEN=*),          INTENT(IN)  :: fname       ! filename
242
    CHARACTER(LEN=*),          INTENT(IN)  :: varname         ! variable name
243
    REAL, DIMENSION(:,:,:,:), INTENT(OUT):: data_file ! INTENT(OUT)
244
 
245
    ! LOCAL
246
    INTEGER                     :: status
247
    CHARACTER(LEN=*), PARAMETER :: substr = 'read_file_4d'
248
    INTEGER,SAVE                :: ncid   ! netCDF-ID
249
    INTEGER                     :: dimid, varid
250
 
251
    ! OPEN FILE
252
    CALL NFERR( status, &
253
         nf90_open(TRIM(fname), NF90_NOWRITE, ncid) &
254
         ,21)
255
 
256
    CALL  NFERR( status, &
257
         nf90_inq_varid(ncid, TRIM(varname), varid ) &
258
         ,22)
259
    IF (status.ne.0) then
260
        write(*,*) "variable not found in NetCDF file, skipping"
261
        STOP
262
    ENDIF
263
 
264
    CALL  NFERR( status, &
265
         nf90_get_var(ncid, varid, data_file ) &
266
         ,23)
267
 
268
    !CLOSE FILE
269
    CALL NFERR( status, &
270
         nf90_close(ncid) &
271
         ,24)
272
 
273
    ! RETURN
274
    status = 0
275
 
276
  END SUBROUTINE read_file_4D
277
  ! ------------------------------------------------------------------------
278
 
279
  ! ------------------------------------------------------------------------
280
  SUBROUTINE read_att(fname, varname, attname, str)
281
 
282
 
283
    IMPLICIT NONE
284
 
285
    ! I/O
286
    CHARACTER(LEN=*),          INTENT(IN)  :: fname      ! filename
287
    CHARACTER(LEN=*),          INTENT(IN)  :: varname    ! variable name
288
    CHARACTER(LEN=*),          INTENT(IN)  :: attname    ! attribute name
289
    CHARACTER(LEN=*),          INTENT(OUT) :: str        ! date string
290
 
291
    ! LOCAL
292
    INTEGER                     :: status
293
    CHARACTER(LEN=*), PARAMETER :: substr = 'read_startdate'
294
    INTEGER,SAVE                :: ncid   ! netCDF-ID
295
    INTEGER                     :: dimid, varid
296
 
297
    ! OPEN FILE
298
    CALL NFERR( status, &
299
         nf90_open(TRIM(fname), NF90_NOWRITE, ncid) &
300
         ,21)
301
 
302
    CALL  NFERR( status, &
303
         nf90_inq_varid(ncid, TRIM(varname), varid ) &
304
         ,22)
305
    IF (status.ne.0) then
306
        write(*,*) "variable not found in NetCDF file, skipping"
307
        STOP
308
    ENDIF
309
 
310
    CALL  NFERR( status, &
311
         nf90_get_att(ncid, varid, TRIM(attname), str ) &
312
         ,23)
313
 
314
    !CLOSE FILE
315
    CALL NFERR( status, &
316
         nf90_close(ncid) &
317
         ,24)
318
 
319
    ! RETURN
320
    status = 0
321
 
322
  END SUBROUTINE read_att
323
  ! ------------------------------------------------------------------------
324
 
325
  ! ------------------------------------------------------------------------
326
  SUBROUTINE nc_dump(fname, x_data, y_data,z_data,t_data,  &
327
                     x_units,y_units, z_units,t_units, aps, ak, bk, &
328
                     label,fold, sfold, mfold, dfold, tp, dp, pmin, pmax)
329
 
330
    USE netcdf
331
 
332
    IMPLICIT NONE
333
 
334
    CHARACTER(LEN=*), INTENT(IN) :: fname
335
    REAL, DIMENSION(:), INTENT(IN) :: x_data, y_data,z_data,t_data
336
    CHARACTER(LEN=*), INTENT(IN) :: x_units,y_units, z_units,t_units
337
    REAL, DIMENSION(:,:,:), INTENT(IN) :: aps
338
    REAL, DIMENSION(:), INTENT(IN) :: ak,bk
339
    REAL, DIMENSION(:,:,:,:), INTENT(IN) :: label
340
    REAL, DIMENSION(:,:,:), INTENT(IN) :: fold
341
    REAL, DIMENSION(:,:,:), INTENT(IN) :: sfold
342
    REAL, DIMENSION(:,:,:), INTENT(IN) :: mfold
343
    REAL, DIMENSION(:,:,:), INTENT(IN) :: dfold
344
    REAL, DIMENSION(:,:,:), INTENT(IN) :: dp
345
    REAL, DIMENSION(:,:,:), INTENT(IN) :: tp
346
    REAL, DIMENSION(:,:,:), INTENT(IN) :: pmin
347
    REAL, DIMENSION(:,:,:), INTENT(IN) :: pmax
348
 
349
    ! LOCAL
350
    INTEGER        :: status
351
    CHARACTER(LEN=*), PARAMETER :: substr = 'nc_dump'
352
    INTEGER :: ncid      ! netCDF-ID
353
    INTEGER :: nlon,nlat,nlev,ntime
354
    INTEGER :: dimid_lat, dimid_lon, dimid_lev, dimid_ilev, dimid_time
355
    INTEGER :: varid_lat, varid_lon, varid_lev, varid_ilev, varid_time
356
    INTEGER :: varid_aps, varid_hyam, varid_hybm, varid_label
357
    INTEGER :: varid_fold, varid_sfold, varid_mfold, varid_dfold, varid_tp, varid_dp
358
    INTEGER :: varid_pmin, varid_pmax
359
 
360
    !
361
    CHARACTER(LEN=8)         :: date
362
    CHARACTER(LEN=10)        :: time
363
    CHARACTER(LEN=5)         :: zone
364
 
365
 
366
     nlon  = SIZE(x_data)
367
     nlat  = SIZE(y_data)
368
     nlev  = SIZE(z_data)
369
     ntime = SIZE(t_data)
370
 
371
    ! CREATE NEW FILE
372
    CALL NFERR(status, &
373
         nf90_create(TRIM(fname), NF90_CLOBBER, ncid) &
374
         ,51)
375
 
376
    ! ADD GLOBALE ATTRIBUTES
377
    ! - VERSION
378
    CALL NFERR(status, &
379
         nf90_put_att(ncid, NF90_GLOBAL, 'contact:',     &
380
         'Andrea Pozzer, MPIC, Mainz') &
381
         ,52)
382
    ! - DATE AND TIME
383
    CALL DATE_AND_TIME(date, time, zone)
384
    CALL NFERR(status, &
385
         nf90_put_att(ncid, NF90_GLOBAL, 'date', date) &
386
         ,53)
387
    CALL NFERR(status, &
388
         nf90_put_att(ncid, NF90_GLOBAL, 'time', TRIM(time)//TRIM(zone)) &
389
         ,54)
390
    ! DEFINE DIMENSIONS
391
    CALL NFERR(status, &
392
         nf90_def_dim(ncid, 'lon', nlon, dimid_lon) &
393
         ,56)
394
    CALL NFERR(status, &
395
         nf90_def_dim(ncid, 'lat', nlat, dimid_lat) &
396
         ,57)
397
    CALL NFERR(status, &
398
         nf90_def_dim(ncid, 'lev', nlev, dimid_lev) &
399
         ,58)
400
    CALL NFERR(status, &
401
         nf90_def_dim(ncid, 'time', NF90_UNLIMITED, dimid_time) &
402
         ,59)
403
 
404
    ! DEFINE COORDINATE VARIABLES WITH ATTRIBUTES
405
    CALL NFERR(status, &
406
         nf90_def_var(ncid, 'lon', NF90_FLOAT, (/ dimid_lon /), varid_lon) &
407
         ,60)
408
    CALL NFERR(status, &
409
         nf90_put_att(ncid, varid_lon, 'long_name', 'longitude') &
410
         ,61)
411
    CALL NFERR(status, &
412
         nf90_put_att(ncid, varid_lon, 'units', x_units) &
413
         ,62)
414
 
415
    CALL NFERR(status, &
416
         nf90_def_var(ncid, 'lat', NF90_FLOAT, (/ dimid_lat /), varid_lat) &
417
         ,63)
418
    CALL NFERR(status, &
419
         nf90_put_att(ncid, varid_lat, 'long_name', 'latitude') &
420
         ,63)
421
    CALL NFERR(status, &
422
         nf90_put_att(ncid, varid_lat, 'units', y_units) &
423
         ,64)
424
 
425
    CALL NFERR(status, &
426
         nf90_def_var(ncid, 'lev', NF90_FLOAT, (/ dimid_lev /), varid_lev) &
427
         ,65)
428
    CALL NFERR(status, &
429
         nf90_put_att(ncid, varid_lev, 'long_name', 'level index') &
430
         ,66)
431
    CALL NFERR(status, &
432
         nf90_put_att(ncid, varid_lev, 'units', z_units) &
433
         ,67)
434
    CALL NFERR(status, &
435
         nf90_def_var(ncid, 'time', NF90_FLOAT, (/ dimid_time /), varid_time) &
436
         ,68)
437
    CALL NFERR(status, &
438
         nf90_put_att(ncid, varid_time, 'long_name', 'time') &
439
         ,69)
440
    CALL NFERR(status, &
441
         nf90_put_att(ncid, varid_time, 'units', t_units) &
442
         ,70)
443
 
444
    ! DEFINE VARIABLES
445
    ! - aps
446
    CALL NFERR(status, &
447
         nf90_def_var(ncid, 'APS', NF90_FLOAT  &
448
         , (/ dimid_lon, dimid_lat, dimid_time /), varid_aps) &
449
         ,78)
450
    CALL NFERR(status, &
451
         nf90_put_att(ncid, varid_aps, 'long_name' &
452
         ,'Surface Pressure') &
453
         ,79)
454
    ! - ak
455
    CALL NFERR(status, &
456
         nf90_def_var(ncid, 'hyam', NF90_FLOAT  &
457
         , (/ dimid_lev /), varid_hyam) &
458
         ,78)
459
    CALL NFERR(status, &
460
         nf90_put_att(ncid, varid_hyam, 'long_name' &
461
         ,'"hybrid A coefficient at layer midpoints') &
462
         ,79)
463
    ! - bk
464
    CALL NFERR(status, &
465
         nf90_def_var(ncid, 'hybm', NF90_FLOAT  &
466
         , (/ dimid_lev /), varid_hybm) &
467
         ,78)
468
    CALL NFERR(status, &
469
         nf90_put_att(ncid, varid_hybm, 'long_name' &
470
         ,'hybrid B coefficient at layer midpoints') &
471
         ,79)
472
    ! - label
473
    CALL NFERR(status, &
474
         nf90_def_var(ncid, 'label', NF90_FLOAT  &
475
         , (/ dimid_lon, dimid_lat, dimid_lev, dimid_time /), varid_label) &
476
         ,78)
477
    CALL NFERR(status, &
478
         nf90_put_att(ncid, varid_label, 'long_name' &
479
         ,'TR=1,ST=2,ST.CUTOFF=3,TR.TR.CUTOFF=4,PV.BLOB=5') &
480
         ,79)
481
    ! - fold
482
    CALL NFERR(status, &
483
         nf90_def_var(ncid, 'fold', NF90_FLOAT  &
484
         , (/ dimid_lon, dimid_lat, dimid_time /), varid_fold) &
485
         ,78)
486
    CALL NFERR(status, &
487
         nf90_put_att(ncid, varid_fold, 'long_name' &
488
         , 'Tropopause folding: 1=YES, 2=NO') &
489
         ,79)
490
    ! - sfold
491
    CALL NFERR(status, &
492
         nf90_def_var(ncid, 'sfold', NF90_FLOAT  &
493
         , (/ dimid_lon, dimid_lat, dimid_time /), varid_sfold) &
494
         ,78)
495
    CALL NFERR(status, &
496
         nf90_put_att(ncid, varid_sfold, 'long_name' &
497
         , 'Shallow tropopause folding (50< >200 hPa) : 1=YES, 2=NO') &
498
         ,79)
499
    ! - mfold
500
    CALL NFERR(status, &
501
         nf90_def_var(ncid, 'mfold', NF90_FLOAT  &
502
         , (/ dimid_lon, dimid_lat, dimid_time /), varid_mfold) &
503
         ,78)
504
    CALL NFERR(status, &
505
         nf90_put_att(ncid, varid_mfold, 'long_name' &
506
         , 'Medium Tropopause folding (200< >350 hPa): 1=YES, 2=NO') &
507
         ,79)
508
    ! - dfold
509
    CALL NFERR(status, &
510
         nf90_def_var(ncid, 'dfold', NF90_FLOAT  &
511
         , (/ dimid_lon, dimid_lat, dimid_time /), varid_dfold) &
512
         ,78)
513
    CALL NFERR(status, &
514
         nf90_put_att(ncid, varid_dfold, 'long_name' &
515
         , 'Deep Tropopause folding (>=350 hPa): 1=YES, 2=NO') &
516
         ,79)
517
    ! - tropopause
518
    CALL NFERR(status, &
519
         nf90_def_var(ncid, 'tp', NF90_FLOAT  &
520
         , (/ dimid_lon, dimid_lat, dimid_time /), varid_tp) &
521
         ,78)
522
    CALL NFERR(status, &
523
         nf90_put_att(ncid, varid_tp, 'long_name' &
524
         , 'Tropopause') &
525
         ,79)
526
    CALL NFERR(status, &
527
         nf90_put_att(ncid, varid_tp, 'units' &
528
         , 'hPa') &
529
         ,79)
530
    ! - folding extension
531
    CALL NFERR(status, &
532
         nf90_def_var(ncid, 'dp', NF90_FLOAT  &
533
         , (/ dimid_lon, dimid_lat, dimid_time /), varid_dp) &
534
         ,78)
535
    CALL NFERR(status, &
536
         nf90_put_att(ncid, varid_dp, 'long_name' &
537
         , 'vertical extend of folding') &
538
         ,79)
539
    CALL NFERR(status, &
540
         nf90_put_att(ncid, varid_dp, 'units' &
541
         , 'hPa') &
542
         ,79)
543
    CALL NFERR(status, &
544
         nf90_def_var(ncid, 'pmin', NF90_FLOAT  &
545
         , (/ dimid_lon, dimid_lat, dimid_time /), varid_pmin) &
546
         ,78)
547
    CALL NFERR(status, &
548
         nf90_put_att(ncid, varid_pmin, 'long_name' &
549
         , 'vertical extend of folding') &
550
         ,79)
551
    CALL NFERR(status, &
552
         nf90_put_att(ncid, varid_pmin, 'units' &
553
         , 'hPa') &
554
         ,79)
555
    CALL NFERR(status, &
556
         nf90_def_var(ncid, 'pmax', NF90_FLOAT  &
557
         , (/ dimid_lon, dimid_lat, dimid_time /), varid_pmax) &
558
         ,78)
559
    CALL NFERR(status, &
560
         nf90_put_att(ncid, varid_pmax, 'long_name' &
561
         , 'vertical extend of folding') &
562
         ,79)
563
    CALL NFERR(status, &
564
         nf90_put_att(ncid, varid_pmax, 'units' &
565
         , 'hPa') &
566
         ,79)
567
 
568
    ! SWITCH MODUS
569
    CALL NFERR(status, &
570
         nf90_enddef(ncid) &
571
         ,82)
572
 
573
    ! SAVE COORDINATE VARIBLES
574
    CALL NFERR(status, &
575
         nf90_put_var(ncid, varid_lon, x_data)  &
576
         ,83)
577
    CALL NFERR(status, &
578
         nf90_put_var(ncid, varid_lat, y_data)  &
579
         ,84)
580
    CALL NFERR(status, &
581
         nf90_put_var(ncid, varid_lev, z_data)  &
582
         ,85)
583
    CALL NFERR(status, &
584
         nf90_put_var(ncid, varid_time, t_data) &
585
         ,86)
586
 
587
    CALL NFERR(status, &
588
         nf90_put_var(ncid, varid_aps, aps) &
589
         ,89)
590
    CALL NFERR(status, &
591
         nf90_put_var(ncid, varid_hyam, ak) &
592
         ,89)
593
    CALL NFERR(status, &
594
         nf90_put_var(ncid, varid_hybm, bk) &
595
         ,89)
596
    CALL NFERR(status, &
597
         nf90_put_var(ncid, varid_label, label) &
598
         ,89)
599
    CALL NFERR(status, &
600
         nf90_put_var(ncid, varid_fold, fold) &
601
         ,89)
602
    CALL NFERR(status, &
603
         nf90_put_var(ncid, varid_sfold, sfold) &
604
         ,89)
605
    CALL NFERR(status, &
606
         nf90_put_var(ncid, varid_mfold, mfold) &
607
         ,89)
608
    CALL NFERR(status, &
609
         nf90_put_var(ncid, varid_dfold, dfold) &
610
         ,89)
611
    CALL NFERR(status, &
612
         nf90_put_var(ncid, varid_tp, tp) &
613
         ,89)
614
    CALL NFERR(status, &
615
         nf90_put_var(ncid, varid_dp, dp) &
616
         ,89)
617
    CALL NFERR(status, &
618
         nf90_put_var(ncid, varid_pmin, pmin) &
619
         ,89)
620
    CALL NFERR(status, &
621
         nf90_put_var(ncid, varid_pmax, pmax) &
622
         ,89)
623
 
624
    ! CLOSE FILE
625
    CALL NFERR(status, &
626
         nf90_close(ncid) &
627
         ,90)
628
 
629
  END SUBROUTINE nc_dump
630
 
631
 
632
  ! ------------------------------------------------------------------
633
  SUBROUTINE NFERR(status,command,pos)
634
 
635
    IMPLICIT NONE
636
 
637
    ! I/O
638
    INTEGER,          INTENT(OUT) :: status
639
    INTEGER,          INTENT(IN) :: command
640
    INTEGER,          INTENT(IN) :: pos
641
 
642
    status=command
643
    IF (status /= NF90_NOERR) THEN
644
       WRITE(*,*) 'netCDF ERROR at position: ', pos
645
       WRITE(*,*) 'netCDF ERROR status     : ',status
646
       WRITE(*,*) 'netCDF ERROR            : ',nf90_strerror(status)
647
    END IF
648
 
649
  END SUBROUTINE NFERR
650
  ! ------------------------------------------------------------------
651
 
652
 
653
 
654
end module netcdf_tools