source: palm/trunk/UTIL/inifor/src/io.f90 @ 2742

Last change on this file since 2742 was 2718, checked in by maronga, 7 years ago

deleting of deprecated files; headers updated where needed

  • Property svn:keywords set to Id
File size: 31.4 KB
RevLine 
[2696]1!> @file src/io.f90
2!------------------------------------------------------------------------------!
[2718]3! This file is part of the PALM model system.
[2696]4!
[2718]5! PALM is free software: you can redistribute it and/or modify it under the
6! terms of the GNU General Public License as published by the Free Software
7! Foundation, either version 3 of the License, or (at your option) any later
[2696]8! version.
9!
[2718]10! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
11! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
12! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
[2696]13!
14! You should have received a copy of the GNU General Public License along with
15! PALM. If not, see <http://www.gnu.org/licenses/>.
16!
[2718]17! Copyright 2017-2018 Leibniz Universitaet Hannover
18! Copyright 2017-2018 Deutscher Wetterdienst Offenbach
[2696]19!------------------------------------------------------------------------------!
20!
21! Current revisions:
22! -----------------
23!
24!
25! Former revisions:
26! -----------------
27! $Id: io.f90 2718 2018-01-02 08:49:38Z suehring $
28! Initial revision
29!
30!
31!
32! Authors:
33! --------
34! @author Eckhard Kadasch
35!
36! Description:
37! ------------
38!> The io module contains the functions needed to read and write netCDF data in
39!> INIFOR.
40!------------------------------------------------------------------------------!
41 MODULE io
42
43    USE control
44    USE defs,                                                                  &
45        ONLY:  DATE, SNAME, PATH, PI, dp, TO_RADIANS, TO_DEGREES, VERSION
46    USE netcdf
47    USE types
48    USE util,                                                                  &
49        ONLY:  reverse, str
50
51    IMPLICIT NONE
52
53 CONTAINS
54
55    SUBROUTINE netcdf_define_variable(var, ncid)
56
57        TYPE(nc_var), INTENT(INOUT) ::  var
58        INTEGER, INTENT(IN)         ::  ncid
59
60        CALL check(nf90_def_var(ncid, var % name, NF90_FLOAT,       var % dimids(1:var % ndim), var % varid))
61        CALL check(nf90_put_att(ncid, var % varid, "standard_name", var % standard_name))
62        CALL check(nf90_put_att(ncid, var % varid, "long_name",     var % long_name))
63        CALL check(nf90_put_att(ncid, var % varid, "units",         var % units))
64        CALL check(nf90_put_att(ncid, var % varid, "lod",           var % lod))
65        CALL check(nf90_put_att(ncid, var % varid, "source",        var % source))
66        CALL check(nf90_put_att(ncid, var % varid, "_FillValue",    NF90_FILL_REAL))
67
68    END SUBROUTINE netcdf_define_variable
69   
70
71    SUBROUTINE netcdf_get_dimensions(var, ncid)
72
73        TYPE(nc_var), INTENT(INOUT) ::  var
74        INTEGER, INTENT(IN)         ::  ncid
75        INTEGER                     ::  i
76        CHARACTER(SNAME)            ::  null
77
78        ! Remember dimension lenghts for NetCDF writing routine
79        DO i = 1, var % ndim
80           CALL check(nf90_inquire_dimension(ncid, var % dimids(i), &
81                                             name = null, &
82                                             len  = var % dimlen(i)  ) )
83        END DO
84
85    END SUBROUTINE netcdf_get_dimensions
86
87
88!------------------------------------------------------------------------------!
89! Description:
90! ------------
91!> This routine initializes Inifor. This includes parsing command-line
92!> arguments, setting the names of the input and output file names as well as
93!> the name of the input namelist and, subsequently, reading in and setting grid
94!> parameters for the PALM-4U computational grid.
95!------------------------------------------------------------------------------!
96    SUBROUTINE parse_command_line_arguments( start_date, hhl_file,             &
97       soiltyp_file, static_driver_file, input_path, output_file,              &
98       namelist_file, ug, vg, p0, z0, mode )
99
100       CHARACTER(LEN=PATH), INTENT(INOUT)  ::  hhl_file, soiltyp_file,         &
101           static_driver_file, input_path, output_file, namelist_file
102       CHARACTER(LEN=SNAME), INTENT(INOUT) ::  mode
103       REAL(dp), INTENT(INOUT)             ::  ug, vg, p0, z0
104       CHARACTER(LEN=DATE), INTENT(INOUT)  ::  start_date
105
106       CHARACTER(LEN=PATH)     ::  option, arg
107       INTEGER                 ::  arg_count, i
108
109       arg_count = COMMAND_ARGUMENT_COUNT()
110       IF (arg_count .GT. 0)  THEN
111
112          ! Every option should have an argument.
113          IF ( MOD(arg_count, 2) .NE. 0 )  THEN
114             message = "Syntax error in command line."
115             CALL abort('parse_command_line_arguments', message)
116          END IF
117         
118          message = "The -clon and -clat command line options are depricated. " // &
119             "Please remove them form your inifor command and specify the " // &
120             "location of the PALM-4U origin either" // NEW_LINE(' ') // &
121             "   - by setting the namelist parameters 'origin_lon' and 'origin_lat, or'" // NEW_LINE(' ') // &
122             "   - by providing a static driver netCDF file via the -static command-line option."
123
124          ! Loop through option/argument pairs.
125          DO i = 1, arg_count, 2
126
127             CALL GET_COMMAND_ARGUMENT( i, option )
128             CALL GET_COMMAND_ARGUMENT( i+1, arg )
129
130             SELECT CASE( TRIM(option) )
131
132             CASE( '-date' )
133                start_date = TRIM(arg)
134
135             ! Elevation of the PALM-4U domain above sea level
136             CASE( '-z0' )
137                READ(arg, *) z0
138
139             ! surface pressure, at z0
140             CASE( '-p0' )
141                READ(arg, *) p0
142
143             ! surface pressure, at z0
144             CASE( '-ug' )
145                READ(arg, *) ug
146
147             ! surface pressure, at z0
148             CASE( '-vg' )
149                READ(arg, *) vg
150
151             ! Domain centre geographical longitude
152             CASE( '-clon' )
153                CALL abort('parse_command_line_arguments', message)         
154                !READ(arg, *) lambda_cg
155                !lambda_cg = lambda_cg * TO_RADIANS
156
157             ! Domain centre geographical latitude
158             CASE( '-clat' )
159                CALL abort('parse_command_line_arguments', message)         
160                !READ(arg, *) phi_cg
161                !phi_cg = phi_cg * TO_RADIANS
162
163             CASE( '-path' )
164                 input_path = TRIM(arg)
165
166             CASE( '-hhl' )
167                 hhl_file = TRIM(arg)
168
169             CASE( '-static' )
170                 static_driver_file = TRIM(arg)
171
172             CASE( '-soil' )
173                 soiltyp_file = TRIM(arg)
174
175             CASE( '-o' )
176                output_file = TRIM(arg)
177
178             CASE( '-n' )
179                namelist_file = TRIM(arg)
180
181             ! Initialization mode: 'profile' / 'volume'
182             CASE( '-mode' )
183                mode = TRIM(arg)
184
185                SELECT CASE( TRIM(mode) )
186
187                CASE( 'profile' )
188
189                CASE DEFAULT
190                   message = "Mode '" // TRIM(mode) // "' is not supported. " //&
191                             "Currently, '-mode profile' is the only supported option. " //&
192                             "Select this one or omit the -mode option entirely."
193                   CALL abort( 'parse_command_line_arguments', message ) 
194                END SELECT
195
196             CASE DEFAULT
197                message = "unknown option '" // TRIM(option(2:)) // "'."
198                CALL abort('parse_command_line_arguments', message)
199
200             END SELECT
201
202          END DO
203
204       ELSE
205           
206          message = "No arguments present, using default input and output files"
207          CALL report('parse_command_line_arguments', message)
208
209       END IF
210
211   END SUBROUTINE parse_command_line_arguments
212
213
214!------------------------------------------------------------------------------!
215! Description:
216! ------------
217!> This routine initializes the Inifor output file, i.e. the PALM-4U
218!> initializing and forcing data as a NetCDF file.
219!>
220!> Besides writing metadata, such as global attributes, coordinates, variables,
221!> in the NetCDF file, various NetCDF IDs are saved for later, when Inifor
222!> writes the actual data.
223!------------------------------------------------------------------------------!
224    SUBROUTINE setup_netcdf_dimensions(output_file, palm_grid)
225
226       TYPE(nc_file), INTENT(INOUT)      ::  output_file
227       TYPE(grid_definition), INTENT(IN) ::  palm_grid
228
229       CHARACTER (LEN=SNAME) ::  date
230       INTEGER               ::  ncid, nx, ny, nz, nt, dimids(3), dimvarids(3)
231       REAL(dp)              ::  z0
232
233       ! Create the NetCDF file. NF90_CLOBBER selects overwrite mode.
234       CALL check(nf90_create(TRIM(output_file % name), OR(NF90_CLOBBER, NF90_HDF5), ncid))
235
236!
237!------------------------------------------------------------------------------
238!- Section 1: Write global NetCDF attributes
239!------------------------------------------------------------------------------
240       CALL date_and_time(date)
241       CALL check(nf90_put_att(ncid, NF90_GLOBAL, 'title',          'PALM input file for scenario ...'))
242       CALL check(nf90_put_att(ncid, NF90_GLOBAL, 'institution',    'Deutscher Wetterdienst, Offenbach'))
243       CALL check(nf90_put_att(ncid, NF90_GLOBAL, 'author',         'Eckhard Kadasch, eckhard.kadasch@dwd.de'))
244       CALL check(nf90_put_att(ncid, NF90_GLOBAL, 'history',        'Created on '//date))
245       CALL check(nf90_put_att(ncid, NF90_GLOBAL, 'references',     '--'))
246       CALL check(nf90_put_att(ncid, NF90_GLOBAL, 'comment',        '--'))
247       CALL check(nf90_put_att(ncid, NF90_GLOBAL, 'origin_lat',     '--'))
248       CALL check(nf90_put_att(ncid, NF90_GLOBAL, 'origin_lon',     '--'))
249       CALL check(nf90_put_att(ncid, NF90_GLOBAL, 'inifor_version', TRIM(VERSION)))
250       CALL check(nf90_put_att(ncid, NF90_GLOBAL, 'palm_version',   '--'))
251
252!
253!------------------------------------------------------------------------------
254!- Section 2: Define NetCDF dimensions and coordinates
255!------------------------------------------------------------------------------
256       nt = SIZE(output_file % time)
257       nx = palm_grid % nx
258       ny = palm_grid % ny
259       nz = palm_grid % nz
260       z0 = palm_grid % z0
261
262!
263!------------------------------------------------------------------------------
264!- Section 2a: Define dimensions for cell centers (scalars in soil and atmosph.)
265!------------------------------------------------------------------------------
266       dimids = (/0, 0, 0/) ! reset dimids
267          CALL check( nf90_def_dim(ncid, "x", nx+1, dimids(1)) )
268          CALL check( nf90_def_dim(ncid, "y", ny+1, dimids(2)) )
269          CALL check( nf90_def_dim(ncid, "z", nz+1, dimids(3)) )
270          output_file % dimids_scl = dimids ! save dimids for later
271
272       dimvarids = (/0, 0, 0/) ! reset dimvarids
273          CALL check(nf90_def_var(ncid, "x", NF90_FLOAT, dimids(1), dimvarids(1)))
274          CALL check(nf90_put_att(ncid, dimvarids(1), "standard_name", "x coordinate of cell centers"))
275          CALL check(nf90_put_att(ncid, dimvarids(1), "units", "m"))
276
277          CALL check(nf90_def_var(ncid, "y", NF90_FLOAT, dimids(2), dimvarids(2)))
278          CALL check(nf90_put_att(ncid, dimvarids(2), "standard_name", "y coordinate of cell centers"))
279          CALL check(nf90_put_att(ncid, dimvarids(2), "units", "m"))
280
281          CALL check(nf90_def_var(ncid, "z", NF90_FLOAT, dimids(3), dimvarids(3)))
282          CALL check(nf90_put_att(ncid, dimvarids(3), "standard_name", "z coordinate of cell centers"))
283          CALL check(nf90_put_att(ncid, dimvarids(3), "units", "m"))
284       output_file % dimvarids_scl = dimvarids ! save dimvarids for later
285
286       ! overwrite third dimid with the one of depth
287       CALL check(nf90_def_dim(ncid, "depth", SIZE(palm_grid % depths), dimids(3)) )
288       output_file % dimids_soil = dimids ! save dimids for later
289
290       ! overwrite third dimvarid with the one of depth
291       CALL check(nf90_def_var(ncid, "depth", NF90_FLOAT, output_file % dimids_soil(3), dimvarids(3)))
292       CALL check(nf90_put_att(ncid, dimvarids(3), "standard_name", "depth_below_land"))
293       CALL check(nf90_put_att(ncid, dimvarids(3), "positive", "down"))
294       CALL check(nf90_put_att(ncid, dimvarids(3), "units", "m"))
295       output_file % dimvarids_soil = dimvarids ! save dimvarids for later
296!
297!------------------------------------------------------------------------------
298!- Section 2b: Define dimensions for cell faces/velocities
299!------------------------------------------------------------------------------
300       dimids = (/0, 0, 0/) ! reset dimids
301          CALL check(nf90_def_dim(ncid, "xu", nx, dimids(1)) )
302          CALL check(nf90_def_dim(ncid, "yv", ny, dimids(2)) )
303          CALL check(nf90_def_dim(ncid, "zw", nz, dimids(3)) )
304       output_file % dimids_vel = dimids ! save dimids for later
305
306       dimvarids = (/0, 0, 0/) ! reset dimvarids
307          CALL check(nf90_def_var(ncid, "xu", NF90_FLOAT, dimids(1), dimvarids(1)))
308          CALL check(nf90_put_att(ncid, dimvarids(1), "standard_name", "x coordinate of cell faces"))
309          CALL check(nf90_put_att(ncid, dimvarids(1), "units", "m"))
310
311          CALL check(nf90_def_var(ncid, "yv", NF90_FLOAT, dimids(2), dimvarids(2)))
312          CALL check(nf90_put_att(ncid, dimvarids(2), "standard_name", "y coordinate of cell faces"))
313          CALL check(nf90_put_att(ncid, dimvarids(2), "units", "m"))
314
315          CALL check(nf90_def_var(ncid, "zw", NF90_FLOAT, dimids(3), dimvarids(3)))
316          CALL check(nf90_put_att(ncid, dimvarids(3), "standard_name", "z coordinate of cell faces"))
317          CALL check(nf90_put_att(ncid, dimvarids(3), "units", "m"))
318       output_file % dimvarids_vel = dimvarids ! save dimvarids for later
319
320!
321!------------------------------------------------------------------------------
322!- Section 2c: Define time dimension
323!------------------------------------------------------------------------------
324       CALL check(nf90_def_dim(ncid, "time", nt, output_file % dimid_time) )
325       CALL check(nf90_def_var(ncid, "time", NF90_FLOAT, &
326                                             output_file % dimid_time, &
327                                             output_file % dimvarid_time))
328       CALL check(nf90_put_att(ncid, output_file % dimvarid_time, "standard_name", "time"))
329       CALL check(nf90_put_att(ncid, output_file % dimvarid_time, "long_name", "time"))
330       CALL check(nf90_put_att(ncid, output_file % dimvarid_time, "units", "seconds since..."))
331
332       CALL check(nf90_enddef(ncid))
333
334!
335!------------------------------------------------------------------------------
336!- Section 3: Write grid coordinates
337!------------------------------------------------------------------------------
338       CALL check(nf90_put_var(ncid, output_file % dimvarids_scl(1), palm_grid%x))
339       CALL check(nf90_put_var(ncid, output_file % dimvarids_scl(2), palm_grid%y))
340       CALL check(nf90_put_var(ncid, output_file % dimvarids_scl(3), palm_grid%z))
341
342       CALL check(nf90_put_var(ncid, output_file % dimvarids_vel(1), palm_grid%xu))
343       CALL check(nf90_put_var(ncid, output_file % dimvarids_vel(2), palm_grid%yv))
344       CALL check(nf90_put_var(ncid, output_file % dimvarids_vel(3), palm_grid%zw))
345
346       ! TODO Read in soil depths from input file before this.
347       CALL check(nf90_put_var(ncid, output_file % dimvarids_soil(3), palm_grid%depths))
348
349       ! Write time vector
350       CALL check(nf90_put_var(ncid, output_file % dimvarid_time, output_file % time))
351
352       ! Close the file
353       CALL check(nf90_close(ncid))
354
355    END SUBROUTINE setup_netcdf_dimensions
356
357
358    SUBROUTINE setup_netcdf_variables(filename, output_variable_table)
359
360       CHARACTER (LEN=*), INTENT(IN)        ::  filename
361       TYPE(nc_var), INTENT(INOUT), TARGET  ::  output_variable_table(:)
362       TYPE(nc_var), POINTER                ::  var
363       INTEGER                              ::  i, ncid
364
365       message = "Initializing PALM-4U forcing file '" // TRIM(filename) // "'."
366       CALL report('setup_netcdf_variables', message)
367
368       CALL check(nf90_open(TRIM(filename), NF90_WRITE, ncid))
369       CALL check(nf90_redef(ncid))
370
371       DO i = 1, SIZE(output_variable_table)
372
373          var => output_variable_table(i)
374
375          IF ( var % to_be_processed )  THEN
376             message = "Defining variable #" // TRIM(str(i)) // " '" // TRIM(var%name) // "'."
377             CALL report('setup_netcdf_variables', message)
378
379             CALL netcdf_define_variable(var, ncid)
380             CALL netcdf_get_dimensions(var, ncid)
381          END IF
382           
383       END DO
384
385       CALL check(nf90_enddef(ncid))
386       CALL check(nf90_close(ncid))
387
388       message = "Forcing file '" // TRIM(filename) // "' initialized successfully."
389       CALL report('setup_netcdf_variables', message)
390
391    END SUBROUTINE setup_netcdf_variables
392
393
394!------------------------------------------------------------------------------!
395! Description:
396! ------------
397!> This routine reads and returns all input variables of the given IO group
398!> It accomodates the data by allocating a container variable that represents a
399!> list of arrays of the same length as the groups variable list. (This list
400!> will typically contain one or two items.) After the container, its members
401!> are allocated one by one with the appropriate, possibly different,
402!> dimensions.
403!>
404!> The 'group' is an INTENT(INOUT) variable so that 'get_netcdf_variable()' can
405!> record netCDF IDs in the 'in_var_list()' member variable.
406!------------------------------------------------------------------------------!
407    SUBROUTINE read_input_variables(group, iter, buffer)
408       TYPE(io_group), INTENT(INOUT), TARGET       ::  group
409       INTEGER, INTENT(IN)                         ::  iter
410       TYPE(container), ALLOCATABLE, INTENT(INOUT) ::  buffer(:)
411       INTEGER                                     ::  hour, buf_id
412       TYPE(nc_var), POINTER                       ::  input_var
413       CHARACTER(LEN=PATH), POINTER                ::  input_file
414       INTEGER                                     ::  ivar, nbuffers
415
416       message = "Reading data for I/O group '" // TRIM(group % in_var_list(1) % name) // "'."
417       CALL report('read_input_variables', message)
418
419       input_file => group % in_files(iter)
420
421!
422!------------------------------------------------------------------------------
423!- Section 1: Load input buffers for accumulated variables
424!------------------------------------------------------------------------------
425       IF (group % kind == 'running average' .OR.                              &
426           group % kind == 'accumulated')  THEN ! radiation budgets, precipitation
427
428          IF (SIZE(group % in_var_list) .GT. 1 ) THEN
429             message = "I/O groups may not contain more than one " // & 
430                       "accumulated variable. Group '" // TRIM(group % kind) //&
431                       "' contains " //                                        &
432                       TRIM( str(SIZE(group % in_var_list)) ) // "."
433             CALL abort('read_input_variables | accumulation', message)
434          END IF
435
436          ! use two buffer arrays
437          nbuffers = 2
438          IF ( .NOT. ALLOCATED( buffer ) )  ALLOCATE( buffer(nbuffers) )
439
440          ! chose correct buffer array
441          hour = iter - 1! hour of the day
442          buf_id = select_buffer(hour)
443
444 CALL run_control('time', 'read')
445          IF ( ALLOCATED(buffer(buf_id) % array) )  DEALLOCATE(buffer(buf_id) % array)
446 CALL run_control('time', 'alloc')
447
448          input_var => group % in_var_list(1)
449          buffer(buf_id) % array = get_netcdf_variable( input_file, input_var ) 
450          CALL report('read_input_variables', "Read accumulated " // TRIM(group % in_var_list(1) % name)) 
451
452          IF ( input_var % is_upside_down )  CALL reverse(buffer(buf_id) % array)
453 CALL run_control('time', 'comp')
454         
455!------------------------------------------------------------------------------
456!- Section 2: Load input buffers for normal I/O groups
457!------------------------------------------------------------------------------
458       ELSE
459
460          nbuffers = SIZE( group % in_var_list )
461          ALLOCATE( buffer(nbuffers) )
462 CALL run_control('time', 'alloc')
463         
464          DO ivar = 1, nbuffers
465
466             input_var => group % in_var_list(ivar)
467
468             ! Check wheather P or PP is present in input file
469             IF (input_var % name == 'P')  THEN
470                input_var % name = TRIM( get_pressure_var(input_file) )
471 CALL run_control('time', 'read')
472             END IF
473
474             buffer(ivar) % array = get_netcdf_variable( input_file, input_var ) 
475
476             IF ( input_var % is_upside_down )  CALL reverse(buffer(ivar) % array)
477 CALL run_control('time', 'comp')
478
479          END DO
480       END IF
481
482    END SUBROUTINE read_input_variables
483
484
485    INTEGER FUNCTION select_buffer(hour)
486       INTEGER, INTENT(IN) ::  hour
487       INTEGER             ::  step
488
489       select_buffer = 0
490       step = MODULO(hour, 3) + 1
491
492       SELECT CASE(step)
493       CASE(1, 3)
494           select_buffer = 1
495       CASE(2)
496           select_buffer = 2
497       CASE DEFAULT
498           message = "Invalid step '" // TRIM(str(step))
499           CALL abort('select_buffer', message)
500       END SELECT
501    END FUNCTION select_buffer
502
503
504!------------------------------------------------------------------------------!
505! Description:
506! ------------
507!> Checks if the input_file contains the absolute pressure, 'P', or the pressure
508!> perturbation, 'PP', and returns the appropriate string.
509!------------------------------------------------------------------------------!
510    CHARACTER(LEN=2) FUNCTION get_pressure_var(input_file) RESULT(var)
511       CHARACTER(LEN=*) ::  input_file
512       INTEGER          ::  ncid, varid
513
514       CALL check(nf90_open( TRIM(input_file), NF90_NOWRITE, ncid ))
515       IF ( nf90_inq_varid( ncid, 'P', varid ) .EQ. NF90_NOERR )  THEN
516
517          var = 'P'
518
519       ELSE IF ( nf90_inq_varid( ncid, 'PP', varid ) .EQ. NF90_NOERR )  THEN
520
521          var = 'PP'
522          CALL report('get_pressure_var', 'Using PP instead of P')
523
524       ELSE
525
526          message = "Failed to read '" // TRIM(var) // &
527                    "' from file '" // TRIM(input_file) // "'."
528          CALL abort('get_pressure_var', message)
529
530       END IF
531
532       CALL check(nf90_close(ncid))
533
534    END FUNCTION get_pressure_var
535
536
537    FUNCTION get_netcdf_attribute(filename, attribute) RESULT(attribute_value)
538
539       CHARACTER(LEN=*), INTENT(IN) ::  filename, attribute
540       REAL(dp)                     ::  attribute_value
541
542       INTEGER                      :: ncid
543
544       IF ( nf90_open( TRIM(filename), NF90_NOWRITE, ncid ) == NF90_NOERR )  THEN
545
546          CALL check(nf90_get_att(ncid, NF90_GLOBAL, TRIM(attribute), attribute_value))
547
548       ELSE
549
550          message = "Failed to read '" // TRIM(attribute) // &
551                    "' from file '" // TRIM(filename) // "'."
552          CALL abort('get_netcdf_attribute', message)
553
554       END IF
555
556    END FUNCTION get_netcdf_attribute
557
558
559
560    FUNCTION get_netcdf_variable(in_file, in_var) RESULT(buffer)
561
562       CHARACTER(LEN=PATH), INTENT(IN)      ::  in_file
563       TYPE(nc_var), INTENT(INOUT)          ::  in_var
564       REAL(dp), ALLOCATABLE                ::  buffer(:,:,:)
565       INTEGER                              ::  i, ncid, start(3)
566
567
568       ! Read in_var NetCDF attributes
569       IF ( nf90_open( TRIM(in_file), NF90_NOWRITE, ncid ) .EQ. NF90_NOERR .AND. &
570            nf90_inq_varid( ncid, in_var % name, in_var % varid ) .EQ. NF90_NOERR )  THEN
571
572          CALL check(nf90_get_att(ncid, in_var % varid, "long_name", in_var % long_name))
573          CALL check(nf90_get_att(ncid, in_var % varid, "units", in_var % units))
574
575          ! Read in_var NetCDF dimensions
576          CALL check(nf90_inquire_variable( ncid, in_var % varid,              &
577                                            ndims  = in_var % ndim,            &
578                                            dimids = in_var % dimids ))
579                                       
580          DO i = 1, in_var % ndim
581             CALL check(nf90_inquire_dimension( ncid, in_var % dimids(i),      &
582                                                name = in_var % dimname(i),    &
583                                                len  = in_var % dimlen(i) ))
584          END DO
585
586          start = (/ 1, 1, 1 /)
587          IF ( TRIM(in_var % name) .EQ. 'T_SO' )  THEN
588             ! Skip depth = 0.0 for T_SO and reduce number of depths from 9 to 8
589             in_var % dimlen(3) = in_var % dimlen(3) - 1
590
591             ! Start reading from second level, e.g. depth = 0.005 instead of 0.0
592             start(3) = 2
593          END IF
594
595          SELECT CASE(in_var % ndim)
596
597          CASE (2)
598
599             ALLOCATE( buffer( in_var % dimlen(1),                             &
600                               in_var % dimlen(2),                             &
601                               1 ) )
602
603          CASE (3)
604
605             ALLOCATE( buffer( in_var % dimlen(1),                             &
606                               in_var % dimlen(2),                             &
607                               in_var % dimlen(3) ) )
608          CASE (4)
609
610             ALLOCATE( buffer( in_var % dimlen(1),                             &
611                               in_var % dimlen(2),                             &
612                               in_var % dimlen(3) ) )
613          CASE DEFAULT
614
615             message = "Failed reading NetCDF variable " //                    &
616                TRIM(in_var % name) // " with " // TRIM(str(in_var%ndim)) //   &
617                " dimensions because only two- and and three-dimensional" //   &
618                " variables are supported."
619             CALL abort('get_netcdf_variable', message) 
620
621          END SELECT
622 CALL run_control('time', 'alloc')
623         
624          ! TODO: Check for matching dimensions of buffer and var
625          CALL check(nf90_get_var( ncid, in_var % varid, buffer,               &
626                                   start = start,                              &
627                                   count = in_var % dimlen(1:3) ) )
628
629 CALL run_control('time', 'read')
630       ELSE
631
632          message = "Failed to read '" // TRIM(in_var % name) // &
633             "' from file '" // TRIM(in_file) // "'."
634          CALL report('get_netcdf_variable', message)
635
636       END IF
637
638       CALL check(nf90_close(ncid))
639
640 CALL run_control('time', 'read')
641
642    END FUNCTION get_netcdf_variable
643
644
645    SUBROUTINE update_output(var, array, iter, output_file)
646       TYPE(nc_var), INTENT(IN)  ::  var
647       REAL(dp), INTENT(IN)      ::  array(:,:,:)
648       INTEGER, INTENT(IN)       ::  iter
649       TYPE(nc_file), INTENT(IN) ::  output_file
650
651       INTEGER ::  ncid, ndim, start(4), count(4)
652       LOGICAL ::  var_is_time_dependent
653
654       var_is_time_dependent = (                                               &
655          var % dimids( var % ndim ) == output_file % dimid_time               &
656       )
657
658       ! Skip time dimension for output
659       IF ( var_is_time_dependent )  THEN
660          ndim = var % ndim - 1
661       ELSE
662          ndim = var % ndim
663       END IF
664
665       start(:)      = (/1,1,1,1/)
666       start(ndim+1) = iter
667       count(1:ndim) = var%dimlen(1:ndim)
668
669       CALL check(nf90_open(output_file % name, NF90_WRITE, ncid))
670
671       ! Reduce dimension of output array according to variable kind
672       SELECT CASE (TRIM(var % kind))
673       
674       CASE ( 'init scalar profile', 'init u profile', 'init v profile',       &
675              'init w profile' )
676
677          CALL check(nf90_put_var( ncid, var%varid, array(1,1,:) ) )
678
679       CASE ( 'init soil', 'init scalar', 'init u', 'init v', 'init w' )
680
681          CALL check(nf90_put_var( ncid, var%varid, array(:,:,:) ) )
682
683       CASE( 'left scalar', 'right scalar', 'left w', 'right w' ) 
684
685          CALL check(nf90_put_var( ncid, var%varid, array(1,:,:),              &
686                                   start=start(1:ndim+1),                      &
687                                   count=count(1:ndim) ) )
688         
689
690          IF (.NOT. SIZE(array, 2) .EQ. var % dimlen(1))  THEN
691             PRINT *, "inifor: update_output: Dimension ", 1, " of variable ", &
692                 TRIM(var % name), " (", var % dimlen(1),                      &
693                 ") does not match the dimension of the output array (",       &
694                 SIZE(array, 2), ")."
695             STOP
696          END IF
697         
698
699       CASE( 'north scalar', 'south scalar', 'north w', 'south w' )
700
701          CALL check(nf90_put_var( ncid, var%varid, array(:,1,:),              &
702                                   start=start(1:ndim+1),                      &
703                                   count=count(1:ndim) ) )
704         
705
706       CASE( 'surface forcing', 'top scalar', 'top w' )
707
708          CALL check(nf90_put_var( ncid, var%varid, array(:,:,1),              &
709                                   start=start(1:ndim+1),                      &
710                                   count=count(1:ndim) ) )
711         
712       CASE ( 'left u', 'right u', 'left v', 'right v' )
713
714          CALL check(nf90_put_var( ncid, var%varid, array(1,:,:),              &
715                                   start=start(1:ndim+1),                      &
716                                   count=count(1:ndim) ) )
717
718       CASE ( 'north u', 'south u', 'north v', 'south v' )
719
720          CALL check(nf90_put_var( ncid, var%varid, array(:,1,:),              &
721                                   start=start(1:ndim+1),                      &
722                                   count=count(1:ndim) ) )
723
724       CASE ( 'top u', 'top v' )
725
726          CALL check(nf90_put_var( ncid, var%varid, array(:,:,1),              &
727                                   start=start(1:ndim+1),                      &
728                                   count=count(1:ndim) ) )
729
730       CASE ( 'time series' )
731
732          CALL check(nf90_put_var( ncid, var%varid, array(1,1,1),              &
733                                   start=start(1:ndim+1) ) )
734
735       CASE ( 'profile' )
736
737          CALL check(nf90_put_var( ncid, var%varid, array(1,1,:),              &
738                                   start=start(1:ndim+1),                      &
739                                   count=count(1:ndim) ) )
740
741       CASE DEFAULT
742
743           message = "Variable kind '" // TRIM(var % kind) //                  &
744                    "' not recognized."
745           CALL abort('update_output', message)
746
747       END SELECT
748
749       CALL check(nf90_close(ncid))
750
751    END SUBROUTINE update_output
752
753
754    SUBROUTINE write_netcdf_variable_2d(var, iter, output_file, buffer)
755       TYPE(nc_var), INTENT(IN)          ::  var
756       INTEGER, INTENT(IN)               ::  iter
757       TYPE(nc_file), INTENT(IN)         ::  output_file
758       REAL(dp), INTENT(IN)              ::  buffer(:,:,:)
759
760       INTEGER ::  ncid, ndim_out, start(4), count(4)
761       LOGICAL ::  last_dimension_is_time
762
763       ndim_out = var % ndim
764
765       last_dimension_is_time = var % dimids( var % ndim ) == output_file % dimid_time
766       IF ( last_dimension_is_time )  THEN
767          ndim_out = ndim_out - 1
768       END IF
769
770       start(:)      = (/1,1,1,iter/)
771       count(1:ndim_out) = var%dimlen(1:ndim_out)
772
773       CALL check(nf90_open(output_file % name, NF90_WRITE, ncid))
774
775       IF (TRIM(var % kind) .EQ. 'left/right scalar')  THEN
776
777          CALL check(nf90_put_var( ncid, var%varid, buffer(1,:,:) ) )
778
779       ELSE IF (TRIM(var % kind) .EQ. 'north/south scalar')  THEN
780
781          CALL check(nf90_put_var( ncid, var%varid, buffer(:,1,:) ) )
782
783       ELSE IF (TRIM(var % kind) .EQ. 'top scalar')  THEN
784
785          CALL check(nf90_put_var( ncid, var%varid, buffer(:,:,1) ) )
786       ELSE
787
788          CALL check(nf90_put_var( ncid, var%varid, buffer ) )
789
790       END IF
791       CALL check(nf90_close(ncid))
792
793    END SUBROUTINE write_netcdf_variable_2d
794
795
796    SUBROUTINE check(status)
797
798       INTEGER, INTENT(IN) ::  status
799
800       IF (status /= nf90_noerr)  THEN
801          message = "NetCDF API call failed with error: " //                     &
802                    TRIM( nf90_strerror(status) )
803          CALL abort('io.check', message) 
804       END IF
805
806    END SUBROUTINE check
807
808 END MODULE io
Note: See TracBrowser for help on using the repository browser.