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

Last change on this file since 2921 was 2921, checked in by Giersch, 6 years ago

further inipar parameter has been added to restart data, bugfix in spinup mechanism

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