source: palm/trunk/SOURCE/parin.f90 @ 4855

Last change on this file since 4855 was 4845, checked in by raasch, 3 years ago

maximum phase velocities are alwasy used for radiation boundary conditions, parameter use_cmax removed

  • Property svn:keywords set to Id
File size: 46.9 KB
RevLine 
[1682]1!> @file parin.f90
[4797]2!--------------------------------------------------------------------------------------------------!
[2696]3! This file is part of the PALM model system.
[1036]4!
[4797]5! PALM is free software: you can redistribute it and/or modify it under the terms of the GNU General
6! Public License as published by the Free Software Foundation, either version 3 of the License, or
7! (at your option) any later version.
[1036]8!
[4797]9! PALM is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the
10! implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General
11! Public License for more details.
[1036]12!
[4797]13! You should have received a copy of the GNU General Public License along with PALM. If not, see
14! <http://www.gnu.org/licenses/>.
[1036]15!
[4828]16! Copyright 1997-2021 Leibniz Universitaet Hannover
[4797]17!--------------------------------------------------------------------------------------------------!
[1036]18!
[257]19! Current revisions:
[1484]20! -----------------
[4680]21!
22!
[2233]23! Former revisions:
24! -----------------
25! $Id: parin.f90 4845 2021-01-18 11:15:37Z raasch $
[4845]26! -use_cmax
27!
28! 4842 2021-01-14 10:42:28Z raasch
[4842]29! reading of namelist file and actions in case of namelist errors revised so that statement labels
30! and goto statements are not required any more,
31! deprecated namelists removed
32!
33! 4828 2021-01-05 11:21:41Z Giersch
[4797]34! file re-formatted to follow the PALM coding standard
35!
36! 4783 2020-11-13 13:58:45Z raasch
[4783]37! bugfix for reading restart data with MPI-I/O (does not work with blockwise I/O)
38!
39! 4680 2020-09-16 10:20:34Z gronemeier
[4680]40! Add option to fix date or time of the simulation
41!
42! 4565 2020-06-15 08:30:38Z oliver.maas
[4565]43! added pt_surface_heating_rate
[4680]44!
[4565]45! 4564 2020-06-12 14:03:36Z raasch
[4564]46! Vertical nesting method of Huq et al. (2019) removed
[4680]47!
[4564]48! 4536 2020-05-17 17:24:13Z raasch
[4535]49! bugfix for restart data format query
[4680]50!
[4535]51! 4505 2020-04-20 15:37:15Z schwenkel
[4505]52! Add flag for saturation check
[4680]53!
[4505]54! 4495 2020-04-13 20:11:20Z raasch
[4495]55! restart data handling with MPI-IO added
[4680]56!
[4495]57! 4360 2020-01-07 11:25:50Z suehring
[4301]58! removed recycling_yshift
[4680]59!
[4301]60! 4227 2019-09-10 18:04:34Z gronemeier
[4227]61! implement new palm_date_time_mod
[4680]62!
[4227]63! 4146 2019-08-07 07:47:36Z gronemeier
[4196]64! added rotation_angle to initialization_parameters
[4680]65!
[4196]66! 4191 2019-08-27 15:45:07Z gronemeier
[4191]67! bugfix: add recycling_method_for_thermodynamic_quantities to inipar namelist
[4680]68!
[4191]69! 4183 2019-08-23 07:33:16Z oliver.maas
[4183]70! replaced recycle_absolute_quantities by recycling_method_for_thermodynamic_quantities
[4680]71!
[4183]72! 4182 2019-08-22 15:20:23Z scharf
[4182]73! Corrected "Former revisions" section
[4680]74!
[4182]75! 4176 2019-08-20 14:10:41Z oliver.maas
[4176]76! added recycle_absolute_quantities to initialization_parameters namelist
[4680]77!
[4176]78! 4173 2019-08-20 12:04:06Z gronemeier
[4173]79! add vdi_internal_controls
[4680]80!
[4173]81! 4131 2019-08-02 11:06:18Z monakurppa
[4131]82! Allocate hom and hom_sum to allow profile output for salsa variables.
[4680]83!
[4131]84! 4079 2019-07-09 18:04:41Z suehring
[4079]85! +monotonic_limiter_z
[4680]86!
[4079]87! 4022 2019-06-12 11:52:39Z suehring
[4797]88! Change default top boundary condition for pressure to Neumann in offline nesting case
[4680]89!
[4022]90! 4017 2019-06-06 12:16:46Z schwenkel
[3987]91! Introduce alternative switch for debug output during timestepping
[4680]92!
[3987]93! 3885 2019-04-11 11:29:34Z kanani
[4797]94! Changes related to global restructuring of location messages and introduction of additional debug
95! messages
[4680]96!
[3885]97! 3806 2019-03-21 12:45:50Z raasch
[3806]98! additional check for lateral boundary conditions added
[4680]99!
[3806]100! 3747 2019-02-16 15:15:23Z gronemeier
[3747]101! removed setting of parameter region
[4680]102!
[3747]103! 3746 2019-02-16 12:41:27Z gronemeier
[3668]104! Removed most_method
[4680]105!
[3668]106! 3649 2019-01-02 16:52:21Z suehring
[3649]107! Delete debug-print statements
[3569]108!
[4182]109! Revision 1.1  1997/07/24 11:22:50  raasch
110! Initial revision
111!
112!
[1]113! Description:
114! ------------
[1682]115!> This subroutine reads variables controling the run from the NAMELIST files
[3298]116!>
117!> @todo: Revise max_pr_cs (profiles for chemistry)
[4797]118!--------------------------------------------------------------------------------------------------!
[1682]119 SUBROUTINE parin
[1]120
[4680]121
[4797]122    USE arrays_3d,                                                                                 &
123        ONLY:  pt_init, q_init, ref_state, s_init, sa_init, ug, u_init, v_init, vg
[1320]124
[2696]125    USE chem_modules
[1320]126
[1762]127    USE control_parameters
[1320]128
[4797]129    USE cpulog,                                                                                    &
[1320]130        ONLY:  cpu_log_barrierwait
131
[4797]132    USE grid_variables,                                                                            &
[1320]133        ONLY:  dx, dy
134
[4797]135    USE indices,                                                                                   &
[1320]136        ONLY:  nx, ny, nz
137
[1764]138    USE kinds
139
[4797]140    USE model_1d_mod,                                                                              &
[1320]141        ONLY:  damp_level_1d, dt_pr_1d, dt_run_control_1d, end_time_1d
142
[4797]143    USE module_interface,                                                                          &
[3637]144        ONLY:  module_interface_parin
[3159]145
[4797]146    USE netcdf_interface,                                                                          &
[1783]147        ONLY:  netcdf_data_format, netcdf_deflate, netcdf_precision
148
[2696]149    USE pegrid
[3448]150
[4797]151    USE pmc_interface,                                                                             &
[1933]152        ONLY:  nested_run, nesting_mode
[1764]153
[4797]154    USE profil_parameter,                                                                          &
[2298]155        ONLY:  cross_profiles, profile_columns, profile_rows
[1320]156
[4797]157    USE progress_bar,                                                                              &
[3313]158        ONLY :  progress_bar_disabled
[1402]159
[4797]160    USE read_restart_data_mod,                                                                     &
[3637]161        ONLY:  rrd_global
[2894]162
[4797]163    USE statistics,                                                                                &
[3746]164        ONLY:  hom, hom_sum, pr_palm, statistic_regions
[1320]165
[4797]166    USE turbulence_closure_mod,                                                                    &
[3083]167        ONLY:  rans_const_c, rans_const_sigma
168
[1914]169
[1]170    IMPLICIT NONE
171
[4842]172    CHARACTER(LEN=100) ::  line  !< dummy string that contains the current line of the parameter
173                                 !< file
[2849]174
[2600]175    INTEGER(iwp) ::  global_id      !< process id with respect to MPI_COMM_WORLD
176    INTEGER(iwp) ::  global_procs   !< # of procs with respect to MPI_COMM_WORLD
177    INTEGER(iwp) ::  i              !<
[4842]178    INTEGER(iwp) ::  io_status      !< status after reading the namelist files
[1]179
[2932]180
[4797]181    NAMELIST /initialization_parameters/  alpha_surface,                                           &
182                                          approximation,                                           &
183                                          bc_e_b,                                                  &
184                                          bc_lr,                                                   &
185                                          bc_ns,                                                   &
186                                          bc_p_b,                                                  &
187                                          bc_p_t,                                                  &
188                                          bc_pt_b,                                                 &
189                                          bc_pt_t,                                                 &
190                                          bc_q_b,                                                  &
191                                          bc_q_t,                                                  &
192                                          bc_s_b,                                                  &
193                                          bc_s_t,                                                  &
194                                          bc_uv_b,                                                 &
195                                          bc_uv_t,                                                 &
196                                          building_height,                                         &
197                                          building_length_x,                                       &
198                                          building_length_y,                                       &
199                                          building_wall_left,                                      &
200                                          building_wall_south,                                     &
201                                          calc_soil_moisture_during_spinup,                        &
202                                          call_psolver_at_all_substeps,                            &
203                                          canyon_height,                                           &
204                                          canyon_wall_left,                                        &
205                                          canyon_wall_south,                                       &
206                                          canyon_width_x,                                          &
207                                          canyon_width_y,                                          &
208                                          cfl_factor,                                              &
209                                          check_realistic_q,                                       &
210                                          cloud_droplets,                                          &
211                                          collective_wait,                                         &
212                                          complex_terrain,                                         &
213                                          conserve_volume_flow,                                    &
214                                          conserve_volume_flow_mode,                               &
215                                          constant_flux_layer,                                     &
216                                          coupling_start_time,                                     &
217                                          cycle_mg,                                                &
218                                          damp_level_1d,                                           &
219                                          data_output_during_spinup,                               &
220                                          dissipation_1d,                                          &
221                                          dp_external,                                             &
222                                          dp_level_b,                                              &
223                                          dp_smooth, dpdxy,                                        &
224                                          dt,                                                      &
225                                          dt_pr_1d,                                                &
226                                          dt_run_control_1d,                                       &
227                                          dt_spinup,                                               &
228                                          dx,                                                      &
229                                          dy,                                                      &
230                                          dz,                                                      &
231                                          dz_max,                                                  &
232                                          dz_stretch_factor,                                       &
233                                          dz_stretch_level,                                        &
234                                          dz_stretch_level_start,                                  &
235                                          dz_stretch_level_end,                                    &
236                                          e_init,                                                  &
237                                          e_min,                                                   &
238                                          end_time_1d,                                             &
239                                          ensemble_member_nr,                                      &
240                                          fft_method,                                              &
241                                          flux_input_mode,                                         &
242                                          flux_output_mode,                                        &
243                                          galilei_transformation,                                  &
244                                          humidity,                                                &
245                                          inflow_damping_height,                                   &
246                                          inflow_damping_width,                                    &
247                                          inflow_disturbance_begin,                                &
248                                          inflow_disturbance_end,                                  &
249                                          initializing_actions,                                    &
250                                          km_constant,                                             &
251                                          large_scale_forcing,                                     &
252                                          large_scale_subsidence,                                  &
253                                          latitude,                                                &
254                                          longitude,                                               &
255                                          loop_optimization,                                       &
256                                          lsf_exception,                                           &
257                                          masking_method,                                          &
258                                          mg_cycles,                                               &
259                                          mg_switch_to_pe0_level,                                  &
260                                          mixing_length_1d,                                        &
261                                          momentum_advec,                                          &
262                                          monotonic_limiter_z,                                     &
263                                          netcdf_precision,                                        &
264                                          neutral,                                                 &
265                                          ngsrb,                                                   &
266                                          nsor,                                                    &
267                                          nsor_ini,                                                &
268                                          nudging,                                                 &
269                                          nx,                                                      &
270                                          ny,                                                      &
271                                          nz,                                                      &
272                                          ocean_mode,                                              &
273                                          omega,                                                   &
274                                          omega_sor,                                               &
275                                          origin_date_time,                                        &
276                                          outflow_source_plane,                                    &
277                                          passive_scalar,                                          &
278                                          prandtl_number,                                          &
279                                          psolver,                                                 &
280                                          pt_damping_factor,                                       &
281                                          pt_damping_width,                                        &
282                                          pt_reference, pt_surface,                                &
283                                          pt_surface_heating_rate,                                 &
284                                          pt_surface_initial_change,                               &
285                                          pt_vertical_gradient,                                    &
286                                          pt_vertical_gradient_level,                              &
287                                          q_surface,                                               &
288                                          q_surface_initial_change,                                &
289                                          q_vertical_gradient,                                     &
290                                          q_vertical_gradient_level,                               &
291                                          random_generator,                                        &
292                                          random_heatflux,                                         &
293                                          rans_const_c,                                            &
294                                          rans_const_sigma,                                        &
295                                          rayleigh_damping_factor,                                 &
296                                          rayleigh_damping_height,                                 &
297                                          recycling_method_for_thermodynamic_quantities,           &
298                                          recycling_width,                                         &
299                                          reference_state,                                         &
300                                          residual_limit,                                          &
301                                          restart_data_format,                                     &
302                                          restart_data_format_input,                               &
303                                          restart_data_format_output,                              &
304                                          rotation_angle,                                          &
305                                          roughness_length,                                        &
306                                          s_surface,                                               &
307                                          s_surface_initial_change,                                &
308                                          s_vertical_gradient,                                     &
309                                          s_vertical_gradient_level,                               &
310                                          scalar_advec,                                            &
311                                          scalar_rayleigh_damping,                                 &
312                                          spinup_time,                                             &
313                                          spinup_pt_amplitude,                                     &
314                                          spinup_pt_mean,                                          &
315                                          statistic_regions,                                       &
316                                          subs_vertical_gradient,                                  &
317                                          subs_vertical_gradient_level,                            &
318                                          surface_heatflux,                                        &
319                                          surface_pressure,                                        &
320                                          surface_scalarflux,                                      &
321                                          surface_waterflux,                                       &
322                                          timestep_scheme,                                         &
323                                          topography,                                              &
324                                          topography_grid_convention,                              &
325                                          top_heatflux,                                            &
326                                          top_momentumflux_u,                                      &
327                                          top_momentumflux_v,                                      &
328                                          top_scalarflux,                                          &
329                                          transpose_compute_overlap,                               &
330                                          tunnel_height,                                           &
331                                          tunnel_length,                                           &
332                                          tunnel_wall_depth,                                       &
333                                          tunnel_width_x,                                          &
334                                          tunnel_width_y,                                          &
335                                          turbulence_closure,                                      &
336                                          turbulent_inflow,                                        &
337                                          turbulent_outflow,                                       &
338                                          u_bulk,                                                  &
339                                          u_profile,                                               &
340                                          ug_surface,                                              &
341                                          ug_vertical_gradient,                                    &
342                                          ug_vertical_gradient_level,                              &
343                                          use_fixed_date,                                          &
344                                          use_fixed_time,                                          &
345                                          use_free_convection_scaling,                             &
346                                          use_ug_for_galilei_tr,                                   &
347                                          use_subsidence_tendencies,                               &
[4845]348                                          use_surface_fluxes,                                      &
[4797]349                                          use_top_fluxes,                                          &
350                                          use_upstream_for_tke,                                    &
351                                          uv_heights,                                              &
352                                          v_bulk,                                                  &
353                                          v_profile,                                               &
354                                          vdi_checks,                                              &
355                                          vg_surface,                                              &
356                                          vg_vertical_gradient,                                    &
357                                          vg_vertical_gradient_level,                              &
358                                          wall_adjustment,                                         &
359                                          wall_heatflux,                                           &
360                                          wall_humidityflux,                                       &
361                                          wall_scalarflux,                                         &
362                                          y_shift,                                                 &
363                                          zeta_max,                                                &
364                                          zeta_min,                                                &
365                                          z0h_factor
[4680]366
[4797]367    NAMELIST /runtime_parameters/  averaging_interval,                                             &
368                                   averaging_interval_pr,                                          &
369                                   cpu_log_barrierwait,                                            &
370                                   create_disturbances,                                            &
371                                   cross_profiles,                                                 &
372                                   data_output,                                                    &
373                                   data_output_2d_on_each_pe,                                      &
374                                   data_output_masks,                                              &
375                                   data_output_pr,                                                 &
376                                   debug_output,                                                   &
377                                   debug_output_timestep,                                          &
378                                   disturbance_amplitude,                                          &
379                                   disturbance_energy_limit,                                       &
380                                   disturbance_level_b,                                            &
381                                   disturbance_level_t,                                            &
382                                   do2d_at_begin,                                                  &
383                                   do3d_at_begin,                                                  &
384                                   dt,                                                             &
385                                   dt_averaging_input,                                             &
386                                   dt_averaging_input_pr,                                          &
387                                   dt_coupling,                                                    &
388                                   dt_data_output,                                                 &
389                                   dt_data_output_av,                                              &
390                                   dt_disturb,                                                     &
391                                   dt_domask,                                                      &
392                                   dt_dopr,                                                        &
393                                   dt_dopr_listing,                                                &
394                                   dt_dots,                                                        &
395                                   dt_do2d_xy,                                                     &
396                                   dt_do2d_xz,                                                     &
397                                   dt_do2d_yz,                                                     &
398                                   dt_do3d,                                                        &
399                                   dt_max,                                                         &
400                                   dt_restart,                                                     &
401                                   dt_run_control,                                                 &
402                                   end_time,                                                       &
403                                   force_print_header,                                             &
404                                   mask_k_over_surface,                                            &
405                                   mask_scale_x,                                                   &
406                                   mask_scale_y,                                                   &
407                                   mask_scale_z,                                                   &
408                                   mask_x,                                                         &
409                                   mask_y,                                                         &
410                                   mask_z,                                                         &
411                                   mask_x_loop,                                                    &
412                                   mask_y_loop,                                                    &
413                                   mask_z_loop,                                                    &
414                                   netcdf_data_format,                                             &
415                                   netcdf_deflate,                                                 &
416                                   normalizing_region,                                             &
417                                   npex,                                                           &
418                                   npey,                                                           &
419                                   nz_do3d,                                                        &
420                                   profile_columns,                                                &
421                                   profile_rows,                                                   &
422                                   restart_time,                                                   &
423                                   section_xy,                                                     &
424                                   section_xz,                                                     &
425                                   section_yz,                                                     &
426                                   restart_data_format,                                            &
427                                   restart_data_format_input,                                      &
428                                   restart_data_format_output,                                     &
429                                   skip_time_data_output,                                          &
430                                   skip_time_data_output_av,                                       &
431                                   skip_time_dopr,                                                 &
432                                   skip_time_do2d_xy,                                              &
433                                   skip_time_do2d_xz,                                              &
434                                   skip_time_do2d_yz,                                              &
435                                   skip_time_do3d,                                                 &
436                                   skip_time_domask,                                               &
437                                   synchronous_exchange,                                           &
438                                   termination_time_needed
[1]439
[4797]440    NAMELIST /envpar/  host,                                                                       &
441                       maximum_cpu_time_allowed,                                                   &
442                       maximum_parallel_io_streams,                                                &
443                       progress_bar_disabled,                                                      &
444                       read_svf,                                                                   &
445                       revision,                                                                   &
446                       run_identifier,                                                             &
447                       tasks_per_node,                                                             &
448                       write_binary,                                                               &
449                       write_svf
[1]450
451!
[4797]452!-- First read values of environment variables (this NAMELIST file is generated by palmrun)
[3885]453    CALL location_message( 'reading environment parameters from ENVPAR', 'start' )
[2310]454
[4842]455    OPEN( 90, FILE='ENVPAR', STATUS='OLD', FORM='FORMATTED', IOSTAT=io_status )
[2310]456
[4842]457    IF ( io_status /= 0 )  THEN
[4797]458       message_string = 'local file ENVPAR not found' //                                           &
[3046]459                        '&some variables for steering may not be properly set'
[2310]460       CALL message( 'parin', 'PA0276', 0, 1, 0, 6, 0 )
461    ELSE
[4842]462       READ( 90, envpar, IOSTAT=io_status )
463       IF ( io_status < 0 )  THEN
[4797]464          message_string = 'no envpar-NAMELIST found in local file '  //                           &
[4842]465                           'ENVPAR& or some variables for steering may not be properly set'
[2310]466          CALL message( 'parin', 'PA0278', 0, 1, 0, 6, 0 )
[4842]467       ELSEIF ( io_status > 0 )  THEN
[4797]468          message_string = 'errors in local file ENVPAR' //                                        &
[3046]469                           '&some variables for steering may not be properly set'
[2310]470          CALL message( 'parin', 'PA0277', 0, 1, 0, 6, 0 )
471       ENDIF
[4797]472       CLOSE( 90 )
[2310]473    ENDIF
474
[3885]475    CALL location_message( 'reading environment parameters from ENVPAR', 'finished' )
[4783]476
477    CALL location_message( 'reading NAMELIST parameters from PARIN', 'start' )
478
[1]479!
[4783]480!-- Open the NAMELIST-file which is send with this job
481    CALL check_open( 11 )
482
483!
484!-- Read the control parameters for initialization.
[4842]485!-- The namelist "initialisation_parameters" must be provided in the NAMELIST-file.
486    READ( 11, initialization_parameters, IOSTAT=io_status )
[4783]487
[4842]488!
489!-- Action depending on the READ status
490    IF ( io_status > 0 )  THEN
491!
492!--    initialisation_parameters namelist was found but countained errors. Print an error message
493!-     including the line that caused the problem.
494       BACKSPACE( 11 )
495       READ( 11 , '(A)' ) line
496       CALL parin_fail_message( 'initialization_parameters', line )
[4783]497
[4842]498    ELSEIF ( io_status < 0 )  THEN
499!
500!--    initialisation_parametes namelist was not found. Return a message.
501       message_string = 'no initialization_parameters-namelist found'
502       CALL message( 'parin', 'PA0272', 1, 2, 0, 6, 0 )
[4783]503
[4842]504    ENDIF
[4783]505!
[4797]506!-- Try to read runtime parameters given by the user for this run (namelist "runtime_parameters").
507!-- The namelist "runtime_parmeters" can be omitted. In that case default values are used for the
[4783]508!-- parameters.
[4797]509    REWIND( 11 )
[4842]510    READ( 11, runtime_parameters, IOSTAT=io_status )
[4783]511
[4842]512    IF ( io_status > 0 )  THEN
[4783]513!
[4842]514!--    Namelist runtime_parameters was found but contained errors. Print an error message including
515!--    the line that caused the problem.
516       BACKSPACE( 11 )
517       READ( 11 , '(A)') line
518       CALL parin_fail_message( 'runtime_parameters', line )
[4783]519
[4842]520    ENDIF
[4783]521
522!
523!-- Check for module namelists and read them
524    CALL module_interface_parin
525
526!
[759]527!-- Calculate the number of groups into which parallel I/O is split.
[4797]528!-- The default for files which are opened by all PEs (or where each PE opens its own independent
529!-- file) is, that all PEs are doing input/output in parallel at the same time. This might cause
530!-- performance or even more severe problems depending on the configuration of the underlying file
[759]531!-- system.
[4797]532!-- Calculation of the number of blocks and the I/O group must be based on all PEs involved in this
533!-- run. Since myid and numprocs are related to the comm2d communicator, which gives only a subset
534!-- of all PEs in case of nested runs, that information must be inquired again from the global
[2600]535!-- communicator.
[759]536!-- First, set the default:
[2967]537#if defined( __parallel )
[2600]538    CALL MPI_COMM_RANK( MPI_COMM_WORLD, global_id, ierr )
539    CALL MPI_COMM_SIZE( MPI_COMM_WORLD, global_procs, ierr )
[2967]540#else
541    global_id    = 0
542    global_procs = 1
543#endif
[4797]544    IF ( maximum_parallel_io_streams == -1  .OR.  maximum_parallel_io_streams > global_procs )  THEN
[2600]545       maximum_parallel_io_streams = global_procs
[759]546    ENDIF
547!
[4797]548!-- Now calculate the number of io_blocks and the io_group to which the respective PE belongs. I/O
549!-- of the groups is done in serial, but in parallel for all PEs belonging to the same group.
[2600]550    io_blocks = global_procs / maximum_parallel_io_streams
551    io_group  = MOD( global_id+1, io_blocks )
[4680]552
[759]553!
[4797]554!-- If required, read control parameters from restart file (produced by a prior run). All PEs are
555!-- reading from file created by PE0 (see check_open)
[4783]556    IF ( TRIM( initializing_actions ) == 'read_restart_data' )  THEN
[559]557
[1]558!
[4797]559!--    If not set by the user in the runtime parameters, the data format for restart input needs to
560!--    be set now! This is normally done later in check parameters.
561       IF ( TRIM( restart_data_format ) /= 'fortran_binary'  .AND.                                 &
562            TRIM( restart_data_format ) /= 'mpi'             .AND.                                 &
[4783]563            TRIM( restart_data_format ) /= 'mpi_shared_memory' )  THEN
564          message_string = 'illegal restart data format "' // TRIM( restart_data_format ) // '"'
565          CALL message( 'parin', 'PA0724', 1, 2, 0, 6, 0 )
566       ENDIF
567       IF ( TRIM( restart_data_format_input ) == 'undefined' )  THEN
568          restart_data_format_input = restart_data_format
569       ENDIF
[559]570
[1]571!
[4783]572!--    Blockwise I/O does not work together with MPI-I/O
573       IF ( restart_data_format_input(1:3) == 'mpi'  .AND.  io_blocks /= 1 )  THEN
574          CALL rrd_global
575       ELSE
[146]576!
[4783]577!--       Data is read in parallel by groups of PEs
578          DO  i = 0, io_blocks-1
579             IF ( i == io_group )  THEN
580                CALL rrd_global
581             ENDIF
582#if defined( __parallel )
583             CALL MPI_BARRIER( comm2d, ierr )
584#endif
[2894]585          ENDDO
[4783]586       ENDIF
[2894]587
[2932]588
[3246]589!
[4783]590!--    Increment the run count
591       runnr = runnr + 1
[2894]592!
[4797]593!--    In case of a restart run, the number of user-defined profiles on the restart file (already
594!--    stored in max_pr_user) has to match the one given for the current run. max_pr_user_tmp is
595!--    calculated in user_parin and max_pr_user is read in via rrd_global.
[4783]596       IF ( max_pr_user /= max_pr_user_tmp )  THEN
[4797]597          WRITE( message_string, * ) 'the number of user-defined ',                                &
598                                     'profiles given in data_output_pr (', max_pr_user_tmp,        &
599                                     ') does not match the one ',                                  &
600                                     'found in the restart file (', max_pr_user, ')'
[4783]601          CALL message( 'user_parin', 'UI0009', 1, 2, 0, 6, 0 )
602       ENDIF
603    ELSE
604       max_pr_user = max_pr_user_tmp
605    ENDIF
[2894]606!
[4783]607!-- Activate spinup
[4797]608    IF ( land_surface  .OR.  urban_surface )  THEN
[4783]609       IF ( spinup_time > 0.0_wp )  THEN
610          coupling_start_time = spinup_time
611          time_since_reference_point = simulated_time - coupling_start_time
612          IF ( spinup_pt_mean == 9999999.9_wp )  THEN
613             spinup_pt_mean = pt_surface
[759]614          ENDIF
[4783]615          end_time = end_time + spinup_time
616          IF ( TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
617             spinup = .TRUE.
[4797]618          ELSEIF ( TRIM( initializing_actions ) == 'read_restart_data'  .AND.                      &
619                   time_since_reference_point > 0.0_wp )  THEN
[4783]620             data_output_during_spinup = .FALSE.  !< required for correct ntdim calculation
621                                                  !< in check_parameters for restart run
[2921]622          ENDIF
[4783]623       ENDIF
624    ENDIF
[2921]625
626!
[4783]627!-- In case of nested runs, explicitly set nesting boundary conditions.
628!-- This will overwrite the user settings and basic defaults.
629!-- bc_lr and bc_ns always need to be cyclic for vertical nesting.
630    IF ( nested_run )  THEN
631       IF ( nesting_mode == 'vertical' )  THEN
[4797]632          IF ( bc_lr /= 'cyclic'  .OR.  bc_ns /= 'cyclic' )  THEN
633             WRITE ( message_string, *) 'bc_lr and bc_ns were set to ,',                           &
634                                        'cyclic for vertical nesting'
[4783]635             CALL message( 'parin', 'PA0428', 0, 0, 0, 6, 0 )
636             bc_lr   = 'cyclic'
637             bc_ns   = 'cyclic'
[1764]638          ENDIF
[4783]639          IF ( child_domain )  THEN
640             bc_uv_t  = 'nested'
641             bc_pt_t  = 'nested'
642             bc_q_t   = 'nested'
643             bc_s_t   = 'nested'
644             bc_cs_t  = 'nested'
645             bc_p_t   = 'neumann'
646          ENDIF
[3182]647!
[4797]648!--    For other nesting modes only set boundary conditions for nested domains.
[4783]649       ELSE
650          IF ( child_domain )  THEN
651             bc_lr    = 'nested'
652             bc_ns    = 'nested'
653             bc_uv_t  = 'nested'
654             bc_pt_t  = 'nested'
655             bc_q_t   = 'nested'
656             bc_s_t   = 'nested'
657             bc_cs_t  = 'nested'
[4022]658             bc_p_t   = 'neumann'
[2696]659          ENDIF
[4783]660       ENDIF
661    ENDIF
662!
[4797]663!-- Set boundary conditions also in case the model is offline-nested in larger-scale models.
[4783]664    IF ( nesting_offline )  THEN
665       bc_lr    = 'nesting_offline'
666       bc_ns    = 'nesting_offline'
667       bc_uv_t  = 'nesting_offline'
668       bc_pt_t  = 'nesting_offline'
669       bc_q_t   = 'nesting_offline'
670     !  bc_s_t   = 'nesting_offline'  ! scalar boundary condition is not clear yet
671     !  bc_cs_t  = 'nesting_offline'  ! same for chemical species
672       bc_p_t   = 'neumann'
673    ENDIF
[2696]674
[4680]675!
[4797]676!-- In case of nested runs, make sure that initializing_actions = 'set_constant_profiles' even
677!-- though the constant-profiles initializations for the prognostic variables will be overwritten by
678!-- pmci_child_initialize and pmci_parent_initialize. This is, however, important e.g. to make sure
679!-- that diagnostic variables are set properly. An exception is made in case of restart runs and if
680!-- user decides to do everything by its own.
681    IF ( child_domain  .AND.  .NOT. ( TRIM( initializing_actions ) == 'read_restart_data'     .OR. &
682                                      TRIM( initializing_actions ) == 'set_constant_profiles' .OR. &
683                                      TRIM( initializing_actions ) == 'by_user' )                  &
684                                    )  THEN
685       message_string = 'initializing_actions = ' // TRIM( initializing_actions ) //               &
686                        ' has been ' // 'changed to set_constant_profiles in child ' // 'domain.'
[4783]687       CALL message( 'parin', 'PA0492', 0, 0, 0, 6, 0 )
[2975]688
[4783]689       initializing_actions = 'set_constant_profiles'
690    ENDIF
[1764]691!
[4797]692!-- Check validity of lateral boundary conditions. This has to be done here because they are already
693!-- used in init_pegrid and init_grid and therefore cannot be check in check_parameters
694    IF ( bc_lr /= 'cyclic'               .AND.  bc_lr /= 'dirichlet/radiation'  .AND.              &
695         bc_lr /= 'radiation/dirichlet'  .AND.  bc_lr /= 'nested'               .AND.              &
[4783]696         bc_lr /= 'nesting_offline' )  THEN
[4797]697       message_string = 'unknown boundary condition: bc_lr = "' // TRIM( bc_lr ) // '"'
[4783]698       CALL message( 'parin', 'PA0049', 1, 2, 0, 6, 0 )
699    ENDIF
[4797]700    IF ( bc_ns /= 'cyclic'               .AND.  bc_ns /= 'dirichlet/radiation'  .AND.              &
701         bc_ns /= 'radiation/dirichlet'  .AND.  bc_ns /= 'nested'               .AND.              &
[4783]702         bc_ns /= 'nesting_offline' )  THEN
[4797]703       message_string = 'unknown boundary condition: bc_ns = "' // TRIM( bc_ns ) // '"'
[4783]704       CALL message( 'parin', 'PA0050', 1, 2, 0, 6, 0 )
705    ENDIF
[1764]706!
[4783]707!-- Set internal variables used for speed optimization in if clauses
708    IF ( bc_lr /= 'cyclic' )               bc_lr_cyc    = .FALSE.
709    IF ( bc_lr == 'dirichlet/radiation' )  bc_lr_dirrad = .TRUE.
710    IF ( bc_lr == 'radiation/dirichlet' )  bc_lr_raddir = .TRUE.
711    IF ( bc_ns /= 'cyclic' )               bc_ns_cyc    = .FALSE.
712    IF ( bc_ns == 'dirichlet/radiation' )  bc_ns_dirrad = .TRUE.
713    IF ( bc_ns == 'radiation/dirichlet' )  bc_ns_raddir = .TRUE.
[1764]714!
[4783]715!-- Radiation-Dirichlet conditions are allowed along one of the horizontal directions only.
[4797]716!-- In general, such conditions along x and y may work, but require a) some code changes (e.g.
717!-- concerning allocation of c_u, c_v ... arrays), and b) a careful model setup by the user, in
718!-- order to guarantee that there is a clearly defined outflow at two sides.
[4783]719!-- Otherwise, the radiation condition may produce wrong results.
720    IF ( ( bc_lr_dirrad .OR. bc_lr_raddir )  .AND.  ( bc_ns_dirrad .OR. bc_ns_raddir ) )  THEN
[4797]721       message_string = 'bc_lr = "' // TRIM( bc_lr ) // '" and bc_ns = "' // TRIM( bc_ns ) //      &
722                        '" are not allowed to be set at the same time'
[4783]723       CALL message( 'parin', 'PA0589', 1, 2, 0, 6, 0 )
724    ENDIF
[3806]725!
[4797]726!-- Check in case of initial run, if the grid point numbers are well defined and allocate some
727!-- arrays which are already needed in init_pegrid or check_parameters. During restart jobs, these
728!-- arrays will be allocated in rrd_global. All other arrays are allocated in init_3d_model.
[4783]729    IF ( TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
[667]730
[4783]731       IF ( nx <= 0 )  THEN
[4797]732          WRITE( message_string, * ) 'no value or wrong value given', ' for nx: nx=', nx
[4783]733          CALL message( 'parin', 'PA0273', 1, 2, 0, 6, 0 )
734       ENDIF
735       IF ( ny <= 0 )  THEN
[4797]736          WRITE( message_string, * ) 'no value or wrong value given', ' for ny: ny=', ny
[4783]737          CALL message( 'parin', 'PA0274', 1, 2, 0, 6, 0 )
738       ENDIF
739       IF ( nz <= 0 )  THEN
[4797]740          WRITE( message_string, * ) 'no value or wrong value given', ' for nz: nz=', nz
[4783]741          CALL message( 'parin', 'PA0275', 1, 2, 0, 6, 0 )
742       ENDIF
[3204]743
[759]744!
[4797]745!--    As a side condition, routine flow_statistics require at least 14 vertical grid levels (they
746!--    are used to store time-series data)
[4783]747!>     @todo   Remove this restriction
748       IF ( nz < 14 )  THEN
749          WRITE( message_string, * ) 'nz >= 14 is required'
750          CALL message( 'parin', 'PA0362', 1, 2, 0, 6, 0 )
751       ENDIF
[3204]752
753!
[4797]754!--    ATTENTION: in case of changes to the following statement please also check the allocate
755!--               statement in routine rrd_global
756       ALLOCATE( pt_init(0:nz+1), q_init(0:nz+1), s_init(0:nz+1),                                  &
757                 ref_state(0:nz+1), sa_init(0:nz+1), ug(0:nz+1),                                   &
758                 u_init(0:nz+1), v_init(0:nz+1), vg(0:nz+1),                                       &
759                 hom(0:nz+1,2,pr_palm+max_pr_user+max_pr_cs+max_pr_salsa,0:statistic_regions),     &
[4783]760                 hom_sum(0:nz+1,pr_palm+max_pr_user+max_pr_cs+max_pr_salsa,0:statistic_regions) )
[1]761
[4783]762       hom = 0.0_wp
[1]763
[4783]764    ENDIF
[759]765
[1]766!
[4783]767!-- NAMELIST-file is not needed anymore
768    CALL close_file( 11 )
[1]769
[3885]770    CALL location_message( 'reading NAMELIST parameters from PARIN', 'finished' )
[1384]771
[1]772 END SUBROUTINE parin
Note: See TracBrowser for help on using the repository browser.