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

Last change on this file since 2599 was 2599, checked in by hellstea, 6 years ago

i/o grouping update for nested runs

  • Property svn:keywords set to Id
File size: 30.5 KB
Line 
1!> @file parin.f90
2!------------------------------------------------------------------------------!
3! This file is part of PALM.
4!
5! PALM is free software: you can redistribute it and/or modify it under the
6! terms of the GNU General Public License as published by the Free Software
7! Foundation, either version 3 of the License, or (at your option) any later
8! version.
9!
10! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
11! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
12! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
13!
14! You should have received a copy of the GNU General Public License along with
15! PALM. If not, see <http://www.gnu.org/licenses/>.
16!
17! Copyright 1997-2017 Leibniz Universitaet Hannover
18!------------------------------------------------------------------------------!
19!
20! Current revisions:
21! -----------------
22!
23!
24! Former revisions:
25! -----------------
26! $Id: parin.f90 2599 2017-11-01 13:18:45Z hellstea $
27! The i/o grouping is updated to work correctly also in nested runs.
28!
29! 2575 2017-10-24 09:57:58Z maronga
30! Renamed phi -> latitude, added longitude
31!
32! 2563 2017-10-19 15:36:10Z Giersch
33! Changed position where restart files are closed.
34!
35! 2550 2017-10-16 17:12:01Z boeske
36! Added complex_terrain
37!
38! 2544 2017-10-13 18:09:32Z maronga
39! Moved day_of_year_init and time_utc_init to inipar.
40!
41! 2397 2017-09-04 16:22:48Z suehring
42! Enable initialization of 3d model by user in the child domain.
43!
44! 2375 2017-08-29 14:10:28Z schwenkel
45! Added aerosol initialization for bulk microphysics
46!
47! 2372 2017-08-25 12:37:32Z sward
48! y_shift added to namelist
49!
50! 2365 2017-08-21 14:59:59Z kanani
51! Vertical grid nesting: add vnest_start_time to d3par (SadiqHuq)
52!
53! 2339 2017-08-07 13:55:26Z gronemeier
54! corrected timestamp in header
55!
56! 2338 2017-08-07 12:15:38Z gronemeier
57! Modularize 1D model
58!
59! 2310 2017-07-11 09:37:02Z gronemeier
60! Bugfix: re-arranged call for error messages for ENVPAR file
61!
62! 2304 2017-07-04 14:35:55Z suehring
63! Bugfix, enable restarts for child domain.
64!
65! 2298 2017-06-29 09:28:18Z raasch
66! -return_addres, return_username in ENVPAR, -cross_ts_uymax, cross_ts_uymin in
67! d3par
68!
69! 2296 2017-06-28 07:53:56Z maronga
70! Added parameters for model spinup
71!
72! 2292 2017-06-20 09:51:42Z schwenkel
73! Implementation of new microphysic scheme: cloud_scheme = 'morrison'
74! includes two more prognostic equations for cloud drop concentration (nc) 
75! and cloud water content (qc).
76!
77! 2267 2017-06-09 09:33:25Z gronemeier
78! Bugfix: removed skipping of reading namelists in case of omitted d3par
79!
80! 2259 2017-06-08 09:09:11Z gronemeier
81! Implemented synthetic turbulence generator
82!
83! 2233 2017-05-30 18:08:54Z suehring
84!
85! 2232 2017-05-30 17:47:52Z suehring
86! typo corrected
87! +wall_salinityflux
88! +tunnel_height, tunnel_lenght, tunnel_width_x, tunnel_width_y,
89!  tunnel_wall_depth
90!
91! 2118 2017-01-17 16:38:49Z raasch
92! -background_communication from inipar
93!
94! 2050 2016-11-08 15:00:55Z gronemeier
95! Implement turbulent outflow condition
96!
97! 2037 2016-10-26 11:15:40Z knoop
98! Anelastic approximation implemented
99!
100! 2035 2016-10-24 15:06:17Z suehring
101! Remove check for npex and npey in nesting case
102!
103! 2011 2016-09-19 17:29:57Z kanani
104! Added flag lsf_exception to allow explicit enabling of large scale forcing
105! together with buildings on flat terrain.
106!
107! 2007 2016-08-24 15:47:17Z kanani
108! Added call to urban surface model for reading of &urban_surface_par
109!
110! 2004 2016-08-24 10:25:59Z suehring
111! Humidity and passive scalar treated separately in nesting mode
112!
113! 2000 2016-08-20 18:09:15Z knoop
114! Forced header and separation lines into 80 columns
115!
116! 1992 2016-08-12 15:14:59Z suehring
117! +top_scalarflux
118!
119! 1960 2016-07-12 16:34:24Z suehring
120! Allocation of s_init
121!
122! 1957 2016-07-07 10:43:48Z suehring
123! flight module added
124!
125! 1955 2016-07-01 12:38:59Z hellstea
126! The parameter intializating_actions is set to 'set_constant_profiles for
127! all nest domains in order to make sure that diagnostic variables are properly
128! initialized for nest domains. Prognostic variables are later initialized by
129! interpolation from the parent domain.
130!
131! 1917 2016-05-27 14:28:12Z witha
132! Initial version of purely vertical nesting introduced.
133!
134! 1914 2016-05-26 14:44:07Z witha
135! Added call to wind turbine model for reading of &wind_turbine_par
136!
137! 1849 2016-04-08 11:33:18Z hoffmann
138! Adapted for modularization of microphysics
139!
140! 1833 2016-04-07 14:23:03Z raasch
141! call of spectra_parin
142!
143! 1831 2016-04-07 13:15:51Z hoffmann
144! turbulence renamed collision_turbulence, drizzle renamed
145! cloud_water_sedimentation
146! curvature_solution_effects removed
147!
148! 1826 2016-04-07 12:01:39Z maronga
149! Added call to radiation model for reading of &radiation_par.
150! Added call to plant canopy model for reading of &canopy_par.
151!
152! 1817 2016-04-06 15:44:20Z maronga
153! Added call to land surface model for reading of &lsm_par
154!
155! 1804 2016-04-05 16:30:18Z maronga
156! Removed code for parameter file check (__check)
157!
158! 1783 2016-03-06 18:36:17Z raasch
159! +netcdf_deflate in d3par, netcdf module and variable names changed
160!
161! 1764 2016-02-28 12:45:19Z raasch
162! cpp-statements for nesting removed, explicit settings of boundary conditions
163! in nest domains,
164! bugfix: npex/npey message moved from inipar to d3par
165! bugfix: check of lateral boundary conditions from check_parameters to here,
166! because they will be already used in init_pegrid and init_grid
167!
168! 1762 2016-02-25 12:31:13Z hellstea
169! Introduction of nested domain feature
170!
171! 1691 2015-10-26 16:17:44Z maronga
172! Added parameter most_method. Renamed prandtl_layer to constant_flux_layer.
173!
174! 1682 2015-10-07 23:56:08Z knoop
175! Code annotations made doxygen readable
176!
177! 1560 2015-03-06 10:48:54Z keck
178! +recycling_yshift
179!
180! 1496 2014-12-02 17:25:50Z maronga
181! Renamed: "radiation -> "cloud_top_radiation"
182!
183! 1484 2014-10-21 10:53:05Z kanani
184! Changes due to new module structure of the plant canopy model:
185!   canopy-model related parameters moved to new package canopy_par in
186!   subroutine package_parin
187!
188! 1429 2014-07-15 12:53:45Z knoop
189! +ensemble_member_nr to prepare the random_generator for ensemble runs
190!
191! 1402 2014-05-09 14:25:13Z raasch
192! location messages modified, batch_job included in envpar-NAMELIST
193!
194! 1384 2014-05-02 14:31:06Z raasch
195! location messages added
196!
197! 1365 2014-04-22 15:03:56Z boeske
198! Usage of large scale forcing enabled:
199! +use_subsidence_tendencies
200!
201! 1361 2014-04-16 15:17:48Z hoffmann
202! +call_microphysics_at_all_substeps
203!
204! 1359 2014-04-11 17:15:14Z hoffmann
205! REAL constants provided with KIND-attribute
206!
207! 1353 2014-04-08 15:21:23Z heinze
208! REAL constants provided with KIND-attribute
209!
210! 1327 2014-03-21 11:00:16Z raasch
211! -data_output_format, do3d_compress, do3d_comp_prec
212!
213! 1320 2014-03-20 08:40:49Z raasch
214! ONLY-attribute added to USE-statements,
215! kind-parameters added to all INTEGER and REAL declaration statements,
216! kinds are defined in new module kinds,
217! old module precision_kind is removed,
218! revision history before 2012 removed,
219! comment fields (!:) to be used for variable explanations added to
220! all variable declaration statements
221!
222! 1318 2014-03-17 13:35:16Z raasch
223! +cpu_log_barrierwait in d3par
224!
225! 1301 2014-03-06 13:29:46Z heinze
226! +large_scale_subsidence
227!
228! 1241 2013-10-30 11:36:58Z heinze
229! +nudging
230! +large_scale_forcing
231!
232! 1216 2013-08-26 09:31:42Z raasch
233! +transpose_compute_overlap in inipar
234!
235! 1195 2013-07-01 12:27:57Z heinze
236! Bugfix: allocate ref_state
237!
238! 1179 2013-06-14 05:57:58Z raasch
239! +reference_state in inipar
240!
241! 1159 2013-05-21 11:58:22Z fricke
242! +use_cmax
243!
244! 1128 2013-04-12 06:19:32Z raasch
245! +background_communication in inipar
246!
247! 1115 2013-03-26 18:16:16Z hoffmann
248! unused variables removed
249!
250! 1092 2013-02-02 11:24:22Z raasch
251! unused variables removed
252!
253! 1065 2012-11-22 17:42:36Z hoffmann
254! +nc, c_sedimentation, limiter_sedimentation, turbulence
255! -mu_constant, mu_constant_value
256!
257! 1053 2012-11-13 17:11:03Z hoffmann
258! necessary expansions according to the two new prognostic equations (nr, qr)
259! of the two-moment cloud physics scheme and steering parameters:
260! +*_init, *_surface, *_surface_initial_change, *_vertical_gradient,
261! +*_vertical_gradient_level, surface_waterflux_*,
262! +cloud_scheme, drizzle, mu_constant, mu_constant_value, ventilation_effect
263!
264! 1036 2012-10-22 13:43:42Z raasch
265! code put under GPL (PALM 3.9)
266!
267! 1015 2012-09-27 09:23:24Z raasch
268! -adjust_mixing_length
269!
270! 1003 2012-09-14 14:35:53Z raasch
271! -grid_matching
272!
273! 1001 2012-09-13 14:08:46Z raasch
274! -cut_spline_overshoot, long_filter_factor, overshoot_limit_*, ups_limit_*
275!
276! 996 2012-09-07 10:41:47Z raasch
277! -use_prior_plot1d_parameters
278!
279! 978 2012-08-09 08:28:32Z fricke
280! -km_damp_max, outflow_damping_width
281! +pt_damping_factor, pt_damping_width
282! +z0h_factor
283!
284! 964 2012-07-26 09:14:24Z raasch
285! -cross_normalized_x, cross_normalized_y, cross_xtext, z_max_do1d,
286! z_max_do1d_normalized
287!
288! 940 2012-07-09 14:31:00Z raasch
289! +neutral in inipar
290!
291! 927 2012-06-06 19:15:04Z raasch
292! +masking_method in inipar
293!
294! 824 2012-02-17 09:09:57Z raasch
295! +curvature_solution_effects in inipar
296!
297! 809 2012-01-30 13:32:58Z maronga
298! Bugfix: replaced .AND. and .NOT. with && and ! in the preprocessor directives
299!
300! 807 2012-01-25 11:53:51Z maronga
301! New cpp directive "__check" implemented which is used by check_namelist_files
302!
303! Revision 1.1  1997/07/24 11:22:50  raasch
304! Initial revision
305!
306!
307! Description:
308! ------------
309!> This subroutine reads variables controling the run from the NAMELIST files
310!------------------------------------------------------------------------------!
311 SUBROUTINE parin
312 
313
314    USE arrays_3d,                                                             &
315        ONLY:  pt_init, q_init, ref_state, s_init, sa_init, ug, u_init, v_init,&
316               vg
317
318    USE date_and_time_mod,                                                     &
319        ONLY:  day_of_year_init, time_utc_init
320               
321    USE plant_canopy_model_mod,                                                &
322         ONLY: pcm_parin 
323
324    USE control_parameters
325
326    USE cpulog,                                                                &
327        ONLY:  cpu_log_barrierwait
328
329    USE dvrp_variables,                                                        &
330        ONLY:  local_dvrserver_running
331
332    USE flight_mod,                                                            &
333        ONLY:  flight_parin
334
335    USE grid_variables,                                                        &
336        ONLY:  dx, dy
337
338    USE indices,                                                               &
339        ONLY:  nx, ny, nz
340
341    USE kinds
342
343    USE land_surface_model_mod,                                                &
344        ONLY: lsm_parin
345
346    USE microphysics_mod,                                                      &
347        ONLY:  c_sedimentation, cloud_water_sedimentation,                     &
348               collision_turbulence, curvature_solution_effects_bulk,          &
349               dry_aerosol_radius, limiter_sedimentation, na_init, nc_const,   &
350               sigma_bulk, ventilation_effect
351
352
353    USE model_1d_mod,                                                          &
354        ONLY:  damp_level_1d, dt_pr_1d, dt_run_control_1d, end_time_1d
355
356    USE pegrid
357
358    USE netcdf_interface,                                                      &
359        ONLY:  netcdf_data_format, netcdf_deflate, netcdf_precision
360
361    USE pmc_interface,                                                         &
362        ONLY:  nested_run, nesting_mode
363
364    USE profil_parameter,                                                      &
365        ONLY:  cross_profiles, profile_columns, profile_rows
366
367    USE progress_bar,                                                          &
368        ONLY :  batch_job
369
370    USE radiation_model_mod,                                                   &
371        ONLY: radiation_parin 
372
373    USE spectra_mod,                                                           &
374        ONLY :  spectra_parin
375
376    USE statistics,                                                            &
377        ONLY:  hom, hom_sum, pr_palm, region, statistic_regions
378
379    USE synthetic_turbulence_generator_mod,                                    &
380        ONLY:  stg_parin
381
382    USE urban_surface_mod,                                                     &
383        ONLY: usm_parin
384
385    USE wind_turbine_model_mod,                                                &
386        ONLY:  wtm_parin
387
388    USE vertical_nesting_mod,                                                  &
389        ONLY:  vnest_start_time
390
391    IMPLICIT NONE
392
393    INTEGER(iwp) ::  i      !<
394    INTEGER(iwp) ::  ioerr  !< error flag for open/read/write
395    INTEGER(iwp) ::  myworldid       !<
396    INTEGER(iwp) ::  numworldprocs   !<
397
398    NAMELIST /inipar/  aerosol_bulk, alpha_surface, approximation, bc_e_b,     &
399                       bc_lr, bc_ns, bc_p_b, bc_p_t, bc_pt_b, bc_pt_t, bc_q_b, &
400             bc_q_t,bc_s_b, bc_s_t, bc_sa_t, bc_uv_b, bc_uv_t,                 &
401             bottom_salinityflux, building_height, building_length_x,          &
402             building_length_y, building_wall_left, building_wall_south,       &
403             call_psolver_at_all_substeps, call_microphysics_at_all_substeps,  &
404             canyon_height,                                                    &
405             canyon_width_x, canyon_width_y, canyon_wall_left,                 &
406             canyon_wall_south, c_sedimentation, cfl_factor, cloud_droplets,   &
407             cloud_physics, cloud_scheme, cloud_top_radiation,                 &
408             cloud_water_sedimentation,                                        &
409             collective_wait, collision_turbulence, complex_terrain,           &
410             conserve_volume_flow,                                             &
411             conserve_volume_flow_mode, constant_flux_layer,                   &
412             coupling_start_time, curvature_solution_effects_bulk,             &
413             cycle_mg, damp_level_1d,                                          &
414             data_output_during_spinup,                                        &
415             day_of_year_init,                                                 &
416             dissipation_1d,                                                   &
417             dp_external, dp_level_b, dp_smooth, dpdxy, dry_aerosol_radius,    &
418             dt, dt_pr_1d, dt_run_control_1d, dt_spinup, dx, dy, dz, dz_max,   &
419             dz_stretch_factor, dz_stretch_level, end_time_1d,                 &
420             ensemble_member_nr, e_init, e_min, fft_method,                    &
421             flux_input_mode, flux_output_mode,                                &
422             galilei_transformation, humidity,                                 &
423             inflow_damping_height, inflow_damping_width,                      &
424             inflow_disturbance_begin, inflow_disturbance_end,                 &
425             initializing_actions, km_constant,                                &
426             large_scale_forcing, large_scale_subsidence, latitude,            &
427             limiter_sedimentation, longitude,                                 &
428             loop_optimization, lsf_exception, masking_method, mg_cycles,      &
429             mg_switch_to_pe0_level, mixing_length_1d, momentum_advec,         &
430             most_method, na_init, nc_const, netcdf_precision, neutral, ngsrb, &
431             nsor, nsor_ini, nudging, nx, ny, nz, ocean, omega, omega_sor,     &
432             outflow_source_plane, passive_scalar,                             &
433             prandtl_number, precipitation, psolver, pt_damping_factor,        &
434             pt_damping_width, pt_reference, pt_surface,                       &
435             pt_surface_initial_change, pt_vertical_gradient,                  &
436             pt_vertical_gradient_level, q_surface, q_surface_initial_change,  &
437             q_vertical_gradient, q_vertical_gradient_level,                   &
438             random_generator, random_heatflux,                                &
439             rayleigh_damping_factor, rayleigh_damping_height,                 &
440             recycling_width, recycling_yshift,                                &
441             reference_state, residual_limit,                                  &
442             roughness_length, sa_surface,                                     &
443             sa_vertical_gradient, sa_vertical_gradient_level, scalar_advec,   &
444             scalar_rayleigh_damping, sigma_bulk,                              &
445             spinup_time, spinup_pt_amplitude, spinup_pt_mean,                 &
446             statistic_regions, subs_vertical_gradient,                        &
447             subs_vertical_gradient_level, surface_heatflux, surface_pressure, &
448             surface_scalarflux, surface_waterflux,                            &
449             s_surface, s_surface_initial_change, s_vertical_gradient,         &
450             s_vertical_gradient_level, time_utc_init, timestep_scheme,        &
451             topography, topography_grid_convention, top_heatflux,             &
452             top_momentumflux_u, top_momentumflux_v, top_salinityflux,         &
453             top_scalarflux, transpose_compute_overlap,                        &
454             tunnel_height, tunnel_length, tunnel_width_x, tunnel_width_y,     &
455             tunnel_wall_depth,                                                &
456             turbulent_inflow, turbulent_outflow,                              &
457             use_subsidence_tendencies, ug_surface, ug_vertical_gradient,      &
458             ug_vertical_gradient_level, use_surface_fluxes, use_cmax,         &
459             use_top_fluxes, use_ug_for_galilei_tr, use_upstream_for_tke,      &
460             uv_heights, u_bulk, u_profile, vg_surface, vg_vertical_gradient,  &
461             vg_vertical_gradient_level, v_bulk, v_profile, ventilation_effect,&
462             wall_adjustment, wall_heatflux, wall_humidityflux,                &
463             wall_salinityflux, wall_scalarflux, y_shift, zeta_max, zeta_min,  &
464             z0h_factor
465     
466    NAMELIST /d3par/  averaging_interval, averaging_interval_pr,               &
467             cpu_log_barrierwait, create_disturbances,                         &
468             cross_profiles, data_output, data_output_masks,                   &
469             data_output_pr, data_output_2d_on_each_pe, disturbance_amplitude, &
470             disturbance_energy_limit, disturbance_level_b,                    &
471             disturbance_level_t, do2d_at_begin, do3d_at_begin,                &
472             dt, dt_averaging_input, dt_averaging_input_pr,                    &
473             dt_coupling, dt_data_output, dt_data_output_av, dt_disturb,       &
474             dt_domask, dt_dopr, dt_dopr_listing, dt_dots, dt_do2d_xy,         &
475             dt_do2d_xz, dt_do2d_yz, dt_do3d, dt_max, dt_restart,              &
476             dt_run_control,end_time, force_print_header, mask_scale_x,        &
477             mask_scale_y, mask_scale_z, mask_x, mask_y, mask_z, mask_x_loop,  &
478             mask_y_loop, mask_z_loop, netcdf_data_format, netcdf_deflate,     &
479             normalizing_region, npex, npey, nz_do3d,                          &
480             precipitation_amount_interval, profile_columns, profile_rows,     &
481             restart_time, section_xy, section_xz, section_yz,                 &
482             skip_time_data_output, skip_time_data_output_av, skip_time_dopr,  &
483             skip_time_do2d_xy, skip_time_do2d_xz, skip_time_do2d_yz,          &
484             skip_time_do3d, skip_time_domask, synchronous_exchange,           &
485             termination_time_needed, vnest_start_time, z_max_do2d
486
487
488    NAMELIST /envpar/  batch_job, host, local_dvrserver_running,               &
489                       maximum_cpu_time_allowed, maximum_parallel_io_streams,  &
490                       revision, run_identifier, tasks_per_node, write_binary
491
492!
493!-- First read values of environment variables (this NAMELIST file is
494!-- generated by mrun)
495    CALL location_message( 'reading environment parameters from ENVPAR', .FALSE. )
496
497    OPEN ( 90, FILE='ENVPAR', STATUS='OLD', FORM='FORMATTED', IOSTAT=ioerr )
498
499    IF ( ioerr /= 0 )  THEN
500       message_string = 'local file ENVPAR not found' //                       &
501                        '&some variables for steering may not be properly set'
502       CALL message( 'parin', 'PA0276', 0, 1, 0, 6, 0 )
503    ELSE
504       READ ( 90, envpar, IOSTAT=ioerr )
505       IF ( ioerr < 0 )  THEN
506          message_string = 'no envpar-NAMELIST found in local file '  //       &
507                           'ENVPAR&some variables for steering may '  //       &
508                           'not be properly set'
509          CALL message( 'parin', 'PA0278', 0, 1, 0, 6, 0 )
510       ELSEIF ( ioerr > 0 )  THEN
511          message_string = 'errors in local file ENVPAR' //                    &
512                           '&some variables for steering may not be properly set'
513          CALL message( 'parin', 'PA0277', 0, 1, 0, 6, 0 )
514       ENDIF
515       CLOSE ( 90 )
516    ENDIF
517
518    CALL location_message( 'finished', .TRUE. )
519!
520!-- Calculate the number of groups into which parallel I/O is split.
521!-- The default for files which are opened by all PEs (or where each
522!-- PE opens his own independent file) is, that all PEs are doing input/output
523!-- in parallel at the same time. This might cause performance or even more
524!-- severe problems depending on the configuration of the underlying file
525!-- system.
526!-- First, set the default:
527    CALL MPI_COMM_RANK( MPI_COMM_WORLD, myworldid, ierr ) 
528    CALL MPI_COMM_SIZE( MPI_COMM_WORLD, numworldprocs, ierr )
529    IF ( maximum_parallel_io_streams == -1  .OR.                               &
530         maximum_parallel_io_streams > numworldprocs )  THEN
531       maximum_parallel_io_streams = numworldprocs
532    ENDIF
533!
534!-- Now calculate the number of io_blocks and the io_group to which the
535!-- respective PE belongs. I/O of the groups is done in serial, but in parallel
536!-- for all PEs belonging to the same group.
537!-- These settings are repeated in init_pegrid for the communicator comm2d,
538!-- which is not available here
539    !io_blocks = numprocs / maximum_parallel_io_streams
540    io_blocks = numworldprocs / maximum_parallel_io_streams
541    !io_group  = MOD( myid+1, io_blocks )
542    io_group  = MOD( myworldid+1, io_blocks )
543   
544    CALL location_message( 'reading NAMELIST parameters from PARIN', .FALSE. )
545!
546!-- Data is read in parallel by groups of PEs
547    DO  i = 0, io_blocks-1
548       IF ( i == io_group )  THEN
549
550!
551!--       Open the NAMELIST-file which is send with this job
552          CALL check_open( 11 )
553
554!
555!--       Read the control parameters for initialization.
556!--       The namelist "inipar" must be provided in the NAMELIST-file.
557          READ ( 11, inipar, ERR=10, END=11 )
558
559          GOTO 12
560
561 10       message_string = 'errors in \$inipar &or no \$inipar-namelist ' //   &
562                           'found (CRAY-machines only)'
563          CALL message( 'parin', 'PA0271', 1, 2, 0, 6, 0 )
564
565 11       message_string = 'no \$inipar-namelist found'
566          CALL message( 'parin', 'PA0272', 1, 2, 0, 6, 0 )
567
568!
569!--       If required, read control parameters from restart file (produced by
570!--       a prior run). All PEs are reading from file created by PE0 (see
571!--       check_open)
572 12       IF ( TRIM( initializing_actions ) == 'read_restart_data' )  THEN
573             CALL read_var_list
574
575!
576!--          Increment the run count
577             runnr = runnr + 1
578          ENDIF
579
580!
581!--       In case of nested runs, explicitly set nesting boundary conditions.
582!--       This will overwrite the user settings. bc_lr and bc_ns always need
583!--       to be cyclic for vertical nesting.
584          IF ( nesting_mode == 'vertical' )  THEN
585             IF (bc_lr /= 'cyclic' .OR. bc_ns /= 'cyclic' )  THEN
586                WRITE ( message_string, *)  'bc_lr and bc_ns were set to ,',   &
587                                      'cyclic for vertical nesting'
588                CALL message( 'parin', 'PA0428', 0, 0, 0, 6, 0 )
589                bc_lr   = 'cyclic'
590                bc_ns   = 'cyclic'
591             ENDIF
592             IF ( nest_domain )  THEN
593                bc_uv_t = 'nested'
594                bc_pt_t = 'nested'
595                bc_q_t  = 'nested'
596                bc_s_t  = 'nested'
597                bc_p_t  = 'neumann'
598             ENDIF
599          ELSE
600         
601!
602!--       For other nesting modes only set boundary conditions for
603!--       nested domains.
604             IF ( nest_domain )  THEN
605                bc_lr   = 'nested'
606                bc_ns   = 'nested'
607                bc_uv_t = 'nested'
608                bc_pt_t = 'nested'
609                bc_q_t  = 'nested'
610                bc_s_t  = 'nested'
611                bc_p_t  = 'neumann'
612             ENDIF
613          ENDIF
614
615!         
616!--       In case of nested runs, make sure that initializing_actions =
617!--       'set_constant_profiles' even though the constant-profiles
618!--       initializations for the prognostic variables will be overwritten
619!--       by pmci_child_initialize and pmci_parent_initialize. This is,
620!--       however, important e.g. to make sure that diagnostic variables
621!--       are set properly. An exception is made in case of restart runs and
622!--       if user decides to do everything by its own.
623!--       MS: is this setting really necessary?
624          IF ( nest_domain  .AND.  (                                           &
625               TRIM( initializing_actions ) /= 'read_restart_data'  .OR.       &
626               TRIM( initializing_actions ) /= 'by_user' ) )  THEN
627             initializing_actions = 'set_constant_profiles'
628          ENDIF
629           
630!
631!--       Check validity of lateral boundary conditions. This has to be done
632!--       here because they are already used in init_pegrid and init_grid and
633!--       therefore cannot be check in check_parameters
634          IF ( bc_lr /= 'cyclic'  .AND.  bc_lr /= 'dirichlet/radiation'  .AND. &
635               bc_lr /= 'radiation/dirichlet'  .AND.  bc_lr /= 'nested' )  THEN
636             message_string = 'unknown boundary condition: bc_lr = "' // &
637                              TRIM( bc_lr ) // '"'
638             CALL message( 'check_parameters', 'PA0049', 1, 2, 0, 6, 0 )
639          ENDIF
640          IF ( bc_ns /= 'cyclic'  .AND.  bc_ns /= 'dirichlet/radiation'  .AND. &
641               bc_ns /= 'radiation/dirichlet'  .AND.  bc_ns /= 'nested' )  THEN
642             message_string = 'unknown boundary condition: bc_ns = "' // &
643                              TRIM( bc_ns ) // '"'
644             CALL message( 'check_parameters', 'PA0050', 1, 2, 0, 6, 0 )
645          ENDIF
646
647!
648!--       Set internal variables used for speed optimization in if clauses
649          IF ( bc_lr /= 'cyclic' )               bc_lr_cyc    = .FALSE.
650          IF ( bc_lr == 'dirichlet/radiation' )  bc_lr_dirrad = .TRUE.
651          IF ( bc_lr == 'radiation/dirichlet' )  bc_lr_raddir = .TRUE.
652          IF ( bc_ns /= 'cyclic' )               bc_ns_cyc    = .FALSE.
653          IF ( bc_ns == 'dirichlet/radiation' )  bc_ns_dirrad = .TRUE.
654          IF ( bc_ns == 'radiation/dirichlet' )  bc_ns_raddir = .TRUE.
655
656!
657!--       Definition of names of areas used for computing statistics. They must
658!--       be defined at this place, because they are allowed to be redefined by
659!--       the user in user_parin.
660          region = 'total domain'
661
662!
663!--       Read runtime parameters given by the user for this run (namelist
664!--       "d3par"). The namelist "d3par" can be omitted. In that case, default
665!--       values are used for the parameters.
666          READ ( 11, d3par, END=20 )
667 20       CONTINUE
668
669!
670!--       Check if land surface model is used and read &lsm_par if required
671          CALL lsm_parin
672
673!
674!--       Check if urban surface model is used and read &urban_surface_par if required
675          CALL usm_parin
676
677!
678!--       Check if spectra shall be calculated and read spectra_par if required
679          CALL spectra_parin
680
681!
682!--       Check if radiation model is used and read &radiation_par if required
683          CALL radiation_parin
684 
685 
686!--       Check if plant canopy model is used and read &canopy_par if required
687          CALL pcm_parin
688 
689!
690!--       Read control parameters for optionally used model software packages
691          CALL package_parin
692
693!
694!--       Check if wind turbine model is used and read &wind_turbine_par if
695!--       required
696          CALL wtm_parin
697!
698!--       Check if virtual flights should be carried out and read &flight_par
699!--       if required
700          CALL flight_parin
701!
702!--       Check if synthetic turbulence generator is used and read stg_par if
703!--       required
704          CALL stg_parin
705!
706!--       Read user-defined variables
707          CALL user_parin
708
709!
710!--       The restart file will be reopened when reading the subdomain data
711          IF ( TRIM( initializing_actions ) == 'read_restart_data' )  THEN
712             CALL close_file( 13 )
713          ENDIF
714
715!
716!--       Check in case of initial run, if the grid point numbers are well
717!--       defined and allocate some arrays which are already needed in
718!--       init_pegrid or check_parameters. During restart jobs, these arrays
719!--       will be allocated in read_var_list. All other arrays are allocated
720!--       in init_3d_model.
721          IF ( TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
722
723             IF ( nx <= 0 )  THEN
724                WRITE( message_string, * ) 'no value or wrong value given',    &
725                                           ' for nx: nx=', nx
726                CALL message( 'parin', 'PA0273', 1, 2, 0, 6, 0 )
727             ENDIF
728             IF ( ny <= 0 )  THEN
729                WRITE( message_string, * ) 'no value or wrong value given',    &
730                                           ' for ny: ny=', ny
731                CALL message( 'parin', 'PA0274', 1, 2, 0, 6, 0 )
732             ENDIF
733             IF ( nz <= 0 )  THEN
734                WRITE( message_string, * ) 'no value or wrong value given',    &
735                                           ' for nz: nz=', nz
736                CALL message( 'parin', 'PA0275', 1, 2, 0, 6, 0 )
737             ENDIF
738!
739!--          ATTENTION: in case of changes to the following statement please
740!--                  also check the allocate statement in routine read_var_list
741             ALLOCATE( pt_init(0:nz+1), q_init(0:nz+1), s_init(0:nz+1),        &
742                       ref_state(0:nz+1), sa_init(0:nz+1), ug(0:nz+1),         &
743                       u_init(0:nz+1), v_init(0:nz+1), vg(0:nz+1),             &
744                       hom(0:nz+1,2,pr_palm+max_pr_user,0:statistic_regions),  &
745                       hom_sum(0:nz+1,pr_palm+max_pr_user,0:statistic_regions) )
746
747             hom = 0.0_wp
748
749          ENDIF
750
751!
752!--       NAMELIST-file is not needed anymore
753          CALL close_file( 11 )
754
755       ENDIF
756#if defined( __parallel )
757       CALL MPI_BARRIER( MPI_COMM_WORLD, ierr )
758#endif
759    ENDDO
760
761    CALL location_message( 'finished', .TRUE. )
762
763 END SUBROUTINE parin
Note: See TracBrowser for help on using the repository browser.