SUBROUTINE check_parameters !--------------------------------------------------------------------------------! ! This file is part of PALM. ! ! 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 1997-2012 Leibniz University Hannover !--------------------------------------------------------------------------------! ! ! Current revisions: ! ----------------- ! ! ! Former revisions: ! ----------------- ! $Id: check_parameters.f90 1217 2013-08-26 11:30:05Z raasch $ ! ! 1216 2013-08-26 09:31:42Z raasch ! check for transpose_compute_overlap (temporary) ! ! 1214 2013-08-21 12:29:17Z kanani ! additional check for simultaneous use of vertical grid stretching ! and particle advection ! ! 1212 2013-08-15 08:46:27Z raasch ! checks for poisfft_hybrid removed ! ! 1210 2013-08-14 10:58:20Z raasch ! check for fftw ! ! 1179 2013-06-14 05:57:58Z raasch ! checks and settings of buoyancy parameters and switches revised, ! initial profile for rho added to hom (id=77) ! ! 1174 2013-05-31 10:28:08Z gryschka ! Bugfix in computing initial profiles for ug, vg, lad, q in case of Atmosphere ! ! 1159 2013-05-21 11:58:22Z fricke ! bc_lr/ns_dirneu/neudir removed ! ! 1115 2013-03-26 18:16:16Z hoffmann ! unused variables removed ! drizzle can be used without precipitation ! ! 1111 2013-03-08 23:54:10Z raasch ! ibc_p_b = 2 removed ! ! 1103 2013-02-20 02:15:53Z raasch ! Bugfix: turbulent inflow must not require cyclic fill in restart runs ! ! 1092 2013-02-02 11:24:22Z raasch ! unused variables removed ! ! 1069 2012-11-28 16:18:43Z maronga ! allow usage of topography in combination with cloud physics ! ! 1065 2012-11-22 17:42:36Z hoffmann ! Bugfix: It is not allowed to use cloud_scheme = seifert_beheng without ! precipitation in order to save computational resources. ! ! 1060 2012-11-21 07:19:51Z raasch ! additional check for parameter turbulent_inflow ! ! 1053 2012-11-13 17:11:03Z hoffmann ! necessary changes for the new two-moment cloud physics scheme added: ! - check cloud physics scheme (Kessler or Seifert and Beheng) ! - plant_canopy is not allowed ! - currently, only cache loop_optimization is allowed ! - initial profiles of nr, qr ! - boundary condition of nr, qr ! - check output quantities (qr, nr, prr) ! ! 1036 2012-10-22 13:43:42Z raasch ! code put under GPL (PALM 3.9) ! ! 1031/1034 2012-10-22 11:32:49Z raasch ! check of netcdf4 parallel file support ! ! 1019 2012-09-28 06:46:45Z raasch ! non-optimized version of prognostic_equations not allowed any more ! ! 1015 2012-09-27 09:23:24Z raasch ! acc allowed for loop optimization, ! checks for adjustment of mixing length to the Prandtl mixing length removed ! ! 1003 2012-09-14 14:35:53Z raasch ! checks for cases with unequal subdomain sizes removed ! ! 1001 2012-09-13 14:08:46Z raasch ! all actions concerning leapfrog- and upstream-spline-scheme removed ! ! 996 2012-09-07 10:41:47Z raasch ! little reformatting ! 978 2012-08-09 08:28:32Z fricke ! setting of bc_lr/ns_dirneu/neudir ! outflow damping layer removed ! check for z0h* ! check for pt_damping_width ! ! 964 2012-07-26 09:14:24Z raasch ! check of old profil-parameters removed ! ! 940 2012-07-09 14:31:00Z raasch ! checks for parameter neutral ! ! 924 2012-06-06 07:44:41Z maronga ! Bugfix: preprocessor directives caused error during compilation ! ! 892 2012-05-02 13:51:44Z maronga ! Bugfix for parameter file check ( excluding __netcdf4 ) ! ! 866 2012-03-28 06:44:41Z raasch ! use only 60% of the geostrophic wind as translation speed in case of Galilean ! transformation and use_ug_for_galilei_tr = .T. in order to mimimize the ! timestep ! ! 861 2012-03-26 14:18:34Z suehring ! Check for topography and ws-scheme removed. ! Check for loop_optimization = 'vector' and ws-scheme removed. ! ! 845 2012-03-07 10:23:05Z maronga ! Bugfix: exclude __netcdf4 directive part from namelist file check compilation ! ! 828 2012-02-21 12:00:36Z raasch ! check of collision_kernel extended ! ! 825 2012-02-19 03:03:44Z raasch ! check for collision_kernel and curvature_solution_effects ! ! 809 2012-01-30 13:32:58Z maronga ! Bugfix: replaced .AND. and .NOT. with && and ! in the preprocessor directives ! ! 807 2012-01-25 11:53:51Z maronga ! New cpp directive "__check" implemented which is used by check_namelist_files ! ! 774 2011-10-27 13:34:16Z letzel ! bugfix for prescribed u,v-profiles ! ! 767 2011-10-14 06:39:12Z raasch ! Calculating u,v-profiles from given profiles by linear interpolation. ! bugfix: dirichlet_0 conditions for ug/vg moved from init_3d_model to here ! ! 707 2011-03-29 11:39:40Z raasch ! setting of bc_lr/ns_dirrad/raddir ! ! 689 2011-02-20 19:31:12z gryschka ! Bugfix for some logical expressions ! (syntax was not compatible with all compilers) ! ! 680 2011-02-04 23:16:06Z gryschka ! init_vortex is not allowed with volume_flow_control ! ! 673 2011-01-18 16:19:48Z suehring ! Declaration of ws_scheme_sca and ws_scheme_mom added (moved from advec_ws). ! ! 667 2010-12-23 12:06:00Z suehring/gryschka ! Exchange of parameters between ocean and atmosphere via PE0 ! Check for illegal combination of ws-scheme and timestep scheme. ! Check for topography and ws-scheme. ! Check for not cyclic boundary conditions in combination with ws-scheme and ! loop_optimization = 'vector'. ! Check for call_psolver_at_all_substeps and ws-scheme for momentum_advec. ! Different processor/grid topology in atmosphere and ocean is now allowed! ! Bugfixes in checking for conserve_volume_flow_mode ! 600 2010-11-24 16:10:51Z raasch ! change due to new default value of surface_waterflux ! 580 2010-10-05 13:59:11Z heinze ! renaming of ws_vertical_gradient_level to subs_vertical_gradient_level ! ! 567 2010-10-01 10:46:30Z helmke ! calculating masks changed ! ! 564 2010-09-30 13:18:59Z helmke ! palm message identifiers of masked output changed, 20 replaced by max_masks ! ! 553 2010-09-01 14:09:06Z weinreis ! masks is calculated and removed from inipar ! ! 531 2010-04-21 06:47:21Z heinze ! Bugfix: unit of hyp changed to dbar ! ! 524 2010-03-30 02:04:51Z raasch ! Bugfix: "/" in netcdf profile variable names replaced by ":" ! ! 493 2010-03-01 08:30:24Z raasch ! netcdf_data_format is checked ! ! 411 2009-12-11 14:15:58Z heinze ! Enabled passive scalar/humidity wall fluxes for non-flat topography ! Initialization of large scale vertical motion (subsidence/ascent) ! ! 410 2009-12-04 17:05:40Z letzel ! masked data output ! ! 388 2009-09-23 09:40:33Z raasch ! Check profiles fpr prho and hyp. ! Bugfix: output of averaged 2d/3d quantities requires that an avaraging ! interval has been set, respective error message is included ! bc_lr_cyc and bc_ns_cyc are set, ! initializing_actions='read_data_for_recycling' renamed to 'cyclic_fill' ! Check for illegal entries in section_xy|xz|yz that exceed nz+1|ny+1|nx+1 ! Coupling with independent precursor runs. ! Check particle_color, particle_dvrpsize, color_interval, dvrpsize_interval ! Bugfix: pressure included for profile output ! Check pressure gradient conditions ! topography_grid_convention moved from user_check_parameters ! 'single_street_canyon' ! Added shf* and qsws* to the list of available output data ! ! 222 2009-01-12 16:04:16Z letzel ! +user_check_parameters ! Output of messages replaced by message handling routine. ! Implementation of an MPI-1 coupling: replaced myid with target_id, ! deleted __mpi2 directives ! Check that PALM is called with mrun -K parallel for coupling ! ! 197 2008-09-16 15:29:03Z raasch ! Bug fix: Construction of vertical profiles when 10 gradients have been ! specified in the parameter list (ug, vg, pt, q, sa, lad) ! ! Strict grid matching along z is not needed for mg-solver. ! Leaf area density (LAD) explicitly set to its surface value at k=0 ! Case of reading data for recycling included in initializing_actions, ! check of turbulent_inflow and calculation of recycling_plane. ! q*2 profile added ! ! 138 2007-11-28 10:03:58Z letzel ! Plant canopy added ! Allow new case bc_uv_t = 'dirichlet_0' for channel flow. ! Multigrid solver allows topography, checking of dt_sort_particles ! Bugfix: initializing u_init and v_init in case of ocean runs ! ! 109 2007-08-28 15:26:47Z letzel ! Check coupling_mode and set default (obligatory) values (like boundary ! conditions for temperature and fluxes) in case of coupled runs. ! +profiles for w*p* and w"e ! Bugfix: Error message concerning output of particle concentration (pc) ! modified ! More checks and more default values for coupled runs ! allow data_output_pr= q, wq, w"q", w*q* for humidity = .T. (instead of ! cloud_physics = .T.) ! Rayleigh damping for ocean fixed. ! Check and, if necessary, set default value for dt_coupling ! ! 97 2007-06-21 08:23:15Z raasch ! Initial salinity profile is calculated, salinity boundary conditions are ! checked, ! z_max_do1d is checked only in case of ocean = .f., ! +initial temperature and geostrophic velocity profiles for the ocean version, ! use_pt_reference renamed use_reference ! ! 89 2007-05-25 12:08:31Z raasch ! Check for user-defined profiles ! ! 75 2007-03-22 09:54:05Z raasch ! "by_user" allowed as initializing action, -data_output_ts, ! leapfrog with non-flat topography not allowed any more, loop_optimization ! and pt_reference are checked, moisture renamed humidity, ! output of precipitation amount/rate and roughnes length + check ! possible negative humidities are avoided in initial profile, ! dirichlet/neumann changed to dirichlet/radiation, etc., ! revision added to run_description_header ! ! 20 2007-02-26 00:12:32Z raasch ! Temperature and humidity gradients at top are now calculated for nzt+1, ! top_heatflux and respective boundary condition bc_pt_t is checked ! ! RCS Log replace by Id keyword, revision history cleaned up ! ! Revision 1.61 2006/08/04 14:20:25 raasch ! do2d_unit and do3d_unit now defined as 2d-arrays, check of ! use_upstream_for_tke, default value for dt_dopts, ! generation of file header moved from routines palm and header to here ! ! Revision 1.1 1997/08/26 06:29:23 raasch ! Initial revision ! ! ! Description: ! ------------ ! Check control parameters and deduce further quantities. !------------------------------------------------------------------------------! USE arrays_3d USE cloud_parameters USE constants USE control_parameters USE dvrp_variables USE grid_variables USE indices USE model_1d USE netcdf_control USE particle_attributes USE pegrid USE profil_parameter USE subsidence_mod USE statistics USE transpose_indices IMPLICIT NONE CHARACTER (LEN=1) :: sq CHARACTER (LEN=6) :: var CHARACTER (LEN=7) :: unit CHARACTER (LEN=8) :: date CHARACTER (LEN=10) :: time CHARACTER (LEN=40) :: coupling_string CHARACTER (LEN=100) :: action INTEGER :: i, ilen, iremote = 0, j, k, kk, position, prec LOGICAL :: found, ldum REAL :: gradient, remote = 0.0, simulation_time_since_reference ! !-- Check for overlap combinations, which are not realized yet IF ( transpose_compute_overlap ) THEN IF ( numprocs == 1 ) STOP '+++ transpose-compute-overlap not implemented for single PE runs' #if defined( __openacc ) STOP '+++ transpose-compute-overlap not implemented for GPU usage' #endif ENDIF ! !-- Warning, if host is not set IF ( host(1:1) == ' ' ) THEN message_string = '"host" is not set. Please check that environment ' // & 'variable "localhost" & is set before running PALM' CALL message( 'check_parameters', 'PA0001', 0, 0, 0, 6, 0 ) ENDIF ! !-- Check the coupling mode IF ( coupling_mode /= 'uncoupled' .AND. & coupling_mode /= 'atmosphere_to_ocean' .AND. & coupling_mode /= 'ocean_to_atmosphere' ) THEN message_string = 'illegal coupling mode: ' // TRIM( coupling_mode ) CALL message( 'check_parameters', 'PA0002', 1, 2, 0, 6, 0 ) ENDIF ! !-- Check dt_coupling, restart_time, dt_restart, end_time, dx, dy, nx and ny IF ( coupling_mode /= 'uncoupled') THEN IF ( dt_coupling == 9999999.9 ) THEN message_string = 'dt_coupling is not set but required for coup' // & 'ling mode "' // TRIM( coupling_mode ) // '"' CALL message( 'check_parameters', 'PA0003', 1, 2, 0, 6, 0 ) ENDIF #if defined( __parallel ) ! !-- NOTE: coupled runs have not been implemented in the check_namelist_files !-- program. !-- check_namelist_files will need the following information of the other !-- model (atmosphere/ocean). ! dt_coupling = remote ! dt_max = remote ! restart_time = remote ! dt_restart= remote ! simulation_time_since_reference = remote ! dx = remote #if ! defined( __check ) IF ( myid == 0 ) THEN CALL MPI_SEND( dt_coupling, 1, MPI_REAL, target_id, 11, comm_inter, & ierr ) CALL MPI_RECV( remote, 1, MPI_REAL, target_id, 11, comm_inter, & status, ierr ) ENDIF CALL MPI_BCAST( remote, 1, MPI_REAL, 0, comm2d, ierr) #endif IF ( dt_coupling /= remote ) THEN WRITE( message_string, * ) 'coupling mode "', TRIM( coupling_mode ), & '": dt_coupling = ', dt_coupling, '& is not equal to ', & 'dt_coupling_remote = ', remote CALL message( 'check_parameters', 'PA0004', 1, 2, 0, 6, 0 ) ENDIF IF ( dt_coupling <= 0.0 ) THEN #if ! defined( __check ) IF ( myid == 0 ) THEN CALL MPI_SEND( dt_max, 1, MPI_REAL, target_id, 19, comm_inter, ierr ) CALL MPI_RECV( remote, 1, MPI_REAL, target_id, 19, comm_inter, & status, ierr ) ENDIF CALL MPI_BCAST( remote, 1, MPI_REAL, 0, comm2d, ierr) #endif dt_coupling = MAX( dt_max, remote ) WRITE( message_string, * ) 'coupling mode "', TRIM( coupling_mode ), & '": dt_coupling <= 0.0 & is not allowed and is reset to ', & 'MAX(dt_max(A,O)) = ', dt_coupling CALL message( 'check_parameters', 'PA0005', 0, 1, 0, 6, 0 ) ENDIF #if ! defined( __check ) IF ( myid == 0 ) THEN CALL MPI_SEND( restart_time, 1, MPI_REAL, target_id, 12, comm_inter, & ierr ) CALL MPI_RECV( remote, 1, MPI_REAL, target_id, 12, comm_inter, & status, ierr ) ENDIF CALL MPI_BCAST( remote, 1, MPI_REAL, 0, comm2d, ierr) #endif IF ( restart_time /= remote ) THEN WRITE( message_string, * ) 'coupling mode "', TRIM( coupling_mode ), & '": restart_time = ', restart_time, '& is not equal to ', & 'restart_time_remote = ', remote CALL message( 'check_parameters', 'PA0006', 1, 2, 0, 6, 0 ) ENDIF #if ! defined( __check ) IF ( myid == 0 ) THEN CALL MPI_SEND( dt_restart, 1, MPI_REAL, target_id, 13, comm_inter, & ierr ) CALL MPI_RECV( remote, 1, MPI_REAL, target_id, 13, comm_inter, & status, ierr ) ENDIF CALL MPI_BCAST( remote, 1, MPI_REAL, 0, comm2d, ierr) #endif IF ( dt_restart /= remote ) THEN WRITE( message_string, * ) 'coupling mode "', TRIM( coupling_mode ), & '": dt_restart = ', dt_restart, '& is not equal to ', & 'dt_restart_remote = ', remote CALL message( 'check_parameters', 'PA0007', 1, 2, 0, 6, 0 ) ENDIF simulation_time_since_reference = end_time - coupling_start_time #if ! defined( __check ) IF ( myid == 0 ) THEN CALL MPI_SEND( simulation_time_since_reference, 1, MPI_REAL, target_id, & 14, comm_inter, ierr ) CALL MPI_RECV( remote, 1, MPI_REAL, target_id, 14, comm_inter, & status, ierr ) ENDIF CALL MPI_BCAST( remote, 1, MPI_REAL, 0, comm2d, ierr) #endif IF ( simulation_time_since_reference /= remote ) THEN WRITE( message_string, * ) 'coupling mode "', TRIM( coupling_mode ), & '": simulation_time_since_reference = ', & simulation_time_since_reference, '& is not equal to ', & 'simulation_time_since_reference_remote = ', remote CALL message( 'check_parameters', 'PA0008', 1, 2, 0, 6, 0 ) ENDIF #if ! defined( __check ) IF ( myid == 0 ) THEN CALL MPI_SEND( dx, 1, MPI_REAL, target_id, 15, comm_inter, ierr ) CALL MPI_RECV( remote, 1, MPI_REAL, target_id, 15, comm_inter, & status, ierr ) ENDIF CALL MPI_BCAST( remote, 1, MPI_REAL, 0, comm2d, ierr) #endif IF ( coupling_mode == 'atmosphere_to_ocean') THEN IF ( dx < remote ) THEN WRITE( message_string, * ) 'coupling mode "', & TRIM( coupling_mode ), & '": dx in Atmosphere is not equal to or not larger then dx in ocean' CALL message( 'check_parameters', 'PA0009', 1, 2, 0, 6, 0 ) ENDIF IF ( (nx_a+1)*dx /= (nx_o+1)*remote ) THEN WRITE( message_string, * ) 'coupling mode "', & TRIM( coupling_mode ), & '": Domain size in x-direction is not equal in ocean and atmosphere' CALL message( 'check_parameters', 'PA0010', 1, 2, 0, 6, 0 ) ENDIF ENDIF #if ! defined( __check ) IF ( myid == 0) THEN CALL MPI_SEND( dy, 1, MPI_REAL, target_id, 16, comm_inter, ierr ) CALL MPI_RECV( remote, 1, MPI_REAL, target_id, 16, comm_inter, & status, ierr ) ENDIF CALL MPI_BCAST( remote, 1, MPI_REAL, 0, comm2d, ierr) #endif IF ( coupling_mode == 'atmosphere_to_ocean') THEN IF ( dy < remote ) THEN WRITE( message_string, * ) 'coupling mode "', & TRIM( coupling_mode ), & '": dy in Atmosphere is not equal to or not larger then dy in ocean' CALL message( 'check_parameters', 'PA0011', 1, 2, 0, 6, 0 ) ENDIF IF ( (ny_a+1)*dy /= (ny_o+1)*remote ) THEN WRITE( message_string, * ) 'coupling mode "', & TRIM( coupling_mode ), & '": Domain size in y-direction is not equal in ocean and atmosphere' CALL message( 'check_parameters', 'PA0012', 1, 2, 0, 6, 0 ) ENDIF IF ( MOD(nx_o+1,nx_a+1) /= 0 ) THEN WRITE( message_string, * ) 'coupling mode "', & TRIM( coupling_mode ), & '": nx+1 in ocean is not divisible without remainder with nx+1 in', & ' atmosphere' CALL message( 'check_parameters', 'PA0339', 1, 2, 0, 6, 0 ) ENDIF IF ( MOD(ny_o+1,ny_a+1) /= 0 ) THEN WRITE( message_string, * ) 'coupling mode "', & TRIM( coupling_mode ), & '": ny+1 in ocean is not divisible without remainder with ny+1 in', & ' atmosphere' CALL message( 'check_parameters', 'PA0340', 1, 2, 0, 6, 0 ) ENDIF ENDIF #else WRITE( message_string, * ) 'coupling requires PALM to be called with', & ' ''mrun -K parallel''' CALL message( 'check_parameters', 'PA0141', 1, 2, 0, 6, 0 ) #endif ENDIF #if defined( __parallel ) && ! defined ( __check ) ! !-- Exchange via intercommunicator IF ( coupling_mode == 'atmosphere_to_ocean' .AND. myid == 0 ) THEN CALL MPI_SEND( humidity, 1, MPI_LOGICAL, target_id, 19, comm_inter, & ierr ) ELSEIF ( coupling_mode == 'ocean_to_atmosphere' .AND. myid == 0) THEN CALL MPI_RECV( humidity_remote, 1, MPI_LOGICAL, target_id, 19, & comm_inter, status, ierr ) ENDIF CALL MPI_BCAST( humidity_remote, 1, MPI_LOGICAL, 0, comm2d, ierr) #endif ! !-- Generate the file header which is used as a header for most of PALM's !-- output files CALL DATE_AND_TIME( date, time ) run_date = date(7:8)//'-'//date(5:6)//'-'//date(3:4) run_time = time(1:2)//':'//time(3:4)//':'//time(5:6) IF ( coupling_mode == 'uncoupled' ) THEN coupling_string = '' ELSEIF ( coupling_mode == 'atmosphere_to_ocean' ) THEN coupling_string = ' coupled (atmosphere)' ELSEIF ( coupling_mode == 'ocean_to_atmosphere' ) THEN coupling_string = ' coupled (ocean)' ENDIF WRITE ( run_description_header, & '(A,2X,A,2X,A,A,A,I2.2,A,2X,A,A,2X,A,1X,A)' ) & TRIM( version ), TRIM( revision ), 'run: ', & TRIM( run_identifier ), '.', runnr, TRIM( coupling_string ), & 'host: ', TRIM( host ), run_date, run_time ! !-- Check the general loop optimization method IF ( loop_optimization == 'default' ) THEN IF ( host(1:3) == 'nec' ) THEN loop_optimization = 'vector' ELSE loop_optimization = 'cache' ENDIF ENDIF SELECT CASE ( TRIM( loop_optimization ) ) CASE ( 'acc', 'cache', 'vector' ) CONTINUE CASE DEFAULT message_string = 'illegal value given for loop_optimization: "' // & TRIM( loop_optimization ) // '"' CALL message( 'check_parameters', 'PA0013', 1, 2, 0, 6, 0 ) END SELECT ! !-- Check if vertical grid stretching is used together with particles IF ( dz_stretch_level < 100000.0 .AND. particle_advection ) THEN message_string = 'Vertical grid stretching is not allowed together ' // & 'with particle advection.' CALL message( 'check_parameters', 'PA0017', 1, 2, 0, 6, 0 ) ENDIF ! !-- Check topography setting (check for illegal parameter combinations) IF ( topography /= 'flat' ) THEN action = ' ' IF ( scalar_advec /= 'pw-scheme' .AND. scalar_advec /= 'ws-scheme') THEN WRITE( action, '(A,A)' ) 'scalar_advec = ', scalar_advec ENDIF IF ( momentum_advec /= 'pw-scheme' .AND. momentum_advec /= 'ws-scheme' ) & THEN WRITE( action, '(A,A)' ) 'momentum_advec = ', momentum_advec ENDIF IF ( psolver == 'sor' ) THEN WRITE( action, '(A,A)' ) 'psolver = ', psolver ENDIF IF ( sloping_surface ) THEN WRITE( action, '(A)' ) 'sloping surface = .TRUE.' ENDIF IF ( galilei_transformation ) THEN WRITE( action, '(A)' ) 'galilei_transformation = .TRUE.' ENDIF IF ( cloud_physics ) THEN WRITE( action, '(A)' ) 'cloud_physics = .TRUE.' ENDIF IF ( cloud_droplets ) THEN WRITE( action, '(A)' ) 'cloud_droplets = .TRUE.' ENDIF IF ( .NOT. prandtl_layer ) THEN WRITE( action, '(A)' ) 'prandtl_layer = .FALSE.' ENDIF IF ( action /= ' ' ) THEN message_string = 'a non-flat topography does not allow ' // & TRIM( action ) CALL message( 'check_parameters', 'PA0014', 1, 2, 0, 6, 0 ) ENDIF ! !-- In case of non-flat topography, check whether the convention how to !-- define the topography grid has been set correctly, or whether the default !-- is applicable. If this is not possible, abort. IF ( TRIM( topography_grid_convention ) == ' ' ) THEN IF ( TRIM( topography ) /= 'single_building' .AND. & TRIM( topography ) /= 'single_street_canyon' .AND. & TRIM( topography ) /= 'read_from_file' ) THEN !-- The default value is not applicable here, because it is only valid !-- for the two standard cases 'single_building' and 'read_from_file' !-- defined in init_grid. WRITE( message_string, * ) & 'The value for "topography_grid_convention" ', & 'is not set. Its default value is & only valid for ', & '"topography" = ''single_building'', ', & '''single_street_canyon'' & or ''read_from_file''.', & ' & Choose ''cell_edge'' or ''cell_center''.' CALL message( 'user_check_parameters', 'PA0239', 1, 2, 0, 6, 0 ) ELSE !-- The default value is applicable here. !-- Set convention according to topography. IF ( TRIM( topography ) == 'single_building' .OR. & TRIM( topography ) == 'single_street_canyon' ) THEN topography_grid_convention = 'cell_edge' ELSEIF ( TRIM( topography ) == 'read_from_file' ) THEN topography_grid_convention = 'cell_center' ENDIF ENDIF ELSEIF ( TRIM( topography_grid_convention ) /= 'cell_edge' .AND. & TRIM( topography_grid_convention ) /= 'cell_center' ) THEN WRITE( message_string, * ) & 'The value for "topography_grid_convention" is ', & 'not recognized. & Choose ''cell_edge'' or ''cell_center''.' CALL message( 'user_check_parameters', 'PA0240', 1, 2, 0, 6, 0 ) ENDIF ENDIF ! !-- Check ocean setting IF ( ocean ) THEN action = ' ' IF ( action /= ' ' ) THEN message_string = 'ocean = .T. does not allow ' // TRIM( action ) CALL message( 'check_parameters', 'PA0015', 1, 2, 0, 6, 0 ) ENDIF ELSEIF ( TRIM( coupling_mode ) == 'uncoupled' .AND. & TRIM( coupling_char ) == '_O' ) THEN ! !-- Check whether an (uncoupled) atmospheric run has been declared as an !-- ocean run (this setting is done via mrun-option -y) message_string = 'ocean = .F. does not allow coupling_char = "' // & TRIM( coupling_char ) // '" set by mrun-option "-y"' CALL message( 'check_parameters', 'PA0317', 1, 2, 0, 6, 0 ) ENDIF ! !-- Check cloud scheme IF ( cloud_scheme == 'seifert_beheng' ) THEN icloud_scheme = 0 ELSEIF ( cloud_scheme == 'kessler' ) THEN icloud_scheme = 1 ELSE message_string = 'unknown cloud microphysics scheme cloud_scheme ="' // & TRIM( cloud_scheme ) // '"' CALL message( 'check_parameters', 'PA0357', 1, 2, 0, 6, 0 ) ENDIF ! !-- Check whether there are any illegal values !-- Pressure solver: IF ( psolver /= 'poisfft' .AND. & psolver /= 'sor' .AND. psolver /= 'multigrid' ) THEN message_string = 'unknown solver for perturbation pressure: psolver' // & ' = "' // TRIM( psolver ) // '"' CALL message( 'check_parameters', 'PA0016', 1, 2, 0, 6, 0 ) ENDIF IF ( psolver == 'multigrid' ) THEN IF ( cycle_mg == 'w' ) THEN gamma_mg = 2 ELSEIF ( cycle_mg == 'v' ) THEN gamma_mg = 1 ELSE message_string = 'unknown multigrid cycle: cycle_mg = "' // & TRIM( cycle_mg ) // '"' CALL message( 'check_parameters', 'PA0020', 1, 2, 0, 6, 0 ) ENDIF ENDIF IF ( fft_method /= 'singleton-algorithm' .AND. & fft_method /= 'temperton-algorithm' .AND. & fft_method /= 'fftw' .AND. & fft_method /= 'system-specific' ) THEN message_string = 'unknown fft-algorithm: fft_method = "' // & TRIM( fft_method ) // '"' CALL message( 'check_parameters', 'PA0021', 1, 2, 0, 6, 0 ) ENDIF IF( momentum_advec == 'ws-scheme' .AND. & .NOT. call_psolver_at_all_substeps ) THEN message_string = 'psolver must be called at each RK3 substep when "'//& TRIM(momentum_advec) // ' "is used for momentum_advec' CALL message( 'check_parameters', 'PA0344', 1, 2, 0, 6, 0 ) END IF ! !-- Advection schemes: IF ( momentum_advec /= 'pw-scheme' .AND. momentum_advec /= 'ws-scheme' ) & THEN message_string = 'unknown advection scheme: momentum_advec = "' // & TRIM( momentum_advec ) // '"' CALL message( 'check_parameters', 'PA0022', 1, 2, 0, 6, 0 ) ENDIF IF ( ( momentum_advec == 'ws-scheme' .OR. scalar_advec == 'ws-scheme' ) & .AND. ( timestep_scheme == 'euler' .OR. & timestep_scheme == 'runge-kutta-2' ) ) & THEN message_string = 'momentum_advec or scalar_advec = "' & // TRIM( momentum_advec ) // '" is not allowed with timestep_scheme = "' // & TRIM( timestep_scheme ) // '"' CALL message( 'check_parameters', 'PA0023', 1, 2, 0, 6, 0 ) ENDIF IF ( scalar_advec /= 'pw-scheme' .AND. scalar_advec /= 'ws-scheme' .AND. & scalar_advec /= 'bc-scheme' ) & THEN message_string = 'unknown advection scheme: scalar_advec = "' // & TRIM( scalar_advec ) // '"' CALL message( 'check_parameters', 'PA0024', 1, 2, 0, 6, 0 ) ENDIF IF ( scalar_advec == 'bc-scheme' .AND. loop_optimization == 'cache' ) & THEN message_string = 'advection_scheme scalar_advec = "' & // TRIM( scalar_advec ) // '" not implemented for & loop_optimization = "' // & TRIM( loop_optimization ) // '"' CALL message( 'check_parameters', 'PA0026', 1, 2, 0, 6, 0 ) ENDIF IF ( use_sgs_for_particles .AND. .NOT. use_upstream_for_tke ) THEN use_upstream_for_tke = .TRUE. message_string = 'use_upstream_for_tke set .TRUE. because ' // & 'use_sgs_for_particles = .TRUE.' CALL message( 'check_parameters', 'PA0025', 0, 1, 0, 6, 0 ) ENDIF IF ( use_sgs_for_particles .AND. curvature_solution_effects ) THEN message_string = 'use_sgs_for_particles = .TRUE. not allowed with ' // & 'curvature_solution_effects = .TRUE.' CALL message( 'check_parameters', 'PA0349', 1, 2, 0, 6, 0 ) ENDIF ! !-- Set LOGICAL switches to enhance performance IF ( momentum_advec == 'ws-scheme' ) ws_scheme_mom = .TRUE. IF ( scalar_advec == 'ws-scheme' ) ws_scheme_sca = .TRUE. ! !-- Timestep schemes: SELECT CASE ( TRIM( timestep_scheme ) ) CASE ( 'euler' ) intermediate_timestep_count_max = 1 CASE ( 'runge-kutta-2' ) intermediate_timestep_count_max = 2 CASE ( 'runge-kutta-3' ) intermediate_timestep_count_max = 3 CASE DEFAULT message_string = 'unknown timestep scheme: timestep_scheme = "' // & TRIM( timestep_scheme ) // '"' CALL message( 'check_parameters', 'PA0027', 1, 2, 0, 6, 0 ) END SELECT IF ( (momentum_advec /= 'pw-scheme' .AND. momentum_advec /= 'ws-scheme') & .AND. timestep_scheme(1:5) == 'runge' ) THEN message_string = 'momentum advection scheme "' // & TRIM( momentum_advec ) // '" & does not work with ' // & 'timestep_scheme "' // TRIM( timestep_scheme ) // '"' CALL message( 'check_parameters', 'PA0029', 1, 2, 0, 6, 0 ) ENDIF ! !-- Collision kernels: SELECT CASE ( TRIM( collision_kernel ) ) CASE ( 'hall', 'hall_fast' ) hall_kernel = .TRUE. CASE ( 'palm' ) palm_kernel = .TRUE. CASE ( 'wang', 'wang_fast' ) wang_kernel = .TRUE. CASE ( 'none' ) CASE DEFAULT message_string = 'unknown collision kernel: collision_kernel = "' // & TRIM( collision_kernel ) // '"' CALL message( 'check_parameters', 'PA0350', 1, 2, 0, 6, 0 ) END SELECT IF ( collision_kernel(6:9) == 'fast' ) use_kernel_tables = .TRUE. IF ( TRIM( initializing_actions ) /= 'read_restart_data' .AND. & TRIM( initializing_actions ) /= 'cyclic_fill' ) THEN ! !-- No restart run: several initialising actions are possible action = initializing_actions DO WHILE ( TRIM( action ) /= '' ) position = INDEX( action, ' ' ) SELECT CASE ( action(1:position-1) ) CASE ( 'set_constant_profiles', 'set_1d-model_profiles', & 'by_user', 'initialize_vortex', 'initialize_ptanom' ) action = action(position+1:) CASE DEFAULT message_string = 'initializing_action = "' // & TRIM( action ) // '" unkown or not allowed' CALL message( 'check_parameters', 'PA0030', 1, 2, 0, 6, 0 ) END SELECT ENDDO ENDIF IF ( TRIM( initializing_actions ) == 'initialize_vortex' .AND. & conserve_volume_flow ) THEN message_string = 'initializing_actions = "initialize_vortex"' // & ' ist not allowed with conserve_volume_flow = .T.' CALL message( 'check_parameters', 'PA0343', 1, 2, 0, 6, 0 ) ENDIF IF ( INDEX( initializing_actions, 'set_constant_profiles' ) /= 0 .AND. & INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 ) THEN message_string = 'initializing_actions = "set_constant_profiles"' // & ' and "set_1d-model_profiles" are not allowed ' // & 'simultaneously' CALL message( 'check_parameters', 'PA0031', 1, 2, 0, 6, 0 ) ENDIF IF ( INDEX( initializing_actions, 'set_constant_profiles' ) /= 0 .AND. & INDEX( initializing_actions, 'by_user' ) /= 0 ) THEN message_string = 'initializing_actions = "set_constant_profiles"' // & ' and "by_user" are not allowed simultaneously' CALL message( 'check_parameters', 'PA0032', 1, 2, 0, 6, 0 ) ENDIF IF ( INDEX( initializing_actions, 'by_user' ) /= 0 .AND. & INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 ) THEN message_string = 'initializing_actions = "by_user" and ' // & '"set_1d-model_profiles" are not allowed simultaneously' CALL message( 'check_parameters', 'PA0033', 1, 2, 0, 6, 0 ) ENDIF IF ( cloud_physics .AND. .NOT. humidity ) THEN WRITE( message_string, * ) 'cloud_physics = ', cloud_physics, ' is ', & 'not allowed with humidity = ', humidity CALL message( 'check_parameters', 'PA0034', 1, 2, 0, 6, 0 ) ENDIF IF ( precipitation .AND. .NOT. cloud_physics ) THEN WRITE( message_string, * ) 'precipitation = ', precipitation, ' is ', & 'not allowed with cloud_physics = ', cloud_physics CALL message( 'check_parameters', 'PA0035', 1, 2, 0, 6, 0 ) ENDIF IF ( humidity .AND. sloping_surface ) THEN message_string = 'humidity = .TRUE. and sloping_surface = .TRUE. ' // & 'are not allowed simultaneously' CALL message( 'check_parameters', 'PA0036', 1, 2, 0, 6, 0 ) ENDIF IF ( passive_scalar .AND. humidity ) THEN message_string = 'humidity = .TRUE. and passive_scalar = .TRUE. ' // & 'is not allowed simultaneously' CALL message( 'check_parameters', 'PA0038', 1, 2, 0, 6, 0 ) ENDIF IF ( plant_canopy .AND. ( drag_coefficient == 0.0 ) ) THEN message_string = 'plant_canopy = .TRUE. requires a non-zero drag ' // & 'coefficient & given value is drag_coefficient = 0.0' CALL message( 'check_parameters', 'PA0041', 1, 2, 0, 6, 0 ) ENDIF IF ( plant_canopy .AND. cloud_physics .AND. icloud_scheme == 0 ) THEN message_string = 'plant_canopy = .TRUE. requires cloud_scheme /=' // & ' seifert_beheng' CALL message( 'check_parameters', 'PA0360', 1, 2, 0, 6, 0 ) ENDIF IF ( loop_optimization /= 'cache' .AND. cloud_physics .AND. & icloud_scheme == 0 ) THEN message_string = 'cloud_scheme = seifert_beheng requires ' // & 'loop_optimization = cache' CALL message( 'check_parameters', 'PA0362', 1, 2, 0, 6, 0 ) ENDIF ! IF ( cloud_physics .AND. icloud_scheme == 0 .AND. & ! .NOT. precipitation .AND. .NOT. drizzle ) THEN ! message_string = 'cloud_scheme = seifert_beheng requires ' // & ! 'precipitation = .TRUE. or drizzle = .TRUE.' ! CALL message( 'check_parameters', 'PA0363', 1, 2, 0, 6, 0 ) ! ENDIF ! !-- In case of no model continuation run, check initialising parameters and !-- deduce further quantities IF ( TRIM( initializing_actions ) /= 'read_restart_data' ) THEN ! !-- Initial profiles for 1D and 3D model, respectively (u,v further below) pt_init = pt_surface IF ( humidity ) THEN q_init = q_surface ENDIF IF ( ocean ) sa_init = sa_surface IF ( passive_scalar ) q_init = s_surface IF ( plant_canopy ) lad = 0.0 ! !-- !-- If required, compute initial profile of the geostrophic wind !-- (component ug) i = 1 gradient = 0.0 IF ( .NOT. ocean ) THEN ug_vertical_gradient_level_ind(1) = 0 ug(0) = ug_surface DO k = 1, nzt+1 IF ( i < 11 ) THEN IF ( ug_vertical_gradient_level(i) < zu(k) .AND. & ug_vertical_gradient_level(i) >= 0.0 ) THEN gradient = ug_vertical_gradient(i) / 100.0 ug_vertical_gradient_level_ind(i) = k - 1 i = i + 1 ENDIF ENDIF IF ( gradient /= 0.0 ) THEN IF ( k /= 1 ) THEN ug(k) = ug(k-1) + dzu(k) * gradient ELSE ug(k) = ug_surface + dzu(k) * gradient ENDIF ELSE ug(k) = ug(k-1) ENDIF ENDDO ELSE ug_vertical_gradient_level_ind(1) = nzt+1 ug(nzt+1) = ug_surface DO k = nzt, nzb, -1 IF ( i < 11 ) THEN IF ( ug_vertical_gradient_level(i) > zu(k) .AND. & ug_vertical_gradient_level(i) <= 0.0 ) THEN gradient = ug_vertical_gradient(i) / 100.0 ug_vertical_gradient_level_ind(i) = k + 1 i = i + 1 ENDIF ENDIF IF ( gradient /= 0.0 ) THEN IF ( k /= nzt ) THEN ug(k) = ug(k+1) - dzu(k+1) * gradient ELSE ug(k) = ug_surface - 0.5 * dzu(k+1) * gradient ug(k+1) = ug_surface + 0.5 * dzu(k+1) * gradient ENDIF ELSE ug(k) = ug(k+1) ENDIF ENDDO ENDIF ! !-- In case of no given gradients for ug, choose a zero gradient IF ( ug_vertical_gradient_level(1) == -9999999.9 ) THEN ug_vertical_gradient_level(1) = 0.0 ENDIF ! !-- !-- If required, compute initial profile of the geostrophic wind !-- (component vg) i = 1 gradient = 0.0 IF ( .NOT. ocean ) THEN vg_vertical_gradient_level_ind(1) = 0 vg(0) = vg_surface DO k = 1, nzt+1 IF ( i < 11 ) THEN IF ( vg_vertical_gradient_level(i) < zu(k) .AND. & vg_vertical_gradient_level(i) >= 0.0 ) THEN gradient = vg_vertical_gradient(i) / 100.0 vg_vertical_gradient_level_ind(i) = k - 1 i = i + 1 ENDIF ENDIF IF ( gradient /= 0.0 ) THEN IF ( k /= 1 ) THEN vg(k) = vg(k-1) + dzu(k) * gradient ELSE vg(k) = vg_surface + dzu(k) * gradient ENDIF ELSE vg(k) = vg(k-1) ENDIF ENDDO ELSE vg_vertical_gradient_level_ind(1) = nzt+1 vg(nzt+1) = vg_surface DO k = nzt, nzb, -1 IF ( i < 11 ) THEN IF ( vg_vertical_gradient_level(i) > zu(k) .AND. & vg_vertical_gradient_level(i) <= 0.0 ) THEN gradient = vg_vertical_gradient(i) / 100.0 vg_vertical_gradient_level_ind(i) = k + 1 i = i + 1 ENDIF ENDIF IF ( gradient /= 0.0 ) THEN IF ( k /= nzt ) THEN vg(k) = vg(k+1) - dzu(k+1) * gradient ELSE vg(k) = vg_surface - 0.5 * dzu(k+1) * gradient vg(k+1) = vg_surface + 0.5 * dzu(k+1) * gradient ENDIF ELSE vg(k) = vg(k+1) ENDIF ENDDO ENDIF ! !-- In case of no given gradients for vg, choose a zero gradient IF ( vg_vertical_gradient_level(1) == -9999999.9 ) THEN vg_vertical_gradient_level(1) = 0.0 ENDIF ! !-- Let the initial wind profiles be the calculated ug/vg profiles or !-- interpolate them from wind profile data (if given) IF ( u_profile(1) == 9999999.9 .AND. v_profile(1) == 9999999.9 ) THEN u_init = ug v_init = vg ELSEIF ( u_profile(1) == 0.0 .AND. v_profile(1) == 0.0 ) THEN IF ( uv_heights(1) /= 0.0 ) THEN message_string = 'uv_heights(1) must be 0.0' CALL message( 'check_parameters', 'PA0345', 1, 2, 0, 6, 0 ) ENDIF use_prescribed_profile_data = .TRUE. kk = 1 u_init(0) = 0.0 v_init(0) = 0.0 DO k = 1, nz+1 IF ( kk < 100 ) THEN DO WHILE ( uv_heights(kk+1) <= zu(k) ) kk = kk + 1 IF ( kk == 100 ) EXIT ENDDO ENDIF IF ( kk < 100 .AND. uv_heights(kk+1) /= 9999999.9 ) THEN u_init(k) = u_profile(kk) + ( zu(k) - uv_heights(kk) ) / & ( uv_heights(kk+1) - uv_heights(kk) ) * & ( u_profile(kk+1) - u_profile(kk) ) v_init(k) = v_profile(kk) + ( zu(k) - uv_heights(kk) ) / & ( uv_heights(kk+1) - uv_heights(kk) ) * & ( v_profile(kk+1) - v_profile(kk) ) ELSE u_init(k) = u_profile(kk) v_init(k) = v_profile(kk) ENDIF ENDDO ELSE message_string = 'u_profile(1) and v_profile(1) must be 0.0' CALL message( 'check_parameters', 'PA0346', 1, 2, 0, 6, 0 ) ENDIF ! !-- Compute initial temperature profile using the given temperature gradients IF ( .NOT. neutral ) THEN i = 1 gradient = 0.0 IF ( .NOT. ocean ) THEN pt_vertical_gradient_level_ind(1) = 0 DO k = 1, nzt+1 IF ( i < 11 ) THEN IF ( pt_vertical_gradient_level(i) < zu(k) .AND. & pt_vertical_gradient_level(i) >= 0.0 ) THEN gradient = pt_vertical_gradient(i) / 100.0 pt_vertical_gradient_level_ind(i) = k - 1 i = i + 1 ENDIF ENDIF IF ( gradient /= 0.0 ) THEN IF ( k /= 1 ) THEN pt_init(k) = pt_init(k-1) + dzu(k) * gradient ELSE pt_init(k) = pt_surface + dzu(k) * gradient ENDIF ELSE pt_init(k) = pt_init(k-1) ENDIF ENDDO ELSE pt_vertical_gradient_level_ind(1) = nzt+1 DO k = nzt, 0, -1 IF ( i < 11 ) THEN IF ( pt_vertical_gradient_level(i) > zu(k) .AND. & pt_vertical_gradient_level(i) <= 0.0 ) THEN gradient = pt_vertical_gradient(i) / 100.0 pt_vertical_gradient_level_ind(i) = k + 1 i = i + 1 ENDIF ENDIF IF ( gradient /= 0.0 ) THEN IF ( k /= nzt ) THEN pt_init(k) = pt_init(k+1) - dzu(k+1) * gradient ELSE pt_init(k) = pt_surface - 0.5 * dzu(k+1) * gradient pt_init(k+1) = pt_surface + 0.5 * dzu(k+1) * gradient ENDIF ELSE pt_init(k) = pt_init(k+1) ENDIF ENDDO ENDIF ENDIF ! !-- In case of no given temperature gradients, choose gradient of neutral !-- stratification IF ( pt_vertical_gradient_level(1) == -9999999.9 ) THEN pt_vertical_gradient_level(1) = 0.0 ENDIF ! !-- Store temperature gradient at the top boundary for possible Neumann !-- boundary condition bc_pt_t_val = ( pt_init(nzt+1) - pt_init(nzt) ) / dzu(nzt+1) ! !-- If required, compute initial humidity or scalar profile using the given !-- humidity/scalar gradient. In case of scalar transport, initially store !-- values of the scalar parameters on humidity parameters IF ( passive_scalar ) THEN bc_q_b = bc_s_b bc_q_t = bc_s_t q_surface = s_surface q_surface_initial_change = s_surface_initial_change q_vertical_gradient = s_vertical_gradient q_vertical_gradient_level = s_vertical_gradient_level surface_waterflux = surface_scalarflux wall_humidityflux = wall_scalarflux ENDIF IF ( humidity .OR. passive_scalar ) THEN i = 1 gradient = 0.0 q_vertical_gradient_level_ind(1) = 0 DO k = 1, nzt+1 IF ( i < 11 ) THEN IF ( q_vertical_gradient_level(i) < zu(k) .AND. & q_vertical_gradient_level(i) >= 0.0 ) THEN gradient = q_vertical_gradient(i) / 100.0 q_vertical_gradient_level_ind(i) = k - 1 i = i + 1 ENDIF ENDIF IF ( gradient /= 0.0 ) THEN IF ( k /= 1 ) THEN q_init(k) = q_init(k-1) + dzu(k) * gradient ELSE q_init(k) = q_init(k-1) + dzu(k) * gradient ENDIF ELSE q_init(k) = q_init(k-1) ENDIF ! !-- Avoid negative humidities IF ( q_init(k) < 0.0 ) THEN q_init(k) = 0.0 ENDIF ENDDO ! !-- In case of no given humidity gradients, choose zero gradient !-- conditions IF ( q_vertical_gradient_level(1) == -1.0 ) THEN q_vertical_gradient_level(1) = 0.0 ENDIF ! !-- Store humidity, rain water content and rain drop concentration !-- gradient at the top boundary for possile Neumann boundary condition bc_q_t_val = ( q_init(nzt+1) - q_init(nzt) ) / dzu(nzt+1) ENDIF ! !-- If required, compute initial salinity profile using the given salinity !-- gradients IF ( ocean ) THEN i = 1 gradient = 0.0 sa_vertical_gradient_level_ind(1) = nzt+1 DO k = nzt, 0, -1 IF ( i < 11 ) THEN IF ( sa_vertical_gradient_level(i) > zu(k) .AND. & sa_vertical_gradient_level(i) <= 0.0 ) THEN gradient = sa_vertical_gradient(i) / 100.0 sa_vertical_gradient_level_ind(i) = k + 1 i = i + 1 ENDIF ENDIF IF ( gradient /= 0.0 ) THEN IF ( k /= nzt ) THEN sa_init(k) = sa_init(k+1) - dzu(k+1) * gradient ELSE sa_init(k) = sa_surface - 0.5 * dzu(k+1) * gradient sa_init(k+1) = sa_surface + 0.5 * dzu(k+1) * gradient ENDIF ELSE sa_init(k) = sa_init(k+1) ENDIF ENDDO ENDIF ! !-- If required compute the profile of leaf area density used in the plant !-- canopy model IF ( plant_canopy ) THEN i = 1 gradient = 0.0 IF ( .NOT. ocean ) THEN lad(0) = lad_surface lad_vertical_gradient_level_ind(1) = 0 DO k = 1, pch_index IF ( i < 11 ) THEN IF ( lad_vertical_gradient_level(i) < zu(k) .AND. & lad_vertical_gradient_level(i) >= 0.0 ) THEN gradient = lad_vertical_gradient(i) lad_vertical_gradient_level_ind(i) = k - 1 i = i + 1 ENDIF ENDIF IF ( gradient /= 0.0 ) THEN IF ( k /= 1 ) THEN lad(k) = lad(k-1) + dzu(k) * gradient ELSE lad(k) = lad_surface + dzu(k) *gradient ENDIF ELSE lad(k) = lad(k-1) ENDIF ENDDO ENDIF ! !-- In case of no given leaf area density gradients, choose a vanishing !-- gradient IF ( lad_vertical_gradient_level(1) == -9999999.9 ) THEN lad_vertical_gradient_level(1) = 0.0 ENDIF ENDIF ENDIF ! !-- Initialize large scale subsidence if required IF ( subs_vertical_gradient_level(1) /= -9999999.9 ) THEN large_scale_subsidence = .TRUE. CALL init_w_subsidence ENDIF ! !-- Compute Coriolis parameter f = 2.0 * omega * SIN( phi / 180.0 * pi ) fs = 2.0 * omega * COS( phi / 180.0 * pi ) ! !-- Check and set buoyancy related parameters and switches IF ( reference_state == 'horizontal_average' ) THEN CONTINUE ELSEIF ( reference_state == 'initial_profile' ) THEN use_initial_profile_as_reference = .TRUE. ELSEIF ( reference_state == 'single_value' ) THEN use_single_reference_value = .TRUE. IF ( pt_reference == 9999999.9 ) pt_reference = pt_surface vpt_reference = pt_reference * ( 1.0 + 0.61 * q_surface ) ELSE message_string = 'illegal value for reference_state: "' // & TRIM( reference_state ) // '"' CALL message( 'check_parameters', 'PA0056', 1, 2, 0, 6, 0 ) ENDIF ! !-- Ocean runs always use reference values in the buoyancy term IF ( ocean ) THEN reference_state = 'single_value' use_single_reference_value = .TRUE. ENDIF ! !-- Sign of buoyancy/stability terms IF ( ocean ) atmos_ocean_sign = -1.0 ! !-- Ocean version must use flux boundary conditions at the top IF ( ocean .AND. .NOT. use_top_fluxes ) THEN message_string = 'use_top_fluxes must be .TRUE. in ocean version' CALL message( 'check_parameters', 'PA0042', 1, 2, 0, 6, 0 ) ENDIF ! !-- In case of a given slope, compute the relevant quantities IF ( alpha_surface /= 0.0 ) THEN IF ( ABS( alpha_surface ) > 90.0 ) THEN WRITE( message_string, * ) 'ABS( alpha_surface = ', alpha_surface, & ' ) must be < 90.0' CALL message( 'check_parameters', 'PA0043', 1, 2, 0, 6, 0 ) ENDIF sloping_surface = .TRUE. cos_alpha_surface = COS( alpha_surface / 180.0 * pi ) sin_alpha_surface = SIN( alpha_surface / 180.0 * pi ) ENDIF ! !-- Check time step and cfl_factor IF ( dt /= -1.0 ) THEN IF ( dt <= 0.0 .AND. dt /= -1.0 ) THEN WRITE( message_string, * ) 'dt = ', dt , ' <= 0.0' CALL message( 'check_parameters', 'PA0044', 1, 2, 0, 6, 0 ) ENDIF dt_3d = dt dt_fixed = .TRUE. ENDIF IF ( cfl_factor <= 0.0 .OR. cfl_factor > 1.0 ) THEN IF ( cfl_factor == -1.0 ) THEN IF ( timestep_scheme == 'runge-kutta-2' ) THEN cfl_factor = 0.8 ELSEIF ( timestep_scheme == 'runge-kutta-3' ) THEN cfl_factor = 0.9 ELSE cfl_factor = 0.9 ENDIF ELSE WRITE( message_string, * ) 'cfl_factor = ', cfl_factor, & ' out of range & 0.0 < cfl_factor <= 1.0 is required' CALL message( 'check_parameters', 'PA0045', 1, 2, 0, 6, 0 ) ENDIF ENDIF ! !-- Store simulated time at begin simulated_time_at_begin = simulated_time ! !-- Store reference time for coupled runs and change the coupling flag, !-- if ... IF ( simulated_time == 0.0 ) THEN IF ( coupling_start_time == 0.0 ) THEN time_since_reference_point = 0.0 ELSEIF ( time_since_reference_point < 0.0 ) THEN run_coupled = .FALSE. ENDIF ENDIF ! !-- Set wind speed in the Galilei-transformed system IF ( galilei_transformation ) THEN IF ( use_ug_for_galilei_tr .AND. & ug_vertical_gradient_level(1) == 0.0 .AND. & ug_vertical_gradient(1) == 0.0 .AND. & vg_vertical_gradient_level(1) == 0.0 .AND. & vg_vertical_gradient(1) == 0.0 ) THEN u_gtrans = ug_surface * 0.6 v_gtrans = vg_surface * 0.6 ELSEIF ( use_ug_for_galilei_tr .AND. & ( ug_vertical_gradient_level(1) /= 0.0 .OR. & ug_vertical_gradient(1) /= 0.0 ) ) THEN message_string = 'baroclinicity (ug) not allowed simultaneously' // & ' with galilei transformation' CALL message( 'check_parameters', 'PA0046', 1, 2, 0, 6, 0 ) ELSEIF ( use_ug_for_galilei_tr .AND. & ( vg_vertical_gradient_level(1) /= 0.0 .OR. & vg_vertical_gradient(1) /= 0.0 ) ) THEN message_string = 'baroclinicity (vg) not allowed simultaneously' // & ' with galilei transformation' CALL message( 'check_parameters', 'PA0047', 1, 2, 0, 6, 0 ) ELSE message_string = 'variable translation speed used for galilei-' // & 'transformation, which may cause & instabilities in stably ' // & 'stratified regions' CALL message( 'check_parameters', 'PA0048', 0, 1, 0, 6, 0 ) ENDIF ENDIF ! !-- In case of using a prandtl-layer, calculated (or prescribed) surface !-- fluxes have to be used in the diffusion-terms IF ( prandtl_layer ) use_surface_fluxes = .TRUE. ! !-- Check boundary conditions and set internal variables: !-- Lateral boundary conditions IF ( bc_lr /= 'cyclic' .AND. bc_lr /= 'dirichlet/radiation' .AND. & bc_lr /= 'radiation/dirichlet' ) THEN message_string = 'unknown boundary condition: bc_lr = "' // & TRIM( bc_lr ) // '"' CALL message( 'check_parameters', 'PA0049', 1, 2, 0, 6, 0 ) ENDIF IF ( bc_ns /= 'cyclic' .AND. bc_ns /= 'dirichlet/radiation' .AND. & bc_ns /= 'radiation/dirichlet' ) THEN message_string = 'unknown boundary condition: bc_ns = "' // & TRIM( bc_ns ) // '"' CALL message( 'check_parameters', 'PA0050', 1, 2, 0, 6, 0 ) ENDIF ! !-- Internal variables used for speed optimization in if clauses IF ( bc_lr /= 'cyclic' ) bc_lr_cyc = .FALSE. IF ( bc_lr == 'dirichlet/radiation' ) bc_lr_dirrad = .TRUE. IF ( bc_lr == 'radiation/dirichlet' ) bc_lr_raddir = .TRUE. IF ( bc_ns /= 'cyclic' ) bc_ns_cyc = .FALSE. IF ( bc_ns == 'dirichlet/radiation' ) bc_ns_dirrad = .TRUE. IF ( bc_ns == 'radiation/dirichlet' ) bc_ns_raddir = .TRUE. ! !-- Non-cyclic lateral boundaries require the multigrid method and Piascek- !-- Willimas or Wicker - Skamarock advection scheme. Several schemes !-- and tools do not work with non-cyclic boundary conditions. IF ( bc_lr /= 'cyclic' .OR. bc_ns /= 'cyclic' ) THEN IF ( psolver /= 'multigrid' ) THEN message_string = 'non-cyclic lateral boundaries do not allow ' // & 'psolver = "' // TRIM( psolver ) // '"' CALL message( 'check_parameters', 'PA0051', 1, 2, 0, 6, 0 ) ENDIF IF ( momentum_advec /= 'pw-scheme' .AND. & momentum_advec /= 'ws-scheme') THEN message_string = 'non-cyclic lateral boundaries do not allow ' // & 'momentum_advec = "' // TRIM( momentum_advec ) // '"' CALL message( 'check_parameters', 'PA0052', 1, 2, 0, 6, 0 ) ENDIF IF ( scalar_advec /= 'pw-scheme' .AND. & scalar_advec /= 'ws-scheme' ) THEN message_string = 'non-cyclic lateral boundaries do not allow ' // & 'scalar_advec = "' // TRIM( scalar_advec ) // '"' CALL message( 'check_parameters', 'PA0053', 1, 2, 0, 6, 0 ) ENDIF IF ( galilei_transformation ) THEN message_string = 'non-cyclic lateral boundaries do not allow ' // & 'galilei_transformation = .T.' CALL message( 'check_parameters', 'PA0054', 1, 2, 0, 6, 0 ) ENDIF ENDIF ! !-- Bottom boundary condition for the turbulent Kinetic energy IF ( bc_e_b == 'neumann' ) THEN ibc_e_b = 1 ELSEIF ( bc_e_b == '(u*)**2+neumann' ) THEN ibc_e_b = 2 IF ( prandtl_layer ) THEN message_string = 'adjust mixing length = FALSE and bc_e_b = "' // & TRIM( bc_e_b ) // '"' CALL message( 'check_parameters', 'PA0056', 0, 1, 0, 6, 0 ) ENDIF IF ( .NOT. prandtl_layer ) THEN bc_e_b = 'neumann' ibc_e_b = 1 message_string = 'boundary condition bc_e_b changed to "' // & TRIM( bc_e_b ) // '"' CALL message( 'check_parameters', 'PA0057', 0, 1, 0, 6, 0 ) ENDIF ELSE message_string = 'unknown boundary condition: bc_e_b = "' // & TRIM( bc_e_b ) // '"' CALL message( 'check_parameters', 'PA0058', 1, 2, 0, 6, 0 ) ENDIF ! !-- Boundary conditions for perturbation pressure IF ( bc_p_b == 'dirichlet' ) THEN ibc_p_b = 0 ELSEIF ( bc_p_b == 'neumann' ) THEN ibc_p_b = 1 ELSE message_string = 'unknown boundary condition: bc_p_b = "' // & TRIM( bc_p_b ) // '"' CALL message( 'check_parameters', 'PA0059', 1, 2, 0, 6, 0 ) ENDIF IF ( bc_p_t == 'dirichlet' ) THEN ibc_p_t = 0 ELSEIF ( bc_p_t == 'neumann' ) THEN ibc_p_t = 1 ELSE message_string = 'unknown boundary condition: bc_p_t = "' // & TRIM( bc_p_t ) // '"' CALL message( 'check_parameters', 'PA0061', 1, 2, 0, 6, 0 ) ENDIF ! !-- Boundary conditions for potential temperature IF ( coupling_mode == 'atmosphere_to_ocean' ) THEN ibc_pt_b = 2 ELSE IF ( bc_pt_b == 'dirichlet' ) THEN ibc_pt_b = 0 ELSEIF ( bc_pt_b == 'neumann' ) THEN ibc_pt_b = 1 ELSE message_string = 'unknown boundary condition: bc_pt_b = "' // & TRIM( bc_pt_b ) // '"' CALL message( 'check_parameters', 'PA0062', 1, 2, 0, 6, 0 ) ENDIF ENDIF IF ( bc_pt_t == 'dirichlet' ) THEN ibc_pt_t = 0 ELSEIF ( bc_pt_t == 'neumann' ) THEN ibc_pt_t = 1 ELSEIF ( bc_pt_t == 'initial_gradient' ) THEN ibc_pt_t = 2 ELSE message_string = 'unknown boundary condition: bc_pt_t = "' // & TRIM( bc_pt_t ) // '"' CALL message( 'check_parameters', 'PA0063', 1, 2, 0, 6, 0 ) ENDIF IF ( surface_heatflux == 9999999.9 ) constant_heatflux = .FALSE. IF ( top_heatflux == 9999999.9 ) constant_top_heatflux = .FALSE. IF ( neutral ) THEN IF ( surface_heatflux /= 0.0 .AND. surface_heatflux /= 9999999.9 ) & THEN message_string = 'heatflux must not be set for pure neutral flow' CALL message( 'check_parameters', 'PA0351', 1, 2, 0, 6, 0 ) ENDIF IF ( top_heatflux /= 0.0 .AND. top_heatflux /= 9999999.9 ) & THEN message_string = 'heatflux must not be set for pure neutral flow' CALL message( 'check_parameters', 'PA0351', 1, 2, 0, 6, 0 ) ENDIF ENDIF IF ( top_momentumflux_u /= 9999999.9 .AND. & top_momentumflux_v /= 9999999.9 ) THEN constant_top_momentumflux = .TRUE. ELSEIF ( .NOT. ( top_momentumflux_u == 9999999.9 .AND. & top_momentumflux_v == 9999999.9 ) ) THEN message_string = 'both, top_momentumflux_u AND top_momentumflux_v ' // & 'must be set' CALL message( 'check_parameters', 'PA0064', 1, 2, 0, 6, 0 ) ENDIF ! !-- A given surface temperature implies Dirichlet boundary condition for !-- temperature. In this case specification of a constant heat flux is !-- forbidden. IF ( ibc_pt_b == 0 .AND. constant_heatflux .AND. & surface_heatflux /= 0.0 ) THEN message_string = 'boundary_condition: bc_pt_b = "' // TRIM( bc_pt_b ) //& '& is not allowed with constant_heatflux = .TRUE.' CALL message( 'check_parameters', 'PA0065', 1, 2, 0, 6, 0 ) ENDIF IF ( constant_heatflux .AND. pt_surface_initial_change /= 0.0 ) THEN WRITE ( message_string, * ) 'constant_heatflux = .TRUE. is not allo', & 'wed with pt_surface_initial_change (/=0) = ', & pt_surface_initial_change CALL message( 'check_parameters', 'PA0066', 1, 2, 0, 6, 0 ) ENDIF ! !-- A given temperature at the top implies Dirichlet boundary condition for !-- temperature. In this case specification of a constant heat flux is !-- forbidden. IF ( ibc_pt_t == 0 .AND. constant_top_heatflux .AND. & top_heatflux /= 0.0 ) THEN message_string = 'boundary_condition: bc_pt_t = "' // TRIM( bc_pt_t ) //& '" is not allowed with constant_top_heatflux = .TRUE.' CALL message( 'check_parameters', 'PA0067', 1, 2, 0, 6, 0 ) ENDIF ! !-- Boundary conditions for salinity IF ( ocean ) THEN IF ( bc_sa_t == 'dirichlet' ) THEN ibc_sa_t = 0 ELSEIF ( bc_sa_t == 'neumann' ) THEN ibc_sa_t = 1 ELSE message_string = 'unknown boundary condition: bc_sa_t = "' // & TRIM( bc_sa_t ) // '"' CALL message( 'check_parameters', 'PA0068', 1, 2, 0, 6, 0 ) ENDIF IF ( top_salinityflux == 9999999.9 ) constant_top_salinityflux = .FALSE. IF ( ibc_sa_t == 1 .AND. top_salinityflux == 9999999.9 ) THEN message_string = 'boundary condition: bc_sa_t = "' // & TRIM( bc_sa_t ) // '" requires to set ' // & 'top_salinityflux' CALL message( 'check_parameters', 'PA0069', 1, 2, 0, 6, 0 ) ENDIF ! !-- A fixed salinity at the top implies Dirichlet boundary condition for !-- salinity. In this case specification of a constant salinity flux is !-- forbidden. IF ( ibc_sa_t == 0 .AND. constant_top_salinityflux .AND. & top_salinityflux /= 0.0 ) THEN message_string = 'boundary condition: bc_sa_t = "' // & TRIM( bc_sa_t ) // '" is not allowed with ' // & 'constant_top_salinityflux = .TRUE.' CALL message( 'check_parameters', 'PA0070', 1, 2, 0, 6, 0 ) ENDIF ENDIF ! !-- In case of humidity or passive scalar, set boundary conditions for total !-- water content / scalar IF ( humidity .OR. passive_scalar ) THEN IF ( humidity ) THEN sq = 'q' ELSE sq = 's' ENDIF IF ( bc_q_b == 'dirichlet' ) THEN ibc_q_b = 0 ELSEIF ( bc_q_b == 'neumann' ) THEN ibc_q_b = 1 ELSE message_string = 'unknown boundary condition: bc_' // TRIM( sq ) // & '_b ="' // TRIM( bc_q_b ) // '"' CALL message( 'check_parameters', 'PA0071', 1, 2, 0, 6, 0 ) ENDIF IF ( bc_q_t == 'dirichlet' ) THEN ibc_q_t = 0 ELSEIF ( bc_q_t == 'neumann' ) THEN ibc_q_t = 1 ELSE message_string = 'unknown boundary condition: bc_' // TRIM( sq ) // & '_t ="' // TRIM( bc_q_t ) // '"' CALL message( 'check_parameters', 'PA0072', 1, 2, 0, 6, 0 ) ENDIF IF ( surface_waterflux == 9999999.9 ) constant_waterflux = .FALSE. ! !-- A given surface humidity implies Dirichlet boundary condition for !-- humidity. In this case specification of a constant water flux is !-- forbidden. IF ( ibc_q_b == 0 .AND. constant_waterflux ) THEN message_string = 'boundary condition: bc_' // TRIM( sq ) // '_b ' // & '= "' // TRIM( bc_q_b ) // '" is not allowed wi' // & 'th prescribed surface flux' CALL message( 'check_parameters', 'PA0073', 1, 2, 0, 6, 0 ) ENDIF IF ( constant_waterflux .AND. q_surface_initial_change /= 0.0 ) THEN WRITE( message_string, * ) 'a prescribed surface flux is not allo', & 'wed with ', sq, '_surface_initial_change (/=0) = ', & q_surface_initial_change CALL message( 'check_parameters', 'PA0074', 1, 2, 0, 6, 0 ) ENDIF ENDIF ! !-- Boundary conditions for horizontal components of wind speed IF ( bc_uv_b == 'dirichlet' ) THEN ibc_uv_b = 0 ELSEIF ( bc_uv_b == 'neumann' ) THEN ibc_uv_b = 1 IF ( prandtl_layer ) THEN message_string = 'boundary condition: bc_uv_b = "' // & TRIM( bc_uv_b ) // '" is not allowed with prandtl_layer = .TRUE.' CALL message( 'check_parameters', 'PA0075', 1, 2, 0, 6, 0 ) ENDIF ELSE message_string = 'unknown boundary condition: bc_uv_b = "' // & TRIM( bc_uv_b ) // '"' CALL message( 'check_parameters', 'PA0076', 1, 2, 0, 6, 0 ) ENDIF ! !-- In case of coupled simulations u and v at the ground in atmosphere will be !-- assigned with the u and v values of the ocean surface IF ( coupling_mode == 'atmosphere_to_ocean' ) THEN ibc_uv_b = 2 ENDIF IF ( coupling_mode == 'ocean_to_atmosphere' ) THEN bc_uv_t = 'neumann' ibc_uv_t = 1 ELSE IF ( bc_uv_t == 'dirichlet' .OR. bc_uv_t == 'dirichlet_0' ) THEN ibc_uv_t = 0 IF ( bc_uv_t == 'dirichlet_0' ) THEN ! !-- Velocities for the initial u,v-profiles are set zero at the top !-- in case of dirichlet_0 conditions u_init(nzt+1) = 0.0 v_init(nzt+1) = 0.0 ENDIF ELSEIF ( bc_uv_t == 'neumann' ) THEN ibc_uv_t = 1 ELSE message_string = 'unknown boundary condition: bc_uv_t = "' // & TRIM( bc_uv_t ) // '"' CALL message( 'check_parameters', 'PA0077', 1, 2, 0, 6, 0 ) ENDIF ENDIF ! !-- Compute and check, respectively, the Rayleigh Damping parameter IF ( rayleigh_damping_factor == -1.0 ) THEN rayleigh_damping_factor = 0.0 ELSE IF ( rayleigh_damping_factor < 0.0 .OR. rayleigh_damping_factor > 1.0 ) & THEN WRITE( message_string, * ) 'rayleigh_damping_factor = ', & rayleigh_damping_factor, ' out of range [0.0,1.0]' CALL message( 'check_parameters', 'PA0078', 1, 2, 0, 6, 0 ) ENDIF ENDIF IF ( rayleigh_damping_height == -1.0 ) THEN IF ( .NOT. ocean ) THEN rayleigh_damping_height = 0.66666666666 * zu(nzt) ELSE rayleigh_damping_height = 0.66666666666 * zu(nzb) ENDIF ELSE IF ( .NOT. ocean ) THEN IF ( rayleigh_damping_height < 0.0 .OR. & rayleigh_damping_height > zu(nzt) ) THEN WRITE( message_string, * ) 'rayleigh_damping_height = ', & rayleigh_damping_height, ' out of range [0.0,', zu(nzt), ']' CALL message( 'check_parameters', 'PA0079', 1, 2, 0, 6, 0 ) ENDIF ELSE IF ( rayleigh_damping_height > 0.0 .OR. & rayleigh_damping_height < zu(nzb) ) THEN WRITE( message_string, * ) 'rayleigh_damping_height = ', & rayleigh_damping_height, ' out of range [0.0,', zu(nzb), ']' CALL message( 'check_parameters', 'PA0079', 1, 2, 0, 6, 0 ) ENDIF ENDIF ENDIF ! !-- Check number of chosen statistic regions. More than 10 regions are not !-- allowed, because so far no more than 10 corresponding output files can !-- be opened (cf. check_open) IF ( statistic_regions > 9 .OR. statistic_regions < 0 ) THEN WRITE ( message_string, * ) 'number of statistic_regions = ', & statistic_regions+1, ' but only 10 regions are allowed' CALL message( 'check_parameters', 'PA0082', 1, 2, 0, 6, 0 ) ENDIF IF ( normalizing_region > statistic_regions .OR. & normalizing_region < 0) THEN WRITE ( message_string, * ) 'normalizing_region = ', & normalizing_region, ' must be >= 0 and <= ',statistic_regions, & ' (value of statistic_regions)' CALL message( 'check_parameters', 'PA0083', 1, 2, 0, 6, 0 ) ENDIF ! !-- Check the interval for sorting particles. !-- Using particles as cloud droplets requires sorting after each timestep. IF ( dt_sort_particles /= 0.0 .AND. cloud_droplets ) THEN dt_sort_particles = 0.0 message_string = 'dt_sort_particles is reset to 0.0 because of cloud' //& '_droplets = .TRUE.' CALL message( 'check_parameters', 'PA0084', 0, 1, 0, 6, 0 ) ENDIF ! !-- Set the default intervals for data output, if necessary !-- NOTE: dt_dosp has already been set in package_parin IF ( dt_data_output /= 9999999.9 ) THEN IF ( dt_dopr == 9999999.9 ) dt_dopr = dt_data_output IF ( dt_dopts == 9999999.9 ) dt_dopts = dt_data_output IF ( dt_do2d_xy == 9999999.9 ) dt_do2d_xy = dt_data_output IF ( dt_do2d_xz == 9999999.9 ) dt_do2d_xz = dt_data_output IF ( dt_do2d_yz == 9999999.9 ) dt_do2d_yz = dt_data_output IF ( dt_do3d == 9999999.9 ) dt_do3d = dt_data_output IF ( dt_data_output_av == 9999999.9 ) dt_data_output_av = dt_data_output DO mid = 1, max_masks IF ( dt_domask(mid) == 9999999.9 ) dt_domask(mid) = dt_data_output ENDDO ENDIF ! !-- Set the default skip time intervals for data output, if necessary IF ( skip_time_dopr == 9999999.9 ) & skip_time_dopr = skip_time_data_output IF ( skip_time_dosp == 9999999.9 ) & skip_time_dosp = skip_time_data_output IF ( skip_time_do2d_xy == 9999999.9 ) & skip_time_do2d_xy = skip_time_data_output IF ( skip_time_do2d_xz == 9999999.9 ) & skip_time_do2d_xz = skip_time_data_output IF ( skip_time_do2d_yz == 9999999.9 ) & skip_time_do2d_yz = skip_time_data_output IF ( skip_time_do3d == 9999999.9 ) & skip_time_do3d = skip_time_data_output IF ( skip_time_data_output_av == 9999999.9 ) & skip_time_data_output_av = skip_time_data_output DO mid = 1, max_masks IF ( skip_time_domask(mid) == 9999999.9 ) & skip_time_domask(mid) = skip_time_data_output ENDDO ! !-- Check the average intervals (first for 3d-data, then for profiles and !-- spectra) IF ( averaging_interval > dt_data_output_av ) THEN WRITE( message_string, * ) 'averaging_interval = ', & averaging_interval, ' must be <= dt_data_output = ', dt_data_output CALL message( 'check_parameters', 'PA0085', 1, 2, 0, 6, 0 ) ENDIF IF ( averaging_interval_pr == 9999999.9 ) THEN averaging_interval_pr = averaging_interval ENDIF IF ( averaging_interval_pr > dt_dopr ) THEN WRITE( message_string, * ) 'averaging_interval_pr = ', & averaging_interval_pr, ' must be <= dt_dopr = ', dt_dopr CALL message( 'check_parameters', 'PA0086', 1, 2, 0, 6, 0 ) ENDIF IF ( averaging_interval_sp == 9999999.9 ) THEN averaging_interval_sp = averaging_interval ENDIF IF ( averaging_interval_sp > dt_dosp ) THEN WRITE( message_string, * ) 'averaging_interval_sp = ', & averaging_interval_sp, ' must be <= dt_dosp = ', dt_dosp CALL message( 'check_parameters', 'PA0087', 1, 2, 0, 6, 0 ) ENDIF ! !-- Set the default interval for profiles entering the temporal average IF ( dt_averaging_input_pr == 9999999.9 ) THEN dt_averaging_input_pr = dt_averaging_input ENDIF ! !-- Set the default interval for the output of timeseries to a reasonable !-- value (tries to minimize the number of calls of flow_statistics) IF ( dt_dots == 9999999.9 ) THEN IF ( averaging_interval_pr == 0.0 ) THEN dt_dots = MIN( dt_run_control, dt_dopr ) ELSE dt_dots = MIN( dt_run_control, dt_averaging_input_pr ) ENDIF ENDIF ! !-- Check the sample rate for averaging (first for 3d-data, then for profiles) IF ( dt_averaging_input > averaging_interval ) THEN WRITE( message_string, * ) 'dt_averaging_input = ', & dt_averaging_input, ' must be <= averaging_interval = ', & averaging_interval CALL message( 'check_parameters', 'PA0088', 1, 2, 0, 6, 0 ) ENDIF IF ( dt_averaging_input_pr > averaging_interval_pr ) THEN WRITE( message_string, * ) 'dt_averaging_input_pr = ', & dt_averaging_input_pr, ' must be <= averaging_interval_pr = ', & averaging_interval_pr CALL message( 'check_parameters', 'PA0089', 1, 2, 0, 6, 0 ) ENDIF ! !-- Set the default value for the integration interval of precipitation amount IF ( precipitation ) THEN IF ( precipitation_amount_interval == 9999999.9 ) THEN precipitation_amount_interval = dt_do2d_xy ELSE IF ( precipitation_amount_interval > dt_do2d_xy ) THEN WRITE( message_string, * ) 'precipitation_amount_interval = ', & precipitation_amount_interval, ' must not be larger than ', & 'dt_do2d_xy = ', dt_do2d_xy CALL message( 'check_parameters', 'PA0090', 1, 2, 0, 6, 0 ) ENDIF ENDIF ENDIF ! !-- Determine the number of output profiles and check whether they are !-- permissible DO WHILE ( data_output_pr(dopr_n+1) /= ' ' ) dopr_n = dopr_n + 1 i = dopr_n ! !-- Determine internal profile number (for hom, homs) !-- and store height levels SELECT CASE ( TRIM( data_output_pr(i) ) ) CASE ( 'u', '#u' ) dopr_index(i) = 1 dopr_unit(i) = 'm/s' hom(:,2,1,:) = SPREAD( zu, 2, statistic_regions+1 ) IF ( data_output_pr(i)(1:1) == '#' ) THEN dopr_initial_index(i) = 5 hom(:,2,5,:) = SPREAD( zu, 2, statistic_regions+1 ) data_output_pr(i) = data_output_pr(i)(2:) ENDIF CASE ( 'v', '#v' ) dopr_index(i) = 2 dopr_unit(i) = 'm/s' hom(:,2,2,:) = SPREAD( zu, 2, statistic_regions+1 ) IF ( data_output_pr(i)(1:1) == '#' ) THEN dopr_initial_index(i) = 6 hom(:,2,6,:) = SPREAD( zu, 2, statistic_regions+1 ) data_output_pr(i) = data_output_pr(i)(2:) ENDIF CASE ( 'w' ) dopr_index(i) = 3 dopr_unit(i) = 'm/s' hom(:,2,3,:) = SPREAD( zw, 2, statistic_regions+1 ) CASE ( 'pt', '#pt' ) IF ( .NOT. cloud_physics ) THEN dopr_index(i) = 4 dopr_unit(i) = 'K' hom(:,2,4,:) = SPREAD( zu, 2, statistic_regions+1 ) IF ( data_output_pr(i)(1:1) == '#' ) THEN dopr_initial_index(i) = 7 hom(:,2,7,:) = SPREAD( zu, 2, statistic_regions+1 ) hom(nzb,2,7,:) = 0.0 ! because zu(nzb) is negative data_output_pr(i) = data_output_pr(i)(2:) ENDIF ELSE dopr_index(i) = 43 dopr_unit(i) = 'K' hom(:,2,43,:) = SPREAD( zu, 2, statistic_regions+1 ) IF ( data_output_pr(i)(1:1) == '#' ) THEN dopr_initial_index(i) = 28 hom(:,2,28,:) = SPREAD( zu, 2, statistic_regions+1 ) hom(nzb,2,28,:) = 0.0 ! because zu(nzb) is negative data_output_pr(i) = data_output_pr(i)(2:) ENDIF ENDIF CASE ( 'e' ) dopr_index(i) = 8 dopr_unit(i) = 'm2/s2' hom(:,2,8,:) = SPREAD( zu, 2, statistic_regions+1 ) hom(nzb,2,8,:) = 0.0 CASE ( 'km', '#km' ) dopr_index(i) = 9 dopr_unit(i) = 'm2/s' hom(:,2,9,:) = SPREAD( zu, 2, statistic_regions+1 ) hom(nzb,2,9,:) = 0.0 IF ( data_output_pr(i)(1:1) == '#' ) THEN dopr_initial_index(i) = 23 hom(:,2,23,:) = hom(:,2,9,:) data_output_pr(i) = data_output_pr(i)(2:) ENDIF CASE ( 'kh', '#kh' ) dopr_index(i) = 10 dopr_unit(i) = 'm2/s' hom(:,2,10,:) = SPREAD( zu, 2, statistic_regions+1 ) hom(nzb,2,10,:) = 0.0 IF ( data_output_pr(i)(1:1) == '#' ) THEN dopr_initial_index(i) = 24 hom(:,2,24,:) = hom(:,2,10,:) data_output_pr(i) = data_output_pr(i)(2:) ENDIF CASE ( 'l', '#l' ) dopr_index(i) = 11 dopr_unit(i) = 'm' hom(:,2,11,:) = SPREAD( zu, 2, statistic_regions+1 ) hom(nzb,2,11,:) = 0.0 IF ( data_output_pr(i)(1:1) == '#' ) THEN dopr_initial_index(i) = 25 hom(:,2,25,:) = hom(:,2,11,:) data_output_pr(i) = data_output_pr(i)(2:) ENDIF CASE ( 'w"u"' ) dopr_index(i) = 12 dopr_unit(i) = 'm2/s2' hom(:,2,12,:) = SPREAD( zw, 2, statistic_regions+1 ) IF ( prandtl_layer ) hom(nzb,2,12,:) = zu(1) CASE ( 'w*u*' ) dopr_index(i) = 13 dopr_unit(i) = 'm2/s2' hom(:,2,13,:) = SPREAD( zw, 2, statistic_regions+1 ) CASE ( 'w"v"' ) dopr_index(i) = 14 dopr_unit(i) = 'm2/s2' hom(:,2,14,:) = SPREAD( zw, 2, statistic_regions+1 ) IF ( prandtl_layer ) hom(nzb,2,14,:) = zu(1) CASE ( 'w*v*' ) dopr_index(i) = 15 dopr_unit(i) = 'm2/s2' hom(:,2,15,:) = SPREAD( zw, 2, statistic_regions+1 ) CASE ( 'w"pt"' ) dopr_index(i) = 16 dopr_unit(i) = 'K m/s' hom(:,2,16,:) = SPREAD( zw, 2, statistic_regions+1 ) CASE ( 'w*pt*' ) dopr_index(i) = 17 dopr_unit(i) = 'K m/s' hom(:,2,17,:) = SPREAD( zw, 2, statistic_regions+1 ) CASE ( 'wpt' ) dopr_index(i) = 18 dopr_unit(i) = 'K m/s' hom(:,2,18,:) = SPREAD( zw, 2, statistic_regions+1 ) CASE ( 'wu' ) dopr_index(i) = 19 dopr_unit(i) = 'm2/s2' hom(:,2,19,:) = SPREAD( zw, 2, statistic_regions+1 ) IF ( prandtl_layer ) hom(nzb,2,19,:) = zu(1) CASE ( 'wv' ) dopr_index(i) = 20 dopr_unit(i) = 'm2/s2' hom(:,2,20,:) = SPREAD( zw, 2, statistic_regions+1 ) IF ( prandtl_layer ) hom(nzb,2,20,:) = zu(1) CASE ( 'w*pt*BC' ) dopr_index(i) = 21 dopr_unit(i) = 'K m/s' hom(:,2,21,:) = SPREAD( zw, 2, statistic_regions+1 ) CASE ( 'wptBC' ) dopr_index(i) = 22 dopr_unit(i) = 'K m/s' hom(:,2,22,:) = SPREAD( zw, 2, statistic_regions+1 ) CASE ( 'sa', '#sa' ) IF ( .NOT. ocean ) THEN message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for ocean = .FALSE.' CALL message( 'check_parameters', 'PA0091', 1, 2, 0, 6, 0 ) ELSE dopr_index(i) = 23 dopr_unit(i) = 'psu' hom(:,2,23,:) = SPREAD( zu, 2, statistic_regions+1 ) IF ( data_output_pr(i)(1:1) == '#' ) THEN dopr_initial_index(i) = 26 hom(:,2,26,:) = SPREAD( zu, 2, statistic_regions+1 ) hom(nzb,2,26,:) = 0.0 ! weil zu(nzb) negativ ist data_output_pr(i) = data_output_pr(i)(2:) ENDIF ENDIF CASE ( 'u*2' ) dopr_index(i) = 30 dopr_unit(i) = 'm2/s2' hom(:,2,30,:) = SPREAD( zu, 2, statistic_regions+1 ) CASE ( 'v*2' ) dopr_index(i) = 31 dopr_unit(i) = 'm2/s2' hom(:,2,31,:) = SPREAD( zu, 2, statistic_regions+1 ) CASE ( 'w*2' ) dopr_index(i) = 32 dopr_unit(i) = 'm2/s2' hom(:,2,32,:) = SPREAD( zw, 2, statistic_regions+1 ) CASE ( 'pt*2' ) dopr_index(i) = 33 dopr_unit(i) = 'K2' hom(:,2,33,:) = SPREAD( zu, 2, statistic_regions+1 ) CASE ( 'e*' ) dopr_index(i) = 34 dopr_unit(i) = 'm2/s2' hom(:,2,34,:) = SPREAD( zu, 2, statistic_regions+1 ) CASE ( 'w*2pt*' ) dopr_index(i) = 35 dopr_unit(i) = 'K m2/s2' hom(:,2,35,:) = SPREAD( zw, 2, statistic_regions+1 ) CASE ( 'w*pt*2' ) dopr_index(i) = 36 dopr_unit(i) = 'K2 m/s' hom(:,2,36,:) = SPREAD( zw, 2, statistic_regions+1 ) CASE ( 'w*e*' ) dopr_index(i) = 37 dopr_unit(i) = 'm3/s3' hom(:,2,37,:) = SPREAD( zw, 2, statistic_regions+1 ) CASE ( 'w*3' ) dopr_index(i) = 38 dopr_unit(i) = 'm3/s3' hom(:,2,38,:) = SPREAD( zw, 2, statistic_regions+1 ) CASE ( 'Sw' ) dopr_index(i) = 39 dopr_unit(i) = 'none' hom(:,2,39,:) = SPREAD( zw, 2, statistic_regions+1 ) CASE ( 'p' ) dopr_index(i) = 40 dopr_unit(i) = 'Pa' hom(:,2,40,:) = SPREAD( zu, 2, statistic_regions+1 ) CASE ( 'q', '#q' ) IF ( .NOT. humidity ) THEN message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for humidity = .FALSE.' CALL message( 'check_parameters', 'PA0092', 1, 2, 0, 6, 0 ) ELSE dopr_index(i) = 41 dopr_unit(i) = 'kg/kg' hom(:,2,41,:) = SPREAD( zu, 2, statistic_regions+1 ) IF ( data_output_pr(i)(1:1) == '#' ) THEN dopr_initial_index(i) = 26 hom(:,2,26,:) = SPREAD( zu, 2, statistic_regions+1 ) hom(nzb,2,26,:) = 0.0 ! weil zu(nzb) negativ ist data_output_pr(i) = data_output_pr(i)(2:) ENDIF ENDIF CASE ( 's', '#s' ) IF ( .NOT. passive_scalar ) THEN message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for passive_scalar = .FALSE.' CALL message( 'check_parameters', 'PA0093', 1, 2, 0, 6, 0 ) ELSE dopr_index(i) = 41 dopr_unit(i) = 'kg/m3' hom(:,2,41,:) = SPREAD( zu, 2, statistic_regions+1 ) IF ( data_output_pr(i)(1:1) == '#' ) THEN dopr_initial_index(i) = 26 hom(:,2,26,:) = SPREAD( zu, 2, statistic_regions+1 ) hom(nzb,2,26,:) = 0.0 ! weil zu(nzb) negativ ist data_output_pr(i) = data_output_pr(i)(2:) ENDIF ENDIF CASE ( 'qv', '#qv' ) IF ( .NOT. cloud_physics ) THEN dopr_index(i) = 41 dopr_unit(i) = 'kg/kg' hom(:,2,41,:) = SPREAD( zu, 2, statistic_regions+1 ) IF ( data_output_pr(i)(1:1) == '#' ) THEN dopr_initial_index(i) = 26 hom(:,2,26,:) = SPREAD( zu, 2, statistic_regions+1 ) hom(nzb,2,26,:) = 0.0 ! weil zu(nzb) negativ ist data_output_pr(i) = data_output_pr(i)(2:) ENDIF ELSE dopr_index(i) = 42 dopr_unit(i) = 'kg/kg' hom(:,2,42,:) = SPREAD( zu, 2, statistic_regions+1 ) IF ( data_output_pr(i)(1:1) == '#' ) THEN dopr_initial_index(i) = 27 hom(:,2,27,:) = SPREAD( zu, 2, statistic_regions+1 ) hom(nzb,2,27,:) = 0.0 ! weil zu(nzb) negativ ist data_output_pr(i) = data_output_pr(i)(2:) ENDIF ENDIF CASE ( 'lpt', '#lpt' ) IF ( .NOT. cloud_physics ) THEN message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for cloud_physics = .FALSE.' CALL message( 'check_parameters', 'PA0094', 1, 2, 0, 6, 0 ) ELSE dopr_index(i) = 4 dopr_unit(i) = 'K' hom(:,2,4,:) = SPREAD( zu, 2, statistic_regions+1 ) IF ( data_output_pr(i)(1:1) == '#' ) THEN dopr_initial_index(i) = 7 hom(:,2,7,:) = SPREAD( zu, 2, statistic_regions+1 ) hom(nzb,2,7,:) = 0.0 ! weil zu(nzb) negativ ist data_output_pr(i) = data_output_pr(i)(2:) ENDIF ENDIF CASE ( 'vpt', '#vpt' ) dopr_index(i) = 44 dopr_unit(i) = 'K' hom(:,2,44,:) = SPREAD( zu, 2, statistic_regions+1 ) IF ( data_output_pr(i)(1:1) == '#' ) THEN dopr_initial_index(i) = 29 hom(:,2,29,:) = SPREAD( zu, 2, statistic_regions+1 ) hom(nzb,2,29,:) = 0.0 ! weil zu(nzb) negativ ist data_output_pr(i) = data_output_pr(i)(2:) ENDIF CASE ( 'w"vpt"' ) dopr_index(i) = 45 dopr_unit(i) = 'K m/s' hom(:,2,45,:) = SPREAD( zw, 2, statistic_regions+1 ) CASE ( 'w*vpt*' ) dopr_index(i) = 46 dopr_unit(i) = 'K m/s' hom(:,2,46,:) = SPREAD( zw, 2, statistic_regions+1 ) CASE ( 'wvpt' ) dopr_index(i) = 47 dopr_unit(i) = 'K m/s' hom(:,2,47,:) = SPREAD( zw, 2, statistic_regions+1 ) CASE ( 'w"q"' ) IF ( .NOT. humidity ) THEN message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for humidity = .FALSE.' CALL message( 'check_parameters', 'PA0092', 1, 2, 0, 6, 0 ) ELSE dopr_index(i) = 48 dopr_unit(i) = 'kg/kg m/s' hom(:,2,48,:) = SPREAD( zw, 2, statistic_regions+1 ) ENDIF CASE ( 'w*q*' ) IF ( .NOT. humidity ) THEN message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for humidity = .FALSE.' CALL message( 'check_parameters', 'PA0092', 1, 2, 0, 6, 0 ) ELSE dopr_index(i) = 49 dopr_unit(i) = 'kg/kg m/s' hom(:,2,49,:) = SPREAD( zw, 2, statistic_regions+1 ) ENDIF CASE ( 'wq' ) IF ( .NOT. humidity ) THEN message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for humidity = .FALSE.' CALL message( 'check_parameters', 'PA0092', 1, 2, 0, 6, 0 ) ELSE dopr_index(i) = 50 dopr_unit(i) = 'kg/kg m/s' hom(:,2,50,:) = SPREAD( zw, 2, statistic_regions+1 ) ENDIF CASE ( 'w"s"' ) IF ( .NOT. passive_scalar ) THEN message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for passive_scalar = .FALSE.' CALL message( 'check_parameters', 'PA0093', 1, 2, 0, 6, 0 ) ELSE dopr_index(i) = 48 dopr_unit(i) = 'kg/m3 m/s' hom(:,2,48,:) = SPREAD( zw, 2, statistic_regions+1 ) ENDIF CASE ( 'w*s*' ) IF ( .NOT. passive_scalar ) THEN message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for passive_scalar = .FALSE.' CALL message( 'check_parameters', 'PA0093', 1, 2, 0, 6, 0 ) ELSE dopr_index(i) = 49 dopr_unit(i) = 'kg/m3 m/s' hom(:,2,49,:) = SPREAD( zw, 2, statistic_regions+1 ) ENDIF CASE ( 'ws' ) IF ( .NOT. passive_scalar ) THEN message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for passive_scalar = .FALSE.' CALL message( 'check_parameters', 'PA0093', 1, 2, 0, 6, 0 ) ELSE dopr_index(i) = 50 dopr_unit(i) = 'kg/m3 m/s' hom(:,2,50,:) = SPREAD( zw, 2, statistic_regions+1 ) ENDIF CASE ( 'w"qv"' ) IF ( humidity .AND. .NOT. cloud_physics ) & THEN dopr_index(i) = 48 dopr_unit(i) = 'kg/kg m/s' hom(:,2,48,:) = SPREAD( zw, 2, statistic_regions+1 ) ELSEIF( humidity .AND. cloud_physics ) THEN dopr_index(i) = 51 dopr_unit(i) = 'kg/kg m/s' hom(:,2,51,:) = SPREAD( zw, 2, statistic_regions+1 ) ELSE message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for cloud_physics = .FALSE. an&' // & 'd humidity = .FALSE.' CALL message( 'check_parameters', 'PA0095', 1, 2, 0, 6, 0 ) ENDIF CASE ( 'w*qv*' ) IF ( humidity .AND. .NOT. cloud_physics ) & THEN dopr_index(i) = 49 dopr_unit(i) = 'kg/kg m/s' hom(:,2,49,:) = SPREAD( zw, 2, statistic_regions+1 ) ELSEIF( humidity .AND. cloud_physics ) THEN dopr_index(i) = 52 dopr_unit(i) = 'kg/kg m/s' hom(:,2,52,:) = SPREAD( zw, 2, statistic_regions+1 ) ELSE message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for cloud_physics = .FALSE. an&' // & 'd humidity = .FALSE.' CALL message( 'check_parameters', 'PA0095', 1, 2, 0, 6, 0 ) ENDIF CASE ( 'wqv' ) IF ( humidity .AND. .NOT. cloud_physics ) & THEN dopr_index(i) = 50 dopr_unit(i) = 'kg/kg m/s' hom(:,2,50,:) = SPREAD( zw, 2, statistic_regions+1 ) ELSEIF( humidity .AND. cloud_physics ) THEN dopr_index(i) = 53 dopr_unit(i) = 'kg/kg m/s' hom(:,2,53,:) = SPREAD( zw, 2, statistic_regions+1 ) ELSE message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for cloud_physics = .FALSE. an&' // & 'd humidity = .FALSE.' CALL message( 'check_parameters', 'PA0095', 1, 2, 0, 6, 0 ) ENDIF CASE ( 'ql' ) IF ( .NOT. cloud_physics .AND. .NOT. cloud_droplets ) THEN message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for cloud_physics = .FALSE. or' // & '&cloud_droplets = .FALSE.' CALL message( 'check_parameters', 'PA0096', 1, 2, 0, 6, 0 ) ELSE dopr_index(i) = 54 dopr_unit(i) = 'kg/kg' hom(:,2,54,:) = SPREAD( zu, 2, statistic_regions+1 ) ENDIF CASE ( 'w*u*u*:dz' ) dopr_index(i) = 55 dopr_unit(i) = 'm2/s3' hom(:,2,55,:) = SPREAD( zu, 2, statistic_regions+1 ) CASE ( 'w*p*:dz' ) dopr_index(i) = 56 dopr_unit(i) = 'm2/s3' hom(:,2,56,:) = SPREAD( zw, 2, statistic_regions+1 ) CASE ( 'w"e:dz' ) dopr_index(i) = 57 dopr_unit(i) = 'm2/s3' hom(:,2,57,:) = SPREAD( zu, 2, statistic_regions+1 ) CASE ( 'u"pt"' ) dopr_index(i) = 58 dopr_unit(i) = 'K m/s' hom(:,2,58,:) = SPREAD( zu, 2, statistic_regions+1 ) CASE ( 'u*pt*' ) dopr_index(i) = 59 dopr_unit(i) = 'K m/s' hom(:,2,59,:) = SPREAD( zu, 2, statistic_regions+1 ) CASE ( 'upt_t' ) dopr_index(i) = 60 dopr_unit(i) = 'K m/s' hom(:,2,60,:) = SPREAD( zu, 2, statistic_regions+1 ) CASE ( 'v"pt"' ) dopr_index(i) = 61 dopr_unit(i) = 'K m/s' hom(:,2,61,:) = SPREAD( zu, 2, statistic_regions+1 ) CASE ( 'v*pt*' ) dopr_index(i) = 62 dopr_unit(i) = 'K m/s' hom(:,2,62,:) = SPREAD( zu, 2, statistic_regions+1 ) CASE ( 'vpt_t' ) dopr_index(i) = 63 dopr_unit(i) = 'K m/s' hom(:,2,63,:) = SPREAD( zu, 2, statistic_regions+1 ) CASE ( 'rho' ) IF ( .NOT. ocean ) THEN message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for ocean = .FALSE.' CALL message( 'check_parameters', 'PA0091', 1, 2, 0, 6, 0 ) ELSE dopr_index(i) = 64 dopr_unit(i) = 'kg/m3' hom(:,2,64,:) = SPREAD( zu, 2, statistic_regions+1 ) IF ( data_output_pr(i)(1:1) == '#' ) THEN dopr_initial_index(i) = 77 hom(:,2,77,:) = SPREAD( zu, 2, statistic_regions+1 ) hom(nzb,2,77,:) = 0.0 ! because zu(nzb) is negative data_output_pr(i) = data_output_pr(i)(2:) ENDIF ENDIF CASE ( 'w"sa"' ) IF ( .NOT. ocean ) THEN message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for ocean = .FALSE.' CALL message( 'check_parameters', 'PA0091', 1, 2, 0, 6, 0 ) ELSE dopr_index(i) = 65 dopr_unit(i) = 'psu m/s' hom(:,2,65,:) = SPREAD( zw, 2, statistic_regions+1 ) ENDIF CASE ( 'w*sa*' ) IF ( .NOT. ocean ) THEN message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for ocean = .FALSE.' CALL message( 'check_parameters', 'PA0091', 1, 2, 0, 6, 0 ) ELSE dopr_index(i) = 66 dopr_unit(i) = 'psu m/s' hom(:,2,66,:) = SPREAD( zw, 2, statistic_regions+1 ) ENDIF CASE ( 'wsa' ) IF ( .NOT. ocean ) THEN message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for ocean = .FALSE.' CALL message( 'check_parameters', 'PA0091', 1, 2, 0, 6, 0 ) ELSE dopr_index(i) = 67 dopr_unit(i) = 'psu m/s' hom(:,2,67,:) = SPREAD( zw, 2, statistic_regions+1 ) ENDIF CASE ( 'w*p*' ) dopr_index(i) = 68 dopr_unit(i) = 'm3/s3' hom(:,2,68,:) = SPREAD( zu, 2, statistic_regions+1 ) CASE ( 'w"e' ) dopr_index(i) = 69 dopr_unit(i) = 'm3/s3' hom(:,2,69,:) = SPREAD( zu, 2, statistic_regions+1 ) CASE ( 'q*2' ) IF ( .NOT. humidity ) THEN message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for humidity = .FALSE.' CALL message( 'check_parameters', 'PA0092', 1, 2, 0, 6, 0 ) ELSE dopr_index(i) = 70 dopr_unit(i) = 'kg2/kg2' hom(:,2,70,:) = SPREAD( zu, 2, statistic_regions+1 ) ENDIF CASE ( 'prho' ) IF ( .NOT. ocean ) THEN message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for ocean = .FALSE.' CALL message( 'check_parameters', 'PA0091', 1, 2, 0, 6, 0 ) ELSE dopr_index(i) = 71 dopr_unit(i) = 'kg/m3' hom(:,2,71,:) = SPREAD( zu, 2, statistic_regions+1 ) ENDIF CASE ( 'hyp' ) dopr_index(i) = 72 dopr_unit(i) = 'dbar' hom(:,2,72,:) = SPREAD( zu, 2, statistic_regions+1 ) CASE ( 'nr' ) IF ( .NOT. cloud_physics ) THEN message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for cloud_physics = .FALSE.' CALL message( 'check_parameters', 'PA0094', 1, 2, 0, 6, 0 ) ELSEIF ( icloud_scheme /= 0 ) THEN message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for cloud_scheme /= seifert_beheng' CALL message( 'check_parameters', 'PA0358', 1, 2, 0, 6, 0 ) ELSEIF ( .NOT. precipitation ) THEN message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for precipitation = .FALSE.' CALL message( 'check_parameters', 'PA0361', 1, 2, 0, 6, 0 ) ELSE dopr_index(i) = 73 dopr_unit(i) = '1/m3' hom(:,2,73,:) = SPREAD( zu, 2, statistic_regions+1 ) ENDIF CASE ( 'qr' ) IF ( .NOT. cloud_physics ) THEN message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for cloud_physics = .FALSE.' CALL message( 'check_parameters', 'PA0094', 1, 2, 0, 6, 0 ) ELSEIF ( icloud_scheme /= 0 ) THEN message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for cloud_scheme /= seifert_beheng' CALL message( 'check_parameters', 'PA0358', 1, 2, 0, 6, 0 ) ELSEIF ( .NOT. precipitation ) THEN message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for precipitation = .FALSE.' CALL message( 'check_parameters', 'PA0361', 1, 2, 0, 6, 0 ) ELSE dopr_index(i) = 74 dopr_unit(i) = 'kg/kg' hom(:,2,74,:) = SPREAD( zu, 2, statistic_regions+1 ) ENDIF CASE ( 'qc' ) IF ( .NOT. cloud_physics ) THEN message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for cloud_physics = .FALSE.' CALL message( 'check_parameters', 'PA0094', 1, 2, 0, 6, 0 ) ELSEIF ( icloud_scheme /= 0 ) THEN message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for cloud_scheme /= seifert_beheng' CALL message( 'check_parameters', 'PA0358', 1, 2, 0, 6, 0 ) ELSE dopr_index(i) = 75 dopr_unit(i) = 'kg/kg' hom(:,2,75,:) = SPREAD( zu, 2, statistic_regions+1 ) ENDIF CASE ( 'prr' ) IF ( .NOT. cloud_physics ) THEN message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for cloud_physics = .FALSE.' CALL message( 'check_parameters', 'PA0094', 1, 2, 0, 6, 0 ) ELSEIF ( icloud_scheme /= 0 ) THEN message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for cloud_scheme /= seifert_beheng' CALL message( 'check_parameters', 'PA0358', 1, 2, 0, 6, 0 ) ELSEIF ( .NOT. precipitation ) THEN message_string = 'data_output_pr = ' // & TRIM( data_output_pr(i) ) // ' is not imp' // & 'lemented for precipitation = .FALSE.' CALL message( 'check_parameters', 'PA0361', 1, 2, 0, 6, 0 ) ELSE dopr_index(i) = 76 dopr_unit(i) = 'kg/kg m/s' hom(:,2,76,:) = SPREAD( zu, 2, statistic_regions+1 ) ENDIF CASE DEFAULT CALL user_check_data_output_pr( data_output_pr(i), i, unit ) IF ( unit == 'illegal' ) THEN IF ( data_output_pr_user(1) /= ' ' ) THEN message_string = 'illegal value for data_output_pr or ' // & 'data_output_pr_user = "' // & TRIM( data_output_pr(i) ) // '"' CALL message( 'check_parameters', 'PA0097', 1, 2, 0, 6, 0 ) ELSE message_string = 'illegal value for data_output_pr = "' // & TRIM( data_output_pr(i) ) // '"' CALL message( 'check_parameters', 'PA0098', 1, 2, 0, 6, 0 ) ENDIF ENDIF END SELECT ENDDO ! !-- Append user-defined data output variables to the standard data output IF ( data_output_user(1) /= ' ' ) THEN i = 1 DO WHILE ( data_output(i) /= ' ' .AND. i <= 100 ) i = i + 1 ENDDO j = 1 DO WHILE ( data_output_user(j) /= ' ' .AND. j <= 100 ) IF ( i > 100 ) THEN message_string = 'number of output quantitities given by data' // & '_output and data_output_user exceeds the limit of 100' CALL message( 'check_parameters', 'PA0102', 1, 2, 0, 6, 0 ) ENDIF data_output(i) = data_output_user(j) i = i + 1 j = j + 1 ENDDO ENDIF ! !-- Check and set steering parameters for 2d/3d data output and averaging i = 1 DO WHILE ( data_output(i) /= ' ' .AND. i <= 100 ) ! !-- Check for data averaging ilen = LEN_TRIM( data_output(i) ) j = 0 ! no data averaging IF ( ilen > 3 ) THEN IF ( data_output(i)(ilen-2:ilen) == '_av' ) THEN j = 1 ! data averaging data_output(i) = data_output(i)(1:ilen-3) ENDIF ENDIF ! !-- Check for cross section or volume data ilen = LEN_TRIM( data_output(i) ) k = 0 ! 3d data var = data_output(i)(1:ilen) IF ( ilen > 3 ) THEN IF ( data_output(i)(ilen-2:ilen) == '_xy' .OR. & data_output(i)(ilen-2:ilen) == '_xz' .OR. & data_output(i)(ilen-2:ilen) == '_yz' ) THEN k = 1 ! 2d data var = data_output(i)(1:ilen-3) ENDIF ENDIF ! !-- Check for allowed value and set units SELECT CASE ( TRIM( var ) ) CASE ( 'e' ) IF ( constant_diffusion ) THEN message_string = 'output of "' // TRIM( var ) // '" requi' // & 'res constant_diffusion = .FALSE.' CALL message( 'check_parameters', 'PA0103', 1, 2, 0, 6, 0 ) ENDIF unit = 'm2/s2' CASE ( 'lpt' ) IF ( .NOT. cloud_physics ) THEN message_string = 'output of "' // TRIM( var ) // '" requi' // & 'res cloud_physics = .TRUE.' CALL message( 'check_parameters', 'PA0108', 1, 2, 0, 6, 0 ) ENDIF unit = 'K' CASE ( 'nr' ) IF ( .NOT. cloud_physics ) THEN message_string = 'output of "' // TRIM( var ) // '" requi' // & 'res cloud_physics = .TRUE.' CALL message( 'check_parameters', 'PA0108', 1, 2, 0, 6, 0 ) ELSEIF ( icloud_scheme /= 0 ) THEN message_string = 'output of "' // TRIM( var ) // '" requi' // & 'res cloud_scheme = seifert_beheng' CALL message( 'check_parameters', 'PA0359', 1, 2, 0, 6, 0 ) ENDIF unit = '1/m3' CASE ( 'pc', 'pr' ) IF ( .NOT. particle_advection ) THEN message_string = 'output of "' // TRIM( var ) // '" requir' // & 'es a "particles_par"-NAMELIST in the parameter file (PARIN)' CALL message( 'check_parameters', 'PA0104', 1, 2, 0, 6, 0 ) ENDIF IF ( TRIM( var ) == 'pc' ) unit = 'number' IF ( TRIM( var ) == 'pr' ) unit = 'm' CASE ( 'prr' ) IF ( .NOT. cloud_physics ) THEN message_string = 'output of "' // TRIM( var ) // '" requi' // & 'res cloud_physics = .TRUE.' CALL message( 'check_parameters', 'PA0108', 1, 2, 0, 6, 0 ) ELSEIF ( icloud_scheme /= 0 ) THEN message_string = 'output of "' // TRIM( var ) // '" requi' // & 'res cloud_scheme = seifert_beheng' CALL message( 'check_parameters', 'PA0359', 1, 2, 0, 6, 0 ) ELSEIF ( .NOT. precipitation ) THEN message_string = 'output of "' // TRIM( var ) // '" requi' // & 'res precipitation = .TRUE.' CALL message( 'check_parameters', 'PA0112', 1, 2, 0, 6, 0 ) ENDIF unit = 'kg/kg m/s' CASE ( 'q', 'vpt' ) IF ( .NOT. humidity ) THEN message_string = 'output of "' // TRIM( var ) // '" requi' // & 'res humidity = .TRUE.' CALL message( 'check_parameters', 'PA0105', 1, 2, 0, 6, 0 ) ENDIF IF ( TRIM( var ) == 'q' ) unit = 'kg/kg' IF ( TRIM( var ) == 'vpt' ) unit = 'K' CASE ( 'qc' ) IF ( .NOT. cloud_physics ) THEN message_string = 'output of "' // TRIM( var ) // '" requi' // & 'res cloud_physics = .TRUE.' CALL message( 'check_parameters', 'PA0108', 1, 2, 0, 6, 0 ) ELSEIF ( icloud_scheme /= 0 ) THEN message_string = 'output of "' // TRIM( var ) // '" requi' // & 'res cloud_scheme = seifert_beheng' CALL message( 'check_parameters', 'PA0359', 1, 2, 0, 6, 0 ) ENDIF unit = 'kg/kg' CASE ( 'ql' ) IF ( .NOT. ( cloud_physics .OR. cloud_droplets ) ) THEN message_string = 'output of "' // TRIM( var ) // '" requi' // & 'res cloud_physics = .TRUE. or cloud_droplets = .TRUE.' CALL message( 'check_parameters', 'PA0106', 1, 2, 0, 6, 0 ) ENDIF unit = 'kg/kg' CASE ( 'ql_c', 'ql_v', 'ql_vp' ) IF ( .NOT. cloud_droplets ) THEN message_string = 'output of "' // TRIM( var ) // '" requi' // & 'res cloud_droplets = .TRUE.' CALL message( 'check_parameters', 'PA0107', 1, 2, 0, 6, 0 ) ENDIF IF ( TRIM( var ) == 'ql_c' ) unit = 'kg/kg' IF ( TRIM( var ) == 'ql_v' ) unit = 'm3' IF ( TRIM( var ) == 'ql_vp' ) unit = 'none' CASE ( 'qr' ) IF ( .NOT. cloud_physics ) THEN message_string = 'output of "' // TRIM( var ) // '" requi' // & 'res cloud_physics = .TRUE.' CALL message( 'check_parameters', 'PA0108', 1, 2, 0, 6, 0 ) ELSEIF ( icloud_scheme /= 0 ) THEN message_string = 'output of "' // TRIM( var ) // '" requi' // & 'res cloud_scheme = seifert_beheng' CALL message( 'check_parameters', 'PA0359', 1, 2, 0, 6, 0 ) ELSEIF ( .NOT. precipitation ) THEN message_string = 'output of "' // TRIM( var ) // '" requi' // & 'res precipitation = .TRUE.' CALL message( 'check_parameters', 'PA0112', 1, 2, 0, 6, 0 ) ENDIF unit = 'kg/kg' CASE ( 'qv' ) IF ( .NOT. cloud_physics ) THEN message_string = 'output of "' // TRIM( var ) // '" requi' // & 'res cloud_physics = .TRUE.' CALL message( 'check_parameters', 'PA0108', 1, 2, 0, 6, 0 ) ENDIF unit = 'kg/kg' CASE ( 'rho' ) IF ( .NOT. ocean ) THEN message_string = 'output of "' // TRIM( var ) // '" requi' // & 'res ocean = .TRUE.' CALL message( 'check_parameters', 'PA0109', 1, 2, 0, 6, 0 ) ENDIF unit = 'kg/m3' CASE ( 's' ) IF ( .NOT. passive_scalar ) THEN message_string = 'output of "' // TRIM( var ) // '" requi' // & 'res passive_scalar = .TRUE.' CALL message( 'check_parameters', 'PA0110', 1, 2, 0, 6, 0 ) ENDIF unit = 'conc' CASE ( 'sa' ) IF ( .NOT. ocean ) THEN message_string = 'output of "' // TRIM( var ) // '" requi' // & 'res ocean = .TRUE.' CALL message( 'check_parameters', 'PA0109', 1, 2, 0, 6, 0 ) ENDIF unit = 'psu' CASE ( 'u*', 't*', 'lwp*', 'pra*', 'prr*', 'qsws*', 'shf*', 'z0*', 'z0h*' ) IF ( k == 0 .OR. data_output(i)(ilen-2:ilen) /= '_xy' ) THEN message_string = 'illegal value for data_output: "' // & TRIM( var ) // '" & only 2d-horizontal ' // & 'cross sections are allowed for this value' CALL message( 'check_parameters', 'PA0111', 1, 2, 0, 6, 0 ) ENDIF IF ( TRIM( var ) == 'lwp*' .AND. .NOT. cloud_physics ) THEN message_string = 'output of "' // TRIM( var ) // '" requi' // & 'res cloud_physics = .TRUE.' CALL message( 'check_parameters', 'PA0108', 1, 2, 0, 6, 0 ) ENDIF IF ( TRIM( var ) == 'pra*' .AND. .NOT. precipitation ) THEN message_string = 'output of "' // TRIM( var ) // '" requi' // & 'res precipitation = .TRUE.' CALL message( 'check_parameters', 'PA0112', 1, 2, 0, 6, 0 ) ENDIF IF ( TRIM( var ) == 'pra*' .AND. j == 1 ) THEN message_string = 'temporal averaging of precipitation ' // & 'amount "' // TRIM( var ) // '" is not possible' CALL message( 'check_parameters', 'PA0113', 1, 2, 0, 6, 0 ) ENDIF IF ( TRIM( var ) == 'prr*' .AND. .NOT. precipitation ) THEN message_string = 'output of "' // TRIM( var ) // '" requi' // & 'res precipitation = .TRUE.' CALL message( 'check_parameters', 'PA0112', 1, 2, 0, 6, 0 ) ENDIF IF ( TRIM( var ) == 'qsws*' .AND. .NOT. humidity ) THEN message_string = 'output of "' // TRIM( var ) // '" requi' // & 'res humidity = .TRUE.' CALL message( 'check_parameters', 'PA0322', 1, 2, 0, 6, 0 ) ENDIF IF ( TRIM( var ) == 'lwp*' ) unit = 'kg/kg*m' IF ( TRIM( var ) == 'pra*' ) unit = 'mm' IF ( TRIM( var ) == 'prr*' ) unit = 'mm/s' IF ( TRIM( var ) == 'qsws*' ) unit = 'kgm/kgs' IF ( TRIM( var ) == 'shf*' ) unit = 'K*m/s' IF ( TRIM( var ) == 't*' ) unit = 'K' IF ( TRIM( var ) == 'u*' ) unit = 'm/s' IF ( TRIM( var ) == 'z0*' ) unit = 'm' IF ( TRIM( var ) == 'z0h*' ) unit = 'm' CASE ( 'p', 'pt', 'u', 'v', 'w' ) IF ( TRIM( var ) == 'p' ) unit = 'Pa' IF ( TRIM( var ) == 'pt' ) unit = 'K' IF ( TRIM( var ) == 'u' ) unit = 'm/s' IF ( TRIM( var ) == 'v' ) unit = 'm/s' IF ( TRIM( var ) == 'w' ) unit = 'm/s' CONTINUE CASE DEFAULT CALL user_check_data_output( var, unit ) IF ( unit == 'illegal' ) THEN IF ( data_output_user(1) /= ' ' ) THEN message_string = 'illegal value for data_output or ' // & 'data_output_user = "' // TRIM( data_output(i) ) // '"' CALL message( 'check_parameters', 'PA0114', 1, 2, 0, 6, 0 ) ELSE message_string = 'illegal value for data_output =' // & TRIM( data_output(i) ) // '"' CALL message( 'check_parameters', 'PA0115', 1, 2, 0, 6, 0 ) ENDIF ENDIF END SELECT ! !-- Set the internal steering parameters appropriately IF ( k == 0 ) THEN do3d_no(j) = do3d_no(j) + 1 do3d(j,do3d_no(j)) = data_output(i) do3d_unit(j,do3d_no(j)) = unit ELSE do2d_no(j) = do2d_no(j) + 1 do2d(j,do2d_no(j)) = data_output(i) do2d_unit(j,do2d_no(j)) = unit IF ( data_output(i)(ilen-2:ilen) == '_xy' ) THEN data_output_xy(j) = .TRUE. ENDIF IF ( data_output(i)(ilen-2:ilen) == '_xz' ) THEN data_output_xz(j) = .TRUE. ENDIF IF ( data_output(i)(ilen-2:ilen) == '_yz' ) THEN data_output_yz(j) = .TRUE. ENDIF ENDIF IF ( j == 1 ) THEN ! !-- Check, if variable is already subject to averaging found = .FALSE. DO k = 1, doav_n IF ( TRIM( doav(k) ) == TRIM( var ) ) found = .TRUE. ENDDO IF ( .NOT. found ) THEN doav_n = doav_n + 1 doav(doav_n) = var ENDIF ENDIF i = i + 1 ENDDO ! !-- Averaged 2d or 3d output requires that an averaging interval has been set IF ( doav_n > 0 .AND. averaging_interval == 0.0 ) THEN WRITE( message_string, * ) 'output of averaged quantity "', & TRIM( doav(1) ), '_av" requires to set a ', & 'non-zero & averaging interval' CALL message( 'check_parameters', 'PA0323', 1, 2, 0, 6, 0 ) ENDIF ! !-- Check sectional planes and store them in one shared array IF ( ANY( section_xy > nz + 1 ) ) THEN WRITE( message_string, * ) 'section_xy must be <= nz + 1 = ', nz + 1 CALL message( 'check_parameters', 'PA0319', 1, 2, 0, 6, 0 ) ENDIF IF ( ANY( section_xz > ny + 1 ) ) THEN WRITE( message_string, * ) 'section_xz must be <= ny + 1 = ', ny + 1 CALL message( 'check_parameters', 'PA0320', 1, 2, 0, 6, 0 ) ENDIF IF ( ANY( section_yz > nx + 1 ) ) THEN WRITE( message_string, * ) 'section_yz must be <= nx + 1 = ', nx + 1 CALL message( 'check_parameters', 'PA0321', 1, 2, 0, 6, 0 ) ENDIF section(:,1) = section_xy section(:,2) = section_xz section(:,3) = section_yz ! !-- Upper plot limit for 2D vertical sections IF ( z_max_do2d == -1.0 ) z_max_do2d = zu(nzt) IF ( z_max_do2d < zu(nzb+1) .OR. z_max_do2d > zu(nzt) ) THEN WRITE( message_string, * ) 'z_max_do2d = ', z_max_do2d, & ' must be >= ', zu(nzb+1), '(zu(nzb+1)) and <= ', zu(nzt), & ' (zu(nzt))' CALL message( 'check_parameters', 'PA0116', 1, 2, 0, 6, 0 ) ENDIF ! !-- Upper plot limit for 3D arrays IF ( nz_do3d == -9999 ) nz_do3d = nzt + 1 ! !-- Determine and check accuracy for compressed 3D plot output IF ( do3d_compress ) THEN ! !-- Compression only permissible on T3E machines IF ( host(1:3) /= 't3e' ) THEN message_string = 'do3d_compress = .TRUE. not allowed on host "' // & TRIM( host ) // '"' CALL message( 'check_parameters', 'PA0117', 1, 2, 0, 6, 0 ) ENDIF i = 1 DO WHILE ( do3d_comp_prec(i) /= ' ' ) ilen = LEN_TRIM( do3d_comp_prec(i) ) IF ( LLT( do3d_comp_prec(i)(ilen:ilen), '0' ) .OR. & LGT( do3d_comp_prec(i)(ilen:ilen), '9' ) ) THEN WRITE( message_string, * ) 'illegal precision: do3d_comp_prec', & '(', i, ') = "', TRIM(do3d_comp_prec(i)),'"' CALL message( 'check_parameters', 'PA0118', 1, 2, 0, 6, 0 ) ENDIF prec = IACHAR( do3d_comp_prec(i)(ilen:ilen) ) - IACHAR( '0' ) var = do3d_comp_prec(i)(1:ilen-1) SELECT CASE ( var ) CASE ( 'u' ) j = 1 CASE ( 'v' ) j = 2 CASE ( 'w' ) j = 3 CASE ( 'p' ) j = 4 CASE ( 'pt' ) j = 5 CASE DEFAULT WRITE( message_string, * ) 'unknown variable "', & TRIM( do3d_comp_prec(i) ), '" given for do3d_comp_prec(', & i, ')' CALL message( 'check_parameters', 'PA0119', 1, 2, 0, 6, 0 ) END SELECT plot_3d_precision(j)%precision = prec i = i + 1 ENDDO ENDIF ! !-- Check the data output format(s) IF ( data_output_format(1) == ' ' ) THEN ! !-- Default value netcdf_output = .TRUE. ELSE i = 1 DO WHILE ( data_output_format(i) /= ' ' ) SELECT CASE ( data_output_format(i) ) CASE ( 'netcdf' ) netcdf_output = .TRUE. CASE ( 'iso2d' ) iso2d_output = .TRUE. CASE ( 'avs' ) avs_output = .TRUE. CASE DEFAULT message_string = 'unknown value for data_output_format "' // & TRIM( data_output_format(i) ) // '"' CALL message( 'check_parameters', 'PA0120', 1, 2, 0, 6, 0 ) END SELECT i = i + 1 IF ( i > 10 ) EXIT ENDDO ENDIF ! !-- Set output format string (used in header) IF ( netcdf_output ) THEN SELECT CASE ( netcdf_data_format ) CASE ( 1 ) output_format_netcdf = 'netCDF classic' CASE ( 2 ) output_format_netcdf = 'netCDF 64bit offset' CASE ( 3 ) output_format_netcdf = 'netCDF4/HDF5' CASE ( 4 ) output_format_netcdf = 'netCDF4/HDF5 classic' CASE ( 5 ) output_format_netcdf = 'parallel netCDF4/HDF5' CASE ( 6 ) output_format_netcdf = 'parallel netCDF4/HDF5 classic' END SELECT ENDIF ! !-- Check mask conditions DO mid = 1, max_masks IF ( data_output_masks(mid,1) /= ' ' .OR. & data_output_masks_user(mid,1) /= ' ' ) THEN masks = masks + 1 ENDIF ENDDO IF ( masks < 0 .OR. masks > max_masks ) THEN WRITE( message_string, * ) 'illegal value: masks must be >= 0 and ', & '<= ', max_masks, ' (=max_masks)' CALL message( 'check_parameters', 'PA0325', 1, 2, 0, 6, 0 ) ENDIF IF ( masks > 0 ) THEN mask_scale(1) = mask_scale_x mask_scale(2) = mask_scale_y mask_scale(3) = mask_scale_z IF ( ANY( mask_scale <= 0.0 ) ) THEN WRITE( message_string, * ) & 'illegal value: mask_scale_x, mask_scale_y and mask_scale_z', & 'must be > 0.0' CALL message( 'check_parameters', 'PA0326', 1, 2, 0, 6, 0 ) ENDIF ! !-- Generate masks for masked data output CALL init_masks ENDIF ! !-- Check the NetCDF data format #if ! defined ( __check ) IF ( netcdf_data_format > 2 ) THEN #if defined( __netcdf4 ) CONTINUE #else message_string = 'netCDF: netCDF4 format requested but no ' // & 'cpp-directive __netcdf4 given & switch ' // & 'back to 64-bit offset format' CALL message( 'check_parameters', 'PA0171', 0, 1, 0, 6, 0 ) netcdf_data_format = 2 #endif ENDIF IF ( netcdf_data_format > 4 ) THEN #if defined( __netcdf4 ) && defined( __netcdf4_parallel ) CONTINUE #else message_string = 'netCDF: netCDF4 parallel output requested but no ' // & 'cpp-directive __netcdf4_parallel given & switch ' // & 'back to netCDF4 non-parallel output' CALL message( 'check_parameters', 'PA0099', 0, 1, 0, 6, 0 ) netcdf_data_format = netcdf_data_format - 2 #endif ENDIF #endif #if ! defined( __check ) ! !-- Check netcdf precison ldum = .FALSE. CALL define_netcdf_header( 'ch', ldum, 0 ) #endif ! !-- Check, whether a constant diffusion coefficient shall be used IF ( km_constant /= -1.0 ) THEN IF ( km_constant < 0.0 ) THEN WRITE( message_string, * ) 'km_constant = ', km_constant, ' < 0.0' CALL message( 'check_parameters', 'PA0121', 1, 2, 0, 6, 0 ) ELSE IF ( prandtl_number < 0.0 ) THEN WRITE( message_string, * ) 'prandtl_number = ', prandtl_number, & ' < 0.0' CALL message( 'check_parameters', 'PA0122', 1, 2, 0, 6, 0 ) ENDIF constant_diffusion = .TRUE. IF ( prandtl_layer ) THEN message_string = 'prandtl_layer is not allowed with fixed ' // & 'value of km' CALL message( 'check_parameters', 'PA0123', 1, 2, 0, 6, 0 ) ENDIF ENDIF ENDIF ! !-- In case of non-cyclic lateral boundaries and a damping layer for the !-- potential temperature, check the width of the damping layer IF ( bc_lr /= 'cyclic' ) THEN IF ( pt_damping_width < 0.0 .OR. pt_damping_width > REAL( nx * dx ) ) THEN message_string = 'pt_damping_width out of range' CALL message( 'check_parameters', 'PA0124', 1, 2, 0, 6, 0 ) ENDIF ENDIF IF ( bc_ns /= 'cyclic' ) THEN IF ( pt_damping_width < 0.0 .OR. pt_damping_width > REAL( ny * dy ) ) THEN message_string = 'pt_damping_width out of range' CALL message( 'check_parameters', 'PA0124', 1, 2, 0, 6, 0 ) ENDIF ENDIF ! !-- Check value range for rif IF ( rif_min >= rif_max ) THEN WRITE( message_string, * ) 'rif_min = ', rif_min, ' must be less ', & 'than rif_max = ', rif_max CALL message( 'check_parameters', 'PA0125', 1, 2, 0, 6, 0 ) ENDIF ! !-- Determine upper and lower hight level indices for random perturbations IF ( disturbance_level_b == -9999999.9 ) THEN IF ( ocean ) THEN disturbance_level_b = zu((nzt*2)/3) disturbance_level_ind_b = ( nzt * 2 ) / 3 ELSE disturbance_level_b = zu(nzb+3) disturbance_level_ind_b = nzb + 3 ENDIF ELSEIF ( disturbance_level_b < zu(3) ) THEN WRITE( message_string, * ) 'disturbance_level_b = ', & disturbance_level_b, ' must be >= ', zu(3), '(zu(3))' CALL message( 'check_parameters', 'PA0126', 1, 2, 0, 6, 0 ) ELSEIF ( disturbance_level_b > zu(nzt-2) ) THEN WRITE( message_string, * ) 'disturbance_level_b = ', & disturbance_level_b, ' must be <= ', zu(nzt-2), '(zu(nzt-2))' CALL message( 'check_parameters', 'PA0127', 1, 2, 0, 6, 0 ) ELSE DO k = 3, nzt-2 IF ( disturbance_level_b <= zu(k) ) THEN disturbance_level_ind_b = k EXIT ENDIF ENDDO ENDIF IF ( disturbance_level_t == -9999999.9 ) THEN IF ( ocean ) THEN disturbance_level_t = zu(nzt-3) disturbance_level_ind_t = nzt - 3 ELSE disturbance_level_t = zu(nzt/3) disturbance_level_ind_t = nzt / 3 ENDIF ELSEIF ( disturbance_level_t > zu(nzt-2) ) THEN WRITE( message_string, * ) 'disturbance_level_t = ', & disturbance_level_t, ' must be <= ', zu(nzt-2), '(zu(nzt-2))' CALL message( 'check_parameters', 'PA0128', 1, 2, 0, 6, 0 ) ELSEIF ( disturbance_level_t < disturbance_level_b ) THEN WRITE( message_string, * ) 'disturbance_level_t = ', & disturbance_level_t, ' must be >= disturbance_level_b = ', & disturbance_level_b CALL message( 'check_parameters', 'PA0129', 1, 2, 0, 6, 0 ) ELSE DO k = 3, nzt-2 IF ( disturbance_level_t <= zu(k) ) THEN disturbance_level_ind_t = k EXIT ENDIF ENDDO ENDIF ! !-- Check again whether the levels determined this way are ok. !-- Error may occur at automatic determination and too few grid points in !-- z-direction. IF ( disturbance_level_ind_t < disturbance_level_ind_b ) THEN WRITE( message_string, * ) 'disturbance_level_ind_t = ', & disturbance_level_ind_t, ' must be >= disturbance_level_b = ', & disturbance_level_b CALL message( 'check_parameters', 'PA0130', 1, 2, 0, 6, 0 ) ENDIF ! !-- Determine the horizontal index range for random perturbations. !-- In case of non-cyclic horizontal boundaries, no perturbations are imposed !-- near the inflow and the perturbation area is further limited to ...(1) !-- after the initial phase of the flow. dist_nxl = 0; dist_nxr = nx dist_nys = 0; dist_nyn = ny IF ( bc_lr /= 'cyclic' ) THEN IF ( inflow_disturbance_begin == -1 ) THEN inflow_disturbance_begin = MIN( 10, nx/2 ) ENDIF IF ( inflow_disturbance_begin < 0 .OR. inflow_disturbance_begin > nx )& THEN message_string = 'inflow_disturbance_begin out of range' CALL message( 'check_parameters', 'PA0131', 1, 2, 0, 6, 0 ) ENDIF IF ( inflow_disturbance_end == -1 ) THEN inflow_disturbance_end = MIN( 100, 3*nx/4 ) ENDIF IF ( inflow_disturbance_end < 0 .OR. inflow_disturbance_end > nx ) & THEN message_string = 'inflow_disturbance_end out of range' CALL message( 'check_parameters', 'PA0132', 1, 2, 0, 6, 0 ) ENDIF ELSEIF ( bc_ns /= 'cyclic' ) THEN IF ( inflow_disturbance_begin == -1 ) THEN inflow_disturbance_begin = MIN( 10, ny/2 ) ENDIF IF ( inflow_disturbance_begin < 0 .OR. inflow_disturbance_begin > ny )& THEN message_string = 'inflow_disturbance_begin out of range' CALL message( 'check_parameters', 'PA0131', 1, 2, 0, 6, 0 ) ENDIF IF ( inflow_disturbance_end == -1 ) THEN inflow_disturbance_end = MIN( 100, 3*ny/4 ) ENDIF IF ( inflow_disturbance_end < 0 .OR. inflow_disturbance_end > ny ) & THEN message_string = 'inflow_disturbance_end out of range' CALL message( 'check_parameters', 'PA0132', 1, 2, 0, 6, 0 ) ENDIF ENDIF IF ( bc_lr == 'radiation/dirichlet' ) THEN dist_nxr = nx - inflow_disturbance_begin dist_nxl(1) = nx - inflow_disturbance_end ELSEIF ( bc_lr == 'dirichlet/radiation' ) THEN dist_nxl = inflow_disturbance_begin dist_nxr(1) = inflow_disturbance_end ENDIF IF ( bc_ns == 'dirichlet/radiation' ) THEN dist_nyn = ny - inflow_disturbance_begin dist_nys(1) = ny - inflow_disturbance_end ELSEIF ( bc_ns == 'radiation/dirichlet' ) THEN dist_nys = inflow_disturbance_begin dist_nyn(1) = inflow_disturbance_end ENDIF ! !-- A turbulent inflow requires Dirichlet conditions at the respective inflow !-- boundary (so far, a turbulent inflow is realized from the left side only) IF ( turbulent_inflow .AND. bc_lr /= 'dirichlet/radiation' ) THEN message_string = 'turbulent_inflow = .T. requires a Dirichlet ' // & 'condition at the inflow boundary' CALL message( 'check_parameters', 'PA0133', 1, 2, 0, 6, 0 ) ENDIF ! !-- Turbulent inflow requires that 3d arrays have been cyclically filled with !-- data from prerun in the first main run IF ( turbulent_inflow .AND. initializing_actions /= 'cyclic_fill' .AND. & initializing_actions /= 'read_restart_data' ) THEN message_string = 'turbulent_inflow = .T. requires ' // & 'initializing_actions = ''cyclic_fill'' ' CALL message( 'check_parameters', 'PA0055', 1, 2, 0, 6, 0 ) ENDIF ! !-- In case of turbulent inflow calculate the index of the recycling plane IF ( turbulent_inflow ) THEN IF ( recycling_width == 9999999.9 ) THEN ! !-- Set the default value for the width of the recycling domain recycling_width = 0.1 * nx * dx ELSE IF ( recycling_width < dx .OR. recycling_width > nx * dx ) THEN WRITE( message_string, * ) 'illegal value for recycling_width:', & ' ', recycling_width CALL message( 'check_parameters', 'PA0134', 1, 2, 0, 6, 0 ) ENDIF ENDIF ! !-- Calculate the index recycling_plane = recycling_width / dx ENDIF ! !-- Check random generator IF ( random_generator /= 'system-specific' .AND. & random_generator /= 'numerical-recipes' ) THEN message_string = 'unknown random generator: random_generator = "' // & TRIM( random_generator ) // '"' CALL message( 'check_parameters', 'PA0135', 1, 2, 0, 6, 0 ) ENDIF ! !-- Determine damping level index for 1D model IF ( INDEX( initializing_actions, 'set_1d-model_profiles' ) /= 0 ) THEN IF ( damp_level_1d == -1.0 ) THEN damp_level_1d = zu(nzt+1) damp_level_ind_1d = nzt + 1 ELSEIF ( damp_level_1d < 0.0 .OR. damp_level_1d > zu(nzt+1) ) THEN WRITE( message_string, * ) 'damp_level_1d = ', damp_level_1d, & ' must be > 0.0 and < ', zu(nzt+1), '(zu(nzt+1))' CALL message( 'check_parameters', 'PA0136', 1, 2, 0, 6, 0 ) ELSE DO k = 1, nzt+1 IF ( damp_level_1d <= zu(k) ) THEN damp_level_ind_1d = k EXIT ENDIF ENDDO ENDIF ENDIF ! !-- Check some other 1d-model parameters IF ( TRIM( mixing_length_1d ) /= 'as_in_3d_model' .AND. & TRIM( mixing_length_1d ) /= 'blackadar' ) THEN message_string = 'mixing_length_1d = "' // TRIM( mixing_length_1d ) // & '" is unknown' CALL message( 'check_parameters', 'PA0137', 1, 2, 0, 6, 0 ) ENDIF IF ( TRIM( dissipation_1d ) /= 'as_in_3d_model' .AND. & TRIM( dissipation_1d ) /= 'detering' ) THEN message_string = 'dissipation_1d = "' // TRIM( dissipation_1d ) // & '" is unknown' CALL message( 'check_parameters', 'PA0138', 1, 2, 0, 6, 0 ) ENDIF ! !-- Set time for the next user defined restart (time_restart is the !-- internal parameter for steering restart events) IF ( restart_time /= 9999999.9 ) THEN IF ( restart_time > time_since_reference_point ) THEN time_restart = restart_time ENDIF ELSE ! !-- In case of a restart run, set internal parameter to default (no restart) !-- if the NAMELIST-parameter restart_time is omitted time_restart = 9999999.9 ENDIF ! !-- Set default value of the time needed to terminate a model run IF ( termination_time_needed == -1.0 ) THEN IF ( host(1:3) == 'ibm' ) THEN termination_time_needed = 300.0 ELSE termination_time_needed = 35.0 ENDIF ENDIF ! !-- Check the time needed to terminate a model run IF ( host(1:3) == 't3e' ) THEN ! !-- Time needed must be at least 30 seconds on all CRAY machines, because !-- MPP_TREMAIN gives the remaining CPU time only in steps of 30 seconds IF ( termination_time_needed <= 30.0 ) THEN WRITE( message_string, * ) 'termination_time_needed = ', & termination_time_needed, ' must be > 30.0 on host "', & TRIM( host ), '"' CALL message( 'check_parameters', 'PA0139', 1, 2, 0, 6, 0 ) ENDIF ELSEIF ( host(1:3) == 'ibm' ) THEN ! !-- On IBM-regatta machines the time should be at least 300 seconds, !-- because the job time consumed before executing palm (for compiling, !-- copying of files, etc.) has to be regarded IF ( termination_time_needed < 300.0 ) THEN WRITE( message_string, * ) 'termination_time_needed = ', & termination_time_needed, ' should be >= 300.0 on host "', & TRIM( host ), '"' CALL message( 'check_parameters', 'PA0140', 1, 2, 0, 6, 0 ) ENDIF ENDIF ! !-- Check pressure gradient conditions IF ( dp_external .AND. conserve_volume_flow ) THEN WRITE( message_string, * ) 'Both dp_external and conserve_volume_flo', & 'w are .TRUE. but one of them must be .FALSE.' CALL message( 'check_parameters', 'PA0150', 1, 2, 0, 6, 0 ) ENDIF IF ( dp_external ) THEN IF ( dp_level_b < zu(nzb) .OR. dp_level_b > zu(nzt) ) THEN WRITE( message_string, * ) 'dp_level_b = ', dp_level_b, ' is out ', & ' of range' CALL message( 'check_parameters', 'PA0151', 1, 2, 0, 6, 0 ) ENDIF IF ( .NOT. ANY( dpdxy /= 0.0 ) ) THEN WRITE( message_string, * ) 'dp_external is .TRUE. but dpdxy is ze', & 'ro, i.e. the external pressure gradient & will not be applied' CALL message( 'check_parameters', 'PA0152', 0, 1, 0, 6, 0 ) ENDIF ENDIF IF ( ANY( dpdxy /= 0.0 ) .AND. .NOT. dp_external ) THEN WRITE( message_string, * ) 'dpdxy is nonzero but dp_external is ', & '.FALSE., i.e. the external pressure gradient & will not be applied' CALL message( 'check_parameters', 'PA0153', 0, 1, 0, 6, 0 ) ENDIF IF ( conserve_volume_flow ) THEN IF ( TRIM( conserve_volume_flow_mode ) == 'default' ) THEN conserve_volume_flow_mode = 'initial_profiles' ELSEIF ( TRIM( conserve_volume_flow_mode ) /= 'initial_profiles' .AND. & TRIM( conserve_volume_flow_mode ) /= 'inflow_profile' .AND. & TRIM( conserve_volume_flow_mode ) /= 'bulk_velocity' ) THEN WRITE( message_string, * ) 'unknown conserve_volume_flow_mode: ', & conserve_volume_flow_mode CALL message( 'check_parameters', 'PA0154', 1, 2, 0, 6, 0 ) ENDIF IF ( (bc_lr /= 'cyclic' .OR. bc_ns /= 'cyclic') .AND. & TRIM( conserve_volume_flow_mode ) == 'bulk_velocity' ) THEN WRITE( message_string, * ) 'non-cyclic boundary conditions ', & 'require conserve_volume_flow_mode = ''initial_profiles''' CALL message( 'check_parameters', 'PA0155', 1, 2, 0, 6, 0 ) ENDIF IF ( bc_lr == 'cyclic' .AND. bc_ns == 'cyclic' .AND. & TRIM( conserve_volume_flow_mode ) == 'inflow_profile' ) THEN WRITE( message_string, * ) 'cyclic boundary conditions ', & 'require conserve_volume_flow_mode = ''initial_profiles''', & ' or ''bulk_velocity''' CALL message( 'check_parameters', 'PA0156', 1, 2, 0, 6, 0 ) ENDIF ENDIF IF ( ( u_bulk /= 0.0 .OR. v_bulk /= 0.0 ) .AND. & ( .NOT. conserve_volume_flow .OR. & TRIM( conserve_volume_flow_mode ) /= 'bulk_velocity' ) ) THEN WRITE( message_string, * ) 'nonzero bulk velocity requires ', & 'conserve_volume_flow = .T. and ', & 'conserve_volume_flow_mode = ''bulk_velocity''' CALL message( 'check_parameters', 'PA0157', 1, 2, 0, 6, 0 ) ENDIF ! !-- Check particle attributes IF ( particle_color /= 'none' ) THEN IF ( particle_color /= 'absuv' .AND. particle_color /= 'pt*' .AND. & particle_color /= 'z' ) THEN message_string = 'illegal value for parameter particle_color: ' // & TRIM( particle_color) CALL message( 'check_parameters', 'PA0313', 1, 2, 0, 6, 0 ) ELSE IF ( color_interval(2) <= color_interval(1) ) THEN message_string = 'color_interval(2) <= color_interval(1)' CALL message( 'check_parameters', 'PA0315', 1, 2, 0, 6, 0 ) ENDIF ENDIF ENDIF IF ( particle_dvrpsize /= 'none' ) THEN IF ( particle_dvrpsize /= 'absw' ) THEN message_string = 'illegal value for parameter particle_dvrpsize:' // & ' ' // TRIM( particle_color) CALL message( 'check_parameters', 'PA0314', 1, 2, 0, 6, 0 ) ELSE IF ( dvrpsize_interval(2) <= dvrpsize_interval(1) ) THEN message_string = 'dvrpsize_interval(2) <= dvrpsize_interval(1)' CALL message( 'check_parameters', 'PA0316', 1, 2, 0, 6, 0 ) ENDIF ENDIF ENDIF ! !-- Check &userpar parameters CALL user_check_parameters END SUBROUTINE check_parameters