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

Last change on this file since 2108 was 2108, checked in by kanani, 7 years ago

last commit documented

  • Property svn:keywords set to Id
File size: 100.2 KB
Line 
1!> @file modules.f90
2!------------------------------------------------------------------------------!
3! This file is part of PALM.
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-2017 Leibniz Universitaet Hannover
18!------------------------------------------------------------------------------!
19!
20! Current revisions:
21! ------------------
22!
23!
24! Former revisions:
25! -----------------
26! $Id: modules.f90 2108 2017-01-09 12:23:20Z kanani $
27!
28! 2107 2017-01-09 12:21:49Z kanani
29! Preparation for doxygen comments (Giersch)
30!
31! 2050 2016-11-08 15:00:55Z gronemeier
32! Implement turbulent outflow condition
33!
34! 2037 2016-10-26 11:15:40Z knoop
35! Anelastic approximation implemented
36!
37! 2031 2016-10-21 15:11:58Z knoop
38! renamed variable rho to rho_ocean and rho_av to rho_ocean_av
39!
40! 2011 2016-09-19 17:29:57Z kanani
41! +urban_surface, +lsf_exception, +varnamelength
42!
43! 2007 2016-08-24 15:47:17Z kanani
44! Increased DIMENSION of data_output, data_output_user, do2d, do3d
45!
46! 2000 2016-08-20 18:09:15Z knoop
47! Forced header and separation lines into 80 columns
48!
49! 1992 2016-08-12 15:14:59Z suehring
50! +constant_top_scalarflux, top_scalarflux
51! default of bc_s_t adjusted
52!
53! 1968 2016-07-18 12:01:49Z suehring
54! Changed dimension for MPI-datatypes type_x_int and type_y_int
55!
56! 1960 2016-07-12 16:34:24Z suehring
57! Separate humidity and passive scalar
58! +bc_s_t_val, diss_s_s, diss_l_s, flux_s_s, flux_l_s, s, sp, s1, s2, s3, ssws_av,
59!  s_init, s_surf, sums_wsss_ws_l, ss, ssws, sswst, ts_m, wall_sflux
60! +constant_scalarflux, ibc_s_b, ibc_s_t, s_vertical_gradient_level_ind
61!
62! Unused variables removed
63! -gamma_x, gamma_y, gamma_z, var_x, var_y, var_z
64!
65! Change initial values (in order to unify gradient calculation):
66! pt_vertical_gradient_level, sa_vertical_gradient_level
67!
68! 1957 2016-07-07 10:43:48Z suehring
69! +fl_max, num_leg, num_var_fl, num_var_fl_user, var_fl_max, virtual_flight
70!
71! 1918 2016-05-27 14:35:57Z raasch
72! default timestep switched from -1.0 to +1.0 in order to avoid wrong sign of
73! initially calculated divergence
74!
75! 1906 2016-05-24 14:38:08Z suehring
76! default value of mg_switch_to_pe0_level changed to -1
77!
78! 1849 2016-04-08 11:33:18Z hoffmann
79! bfactor, mass_of_solute, molecular_weight_of_solute, molecular_weight_of_water,
80! vanthoff moved to mod_particle_attributes.
81! dt_micro and several cloud_parameters moved to microphysics_mod.
82! 1d-microphysics profiles moved to microphysics_mod.
83!
84! 1845 2016-04-08 08:29:13Z raasch
85! -nzb_2d
86!
87! 1833 2016-04-07 14:23:03Z raasch
88! spectra parameter moved to spectra module
89!
90! 1831 2016-04-07 13:15:51Z hoffmann
91! curvature_solution_effects removed
92! turbulence renamed collision_turbulence, drizzle renamed
93! cloud_water_sedimentation
94!
95! 1822 2016-04-07 07:49:42Z hoffmann
96! icloud_scheme removed. microphysics_sat_adjust, microphysics_kessler,
97! microphysics_seifert added.
98!
99! 1815 2016-04-06 13:49:59Z raasch
100! cpp-directive for decalpha removed
101!
102! 1808 2016-04-05 19:44:00Z raasch
103! MPI module used by default on all machines
104!
105! 1804 2016-04-05 16:30:18Z maronga
106! Removed code for parameter file check (__check)
107!
108! 1788 2016-03-10 11:01:04Z maronga
109! Added roughness length for moisture (z0q)
110!
111! 1786 2016-03-08 05:49:27Z raasch
112! module spectrum moved to new separate module
113!
114! 1783 2016-03-06 18:36:17Z raasch
115! netcdf variables moved to the netcdf-interface module
116!
117! 1779 2016-03-03 08:01:28Z raasch
118! coupling_char extended to LEN=3
119!
120! 1764 2016-02-28 12:45:19Z raasch
121! some reformatting
122!
123! 1762 2016-02-25 12:31:13Z hellstea
124! +nest_* variables, size of volume_flow arrays increased by one element
125!
126! 1738 2015-12-18 13:56:05Z raasch
127! +mean_surface_level_height
128!
129! 1695 2015-10-27 10:03:11Z maronga
130! Removed rif (forgotten in last revision)
131!
132! 1693 2015-10-27 08:35:45Z maronga
133! Renamed zp -> z_mo
134!
135! 1691 2015-10-26 16:17:44Z maronga
136! Renamed Obukhov length. Added ol, removed rif. Increased number of profiles
137! (pr_palm). Changed default values for dissipation_1d to 'detering' and
138! (mixing_length_1d to 'blackadar'. Added most_method. rif_min and rif_max
139! renamed to zeta_min and zeta_max and new values assigned.
140!
141! 1682 2015-10-07 23:56:08Z knoop
142! Code annotations made doxygen readable
143!
144! 1677 2015-10-02 13:25:23Z boeske
145! +ngp_yz_int, type_xz_int, type_yz_int
146!
147! 1666 2015-09-23 07:31:10Z raasch
148! +user_interface_current_revision, user_interface_required_revision in
149! control_parameters
150!
151! 1639 2015-08-31 14:46:48Z knoop
152! Bugfix: string 'unknown' extended to match LEN=13
153!
154! 1575 2015-03-27 09:56:27Z raasch
155! +ngp_xz
156!
157! 1560 2015-03-06 10:48:54Z keck
158! +recycling_yshift
159!
160! 1557 2015-03-05 16:43:04Z suehring
161! +monotonic_adjustment
162!
163! 1551 2015-03-03 14:18:16Z maronga
164! Increased pr_palm to 120. Increased length of dots_unit and dots_label to 13
165! digits. Increased length of domask, do2d, and do3d to 20 digits.
166!
167! 1496 2014-12-02 17:25:50Z maronga
168! Renamed "radiation" -> "cloud_top_radiation"
169!
170! 1484 2014-10-21 10:53:05Z kanani
171! Changes due to new module structure of the plant canopy model:
172!   canopy-model related parameters/variables moved to module
173!   plant_canopy_model_mod
174!
175! 1468 2014-09-24 14:06:57Z maronga
176! Adapted for use on up to 6-digit processor cores.
177! Increased identifier string length for user-defined quantities to 20.
178!
179! 1450 2014-08-21 07:31:51Z heinze
180! ntnudge from 100 to 1000 increased to allow longer simulations
181!
182! 1431 2014-07-15 14:47:17Z suehring
183! +var_d
184!
185! 1429 2014-07-15 12:53:45Z knoop
186! +ensemble_member_nr to prepare the random_generator for ensemble runs
187!
188! 1382 2014-04-30 12:15:41Z boeske
189! Renamed variables which store large scale forcing tendencies
190! pt_lsa -> td_lsa_lpt, pt_subs -> td_sub_lpt,
191! q_lsa  -> td_lsa_q,   q_subs  -> td_sub_q
192!
193! 1365 2014-04-22 15:03:56Z boeske
194! Usage of large scale forcing enabled:
195! increase pr_palm from 90 to 100 to allow for more standard profiles
196! + ngp_sums_ls, pt_lsa, pt_subs, q_lsa, q_subs, tmp_tnudge, sums_ls_l,
197! use_subsidence_tendencies
198!
199! 1361 2014-04-16 15:17:48Z hoffmann
200! tend_* removed
201! call_microphysics_at_all_substeps added
202! default of drizzle set to true
203!
204! 1359 2014-04-11 17:15:14Z hoffmann
205! particle_attributes moved to mod_particle_attributes.f90
206!
207! 1353 2014-04-08 15:21:23Z heinze
208! REAL constants provided with KIND-attribute
209!
210! 1327 2014-03-21 11:00:16Z raasch
211! REAL constants defined as wp-kind
212! -avs_output, data_output_format, do3d_compress, iso2d_output, netcdf_output
213!
214! 1320 2014-03-20 08:40:49Z raasch
215! ONLY-attribute added to USE-statements,
216! kind-parameters added to all INTEGER and REAL declaration statements,
217! kinds are defined in new module kinds,
218! old module precision_kind is removed,
219! revision history before 2012 removed,
220! comment fields (!:) to be used for variable explanations added to
221! all variable declaration statements
222!
223! 1318 2014-03-17 13:35:16Z raasch
224! module cpulog moved to new separate module-file
225! interface for cpu_log removed
226!
227! 1314 2014-03-14 18:25:17Z suehring
228! + log_z_z0, number_of_sublayers, z0_av_global
229! 1308 2014-03-13 14:58:42Z fricke
230! +ntdim_2d_xy, ntdim_2d_xz, ntdim_2d_yz, ntdim_3d
231!
232! 1257 2013-11-08 15:18:40Z raasch
233! set default values for grid indices of maximum velocity components
234! u|v|w_max_ijk
235!
236! 1241 2013-10-30 11:36:58Z heinze
237! Usage of nudging enabled
238! +nudging, ntnudge, ptnudge, qnudge, tnudge, unudge, vnudge, wnudge
239! increase pr_palm from 80 to 90 to allow for more standard profiles
240!
241! Enable prescribed time depenend surface fluxes and geostrophic wind read in
242! from external file LSF_DATA
243! +large_scale_forcing, lsf_surf, lsf_vert, nlsf, time_surf, shf_surf, qsws_surf,
244!  pt_surf, q_surf, p_surf, time_vert, ug_vert, vg_vert, wsubs_vert
245!
246! 1221 2013-09-10 08:59:13Z raasch
247! wall_flags_0 changed to 32bit int, +wall_flags_00,
248! +rflags_s_inner, rflags_invers
249!
250! 1216 2013-08-26 09:31:42Z raasch
251! +transpose_compute_overlap,
252! several variables are now defined in the serial (non-parallel) case also
253!
254! 1212 2013-08-15 08:46:27Z raasch
255! +tri
256!
257! 1179 2013-06-14 05:57:58Z raasch
258! +reference_state, ref_state, use_initial_profile_as_reference, vpt_reference,
259! use_reference renamed use_single_reference_value
260!
261! 1159 2013-05-21 11:58:22Z fricke
262! -bc_lr_dirneu, bc_lr_neudir, bc_ns_dirneu, bc_ns_neudir
263! +use_cmax
264!
265! 1128 2013-04-12 06:19:32Z raasch
266! +background_communication, i_left, i_right, j_north, j_south, req, req_count,
267! send_receive, sendrecv_in_background, wait_stat
268!
269! 1115 2013-03-26 18:16:16Z hoffmann
270! unused variables removed
271!
272! 1113 2013-03-10 02:48:14Z raasch
273! +on_device
274!
275! 1111 2013-03-08 23:54:10Z raasch
276! +tric, nr_timesteps_this_run
277!
278! 1106 2013-03-04 05:31:38Z raasch
279! array_kind renamed precision_kind, pdims defined in serial code
280! bugfix: default value assigned to coupling_start_time
281!
282! 1095 2013-02-03 02:21:01Z raasch
283! FORTRAN error in r1092 removed
284!
285! 1092 2013-02-02 11:24:22Z raasch
286! character length in some derived type changed for better alignment
287!
288! 1065 2012-11-22 17:42:36Z hoffmann
289! + c_sedimentation, limiter_sedimentation, turbulence, a_1, a_2, a_3, b_1, b_2,
290! + b_3, c_1, c_2, c_3, beta_cc
291!
292! bottom boundary condition of qr, nr changed from Dirichlet to Neumann
293!
294! 1053 2012-11-13 17:11:03Z hoffmann
295! necessary expansions according to the two new prognostic equations (nr, qr)
296! of the two-moment cloud physics scheme:
297! +*_init, flux_l_*, diss_l_*, flux_s_*, diss_s_*, *sws, *swst, tend_*, *, *_p
298! +t*_m, *_1, *_2, *_3, *_av, bc_*_b, bc_*_t, ibc_*_b, ibc_*_t, bc_*_t_val,
299! +*_vertical_gradient, *_surface_initial_change, *_vertical_gradient_level,
300! +*_vertical_gradient_level_ind, *_surface, constant_waterflux_*, 
301! +cloud_scheme, icloud_scheme, surface_waterflux_*, sums_ws*s_ws_l, wall_*flux
302!
303! constants for the two-moment scheme:
304! +a_vent, a_term, b_vent, b_term, c_evap, c_term, cof, eps_sb, k_cc, k_cr, k_rr,
305! +k_br, kappa_rr, kin_vis_air, mu_constant_value, nc, pirho_l, dpirho_l, rho_1,
306! +schmidt, schmidt_p_1d3, stp, x0, xmin, xmax, dt_precipitation, w_precipitation
307!
308! steering parameters for the two_moment scheme:
309! +mu_constant, ventilation_effect
310!
311! 1036 2012-10-22 13:43:42Z raasch
312! code put under GPL (PALM 3.9)
313!
314! 1031 2012-10-19 14:35:30Z raasch
315! +output_format_netcdf
316!
317! 1015 2012-09-27 09:23:24Z raasch
318! +acc_rank, num_acc_per_node,
319! -adjust_mixing_length
320!
321! 1010 2012-09-20 07:59:54Z raasch
322! pointer free version can be generated with cpp switch __nopointer
323!
324! 1003 2012-09-14 14:35:53Z raasch
325! -grid_matching, nxa, nya, etc., nnx_pe, nny_pe, spl_*
326!
327! 1001 2012-09-13 14:08:46Z raasch
328! -asselin_filter_factor, cut_spline_overshoot, dt_changed, last_dt_change,
329! last_dt_change_1d, long_filter_factor, overshoot_limit_*, ups_limit_*
330! several pointer/target arrays converted to normal ones
331!
332! 996 2012-09-07 10:41:47Z raasch
333! -use_prior_plot1d_parameters
334!
335! 978 2012-08-09 08:28:32Z fricke
336! +c_u_m, c_u_m_l, c_v_m, c_v_m_l, c_w_m, c_w_m_l,
337! +bc_lr_dirneu, bc_lr_neudir, bc_ns_dirneu, bc_ns_neudir
338! -km_damp_x, km_damp_y, km_damp_max, outflow_damping_width
339! +z0h, z0h_av, z0h_factor, z0h1d
340! +ptdf_x, ptdf_y, pt_damping_width, pt_damping_factor
341!
342! 964 2012-07-26 09:14:24Z raasch
343! -cross_linecolors, cross_linestyles, cross_normalized_x, cross_normx_factor,
344! cross_normalized_y, cross_normy_factor, cross_pnc_local,
345! cross_profile_numbers, cross_profile_number_counter, cross_uxmax,
346! cross_uxmax_computed, cross_uxmax_normalized,
347! cross_uxmax_normalized_computed, cross_uxmin, cross_uxmin_computed,
348! cross_uxmin_normalized, cross_uxmin_normalized_computed, cross_uymax,
349! cross_uymin, cross_xtext, dopr_crossindex, dopr_label, linecolors, linestyles,
350! nz_do1d, profil_output, z_max_do1d, z_max_do1d_normalized
351!
352! 951 2012-07-19 14:22:52Z hoffmann
353! changing profile_columns and profile_rows
354!
355! 940 2012-07-09 14:31:00Z raasch
356! +neutral
357!
358! 927 2012-06-06 19:15:04Z raasch
359! +masking_method
360!
361! 880 2012-04-13 06:28:59Z raasch
362! gathered_size, subdomain_size moved to control_parameters
363!
364! 866 2012-03-28 06:44:41Z raasch
365! interface for global_min_max changed
366!
367! 861 2012-03-26 14:18:34Z suehring
368! +wall_flags_0
369! -boundary_flags
370! +nzb_max
371! +adv_sca_1, +adv_mom_1
372!
373! 849 2012-03-15 10:35:09Z raasch
374! +deleted_particles, deleted_tails, tr.._count_sum, tr.._count_recv_sum in
375! particle_attributes,
376! +de_dx, de_dy, de_dz in arrays_3d,
377! first_call_advec_particles renamed first_call_lpm
378!
379! 828 2012-02-21 12:00:36Z raasch
380! +dissipation_classes, radius_classes, use_kernel_tables,
381! particle feature color renamed class
382!
383! 825 2012-02-19 03:03:44Z raasch
384! +bfactor, curvature_solution_effects, eps_ros, molecular_weight_of_water,
385! vanthoff, -b_cond in cloud_parameters,
386! dopts_num increased to 29, particle attributes speed_x|y|z_sgs renamed
387! rvar1|2|3
388! wang_collision_kernel and turbulence_effects_on_collision replaced by
389! collision_kernel, hall_kernel, palm_kernel, wang_kernel
390!
391! 809 2012-01-30 13:32:58Z marongas
392! Bugfix: replaced .AND. and .NOT. with && and ! in the preprocessor directives
393!
394! 807 2012-01-25 11:53:51Z maronga
395! New cpp directive "__check" implemented which is used by check_namelist_files.
396! New parameter check_restart has been defined which is needed by
397! check_namelist_files only.
398!
399! 805 2012-01-17 15:53:28Z franke
400! Bugfix collective_wait must be out of parallel branch for runs in serial mode
401!
402! 801 2012-01-10 17:30:36Z suehring
403! Dimesion of sums_wsus_ws_l, ! sums_wsvs_ws_l, sums_us2_ws_l, sums_vs2_ws_l,
404! sums_ws2_ws_l, sums_wspts_ws_l, sums_wsqs_ws_l, sums_wssas_ws_l increased.
405! for thread-safe summation in advec_ws.
406!
407! RCS Log replace by Id keyword, revision history cleaned up
408!
409! Revision 1.95  2007/02/11 13:18:30  raasch
410! version 3.1b (last under RCS control)
411!
412! Revision 1.1  1997/07/24 11:21:26  raasch
413! Initial revision
414!
415!
416!------------------------------------------------------------------------------!
417! Description:
418! ------------
419!> Definition of all variables
420!>
421!> @todo Add description for each variable
422!------------------------------------------------------------------------------!
423
424
425!------------------------------------------------------------------------------!
426! Description:
427! ------------
428!> Definition of variables for special advection schemes.
429!------------------------------------------------------------------------------!
430 MODULE advection
431 
432    USE kinds
433
434    REAL(wp), DIMENSION(:), ALLOCATABLE ::  aex  !<
435    REAL(wp), DIMENSION(:), ALLOCATABLE ::  bex  !<
436    REAL(wp), DIMENSION(:), ALLOCATABLE ::  dex  !<
437    REAL(wp), DIMENSION(:), ALLOCATABLE ::  eex  !<
438   
439    SAVE
440
441 END MODULE advection
442
443
444!------------------------------------------------------------------------------!
445! Description:
446! ------------
447!> Definition of all arrays defined on the computational grid.
448!------------------------------------------------------------------------------!
449 MODULE arrays_3d
450
451    USE kinds
452
453    REAL(wp), DIMENSION(:), ALLOCATABLE ::  c_u_m                  !<
454    REAL(wp), DIMENSION(:), ALLOCATABLE ::  c_u_m_l                !<
455    REAL(wp), DIMENSION(:), ALLOCATABLE ::  c_v_m                  !<
456    REAL(wp), DIMENSION(:), ALLOCATABLE ::  c_v_m_l                !<
457    REAL(wp), DIMENSION(:), ALLOCATABLE ::  c_w_m                  !<
458    REAL(wp), DIMENSION(:), ALLOCATABLE ::  c_w_m_l                !<   
459    REAL(wp), DIMENSION(:), ALLOCATABLE ::  ddzu                   !<
460    REAL(wp), DIMENSION(:), ALLOCATABLE ::  ddzu_pres              !<
461    REAL(wp), DIMENSION(:), ALLOCATABLE ::  dd2zu                  !<
462    REAL(wp), DIMENSION(:), ALLOCATABLE ::  dzu                    !<
463    REAL(wp), DIMENSION(:), ALLOCATABLE ::  ddzw                   !<
464    REAL(wp), DIMENSION(:), ALLOCATABLE ::  dzw                    !<
465    REAL(wp), DIMENSION(:), ALLOCATABLE ::  hyp                    !<
466    REAL(wp), DIMENSION(:), ALLOCATABLE ::  inflow_damping_factor  !<
467    REAL(wp), DIMENSION(:), ALLOCATABLE ::  l_grid                 !<
468    REAL(wp), DIMENSION(:), ALLOCATABLE ::  ptdf_x                 !<
469    REAL(wp), DIMENSION(:), ALLOCATABLE ::  ptdf_y                 !<
470    REAL(wp), DIMENSION(:), ALLOCATABLE ::  p_surf                 !<
471    REAL(wp), DIMENSION(:), ALLOCATABLE ::  pt_surf                !<
472    REAL(wp), DIMENSION(:), ALLOCATABLE ::  pt_init                !<
473    REAL(wp), DIMENSION(:), ALLOCATABLE ::  qsws_surf              !<
474    REAL(wp), DIMENSION(:), ALLOCATABLE ::  q_init                 !<
475    REAL(wp), DIMENSION(:), ALLOCATABLE ::  q_surf                 !<
476    REAL(wp), DIMENSION(:), ALLOCATABLE ::  rdf                    !<
477    REAL(wp), DIMENSION(:), ALLOCATABLE ::  rdf_sc                 !<
478    REAL(wp), DIMENSION(:), ALLOCATABLE ::  ref_state              !<
479    REAL(wp), DIMENSION(:), ALLOCATABLE ::  s_init                 !<
480    REAL(wp), DIMENSION(:), ALLOCATABLE ::  s_surf                 !<
481    REAL(wp), DIMENSION(:), ALLOCATABLE ::  sa_init                !<
482    REAL(wp), DIMENSION(:), ALLOCATABLE ::  shf_surf               !<
483    REAL(wp), DIMENSION(:), ALLOCATABLE ::  timenudge              !<
484    REAL(wp), DIMENSION(:), ALLOCATABLE ::  time_surf              !<
485    REAL(wp), DIMENSION(:), ALLOCATABLE ::  time_vert              !<
486    REAL(wp), DIMENSION(:), ALLOCATABLE ::  tmp_tnudge             !<
487    REAL(wp), DIMENSION(:), ALLOCATABLE ::  ug                     !<
488    REAL(wp), DIMENSION(:), ALLOCATABLE ::  u_init                 !<
489    REAL(wp), DIMENSION(:), ALLOCATABLE ::  u_nzb_p1_for_vfc       !<
490    REAL(wp), DIMENSION(:), ALLOCATABLE ::  vg                     !<
491    REAL(wp), DIMENSION(:), ALLOCATABLE ::  v_init                 !<
492    REAL(wp), DIMENSION(:), ALLOCATABLE ::  v_nzb_p1_for_vfc       !<
493    REAL(wp), DIMENSION(:), ALLOCATABLE ::  w_subs                 !<
494    REAL(wp), DIMENSION(:), ALLOCATABLE ::  zu                     !<
495    REAL(wp), DIMENSION(:), ALLOCATABLE ::  zw                     !<
496
497    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  c_u                   !<
498    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  c_v                   !<
499    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  c_w                   !<
500    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  diss_s_e              !<
501    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  diss_s_nr             !<     
502    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  diss_s_pt             !<
503    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  diss_s_q              !<
504    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  diss_s_qr             !<
505    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  diss_s_s              !<
506    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  diss_s_sa             !<
507    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  diss_s_u              !<
508    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  diss_s_v              !<
509    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  diss_s_w              !<
510    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  dzu_mg                !<
511    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  dzw_mg                !<
512    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  flux_s_e              !<
513    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  flux_s_nr             !<
514    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  flux_s_pt             !<
515    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  flux_s_q              !<
516    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  flux_s_qr             !<
517    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  flux_s_s              !<
518    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  flux_s_sa             !<
519    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  flux_s_u              !<
520    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  flux_s_v              !<
521    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  flux_s_w              !<
522    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  f1_mg                 !<
523    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  f2_mg                 !<
524    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  f3_mg                 !<
525    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  mean_inflow_profiles  !<
526    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  nrs                   !<
527    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  nrsws                 !<
528    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  nrswst                !<
529    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  ol                    !< Obukhov length
530    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  precipitation_amount  !<
531    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  precipitation_rate    !<
532    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  ptnudge               !<
533    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  pt_slope_ref          !<
534    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  qnudge                !<
535    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  qs                    !<
536    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  qsws                  !<
537    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  qswst                 !<
538    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  qswst_remote          !<
539    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  qrs                   !<
540    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  qrsws                 !<
541    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  qrswst                !<
542    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  saswsb                !<
543    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  saswst                !<
544    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  shf                   !<
545    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  ss                    !<
546    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  ssws                  !<
547    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  sswst                 !<
548    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  tnudge                !<
549    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  td_lsa_lpt            !<
550    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  td_lsa_q              !<
551    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  td_sub_lpt            !<
552    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  td_sub_q              !<
553    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  total_2d_a            !<
554    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  total_2d_o            !<
555    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  ts                    !<
556    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  tswst                 !<
557    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  ug_vert               !<
558    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  unudge                !<
559    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  us                    !<
560    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  usws                  !<
561    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  uswst                 !<
562    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  vnudge                !<
563    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  vg_vert               !<
564    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  vsws                  !<
565    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  vswst                 !<
566    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  wnudge                !<
567    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  wsubs_vert            !<
568    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  z0                    !< roughness length for momentum
569    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  z0h                   !< roughness length for heat
570    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  z0q                   !< roughness length for moisture
571   
572    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  d          !<
573    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  de_dx      !<
574    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  de_dy      !<
575    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  de_dz      !<
576    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  diss       !<
577    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  diss_l_e   !<
578    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  diss_l_nr  !<
579    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  diss_l_pt  !<
580    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  diss_l_q   !<
581    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  diss_l_qr  !<
582    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  diss_l_s   !<
583    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  diss_l_sa  !<
584    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  diss_l_u   !<
585    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  diss_l_v   !<
586    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  diss_l_w   !<
587    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  flux_l_e   !<
588    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  flux_l_nr  !<
589    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  flux_l_pt  !<
590    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  flux_l_q   !<
591    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  flux_l_qr  !<
592    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  flux_l_s   !<
593    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  flux_l_sa  !<
594    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  flux_l_u   !<
595    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  flux_l_v   !<
596    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  flux_l_w   !<
597    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  kh         !<
598    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  km         !<
599    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  l_wall     !<
600    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  prr        !<
601    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  p_loc      !<
602    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  tend       !<
603    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  tric       !<
604    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  u_m_l      !<
605    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  u_m_n      !<
606    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  u_m_r      !<
607    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  u_m_s      !<
608    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  v_m_l      !<
609    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  v_m_n      !<
610    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  v_m_r      !<
611    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  v_m_s      !<
612    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  w_m_l      !<
613    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  w_m_n      !< 
614    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  w_m_r      !<
615    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  w_m_s      !<
616
617#if defined( __nopointer )
618    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  e          !<
619    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  e_p        !<
620    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  nr         !<
621    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  nr_p       !<
622    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  p          !<
623    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  prho       !<
624    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  pt         !<
625    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  pt_p       !<
626    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  q          !<
627    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  q_p        !<
628    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  qc         !<
629    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  ql         !<
630    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  ql_c       !<
631    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  ql_v       !<
632    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  ql_vp      !<
633    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  qr         !<
634    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  qr_p       !<
635    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  rho_ocean  !<
636    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  s          !<
637    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  s_p        !<
638    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  sa         !<
639    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  sa_p       !<
640    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  te_m       !<
641    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  tnr_m      !<
642    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  tpt_m      !<
643    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  tq_m       !<
644    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  tqr_m      !<
645    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  ts_m       !<
646    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  tsa_m      !<
647    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  tu_m       !<
648    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  tv_m       !<
649    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  tw_m       !< 
650    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  u          !<
651    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  u_p        !<
652    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  v          !<
653    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  v_p        !<
654    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  vpt        !<
655    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  w          !<
656    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  w_p        !<
657#else
658    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  e_1     !<
659    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  e_2     !<
660    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  e_3     !<
661    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  p       !<
662    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  prho_1  !<
663    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  nr_1    !<
664    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  nr_2    !<
665    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  nr_3    !<
666    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  pt_1    !<
667    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  pt_2    !<
668    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  pt_3    !<
669    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  q_1     !<
670    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  q_2     !<
671    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  q_3     !<
672    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  qc_1    !<
673    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  ql_v    !<
674    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  ql_vp   !<
675    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  ql_1    !<
676    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  ql_2    !<
677    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  qr_1    !<
678    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  qr_2    !<
679    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  qr_3    !<
680    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  rho_1   !<
681    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  s_1     !<
682    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  s_2     !<
683    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  s_3     !<
684    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  sa_1    !<
685    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  sa_2    !<
686    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  sa_3    !<
687    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  u_1     !<
688    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  u_2     !<
689    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  u_3     !<
690    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  v_1     !<
691    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  v_2     !<
692    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  v_3     !<
693    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  vpt_1   !<
694    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  w_1     !<
695    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  w_2     !<
696    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  w_3     !<
697
698    REAL(wp), DIMENSION(:,:,:), POINTER ::  e          !<
699    REAL(wp), DIMENSION(:,:,:), POINTER ::  e_p        !<
700    REAL(wp), DIMENSION(:,:,:), POINTER ::  nr         !<
701    REAL(wp), DIMENSION(:,:,:), POINTER ::  nr_p       !<
702    REAL(wp), DIMENSION(:,:,:), POINTER ::  prho       !<
703    REAL(wp), DIMENSION(:,:,:), POINTER ::  pt         !<
704    REAL(wp), DIMENSION(:,:,:), POINTER ::  pt_p       !<
705    REAL(wp), DIMENSION(:,:,:), POINTER ::  q          !<
706    REAL(wp), DIMENSION(:,:,:), POINTER ::  q_p        !<
707    REAL(wp), DIMENSION(:,:,:), POINTER ::  qc         !<
708    REAL(wp), DIMENSION(:,:,:), POINTER ::  ql         !<
709    REAL(wp), DIMENSION(:,:,:), POINTER ::  ql_c       !<
710    REAL(wp), DIMENSION(:,:,:), POINTER ::  qr         !<
711    REAL(wp), DIMENSION(:,:,:), POINTER ::  qr_p       !<
712    REAL(wp), DIMENSION(:,:,:), POINTER ::  rho_ocean  !<
713    REAL(wp), DIMENSION(:,:,:), POINTER ::  s          !<
714    REAL(wp), DIMENSION(:,:,:), POINTER ::  s_p        !<
715    REAL(wp), DIMENSION(:,:,:), POINTER ::  sa         !<
716    REAL(wp), DIMENSION(:,:,:), POINTER ::  sa_p       !<
717    REAL(wp), DIMENSION(:,:,:), POINTER ::  te_m       !<
718    REAL(wp), DIMENSION(:,:,:), POINTER ::  tnr_m      !<
719    REAL(wp), DIMENSION(:,:,:), POINTER ::  tpt_m      !<
720    REAL(wp), DIMENSION(:,:,:), POINTER ::  tq_m       !<
721    REAL(wp), DIMENSION(:,:,:), POINTER ::  tqr_m      !<
722    REAL(wp), DIMENSION(:,:,:), POINTER ::  ts_m       !<
723    REAL(wp), DIMENSION(:,:,:), POINTER ::  tsa_m      !<
724    REAL(wp), DIMENSION(:,:,:), POINTER ::  tu_m       !<
725    REAL(wp), DIMENSION(:,:,:), POINTER ::  tv_m       !<
726    REAL(wp), DIMENSION(:,:,:), POINTER ::  tw_m       !<
727    REAL(wp), DIMENSION(:,:,:), POINTER ::  u          !<
728    REAL(wp), DIMENSION(:,:,:), POINTER ::  u_p        !<
729    REAL(wp), DIMENSION(:,:,:), POINTER ::  v          !<
730    REAL(wp), DIMENSION(:,:,:), POINTER ::  v_p        !<
731    REAL(wp), DIMENSION(:,:,:), POINTER ::  vpt        !<
732    REAL(wp), DIMENSION(:,:,:), POINTER ::  w          !<
733    REAL(wp), DIMENSION(:,:,:), POINTER ::  w_p        !<
734#endif
735
736    REAL(wp), DIMENSION(:,:,:,:), ALLOCATABLE ::  rif_wall !<
737    REAL(wp), DIMENSION(:,:,:,:), ALLOCATABLE ::  tri      !<
738
739    REAL(wp), DIMENSION(:), ALLOCATABLE ::  rho_air      !< air density profile on the uv grid
740    REAL(wp), DIMENSION(:), ALLOCATABLE ::  rho_air_zw   !< air density profile on the w grid
741    REAL(wp), DIMENSION(:), ALLOCATABLE ::  drho_air     !< inverse air density profile on the uv grid
742    REAL(wp), DIMENSION(:), ALLOCATABLE ::  drho_air_zw  !< inverse air density profile on the w grid
743
744    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  rho_air_mg     !< air density profiles on the uv grid for multigrid
745    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  rho_air_zw_mg  !< air density profiles on the w grid for multigrid
746
747    REAL(wp), DIMENSION(:), ALLOCATABLE ::  heatflux_input_conversion       !< conversion factor array for heatflux input
748    REAL(wp), DIMENSION(:), ALLOCATABLE ::  waterflux_input_conversion      !< conversion factor array for waterflux input
749    REAL(wp), DIMENSION(:), ALLOCATABLE ::  momentumflux_input_conversion   !< conversion factor array for momentumflux input
750    REAL(wp), DIMENSION(:), ALLOCATABLE ::  heatflux_output_conversion      !< conversion factor array for heatflux output
751    REAL(wp), DIMENSION(:), ALLOCATABLE ::  waterflux_output_conversion     !< conversion factor array for waterflux output
752    REAL(wp), DIMENSION(:), ALLOCATABLE ::  momentumflux_output_conversion  !< conversion factor array for momentumflux output
753
754    SAVE
755
756 END MODULE arrays_3d
757
758
759!------------------------------------------------------------------------------!
760! Description:
761! ------------
762!> Definition of variables needed for time-averaging of 2d/3d data.
763!------------------------------------------------------------------------------!
764 MODULE averaging
765 
766    USE kinds
767
768    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  lwp_av                 !< Avg. liquid water path
769    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  precipitation_rate_av  !< Avg. of precipitation rate
770    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  ol_av                  !< Avg. of Obukhov length
771    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  qsws_av                !< Avg. of surface moisture flux
772    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  ssws_av                !< Avg. of surface scalar flux
773    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  shf_av                 !< Avg. of surface heat flux
774    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  ts_av                  !< Avg. of characteristic temperature scale
775    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  us_av                  !< Avg. of friction velocity
776    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  z0_av                  !< Avg. of roughness length for momentum
777    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  z0h_av                 !< Avg. of roughness length for heat
778    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  z0q_av                 !< Avg. of roughness length for moisture
779
780    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  e_av          !<
781    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  lpt_av        !<
782    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  nr_av         !<
783    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  p_av          !<
784    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  pc_av         !<
785    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  pr_av         !<
786    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  prr_av        !<
787    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  pt_av         !<
788    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  q_av          !<
789    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  qc_av         !<
790    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  ql_av         !<
791    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  ql_c_av       !<
792    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  ql_v_av       !<
793    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  ql_vp_av      !<
794    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  qr_av         !<
795    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  qv_av         !<
796    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  rho_ocean_av  !<
797    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  s_av          !<
798    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  sa_av         !<
799    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  u_av          !<
800    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  v_av          !<
801    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  vpt_av        !<
802    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE, TARGET ::  w_av          !<
803 
804 END MODULE averaging
805
806
807!------------------------------------------------------------------------------!
808! Description:
809! ------------
810!> Definition of variables and constants for cloud physics.
811!------------------------------------------------------------------------------!
812 MODULE cloud_parameters
813 
814    USE kinds
815
816    REAL(wp) ::  cp = 1005.0_wp         !< heat capacity of dry air (J kg-1 K-1)
817    REAL(wp) ::  l_v = 2.5E+06_wp       !< latent heat of vaporization (J kg-1)
818    REAL(wp) ::  l_d_cp                 !< l_v / cp
819    REAL(wp) ::  l_d_r                  !< l_v / r_d
820    REAL(wp) ::  l_d_rv                 !< l_v / r_v   
821    REAL(wp) ::  rho_l = 1.0E3_wp       !< density of water (kg m-3)
822    REAL(wp) ::  r_d = 287.0_wp         !< sp. gas const. dry air (J kg-1 K-1)
823    REAL(wp) ::  r_v = 461.51_wp        !< sp. gas const. water vapor (J kg-1 K-1)
824
825
826    REAL(wp), DIMENSION(:), ALLOCATABLE ::  hyrho   !<
827    REAL(wp), DIMENSION(:), ALLOCATABLE ::  pt_d_t  !<
828    REAL(wp), DIMENSION(:), ALLOCATABLE ::  t_d_pt  !<
829
830    SAVE
831
832 END MODULE cloud_parameters
833
834
835!------------------------------------------------------------------------------!
836! Description:
837! ------------
838!> Definition of general constants.
839!------------------------------------------------------------------------------!
840 MODULE constants
841 
842    USE kinds
843
844    REAL(wp) ::  pi = 3.141592654_wp  !<
845    REAL(wp) ::  adv_mom_1            !<
846    REAL(wp) ::  adv_mom_3            !<
847    REAL(wp) ::  adv_mom_5            !<
848    REAL(wp) ::  adv_sca_1            !<
849    REAL(wp) ::  adv_sca_3            !<
850    REAL(wp) ::  adv_sca_5            !<
851   
852    SAVE
853
854 END MODULE constants
855
856 
857!------------------------------------------------------------------------------!
858! Description:
859! ------------
860!> Definition of parameters for program control
861!------------------------------------------------------------------------------!
862 MODULE control_parameters
863
864    USE kinds
865
866    TYPE plot_precision
867       CHARACTER (LEN=8) ::  variable   !<
868       INTEGER(iwp)      ::  precision  !<
869    END TYPE plot_precision
870
871    TYPE(plot_precision), DIMENSION(100) ::  plot_3d_precision =               &  !<
872                        (/ plot_precision( 'u', 2 ), plot_precision( 'v', 2 ), & 
873                           plot_precision( 'w', 2 ), plot_precision( 'p', 5 ), & 
874                           plot_precision( 'pt', 2 ),                          & 
875                           ( plot_precision( ' ', 1 ), i9 = 1,95 ) /)
876
877    TYPE file_status
878       LOGICAL ::  opened         !<
879       LOGICAL ::  opened_before  !<
880    END TYPE file_status
881   
882    INTEGER, PARAMETER      ::  mask_xyz_dimension = 100  !<
883    INTEGER, PARAMETER      ::  max_masks = 50            !<
884    INTEGER(iwp), PARAMETER ::  varnamelength = 30        !< length of output variable names
885
886    TYPE(file_status), DIMENSION(200+2*max_masks) ::                &  !<
887                             openfile = file_status(.FALSE.,.FALSE.)
888
889    CHARACTER (LEN=1)    ::  cycle_mg = 'w'                               !<
890    CHARACTER (LEN=1)    ::  timestep_reason = ' '                        !<
891    CHARACTER (LEN=3)    ::  coupling_char = ''                           !<
892    CHARACTER (LEN=5)    ::  write_binary = 'false'                       !<
893    CHARACTER (LEN=8)    ::  most_method = 'lookup'                       !< NAMELIST parameter defining method to be used to calculate Okukhov length,
894    CHARACTER (LEN=8)    ::  run_date                                     !<
895    CHARACTER (LEN=8)    ::  run_time                                     !<
896    CHARACTER (LEN=9)    ::  simulated_time_chr                           !<
897    CHARACTER (LEN=11)   ::  topography_grid_convention = ' '             !<
898    CHARACTER (LEN=12)   ::  version = ' '                                !<
899    CHARACTER (LEN=12)   ::  revision = ' '                               !<
900    CHARACTER (LEN=12)   ::  user_interface_current_revision = ' '        !<
901    CHARACTER (LEN=12)   ::  user_interface_required_revision = ' '       !<
902    CHARACTER (LEN=16)   ::  conserve_volume_flow_mode = 'default'        !<
903    CHARACTER (LEN=16)   ::  loop_optimization = 'default'                !<
904    CHARACTER (LEN=16)   ::  momentum_advec = 'ws-scheme'                 !<
905    CHARACTER (LEN=16)   ::  psolver = 'poisfft'                          !< 
906    CHARACTER (LEN=16)   ::  scalar_advec = 'ws-scheme'                   !<
907    CHARACTER (LEN=20)   ::  approximation = 'boussinesq'                 !< 
908    CHARACTER (LEN=40)   ::  flux_input_mode = 'approximation-specific'   !<
909    CHARACTER (LEN=40)   ::  flux_output_mode = 'approximation-specific'  !<
910    CHARACTER (LEN=20)   ::  bc_e_b = 'neumann'                           !<
911    CHARACTER (LEN=20)   ::  bc_lr = 'cyclic'                             !<
912    CHARACTER (LEN=20)   ::  bc_ns = 'cyclic'                             !<
913    CHARACTER (LEN=20)   ::  bc_p_b = 'neumann'                           !<
914    CHARACTER (LEN=20)   ::  bc_p_t = 'dirichlet'                         !<
915    CHARACTER (LEN=20)   ::  bc_pt_b = 'dirichlet'                        !<
916    CHARACTER (LEN=20)   ::  bc_pt_t = 'initial_gradient'                 !<
917    CHARACTER (LEN=20)   ::  bc_q_b = 'dirichlet'                         !<
918    CHARACTER (LEN=20)   ::  bc_q_t = 'neumann'                           !<
919    CHARACTER (LEN=20)   ::  bc_s_b = 'dirichlet'                         !<
920    CHARACTER (LEN=20)   ::  bc_s_t = 'initial_gradient'                  !<
921    CHARACTER (LEN=20)   ::  bc_sa_t = 'neumann'                          !<
922    CHARACTER (LEN=20)   ::  bc_uv_b = 'dirichlet'                        !<
923    CHARACTER (LEN=20)   ::  bc_uv_t = 'dirichlet'                        !<
924    CHARACTER (LEN=20)   ::  cloud_scheme = 'saturation_adjust'           !<
925    CHARACTER (LEN=20)   ::  coupling_mode = 'uncoupled'                  !<   
926    CHARACTER (LEN=20)   ::  coupling_mode_remote = 'uncoupled'           !<
927    CHARACTER (LEN=20)   ::  dissipation_1d = 'detering'                  !<
928    CHARACTER (LEN=20)   ::  fft_method = 'system-specific'               !<
929    CHARACTER (LEN=20)   ::  mixing_length_1d = 'blackadar'               !<
930    CHARACTER (LEN=20)   ::  random_generator = 'numerical-recipes'       !<
931    CHARACTER (LEN=20)   ::  reference_state = 'initial_profile'          !<   
932    CHARACTER (LEN=20)   ::  return_addres                                !<
933    CHARACTER (LEN=20)   ::  return_username                              !<
934    CHARACTER (LEN=20)   ::  timestep_scheme = 'runge-kutta-3'            !<       
935    CHARACTER (LEN=40)   ::  avs_data_file                                !<
936    CHARACTER (LEN=40)   ::  topography = 'flat'                          !< 
937    CHARACTER (LEN=64)   ::  host = ' '                                   !<
938    CHARACTER (LEN=80)   ::  log_message                                  !<
939    CHARACTER (LEN=80)   ::  run_identifier                               !< 
940    CHARACTER (LEN=100)  ::  initializing_actions = ' '                   !<
941    CHARACTER (LEN=110)  ::  run_description_header                       !<
942    CHARACTER (LEN=1000) ::  message_string = ' '                         !<
943
944    CHARACTER (LEN=varnamelength), DIMENSION(500) ::  data_output = ' '       !<
945    CHARACTER (LEN=varnamelength), DIMENSION(500) ::  data_output_user = ' '  !<
946    CHARACTER (LEN=varnamelength), DIMENSION(500) ::  doav = ' '              !<
947                                           
948    CHARACTER (LEN=varnamelength), DIMENSION(max_masks,100) ::  data_output_masks = ' '       !<
949    CHARACTER (LEN=varnamelength), DIMENSION(max_masks,100) ::  data_output_masks_user = ' '  !<
950
951    CHARACTER (LEN=varnamelength), DIMENSION(300) ::  data_output_pr = ' '  !<
952   
953    CHARACTER (LEN=varnamelength), DIMENSION(200) ::  data_output_pr_user = ' '  !<
954   
955    CHARACTER (LEN=varnamelength), DIMENSION(max_masks,0:1,100) ::  domask = ' ' !<
956   
957    CHARACTER (LEN=varnamelength), DIMENSION(0:1,500) ::  do2d = ' '  !<
958    CHARACTER (LEN=varnamelength), DIMENSION(0:1,500) ::  do3d = ' '  !<
959
960    INTEGER(iwp), PARAMETER ::  fl_max = 100     !<
961    INTEGER(iwp), PARAMETER ::  var_fl_max = 20  !<
962   
963    INTEGER(iwp) ::  abort_mode = 1                    !<
964    INTEGER(iwp) ::  average_count_pr = 0              !<
965    INTEGER(iwp) ::  average_count_3d = 0              !<
966    INTEGER(iwp) ::  current_timestep_number = 0       !<
967    INTEGER(iwp) ::  coupling_topology = 0             !<
968    INTEGER(iwp) ::  dist_range = 0                    !<
969    INTEGER(iwp) ::  disturbance_level_ind_b           !<
970    INTEGER(iwp) ::  disturbance_level_ind_t           !<
971    INTEGER(iwp) ::  doav_n = 0                        !<
972    INTEGER(iwp) ::  dopr_n = 0                        !<
973    INTEGER(iwp) ::  dopr_time_count = 0               !<
974    INTEGER(iwp) ::  dopts_time_count = 0              !<
975    INTEGER(iwp) ::  dots_time_count = 0               !<
976    INTEGER(iwp) ::  do2d_xy_n = 0                     !<
977    INTEGER(iwp) ::  do2d_xz_n = 0                     !<
978    INTEGER(iwp) ::  do2d_yz_n = 0                     !<
979    INTEGER(iwp) ::  do3d_avs_n = 0                    !<
980    INTEGER(iwp) ::  dp_level_ind_b = 0                !<
981    INTEGER(iwp) ::  dvrp_filecount = 0                !<
982    INTEGER(iwp) ::  dz_stretch_level_index            !<
983    INTEGER(iwp) ::  ensemble_member_nr = 0            !<
984    INTEGER(iwp) ::  gamma_mg                          !<
985    INTEGER(iwp) ::  gathered_size                     !<
986    INTEGER(iwp) ::  grid_level                        !<
987    INTEGER(iwp) ::  ibc_e_b                           !<
988    INTEGER(iwp) ::  ibc_p_b                           !<
989    INTEGER(iwp) ::  ibc_p_t                           !<
990    INTEGER(iwp) ::  ibc_pt_b                          !<
991    INTEGER(iwp) ::  ibc_pt_t                          !<
992    INTEGER(iwp) ::  ibc_q_b                           !<
993    INTEGER(iwp) ::  ibc_q_t                           !<
994    INTEGER(iwp) ::  ibc_s_b                           !<
995    INTEGER(iwp) ::  ibc_s_t                           !<
996    INTEGER(iwp) ::  ibc_sa_t                          !<
997    INTEGER(iwp) ::  ibc_uv_b                          !<
998    INTEGER(iwp) ::  ibc_uv_t                          !<
999    INTEGER(iwp) ::  inflow_disturbance_begin = -1     !<
1000    INTEGER(iwp) ::  inflow_disturbance_end = -1       !<
1001    INTEGER(iwp) ::  intermediate_timestep_count       !<
1002    INTEGER(iwp) ::  intermediate_timestep_count_max   !<
1003    INTEGER(iwp) ::  io_group = 0                      !<
1004    INTEGER(iwp) ::  io_blocks = 1                     !<
1005    INTEGER(iwp) ::  iran = -1234567                   !<
1006    INTEGER(iwp) ::  masks = 0                         !<
1007    INTEGER(iwp) ::  maximum_grid_level                !<
1008    INTEGER(iwp) ::  maximum_parallel_io_streams = -1  !<
1009    INTEGER(iwp) ::  max_pr_user = 0                   !<
1010    INTEGER(iwp) ::  mgcycles = 0                      !<
1011    INTEGER(iwp) ::  mg_cycles = -1                    !<
1012    INTEGER(iwp) ::  mg_switch_to_pe0_level = -1       !<
1013    INTEGER(iwp) ::  mid                               !<
1014    INTEGER(iwp) ::  nlsf = 1000                       !<
1015    INTEGER(iwp) ::  ntnudge = 1000                    !<
1016    INTEGER(iwp) ::  ngsrb = 2                         !<
1017    INTEGER(iwp) ::  nr_timesteps_this_run = 0         !<
1018    INTEGER(iwp) ::  nsor = 20                         !<
1019    INTEGER(iwp) ::  nsor_ini = 100                    !<
1020    INTEGER(iwp) ::  n_sor                             !<
1021    INTEGER(iwp) ::  normalizing_region = 0            !<
1022    INTEGER(iwp) ::  num_leg=0                         !<
1023    INTEGER(iwp) ::  num_var_fl                        !<
1024    INTEGER(iwp) ::  num_var_fl_user=0                 !<
1025    INTEGER(iwp) ::  nz_do3d = -9999                   !<
1026    INTEGER(iwp) ::  prt_time_count = 0                !<
1027    INTEGER(iwp) ::  recycling_plane                   !<
1028    INTEGER(iwp) ::  runnr = 0                         !<
1029    INTEGER(iwp) ::  skip_do_avs = 0                   !<
1030    INTEGER(iwp) ::  subdomain_size                    !<
1031    INTEGER(iwp) ::  terminate_coupled = 0             !<
1032    INTEGER(iwp) ::  terminate_coupled_remote = 0      !<
1033    INTEGER(iwp) ::  timestep_count = 0                !<
1034
1035    INTEGER(iwp) ::  dist_nxl(0:1)                               !<
1036    INTEGER(iwp) ::  dist_nxr(0:1)                               !<
1037    INTEGER(iwp) ::  dist_nyn(0:1)                               !<
1038    INTEGER(iwp) ::  dist_nys(0:1)                               !<
1039    INTEGER(iwp) ::  do2d_no(0:1) = 0                            !<
1040    INTEGER(iwp) ::  do2d_xy_time_count(0:1)                     !<
1041    INTEGER(iwp) ::  do2d_xz_time_count(0:1)                     !<
1042    INTEGER(iwp) ::  do2d_yz_time_count(0:1)                     !<
1043    INTEGER(iwp) ::  do3d_no(0:1) = 0                            !<
1044    INTEGER(iwp) ::  do3d_time_count(0:1)                        !<
1045    INTEGER(iwp) ::  domask_no(max_masks,0:1) = 0                !<
1046    INTEGER(iwp) ::  domask_time_count(max_masks,0:1)            !<
1047    INTEGER(iwp) ::  mask_size(max_masks,3) = -1                 !<
1048    INTEGER(iwp) ::  mask_size_l(max_masks,3) = -1               !<
1049    INTEGER(iwp) ::  mask_start_l(max_masks,3) = -1              !<
1050    INTEGER(iwp) ::  pt_vertical_gradient_level_ind(10) = -9999  !<
1051    INTEGER(iwp) ::  q_vertical_gradient_level_ind(10) = -9999   !<
1052    INTEGER(iwp) ::  s_vertical_gradient_level_ind(10) = -9999   !<               
1053    INTEGER(iwp) ::  sa_vertical_gradient_level_ind(10) = -9999  !<
1054    INTEGER(iwp) ::  section(100,3)                              !<
1055    INTEGER(iwp) ::  section_xy(100) = -9999                     !<
1056    INTEGER(iwp) ::  section_xz(100) = -9999                     !<
1057    INTEGER(iwp) ::  section_yz(100) = -9999                     !<
1058    INTEGER(iwp) ::  ug_vertical_gradient_level_ind(10) = -9999  !<
1059    INTEGER(iwp) ::  vg_vertical_gradient_level_ind(10) = -9999  !<
1060    INTEGER(iwp) ::  subs_vertical_gradient_level_i(10) = -9999  !<
1061
1062
1063    INTEGER(iwp), DIMENSION(0:1) ::  ntdim_2d_xy  !<
1064    INTEGER(iwp), DIMENSION(0:1) ::  ntdim_2d_xz  !<
1065    INTEGER(iwp), DIMENSION(0:1) ::  ntdim_2d_yz  !<
1066    INTEGER(iwp), DIMENSION(0:1) ::  ntdim_3d     !<
1067
1068    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  grid_level_count  !<
1069
1070    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  mask_i         !<
1071    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  mask_j         !<
1072    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  mask_k         !<
1073    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  mask_i_global  !<   
1074    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  mask_j_global  !<
1075    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  mask_k_global  !<
1076
1077    LOGICAL ::  bc_lr_cyc =.TRUE.                            !<
1078    LOGICAL ::  bc_lr_dirrad = .FALSE.                       !<
1079    LOGICAL ::  bc_lr_raddir = .FALSE.                       !<
1080    LOGICAL ::  bc_ns_cyc = .TRUE.                           !<
1081    LOGICAL ::  bc_ns_dirrad = .FALSE.                       !<
1082    LOGICAL ::  bc_ns_raddir = .FALSE.                       !<
1083    LOGICAL ::  call_microphysics_at_all_substeps = .FALSE.  !<
1084    LOGICAL ::  call_psolver_at_all_substeps = .TRUE.        !<
1085    LOGICAL ::  cloud_droplets = .FALSE.                     !<
1086    LOGICAL ::  cloud_physics = .FALSE.                      !<
1087    LOGICAL ::  cloud_top_radiation = .FALSE.                !<
1088    LOGICAL ::  conserve_volume_flow = .FALSE.               !<
1089    LOGICAL ::  constant_diffusion = .FALSE.                 !<
1090    LOGICAL ::  constant_flux_layer = .TRUE.                 !<
1091    LOGICAL ::  constant_heatflux = .TRUE.                   !< 
1092    LOGICAL ::  constant_top_heatflux = .TRUE.               !<
1093    LOGICAL ::  constant_top_momentumflux = .FALSE.          !<
1094    LOGICAL ::  constant_top_salinityflux = .TRUE.           !<
1095    LOGICAL ::  constant_top_scalarflux = .TRUE.             !<
1096    LOGICAL ::  constant_scalarflux = .TRUE.                 !<
1097    LOGICAL ::  constant_waterflux = .TRUE.                  !<
1098    LOGICAL ::  create_disturbances = .TRUE.                 !<
1099    LOGICAL ::  data_output_2d_on_each_pe = .TRUE.           !<
1100    LOGICAL ::  dissipation_control = .FALSE.                !<
1101    LOGICAL ::  disturbance_created = .FALSE.                !<
1102    LOGICAL ::  do2d_at_begin = .FALSE.                      !<
1103    LOGICAL ::  do3d_at_begin = .FALSE.                      !<
1104    LOGICAL ::  do_sum = .FALSE.                             !<
1105    LOGICAL ::  dp_external = .FALSE.                        !<
1106    LOGICAL ::  dp_smooth = .FALSE.                          !<
1107    LOGICAL ::  dt_fixed = .FALSE.                           !<
1108    LOGICAL ::  dt_3d_reached                                !<
1109    LOGICAL ::  dt_3d_reached_l                              !<
1110    LOGICAL ::  exchange_mg = .FALSE.                        !<
1111    LOGICAL ::  first_call_lpm = .TRUE.                      !<
1112    LOGICAL ::  force_print_header = .FALSE.                 !<
1113    LOGICAL ::  galilei_transformation = .FALSE.             !<
1114    LOGICAL ::  humidity = .FALSE.                           !<
1115    LOGICAL ::  humidity_remote = .FALSE.                    !<
1116    LOGICAL ::  inflow_l = .FALSE.                           !<
1117    LOGICAL ::  inflow_n = .FALSE.                           !<
1118    LOGICAL ::  inflow_r = .FALSE.                           !<
1119    LOGICAL ::  inflow_s = .FALSE.                           !<
1120    LOGICAL ::  large_scale_forcing = .FALSE.                !<
1121    LOGICAL ::  large_scale_subsidence = .FALSE.             !<
1122    LOGICAL ::  lsf_exception = .FALSE.                      !< temporary flag for use of lsf with buildings on flat terrain
1123    LOGICAL ::  lsf_surf = .TRUE.                            !<
1124    LOGICAL ::  lsf_vert = .TRUE.                            !<
1125    LOGICAL ::  lptnudge = .FALSE.                           !<
1126    LOGICAL ::  lqnudge = .FALSE.                            !<
1127    LOGICAL ::  lunudge = .FALSE.                            !<
1128    LOGICAL ::  lvnudge = .FALSE.                            !<
1129    LOGICAL ::  lwnudge = .FALSE.                            !<
1130    LOGICAL ::  masking_method = .FALSE.                     !<
1131    LOGICAL ::  microphysics_sat_adjust = .FALSE.            !<
1132    LOGICAL ::  microphysics_kessler = .FALSE.               !<
1133    LOGICAL ::  microphysics_seifert = .FALSE.               !<
1134    LOGICAL ::  mg_switch_to_pe0 = .FALSE.                   !<
1135    LOGICAL ::  monotonic_adjustment = .FALSE.               !<
1136    LOGICAL ::  nest_bound_l = .FALSE.                       !< nested boundary on left side
1137    LOGICAL ::  nest_bound_n = .FALSE.                       !< nested boundary on north side
1138    LOGICAL ::  nest_bound_r = .FALSE.                       !< nested boundary on right side
1139    LOGICAL ::  nest_bound_s = .FALSE.                       !< nested boundary on south side
1140    LOGICAL ::  nest_domain  = .FALSE.                       !< domain is nested into a parent domain
1141    LOGICAL ::  neutral = .FALSE.                            !<
1142    LOGICAL ::  nudging = .FALSE.                            !<
1143    LOGICAL ::  ocean = .FALSE.                              !<
1144    LOGICAL ::  on_device = .FALSE.                          !<
1145    LOGICAL ::  outflow_l = .FALSE.                          !<
1146    LOGICAL ::  outflow_n = .FALSE.                          !<
1147    LOGICAL ::  outflow_r = .FALSE.                          !<
1148    LOGICAL ::  outflow_s = .FALSE.                          !<
1149    LOGICAL ::  passive_scalar = .FALSE.                     !<
1150    LOGICAL ::  precipitation = .FALSE.                      !<
1151    LOGICAL ::  random_heatflux = .FALSE.                    !<
1152    LOGICAL ::  recycling_yshift = .FALSE.                   !<
1153    LOGICAL ::  run_control_header = .FALSE.                 !<
1154    LOGICAL ::  run_coupled = .TRUE.                         !<
1155    LOGICAL ::  scalar_rayleigh_damping = .TRUE.             !<
1156    LOGICAL ::  sloping_surface = .FALSE.                    !<
1157    LOGICAL ::  stop_dt = .FALSE.                            !<
1158    LOGICAL ::  synchronous_exchange = .FALSE.               !<
1159    LOGICAL ::  terminate_run = .FALSE.                      !<
1160    LOGICAL ::  transpose_compute_overlap = .FALSE.          !<
1161    LOGICAL ::  turbulent_inflow = .FALSE.                   !<
1162    LOGICAL ::  turbulent_outflow = .FALSE.                  !< flag for turbulent outflow condition
1163    LOGICAL ::  urban_surface = .FALSE.                      !< flag for urban surface model
1164    LOGICAL ::  use_cmax = .TRUE.                            !<
1165    LOGICAL ::  use_initial_profile_as_reference = .FALSE.   !<
1166    LOGICAL ::  use_prescribed_profile_data = .FALSE.        !<
1167    LOGICAL ::  use_single_reference_value = .FALSE.         !<
1168    LOGICAL ::  use_subsidence_tendencies = .FALSE.          !<
1169    LOGICAL ::  use_surface_fluxes = .FALSE.                 !<
1170    LOGICAL ::  use_top_fluxes = .FALSE.                     !<
1171    LOGICAL ::  use_ug_for_galilei_tr = .TRUE.               !<
1172    LOGICAL ::  use_upstream_for_tke = .FALSE.               !<
1173    LOGICAL ::  virtual_flight = .FALSE.                     !< flag for virtual flight model
1174    LOGICAL ::  wall_adjustment = .TRUE.                     !<
1175    LOGICAL ::  ws_scheme_sca = .FALSE.                      !<
1176    LOGICAL ::  ws_scheme_mom = .FALSE.                      !<
1177
1178    LOGICAL ::  data_output_xy(0:1) = .FALSE.                !<
1179    LOGICAL ::  data_output_xz(0:1) = .FALSE.                !<
1180    LOGICAL ::  data_output_yz(0:1) = .FALSE.                !<
1181
1182    REAL(wp) ::  advected_distance_x = 0.0_wp                  !<
1183    REAL(wp) ::  advected_distance_y = 0.0_wp                  !<
1184    REAL(wp) ::  alpha_surface = 0.0_wp                        !<
1185    REAL(wp) ::  atmos_ocean_sign = 1.0_wp                     !<
1186    REAL(wp) ::  averaging_interval = 0.0_wp                   !<
1187    REAL(wp) ::  averaging_interval_pr = 9999999.9_wp          !<
1188    REAL(wp) ::  bc_pt_t_val                                   !<
1189    REAL(wp) ::  bc_q_t_val                                    !<
1190    REAL(wp) ::  bc_s_t_val                                    !<
1191    REAL(wp) ::  bottom_salinityflux = 0.0_wp                  !<
1192    REAL(wp) ::  building_height = 50.0_wp                     !<
1193    REAL(wp) ::  building_length_x = 50.0_wp                   !<
1194    REAL(wp) ::  building_length_y = 50.0_wp                   !<
1195    REAL(wp) ::  building_wall_left = 9999999.9_wp             !<
1196    REAL(wp) ::  building_wall_south = 9999999.9_wp            !<
1197    REAL(wp) ::  canyon_height = 50.0_wp                       !<
1198    REAL(wp) ::  canyon_width_x = 9999999.9_wp                 !<
1199    REAL(wp) ::  canyon_width_y = 9999999.9_wp                 !<
1200    REAL(wp) ::  canyon_wall_left = 9999999.9_wp               !<
1201    REAL(wp) ::  canyon_wall_south = 9999999.9_wp              !<
1202    REAL(wp) ::  cfl_factor = -1.0_wp                          !<
1203    REAL(wp) ::  cos_alpha_surface                             !<
1204    REAL(wp) ::  coupling_start_time = 0.0_wp                  !<
1205    REAL(wp) ::  disturbance_amplitude = 0.25_wp               !<
1206    REAL(wp) ::  disturbance_energy_limit = 0.01_wp            !<
1207    REAL(wp) ::  disturbance_level_b = -9999999.9_wp           !<
1208    REAL(wp) ::  disturbance_level_t = -9999999.9_wp           !<
1209    REAL(wp) ::  dp_level_b = 0.0_wp                           !<
1210    REAL(wp) ::  dt = -1.0_wp                                  !<
1211    REAL(wp) ::  dt_averaging_input = 0.0_wp                   !<
1212    REAL(wp) ::  dt_averaging_input_pr = 9999999.9_wp          !<
1213    REAL(wp) ::  dt_coupling = 9999999.9_wp                    !<
1214    REAL(wp) ::  dt_data_output = 9999999.9_wp                 !<
1215    REAL(wp) ::  dt_data_output_av = 9999999.9_wp              !<
1216    REAL(wp) ::  dt_disturb = 9999999.9_wp                     !<
1217    REAL(wp) ::  dt_dopr = 9999999.9_wp                        !<
1218    REAL(wp) ::  dt_dopr_listing = 9999999.9_wp                !<
1219    REAL(wp) ::  dt_dopts = 9999999.9_wp                       !<
1220    REAL(wp) ::  dt_dots = 9999999.9_wp                        !<
1221    REAL(wp) ::  dt_do2d_xy = 9999999.9_wp                     !<
1222    REAL(wp) ::  dt_do2d_xz = 9999999.9_wp                     !<
1223    REAL(wp) ::  dt_do2d_yz = 9999999.9_wp                     !<
1224    REAL(wp) ::  dt_do3d = 9999999.9_wp                        !<
1225    REAL(wp) ::  dt_dvrp = 9999999.9_wp                        !<
1226    REAL(wp) ::  dt_max = 20.0_wp                              !<
1227    REAL(wp) ::  dt_restart = 9999999.9_wp                     !<
1228    REAL(wp) ::  dt_run_control = 60.0_wp                      !<
1229    REAL(wp) ::  dt_3d = 1.0_wp                                !<
1230    REAL(wp) ::  dz = -1.0_wp                                  !<
1231    REAL(wp) ::  dz_max = 9999999.9_wp                         !<
1232    REAL(wp) ::  dz_stretch_factor = 1.08_wp                   !<
1233    REAL(wp) ::  dz_stretch_level = 100000.0_wp                !<
1234    REAL(wp) ::  e_init = 0.0_wp                               !<
1235    REAL(wp) ::  e_min = 0.0_wp                                !<
1236    REAL(wp) ::  end_time = 0.0_wp                             !<
1237    REAL(wp) ::  f = 0.0_wp                                    !<
1238    REAL(wp) ::  fs = 0.0_wp                                   !<
1239    REAL(wp) ::  g = 9.81_wp                                   !<
1240    REAL(wp) ::  inflow_damping_height = 9999999.9_wp          !<
1241    REAL(wp) ::  inflow_damping_width = 9999999.9_wp           !<
1242    REAL(wp) ::  kappa = 0.4_wp                                !<
1243    REAL(wp) ::  km_constant = -1.0_wp                         !<
1244    REAL(wp) ::  mask_scale_x = 1.0_wp                         !<
1245    REAL(wp) ::  mask_scale_y = 1.0_wp                         !<
1246    REAL(wp) ::  mask_scale_z = 1.0_wp                         !<
1247    REAL(wp) ::  maximum_cpu_time_allowed = 0.0_wp             !<
1248    REAL(wp) ::  molecular_viscosity = 1.461E-5_wp             !<
1249    REAL(wp) ::  old_dt = 1.0E-10_wp                           !<
1250    REAL(wp) ::  omega = 7.29212E-5_wp                         !<
1251    REAL(wp) ::  omega_sor = 1.8_wp                            !<
1252    REAL(wp) ::  outflow_source_plane = -9999999.9_wp          !< x-position of outflow-source plane (turbulent outflow method)
1253    REAL(wp) ::  particle_maximum_age = 9999999.9_wp           !<
1254    REAL(wp) ::  phi = 55.0_wp                                 !<
1255    REAL(wp) ::  prandtl_number = 1.0_wp                       !<
1256    REAL(wp) ::  precipitation_amount_interval = 9999999.9_wp  !<
1257    REAL(wp) ::  prho_reference                                !<
1258    REAL(wp) ::  pt_damping_factor = 0.0_wp                    !<
1259    REAL(wp) ::  pt_damping_width = 0.0_wp                     !<
1260    REAL(wp) ::  pt_reference = 9999999.9_wp                   !<
1261    REAL(wp) ::  pt_slope_offset = 0.0_wp                      !<
1262    REAL(wp) ::  pt_surface = 300.0_wp                         !<
1263    REAL(wp) ::  pt_surface_initial_change = 0.0_wp            !<
1264    REAL(wp) ::  q_surface = 0.0_wp                            !<
1265    REAL(wp) ::  q_surface_initial_change = 0.0_wp             !<
1266    REAL(wp) ::  rayleigh_damping_factor = -1.0_wp             !<
1267    REAL(wp) ::  rayleigh_damping_height = -1.0_wp             !<
1268    REAL(wp) ::  recycling_width = 9999999.9_wp                !<
1269    REAL(wp) ::  residual_limit = 1.0E-4_wp                    !<
1270    REAL(wp) ::  restart_time = 9999999.9_wp                   !<
1271    REAL(wp) ::  rho_reference                                 !<
1272    REAL(wp) ::  rho_surface                                   !<
1273    REAL(wp) ::  roughness_length = 0.1_wp                     !<
1274    REAL(wp) ::  sa_surface = 35.0_wp                          !<
1275    REAL(wp) ::  simulated_time = 0.0_wp                       !<
1276    REAL(wp) ::  simulated_time_at_begin                       !<
1277    REAL(wp) ::  sin_alpha_surface                             !<
1278    REAL(wp) ::  skip_time_data_output = 0.0_wp                !<
1279    REAL(wp) ::  skip_time_data_output_av = 9999999.9_wp       !<
1280    REAL(wp) ::  skip_time_dopr = 9999999.9_wp                 !<
1281    REAL(wp) ::  skip_time_do2d_xy = 9999999.9_wp              !<
1282    REAL(wp) ::  skip_time_do2d_xz = 9999999.9_wp              !<
1283    REAL(wp) ::  skip_time_do2d_yz = 9999999.9_wp              !<
1284    REAL(wp) ::  skip_time_do3d = 9999999.9_wp                 !<
1285    REAL(wp) ::  surface_heatflux = 9999999.9_wp               !<
1286    REAL(wp) ::  surface_pressure = 1013.25_wp                 !<
1287    REAL(wp) ::  surface_scalarflux = 9999999.9_wp             !<
1288    REAL(wp) ::  surface_waterflux = 9999999.9_wp              !<
1289    REAL(wp) ::  s_surface = 0.0_wp                            !<
1290    REAL(wp) ::  s_surface_initial_change = 0.0_wp             !<
1291    REAL(wp) ::  termination_time_needed = -1.0_wp             !<
1292    REAL(wp) ::  time_coupling = 0.0_wp                        !<
1293    REAL(wp) ::  time_disturb = 0.0_wp                         !<
1294    REAL(wp) ::  time_dopr = 0.0_wp                            !<
1295    REAL(wp) ::  time_dopr_av = 0.0_wp                         !<
1296    REAL(wp) ::  time_dopr_listing = 0.0_wp                    !<
1297    REAL(wp) ::  time_dopts = 0.0_wp                           !<
1298    REAL(wp) ::  time_dosp = 0.0_wp                            !<
1299    REAL(wp) ::  time_dosp_av = 0.0_wp                         !<
1300    REAL(wp) ::  time_dots = 0.0_wp                            !<
1301    REAL(wp) ::  time_do2d_xy = 0.0_wp                         !<
1302    REAL(wp) ::  time_do2d_xz = 0.0_wp                         !<
1303    REAL(wp) ::  time_do2d_yz = 0.0_wp                         !<
1304    REAL(wp) ::  time_do3d = 0.0_wp                            !<
1305    REAL(wp) ::  time_do_av = 0.0_wp                           !<
1306    REAL(wp) ::  time_do_sla = 0.0_wp                          !<
1307    REAL(wp) ::  time_dvrp = 0.0_wp                            !<
1308    REAL(wp) ::  time_restart = 9999999.9_wp                   !<
1309    REAL(wp) ::  time_run_control = 0.0_wp                     !<
1310    REAL(wp) ::  time_since_reference_point                    !<
1311    REAL(wp) ::  top_heatflux = 9999999.9_wp                   !<
1312    REAL(wp) ::  top_momentumflux_u = 9999999.9_wp             !<
1313    REAL(wp) ::  top_momentumflux_v = 9999999.9_wp             !<
1314    REAL(wp) ::  top_salinityflux = 9999999.9_wp               !<
1315    REAL(wp) ::  top_scalarflux = 9999999.9_wp                 !<
1316    REAL(wp) ::  ug_surface = 0.0_wp                           !<
1317    REAL(wp) ::  u_bulk = 0.0_wp                               !<
1318    REAL(wp) ::  u_gtrans = 0.0_wp                             !<
1319    REAL(wp) ::  vg_surface = 0.0_wp                           !<
1320    REAL(wp) ::  vpt_reference = 9999999.9_wp                  !<
1321    REAL(wp) ::  v_bulk = 0.0_wp                               !<
1322    REAL(wp) ::  v_gtrans = 0.0_wp                             !<
1323    REAL(wp) ::  wall_adjustment_factor = 1.8_wp               !<
1324    REAL(wp) ::  zeta_max = 20.0_wp                            !< Maximum value of zeta = z/L
1325    REAL(wp) ::  zeta_min = -9990.0_wp                         !< Minimum value of zeta = z/L
1326    REAL(wp) ::  z_max_do2d = -1.0_wp                          !<
1327    REAL(wp) ::  z0h_factor = 1.0_wp                           !<
1328
1329    REAL(wp) ::  do2d_xy_last_time(0:1) = -1.0_wp                  !<
1330    REAL(wp) ::  do2d_xz_last_time(0:1) = -1.0_wp                  !<
1331    REAL(wp) ::  do2d_yz_last_time(0:1) = -1.0_wp                  !<
1332    REAL(wp) ::  dpdxy(1:2) = 0.0_wp                               !<
1333    REAL(wp) ::  dt_domask(max_masks) = 9999999.9_wp               !<
1334    REAL(wp) ::  mask_scale(3)                                     !<
1335    REAL(wp) ::  pt_vertical_gradient(10) = 0.0_wp                 !<
1336    REAL(wp) ::  pt_vertical_gradient_level(10) = -1.0_wp          !<
1337    REAL(wp) ::  q_vertical_gradient(10) = 0.0_wp                  !<
1338    REAL(wp) ::  q_vertical_gradient_level(10) = -1.0_wp           !<
1339    REAL(wp) ::  s_vertical_gradient(10) = 0.0_wp                  !<
1340    REAL(wp) ::  s_vertical_gradient_level(10) = -1.0_wp           !<
1341    REAL(wp) ::  sa_vertical_gradient(10) = 0.0_wp                 !<
1342    REAL(wp) ::  sa_vertical_gradient_level(10) = -1.0_wp          !<
1343    REAL(wp) ::  skip_time_domask(max_masks) = 9999999.9_wp        !<
1344    REAL(wp) ::  threshold(20) = 0.0_wp                            !<
1345    REAL(wp) ::  time_domask(max_masks) = 0.0_wp                   !<
1346    REAL(wp) ::  tsc(10) = (/ 1.0_wp, 1.0_wp, 0.0_wp, 0.0_wp, &    !<
1347                 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp /)
1348    REAL(wp) ::  u_profile(100) = 9999999.9_wp                     !<
1349    REAL(wp) ::  uv_heights(100) = 9999999.9_wp                    !<
1350    REAL(wp) ::  v_profile(100) = 9999999.9_wp                     !<
1351    REAL(wp) ::  ug_vertical_gradient(10) = 0.0_wp                 !<
1352    REAL(wp) ::  ug_vertical_gradient_level(10) = -9999999.9_wp    !<
1353    REAL(wp) ::  vg_vertical_gradient(10) = 0.0_wp                 !<
1354    REAL(wp) ::  vg_vertical_gradient_level(10) = -9999999.9_wp    !<
1355    REAL(wp) ::  volume_flow(1:3) = 0.0_wp                         !<
1356    REAL(wp) ::  volume_flow_area(1:3) = 0.0_wp                    !<
1357    REAL(wp) ::  volume_flow_initial(1:3) = 0.0_wp                 !<
1358    REAL(wp) ::  wall_heatflux(0:4) = 0.0_wp                       !<
1359    REAL(wp) ::  wall_humidityflux(0:4) = 0.0_wp                   !<
1360    REAL(wp) ::  wall_nrflux(0:4) = 0.0_wp                         !<
1361    REAL(wp) ::  wall_qflux(0:4) = 0.0_wp                          !<
1362    REAL(wp) ::  wall_qrflux(0:4) = 0.0_wp                         !<
1363    REAL(wp) ::  wall_salinityflux(0:4) = 0.0_wp                   !<
1364    REAL(wp) ::  wall_sflux(0:4) = 0.0_wp                          !<
1365    REAL(wp) ::  wall_scalarflux(0:4) = 0.0_wp                     !<
1366    REAL(wp) ::  subs_vertical_gradient(10) = 0.0_wp               !<
1367    REAL(wp) ::  subs_vertical_gradient_level(10) = -9999999.9_wp  !<
1368
1369    REAL(wp), DIMENSION(:), ALLOCATABLE ::  dp_smooth_factor  !<
1370
1371    REAL(wp), DIMENSION(max_masks,mask_xyz_dimension) ::  mask_x = -1.0_wp  !<
1372    REAL(wp), DIMENSION(max_masks,mask_xyz_dimension) ::  mask_y = -1.0_wp  !<
1373    REAL(wp), DIMENSION(max_masks,mask_xyz_dimension) ::  mask_z = -1.0_wp  !<
1374   
1375    REAL(wp), DIMENSION(max_masks,3) ::  mask_x_loop = -1.0_wp  !<
1376    REAL(wp), DIMENSION(max_masks,3) ::  mask_y_loop = -1.0_wp  !<
1377    REAL(wp), DIMENSION(max_masks,3) ::  mask_z_loop = -1.0_wp  !<
1378   
1379!
1380!--    internal mask arrays ("mask,dimension,selection")
1381       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  mask       !<
1382       REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  mask_loop  !<
1383
1384    SAVE
1385
1386 END MODULE control_parameters
1387
1388
1389!------------------------------------------------------------------------------!
1390! Description:
1391! ------------
1392!> Definition of variables used with dvrp-software.
1393!------------------------------------------------------------------------------!
1394 MODULE dvrp_variables
1395
1396    USE kinds
1397
1398    CHARACTER (LEN=10) ::  dvrp_output = 'rtsp'        !<
1399    CHARACTER (LEN=10) ::  particle_color = 'none'     !<
1400    CHARACTER (LEN=10) ::  particle_dvrpsize = 'none'  !<
1401
1402    CHARACTER (LEN=20), DIMENSION(10) ::  mode_dvrp = &  !<
1403                                     (/ ( '                    ', i9 = 1,10 ) /)
1404
1405    CHARACTER (LEN=80) ::  dvrp_directory = 'default'                     !<
1406    CHARACTER (LEN=80) ::  dvrp_file      = 'default'                     !<
1407    CHARACTER (LEN=80) ::  dvrp_host      = 'origin.rvs.uni-hannover.de'  !<
1408    CHARACTER (LEN=80) ::  dvrp_password  = '********'                    !<
1409    CHARACTER (LEN=80) ::  dvrp_username  = ' '                           !<
1410
1411    INTEGER(iwp) ::  cluster_size = 1                   !<
1412    INTEGER(iwp) ::  dvrp_colortable_entries = 4        !<
1413    INTEGER(iwp) ::  dvrp_colortable_entries_prt = 22   !<
1414    INTEGER(iwp) ::  islice_dvrp                        !<
1415    INTEGER(iwp) ::  nx_dvrp                            !<
1416    INTEGER(iwp) ::  nxl_dvrp                           !<
1417    INTEGER(iwp) ::  nxr_dvrp                           !<
1418    INTEGER(iwp) ::  ny_dvrp                            !<
1419    INTEGER(iwp) ::  nyn_dvrp                           !<
1420    INTEGER(iwp) ::  nys_dvrp                           !<
1421    INTEGER(iwp) ::  nz_dvrp, pathlines_fadeintime = 5  !<
1422    INTEGER(iwp) ::  pathlines_fadeouttime = 5          !<
1423    INTEGER(iwp) ::  pathlines_linecount = 1000         !<
1424    INTEGER(iwp) ::  pathlines_maxhistory = 40          !<
1425    INTEGER(iwp) ::  pathlines_wavecount = 10           !<
1426    INTEGER(iwp) ::  pathlines_wavetime = 50            !<
1427    INTEGER(iwp) ::  vc_gradient_normals = 0            !<
1428    INTEGER(iwp) ::  vc_mode = 0                        !<
1429    INTEGER(iwp) ::  vc_size_x = 2                      !<
1430    INTEGER(iwp) ::  vc_size_y = 2                      !<
1431    INTEGER(iwp) ::  vc_size_z = 2                      !<
1432
1433    INTEGER(iwp), DIMENSION(10) ::  slicer_position_dvrp  !<
1434
1435    LOGICAL ::  cyclic_dvrp = .FALSE.                      !<
1436    LOGICAL ::  dvrp_overlap                               !<
1437    LOGICAL ::  dvrp_total_overlap                         !<
1438    LOGICAL ::  local_dvrserver_running                    !<
1439    LOGICAL ::  lock_steering_update = .FALSE.             !<
1440    LOGICAL ::  use_seperate_pe_for_dvrp_output = .FALSE.  !<
1441
1442    REAL(wp) ::  clip_dvrp_l = 9999999.9_wp  !<
1443    REAL(wp) ::  clip_dvrp_n = 9999999.9_wp  !<
1444    REAL(wp) ::  clip_dvrp_r = 9999999.9_wp  !<
1445    REAL(wp) ::  clip_dvrp_s = 9999999.9_wp  !<
1446    REAL(wp) ::  superelevation = 1.0_wp     !<
1447    REAL(wp) ::  superelevation_x = 1.0_wp   !<
1448    REAL(wp) ::  superelevation_y = 1.0_wp   !<
1449    REAL(wp) ::  vc_alpha = 38.0_wp          !<
1450
1451    REAL(wp), DIMENSION(2) ::  color_interval = (/ 0.0_wp, 1.0_wp /)     !<
1452    REAL(wp), DIMENSION(2) ::  dvrpsize_interval = (/ 0.0_wp, 1.0_wp /)  !<
1453
1454    REAL(wp), DIMENSION(3) ::  groundplate_color = (/ 0.0_wp, 0.6_wp, 0.0_wp /)  !<
1455    REAL(wp), DIMENSION(3) ::  topography_color = (/ 0.8_wp, 0.7_wp, 0.6_wp /)   !<
1456
1457    REAL(wp), DIMENSION(2,10) ::  slicer_range_limits_dvrp  !<
1458
1459    REAL(wp), DIMENSION(3,10) ::  isosurface_color  !<
1460
1461    REAL(sp), DIMENSION(2,100) ::  interval_values_dvrp          !<
1462    REAL(sp), DIMENSION(2,100) ::  interval_values_dvrp_prt      !<
1463    REAL(sp), DIMENSION(2,100) ::  interval_h_dvrp               !<
1464    REAL(sp), DIMENSION(2,100) ::  interval_h_dvrp_prt           !<
1465    REAL(sp), DIMENSION(2,100) ::  interval_l_dvrp = 0.5_sp      !<
1466    REAL(sp), DIMENSION(2,100) ::  interval_l_dvrp_prt = 0.5_sp  !<
1467    REAL(sp), DIMENSION(2,100) ::  interval_s_dvrp = 1.0_sp      !<
1468    REAL(sp), DIMENSION(2,100) ::  interval_s_dvrp_prt = 1.0_sp  !<
1469    REAL(sp), DIMENSION(2,100) ::  interval_a_dvrp = 0.0_sp      !<
1470    REAL(sp), DIMENSION(2,100) ::  interval_a_dvrp_prt = 0.0_sp  !<
1471
1472    DATA  slicer_range_limits_dvrp / -1.0_wp, 1.0_wp, -1.0_wp, 1.0_wp, -1.0_wp, 1.0_wp, &  !<
1473                                     -1.0_wp, 1.0_wp, -1.0_wp, 1.0_wp, -1.0_wp, 1.0_wp, &
1474                                     -1.0_wp, 1.0_wp, -1.0_wp, 1.0_wp, -1.0_wp, 1.0_wp, &
1475                                     -1.0_wp, 1.0_wp /
1476
1477    DATA  isosurface_color / 0.9_wp, 0.9_wp, 0.9_wp,  0.8_wp, 0.1_wp, 0.1_wp,  0.1_wp, 0.1_wp, 0.8_wp, &  !<
1478                             0.1_wp, 0.8_wp, 0.1_wp,  0.6_wp, 0.1_wp, 0.1_wp,  0.1_wp, 0.1_wp, 0.6_wp, &
1479                             0.1_wp, 0.6_wp, 0.1_wp,  0.4_wp, 0.1_wp, 0.1_wp,  0.1_wp, 0.1_wp, 0.4_wp, &
1480                             0.1_wp, 0.4_wp, 0.1_wp /
1481
1482    DATA  interval_h_dvrp / 270.0_wp, 225.0_wp, 225.0_wp, 180.0_wp, 70.0_wp, 25.0_wp, &  !<
1483                            25.0_wp, -25.0_wp, 192 * 0.0_wp /
1484
1485    DATA  interval_h_dvrp_prt / 270.0_wp, 225.0_wp, 225.0_wp, 180.0_wp, 70.0_wp, 25.0_wp, &  !<
1486                                25.0_wp, -25.0_wp, 192 * 0.0_wp /
1487
1488    REAL(sp), DIMENSION(:), ALLOCATABLE ::  xcoor_dvrp  !<
1489    REAL(sp), DIMENSION(:), ALLOCATABLE ::  ycoor_dvrp  !<
1490    REAL(sp), DIMENSION(:), ALLOCATABLE ::  zcoor_dvrp  !<
1491
1492    TYPE steering
1493       CHARACTER (LEN=24) ::  name  !<
1494       REAL(sp)           ::  min   !<
1495       REAL(sp)           ::  max   !<
1496       INTEGER(iwp)       ::  imin  !<
1497       INTEGER(iwp)       ::  imax  !<
1498    END TYPE steering
1499
1500    TYPE(steering), DIMENSION(:), ALLOCATABLE ::  steering_dvrp  !<
1501
1502    SAVE
1503
1504 END MODULE dvrp_variables
1505
1506
1507!------------------------------------------------------------------------------!
1508! Description:
1509! ------------
1510!> Definition of grid spacings.
1511!------------------------------------------------------------------------------!
1512 MODULE grid_variables
1513
1514    USE kinds
1515
1516    REAL(wp) ::  ddx          !<
1517    REAL(wp) ::  ddx2         !<
1518    REAL(wp) ::  dx = 1.0_wp  !<
1519    REAL(wp) ::  dx2          !<
1520    REAL(wp) ::  ddy          !<
1521    REAL(wp) ::  ddy2         !<
1522    REAL(wp) ::  dy = 1.0_wp  !<
1523    REAL(wp) ::  dy2          !<
1524
1525    REAL(wp), DIMENSION(:), ALLOCATABLE ::  ddx2_mg  !<
1526    REAL(wp), DIMENSION(:), ALLOCATABLE ::  ddy2_mg  !<
1527
1528    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  fwxm        !<
1529    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  fwxp        !<
1530    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  fwym        !<
1531    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  fwyp        !<
1532    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  fxm         !<
1533    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  fxp         !<
1534    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  fym         !<
1535    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  fyp         !<
1536    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  wall_e_x    !<
1537    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  wall_e_y    !<
1538    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  wall_u      !<
1539    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  wall_v      !<
1540    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  wall_w_x    !<
1541    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  wall_w_y    !<
1542    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  zu_s_inner  !< 
1543    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  zw_w_inner  !<
1544                                             
1545    SAVE
1546
1547 END MODULE grid_variables
1548
1549
1550!------------------------------------------------------------------------------!
1551! Description:
1552! ------------
1553!> Definition of array bounds, number of gridpoints, and wall flag arrays.
1554!------------------------------------------------------------------------------!
1555 MODULE indices
1556
1557    USE kinds
1558
1559    INTEGER(iwp) ::  i_left       !<
1560    INTEGER(iwp) ::  i_right      !<
1561    INTEGER(iwp) ::  j_north      !<
1562    INTEGER(iwp) ::  j_south      !<
1563    INTEGER(iwp) ::  nbgp = 3     !<
1564    INTEGER(iwp) ::  ngp_sums     !<
1565    INTEGER(iwp) ::  ngp_sums_ls  !<
1566    INTEGER(iwp) ::  nnx          !<
1567    INTEGER(iwp) ::  nx = 0       !<
1568    INTEGER(iwp) ::  nx_a         !<
1569    INTEGER(iwp) ::  nx_o         !<
1570    INTEGER(iwp) ::  nxl          !<
1571    INTEGER(iwp) ::  nxlg         !<
1572    INTEGER(iwp) ::  nxlu         !<
1573    INTEGER(iwp) ::  nxr          !<
1574    INTEGER(iwp) ::  nxrg         !<
1575    INTEGER(iwp) ::  nx_on_file   !<
1576    INTEGER(iwp) ::  nny          !<
1577    INTEGER(iwp) ::  ny = 0       !<
1578    INTEGER(iwp) ::  ny_a         !<
1579    INTEGER(iwp) ::  ny_o         !<
1580    INTEGER(iwp) ::  nyn          !<
1581    INTEGER(iwp) ::  nyng         !<
1582    INTEGER(iwp) ::  nys          !<
1583    INTEGER(iwp) ::  nysg         !<
1584    INTEGER(iwp) ::  nysv         !<
1585    INTEGER(iwp) ::  ny_on_file   !<
1586    INTEGER(iwp) ::  nnz          !<
1587    INTEGER(iwp) ::  nz = 0       !<
1588    INTEGER(iwp) ::  nzb          !<
1589    INTEGER(iwp) ::  nzb_diff     !<
1590    INTEGER(iwp) ::  nzb_max      !<
1591    INTEGER(iwp) ::  nzt          !<
1592    INTEGER(iwp) ::  nzt_diff     !<
1593
1594    INTEGER(idp), DIMENSION(:), ALLOCATABLE ::  ngp_3d        !<
1595    INTEGER(idp), DIMENSION(:), ALLOCATABLE ::  ngp_3d_inner  !< ! need to have 64 bit for grids > 2E9
1596                   
1597    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  ngp_2dh  !<
1598    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  nxl_mg   !<
1599    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  nxr_mg   !<
1600    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  nyn_mg   !<
1601    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  nys_mg   !<
1602    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  nzt_mg   !<
1603
1604
1605    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  ngp_2dh_outer     !<
1606    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  ngp_2dh_s_inner   !<
1607    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  mg_loc_ind        !<
1608    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  nzb_diff_s_inner  !<
1609    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  nzb_diff_s_outer  !<
1610    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  nzb_diff_u        !<
1611    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  nzb_diff_v        !<
1612    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  nzb_inner         !<
1613    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  nzb_outer         !<
1614    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  nzb_s_inner       !<
1615    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  nzb_s_outer       !<
1616    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  nzb_u_inner       !<
1617    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  nzb_u_outer       !<
1618    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  nzb_v_inner       !<
1619    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  nzb_v_outer       !<
1620    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  nzb_w_inner       !<
1621    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  nzb_w_outer       !<
1622
1623    INTEGER(iwp), DIMENSION(:,:,:), POINTER ::  flags  !<
1624
1625    INTEGER(iwp), DIMENSION(:,:,:), ALLOCATABLE ::  wall_flags_0   !<
1626    INTEGER(iwp), DIMENSION(:,:,:), ALLOCATABLE ::  wall_flags_00  !<
1627
1628    INTEGER(iwp), DIMENSION(:,:,:), ALLOCATABLE,  TARGET ::  wall_flags_1   !<
1629    INTEGER(iwp), DIMENSION(:,:,:), ALLOCATABLE,  TARGET ::  wall_flags_2   !<
1630    INTEGER(iwp), DIMENSION(:,:,:), ALLOCATABLE,  TARGET ::  wall_flags_3   !<
1631    INTEGER(iwp), DIMENSION(:,:,:), ALLOCATABLE,  TARGET ::  wall_flags_4   !<
1632    INTEGER(iwp), DIMENSION(:,:,:), ALLOCATABLE,  TARGET ::  wall_flags_5   !<
1633    INTEGER(iwp), DIMENSION(:,:,:), ALLOCATABLE,  TARGET ::  wall_flags_6   !<
1634    INTEGER(iwp), DIMENSION(:,:,:), ALLOCATABLE,  TARGET ::  wall_flags_7   !<
1635    INTEGER(iwp), DIMENSION(:,:,:), ALLOCATABLE,  TARGET ::  wall_flags_8   !<
1636    INTEGER(iwp), DIMENSION(:,:,:), ALLOCATABLE,  TARGET ::  wall_flags_9   !<
1637    INTEGER(iwp), DIMENSION(:,:,:), ALLOCATABLE,  TARGET ::  wall_flags_10  !<
1638
1639    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  rflags_s_inner  !<
1640    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  rflags_invers   !<
1641
1642    SAVE
1643
1644 END MODULE indices
1645
1646
1647!------------------------------------------------------------------------------!
1648! Description:
1649! ------------
1650!> Interfaces for special subroutines which use optional parameters.
1651!------------------------------------------------------------------------------!
1652 MODULE interfaces
1653
1654    INTERFACE
1655
1656!------------------------------------------------------------------------------!
1657! Description:
1658! ------------
1659!> @todo Missing subroutine description.
1660!------------------------------------------------------------------------------!
1661       SUBROUTINE global_min_max ( i1, i2, j1, j2, k1, k2, feld, mode, offset, &
1662                                   wert, wert_ijk, wert1, wert1_ijk )
1663
1664          USE kinds
1665
1666          CHARACTER (LEN=*), INTENT(IN) ::  mode                     !<
1667          INTEGER(iwp), INTENT(IN)      ::  i1                       !<
1668          INTEGER(iwp), INTENT(IN)      ::  i2                       !<
1669          INTEGER(iwp), INTENT(IN)      ::  j1                       !<
1670          INTEGER(iwp), INTENT(IN)      ::  j2                       !<
1671          INTEGER(iwp), INTENT(IN)      ::  k1                       !<
1672          INTEGER(iwp), INTENT(IN)      ::  k2                       !< 
1673          INTEGER(iwp)                  ::  wert_ijk(3)              !<
1674          INTEGER(iwp), OPTIONAL        ::  wert1_ijk(3)             !< 
1675          REAL(wp)                      ::  offset                   !<
1676          REAL(wp)                      ::  wert                     !<
1677          REAL(wp), OPTIONAL            ::  wert1                    !<
1678          REAL(wp), INTENT(IN)          ::  feld(i1:i2,j1:j2,k1:k2)  !<
1679
1680       END SUBROUTINE global_min_max
1681
1682    END INTERFACE
1683
1684    SAVE
1685
1686 END MODULE interfaces
1687
1688
1689!------------------------------------------------------------------------------!
1690! Description:
1691! ------------
1692!> Interfaces for subroutines with pointer arguments called in
1693!> prognostic_equations.
1694!------------------------------------------------------------------------------!
1695 MODULE pointer_interfaces
1696
1697    INTERFACE
1698
1699!------------------------------------------------------------------------------!
1700! Description:
1701! ------------
1702!> @todo Missing subroutine description.
1703!------------------------------------------------------------------------------!
1704       SUBROUTINE advec_s_bc( sk, sk_char )
1705
1706          USE kinds
1707
1708          CHARACTER (LEN=*), INTENT(IN) ::  sk_char  !<
1709#if defined( __nopointer )
1710          REAL(wp), DIMENSION(:,:,:) ::  sk  !<
1711#else
1712          REAL(wp), DIMENSION(:,:,:), POINTER ::  sk  !<
1713#endif
1714       END SUBROUTINE advec_s_bc
1715
1716    END INTERFACE
1717
1718    SAVE
1719
1720 END MODULE pointer_interfaces
1721
1722
1723!------------------------------------------------------------------------------!
1724! Description:
1725! ------------
1726!> Definition of variables for the 1D-model.
1727!------------------------------------------------------------------------------!
1728 MODULE model_1d
1729
1730    USE kinds
1731
1732    INTEGER(iwp) ::  current_timestep_number_1d = 0  !<
1733    INTEGER(iwp) ::  damp_level_ind_1d               !<
1734
1735    LOGICAL ::  run_control_header_1d = .FALSE.  !<
1736    LOGICAL ::  stop_dt_1d = .FALSE.             !<
1737
1738    REAL(wp) ::  damp_level_1d = -1.0_wp       !<
1739    REAL(wp) ::  dt_1d = 60.0_wp               !<
1740    REAL(wp) ::  dt_max_1d = 300.0_wp          !<
1741    REAL(wp) ::  dt_pr_1d = 9999999.9_wp       !<
1742    REAL(wp) ::  dt_run_control_1d = 60.0_wp   !<
1743    REAL(wp) ::  end_time_1d = 864000.0_wp     !<
1744    REAL(wp) ::  old_dt_1d = 1.0E-10_wp        !<
1745    REAL(wp) ::  qs1d                          !<
1746    REAL(wp) ::  simulated_time_1d = 0.0_wp    !<
1747    REAL(wp) ::  time_pr_1d = 0.0_wp           !<
1748    REAL(wp) ::  time_run_control_1d = 0.0_wp  !<
1749    REAL(wp) ::  ts1d                          !<
1750    REAL(wp) ::  us1d                          !<
1751    REAL(wp) ::  usws1d                        !<
1752    REAL(wp) ::  vsws1d                        !<               
1753    REAL(wp) ::  z01d                          !<
1754    REAL(wp) ::  z0h1d                         !<
1755
1756
1757    REAL(wp), DIMENSION(:), ALLOCATABLE ::  e1d      !<
1758    REAL(wp), DIMENSION(:), ALLOCATABLE ::  e1d_p    !<
1759    REAL(wp), DIMENSION(:), ALLOCATABLE ::  kh1d     !<
1760    REAL(wp), DIMENSION(:), ALLOCATABLE ::  km1d     !<
1761    REAL(wp), DIMENSION(:), ALLOCATABLE ::  l_black  !<
1762    REAL(wp), DIMENSION(:), ALLOCATABLE ::  l1d      !<
1763    REAL(wp), DIMENSION(:), ALLOCATABLE ::  rif1d    !<
1764    REAL(wp), DIMENSION(:), ALLOCATABLE ::  te_e     !<
1765    REAL(wp), DIMENSION(:), ALLOCATABLE ::  te_em    !<
1766    REAL(wp), DIMENSION(:), ALLOCATABLE ::  te_u     !<
1767    REAL(wp), DIMENSION(:), ALLOCATABLE ::  te_um    !<
1768    REAL(wp), DIMENSION(:), ALLOCATABLE ::  te_v     !<
1769    REAL(wp), DIMENSION(:), ALLOCATABLE ::  te_vm    !<
1770    REAL(wp), DIMENSION(:), ALLOCATABLE ::  u1d      !<
1771    REAL(wp), DIMENSION(:), ALLOCATABLE ::  u1d_p    !<
1772    REAL(wp), DIMENSION(:), ALLOCATABLE ::  v1d      !<
1773    REAL(wp), DIMENSION(:), ALLOCATABLE ::  v1d_p    !<
1774
1775    SAVE
1776
1777 END MODULE model_1d
1778
1779
1780
1781
1782!------------------------------------------------------------------------------!
1783! Description:
1784! ------------
1785!> Definition of variables which define processor topology and the exchange of
1786!> ghost point layers. This modules must be placed in all routines which contain
1787!> MPI-calls.
1788!------------------------------------------------------------------------------!
1789 MODULE pegrid
1790
1791    USE kinds
1792
1793#if defined( __parallel )
1794#if defined( __mpifh )
1795    INCLUDE "mpif.h"
1796#else
1797    USE MPI
1798#endif
1799#endif
1800    CHARACTER(LEN=2) ::  send_receive = 'al'
1801    CHARACTER(LEN=7) ::  myid_char = ''
1802   
1803    INTEGER(iwp) ::  acc_rank                    !<
1804    INTEGER(iwp) ::  comm1dx                     !<
1805    INTEGER(iwp) ::  comm1dy                     !<
1806    INTEGER(iwp) ::  comm2d                      !<
1807    INTEGER(iwp) ::  comm_inter                  !<
1808    INTEGER(iwp) ::  comm_palm                   !<
1809    INTEGER(iwp) ::  id_inflow = 0               !<
1810    INTEGER(iwp) ::  id_outflow = 0              !< myidx of procs at outflow (turbulent outflow method)
1811    INTEGER(iwp) ::  id_outflow_source = 0       !< myidx of procs including ouflow source plane (turbulent outflow method)
1812    INTEGER(iwp) ::  id_recycling = 0            !<
1813    INTEGER(iwp) ::  ierr                        !<
1814    INTEGER(iwp) ::  myid = 0                    !<
1815    INTEGER(iwp) ::  myidx = 0                   !<
1816    INTEGER(iwp) ::  myidy = 0                   !<
1817    INTEGER(iwp) ::  ndim = 2                    !<
1818    INTEGER(iwp) ::  ngp_a                       !<
1819    INTEGER(iwp) ::  ngp_o                       !<
1820    INTEGER(iwp) ::  ngp_xy                      !<
1821    INTEGER(iwp) ::  ngp_y                       !<
1822    INTEGER(iwp) ::  npex = -1                   !<
1823    INTEGER(iwp) ::  npey = -1                   !<
1824    INTEGER(iwp) ::  numprocs = 1                !<
1825    INTEGER(iwp) ::  numprocs_previous_run = -1  !<
1826    INTEGER(iwp) ::  num_acc_per_node = 0        !<
1827    INTEGER(iwp) ::  pleft                       !<
1828    INTEGER(iwp) ::  pnorth                      !<
1829    INTEGER(iwp) ::  pright                      !<
1830    INTEGER(iwp) ::  psouth                      !<
1831    INTEGER(iwp) ::  req_count = 0               !<
1832    INTEGER(iwp) ::  sendrecvcount_xy            !<
1833    INTEGER(iwp) ::  sendrecvcount_yz            !<
1834    INTEGER(iwp) ::  sendrecvcount_zx            !<
1835    INTEGER(iwp) ::  sendrecvcount_zyd           !<
1836    INTEGER(iwp) ::  sendrecvcount_yxd           !<
1837    INTEGER(iwp) ::  target_id                   !<
1838    INTEGER(iwp) ::  tasks_per_node = -9999      !<
1839    INTEGER(iwp) ::  threads_per_task = 1        !<
1840    INTEGER(iwp) ::  type_x                      !<
1841    INTEGER(iwp) ::  type_xy                     !<
1842    INTEGER(iwp) ::  type_y                      !<
1843
1844    INTEGER(iwp) ::  pdims(2) = 1  !<
1845    INTEGER(iwp) ::  req(100)      !<
1846
1847    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  hor_index_bounds               !<
1848    INTEGER(iwp), DIMENSION(:,:), ALLOCATABLE ::  hor_index_bounds_previous_run  !<
1849
1850    LOGICAL ::  background_communication =.FALSE.  !<
1851    LOGICAL ::  collective_wait = .FALSE.          !<
1852    LOGICAL ::  sendrecv_in_background = .FALSE.   !<
1853
1854#if defined( __parallel )
1855#if defined( __mpi2 )
1856    CHARACTER (LEN=MPI_MAX_PORT_NAME) ::  port_name  !<
1857#endif
1858
1859    INTEGER(iwp) ::  ibuf(12)                 !<
1860    INTEGER(iwp) ::  pcoord(2)                !<
1861    INTEGER(iwp) ::  status(MPI_STATUS_SIZE)  !<
1862   
1863    INTEGER(iwp), DIMENSION(MPI_STATUS_SIZE,100) ::  wait_stat  !<
1864   
1865    INTEGER(iwp) ::  ngp_yz_int   !<
1866    INTEGER(iwp) ::  type_xz_int  !<
1867    INTEGER(iwp) ::  type_yz_int  !<
1868
1869    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  ngp_xz      !<
1870    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  ngp_yz      !<
1871    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  type_x_int  !<
1872    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  type_xz     !<
1873    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  type_y_int  !<
1874    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  type_yz     !<
1875
1876    LOGICAL ::  left_border_pe  = .FALSE.  !<
1877    LOGICAL ::  north_border_pe = .FALSE.  !<
1878    LOGICAL ::  reorder = .TRUE.           !<
1879    LOGICAL ::  right_border_pe = .FALSE.  !<
1880    LOGICAL ::  south_border_pe = .FALSE.  !<
1881
1882    LOGICAL, DIMENSION(2) ::  cyclic = (/ .TRUE. , .TRUE. /)  !<
1883    LOGICAL, DIMENSION(2) ::  remain_dims                     !<
1884#endif
1885
1886    SAVE
1887
1888 END MODULE pegrid
1889
1890
1891!------------------------------------------------------------------------------!
1892! Description:
1893! ------------
1894!> Definition of variables which control PROFIL-output.
1895!------------------------------------------------------------------------------!
1896 MODULE profil_parameter
1897
1898    USE kinds
1899
1900    INTEGER(iwp), PARAMETER ::  crmax = 100  !<
1901
1902    CHARACTER (LEN=20), DIMENSION(20) ::  cross_ts_profiles = &  !<
1903                           (/ ' E E*               ', ' dt                 ', &
1904                              ' u* w*              ', ' th*                ', &
1905                              ' umax vmax wmax     ', ' div_old div_new    ', &
1906                              ' z_i_wpt z_i_pt     ', ' w"pt"0 w"pt" wpt   ', &
1907                              ' pt(0) pt(zp)       ', ' splux spluy spluz  ', &
1908                              ' L                  ',                         &
1909                            ( '                    ', i9 = 1, 9 ) /)
1910
1911    CHARACTER (LEN=100), DIMENSION(crmax) ::  cross_profiles = &  !<
1912                           (/ ' u v                           ', &
1913                              ' pt                            ', &
1914                              ' w"pt" w*pt* w*pt*BC wpt wptBC ', &
1915                              ' w"u" w*u* wu w"v" w*v* wv     ', &
1916                              ' km kh                         ', &
1917                              ' l                             ', &
1918                         ( '                               ', i9 = 1, 94 ) /)
1919
1920    INTEGER(iwp) ::  profile_columns = 2  !<
1921    INTEGER(iwp) ::  profile_rows = 3     !<
1922    INTEGER(iwp) ::  profile_number = 0   !<
1923
1924    INTEGER(iwp) ::  cross_ts_numbers(crmax,crmax) = 0  !<
1925    INTEGER(iwp) ::  cross_ts_number_count(crmax) = 0   !<
1926    INTEGER(iwp) ::  dopr_index(300) = 0                !<
1927    INTEGER(iwp) ::  dopr_initial_index(300) = 0        !<
1928    INTEGER(iwp) ::  dots_crossindex(100) = 0           !<
1929    INTEGER(iwp) ::  dots_index(100) = 0                !<
1930               
1931
1932    REAL(wp) ::  cross_ts_uymax(20) = &                    !<
1933                             (/ 999.999_wp, 999.999_wp, 999.999_wp, 999.999_wp, 999.999_wp,   &
1934                                999.999_wp, 999.999_wp, 999.999_wp, 999.999_wp, 999.999_wp,   &
1935                                999.999_wp, 999.999_wp, 999.999_wp, 999.999_wp, 999.999_wp,   &
1936                                999.999_wp, 999.999_wp, 999.999_wp, 999.999_wp, 999.999_wp /)
1937    REAL(wp) ::  cross_ts_uymax_computed(20) = 999.999_wp  !<
1938    REAL(wp) ::  cross_ts_uymin(20) = &                    !<
1939                             (/ 999.999_wp, 999.999_wp, 999.999_wp,  -5.000_wp, 999.999_wp,   &
1940                                999.999_wp,   0.000_wp, 999.999_wp, 999.999_wp, 999.999_wp,   &
1941                                999.999_wp, 999.999_wp, 999.999_wp, 999.999_wp, 999.999_wp,   &
1942                                999.999_wp, 999.999_wp, 999.999_wp, 999.999_wp, 999.999_wp /)
1943    REAL(wp) ::  cross_ts_uymin_computed(20) = 999.999_wp  !<
1944
1945    SAVE
1946
1947 END MODULE profil_parameter
1948
1949!------------------------------------------------------------------------------!
1950! Description:
1951! ------------
1952!> Definition of statistical quantities, e.g. global sums.
1953!------------------------------------------------------------------------------!
1954 MODULE statistics
1955
1956    USE kinds
1957
1958    CHARACTER (LEN=40) ::  region(0:9)  !<
1959   
1960    INTEGER(iwp) ::  pr_palm = 130          !<
1961    INTEGER(iwp) ::  statistic_regions = 0  !<
1962   
1963    INTEGER(iwp) ::  u_max_ijk(3) = -1  !<
1964    INTEGER(iwp) ::  v_max_ijk(3) = -1  !<
1965    INTEGER(iwp) ::  w_max_ijk(3) = -1  !<
1966   
1967    LOGICAL ::  flow_statistics_called = .FALSE.  !<
1968   
1969    REAL(wp) ::  u_max  !<
1970    REAL(wp) ::  v_max  !<
1971    REAL(wp) ::  w_max  !<
1972   
1973    REAL(wp), DIMENSION(:), ALLOCATABLE ::  mean_surface_level_height  !<
1974    REAL(wp), DIMENSION(:), ALLOCATABLE ::  sums_divnew_l              !<
1975    REAL(wp), DIMENSION(:), ALLOCATABLE ::  sums_divold_l              !<
1976    REAL(wp), DIMENSION(:), ALLOCATABLE ::  weight_substep             !<
1977    REAL(wp), DIMENSION(:), ALLOCATABLE ::  weight_pres                !<
1978   
1979    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  sums             !<
1980    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  sums_wsts_bc_l   !<
1981    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  ts_value         !<
1982    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  sums_wsus_ws_l   !<
1983    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  sums_wsvs_ws_l   !<
1984    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  sums_us2_ws_l    !<
1985    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  sums_vs2_ws_l    !<
1986    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  sums_ws2_ws_l    !<
1987    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  sums_wsnrs_ws_l  !<
1988    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  sums_wspts_ws_l  !<
1989    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  sums_wsqs_ws_l   !<
1990    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  sums_wsqrs_ws_l  !<
1991    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  sums_wssas_ws_l  !<
1992    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  sums_wsss_ws_l   !<
1993    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  sums_ls_l        !<
1994                                             
1995    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  hom_sum             !<
1996    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  rmask               !<
1997    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  sums_l              !<
1998    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  sums_l_l            !<
1999    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  sums_up_fraction_l  !<
2000
2001    REAL(wp), DIMENSION(:,:,:,:), ALLOCATABLE ::  hom  !<
2002
2003    SAVE
2004
2005 END MODULE statistics
2006
2007
2008
2009!------------------------------------------------------------------------------!
2010! Description:
2011! ------------
2012!> Definition of indices for transposed arrays.
2013!------------------------------------------------------------------------------!
2014 MODULE transpose_indices
2015
2016    USE kinds
2017
2018    INTEGER(iwp) ::  nxl_y   !<
2019    INTEGER(iwp) ::  nxl_yd  !<
2020    INTEGER(iwp) ::  nxl_z   !<
2021    INTEGER(iwp) ::  nxr_y   !<
2022    INTEGER(iwp) ::  nxr_yd  !<
2023    INTEGER(iwp) ::  nxr_z   !<
2024    INTEGER(iwp) ::  nyn_x   !<
2025    INTEGER(iwp) ::  nyn_z   !<
2026    INTEGER(iwp) ::  nys_x   !<
2027    INTEGER(iwp) ::  nys_z   !<
2028    INTEGER(iwp) ::  nzb_x   !<
2029    INTEGER(iwp) ::  nzb_y   !<
2030    INTEGER(iwp) ::  nzb_yd  !<
2031    INTEGER(iwp) ::  nzt_x   !<
2032    INTEGER(iwp) ::  nzt_y   !<
2033    INTEGER(iwp) ::  nzt_yd  !<
2034               
2035    SAVE
2036
2037 END MODULE transpose_indices
Note: See TracBrowser for help on using the repository browser.