!> @file src/io.f90 !------------------------------------------------------------------------------! ! This file is part of the PALM model system. ! ! PALM is free software: you can redistribute it and/or modify it under the ! terms of the GNU General Public License as published by the Free Software ! Foundation, either version 3 of the License, or (at your option) any later ! version. ! ! PALM is distributed in the hope that it will be useful, but WITHOUT ANY ! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR ! A PARTICULAR PURPOSE. See the GNU General Public License for more details. ! ! You should have received a copy of the GNU General Public License along with ! PALM. If not, see . ! ! Copyright 2017-2018 Leibniz Universitaet Hannover ! Copyright 2017-2018 Deutscher Wetterdienst Offenbach !------------------------------------------------------------------------------! ! ! Current revisions: ! ----------------- ! ! ! Former revisions: ! ----------------- ! $Id: io.f90 3395 2018-10-22 17:32:49Z gronemeier $ ! Added command-line options for configuring the computation of geostrophic ! winds (--averagin-mode, --averaging-angle) ! Added command-line option --input-prefix for setting input file prefixes all ! at once ! Added --debug option for more verbose terminal output ! Added option-specific *_is_set LOGICALs to indicate invocation from the ! command-line ! Improved error messages in case of empty file-name strings ! Improved routine naming ! ! 3183 2018-07-27 14:25:55Z suehring ! Introduced new PALM grid stretching ! Updated variable names and metadata for PIDS v1.9 compatibility ! Improved handling of the start date string ! Better compatibility with older Intel compilers: ! - avoiding implicit array allocation with new get_netcdf_variable() ! subroutine instead of function ! Improved command line interface: ! - Added configuration validation ! - New options to configure input file prefixes ! - GNU-style short and long option names ! - Added version and copyright output ! ! ! 3182 2018-07-27 13:36:03Z suehring ! Initial revision ! ! ! ! Authors: ! -------- ! @author Eckhard Kadasch ! ! Description: ! ------------ !> The io module contains the functions needed to read and write netCDF data in !> INIFOR. !------------------------------------------------------------------------------! MODULE io USE control USE defs, & ONLY: DATE, SNAME, PATH, PI, dp, hp, TO_RADIANS, TO_DEGREES, VERSION USE netcdf USE types USE util, & ONLY: reverse, str, real_to_str IMPLICIT NONE INTERFACE get_netcdf_variable MODULE PROCEDURE get_netcdf_variable_int MODULE PROCEDURE get_netcdf_variable_real END INTERFACE get_netcdf_variable PRIVATE :: get_netcdf_variable_int, get_netcdf_variable_real CONTAINS SUBROUTINE get_netcdf_variable_int(in_file, in_var, buffer) CHARACTER(LEN=PATH), INTENT(IN) :: in_file TYPE(nc_var), INTENT(INOUT) :: in_var INTEGER(hp), ALLOCATABLE, INTENT(INOUT) :: buffer(:,:,:) INCLUDE 'get_netcdf_variable.inc' END SUBROUTINE get_netcdf_variable_int SUBROUTINE get_netcdf_variable_real(in_file, in_var, buffer) CHARACTER(LEN=PATH), INTENT(IN) :: in_file TYPE(nc_var), INTENT(INOUT) :: in_var REAL(dp), ALLOCATABLE, INTENT(INOUT) :: buffer(:,:,:) INCLUDE 'get_netcdf_variable.inc' END SUBROUTINE get_netcdf_variable_real SUBROUTINE netcdf_define_variable(var, ncid) TYPE(nc_var), INTENT(INOUT) :: var INTEGER, INTENT(IN) :: ncid CALL check(nf90_def_var(ncid, var % name, NF90_FLOAT, var % dimids(1:var % ndim), var % varid)) CALL check(nf90_put_att(ncid, var % varid, "long_name", var % long_name)) CALL check(nf90_put_att(ncid, var % varid, "units", var % units)) IF ( var % lod .GE. 0 ) THEN CALL check(nf90_put_att(ncid, var % varid, "lod", var % lod)) END IF CALL check(nf90_put_att(ncid, var % varid, "source", var % source)) CALL check(nf90_put_att(ncid, var % varid, "_FillValue", NF90_FILL_REAL)) END SUBROUTINE netcdf_define_variable SUBROUTINE netcdf_get_dimensions(var, ncid) TYPE(nc_var), INTENT(INOUT) :: var INTEGER, INTENT(IN) :: ncid INTEGER :: i CHARACTER(SNAME) :: null ! Remember dimension lenghts for NetCDF writing routine DO i = 1, var % ndim CALL check(nf90_inquire_dimension(ncid, var % dimids(i), & name = null, & len = var % dimlen(i) ) ) END DO END SUBROUTINE netcdf_get_dimensions !------------------------------------------------------------------------------! ! Description: ! ------------ !> This routine initializes Inifor. This includes parsing command-line !> arguments, setting the names of the input and output file names as well as !> the name of the input namelist and, subsequently, reading in and setting grid !> parameters for the PALM-4U computational grid. !------------------------------------------------------------------------------! SUBROUTINE parse_command_line_arguments( cfg ) TYPE(inifor_config), INTENT(INOUT) :: cfg CHARACTER(LEN=PATH) :: option, arg INTEGER :: arg_count, i cfg % p0_is_set = .FALSE. cfg % ug_is_set = .FALSE. cfg % vg_is_set = .FALSE. cfg % flow_prefix_is_set = .FALSE. cfg % input_prefix_is_set = .FALSE. cfg % radiation_prefix_is_set = .FALSE. cfg % soil_prefix_is_set = .FALSE. cfg % soilmoisture_prefix_is_set = .FALSE. arg_count = COMMAND_ARGUMENT_COUNT() IF (arg_count .GT. 0) THEN ! Every option should have an argument. !IF ( MOD(arg_count, 2) .NE. 0 ) THEN ! message = "Syntax error in command line." ! CALL abort('parse_command_line_arguments', message) !END IF message = "The -clon and -clat command line options are depricated. " // & "Please remove them form your inifor command and specify the " // & "location of the PALM-4U origin either" // NEW_LINE(' ') // & " - by setting the namelist parameters 'longitude' and 'latitude', or" // NEW_LINE(' ') // & " - by providing a static driver netCDF file via the -static command-line option." i = 1 DO WHILE (i .LE. arg_count) CALL GET_COMMAND_ARGUMENT( i, option ) SELECT CASE( TRIM(option) ) CASE( '--averaging-mode' ) CALL get_option_argument( i, arg ) cfg % averaging_mode = TRIM(arg) CASE( '-date', '-d', '--date' ) CALL get_option_argument( i, arg ) cfg % start_date = TRIM(arg) CASE( '--debug' ) cfg % debug = .TRUE. ! Elevation of the PALM-4U domain above sea level CASE( '-z0', '-z', '--elevation' ) CALL get_option_argument( i, arg ) READ(arg, *) cfg % z0 ! surface pressure, at z0 CASE( '-p0', '-r', '--surface-pressure' ) cfg % p0_is_set = .TRUE. CALL get_option_argument( i, arg ) READ(arg, *) cfg % p0 ! geostrophic wind in x direction CASE( '-ug', '-u', '--geostrophic-u' ) cfg % ug_is_set = .TRUE. CALL get_option_argument( i, arg ) READ(arg, *) cfg % ug ! geostrophic wind in y direction CASE( '-vg', '-v', '--geostrophic-v' ) cfg % vg_is_set = .TRUE. CALL get_option_argument( i, arg ) READ(arg, *) cfg % vg ! domain centre geographical longitude and latitude CASE( '-clon', '-clat' ) CALL abort('parse_command_line_arguments', message) !READ(arg, *) lambda_cg !lambda_cg = lambda_cg * TO_RADIANS !READ(arg, *) phi_cg !phi_cg = phi_cg * TO_RADIANS CASE( '-path', '-p', '--path' ) CALL get_option_argument( i, arg ) cfg % input_path = TRIM(arg) CASE( '-hhl', '-l', '--hhl-file' ) CALL get_option_argument( i, arg ) cfg % hhl_file = TRIM(arg) CASE( '--input-prefix') CALL get_option_argument( i, arg ) cfg % input_prefix = TRIM(arg) cfg % input_prefix_is_set = .TRUE. CASE( '-a', '--averaging-angle' ) CALL get_option_argument( i, arg ) READ(arg, *) cfg % averaging_angle CASE( '-static', '-t', '--static-driver' ) CALL get_option_argument( i, arg ) cfg % static_driver_file = TRIM(arg) CASE( '-soil', '-s', '--soil-file') CALL get_option_argument( i, arg ) cfg % soiltyp_file = TRIM(arg) CASE( '--flow-prefix') CALL get_option_argument( i, arg ) cfg % flow_prefix = TRIM(arg) cfg % flow_prefix_is_set = .TRUE. CASE( '--radiation-prefix') CALL get_option_argument( i, arg ) cfg % radiation_prefix = TRIM(arg) cfg % radiation_prefix_is_set = .TRUE. CASE( '--soil-prefix') CALL get_option_argument( i, arg ) cfg % soil_prefix = TRIM(arg) cfg % soil_prefix_is_set = .TRUE. CASE( '--soilmoisture-prefix') CALL get_option_argument( i, arg ) cfg % soilmoisture_prefix = TRIM(arg) cfg % soilmoisture_prefix_is_set = .TRUE. CASE( '-o', '--output' ) CALL get_option_argument( i, arg ) cfg % output_file = TRIM(arg) CASE( '-n', '--namelist' ) CALL get_option_argument( i, arg ) cfg % namelist_file = TRIM(arg) ! initial condition mode: 'profile' / 'volume' CASE( '-mode', '-i', '--init-mode' ) CALL get_option_argument( i, arg ) cfg % ic_mode = TRIM(arg) ! boundary conditions / forcing mode: 'ideal' / 'real' CASE( '-f', '--forcing-mode' ) CALL get_option_argument( i, arg ) cfg % bc_mode = TRIM(arg) CASE( '--version' ) CALL print_version() STOP CASE( '--help' ) CALL print_version() PRINT *, "" PRINT *, "For a list of command-line options have a look at the README file." STOP CASE DEFAULT message = "unknown option '" // TRIM(option) // "'." CALL abort('parse_command_line_arguments', message) END SELECT i = i + 1 END DO ELSE message = "No arguments present, using default input and output files" CALL report('parse_command_line_arguments', message) END IF END SUBROUTINE parse_command_line_arguments SUBROUTINE get_option_argument(i, arg) CHARACTER(LEN=PATH), INTENT(INOUT) :: arg INTEGER, INTENT(INOUT) :: i i = i + 1 CALL GET_COMMAND_ARGUMENT(i, arg) END SUBROUTINE SUBROUTINE validate_config(cfg) TYPE(inifor_config), INTENT(IN) :: cfg LOGICAL :: all_files_present all_files_present = .TRUE. all_files_present = all_files_present .AND. file_present(cfg % hhl_file, 'HHL') all_files_present = all_files_present .AND. file_present(cfg % namelist_file, 'NAMELIST') all_files_present = all_files_present .AND. file_present(cfg % soiltyp_file, 'SOILTYP') ! Only check optional static driver file name, if it has been given. IF (TRIM(cfg % static_driver_file) .NE. '') THEN all_files_present = all_files_present .AND. file_present(cfg % static_driver_file, 'static driver') END IF IF (.NOT. all_files_present) THEN message = "INIFOR configuration invalid; some input files are missing." CALL abort( 'validate_config', message ) END IF SELECT CASE( TRIM(cfg % ic_mode) ) CASE( 'profile', 'volume') CASE DEFAULT message = "Initialization mode '" // TRIM(cfg % ic_mode) //& "' is not supported. " //& "Please select either 'profile' or 'volume', " //& "or omit the -i/--init-mode/-mode option entirely, which corresponds "//& "to the latter." CALL abort( 'validate_config', message ) END SELECT SELECT CASE( TRIM(cfg % bc_mode) ) CASE( 'real', 'ideal') CASE DEFAULT message = "Forcing mode '" // TRIM(cfg % bc_mode) //& "' is not supported. " //& "Please select either 'real' or 'ideal', " //& "or omit the -f/--forcing-mode option entirely, which corresponds "//& "to the latter." CALL abort( 'validate_config', message ) END SELECT SELECT CASE( TRIM(cfg % averaging_mode) ) CASE( 'level', 'height') CASE DEFAULT message = "Averaging mode '" // TRIM(cfg % averaging_mode) //& "' is not supported. " //& "Please select either 'height' or 'level', " //& "or omit the --averaging-mode option entirely, which corresponds "//& "to the latter." CALL abort( 'validate_config', message ) END SELECT IF ( cfg % ug_is_set .NEQV. cfg % vg_is_set ) THEN message = "You specified only one component of the geostrophic " // & "wind. Please specify either both or none." CALL abort( 'validate_config', message ) END IF END SUBROUTINE validate_config LOGICAL FUNCTION file_present(filename, kind) CHARACTER(LEN=PATH), INTENT(IN) :: filename CHARACTER(LEN=*), INTENT(IN) :: kind IF (LEN(TRIM(filename))==0) THEN file_present = .FALSE. message = "No name was given for the " // TRIM(kind) // " file." CALL report('file_present', message) ELSE INQUIRE(FILE=filename, EXIST=file_present) IF (.NOT. file_present) THEN message = "The given " // TRIM(kind) // " file '" // & TRIM(filename) // "' does not exist." CALL report('file_present', message) END IF END IF END FUNCTION file_present !------------------------------------------------------------------------------! ! Description: ! ------------ !> This routine initializes the Inifor output file, i.e. the PALM-4U !> initializing and forcing data as a NetCDF file. !> !> Besides writing metadata, such as global attributes, coordinates, variables, !> in the NetCDF file, various NetCDF IDs are saved for later, when Inifor !> writes the actual data. !------------------------------------------------------------------------------! SUBROUTINE setup_netcdf_dimensions(output_file, palm_grid, & start_date_string, origin_lon, origin_lat) TYPE(nc_file), INTENT(INOUT) :: output_file TYPE(grid_definition), INTENT(IN) :: palm_grid CHARACTER (LEN=DATE), INTENT(IN) :: start_date_string REAL(dp), INTENT(IN) :: origin_lon, origin_lat CHARACTER (LEN=8) :: date_string CHARACTER (LEN=10) :: time_string CHARACTER (LEN=5) :: zone_string CHARACTER (LEN=SNAME) :: history_string INTEGER :: ncid, nx, ny, nz, nt, dimids(3), dimvarids(3) REAL(dp) :: z0 message = "Initializing PALM-4U dynamic driver file '" // & TRIM(output_file % name) // "' and setting up dimensions." CALL report('setup_netcdf_dimensions', message) ! Create the NetCDF file. NF90_CLOBBER selects overwrite mode. #if defined( __netcdf4 ) CALL check(nf90_create(TRIM(output_file % name), OR(NF90_CLOBBER, NF90_HDF5), ncid)) #else CALL check(nf90_create(TRIM(output_file % name), NF90_CLOBBER, ncid)) #endif !------------------------------------------------------------------------------ !- Section 1: Define NetCDF dimensions and coordinates !------------------------------------------------------------------------------ nt = SIZE(output_file % time) nx = palm_grid % nx ny = palm_grid % ny nz = palm_grid % nz z0 = palm_grid % z0 ! !------------------------------------------------------------------------------ !- Section 2: Write global NetCDF attributes !------------------------------------------------------------------------------ CALL date_and_time(DATE=date_string, TIME=time_string, ZONE=zone_string) history_string = & 'Created on '// date_string // & ' at ' // time_string(1:2) // ':' // time_string(3:4) // & ' (UTC' // zone_string // ')' CALL check(nf90_put_att(ncid, NF90_GLOBAL, 'title', 'PALM input file for scenario ...')) CALL check(nf90_put_att(ncid, NF90_GLOBAL, 'institution', 'Deutscher Wetterdienst, Offenbach')) CALL check(nf90_put_att(ncid, NF90_GLOBAL, 'author', 'Eckhard Kadasch, eckhard.kadasch@dwd.de')) CALL check(nf90_put_att(ncid, NF90_GLOBAL, 'history', TRIM(history_string))) CALL check(nf90_put_att(ncid, NF90_GLOBAL, 'references', '--')) CALL check(nf90_put_att(ncid, NF90_GLOBAL, 'comment', '--')) CALL check(nf90_put_att(ncid, NF90_GLOBAL, 'origin_lat', TRIM(real_to_str(origin_lat*TO_DEGREES, '(F18.13)')))) CALL check(nf90_put_att(ncid, NF90_GLOBAL, 'origin_lon', TRIM(real_to_str(origin_lon*TO_DEGREES, '(F18.13)')))) ! FIXME: This is the elevation relative to COSMO-DE/D2 sea level and does ! FIXME: not necessarily comply with DHHN2016 (c.f. PALM Input Data ! FIXME: Standard v1.9., origin_z) CALL check(nf90_put_att(ncid, NF90_GLOBAL, 'origin_z', TRIM(real_to_str(z0, '(F18.13)')))) CALL check(nf90_put_att(ncid, NF90_GLOBAL, 'inifor_version', TRIM(VERSION))) CALL check(nf90_put_att(ncid, NF90_GLOBAL, 'palm_version', '--')) ! ! !------------------------------------------------------------------------------ !- Section 2a: Define dimensions for cell centers (scalars in soil and atmosph.) !------------------------------------------------------------------------------ dimids = (/0, 0, 0/) ! reset dimids CALL check( nf90_def_dim(ncid, "x", nx+1, dimids(1)) ) CALL check( nf90_def_dim(ncid, "y", ny+1, dimids(2)) ) CALL check( nf90_def_dim(ncid, "z", nz, dimids(3)) ) output_file % dimids_scl = dimids ! save dimids for later dimvarids = (/0, 0, 0/) ! reset dimvarids CALL check(nf90_def_var(ncid, "x", NF90_FLOAT, dimids(1), dimvarids(1))) CALL check(nf90_put_att(ncid, dimvarids(1), "standard_name", "x coordinate of cell centers")) CALL check(nf90_put_att(ncid, dimvarids(1), "units", "m")) CALL check(nf90_def_var(ncid, "y", NF90_FLOAT, dimids(2), dimvarids(2))) CALL check(nf90_put_att(ncid, dimvarids(2), "standard_name", "y coordinate of cell centers")) CALL check(nf90_put_att(ncid, dimvarids(2), "units", "m")) CALL check(nf90_def_var(ncid, "z", NF90_FLOAT, dimids(3), dimvarids(3))) CALL check(nf90_put_att(ncid, dimvarids(3), "standard_name", "z coordinate of cell centers")) CALL check(nf90_put_att(ncid, dimvarids(3), "units", "m")) output_file % dimvarids_scl = dimvarids ! save dimvarids for later ! overwrite third dimid with the one of depth CALL check(nf90_def_dim(ncid, "zsoil", SIZE(palm_grid % depths), dimids(3)) ) output_file % dimids_soil = dimids ! save dimids for later ! overwrite third dimvarid with the one of depth CALL check(nf90_def_var(ncid, "zsoil", NF90_FLOAT, output_file % dimids_soil(3), dimvarids(3))) CALL check(nf90_put_att(ncid, dimvarids(3), "standard_name", "depth_below_land")) CALL check(nf90_put_att(ncid, dimvarids(3), "positive", "down")) CALL check(nf90_put_att(ncid, dimvarids(3), "units", "m")) output_file % dimvarids_soil = dimvarids ! save dimvarids for later ! !------------------------------------------------------------------------------ !- Section 2b: Define dimensions for cell faces/velocities !------------------------------------------------------------------------------ dimids = (/0, 0, 0/) ! reset dimids CALL check(nf90_def_dim(ncid, "xu", nx, dimids(1)) ) CALL check(nf90_def_dim(ncid, "yv", ny, dimids(2)) ) CALL check(nf90_def_dim(ncid, "zw", nz-1, dimids(3)) ) output_file % dimids_vel = dimids ! save dimids for later dimvarids = (/0, 0, 0/) ! reset dimvarids CALL check(nf90_def_var(ncid, "xu", NF90_FLOAT, dimids(1), dimvarids(1))) CALL check(nf90_put_att(ncid, dimvarids(1), "standard_name", "x coordinate of cell faces")) CALL check(nf90_put_att(ncid, dimvarids(1), "units", "m")) CALL check(nf90_def_var(ncid, "yv", NF90_FLOAT, dimids(2), dimvarids(2))) CALL check(nf90_put_att(ncid, dimvarids(2), "standard_name", "y coordinate of cell faces")) CALL check(nf90_put_att(ncid, dimvarids(2), "units", "m")) CALL check(nf90_def_var(ncid, "zw", NF90_FLOAT, dimids(3), dimvarids(3))) CALL check(nf90_put_att(ncid, dimvarids(3), "standard_name", "z coordinate of cell faces")) CALL check(nf90_put_att(ncid, dimvarids(3), "units", "m")) output_file % dimvarids_vel = dimvarids ! save dimvarids for later ! !------------------------------------------------------------------------------ !- Section 2c: Define time dimension !------------------------------------------------------------------------------ CALL check(nf90_def_dim(ncid, "time", nt, output_file % dimid_time) ) CALL check(nf90_def_var(ncid, "time", NF90_FLOAT, & output_file % dimid_time, & output_file % dimvarid_time)) CALL check(nf90_put_att(ncid, output_file % dimvarid_time, "standard_name", "time")) CALL check(nf90_put_att(ncid, output_file % dimvarid_time, "long_name", "time")) CALL check(nf90_put_att(ncid, output_file % dimvarid_time, "units", & "seconds since " // start_date_string // " UTC")) CALL check(nf90_enddef(ncid)) ! !------------------------------------------------------------------------------ !- Section 3: Write grid coordinates !------------------------------------------------------------------------------ CALL check(nf90_put_var(ncid, output_file % dimvarids_scl(1), palm_grid%x)) CALL check(nf90_put_var(ncid, output_file % dimvarids_scl(2), palm_grid%y)) CALL check(nf90_put_var(ncid, output_file % dimvarids_scl(3), palm_grid%z)) CALL check(nf90_put_var(ncid, output_file % dimvarids_vel(1), palm_grid%xu)) CALL check(nf90_put_var(ncid, output_file % dimvarids_vel(2), palm_grid%yv)) CALL check(nf90_put_var(ncid, output_file % dimvarids_vel(3), palm_grid%zw)) ! TODO Read in soil depths from input file before this. CALL check(nf90_put_var(ncid, output_file % dimvarids_soil(3), palm_grid%depths)) ! Write time vector CALL check(nf90_put_var(ncid, output_file % dimvarid_time, output_file % time)) ! Close the file CALL check(nf90_close(ncid)) END SUBROUTINE setup_netcdf_dimensions SUBROUTINE setup_netcdf_variables(filename, output_variable_table) CHARACTER (LEN=*), INTENT(IN) :: filename TYPE(nc_var), INTENT(INOUT), TARGET :: output_variable_table(:) TYPE(nc_var), POINTER :: var INTEGER :: i, ncid message = "Defining variables in dynamic driver '" // TRIM(filename) // "'." CALL report('setup_netcdf_variables', message) CALL check(nf90_open(TRIM(filename), NF90_WRITE, ncid)) CALL check(nf90_redef(ncid)) DO i = 1, SIZE(output_variable_table) var => output_variable_table(i) IF ( var % to_be_processed ) THEN message = " variable #" // TRIM(str(i)) // " '" // TRIM(var%name) // "'." CALL report('setup_netcdf_variables', message) CALL netcdf_define_variable(var, ncid) CALL netcdf_get_dimensions(var, ncid) END IF END DO CALL check(nf90_enddef(ncid)) CALL check(nf90_close(ncid)) message = "Dynamic driver '" // TRIM(filename) // "' initialized successfully." CALL report('setup_netcdf_variables', message) END SUBROUTINE setup_netcdf_variables !------------------------------------------------------------------------------! ! Description: ! ------------ !> This routine reads and returns all input variables of the given IO group !> It accomodates the data by allocating a container variable that represents a !> list of arrays of the same length as the groups variable list. (This list !> will typically contain one or two items.) After the container, its members !> are allocated one by one with the appropriate, possibly different, !> dimensions. !> !> The 'group' is an INTENT(INOUT) variable so that 'get_netcdf_variable()' can !> record netCDF IDs in the 'in_var_list()' member variable. !------------------------------------------------------------------------------! SUBROUTINE read_input_variables(group, iter, buffer) TYPE(io_group), INTENT(INOUT), TARGET :: group INTEGER, INTENT(IN) :: iter TYPE(container), ALLOCATABLE, INTENT(INOUT) :: buffer(:) INTEGER :: hour, buf_id TYPE(nc_var), POINTER :: input_var CHARACTER(LEN=PATH), POINTER :: input_file INTEGER :: ivar, nbuffers message = "Reading data for I/O group '" // TRIM(group % in_var_list(1) % name) // "'." CALL report('read_input_variables', message) input_file => group % in_files(iter) ! !------------------------------------------------------------------------------ !- Section 1: Load input buffers for accumulated variables !------------------------------------------------------------------------------ IF (group % kind == 'running average' .OR. & group % kind == 'accumulated') THEN ! radiation budgets, precipitation IF (SIZE(group % in_var_list) .GT. 1 ) THEN message = "I/O groups may not contain more than one " // & "accumulated variable. Group '" // TRIM(group % kind) //& "' contains " // & TRIM( str(SIZE(group % in_var_list)) ) // "." CALL abort('read_input_variables | accumulation', message) END IF ! use two buffer arrays nbuffers = 2 IF ( .NOT. ALLOCATED( buffer ) ) ALLOCATE( buffer(nbuffers) ) ! chose correct buffer array hour = iter - 1! hour of the day buf_id = select_buffer(hour) CALL run_control('time', 'read') IF ( ALLOCATED(buffer(buf_id) % array) ) DEALLOCATE(buffer(buf_id) % array) CALL run_control('time', 'alloc') input_var => group % in_var_list(1) CALL get_netcdf_variable(input_file, input_var, buffer(buf_id) % array) CALL report('read_input_variables', "Read accumulated " // TRIM(group % in_var_list(1) % name)) IF ( input_var % is_upside_down ) CALL reverse(buffer(buf_id) % array) CALL run_control('time', 'comp') !------------------------------------------------------------------------------ !- Section 2: Load input buffers for normal I/O groups !------------------------------------------------------------------------------ ELSE ! Allocate one input buffer per input_variable. If more quantities ! have to be computed than input variables exist in this group, ! allocate more buffers. For instance, in the thermodynamics group, ! there are three input variabels (PP, T, Qv) and four quantities ! necessart (P, Theta, Rho, qv) for the corresponding output fields ! (p0, Theta, qv, ug, and vg) nbuffers = MAX( group % n_inputs, group % n_output_quantities ) ALLOCATE( buffer(nbuffers) ) CALL run_control('time', 'alloc') ! Read in all input variables, leave extra buffers-if any-untouched. DO ivar = 1, group % n_inputs input_var => group % in_var_list(ivar) ! Check wheather P or PP is present in input file IF (input_var % name == 'P') THEN input_var % name = TRIM( get_pressure_varname(input_file) ) CALL run_control('time', 'read') END IF CALL get_netcdf_variable(input_file, input_var, buffer(ivar) % array) IF ( input_var % is_upside_down ) CALL reverse(buffer(ivar) % array) CALL run_control('time', 'comp') END DO END IF END SUBROUTINE read_input_variables INTEGER FUNCTION select_buffer(hour) INTEGER, INTENT(IN) :: hour INTEGER :: step select_buffer = 0 step = MODULO(hour, 3) + 1 SELECT CASE(step) CASE(1, 3) select_buffer = 1 CASE(2) select_buffer = 2 CASE DEFAULT message = "Invalid step '" // TRIM(str(step)) CALL abort('select_buffer', message) END SELECT END FUNCTION select_buffer !------------------------------------------------------------------------------! ! Description: ! ------------ !> Checks if the input_file contains the absolute pressure, 'P', or the pressure !> perturbation, 'PP', and returns the appropriate string. !------------------------------------------------------------------------------! CHARACTER(LEN=2) FUNCTION get_pressure_varname(input_file) RESULT(var) CHARACTER(LEN=*) :: input_file INTEGER :: ncid, varid CALL check(nf90_open( TRIM(input_file), NF90_NOWRITE, ncid )) IF ( nf90_inq_varid( ncid, 'P', varid ) .EQ. NF90_NOERR ) THEN var = 'P' ELSE IF ( nf90_inq_varid( ncid, 'PP', varid ) .EQ. NF90_NOERR ) THEN var = 'PP' CALL report('get_pressure_var', 'Using PP instead of P') ELSE message = "Failed to read '" // TRIM(var) // & "' from file '" // TRIM(input_file) // "'." CALL abort('get_pressure_var', message) END IF CALL check(nf90_close(ncid)) END FUNCTION get_pressure_varname FUNCTION get_netcdf_attribute(filename, attribute) RESULT(attribute_value) CHARACTER(LEN=*), INTENT(IN) :: filename, attribute REAL(dp) :: attribute_value INTEGER :: ncid IF ( nf90_open( TRIM(filename), NF90_NOWRITE, ncid ) == NF90_NOERR ) THEN CALL check(nf90_get_att(ncid, NF90_GLOBAL, TRIM(attribute), attribute_value)) CALL check(nf90_close(ncid)) ELSE message = "Failed to read '" // TRIM(attribute) // & "' from file '" // TRIM(filename) // "'." CALL abort('get_netcdf_attribute', message) END IF END FUNCTION get_netcdf_attribute SUBROUTINE update_output(var, array, iter, output_file) TYPE(nc_var), INTENT(IN) :: var REAL(dp), INTENT(IN) :: array(:,:,:) INTEGER, INTENT(IN) :: iter TYPE(nc_file), INTENT(IN) :: output_file INTEGER :: ncid, ndim, start(4), count(4) LOGICAL :: var_is_time_dependent var_is_time_dependent = ( & var % dimids( var % ndim ) == output_file % dimid_time & ) ! Skip time dimension for output ndim = var % ndim IF ( var_is_time_dependent ) ndim = var % ndim - 1 start(:) = (/1,1,1,1/) start(ndim+1) = iter count(1:ndim) = var%dimlen(1:ndim) CALL check(nf90_open(output_file % name, NF90_WRITE, ncid)) ! Reduce dimension of output array according to variable kind SELECT CASE (TRIM(var % kind)) CASE ( 'init scalar profile', 'init u profile', 'init v profile', & 'init w profile' ) CALL check(nf90_put_var( ncid, var%varid, array(1,1,:) ) ) CASE ( 'init soil', 'init scalar', 'init u', 'init v', 'init w' ) CALL check(nf90_put_var( ncid, var%varid, array(:,:,:) ) ) CASE( 'left scalar', 'right scalar', 'left w', 'right w' ) CALL check(nf90_put_var( ncid, var%varid, array(1,:,:), & start=start(1:ndim+1), & count=count(1:ndim) ) ) IF (.NOT. SIZE(array, 2) .EQ. var % dimlen(1)) THEN PRINT *, "inifor: update_output: Dimension ", 1, " of variable ", & TRIM(var % name), " (", var % dimlen(1), & ") does not match the dimension of the output array (", & SIZE(array, 2), ")." STOP END IF CASE( 'north scalar', 'south scalar', 'north w', 'south w' ) CALL check(nf90_put_var( ncid, var%varid, array(:,1,:), & start=start(1:ndim+1), & count=count(1:ndim) ) ) CASE( 'surface forcing', 'top scalar', 'top w' ) CALL check(nf90_put_var( ncid, var%varid, array(:,:,1), & start=start(1:ndim+1), & count=count(1:ndim) ) ) CASE ( 'left u', 'right u', 'left v', 'right v' ) CALL check(nf90_put_var( ncid, var%varid, array(1,:,:), & start=start(1:ndim+1), & count=count(1:ndim) ) ) CASE ( 'north u', 'south u', 'north v', 'south v' ) CALL check(nf90_put_var( ncid, var%varid, array(:,1,:), & start=start(1:ndim+1), & count=count(1:ndim) ) ) CASE ( 'top u', 'top v' ) CALL check(nf90_put_var( ncid, var%varid, array(:,:,1), & start=start(1:ndim+1), & count=count(1:ndim) ) ) CASE ( 'time series' ) CALL check(nf90_put_var( ncid, var%varid, array(1,1,1), & start=start(1:ndim+1) ) ) CASE ( 'constant scalar profile', 'geostrophic', 'internal profile' ) CALL check(nf90_put_var( ncid, var%varid, array(1,1,:), & start=start(1:ndim+1), & count=count(1:ndim) ) ) CASE ( 'large-scale scalar forcing', 'large-scale w forcing' ) message = "Doing nothing in terms of writing large-scale forings." CALL report('update_output', message) CASE DEFAULT message = "Variable kind '" // TRIM(var % kind) // & "' not recognized." CALL abort('update_output', message) END SELECT CALL check(nf90_close(ncid)) END SUBROUTINE update_output SUBROUTINE write_netcdf_variable_2d(var, iter, output_file, buffer) TYPE(nc_var), INTENT(IN) :: var INTEGER, INTENT(IN) :: iter TYPE(nc_file), INTENT(IN) :: output_file REAL(dp), INTENT(IN) :: buffer(:,:,:) INTEGER :: ncid, ndim_out, start(4), count(4) LOGICAL :: last_dimension_is_time ndim_out = var % ndim last_dimension_is_time = var % dimids( var % ndim ) == output_file % dimid_time IF ( last_dimension_is_time ) THEN ndim_out = ndim_out - 1 END IF start(:) = (/1,1,1,iter/) count(1:ndim_out) = var%dimlen(1:ndim_out) CALL check(nf90_open(output_file % name, NF90_WRITE, ncid)) IF (TRIM(var % kind) .EQ. 'left/right scalar') THEN CALL check(nf90_put_var( ncid, var%varid, buffer(1,:,:) ) ) ELSE IF (TRIM(var % kind) .EQ. 'north/south scalar') THEN CALL check(nf90_put_var( ncid, var%varid, buffer(:,1,:) ) ) ELSE IF (TRIM(var % kind) .EQ. 'top scalar') THEN CALL check(nf90_put_var( ncid, var%varid, buffer(:,:,1) ) ) ELSE CALL check(nf90_put_var( ncid, var%varid, buffer ) ) END IF CALL check(nf90_close(ncid)) END SUBROUTINE write_netcdf_variable_2d SUBROUTINE check(status) INTEGER, INTENT(IN) :: status IF (status /= nf90_noerr) THEN message = "NetCDF API call failed with error: " // & TRIM( nf90_strerror(status) ) CALL abort('io.check', message) END IF END SUBROUTINE check END MODULE io