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

Last change on this file since 3355 was 3355, checked in by knoop, 5 years ago

removed calc_radiation.f90 and related cloud_top_radiation namelist parameter (functionality replaced by radiation model)

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