source: palm/trunk/SOURCE/modules.f90 @ 349

Last change on this file since 349 was 345, checked in by heinze, 15 years ago

In case of restart runs without extension, initial profiles are not written to NetCDF-file anymore.

  • Property svn:keywords set to Id
File size: 54.9 KB
Line 
1 MODULE advection
2
3
4!------------------------------------------------------------------------------!
5! Current revisions:
6! -----------------
7! +output_for_t0
8! translation error of actual -> current revisions fixed
9! +q* in dots_label, dots_unit. increased dots_num respectively
10! typographical error in dots_unit fixed
11! +clip_dvrp_*, cluster_size, color_interval, dvrpsize_interval, dvrp_overlap,
12! dvrp_total_overlap, groundplate_color, local_dvrserver_running, n*_dvrp,
13! interval_*_dvrp_prt, isosurface_color, particle_color, particle_dvrpsize,
14! topography color in dvrp_variables,
15! vertical_particle_advection is a 1d-array,
16! +canyon_height, canyon_width_x, canyon_width_y, canyon_wall_left,
17! canyon_wall_south, conserve_volume_flow_mode, coupling_start_time,
18! dp_external, dp_level_b, dp_level_ind_b, dp_smooth, dp_smooth_factor, dpdxy,
19! run_coupled, time_since_reference_point, u_bulk, v_bulk in control_parameters,
20! default value of grid_matching changed to strict
21!
22! Former revisions:
23! -----------------
24! $Id: modules.f90 345 2009-07-01 14:37:56Z maronga $
25!
26! 217 2008-12-09 18:00:48Z letzel
27! +topography_grid_convention
28! some dvrp-variables changed to single precision, variables for dvrp-mode
29! pathlines added, +target_id, abort_mode, message_string
30!
31! 197 2008-09-16 15:29:03Z raasch
32! allow 100 spectra levels instead of 10 for consistency with
33! define_netcdf_header, +canopy_heat_flux, cthf, lai,
34! +leaf_surface_concentration, scalar_exchange_coefficient, sec, sls
35! +hor_index_bounds, hor_index_bounds_previous_run, id_inflow, id_recycling,
36! inflow_damping_*, mean_inflow_profiles, numprocs_previous_run, nx_on_file,
37! ny_on_file, offset_ocean_*, recycling_plane, recycling_width, turbulent_inflow
38! -myid_char_14
39!
40! 138 2007-11-28 10:03:58Z letzel
41! +drag_coefficient, pch_index, lad_surface, lad_vertical_gradient,
42! lad_vertical_gradient_level, plant_canopy, lad, lad_s, lad_u, lad_v,
43! lad_w, cdc, lad_vertical_gradient_level_ind, canopy_mode
44! +dt_sort_particles, ngp_2dh_s_inner, time_sort_particles, flags,
45! wall_flags_1..10, wall_humidityflux(0:4), wall_qflux(0:4),
46! wall_salinityflux(0:4), wall_scalarflux(0:4)
47!
48! 108 2007-08-24 15:10:38Z letzel
49! +comm_inter, constant_top_momentumflux, coupling_char, coupling_mode,
50! coupling_mode_remote, c_u, c_v, c_w, dt_coupling, e_init, humidity_remote,
51! ngp_xy, nxlu, nysv, port_name, qswst_remote, terminate_coupled,
52! terminate_coupled_remote, time_coupling, top_momentumflux_u|v, type_xy,
53! uswst*, vswst*
54!
55! 97 2007-06-21 08:23:15Z raasch
56! +atmos_ocean_sign, ocean, r, + salinity variables
57! defaults of .._vertical_gradient_levels changed from -1.0 to -9999999.9
58! hydro_press renamed hyp, use_pt_reference renamed use_reference
59!
60! 89 2007-05-25 12:08:31Z raasch
61! +data_output_pr_user, max_pr_user, size of data_output_pr, dopr_index,
62! dopr_initial_index and dopr_unit enlarged,
63! var_hom and var_sum renamed pr_palm
64!
65! 82 2007-04-16 15:40:52Z raasch
66! +return_addres, return_username
67! Cpp-directive lcmuk renamed lc
68!
69! 75 2007-03-22 09:54:05Z raasch
70! +arrays precipitation_amount, precipitation_rate, precipitation_rate_av,
71! rif_wall, z0_av, +arrays u_m_l, u_m_r, etc. for radiation boundary conditions,
72! +loop_optimization, netcdf_64bit_3d, zu_s_inner, zw_w_inner, id_var_zusi_*,
73! id_var_zwwi_*, ts_value, u_nzb_p1_for_vfc, v_nzb_p1_for_vfc, pt_reference,
74! use_pt_reference, precipitation_amount_interval, revision
75! +age_m in particle_type, moisture renamed humidity,
76! -data_output_ts, dots_n, uvmean_outflow, uxrp, vynp,
77! arrays dots_label and dots_unit now dimensioned with dots_max,
78! setting of palm version moved to main program
79!
80! 37 2007-03-01 08:33:54Z raasch
81! +constant_top_heatflux, top_heatflux, use_top_fluxes, +arrays for top fluxes,
82! +nzt_diff, default of bc_pt_t renamed "initial_gradient"
83! Bugfix: p is not a pointer
84!
85! RCS Log replace by Id keyword, revision history cleaned up
86!
87! Revision 1.95  2007/02/11 13:18:30  raasch
88! version 3.1b (last under RCS control)
89!
90! Revision 1.1  1997/07/24 11:21:26  raasch
91! Initial revision
92!
93!
94! Description:
95! ------------
96! Definition of variables for special advection schemes
97!------------------------------------------------------------------------------!
98
99    REAL ::  spl_gamma_x, spl_gamma_y
100
101    REAL, DIMENSION(:), ALLOCATABLE   ::  aex, bex, dex, eex, spl_z_x, spl_z_y
102    REAL, DIMENSION(:,:), ALLOCATABLE ::  spl_tri_x, spl_tri_y, spl_tri_zu, &
103                                          spl_tri_zw
104
105    SAVE
106
107 END MODULE advection
108
109
110
111
112 MODULE array_kind
113
114!------------------------------------------------------------------------------!
115! Description:
116! ------------
117! Definition of type parameters (used for the definition of single or double
118! precision variables)
119!------------------------------------------------------------------------------!
120
121    INTEGER, PARAMETER ::  dpk = SELECTED_REAL_KIND( 12 ), &
122                           spk = SELECTED_REAL_KIND( 6 )
123
124    SAVE
125
126 END MODULE array_kind
127
128
129
130
131 MODULE arrays_3d
132
133!------------------------------------------------------------------------------!
134! Description:
135! ------------
136! Definition of all arrays defined on the computational grid
137!------------------------------------------------------------------------------!
138
139    USE array_kind
140
141    REAL, DIMENSION(:), ALLOCATABLE ::                                         &
142          ddzu, dd2zu, dzu, ddzw, dzw, hyp, inflow_damping_factor, km_damp_x,  &
143          km_damp_y, lad, l_grid, pt_init, q_init, rdf, sa_init, ug, u_init,   &
144          u_nzb_p1_for_vfc, vg, v_init, v_nzb_p1_for_vfc, zu, zw
145
146    REAL, DIMENSION(:,:), ALLOCATABLE ::                                       &
147          c_u, c_v, c_w, dzu_mg, dzw_mg, f1_mg, f2_mg, f3_mg,                  &
148          mean_inflow_profiles, pt_slope_ref, qs, qswst_remote, ts, us, z0
149
150    REAL, DIMENSION(:,:), ALLOCATABLE, TARGET ::                               &
151          qsws_1, qsws_2, qswst_1, qswst_2, rif_1, rif_2, saswsb_1, saswst_1,  &
152          shf_1, shf_2, tswst_1, tswst_2, usws_1, usws_2, uswst_1, uswst_2,    &
153          vsws_1, vsws_2, vswst_1, vswst_2
154
155    REAL, DIMENSION(:,:), POINTER ::                                           &
156          qsws, qsws_m, qswst, qswst_m, rif, rif_m, saswsb, saswst, shf,       &
157          shf_m, tswst, tswst_m, usws, uswst, usws_m, uswst_m, vsws, vswst,    &
158          vsws_m, vswst_m
159
160    REAL, DIMENSION(:,:,:), ALLOCATABLE ::                                     &
161          canopy_heat_flux, cdc, d, diss, lad_s, lad_u, lad_v, lad_w, lai,     &
162          l_wall, sec, sls, tend, u_m_l, u_m_n, u_m_r, u_m_s, v_m_l, v_m_n,    &
163          v_m_r, v_m_s, w_m_l, w_m_n, w_m_r, w_m_s
164
165    REAL, DIMENSION(:,:,:), ALLOCATABLE, TARGET ::                             &
166          ql_v, ql_vp
167
168    REAL, DIMENSION(:,:,:), ALLOCATABLE, TARGET ::                             &
169          e_1, e_2, e_3, kh_1, kh_2, km_1, km_2, p, pt_1, pt_2, pt_3, q_1,     &
170          q_2, q_3, ql_1, ql_2, rho_1, sa_1, sa_2, sa_3, u_1, u_2, u_3, v_1,   &
171          v_2, v_3, vpt_1, vpt_2, w_1, w_2, w_3
172
173    REAL, DIMENSION(:,:,:), POINTER ::                                         &
174          e, e_m, e_p, kh, kh_m, km, km_m, pt, pt_m, pt_p, q, q_m, q_p, ql,    &
175          ql_c, rho, sa, sa_p, te_m, tpt_m, tq_m, tsa_m, tu_m, tv_m, tw_m, u,  &
176          u_m, u_p, v, v_m, v_p, vpt, vpt_m, w, w_m, w_p
177
178    REAL, DIMENSION(:,:,:,:), ALLOCATABLE ::  rif_wall
179
180    SAVE
181
182 END MODULE arrays_3d
183
184
185
186
187 MODULE averaging
188
189!------------------------------------------------------------------------------!
190! Description:
191! ------------
192! Definition of variables needed for time-averaging of 2d/3d data
193!------------------------------------------------------------------------------!
194
195    REAL, DIMENSION(:,:), ALLOCATABLE ::  lwp_av, precipitation_rate_av, &
196                                          ts_av, us_av, z0_av
197
198    REAL, DIMENSION(:,:,:), ALLOCATABLE, TARGET :: &
199          e_av, p_av, pc_av, pr_av, pt_av, q_av, ql_av, ql_c_av, ql_v_av, &
200          ql_vp_av, qv_av, rho_av, s_av, sa_av, u_av, v_av, vpt_av, w_av
201
202 END MODULE averaging
203
204
205
206
207 MODULE cloud_parameters
208
209!------------------------------------------------------------------------------!
210! Description:
211! ------------
212! Definition of variables and constants for cloud physics
213!------------------------------------------------------------------------------!
214
215    REAL  ::  b_cond, cp = 1005.0, diff_coeff_l = 0.23E-4,                     &
216              effective_coll_efficiency, l_d_cp, l_d_r, l_d_rv, l_v = 2.5E+06, &
217              mass_of_solute, molecular_weight_of_solute,                      &
218              prec_time_const = 0.001, ql_crit = 0.0005, rho_l = 1.0E3,        &
219              r_d = 287.0, r_v = 461.51, thermal_conductivity_l = 2.43E-2
220
221    REAL, DIMENSION(:), ALLOCATABLE   ::  pt_d_t, t_d_pt
222
223    REAL, DIMENSION(:,:), ALLOCATABLE ::  precipitation_amount, &
224                                          precipitation_rate
225
226    SAVE
227
228 END MODULE cloud_parameters
229
230
231
232
233 MODULE constants
234
235!------------------------------------------------------------------------------!
236! Description:
237! ------------
238! Definition of general constants
239!------------------------------------------------------------------------------!
240
241    REAL    ::  pi = 3.141592654
242
243    SAVE
244
245 END MODULE constants
246
247
248
249
250 MODULE control_parameters
251
252!------------------------------------------------------------------------------!
253! Description:
254! ------------
255! Definition of parameters for program control
256!------------------------------------------------------------------------------!
257
258    TYPE plot_precision
259       CHARACTER (LEN=6) ::  variable
260       INTEGER           ::  precision
261    END TYPE plot_precision
262
263    TYPE(plot_precision), DIMENSION(100) ::  plot_3d_precision =               &
264                        (/ plot_precision( 'u', 2 ), plot_precision( 'v', 2 ), &
265                           plot_precision( 'w', 2 ), plot_precision( 'p', 5 ), &
266                           plot_precision( 'pt', 2 ),                          &
267                           ( plot_precision( ' ', 1 ), i9 = 1,95 ) /)
268
269    TYPE file_status
270       LOGICAL ::  opened, opened_before
271    END TYPE file_status
272
273    TYPE(file_status), DIMENSION(200) :: openfile = file_status(.FALSE.,.FALSE.)
274
275
276    CHARACTER (LEN=1)    ::  cycle_mg = 'w', timestep_reason = ' '
277    CHARACTER (LEN=2)    ::  coupling_char = ''
278    CHARACTER (LEN=5)    ::  write_binary = 'false'
279    CHARACTER (LEN=6)    ::  grid_matching = 'strict'
280    CHARACTER (LEN=8)    ::  run_date, run_time
281    CHARACTER (LEN=9)    ::  simulated_time_chr
282    CHARACTER (LEN=11)   ::  topography_grid_convention = ' '
283    CHARACTER (LEN=12)   ::  version = ' ', revision = ' '
284    CHARACTER (LEN=16)   ::  conserve_volume_flow_mode = 'default', &
285                             loop_optimization = 'default', &
286                             momentum_advec = 'pw-scheme', &
287                             psolver = 'poisfft', &
288                             scalar_advec = 'pw-scheme'
289    CHARACTER (LEN=20)   ::  bc_e_b = 'neumann', bc_lr = 'cyclic', &
290                             bc_ns = 'cyclic', bc_p_b = 'neumann', &
291                             bc_p_t = 'dirichlet', bc_pt_b = 'dirichlet', &
292                             bc_pt_t = 'initial_gradient', &
293                             bc_q_b = 'dirichlet', bc_q_t = 'neumann', &
294                             bc_s_b = 'dirichlet', bc_s_t = 'neumann', &
295                             bc_sa_t = 'neumann', &
296                             bc_uv_b = 'dirichlet', bc_uv_t = 'dirichlet', &
297                             canopy_mode = 'block', &
298                             coupling_mode = 'uncoupled', &
299                             coupling_mode_remote = 'uncoupled', &
300                             dissipation_1d = 'as_in_3d_model', &
301                             fft_method = 'system-specific', &
302                             mixing_length_1d = 'as_in_3d_model', &
303                             random_generator = 'numerical-recipes', &
304                             return_addres, return_username, &
305                             timestep_scheme = 'runge-kutta-3'
306    CHARACTER (LEN=40)   ::  avs_data_file, topography = 'flat'
307    CHARACTER (LEN=64)   ::  host = ' '
308    CHARACTER (LEN=80)   ::  log_message, run_identifier
309    CHARACTER (LEN=100)  ::  initializing_actions = ' ', run_description_header
310    CHARACTER (LEN=1000) ::  message_string = ' '
311
312    CHARACTER (LEN=7),  DIMENSION(100) ::  do3d_comp_prec = ' '
313    CHARACTER (LEN=10), DIMENSION(10)  ::  data_output_format = ' '
314    CHARACTER (LEN=10), DIMENSION(100) ::  data_output = ' ',    &
315                                           data_output_user = ' ', doav = ' '
316    CHARACTER (LEN=10), DIMENSION(300) ::  data_output_pr = ' '
317    CHARACTER (LEN=10), DIMENSION(200) ::  data_output_pr_user = ' '
318    CHARACTER (LEN=20), DIMENSION(10)  ::  netcdf_precision = ' '
319
320    CHARACTER (LEN=10), DIMENSION(0:1,100) ::  do2d = ' ', do3d = ' '
321
322    INTEGER ::  abort_mode = 1, average_count_pr = 0, average_count_sp = 0, &
323                average_count_3d = 0, current_timestep_number = 0, &
324                dist_range = 0, disturbance_level_ind_b, &
325                disturbance_level_ind_t, doav_n = 0, dopr_n = 0, &
326                dopr_time_count = 0, dopts_time_count = 0, &
327                dosp_time_count = 0, dots_time_count = 0, &
328                do2d_xy_n = 0, do2d_xz_n = 0, do2d_yz_n = 0, do3d_avs_n = 0, &
329                dp_level_ind_b = 0, &
330                dvrp_filecount = 0, dz_stretch_level_index, gamma_mg, &
331                grid_level, ibc_e_b, ibc_p_b, ibc_p_t, ibc_pt_b, ibc_pt_t, &
332                ibc_q_b, ibc_q_t, ibc_sa_t, ibc_uv_b, ibc_uv_t, &
333                inflow_disturbance_begin = -1, inflow_disturbance_end = -1, &
334                intermediate_timestep_count, intermediate_timestep_count_max, &
335                iran = -1234567, last_dt_change = 0, maximum_grid_level, &
336                max_pr_user = 0, mgcycles = 0, mg_cycles = -1, &
337                mg_switch_to_pe0_level = 0, ngsrb = 2, nsor = 20, &
338                nsor_ini = 100, n_sor, normalizing_region = 0, &
339                nz_do1d, nz_do3d = -9999, outflow_damping_width = -1, &
340                pch_index = 0, prt_time_count = 0, recycling_plane, runnr = 0, &
341                skip_do_avs = 0, terminate_coupled = 0, &
342                terminate_coupled_remote = 0, timestep_count = 0
343
344    INTEGER ::  dist_nxl(0:1), dist_nxr(0:1), dist_nyn(0:1), dist_nys(0:1), &
345                do2d_no(0:1) = 0, do2d_xy_time_count(0:1), &
346                do2d_xz_time_count(0:1), do2d_yz_time_count(0:1), &
347                do3d_no(0:1) = 0, do3d_time_count(0:1), &
348                lad_vertical_gradient_level_ind(10) = -9999, &
349                pt_vertical_gradient_level_ind(10) = -9999, &
350                q_vertical_gradient_level_ind(10) = -9999, &
351                sa_vertical_gradient_level_ind(10) = -9999, &
352                section(100,3), section_xy(100) = -9999, &
353                section_xz(100) = -9999, section_yz(100) = -9999, &
354                ug_vertical_gradient_level_ind(10) = -9999, &
355                vg_vertical_gradient_level_ind(10) = -9999
356
357    INTEGER, DIMENSION(:), ALLOCATABLE ::  grid_level_count
358
359    LOGICAL ::  adjust_mixing_length = .FALSE., avs_output = .FALSE., &
360                call_psolver_at_all_substeps = .TRUE., &
361                cloud_droplets = .FALSE., cloud_physics = .FALSE., &
362                conserve_volume_flow = .FALSE., constant_diffusion = .FALSE., &
363                constant_heatflux = .TRUE., constant_top_heatflux = .TRUE., &
364                constant_top_momentumflux = .FALSE., &
365                constant_top_salinityflux = .TRUE., &
366                constant_waterflux = .TRUE., create_disturbances = .TRUE., &
367                cut_spline_overshoot = .TRUE., &
368                data_output_2d_on_each_pe = .TRUE., do2d_at_begin = .FALSE., &
369                do3d_at_begin = .FALSE., do3d_compress = .FALSE., &
370                do_sum = .FALSE., dp_external = .FALSE., dp_smooth = .FALSE., &
371                dt_changed = .FALSE., dt_fixed = .FALSE., &
372                disturbance_created = .FALSE., &
373                first_call_advec_particles = .TRUE., &
374                force_print_header = .FALSE., galilei_transformation = .FALSE.,&
375                humidity = .FALSE., humidity_remote = .FALSE., &
376                inflow_l = .FALSE., inflow_n = .FALSE., &
377                inflow_r = .FALSE., inflow_s = .FALSE., iso2d_output = .FALSE.,&
378                mg_switch_to_pe0 = .FALSE., &
379                netcdf_output = .FALSE., netcdf_64bit = .FALSE., &
380                netcdf_64bit_3d = .TRUE., ocean = .FALSE., &
381                outflow_l = .FALSE., outflow_n = .FALSE., outflow_r = .FALSE., &
382                outflow_s = .FALSE., passive_scalar = .FALSE., &
383                plant_canopy = .FALSE., &
384                prandtl_layer = .TRUE., precipitation = .FALSE., &
385                profil_output = .FALSE., radiation = .FALSE., &
386                random_heatflux = .FALSE., run_control_header = .FALSE., &
387                run_coupled = .TRUE., sloping_surface = .FALSE., &
388                stop_dt = .FALSE., terminate_run = .FALSE., &
389                turbulent_inflow = .FALSE., &
390                use_prior_plot1d_parameters = .FALSE., use_reference = .FALSE.,&
391                use_surface_fluxes = .FALSE., use_top_fluxes = .FALSE., &
392                use_ug_for_galilei_tr = .TRUE., use_upstream_for_tke = .FALSE.,&
393                wall_adjustment = .TRUE.
394
395    LOGICAL ::  data_output_xy(0:1) = .FALSE., data_output_xz(0:1) = .FALSE., &
396                data_output_yz(0:1) = .FALSE.
397
398    REAL ::  advected_distance_x = 0.0, advected_distance_y = 0.0, &
399             alpha_surface = 0.0, asselin_filter_factor = 0.1, &
400             atmos_ocean_sign = 1.0, &
401             averaging_interval = 0.0, averaging_interval_pr = 9999999.9, &
402             averaging_interval_sp = 9999999.9, bc_pt_t_val, bc_q_t_val, &
403             bottom_salinityflux = 0.0, &
404             building_height = 50.0, building_length_x = 50.0, &
405             building_length_y = 50.0, building_wall_left = 9999999.9, &
406             building_wall_south = 9999999.9, canyon_height = 50.0, &
407             canyon_width_x = 9999999.9, canyon_width_y = 9999999.9, &
408             canyon_wall_left = 9999999.9, canyon_wall_south = 9999999.9, &
409             cthf = 0.0, cfl_factor = -1.0, cos_alpha_surface, &
410             coupling_start_time, disturbance_amplitude = 0.25, &
411             disturbance_energy_limit = 0.01, &
412             disturbance_level_b = -9999999.9, &
413             disturbance_level_t = -9999999.9, &
414             dp_level_b = 0.0, drag_coefficient = 0.0, &
415             dt = -1.0, dt_averaging_input = 0.0, &
416             dt_averaging_input_pr = 9999999.9, dt_coupling = 9999999.9, &
417             dt_data_output = 9999999.9, &
418             dt_data_output_av = 9999999.9, dt_disturb = 9999999.9, &
419             dt_dopr = 9999999.9, dt_dopr_listing = 9999999.9, &
420             dt_dopts = 9999999.9, dt_dosp = 9999999.9, dt_dots = 9999999.9, &
421             dt_do2d_xy = 9999999.9, dt_do2d_xz = 9999999.9, &
422             dt_do2d_yz = 9999999.9, dt_do3d = 9999999.9, dt_dvrp = 9999999.9, &
423             dt_max = 20.0, dt_prel = 9999999.9, dt_restart = 9999999.9, &
424             dt_run_control = 60.0, dt_3d = -1.0, dz = -1.0, &
425             dz_max = 9999999.9, dz_stretch_factor = 1.08, &
426             dz_stretch_level = 100000.0, e_init = 0.0, e_min = 0.0, &
427             end_time = 0.0, &
428             f = 0.0, fs = 0.0, g = 9.81, inflow_damping_height = 9999999.9, &
429             inflow_damping_width = 9999999.9, kappa = 0.4, km_constant = -1.0,&
430             km_damp_max = -1.0, lad_surface = 0.0,  &
431             leaf_surface_concentration = 0.0, long_filter_factor = 0.0, &
432             maximum_cpu_time_allowed = 0.0, molecular_viscosity = 1.461E-5, &
433             old_dt = 1.0E-10, omega = 7.29212E-5, omega_sor = 1.8, &
434             overshoot_limit_e = 0.0, overshoot_limit_pt = 0.0, &
435             overshoot_limit_u = 0.0, overshoot_limit_v = 0.0, &
436             overshoot_limit_w = 0.0, particle_maximum_age = 9999999.9, &
437             phi = 55.0, prandtl_number = 1.0, &
438             precipitation_amount_interval = 9999999.9, prho_reference, &
439             pt_reference = 9999999.9, pt_slope_offset = 0.0, &
440             pt_surface = 300.0, pt_surface_initial_change = 0.0, &
441             q_surface = 0.0, q_surface_initial_change = 0.0, &
442             rayleigh_damping_factor = -1.0, rayleigh_damping_height = -1.0, &
443             recycling_width = 9999999.9, residual_limit = 1.0E-4, &
444             restart_time = 9999999.9, rho_reference, rho_surface, &
445             rif_max = 1.0, rif_min = -5.0, roughness_length = 0.1, &
446             sa_surface = 35.0, scalar_exchange_coefficient = 0.0, &
447             simulated_time = 0.0, simulated_time_at_begin, sin_alpha_surface, &
448             skip_time_data_output = 0.0, skip_time_data_output_av = 9999999.9,&
449             skip_time_dopr = 9999999.9, skip_time_dosp = 9999999.9, &
450             skip_time_do2d_xy = 9999999.9, skip_time_do2d_xz = 9999999.9, &
451             skip_time_do2d_yz = 9999999.9, skip_time_do3d = 9999999.9, &
452             surface_heatflux = 9999999.9, surface_pressure = 1013.25, &
453             surface_scalarflux = 0.0, surface_waterflux = 0.0, &
454             s_surface = 0.0, s_surface_initial_change = 0.0, &
455             termination_time_needed = -1.0, time_coupling = 0.0, &
456             time_disturb = 0.0, time_dopr = 0.0, time_dopr_av = 0.0, &
457             time_dopr_listing = 0.0, time_dopts = 0.0, time_dosp = 0.0, &
458             time_dosp_av = 0.0, time_dots = 0.0, time_do2d_xy = 0.0, &
459             time_do2d_xz = 0.0, time_do2d_yz = 0.0, time_do3d = 0.0, &
460             time_do_av = 0.0, time_do_sla = 0.0, time_dvrp = 0.0, &
461             time_prel = 0.0, time_restart = 9999999.9, time_run_control = 0.0,&
462             time_since_reference_point, top_heatflux = 9999999.9, &
463             top_momentumflux_u = 9999999.9, &
464             top_momentumflux_v = 9999999.9, top_salinityflux = 9999999.9, &
465             ug_surface = 0.0, u_bulk = 0.0, u_gtrans = 0.0, &
466             ups_limit_e = 0.0, ups_limit_pt = 0.0, ups_limit_u = 0.0, &
467             ups_limit_v = 0.0, ups_limit_w = 0.0, vg_surface = 0.0, &
468             v_bulk = 0.0, v_gtrans = 0.0, wall_adjustment_factor = 1.8, &
469             z_max_do1d = -1.0, z_max_do1d_normalized = -1.0, z_max_do2d = -1.0
470
471    REAL ::  do2d_xy_last_time(0:1) = -1.0, do2d_xz_last_time(0:1) = -1.0, &
472             do2d_yz_last_time(0:1) = -1.0, dpdxy(1:2) = 0.0, &
473             lad_vertical_gradient(10) = 0.0, &
474             lad_vertical_gradient_level(10) = -9999999.9, &
475             pt_vertical_gradient(10) = 0.0, &
476             pt_vertical_gradient_level(10) = -9999999.9, &
477             q_vertical_gradient(10) = 0.0, &
478             q_vertical_gradient_level(10) = -1.0, &
479             s_vertical_gradient(10) = 0.0, &
480             s_vertical_gradient_level(10) = -1.0, &
481             sa_vertical_gradient(10) = 0.0, &
482             sa_vertical_gradient_level(10) = -9999999.9, threshold(20) = 0.0, &
483             tsc(10) = (/ 0.0, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 /), &
484             ug_vertical_gradient(10) = 0.0, &
485             ug_vertical_gradient_level(10) = -9999999.9, &
486             vg_vertical_gradient(10) = 0.0, &
487             vg_vertical_gradient_level(10) = -9999999.9, &
488             volume_flow(1:2) = 0.0, volume_flow_area(1:2) = 0.0, &
489             volume_flow_initial(1:2) = 0.0, wall_heatflux(0:4) = 0.0, &
490             wall_humidityflux(0:4) = 0.0, wall_qflux(0:4) = 0.0, &
491             wall_salinityflux(0:4) = 0.0, wall_scalarflux(0:4) = 0.0
492
493    REAL, DIMENSION(:), ALLOCATABLE ::  dp_smooth_factor
494
495
496    SAVE
497
498 END MODULE control_parameters
499
500
501
502
503 MODULE cpulog
504
505!------------------------------------------------------------------------------!
506! Description:
507! ------------
508! Definition of variables for cpu-time measurements
509!------------------------------------------------------------------------------!
510
511    REAL ::  initial_wallclock_time
512
513    TYPE logpoint
514       REAL               ::  isum, ivect, mean, mtime, mtimevec, sum, vector
515       INTEGER            ::  counts
516       CHARACTER (LEN=20) ::  place
517    END TYPE logpoint
518
519    TYPE(logpoint), DIMENSION(100) ::  log_point = logpoint( 0.0, 0.0, 0.0,   &
520                                       0.0, 0.0, 0.0, 0.0, 0, ' ' ),          &
521                                       log_point_s = logpoint( 0.0, 0.0, 0.0, &
522                                       0.0, 0.0, 0.0, 0.0, 0, ' ' )
523
524    SAVE
525
526 END MODULE cpulog
527
528
529
530
531 MODULE dvrp_variables
532
533!------------------------------------------------------------------------------!
534! Description:
535! ------------
536! Definition of variables used with dvrp-software
537!------------------------------------------------------------------------------!
538
539    CHARACTER (LEN=10) ::  dvrp_output = 'rtsp', particle_color = 'none', &
540                           particle_dvrpsize = 'none'
541
542    CHARACTER (LEN=20), DIMENSION(10) ::  mode_dvrp = &
543                                     (/ ( '                    ', i9 = 1,10 ) /)
544
545    CHARACTER (LEN=80) ::  dvrp_directory = 'default',                    &
546                           dvrp_file      = 'default',                    &
547                           dvrp_host      = 'origin.rvs.uni-hannover.de', &
548                           dvrp_password  = '********',                   &
549                           dvrp_username  = ' '
550
551    INTEGER ::  cluster_size = 1, dvrp_colortable_entries = 4,                 &
552                dvrp_colortable_entries_prt = 22, islice_dvrp,                 &
553                nx_dvrp, nxl_dvrp, nxr_dvrp, ny_dvrp, nyn_dvrp, nys_dvrp,      &
554                nz_dvrp, pathlines_fadeintime = 5, pathlines_fadeouttime = 5,  &
555                pathlines_linecount = 1000, pathlines_maxhistory = 40,         &
556                pathlines_wavecount = 10, pathlines_wavetime = 50,             &
557                vc_gradient_normals = 0, vc_mode = 0, vc_size_x = 2,           &
558                vc_size_y = 2, vc_size_z = 2
559
560    INTEGER, DIMENSION(10) ::  slicer_position_dvrp
561
562    LOGICAL ::  cyclic_dvrp = .FALSE., dvrp_overlap, dvrp_total_overlap, &
563                local_dvrserver_running, lock_steering_update = .FALSE., &
564                use_seperate_pe_for_dvrp_output = .FALSE.
565
566    REAL    ::  clip_dvrp_l = 9999999.9, clip_dvrp_n = 9999999.9, &
567                clip_dvrp_r = 9999999.9, clip_dvrp_s = 9999999.9, &
568                superelevation = 1.0, superelevation_x = 1.0,     &
569                superelevation_y = 1.0, vc_alpha = 38.0
570
571    REAL, DIMENSION(2) ::  color_interval = (/ 0.0, 1.0 /), &
572                           dvrpsize_interval = (/ 0.0, 1.0 /)
573
574    REAL, DIMENSION(3) ::  groundplate_color = (/ 0.0, 0.6, 0.0 /), &
575                           topography_color = (/ 0.8, 0.7, 0.6 /)
576
577#if defined( __decalpha )
578    REAL, DIMENSION(2,10)  ::  slicer_range_limits_dvrp = (/ &
579                                -1.0, 1.0, -1.0, 1.0, -1.0, 1.0, -1.0, 1.0, &
580                                -1.0, 1.0, -1.0, 1.0, -1.0, 1.0, -1.0, 1.0, &
581                                -1.0, 1.0, -1.0, 1.0 /)
582
583    REAL, DIMENSION(3,10)  ::  isosurface_color = (/ &
584                                0.9, 0.9, 0.9,  0.8, 0.1, 0.1,  0.1, 0.1, 0.8, &
585                                0.1, 0.8, 0.1,  0.6, 0.1, 0.1,  0.1, 0.1, 0.6, &
586                                0.1, 0.6, 0.1,  0.4, 0.1, 0.1,  0.1, 0.1, 0.4  &
587                                0.1, 0.4, 0.1 /)
588
589    REAL(4), DIMENSION(2,100) ::  interval_values_dvrp, interval_h_dvrp =      &
590                                  (/ 270.0, 225.0, 225.0, 180.0, 70.0, 25.0,   &
591                                     25.0, -25.0, ( 0.0, i9 = 1, 192 ) /),     &
592                                  interval_l_dvrp = 0.5, interval_s_dvrp = 1.0,&
593                                  interval_a_dvrp = 0.0,                       &
594                                  interval_values_dvrp_prt,                    &
595                                  interval_h_dvrp_prt,                         &
596                                  (/ 270.0, 225.0, 225.0, 180.0, 70.0, 25.0,   &
597                                     25.0, -25.0, ( 0.0, i9 = 1, 192 ) /),     &
598                                  interval_l_dvrp_prt = 0.5,                   &
599                                  interval_s_dvrp_prt = 1.0,                   &
600                                  interval_a_dvrp_prt = 0.0
601#else
602    REAL, DIMENSION(2,10)     ::  slicer_range_limits_dvrp
603
604    REAL, DIMENSION(3,10)     ::  isosurface_color
605
606    REAL(4), DIMENSION(2,100) ::  interval_values_dvrp,                       &
607                                  interval_values_dvrp_prt, interval_h_dvrp,  &
608                                  interval_h_dvrp_prt, interval_l_dvrp = 0.5, &
609                                  interval_l_dvrp_prt = 0.5, interval_s_dvrp = 1.0, &
610                                  interval_s_dvrp_prt = 1.0, interval_a_dvrp = 0.0, &
611                                  interval_a_dvrp_prt = 0.0
612
613    DATA  slicer_range_limits_dvrp / -1.0, 1.0, -1.0, 1.0, -1.0, 1.0, &
614                                     -1.0, 1.0, -1.0, 1.0, -1.0, 1.0, &
615                                     -1.0, 1.0, -1.0, 1.0, -1.0, 1.0, &
616                                     -1.0, 1.0 /
617
618    DATA  isosurface_color / 0.9, 0.9, 0.9,  0.8, 0.1, 0.1,  0.1, 0.1, 0.8, &
619                             0.1, 0.8, 0.1,  0.6, 0.1, 0.1,  0.1, 0.1, 0.6, &
620                             0.1, 0.6, 0.1,  0.4, 0.1, 0.1,  0.1, 0.1, 0.4, &
621                             0.1, 0.4, 0.1 /
622
623    DATA  interval_h_dvrp / 270.0, 225.0, 225.0, 180.0, 70.0, 25.0, &
624                            25.0, -25.0, 192 * 0.0 /
625
626    DATA  interval_h_dvrp_prt / 270.0, 225.0, 225.0, 180.0, 70.0, 25.0, &
627                                25.0, -25.0, 192 * 0.0 /
628#endif
629
630    REAL(4), DIMENSION(:), ALLOCATABLE ::  xcoor_dvrp, ycoor_dvrp, zcoor_dvrp
631
632    TYPE steering
633       CHARACTER (LEN=20) ::  name
634       REAL(4)            ::  min, max
635       INTEGER            ::  imin, imax
636    END TYPE steering
637
638    TYPE(steering), DIMENSION(:), ALLOCATABLE ::  steering_dvrp
639
640    SAVE
641
642 END MODULE dvrp_variables
643
644
645
646
647 MODULE grid_variables
648
649!------------------------------------------------------------------------------!
650! Description:
651! ------------
652! Definition of grid spacings
653!------------------------------------------------------------------------------!
654
655    REAL ::  ddx, ddx2, dx = 1.0, dx2, ddy, ddy2, dy = 1.0, dy2
656
657    REAL, DIMENSION(:), ALLOCATABLE ::  ddx2_mg, ddy2_mg
658
659    REAL, DIMENSION(:,:), ALLOCATABLE ::  fwxm, fwxp, fwym, fwyp, fxm, fxp,   &
660                                          fym, fyp, wall_e_x, wall_e_y,       &
661                                          wall_u, wall_v, wall_w_x, wall_w_y, &
662                                          zu_s_inner, zw_w_inner
663
664    SAVE
665
666 END MODULE grid_variables
667
668
669
670
671 MODULE indices
672
673!------------------------------------------------------------------------------!
674! Description:
675! ------------
676! Definition of array bounds and number of gridpoints
677!------------------------------------------------------------------------------!
678
679    INTEGER ::  ngp_sums, nnx, nx = 0, nxa, nxl, nxlu, nxr, nxra, nx_on_file,  &
680                nny, ny = 0, nya, nyn, nyna, nys, nysv, ny_on_file, nnz,       &
681                nz = 0, nza, nzb, nzb_diff, nzt, nzta, nzt_diff
682
683    INTEGER, DIMENSION(:), ALLOCATABLE ::                                      &
684                ngp_2dh, ngp_3d, ngp_3d_inner,                                 &
685                nnx_pe, nny_pe, nxl_mg, nxr_mg, nyn_mg, nys_mg, nzt_mg
686
687    INTEGER, DIMENSION(:,:), ALLOCATABLE ::                                    &
688                ngp_2dh_outer, ngp_2dh_s_inner, mg_loc_ind, nzb_diff_s_inner,  &
689                nzb_diff_s_outer, nzb_diff_u, nzb_diff_v, nzb_inner, nzb_outer,&
690                nzb_s_inner, nzb_s_outer, nzb_u_inner, nzb_u_outer,            &
691                nzb_v_inner, nzb_v_outer, nzb_w_inner, nzb_w_outer, nzb_2d
692
693    INTEGER, DIMENSION(:,:,:), POINTER ::  flags
694
695    INTEGER, DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  wall_flags_1,           &
696                wall_flags_2, wall_flags_3, wall_flags_4, wall_flags_5,        &
697                wall_flags_6, wall_flags_7, wall_flags_8, wall_flags_9,        &
698                wall_flags_10
699
700    SAVE
701
702 END MODULE indices
703
704
705
706
707 MODULE interfaces
708
709!------------------------------------------------------------------------------!
710! Description:
711! ------------
712! Interfaces for special subroutines which use optional parameters
713!------------------------------------------------------------------------------!
714
715    INTERFACE
716
717       SUBROUTINE cpu_log( log_event, place, modus, barrierwait )
718
719          USE cpulog
720
721          CHARACTER (LEN=*)           ::  modus, place
722          CHARACTER (LEN=*), OPTIONAL ::  barrierwait
723          TYPE(logpoint)              ::  log_event
724
725       END SUBROUTINE cpu_log
726
727    END INTERFACE
728
729
730
731    INTERFACE
732
733       SUBROUTINE global_min_max ( i1, i2, j1, j2, k1, k2, feld, mode, wert, &
734                                   wert_ijk, wert1, wert1_ijk )
735          CHARACTER (LEN=*), INTENT(IN) ::  mode
736          INTEGER, INTENT(IN)           ::  i1, i2, j1, j2, k1, k2
737          INTEGER                       ::  wert_ijk(3)
738          INTEGER, OPTIONAL             ::  wert1_ijk(3)
739          REAL                          ::  wert
740          REAL, OPTIONAL                ::  wert1
741          REAL, INTENT(IN)              ::  feld(i1:i2,j1:j2,k1:k2)
742
743       END SUBROUTINE global_min_max
744
745    END INTERFACE
746
747    SAVE
748
749 END MODULE interfaces
750
751
752
753 MODULE pointer_interfaces
754
755!------------------------------------------------------------------------------!
756! Description:
757! ------------
758! Interfaces for subroutines with pointer arguments called in
759! prognostic_equations
760!------------------------------------------------------------------------------!
761
762    INTERFACE
763
764       SUBROUTINE advec_s_bc( sk, sk_char )
765
766          CHARACTER (LEN=*), INTENT(IN)   ::  sk_char
767          REAL, DIMENSION(:,:,:), POINTER ::  sk
768
769       END SUBROUTINE advec_s_bc
770
771    END INTERFACE
772
773
774    SAVE
775
776 END MODULE pointer_interfaces
777
778
779
780
781 MODULE model_1d
782
783!------------------------------------------------------------------------------!
784! Description:
785! ------------
786! Definition of variables for the 1D-model
787!------------------------------------------------------------------------------!
788
789    INTEGER ::  current_timestep_number_1d = 0, damp_level_ind_1d, &
790                last_dt_change_1d = 0 
791
792    LOGICAL ::  run_control_header_1d = .FALSE., stop_dt_1d = .FALSE. 
793
794    REAL ::     damp_level_1d = -1.0, dt_1d = 60.0, dt_max_1d = 300.0, &
795                dt_pr_1d = 9999999.9, dt_run_control_1d = 60.0, &
796                end_time_1d = 864000.0, old_dt_1d = 1.0E-10, &
797                qs1d, simulated_time_1d = 0.0, time_pr_1d = 0.0, &
798                time_run_control_1d = 0.0, ts1d, us1d, usws1d, usws1d_m, &
799                vsws1d, vsws1d_m, z01d
800
801
802    REAL, DIMENSION(:), ALLOCATABLE ::  e1d, e1d_m, e1d_p, kh1d, kh1d_m, km1d, &
803                                        km1d_m, l_black, l1d, l1d_m, rif1d,    &
804                                        te_e, te_em, te_u, te_um, te_v, te_vm, &
805                                        u1d, u1d_m, u1d_p, v1d, v1d_m, v1d_p
806
807    SAVE
808
809 END MODULE model_1d
810
811
812
813
814 MODULE netcdf_control
815
816!------------------------------------------------------------------------------!
817! Description:
818! ------------
819! Definition of parameters and variables for netcdf control.
820!------------------------------------------------------------------------------!
821
822#if defined( __netcdf )
823    USE netcdf
824#endif
825
826    INTEGER, PARAMETER ::  dopr_norm_num = 7, dopts_num = 26, dots_max = 100, &
827                           replace_num = 6
828
829    INTEGER ::  dots_num = 23
830
831    CHARACTER, DIMENSION( replace_num ) :: &
832                           replace_char = (/ '''', '"', '*', '/', '(', ')' /), &
833                           replace_by   = (/ 'p' , 'p', 's', 'o', '_', '_' /)
834
835    CHARACTER (LEN=6), DIMENSION(dopr_norm_num) ::  dopr_norm_names =   &
836         (/ 'wpt0  ', 'ws2   ', 'tsw2  ', 'ws3   ', 'ws2tsw', 'wstsw2', &
837            'z_i   ' /)
838
839    CHARACTER (LEN=6), DIMENSION(dopr_norm_num) ::  dopr_norm_longnames = &
840         (/ 'wpt0  ', 'w*2   ', 't*w2  ', 'w*3   ', 'w*2t*w', 'w*t*w2',   &
841            'z_i   ' /)
842
843    CHARACTER (LEN=7), DIMENSION(dopts_num) :: dopts_label = &
844          (/ 'tnpt   ', 'x_     ', 'y_     ', 'z_     ', 'z_abs  ', 'u      ', &
845             'v      ', 'w      ', 'u"     ', 'v"     ', 'w"     ', 'npt_up ', &
846             'w_up   ', 'w_down ', 'npt_max', 'npt_min', 'x*2    ', 'y*2    ', &
847             'z*2    ', 'u*2    ', 'v*2    ', 'w*2    ', 'u"2    ', 'v"2    ', &
848             'w"2    ', 'npt*2  ' /)
849
850    CHARACTER (LEN=7), DIMENSION(dopts_num) :: dopts_unit = &
851          (/ 'number ', 'm      ', 'm      ', 'm      ', 'm      ', 'm/s    ', &
852             'm/s    ', 'm/s    ', 'm/s    ', 'm/s    ', 'm/s    ', 'number ', &
853             'm/s    ', 'm/s    ', 'number ', 'number ', 'm2     ', 'm2     ', &
854             'm2     ', 'm2/s2  ', 'm2/s2  ', 'm2/s2  ', 'm2/s2  ', 'm2/s2  ', &
855             'm2/s2  ', 'number2' /)
856
857    CHARACTER (LEN=7), DIMENSION(dots_max) :: dots_label = &
858          (/ 'E      ', 'E*     ', 'dt     ', 'u*     ', 'th*    ', 'umax   ', &
859             'vmax   ', 'wmax   ', 'div_new', 'div_old', 'z_i_wpt', 'z_i_pt ', &
860             'w*     ', 'w"pt"0 ', 'w"pt"  ', 'wpt    ', 'pt(0)  ', 'pt(zp) ', &
861             'w"u"0  ', 'w"v"0  ', 'w"q"0  ', 'mo_L   ', 'q*     ',            &
862             ( 'unknown', i9 = 1, 77) /)
863
864    CHARACTER (LEN=7), DIMENSION(dots_max) :: dots_unit = &
865          (/ 'm2/s2  ', 'm2/s2  ', 's      ', 'm/s    ', 'K      ', 'm/s    ', &
866             'm/s    ', 'm/s    ', 's-1    ', 's-1    ', 'm      ', 'm      ', &
867             'm/s    ', 'K m/s  ', 'K m/s  ', 'K m/s  ', 'K      ', 'K      ', &
868             'm2/s2  ', 'm2/s2  ', 'kg m/s ', 'm      ', 'kg/kg  ',            &
869             ( 'unknown', i9 = 1, 77 ) /)
870
871    CHARACTER (LEN=7), DIMENSION(300) ::  dopr_unit = 'unknown'
872
873    CHARACTER (LEN=7), DIMENSION(0:1,100) ::  do2d_unit, do3d_unit
874
875    CHARACTER (LEN=16), DIMENSION(25) ::  prt_var_names = &
876          (/ 'pt_age          ', 'pt_dvrp_size    ', 'pt_origin_x     ', &
877             'pt_origin_y     ', 'pt_origin_z     ', 'pt_radius       ', &
878             'pt_speed_x      ', 'pt_speed_y      ', 'pt_speed_z      ', &
879             'pt_weight_factor', 'pt_x            ', 'pt_y            ', &
880             'pt_z            ', 'pt_color        ', 'pt_group        ', &
881             'pt_tailpoints   ', 'pt_tail_id      ', 'pt_density_ratio', &
882             'pt_exp_arg      ', 'pt_exp_term     ', 'not_used        ', &
883             'not_used        ', 'not_used        ', 'not_used        ', &
884             'not_used        ' /)
885
886    CHARACTER (LEN=16), DIMENSION(25) ::  prt_var_units = &
887          (/ 'seconds         ', 'meters          ', 'meters          ', &
888             'meters          ', 'meters          ', 'meters          ', &
889             'm/s             ', 'm/s             ', 'm/s             ', &
890             'factor          ', 'meters          ', 'meters          ', &
891             'meters          ', 'none            ', 'none            ', &
892             'none            ', 'none            ', 'ratio           ', &
893             'none            ', 'none            ', 'not_used        ', &
894             'not_used        ', 'not_used        ', 'not_used        ', &
895             'not_used        ' /)
896
897    INTEGER ::  id_dim_prtnum, id_dim_time_pr, id_dim_time_prt, &
898                id_dim_time_pts, id_dim_time_sp, id_dim_time_ts, id_dim_x_sp, &
899                id_dim_y_sp, id_dim_zu_sp, id_dim_zw_sp, id_set_pr, &
900                id_set_prt, id_set_pts, id_set_sp, id_set_ts, id_var_prtnum, &
901                id_var_rnop_prt, id_var_time_pr, id_var_time_prt, &
902                id_var_time_pts, id_var_time_sp, id_var_time_ts, id_var_x_sp, &
903                id_var_y_sp, id_var_zu_sp, id_var_zw_sp, nc_stat
904
905    INTEGER, DIMENSION(0:1) ::  id_dim_time_xy, id_dim_time_xz, &
906                id_dim_time_yz, id_dim_time_3d, id_dim_x_xy, id_dim_xu_xy, &
907                id_dim_x_xz, id_dim_xu_xz, id_dim_x_yz, id_dim_xu_yz, &
908                id_dim_x_3d, id_dim_xu_3d, id_dim_y_xy, id_dim_yv_xy, &
909                id_dim_y_xz, id_dim_yv_xz, id_dim_y_yz, id_dim_yv_yz, &
910                id_dim_y_3d, id_dim_yv_3d, id_dim_zu_xy, id_dim_zu1_xy, &
911                id_dim_zu_xz, id_dim_zu_yz, id_dim_zu_3d, id_dim_zw_xy, &
912                id_dim_zw_xz, id_dim_zw_yz, id_dim_zw_3d, id_set_xy, &
913                id_set_xz, id_set_yz, id_set_3d, id_var_ind_x_yz, &
914                id_var_ind_y_xz, id_var_ind_z_xy, id_var_time_xy, &
915                id_var_time_xz, id_var_time_yz, id_var_time_3d, id_var_x_xy, &
916                id_var_xu_xy, id_var_x_xz, id_var_xu_xz, id_var_x_yz, &
917                id_var_xu_yz, id_var_x_3d, id_var_xu_3d, id_var_y_xy, &
918                id_var_yv_xy, id_var_y_xz, id_var_yv_xz, id_var_y_yz, &
919                id_var_yv_yz, id_var_y_3d, id_var_yv_3d, id_var_zusi_xy, &
920                id_var_zusi_3d, id_var_zu_xy, id_var_zu1_xy, id_var_zu_xz, &
921                id_var_zu_yz, id_var_zu_3d, id_var_zwwi_xy, id_var_zwwi_3d, &
922                id_var_zw_xy, id_var_zw_xz, id_var_zw_yz, id_var_zw_3d
923
924    INTEGER, DIMENSION(10)  ::  id_var_dospx, id_var_dospy, nc_precision 
925    INTEGER, DIMENSION(20)  ::  id_var_prt
926    INTEGER, DIMENSION(dopr_norm_num) ::  id_var_norm_dopr
927
928    INTEGER, DIMENSION(dopts_num,0:10) ::  id_var_dopts
929    INTEGER, DIMENSION(0:1,100)        ::  id_var_do2d, id_var_do3d
930    INTEGER, DIMENSION(100,0:9)        ::  id_dim_z_pr, id_var_dopr, &
931                                           id_var_z_pr
932    INTEGER, DIMENSION(dots_max,0:9)   ::  id_var_dots
933
934    LOGICAL ::  output_for_t0 = .FALSE.
935
936    SAVE
937
938 END MODULE netcdf_control
939
940
941
942 MODULE particle_attributes
943
944!------------------------------------------------------------------------------!
945! Description:
946! ------------
947! Definition of variables used to compute particle transport
948!------------------------------------------------------------------------------!
949
950    USE array_kind
951
952    CHARACTER (LEN=15)  ::  bc_par_lr = 'cyclic',  bc_par_ns = 'cyclic', &
953                            bc_par_b  = 'reflect', bc_par_t  = 'absorb'
954
955#if defined( __parallel )
956    INTEGER ::  mpi_particle_type
957#endif
958    INTEGER ::  ibc_par_lr, ibc_par_ns, ibc_par_b, ibc_par_t,                  &
959                iran_part = -1234567, maximum_number_of_particles = 1000,      &
960                maximum_number_of_tailpoints = 100,                            &
961                maximum_number_of_tails = 0,                                   &
962                number_of_initial_particles = 0, number_of_particles = 0,      &
963                number_of_particle_groups = 1, number_of_tails = 0,            &
964                number_of_initial_tails = 0, offset_ocean_nzt = 0,             &
965                offset_ocean_nzt_m1 = 0, particles_per_point = 1,              &
966                particle_file_count = 0, skip_particles_for_tail = 100,        &
967                total_number_of_particles, total_number_of_tails = 0
968
969    INTEGER, PARAMETER ::  max_number_of_particle_groups = 10
970
971    INTEGER, DIMENSION(:), ALLOCATABLE     ::  new_tail_id
972    INTEGER, DIMENSION(:,:,:), ALLOCATABLE ::  prt_count, prt_start_index
973
974    LOGICAL ::  particle_advection = .FALSE., random_start_position = .FALSE., &
975                read_particles_from_restartfile = .TRUE.,                      &
976                uniform_particles = .TRUE., use_particle_tails = .FALSE.,      &
977                use_sgs_for_particles = .FALSE.,                               &
978                write_particle_statistics = .FALSE.
979
980    LOGICAL, DIMENSION(max_number_of_particle_groups) ::                       &
981                vertical_particle_advection = .TRUE.
982
983    LOGICAL, DIMENSION(:), ALLOCATABLE ::  particle_mask, tail_mask
984
985    REAL    ::  c_0 = 3.0, dt_min_part = 0.0002, dt_sort_particles = 0.0,      &
986                dt_write_particle_data = 9999999.9, dvrp_psize = 9999999.9,    &
987                end_time_prel = 9999999.9, initial_weighting_factor = 1.0,     &
988                maximum_tailpoint_age = 100000.0,                              &
989                minimum_tailpoint_distance = 0.0,                              &
990                particle_advection_start = 0.0, sgs_wfu_part = 0.3333333,      &
991                sgs_wfv_part = 0.3333333, sgs_wfw_part = 0.3333333,            &
992                time_sort_particles = 0.0, time_write_particle_data = 0.0
993
994    REAL, DIMENSION(max_number_of_particle_groups) ::  &
995                density_ratio = 9999999.9, pdx = 9999999.9, pdy = 9999999.9, &
996                pdz = 9999999.9, psb = 9999999.9, psl = 9999999.9,           &
997                psn = 9999999.9, psr = 9999999.9, pss = 9999999.9,           &
998                pst = 9999999.9, radius = 9999999.9
999
1000    REAL, DIMENSION(:,:,:), ALLOCATABLE ::  particle_tail_coordinates
1001
1002
1003    TYPE particle_type
1004       SEQUENCE
1005       REAL    ::  age, age_m, dt_sum, dvrp_psize, e_m, origin_x, origin_y, &
1006                   origin_z, radius, speed_x, speed_x_sgs, speed_y,         &
1007                   speed_y_sgs, speed_z, speed_z_sgs, weight_factor, x, y, z
1008       INTEGER ::  color, group, tailpoints, tail_id
1009    END TYPE particle_type
1010
1011    TYPE(particle_type), DIMENSION(:), ALLOCATABLE ::  initial_particles, &
1012                                                       particles
1013
1014    TYPE particle_groups_type
1015       SEQUENCE
1016       REAL    ::  density_ratio, radius, exp_arg, exp_term
1017    END TYPE particle_groups_type
1018
1019    TYPE(particle_groups_type), DIMENSION(max_number_of_particle_groups) ::&
1020                   particle_groups
1021
1022    SAVE
1023
1024 END MODULE particle_attributes
1025
1026
1027
1028
1029
1030 MODULE pegrid
1031
1032!------------------------------------------------------------------------------!
1033! Description:
1034! ------------
1035! Definition of variables which define processor topology and the exchange of
1036! ghost point layers. This modules must be placed in all routines which contain
1037! MPI-calls.
1038!------------------------------------------------------------------------------!
1039
1040#if defined( __parallel )
1041#if defined( __lc )
1042    USE MPI
1043#else
1044    INCLUDE "mpif.h"
1045#endif
1046#endif
1047    CHARACTER(LEN=5)       ::  myid_char = ''
1048    INTEGER                ::  id_inflow = 0, id_recycling = 0, myid = 0,      &
1049                               target_id, npex = -1, npey = -1, numprocs = 1,  &
1050                               numprocs_previous_run = -1,                     &
1051                               tasks_per_node = -9999, threads_per_task = 1
1052
1053    INTEGER, DIMENSION(:,:), ALLOCATABLE ::  hor_index_bounds, &
1054                                             hor_index_bounds_previous_run
1055
1056#if defined( __parallel )
1057#if defined( __mpi2 )
1058    CHARACTER (LEN=MPI_MAX_PORT_NAME) ::  port_name
1059#endif
1060
1061    INTEGER ::  comm1dx, comm1dy, comm2d, comm_inter, comm_palm, ierr, myidx,  &
1062                myidy, ndim = 2, ngp_xy, ngp_y, pleft, pnorth, pright, psouth, &
1063                sendrecvcount_xy, sendrecvcount_yz, sendrecvcount_zx,          &
1064                sendrecvcount_zyd, sendrecvcount_yxd,                          &
1065                type_x, type_x_int, type_xy
1066
1067    INTEGER ::  ibuf(12), pcoord(2), pdims(2), status(MPI_STATUS_SIZE)
1068
1069    INTEGER, DIMENSION(:), ALLOCATABLE ::  ngp_yz, type_xz
1070
1071    LOGICAL ::  reorder = .TRUE.
1072    LOGICAL, DIMENSION(2) ::  cyclic = (/ .TRUE. , .TRUE. /), &
1073                              remain_dims
1074#endif
1075
1076    SAVE
1077
1078 END MODULE pegrid
1079
1080
1081
1082
1083 MODULE profil_parameter
1084
1085!------------------------------------------------------------------------------!
1086! Description:
1087! ------------
1088! Definition of variables which control PROFIL-output
1089!------------------------------------------------------------------------------!
1090
1091    INTEGER, PARAMETER ::  crmax = 100
1092
1093    CHARACTER (LEN=10), DIMENSION(100) ::  dopr_label = ' '
1094
1095    CHARACTER (LEN=10), DIMENSION(crmax) ::  cross_normalized_x = ' ', &
1096                                             cross_normalized_y = ' '
1097
1098    CHARACTER (LEN=20), DIMENSION(20) ::  cross_ts_profiles = &
1099                           (/ ' E E*               ', ' dt                 ', &
1100                              ' u* w*              ', ' th*                ', &
1101                              ' umax vmax wmax     ', ' div_old div_new    ', &
1102                              ' z_i_wpt z_i_pt     ', ' w"pt"0 w"pt" wpt   ', &
1103                              ' pt(0) pt(zp)       ', ' splux spluy spluz  ', &
1104                              ' L                  ',                         &
1105                            ( '                    ', i9 = 1, 9 ) /)
1106
1107    CHARACTER (LEN=40), DIMENSION(crmax) ::  cross_xtext = &
1108                           (/ 'windspeed in ms>->1         ', &
1109                              'pot. temperature in K       ', &
1110                              'heat flux in K ms>->1       ', &
1111                              'momentum flux in m>2s>2     ', &
1112                              'eddy diffusivity in m>2s>->1', &
1113                              'mixing length in m          ', &
1114                            ( '                            ', i9 = 1, 94 ) /)
1115
1116    CHARACTER (LEN=100), DIMENSION(crmax) ::  cross_profiles = &
1117                           (/ ' u v                           ', &
1118                              ' pt                            ', &
1119                              ' w"pt" w*pt* w*pt*BC wpt wptBC ', &
1120                              ' w"u" w*u* wu w"v" w*v* wv     ', &
1121                              ' km kh                         ', &
1122                              ' l                             ', &
1123                         ( '                               ', i9 = 1, 94 ) /)
1124
1125    INTEGER ::  profile_columns = 3, profile_rows = 2, profile_number = 0
1126
1127    INTEGER ::  cross_linecolors(100,crmax) = 1, &
1128                cross_linestyles(100,crmax) = 0, &
1129                cross_profile_numbers(100,crmax) = 0, &
1130                cross_pnc_local(crmax), cross_profile_number_count(crmax) = 0, &
1131                cross_ts_numbers(crmax,crmax) = 0, &
1132                cross_ts_number_count(crmax) = 0, dopr_crossindex(100) = 0, &
1133                dopr_index(300) = 0, dopr_initial_index(300) = 0, &
1134                dots_crossindex(100) = 0, dots_index(100) = 0, &
1135                linecolors(10) = (/ 2, 3, 4,  5, 7, 8, 12, 15, 16, 23 /), &
1136                linestyles(11) = (/ 0, 7, 3, 10, 4, 1,  9,  2,  5,  8, 6 /)
1137               
1138
1139    REAL ::  cross_normx_factor(100,crmax) = 1.0, &
1140             cross_normy_factor(100,crmax) = 1.0, &
1141             cross_ts_uymax(20) = &
1142                             (/ 999.999, 999.999, 999.999, 999.999, 999.999,   &
1143                                999.999, 999.999, 999.999, 999.999, 999.999,   &
1144                                999.999, 999.999, 999.999, 999.999, 999.999,   &
1145                                999.999, 999.999, 999.999, 999.999, 999.999 /),&
1146             cross_ts_uymax_computed(20) = 999.999, &
1147             cross_ts_uymin(20) = &
1148                             (/ 999.999, 999.999, 999.999,  -5.000, 999.999,   &
1149                                999.999,   0.000, 999.999, 999.999, 999.999,   &
1150                                999.999, 999.999, 999.999, 999.999, 999.999,   &
1151                                999.999, 999.999, 999.999, 999.999, 999.999 /),&
1152             cross_ts_uymin_computed(20) = 999.999, &
1153             cross_uxmax(crmax) = 0.0, cross_uxmax_computed(crmax) = -1.0, &
1154             cross_uxmax_normalized(crmax) = 0.0, &
1155             cross_uxmax_normalized_computed(crmax) = -1.0, &
1156             cross_uxmin(crmax) = 0.0, cross_uxmin_computed(crmax) =  1.0, &
1157             cross_uxmin_normalized(crmax) = 0.0, &
1158             cross_uxmin_normalized_computed(crmax) = 1.0, &
1159             cross_uymax(crmax), cross_uymin(crmax)
1160
1161    SAVE
1162
1163 END MODULE profil_parameter
1164
1165
1166
1167
1168 MODULE spectrum
1169
1170!------------------------------------------------------------------------------!
1171! Description:
1172! ------------
1173! Definition of quantities used for computing spectra
1174!------------------------------------------------------------------------------!
1175
1176    CHARACTER (LEN=6),  DIMENSION(1:5) ::  header_char = (/ 'PS(u) ', 'PS(v) ',&
1177                                           'PS(w) ', 'PS(pt)', 'PS(q) ' /)
1178    CHARACTER (LEN=2),  DIMENSION(10)  ::  spectra_direction = 'x'
1179    CHARACTER (LEN=10), DIMENSION(10)  ::  data_output_sp  = ' '
1180    CHARACTER (LEN=25), DIMENSION(1:5) ::  utext_char =                    &
1181                                           (/ '-power spectrum of u     ', &
1182                                              '-power spectrum of v     ', &
1183                                              '-power spectrum of w     ', &
1184                                              '-power spectrum of ^1185 ', &
1185                                              '-power spectrum of q     ' /)
1186    CHARACTER (LEN=39), DIMENSION(1:5) ::  ytext_char =                        &
1187                                 (/ 'k ^2236 ^2566^2569<u(k) in m>2s>->2    ', &
1188                                    'k ^2236 ^2566^2569<v(k) in m>2s>->2    ', &
1189                                    'k ^2236 ^2566^2569<w(k) in m>2s>->2    ', &
1190                                    'k ^2236 ^2566^2569<^1185(k) in m>2s>->2', &
1191                                    'k ^2236 ^2566^2569<q(k) in m>2s>->2    ' /)
1192
1193    INTEGER ::  klist_x = 0, klist_y = 0, n_sp_x = 0, n_sp_y = 0
1194
1195    INTEGER ::  comp_spectra_level(100) = 999999,                   &
1196                lstyles(100) = (/ 0, 7, 3, 10, 1, 4, 9, 2, 6, 8,    &
1197                                  0, 7, 3, 10, 1, 4, 9, 2, 6, 8,    &
1198                                  0, 7, 3, 10, 1, 4, 9, 2, 6, 8,    &
1199                                  0, 7, 3, 10, 1, 4, 9, 2, 6, 8,    &
1200                                  0, 7, 3, 10, 1, 4, 9, 2, 6, 8,    &
1201                                  0, 7, 3, 10, 1, 4, 9, 2, 6, 8,    &
1202                                  0, 7, 3, 10, 1, 4, 9, 2, 6, 8,    &
1203                                  0, 7, 3, 10, 1, 4, 9, 2, 6, 8,    &
1204                                  0, 7, 3, 10, 1, 4, 9, 2, 6, 8,    &
1205                                  0, 7, 3, 10, 1, 4, 9, 2, 6, 8 /), &
1206                plot_spectra_level(100) = 999999
1207
1208    REAL    ::  time_to_start_sp = 0.0
1209
1210    SAVE
1211
1212 END MODULE spectrum
1213
1214
1215
1216
1217 MODULE statistics
1218
1219!------------------------------------------------------------------------------!
1220! Description:
1221! ------------
1222! Definition of statistical quantities, e.g. global sums
1223!------------------------------------------------------------------------------!
1224
1225    CHARACTER (LEN=40) ::  region(0:9)
1226    INTEGER ::  pr_palm = 80, statistic_regions = 0, var_ts = 100
1227    INTEGER ::  u_max_ijk(3), v_max_ijk(3), w_max_ijk(3)
1228    LOGICAL ::  flow_statistics_called = .FALSE.
1229    REAL ::     u_max, v_max, w_max
1230    REAL, DIMENSION(:), ALLOCATABLE       ::  sums_divnew_l, sums_divold_l
1231    REAL, DIMENSION(:,:), ALLOCATABLE     ::  sums, sums_wsts_bc_l, ts_value
1232    REAL, DIMENSION(:,:,:), ALLOCATABLE   ::  hom_sum, rmask, spectrum_x, &
1233                                              spectrum_y, sums_l, sums_l_l, &
1234                                              sums_up_fraction_l
1235    REAL, DIMENSION(:,:,:,:), ALLOCATABLE ::  hom
1236
1237    SAVE
1238
1239 END MODULE statistics
1240
1241
1242
1243
1244 MODULE transpose_indices
1245
1246!------------------------------------------------------------------------------!
1247! Description:
1248! ------------
1249! Definition of indices for transposed arrays
1250!------------------------------------------------------------------------------!
1251
1252    INTEGER ::  nxl_y, nxl_yd, nxl_z, nxr_y, nxr_ya, nxr_yd, nxr_yda, nxr_z, &
1253                nxr_za, nyn_x, nyn_xa, nyn_z, nyn_za, nys_x, nys_z, nzb_x,   &
1254                nzb_y, nzb_yd, nzt_x, nzt_xa, nzt_y, nzt_ya, nzt_yd, nzt_yda
1255               
1256
1257    SAVE
1258
1259 END MODULE transpose_indices
Note: See TracBrowser for help on using the repository browser.