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

Last change on this file since 3294 was 3294, checked in by raasch, 6 years ago

modularization of the ocean code

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