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

Last change on this file since 1179 was 1179, checked in by raasch, 11 years ago

New:
---
Initial profiles can be used as reference state in the buoyancy term. New parameter
reference_state introduced. Calculation and handling of reference state in buoyancy term revised.
binary version for restart files changed from 3.9 to 3.9a (no downward compatibility!),
initial profile for rho added to hom (id=77)

Errors:


small bugfix for background communication (time_integration)

  • Property svn:keywords set to Id
File size: 17.8 KB
Line 
1 SUBROUTINE parin
2
3!--------------------------------------------------------------------------------!
4! This file is part of PALM.
5!
6! PALM is free software: you can redistribute it and/or modify it under the terms
7! of the GNU General Public License as published by the Free Software Foundation,
8! either version 3 of the License, or (at your option) any later 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-2012  Leibniz University Hannover
18!--------------------------------------------------------------------------------!
19!
20! Current revisions:
21! -----------------
22! +reference_state in inipar
23!
24! Former revisions:
25! -----------------
26! $Id: parin.f90 1179 2013-06-14 05:57:58Z raasch $
27!
28! 1159 2013-05-21 11:58:22Z fricke
29! +use_cmax
30!
31! 1128 2013-04-12 06:19:32Z raasch
32! +background_communication in inipar
33!
34! 1115 2013-03-26 18:16:16Z hoffmann
35! unused variables removed
36!
37! 1092 2013-02-02 11:24:22Z raasch
38! unused variables removed
39!
40! 1065 2012-11-22 17:42:36Z hoffmann
41! +nc, c_sedimentation, limiter_sedimentation, turbulence
42! -mu_constant, mu_constant_value
43!
44! 1053 2012-11-13 17:11:03Z hoffmann
45! necessary expansions according to the two new prognostic equations (nr, qr)
46! of the two-moment cloud physics scheme and steering parameters:
47! +*_init, *_surface, *_surface_initial_change, *_vertical_gradient,
48! +*_vertical_gradient_level, surface_waterflux_*,
49! +cloud_scheme, drizzle, mu_constant, mu_constant_value, ventilation_effect
50!
51! 1036 2012-10-22 13:43:42Z raasch
52! code put under GPL (PALM 3.9)
53!
54! 1015 2012-09-27 09:23:24Z raasch
55! -adjust_mixing_length
56!
57! 1003 2012-09-14 14:35:53Z raasch
58! -grid_matching
59!
60! 1001 2012-09-13 14:08:46Z raasch
61! -cut_spline_overshoot, long_filter_factor, overshoot_limit_*, ups_limit_*
62!
63! 996 2012-09-07 10:41:47Z raasch
64! -use_prior_plot1d_parameters
65!
66! 978 2012-08-09 08:28:32Z fricke
67! -km_damp_max, outflow_damping_width
68! +pt_damping_factor, pt_damping_width
69! +z0h_factor
70!
71! 964 2012-07-26 09:14:24Z raasch
72! -cross_normalized_x, cross_normalized_y, cross_xtext, z_max_do1d,
73! z_max_do1d_normalized
74!
75! 940 2012-07-09 14:31:00Z raasch
76! +neutral in inipar
77!
78! 927 2012-06-06 19:15:04Z raasch
79! +masking_method in inipar
80!
81! 824 2012-02-17 09:09:57Z raasch
82! +curvature_solution_effects in inipar
83!
84! 809 2012-01-30 13:32:58Z maronga
85! Bugfix: replaced .AND. and .NOT. with && and ! in the preprocessor directives
86!
87! 807 2012-01-25 11:53:51Z maronga
88! New cpp directive "__check" implemented which is used by check_namelist_files
89!
90! 785 2011-11-28 09:47:19Z raasch
91! +scalar_rayleigh_damping in inipar
92!
93! 767 2011-10-14 06:39:12Z raasch
94! +u_profile, v_profile, uv_heights in inipar
95!
96! 759 2011-09-15 13:58:31Z raasch
97! +maximum_parallel_io_streams in envpar,
98! splitting of parallel I/O in blocks of PEs
99!
100! 683 2011-02-09 14:25:15Z raasch
101! +synchronous_exchange in d3par
102!
103! 667 2010-12-23 12:06:00Z suehring/gryschka
104! Steering parameter dissipation_control added in inipar. (commented out)
105!
106! 622 2010-12-10 08:08:13Z raasch
107! +collective_wait in inipar
108!
109! 600 2010-11-24 16:10:51Z raasch
110! parameters moved from d3par to inipar: call_psolver_at_all_substeps,
111! cfl_factor, cycle_mg, mg_cycles, mg_switch_to_pe0_level, ngsrb, nsor,
112! omega_sor, prandtl_number, psolver, rayleigh_damping_factor,
113! rayleigh_damping_height, residual_limit
114!
115! 580 2010-10-05 13:59:11Z heinze
116! Renaming of ws_vertical_gradient to subs_vertical_gradient and
117! ws_vertical_gradient_level to subs_vertical_gradient_level
118!
119! 553 2010-09-01 14:09:06Z weinreis
120! parameters for masked output are replaced by arrays
121!
122! 493 2010-03-01 08:30:24Z raasch
123! +netcdf_data_format in d3par, -netcdf_64bit, -netcdf_64bit_3d
124!
125! 449 2010-02-02 11:23:59Z raasch
126! +wall_humidityflux, wall_scalarflux
127! +ws_vertical_gradient, ws_vertical_gradient_level
128!
129! 410 2009-12-04 17:05:40Z letzel
130! masked data output: + dt_domask, mask_01~20_x|y|z, mask_01~20_x|y|z_loop,
131! mask_scale_x|y|z, masks, skip_time_domask
132!
133! 291 2009-04-16 12:07:26Z raasch
134! +local_dvrserver_running in envpar
135! Output of messages replaced by message handling routine.
136! +canyon_height, canyon_width_x, canyon_width_y, canyon_wall_left,
137! canyon_wall_south, conserve_volume_flow_mode, coupling_start_time,
138! dp_external, dp_level_b, dp_smooth, dpdxy, u_bulk, v_bulk in inipar
139! topography_grid_convention moved from userpar
140!
141! 197 2008-09-16 15:29:03Z raasch
142! +cthf,leaf_surface_concentration, scalar_exchange_coefficient
143! +inflow_damping_height, inflow_damping_width, recycling_width,
144! turbulent_inflow in inipar, -skip_time_dosp in d3par,
145! allocation of hom_sum moved from init_3d_model to here,
146! npex, npey moved from inipar to d3par, setting of myid_char_14 removed,
147! lad is allways allocated
148!
149! 138 2007-11-28 10:03:58Z letzel
150! +canopy_mode, drag_coefficient, lad_surface, lad_vertical_gradient,
151! lad_vertical_gradient_level, pch_index, plant_canopy,
152! +allocation of leaf area density field
153!
154! 108 2007-08-24 15:10:38Z letzel
155! +e_init, top_momentumflux_u|v in inipar, +dt_coupling in d3par
156!
157! 95 2007-06-02 16:48:38Z raasch
158! +bc_sa_t, bottom_salinityflux, ocean, sa_surface, sa_vertical_gradient,
159! sa_vertical_gradient_level, top_salinityflux in inipar,
160! sa_init is allocated
161!
162! 87 2007-05-22 15:46:47Z raasch
163! Size of hom increased by the maximum number of user-defined profiles,
164! var_hom renamed pr_palm
165!
166! 82 2007-04-16 15:40:52Z raasch
167! +return_addres, return_username in envpar
168!
169! 75 2007-03-22 09:54:05Z raasch
170! +dt_max, netcdf_64bit_3d, precipitation_amount_interval in d3par,
171! +loop_optimization, pt_reference in inipar, -data_output_ts,
172! moisture renamed humidity
173!
174! 20 2007-02-26 00:12:32Z raasch
175! +top_heatflux, use_top_fluxes in inipar
176!
177! 3 2007-02-13 11:30:58Z raasch
178! +netcdf_64bit_3d in d3par,
179! RCS Log replace by Id keyword, revision history cleaned up
180!
181! Revision 1.57  2007/02/11 13:11:22  raasch
182! Values of environment variables are now read from file ENVPAR instead of
183! reading them with a system call, + NAMELIST envpar
184!
185! Revision 1.1  1997/07/24 11:22:50  raasch
186! Initial revision
187!
188!
189! Description:
190! ------------
191! This subroutine reads variables controling the run from the NAMELIST files
192!------------------------------------------------------------------------------!
193
194    USE arrays_3d
195    USE averaging
196    USE cloud_parameters
197    USE control_parameters
198    USE dvrp_variables
199    USE grid_variables
200    USE indices
201    USE model_1d
202    USE pegrid
203    USE profil_parameter
204    USE statistics
205
206    IMPLICIT NONE
207
208    INTEGER ::  i
209
210
211    NAMELIST /inipar/  alpha_surface, background_communication, bc_e_b, bc_lr, &
212                       bc_ns, bc_p_b, bc_p_t, bc_pt_b, bc_pt_t, bc_q_b, &
213             bc_q_t,bc_s_b, bc_s_t, bc_sa_t, bc_uv_b, bc_uv_t, &
214             bottom_salinityflux, building_height, building_length_x, &
215             building_length_y, building_wall_left, building_wall_south, &
216             call_psolver_at_all_substeps, canopy_mode, canyon_height, &
217             canyon_width_x, canyon_width_y, canyon_wall_left, &
218             canyon_wall_south, c_sedimentation, cfl_factor, cloud_droplets, &
219             cloud_physics, cloud_scheme, collective_wait, &
220             conserve_volume_flow, &
221             conserve_volume_flow_mode, coupling_start_time, cthf, &
222             curvature_solution_effects, cycle_mg, damp_level_1d, &
223             dissipation_1d, & !dissipation_control, &
224             dp_external, dp_level_b, dp_smooth, dpdxy, drag_coefficient, &
225             drizzle, dt, dt_pr_1d, dt_run_control_1d, dx, dy, dz, dz_max, & 
226             dz_stretch_factor, dz_stretch_level, e_init, e_min, end_time_1d, &
227             fft_method, galilei_transformation, humidity, &
228             inflow_damping_height, inflow_damping_width, &
229             inflow_disturbance_begin, inflow_disturbance_end, &
230             initializing_actions, km_constant, lad_surface, &
231             lad_vertical_gradient, lad_vertical_gradient_level, &
232             leaf_surface_concentration, limiter_sedimentation, &
233             loop_optimization, masking_method, mg_cycles, &
234             mg_switch_to_pe0_level, mixing_length_1d, momentum_advec, &
235             nc_const, netcdf_precision, neutral, ngsrb, &
236             nsor, nsor_ini, nx, ny, nz, ocean, omega, omega_sor, &
237             passive_scalar, pch_index, phi, plant_canopy, prandtl_layer, &
238             prandtl_number, precipitation, psolver, pt_damping_factor, &
239             pt_damping_width, pt_reference, pt_surface, &
240             pt_surface_initial_change, pt_vertical_gradient, &
241             pt_vertical_gradient_level, q_surface, q_surface_initial_change, &
242             q_vertical_gradient, q_vertical_gradient_level, &
243             radiation, random_generator, random_heatflux, &
244             rayleigh_damping_factor, rayleigh_damping_height, recycling_width,&
245             reference_state, residual_limit, &
246             rif_max, rif_min, roughness_length, sa_surface, &
247             sa_vertical_gradient, sa_vertical_gradient_level, scalar_advec, &
248             scalar_exchange_coefficient, scalar_rayleigh_damping, &
249             statistic_regions, subs_vertical_gradient, &
250             subs_vertical_gradient_level, surface_heatflux, surface_pressure, &
251             surface_scalarflux, surface_waterflux, &
252             s_surface, &
253             s_surface_initial_change, s_vertical_gradient, &
254             s_vertical_gradient_level, timestep_scheme, &
255             topography, topography_grid_convention, top_heatflux, &
256             top_momentumflux_u, top_momentumflux_v, top_salinityflux, &
257             turbulence, turbulent_inflow, ug_surface, ug_vertical_gradient, &
258             ug_vertical_gradient_level, use_surface_fluxes, use_cmax, &
259             use_top_fluxes, use_ug_for_galilei_tr, use_upstream_for_tke, &
260             uv_heights, u_bulk, u_profile, vg_surface, vg_vertical_gradient, &
261             vg_vertical_gradient_level, v_bulk, v_profile, ventilation_effect, &
262             wall_adjustment, wall_heatflux, wall_humidityflux, wall_scalarflux, &
263             z0h_factor
264     
265    NAMELIST /d3par/  averaging_interval, averaging_interval_pr, &
266             create_disturbances, &
267             cross_profiles, cross_ts_uymax, cross_ts_uymin, &
268             data_output, data_output_format, data_output_masks, &
269             data_output_pr, data_output_2d_on_each_pe, disturbance_amplitude, &
270             disturbance_energy_limit, disturbance_level_b, &
271             disturbance_level_t, do2d_at_begin, do3d_at_begin, do3d_compress, &
272             do3d_comp_prec, dt, dt_averaging_input, dt_averaging_input_pr, &
273             dt_coupling, dt_data_output, dt_data_output_av, dt_disturb, &
274             dt_domask, dt_dopr, dt_dopr_listing, dt_dots, dt_do2d_xy, &
275             dt_do2d_xz, dt_do2d_yz, dt_do3d, dt_max, dt_restart, &
276             dt_run_control,end_time, force_print_header, mask_scale_x, &
277             mask_scale_y, mask_scale_z, mask_x, mask_y, mask_z, mask_x_loop, &
278             mask_y_loop, mask_z_loop, netcdf_data_format, normalizing_region, &
279             npex, npey, nz_do3d, precipitation_amount_interval, &
280             profile_columns, profile_rows, restart_time, section_xy, &
281             section_xz, section_yz, skip_time_data_output, &
282             skip_time_data_output_av, skip_time_dopr, skip_time_do2d_xy, &
283             skip_time_do2d_xz, skip_time_do2d_yz, skip_time_do3d, &
284             skip_time_domask, synchronous_exchange, termination_time_needed, &
285             z_max_do2d
286
287
288    NAMELIST /envpar/  host, local_dvrserver_running, maximum_cpu_time_allowed,&
289                       maximum_parallel_io_streams, revision, return_addres, &
290                       return_username, run_identifier, tasks_per_node, &
291                       write_binary
292
293!
294!-- First read values of environment variables (this NAMELIST file is
295!-- generated by mrun)
296    OPEN ( 90, FILE='ENVPAR', STATUS='OLD', FORM='FORMATTED', ERR=30 )
297    READ ( 90, envpar, ERR=31, END=32 )
298    CLOSE ( 90 )
299
300!
301!-- Calculate the number of groups into which parallel I/O is split.
302!-- The default for files which are opened by all PEs (or where each
303!-- PE opens his own independent file) is, that all PEs are doing input/output
304!-- in parallel at the same time. This might cause performance or even more
305!-- severe problems depending on the configuration of the underlying file
306!-- system.
307!-- First, set the default:
308    IF ( maximum_parallel_io_streams == -1  .OR. &
309         maximum_parallel_io_streams > numprocs )  THEN
310       maximum_parallel_io_streams = numprocs
311    ENDIF
312!
313!-- Now calculate the number of io_blocks and the io_group to which the
314!-- respective PE belongs. I/O of the groups is done in serial, but in parallel
315!-- for all PEs belonging to the same group.
316!-- These settings are repeated in init_pegrid for the communicator comm2d,
317!-- which is not available here
318    io_blocks = numprocs / maximum_parallel_io_streams
319    io_group  = MOD( myid+1, io_blocks )
320
321!
322!-- Data is read in parallel by groups of PEs
323    DO  i = 0, io_blocks-1
324       IF ( i == io_group )  THEN
325
326!
327!--       Open the NAMELIST-file which is send with this job
328          CALL check_open( 11 )
329
330!
331!--       Read the control parameters for initialization.
332!--       The namelist "inipar" must be provided in the NAMELIST-file.
333          READ ( 11, inipar, ERR=10, END=11 )
334
335#if defined ( __check )
336!
337!--       In case of a namelist file check, &inipar from the p3d file is
338!--       used. The p3d file here must be closed and the p3df file for reading
339!--       3dpar is opened.
340          IF ( check_restart == 1 )  THEN
341             CALL close_file( 11 )
342             check_restart = 2
343             CALL check_open( 11 )             
344             initializing_actions = 'read_restart_data'
345          ENDIF
346#endif
347          GOTO 12
348
349 10       message_string = 'errors in \$inipar &or no \$inipar-namelist ' // &
350                           'found (CRAY-machines only)'
351          CALL message( 'parin', 'PA0271', 1, 2, 0, 6, 0 )
352
353 11       message_string = 'no \$inipar-namelist found'
354          CALL message( 'parin', 'PA0272', 1, 2, 0, 6, 0 )
355
356!
357!--       If required, read control parameters from restart file (produced by
358!--       a prior run). All PEs are reading from file created by PE0 (see
359!--       check_open)
360 12       IF ( TRIM( initializing_actions ) == 'read_restart_data' )  THEN
361#if ! defined ( __check )
362             CALL read_var_list
363!
364!--          The restart file will be reopened when reading the subdomain data
365             CALL close_file( 13 )
366
367!
368!--          Increment the run count
369             runnr = runnr + 1
370#endif
371          ENDIF
372
373!
374!--       Definition of names of areas used for computing statistics. They must
375!--       be defined at this place, because they are allowed to be redefined by
376!--       the user in user_parin.
377          region = 'total domain'
378
379!
380!--       Read runtime parameters given by the user for this run (namelist
381!--       "d3par"). The namelist "d3par" can be omitted. In that case, default
382!--       values are used for the parameters.
383          READ ( 11, d3par, END=20 )
384
385!
386!--       Read control parameters for optionally used model software packages
387 20       CALL package_parin
388
389!
390!--       Read user-defined variables
391          CALL user_parin
392
393!
394!--       Check in case of initial run, if the grid point numbers are well
395!--       defined and allocate some arrays which are already needed in
396!--       init_pegrid or check_parameters. During restart jobs, these arrays
397!--       will be allocated in read_var_list. All other arrays are allocated
398!--       in init_3d_model.
399          IF ( TRIM( initializing_actions ) /= 'read_restart_data' )  THEN
400
401             IF ( nx <= 0 )  THEN
402                WRITE( message_string, * ) 'no value or wrong value given', &
403                                           ' for nx: nx=', nx
404                CALL message( 'parin', 'PA0273', 1, 2, 0, 6, 0 )
405             ENDIF
406             IF ( ny <= 0 )  THEN
407                WRITE( message_string, * ) 'no value or wrong value given', &
408                                           ' for ny: ny=', ny
409                CALL message( 'parin', 'PA0274', 1, 2, 0, 6, 0 )
410             ENDIF
411             IF ( nz <= 0 )  THEN
412                WRITE( message_string, * ) 'no value or wrong value given', &
413                                           ' for nz: nz=', nz
414                CALL message( 'parin', 'PA0275', 1, 2, 0, 6, 0 )
415             ENDIF
416!
417!--          ATTENTION: in case of changes to the following statement please
418!--                  also check the allocate statement in routine read_var_list
419             ALLOCATE( lad(0:nz+1),pt_init(0:nz+1), q_init(0:nz+1),           &
420                       sa_init(0:nz+1), ug(0:nz+1), u_init(0:nz+1),           &
421                       v_init(0:nz+1), vg(0:nz+1),                            &
422                       hom(0:nz+1,2,pr_palm+max_pr_user,0:statistic_regions), &
423                       hom_sum(0:nz+1,pr_palm+max_pr_user,0:statistic_regions) )
424
425             hom = 0.0
426
427          ENDIF
428
429!
430!--       NAMELIST-file is not needed anymore
431          CALL close_file( 11 )
432
433       ENDIF
434#if defined( __parallel ) && ! ( __check )
435       CALL MPI_BARRIER( MPI_COMM_WORLD, ierr )
436#endif
437    ENDDO
438
439    RETURN
440
441 30 message_string = 'local file ENVPAR not found' // &
442                     '&some variables for steering may not be properly set'
443    CALL message( 'parin', 'PA0276', 0, 1, 0, 6, 0 )
444    RETURN
445
446 31 message_string = 'errors in local file ENVPAR' // &
447                     '&some variables for steering may not be properly set'
448    CALL message( 'parin', 'PA0277', 0, 1, 0, 6, 0 )
449    RETURN
450
451 32 message_string = 'no envpar-NAMELIST found in local file ENVPAR'  // &
452                     '&some variables for steering may not be properly set'
453    CALL message( 'parin', 'PA0278', 0, 1, 0, 6, 0 )
454
455 END SUBROUTINE parin
Note: See TracBrowser for help on using the repository browser.