source: palm/trunk/SOURCE/lagrangian_particle_model_mod.f90 @ 4028

Last change on this file since 4028 was 4028, checked in by schwenkel, 5 years ago

Further modularization of particle code components

  • Property svn:keywords set to Id
File size: 336.8 KB
Line 
1!> @file lagrangian_particle_model_mod.f90
2!------------------------------------------------------------------------------!
3! This file is part of the PALM model system.
4!
5! PALM is free software: you can redistribute it and/or modify it under the
6! terms of the GNU General Public License as published by the Free Software
7! Foundation, either version 3 of the License, or (at your option) any later
8! version.
9!
10! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
11! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
12! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
13!
14! You should have received a copy of the GNU General Public License along with
15! PALM. If not, see <http://www.gnu.org/licenses/>.
16!
17! Copyright 1997-2019 Leibniz Universitaet Hannover
18!------------------------------------------------------------------------------!
19!
20! Current revisions:
21! ------------------
22!
23!
24! Former revisions:
25! -----------------
26! $Id: lagrangian_particle_model_mod.f90 4028 2019-06-13 12:21:37Z schwenkel $
27! Further modularization of particle code components
28!
29! 4020 2019-06-06 14:57:48Z schwenkel
30! Removing submodules
31!
32! 4018 2019-06-06 13:41:50Z eckhard
33! Bugfix for former revision
34!
35! 4017 2019-06-06 12:16:46Z schwenkel
36! Modularization of all lagrangian particle model code components
37!
38! 3655 2019-01-07 16:51:22Z knoop
39! bugfix to guarantee correct particle releases in case that the release
40! interval is smaller than the model timestep
41!
42! 2801 2018-02-14 16:01:55Z thiele
43! Changed lpm from subroutine to module.
44! Introduce particle transfer in nested models.
45!
46! 2718 2018-01-02 08:49:38Z maronga
47! Corrected "Former revisions" section
48!
49! 2701 2017-12-15 15:40:50Z suehring
50! Changes from last commit documented
51!
52! 2698 2017-12-14 18:46:24Z suehring
53! Grid indices passed to lpm_boundary_conds. (responsible Philipp Thiele)
54!
55! 2696 2017-12-14 17:12:51Z kanani
56! Change in file header (GPL part)
57!
58! 2606 2017-11-10 10:36:31Z schwenkel
59! Changed particle box locations: center of particle box now coincides
60! with scalar grid point of same index.
61! Renamed module and subroutines: lpm_pack_arrays_mod -> lpm_pack_and_sort_mod
62! lpm_pack_all_arrays -> lpm_sort_in_subboxes, lpm_pack_arrays -> lpm_pack
63! lpm_sort -> lpm_sort_timeloop_done
64!
65! 2418 2017-09-06 15:24:24Z suehring
66! Major bugfixes in modeling SGS particle speeds (since revision 1359).
67! Particle sorting added to distinguish between already completed and
68! non-completed particles.
69!
70! 2263 2017-06-08 14:59:01Z schwenkel
71! Implemented splitting and merging algorithm
72!
73! 2233 2017-05-30 18:08:54Z suehring
74!
75! 2232 2017-05-30 17:47:52Z suehring
76! Adjustments to new topography concept
77!
78! 2000 2016-08-20 18:09:15Z knoop
79! Forced header and separation lines into 80 columns
80!
81! 1936 2016-06-13 13:37:44Z suehring
82! Call routine for deallocation of unused memory.
83! Formatting adjustments
84!
85! 1929 2016-06-09 16:25:25Z suehring
86! Call wall boundary conditions only if particles are in the vertical range of
87! topography.
88!
89! 1822 2016-04-07 07:49:42Z hoffmann
90! Tails removed.
91!
92! Initialization of sgs model not necessary for the use of cloud_droplets and
93! use_sgs_for_particles.
94!
95! lpm_release_set integrated.
96!
97! Unused variabled removed.
98!
99! 1682 2015-10-07 23:56:08Z knoop
100! Code annotations made doxygen readable
101!
102! 1416 2014-06-04 16:04:03Z suehring
103! user_lpm_advec is called for each gridpoint.
104! Bugfix: in order to prevent an infinite loop, time_loop_done is set .TRUE.
105! at the head of the do-loop. 
106!
107! 1359 2014-04-11 17:15:14Z hoffmann
108! New particle structure integrated.
109! Kind definition added to all floating point numbers.
110!
111! 1320 2014-03-20 08:40:49Z raasch
112! ONLY-attribute added to USE-statements,
113! kind-parameters added to all INTEGER and REAL declaration statements,
114! kinds are defined in new module kinds,
115! revision history before 2012 removed,
116! comment fields (!:) to be used for variable explanations added to
117! all variable declaration statements
118!
119! 1318 2014-03-17 13:35:16Z raasch
120! module interfaces removed
121!
122! 1036 2012-10-22 13:43:42Z raasch
123! code put under GPL (PALM 3.9)
124!
125! 851 2012-03-15 14:32:58Z raasch
126! Bugfix: resetting of particle_mask and tail mask moved from routine
127! lpm_exchange_horiz to here (end of sub-timestep loop)
128!
129! 849 2012-03-15 10:35:09Z raasch
130! original routine advec_particles split into several subroutines and renamed
131! lpm
132!
133! 831 2012-02-22 00:29:39Z raasch
134! thermal_conductivity_l and diff_coeff_l now depend on temperature and
135! pressure
136!
137! 828 2012-02-21 12:00:36Z raasch
138! fast hall/wang kernels with fixed radius/dissipation classes added,
139! particle feature color renamed class, routine colker renamed
140! recalculate_kernel,
141! lower limit for droplet radius changed from 1E-7 to 1E-8
142!
143! Bugfix: transformation factor for dissipation changed from 1E5 to 1E4
144!
145! 825 2012-02-19 03:03:44Z raasch
146! droplet growth by condensation may include curvature and solution effects,
147! initialisation of temporary particle array for resorting removed,
148! particle attributes speed_x|y|z_sgs renamed rvar1|2|3,
149! module wang_kernel_mod renamed lpm_collision_kernels_mod,
150! wang_collision_kernel renamed wang_kernel
151!
152!
153! Revision 1.1  1999/11/25 16:16:06  raasch
154! Initial revision
155!
156!
157! Description:
158! ------------
159!>
160!------------------------------------------------------------------------------!
161 MODULE lagrangian_particle_model_mod
162       
163    USE, INTRINSIC ::  ISO_C_BINDING
164       
165    USE arrays_3d,                                                             &
166        ONLY:  de_dx, de_dy, de_dz, dzw, zu, zw,  ql_c, ql_v, ql_vp, hyp,      &
167               pt, q, exner, ql, diss, e, u, v, w, km, ql_1, ql_2     
168 
169    USE averaging,                                                             &
170        ONLY:  ql_c_av, pr_av, pc_av, ql_vp_av, ql_v_av
171       
172    USE basic_constants_and_equations_mod,                                     &
173        ONLY: molecular_weight_of_solute, molecular_weight_of_water, magnus,   &
174              pi, rd_d_rv, rho_l, r_v, rho_s, vanthoff, l_v, kappa, g               
175       
176    USE control_parameters,                                                    &
177        ONLY:  bc_dirichlet_l, bc_dirichlet_n, bc_dirichlet_r, bc_dirichlet_s, &
178               cloud_droplets, constant_flux_layer, current_timestep_number,   &
179               dt_3d, dt_3d_reached, humidity,                                 &
180               dt_3d_reached_l, dt_dopts, dz, initializing_actions,            &
181               message_string, molecular_viscosity, ocean_mode,                &
182               particle_maximum_age, iran,                                     & 
183               simulated_time, topography, dopts_time_count,                   &
184               time_since_reference_point, rho_surface, u_gtrans, v_gtrans
185
186    USE cpulog,                                                                &
187        ONLY:  cpu_log, log_point, log_point_s
188       
189    USE indices,                                                               &
190        ONLY:  nx, nxl, nxlg, nxrg, nxr, ny, nyn, nys, nyng, nysg, nz, nzb,    &
191               nzb_max, nzt, wall_flags_0,nbgp, ngp_2dh_outer               
192       
193    USE kinds
194           
195    USE pegrid
196   
197    USE particle_attributes
198   
199    USE pmc_particle_interface,                                                &
200        ONLY: pmcp_c_get_particle_from_parent, pmcp_p_fill_particle_win,       &
201              pmcp_c_send_particle_to_parent, pmcp_p_empty_particle_win,       &
202              pmcp_p_delete_particles_in_fine_grid_area, pmcp_g_init,          &
203              pmcp_g_print_number_of_particles
204
205    USE pmc_interface,                                                         &
206        ONLY: nested_run
207
208    USE grid_variables,                                                        &
209        ONLY:  ddx, dx, ddy, dy
210
211    USE netcdf_interface,                                                      &
212        ONLY:  netcdf_data_format, netcdf_deflate, dopts_num, id_set_pts,      &
213               id_var_dopts, id_var_time_pts, nc_stat,                         &
214               netcdf_handle_error
215             
216    USE random_function_mod,                                                   &
217        ONLY:  random_function
218                   
219    USE statistics,                                                            &
220        ONLY:  hom                   
221
222    USE surface_mod,                                                           &
223        ONLY:  get_topography_top_index_ji, surf_def_h, surf_lsm_h, surf_usm_h,&
224               bc_h
225       
226#if defined( __parallel )  &&  !defined( __mpifh )
227    USE MPI
228#endif
229
230#if defined( __parallel )  &&  defined( __mpifh )
231    INCLUDE "mpif.h"
232#endif     
233
234#if defined( __netcdf )
235    USE NETCDF
236#endif
237
238    IMPLICIT NONE
239   
240    CHARACTER(LEN=15) ::  aero_species = 'nacl'                    !< aerosol species
241    CHARACTER(LEN=15) ::  aero_type    = 'maritime'                !< aerosol type
242    CHARACTER(LEN=15) ::  bc_par_lr    = 'cyclic'                  !< left/right boundary condition
243    CHARACTER(LEN=15) ::  bc_par_ns    = 'cyclic'                  !< north/south boundary condition
244    CHARACTER(LEN=15) ::  bc_par_b     = 'reflect'                 !< bottom boundary condition
245    CHARACTER(LEN=15) ::  bc_par_t     = 'absorb'                  !< top boundary condition
246    CHARACTER(LEN=15) ::  collision_kernel   = 'none'              !< collision kernel   
247   
248    CHARACTER(LEN=5)  ::  splitting_function = 'gamma'             !< function for calculation critical weighting factor
249    CHARACTER(LEN=5)  ::  splitting_mode     = 'const'             !< splitting mode
250
251    INTEGER(iwp) ::  deleted_particles = 0                        !< number of deleted particles per time step   
252    INTEGER(iwp) ::  i_splitting_mode                             !< dummy for splitting mode
253    INTEGER(iwp) ::  iran_part = -1234567                         !< number for random generator   
254    INTEGER(iwp) ::  max_number_particles_per_gridbox = 100       !< namelist parameter (see documentation)
255    INTEGER(iwp) ::  isf                                          !< dummy for splitting function
256    INTEGER(iwp) ::  number_particles_per_gridbox = -1            !< namelist parameter (see documentation)
257    INTEGER(iwp) ::  number_of_sublayers = 20                     !< number of sublayers for particle velocities betwenn surface and first grid level
258    INTEGER(iwp) ::  offset_ocean_nzt = 0                         !< in case of oceans runs, the vertical index calculations need an offset
259    INTEGER(iwp) ::  offset_ocean_nzt_m1 = 0                      !< in case of oceans runs, the vertical index calculations need an offset
260    INTEGER(iwp) ::  particles_per_point = 1                      !< namelist parameter (see documentation)
261    INTEGER(iwp) ::  radius_classes = 20                          !< namelist parameter (see documentation)
262   
263    INTEGER(iwp) ::  splitting_factor = 2                         !< namelist parameter (see documentation)
264    INTEGER(iwp) ::  splitting_factor_max = 5                     !< namelist parameter (see documentation)
265    INTEGER(iwp) ::  step_dealloc = 100                           !< namelist parameter (see documentation)
266    INTEGER(iwp) ::  total_number_of_particles                    !< total number of particles in the whole model domain
267    INTEGER(iwp) ::  trlp_count_sum                               !< parameter for particle exchange of PEs
268    INTEGER(iwp) ::  trlp_count_recv_sum                          !< parameter for particle exchange of PEs
269    INTEGER(iwp) ::  trrp_count_sum                               !< parameter for particle exchange of PEs
270    INTEGER(iwp) ::  trrp_count_recv_sum                          !< parameter for particle exchange of PEs
271    INTEGER(iwp) ::  trsp_count_sum                               !< parameter for particle exchange of PEs
272    INTEGER(iwp) ::  trsp_count_recv_sum                          !< parameter for particle exchange of PEs
273    INTEGER(iwp) ::  trnp_count_sum                               !< parameter for particle exchange of PEs
274    INTEGER(iwp) ::  trnp_count_recv_sum                          !< parameter for particle exchange of PEs   
275   
276    LOGICAL ::  lagrangian_particle_model = .FALSE.       !< namelist parameter (see documentation)
277    LOGICAL ::  curvature_solution_effects = .FALSE.      !< namelist parameter (see documentation)
278    LOGICAL ::  deallocate_memory = .TRUE.                !< namelist parameter (see documentation)
279    LOGICAL ::  hall_kernel = .FALSE.                     !< flag for collision kernel
280    LOGICAL ::  merging = .FALSE.                         !< namelist parameter (see documentation)
281    LOGICAL ::  random_start_position = .FALSE.           !< namelist parameter (see documentation)
282    LOGICAL ::  read_particles_from_restartfile = .TRUE.  !< namelist parameter (see documentation)
283    LOGICAL ::  seed_follows_topography = .FALSE.         !< namelist parameter (see documentation)
284    LOGICAL ::  splitting = .FALSE.                       !< namelist parameter (see documentation)
285    LOGICAL ::  use_kernel_tables = .FALSE.               !< parameter, which turns on the use of precalculated collision kernels
286    LOGICAL ::  write_particle_statistics = .FALSE.       !< namelist parameter (see documentation)
287   
288    LOGICAL, DIMENSION(max_number_of_particle_groups) ::   vertical_particle_advection = .TRUE. !< Switch for vertical particle transport
289   
290    REAL(wp) ::  aero_weight = 1.0_wp                      !< namelist parameter (see documentation)
291    REAL(wp) ::  dt_min_part = 0.0002_wp                   !< minimum particle time step when SGS velocities are used (s)
292    REAL(wp) ::  dt_prel = 9999999.9_wp                    !< namelist parameter (see documentation)
293    REAL(wp) ::  dt_write_particle_data = 9999999.9_wp     !< namelist parameter (see documentation)
294    REAL(wp) ::  end_time_prel = 9999999.9_wp              !< namelist parameter (see documentation)
295    REAL(wp) ::  initial_weighting_factor = 1.0_wp         !< namelist parameter (see documentation)
296    REAL(wp) ::  last_particle_release_time = 0.0_wp       !< last time of particle release
297    REAL(wp) ::  log_sigma(3) = 1.0_wp                     !< namelist parameter (see documentation)
298    REAL(wp) ::  na(3) = 0.0_wp                            !< namelist parameter (see documentation)
299    REAL(wp) ::  number_concentration = -1.0_wp            !< namelist parameter (see documentation)
300    REAL(wp) ::  radius_merge = 1.0E-7_wp                  !< namelist parameter (see documentation)
301    REAL(wp) ::  radius_split = 40.0E-6_wp                 !< namelist parameter (see documentation)
302    REAL(wp) ::  rm(3) = 1.0E-6_wp                         !< namelist parameter (see documentation)
303    REAL(wp) ::  sgs_wf_part                               !< parameter for sgs
304    REAL(wp) ::  time_write_particle_data = 0.0_wp         !< write particle data at current time on file
305    REAL(wp) ::  weight_factor_merge = -1.0_wp             !< namelist parameter (see documentation)
306    REAL(wp) ::  weight_factor_split = -1.0_wp             !< namelist parameter (see documentation)
307    REAL(wp) ::  z0_av_global                              !< horizontal mean value of z0
308   
309    REAL(wp) ::  rclass_lbound !<
310    REAL(wp) ::  rclass_ubound !<
311
312    REAL(wp), PARAMETER ::  c_0 = 3.0_wp         !< parameter for lagrangian timescale
313   
314    REAL(wp), DIMENSION(max_number_of_particle_groups) ::  density_ratio = 9999999.9_wp  !< namelist parameter (see documentation)
315    REAL(wp), DIMENSION(max_number_of_particle_groups) ::  pdx = 9999999.9_wp            !< namelist parameter (see documentation)
316    REAL(wp), DIMENSION(max_number_of_particle_groups) ::  pdy = 9999999.9_wp            !< namelist parameter (see documentation)
317    REAL(wp), DIMENSION(max_number_of_particle_groups) ::  pdz = 9999999.9_wp            !< namelist parameter (see documentation)
318    REAL(wp), DIMENSION(max_number_of_particle_groups) ::  psb = 9999999.9_wp            !< namelist parameter (see documentation)
319    REAL(wp), DIMENSION(max_number_of_particle_groups) ::  psl = 9999999.9_wp            !< namelist parameter (see documentation)
320    REAL(wp), DIMENSION(max_number_of_particle_groups) ::  psn = 9999999.9_wp            !< namelist parameter (see documentation)
321    REAL(wp), DIMENSION(max_number_of_particle_groups) ::  psr = 9999999.9_wp            !< namelist parameter (see documentation)
322    REAL(wp), DIMENSION(max_number_of_particle_groups) ::  pss = 9999999.9_wp            !< namelist parameter (see documentation)
323    REAL(wp), DIMENSION(max_number_of_particle_groups) ::  pst = 9999999.9_wp            !< namelist parameter (see documentation).
324    REAL(wp), DIMENSION(max_number_of_particle_groups) ::  radius = 9999999.9_wp         !< namelist parameter (see documentation)
325
326    REAL(wp), DIMENSION(:), ALLOCATABLE     ::  log_z_z0   !< Precalculate LOG(z/z0) 
327   
328    INTEGER(iwp), PARAMETER ::  NR_2_direction_move = 10000 !<
329    INTEGER(iwp)            ::  nr_move_north               !<
330    INTEGER(iwp)            ::  nr_move_south               !<
331
332    TYPE(particle_type), DIMENSION(:), ALLOCATABLE ::  move_also_north
333    TYPE(particle_type), DIMENSION(:), ALLOCATABLE ::  move_also_south
334   
335    REAL(wp) ::  epsilon       !<
336    REAL(wp) ::  urms          !<
337
338    REAL(wp), DIMENSION(:),   ALLOCATABLE ::  epsclass  !< dissipation rate class
339    REAL(wp), DIMENSION(:),   ALLOCATABLE ::  radclass  !< radius class
340    REAL(wp), DIMENSION(:),   ALLOCATABLE ::  winf      !<
341   
342    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  ec        !<
343    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  ecf       !<
344    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  gck       !<
345    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  hkernel   !<
346    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  hwratio   !<
347       
348    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  ckernel !<       
349       
350    INTEGER(iwp), PARAMETER         :: PHASE_INIT    = 1  !<
351    INTEGER(iwp), PARAMETER, PUBLIC :: PHASE_RELEASE = 2  !<
352
353    SAVE
354   
355    PRIVATE
356   
357    PUBLIC lpm_parin,     &
358           lpm_header,    &
359           lpm_init_arrays,&           
360           lpm_init,      &
361           lpm_actions,   &
362           lpm_data_output_ptseries, &           
363           lpm_rrd_local_particles, &
364           lpm_wrd_local, &
365           lpm_rrd_global, &
366           lpm_wrd_global, &
367           lpm_rrd_local, &
368           lpm_check_parameters
369
370    PUBLIC lagrangian_particle_model, & 
371           max_number_particles_per_gridbox, &
372           radius_merge, &
373           radius_split, &
374           splitting_factor, &
375           splitting_factor_max, &
376           weight_factor_merge, &
377           weight_factor_split
378           
379
380    INTERFACE lpm_check_parameters
381       MODULE PROCEDURE lpm_check_parameters
382    END INTERFACE lpm_check_parameters
383
384    INTERFACE lpm_parin
385       MODULE PROCEDURE lpm_parin
386    END INTERFACE lpm_parin   
387   
388    INTERFACE lpm_header
389       MODULE PROCEDURE lpm_header
390    END INTERFACE lpm_header
391   
392    INTERFACE lpm_init_arrays
393       MODULE PROCEDURE lpm_init_arrays
394    END INTERFACE lpm_init_arrays             
395 
396    INTERFACE lpm_init
397       MODULE PROCEDURE lpm_init
398    END INTERFACE lpm_init       
399   
400    INTERFACE lpm_actions
401       MODULE PROCEDURE lpm_actions
402    END INTERFACE lpm_actions
403   
404    INTERFACE lpm_data_output_ptseries
405       MODULE PROCEDURE lpm_data_output_ptseries
406    END INTERFACE
407       
408    INTERFACE lpm_rrd_local_particles
409       MODULE PROCEDURE lpm_rrd_local_particles
410    END INTERFACE lpm_rrd_local_particles   
411       
412    INTERFACE lpm_rrd_global
413       MODULE PROCEDURE lpm_rrd_global
414    END INTERFACE lpm_rrd_global
415   
416    INTERFACE lpm_rrd_local
417       MODULE PROCEDURE lpm_rrd_local
418    END INTERFACE lpm_rrd_local       
419
420    INTERFACE lpm_wrd_local
421       MODULE PROCEDURE lpm_wrd_local
422    END INTERFACE lpm_wrd_local
423   
424    INTERFACE lpm_wrd_global
425       MODULE PROCEDURE lpm_wrd_global
426    END INTERFACE lpm_wrd_global   
427       
428    INTERFACE lpm_advec
429       MODULE PROCEDURE lpm_advec
430    END INTERFACE lpm_advec
431   
432    INTERFACE lpm_calc_liquid_water_content
433       MODULE PROCEDURE lpm_calc_liquid_water_content
434    END INTERFACE
435   
436    INTERFACE lpm_boundary_conds
437       MODULE PROCEDURE lpm_boundary_conds
438    END INTERFACE lpm_boundary_conds
439   
440    INTERFACE lpm_droplet_condensation
441       MODULE PROCEDURE lpm_droplet_condensation
442    END INTERFACE
443   
444    INTERFACE lpm_droplet_collision
445       MODULE PROCEDURE lpm_droplet_collision
446    END INTERFACE lpm_droplet_collision
447   
448    INTERFACE lpm_init_kernels
449       MODULE PROCEDURE lpm_init_kernels
450    END INTERFACE lpm_init_kernels
451   
452    INTERFACE lpm_splitting
453       MODULE PROCEDURE lpm_splitting
454    END INTERFACE lpm_splitting
455   
456    INTERFACE lpm_merging
457       MODULE PROCEDURE lpm_merging
458    END INTERFACE lpm_merging
459   
460    INTERFACE lpm_exchange_horiz
461       MODULE PROCEDURE lpm_exchange_horiz
462    END INTERFACE lpm_exchange_horiz
463
464    INTERFACE lpm_move_particle
465       MODULE PROCEDURE lpm_move_particle
466    END INTERFACE lpm_move_particle
467
468    INTERFACE realloc_particles_array
469       MODULE PROCEDURE realloc_particles_array
470    END INTERFACE realloc_particles_array
471
472    INTERFACE dealloc_particles_array
473       MODULE PROCEDURE dealloc_particles_array
474    END INTERFACE dealloc_particles_array
475   
476    INTERFACE lpm_sort_in_subboxes
477       MODULE PROCEDURE lpm_sort_in_subboxes
478    END INTERFACE lpm_sort_in_subboxes
479
480    INTERFACE lpm_sort_timeloop_done
481       MODULE PROCEDURE lpm_sort_timeloop_done
482    END INTERFACE lpm_sort_timeloop_done
483   
484    INTERFACE lpm_pack
485       MODULE PROCEDURE lpm_pack
486    END INTERFACE lpm_pack
487   
488   
489   
490 CONTAINS
491 
492
493!------------------------------------------------------------------------------!
494! Description:
495! ------------
496!> Parin for &particle_parameters for the Lagrangian particle model
497!------------------------------------------------------------------------------!
498 SUBROUTINE lpm_parin
499 
500    CHARACTER (LEN=80) ::  line  !<
501   
502    NAMELIST /particles_par/ &
503       aero_species, &
504       aero_type, &
505       aero_weight, &       
506       alloc_factor, &
507       bc_par_b, &
508       bc_par_lr, &
509       bc_par_ns, &
510       bc_par_t, &
511       collision_kernel, &
512       curvature_solution_effects, &
513       deallocate_memory, &
514       density_ratio, &
515       dissipation_classes, &
516       dt_dopts, &
517       dt_min_part, &
518       dt_prel, &
519       dt_write_particle_data, &
520       end_time_prel, &
521       initial_weighting_factor, &
522       log_sigma, &
523       max_number_particles_per_gridbox, &
524       merging, &
525       min_nr_particle, &
526       na, &
527       number_concentration, &
528       number_of_particle_groups, &
529       number_particles_per_gridbox, &
530       particles_per_point, &
531       particle_advection_start, &
532       particle_maximum_age, &
533       pdx, &
534       pdy, &
535       pdz, &
536       psb, &
537       psl, &
538       psn, &
539       psr, &
540       pss, &
541       pst, &
542       radius, &
543       radius_classes, &
544       radius_merge, &
545       radius_split, &
546       random_start_position, &
547       read_particles_from_restartfile, &
548       rm, &
549       seed_follows_topography, &
550       splitting, &
551       splitting_factor, &
552       splitting_factor_max, &
553       splitting_function, &
554       splitting_mode, &
555       step_dealloc, &
556       use_sgs_for_particles, &
557       vertical_particle_advection, &
558       weight_factor_merge, &
559       weight_factor_split, &
560       write_particle_statistics
561                                 
562    NAMELIST /particle_parameters/ &
563       aero_species, &
564       aero_type, &
565       aero_weight, &       
566       alloc_factor, &
567       bc_par_b, &
568       bc_par_lr, &
569       bc_par_ns, &
570       bc_par_t, &
571       collision_kernel, &
572       curvature_solution_effects, &
573       deallocate_memory, &
574       density_ratio, &
575       dissipation_classes, &
576       dt_dopts, &
577       dt_min_part, &
578       dt_prel, &
579       dt_write_particle_data, &
580       end_time_prel, &
581       initial_weighting_factor, &
582       log_sigma, &
583       max_number_particles_per_gridbox, &
584       merging, &
585       min_nr_particle, &
586       na, &
587       number_concentration, &
588       number_of_particle_groups, &
589       number_particles_per_gridbox, &
590       particles_per_point, &
591       particle_advection_start, &
592       particle_maximum_age, &
593       pdx, &
594       pdy, &
595       pdz, &
596       psb, &
597       psl, &
598       psn, &
599       psr, &
600       pss, &
601       pst, &
602       radius, &
603       radius_classes, &
604       radius_merge, &
605       radius_split, &
606       random_start_position, &
607       read_particles_from_restartfile, &
608       rm, &
609       seed_follows_topography, &
610       splitting, &
611       splitting_factor, &
612       splitting_factor_max, &
613       splitting_function, &
614       splitting_mode, &
615       step_dealloc, &
616       use_sgs_for_particles, &
617       vertical_particle_advection, &
618       weight_factor_merge, &
619       weight_factor_split, &
620       write_particle_statistics
621       
622!
623!-- Position the namelist-file at the beginning (it was already opened in
624!-- parin), search for the namelist-group of the package and position the
625!-- file at this line. Do the same for each optionally used package.
626    line = ' '
627   
628!
629!-- Try to find particles package
630    REWIND ( 11 )
631    line = ' '
632    DO   WHILE ( INDEX( line, '&particle_parameters' ) == 0 )
633       READ ( 11, '(A)', END=12 )  line
634    ENDDO
635    BACKSPACE ( 11 )
636!
637!-- Read user-defined namelist
638    READ ( 11, particle_parameters, ERR = 10 )
639!
640!-- Set flag that indicates that particles are switched on
641    particle_advection = .TRUE.
642   
643    GOTO 14
644
64510  BACKSPACE( 11 )
646    READ( 11 , '(A)') line
647    CALL parin_fail_message( 'particle_parameters', line )
648!
649!-- Try to find particles package (old namelist)
65012  REWIND ( 11 )
651    line = ' '
652    DO WHILE ( INDEX( line, '&particles_par' ) == 0 )
653       READ ( 11, '(A)', END=14 )  line
654    ENDDO
655    BACKSPACE ( 11 )
656!
657!-- Read user-defined namelist
658    READ ( 11, particles_par, ERR = 13, END = 14 )
659   
660   
661    message_string = 'namelist particles_par is deprecated and will be ' //    &
662                     'removed in near future. Please use namelist ' //         &
663                     'particle_parameters instead'
664    CALL message( 'package_parin', 'PA0487', 0, 1, 0, 6, 0 )
665
666!
667!-- Set flag that indicates that particles are switched on
668    particle_advection = .TRUE.
669
670    GOTO 14
671
67213    BACKSPACE( 11 )
673       READ( 11 , '(A)') line
674       CALL parin_fail_message( 'particles_par', line )
675
67614 CONTINUE
677   
678 END SUBROUTINE lpm_parin
679 
680!------------------------------------------------------------------------------!
681! Description:
682! ------------
683!> Writes used particle attributes in header file.
684!------------------------------------------------------------------------------!
685 SUBROUTINE lpm_header ( io )
686
687    CHARACTER (LEN=40) ::  output_format       !< netcdf format
688 
689    INTEGER(iwp) ::  i               !<
690    INTEGER(iwp), INTENT(IN) ::  io  !< Unit of the output file
691
692 
693     IF ( humidity  .AND.  cloud_droplets )  THEN
694       WRITE ( io, 433 )
695       IF ( curvature_solution_effects )  WRITE ( io, 434 )
696       IF ( collision_kernel /= 'none' )  THEN
697          WRITE ( io, 435 )  TRIM( collision_kernel )
698          IF ( collision_kernel(6:9) == 'fast' )  THEN
699             WRITE ( io, 436 )  radius_classes, dissipation_classes
700          ENDIF
701       ELSE
702          WRITE ( io, 437 )
703       ENDIF
704    ENDIF
705 
706    IF ( particle_advection )  THEN
707!
708!--    Particle attributes
709       WRITE ( io, 480 )  particle_advection_start, dt_prel, bc_par_lr, &
710                          bc_par_ns, bc_par_b, bc_par_t, particle_maximum_age, &
711                          end_time_prel
712       IF ( use_sgs_for_particles )  WRITE ( io, 488 )  dt_min_part
713       IF ( random_start_position )  WRITE ( io, 481 )
714       IF ( seed_follows_topography )  WRITE ( io, 496 )
715       IF ( particles_per_point > 1 )  WRITE ( io, 489 )  particles_per_point
716       WRITE ( io, 495 )  total_number_of_particles
717       IF ( dt_write_particle_data /= 9999999.9_wp )  THEN
718          WRITE ( io, 485 )  dt_write_particle_data
719          IF ( netcdf_data_format > 1 )  THEN
720             output_format = 'netcdf (64 bit offset) and binary'
721          ELSE
722             output_format = 'netcdf and binary'
723          ENDIF
724          IF ( netcdf_deflate == 0 )  THEN
725             WRITE ( io, 344 )  output_format
726          ELSE
727             WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
728          ENDIF
729       ENDIF
730       IF ( dt_dopts /= 9999999.9_wp )  WRITE ( io, 494 )  dt_dopts
731       IF ( write_particle_statistics )  WRITE ( io, 486 )
732
733       WRITE ( io, 487 )  number_of_particle_groups
734
735       DO  i = 1, number_of_particle_groups
736          IF ( i == 1  .AND.  density_ratio(i) == 9999999.9_wp )  THEN
737             WRITE ( io, 490 )  i, 0.0_wp
738             WRITE ( io, 492 )
739          ELSE
740             WRITE ( io, 490 )  i, radius(i)
741             IF ( density_ratio(i) /= 0.0_wp )  THEN
742                WRITE ( io, 491 )  density_ratio(i)
743             ELSE
744                WRITE ( io, 492 )
745             ENDIF
746          ENDIF
747          WRITE ( io, 493 )  psl(i), psr(i), pss(i), psn(i), psb(i), pst(i), &
748                             pdx(i), pdy(i), pdz(i)
749          IF ( .NOT. vertical_particle_advection(i) )  WRITE ( io, 482 )
750       ENDDO
751
752    ENDIF
753   
754344 FORMAT ('       Output format: ',A/)
755354 FORMAT ('       Output format: ',A, '   compressed with level: ',I1/)
756
757433 FORMAT ('    Cloud droplets treated explicitly using the Lagrangian part', &
758                 'icle model')
759434 FORMAT ('    Curvature and solution effecs are considered for growth of', &
760                 ' droplets < 1.0E-6 m')
761435 FORMAT ('    Droplet collision is handled by ',A,'-kernel')
762436 FORMAT ('       Fast kernel with fixed radius- and dissipation classes ', &
763                    'are used'/ &
764            '          number of radius classes:       ',I3,'    interval ', &
765                       '[1.0E-6,2.0E-4] m'/ &
766            '          number of dissipation classes:   ',I2,'    interval ', &
767                       '[0,1000] cm**2/s**3')
768437 FORMAT ('    Droplet collision is switched off')
769
770480 FORMAT ('    Particles:'/ &
771            '    ---------'// &
772            '       Particle advection is active (switched on at t = ', F7.1, &
773                    ' s)'/ &
774            '       Start of new particle generations every  ',F6.1,' s'/ &
775            '       Boundary conditions: left/right: ', A, ' north/south: ', A/&
776            '                            bottom:     ', A, ' top:         ', A/&
777            '       Maximum particle age:                 ',F9.1,' s'/ &
778            '       Advection stopped at t = ',F9.1,' s'/)
779481 FORMAT ('       Particles have random start positions'/)
780482 FORMAT ('          Particles are advected only horizontally'/)
781485 FORMAT ('       Particle data are written on file every ', F9.1, ' s')
782486 FORMAT ('       Particle statistics are written on file'/)
783487 FORMAT ('       Number of particle groups: ',I2/)
784488 FORMAT ('       SGS velocity components are used for particle advection'/ &
785            '          minimum timestep for advection:', F8.5/)
786489 FORMAT ('       Number of particles simultaneously released at each ', &
787                    'point: ', I5/)
788490 FORMAT ('       Particle group ',I2,':'/ &
789            '          Particle radius: ',E10.3, 'm')
790491 FORMAT ('          Particle inertia is activated'/ &
791            '             density_ratio (rho_fluid/rho_particle) =',F6.3/)
792492 FORMAT ('          Particles are advected only passively (no inertia)'/)
793493 FORMAT ('          Boundaries of particle source: x:',F8.1,' - ',F8.1,' m'/&
794            '                                         y:',F8.1,' - ',F8.1,' m'/&
795            '                                         z:',F8.1,' - ',F8.1,' m'/&
796            '          Particle distances:  dx = ',F8.1,' m  dy = ',F8.1, &
797                       ' m  dz = ',F8.1,' m'/)
798494 FORMAT ('       Output of particle time series in NetCDF format every ', &
799                    F8.2,' s'/)
800495 FORMAT ('       Number of particles in total domain: ',I10/)
801496 FORMAT ('       Initial vertical particle positions are interpreted ', &
802                    'as relative to the given topography')
803   
804 END SUBROUTINE lpm_header
805 
806!------------------------------------------------------------------------------!
807! Description:
808! ------------
809!> Writes used particle attributes in header file.
810!------------------------------------------------------------------------------! 
811 SUBROUTINE lpm_check_parameters
812 
813!
814!-- Collision kernels:
815    SELECT CASE ( TRIM( collision_kernel ) )
816
817       CASE ( 'hall', 'hall_fast' )
818          hall_kernel = .TRUE.
819
820       CASE ( 'wang', 'wang_fast' )
821          wang_kernel = .TRUE.
822
823       CASE ( 'none' )
824
825
826       CASE DEFAULT
827          message_string = 'unknown collision kernel: collision_kernel = "' // &
828                           TRIM( collision_kernel ) // '"'
829          CALL message( 'check_parameters', 'PA0350', 1, 2, 0, 6, 0 )
830
831    END SELECT
832    IF ( collision_kernel(6:9) == 'fast' )  use_kernel_tables = .TRUE.
833
834 END SUBROUTINE
835 
836!------------------------------------------------------------------------------!
837! Description:
838! ------------
839!> Initialize arrays for lpm
840!------------------------------------------------------------------------------!   
841 SUBROUTINE lpm_init_arrays
842 
843    IF ( cloud_droplets )  THEN
844!
845!--    Liquid water content, change in liquid water content
846       ALLOCATE ( ql_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                      &
847                  ql_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
848!
849!--    Real volume of particles (with weighting), volume of particles
850       ALLOCATE ( ql_v(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                      &
851                     ql_vp(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
852    ENDIF
853   
854!
855!--    Initial assignment of the pointers   
856    IF ( cloud_droplets )  THEN
857       ql   => ql_1
858       ql_c => ql_2
859    ENDIF
860   
861 END SUBROUTINE lpm_init_arrays
862 
863!------------------------------------------------------------------------------!
864! Description:
865! ------------
866!> Initialize Lagrangian particle model
867!------------------------------------------------------------------------------!
868 SUBROUTINE lpm_init
869
870    INTEGER(iwp) ::  i                           !<
871    INTEGER(iwp) ::  j                           !<
872    INTEGER(iwp) ::  k                           !<
873
874    REAL(wp) ::  div                             !<
875    REAL(wp) ::  height_int                      !<
876    REAL(wp) ::  height_p                        !<
877    REAL(wp) ::  z_p                             !<
878    REAL(wp) ::  z0_av_local                     !<
879
880!
881!-- In case of oceans runs, the vertical index calculations need an offset,
882!-- because otherwise the k indices will become negative
883    IF ( ocean_mode )  THEN
884       offset_ocean_nzt    = nzt
885       offset_ocean_nzt_m1 = nzt - 1
886    ENDIF
887
888!
889!-- Define block offsets for dividing a gridcell in 8 sub cells
890!-- See documentation for List of subgrid boxes
891!-- See pack_and_sort in lpm_pack_arrays.f90 for assignment of the subgrid boxes
892    block_offset(0) = block_offset_def ( 0, 0, 0)
893    block_offset(1) = block_offset_def ( 0, 0,-1)
894    block_offset(2) = block_offset_def ( 0,-1, 0)
895    block_offset(3) = block_offset_def ( 0,-1,-1)
896    block_offset(4) = block_offset_def (-1, 0, 0)
897    block_offset(5) = block_offset_def (-1, 0,-1)
898    block_offset(6) = block_offset_def (-1,-1, 0)
899    block_offset(7) = block_offset_def (-1,-1,-1)
900!
901!-- Check the number of particle groups.
902    IF ( number_of_particle_groups > max_number_of_particle_groups )  THEN
903       WRITE( message_string, * ) 'max_number_of_particle_groups =',           &
904                                  max_number_of_particle_groups ,              &
905                                  '&number_of_particle_groups reset to ',      &
906                                  max_number_of_particle_groups
907       CALL message( 'lpm_init', 'PA0213', 0, 1, 0, 6, 0 )
908       number_of_particle_groups = max_number_of_particle_groups
909    ENDIF
910!
911!-- Check if downward-facing walls exist. This case, reflection boundary
912!-- conditions (as well as subgrid-scale velocities) may do not work
913!-- propably (not realized so far).
914    IF ( surf_def_h(1)%ns >= 1 )  THEN
915       WRITE( message_string, * ) 'Overhanging topography do not work '//      &
916                                  'with particles'
917       CALL message( 'lpm_init', 'PA0212', 0, 1, 0, 6, 0 )
918
919    ENDIF
920
921!
922!-- Set default start positions, if necessary
923    IF ( psl(1) == 9999999.9_wp )  psl(1) = 0.0_wp
924    IF ( psr(1) == 9999999.9_wp )  psr(1) = ( nx +1 ) * dx
925    IF ( pss(1) == 9999999.9_wp )  pss(1) = 0.0_wp
926    IF ( psn(1) == 9999999.9_wp )  psn(1) = ( ny +1 ) * dy
927    IF ( psb(1) == 9999999.9_wp )  psb(1) = zu(nz/2)
928    IF ( pst(1) == 9999999.9_wp )  pst(1) = psb(1)
929
930    IF ( pdx(1) == 9999999.9_wp  .OR.  pdx(1) == 0.0_wp )  pdx(1) = dx
931    IF ( pdy(1) == 9999999.9_wp  .OR.  pdy(1) == 0.0_wp )  pdy(1) = dy
932    IF ( pdz(1) == 9999999.9_wp  .OR.  pdz(1) == 0.0_wp )  pdz(1) = zu(2) - zu(1)
933
934!
935!-- If number_particles_per_gridbox is set, the parametres pdx, pdy and pdz are
936!-- calculated diagnostically. Therfore an isotropic distribution is prescribed.
937    IF ( number_particles_per_gridbox /= -1 .AND.   &
938         number_particles_per_gridbox >= 1 )    THEN
939       pdx(1) = (( dx * dy * ( zu(2) - zu(1) ) ) /  &
940             REAL(number_particles_per_gridbox))**0.3333333_wp
941!
942!--    Ensure a smooth value (two significant digits) of distance between
943!--    particles (pdx, pdy, pdz).
944       div = 1000.0_wp
945       DO  WHILE ( pdx(1) < div )
946          div = div / 10.0_wp
947       ENDDO
948       pdx(1) = NINT( pdx(1) * 100.0_wp / div ) * div / 100.0_wp
949       pdy(1) = pdx(1)
950       pdz(1) = pdx(1)
951
952    ENDIF
953
954    DO  j = 2, number_of_particle_groups
955       IF ( psl(j) == 9999999.9_wp )  psl(j) = psl(j-1)
956       IF ( psr(j) == 9999999.9_wp )  psr(j) = psr(j-1)
957       IF ( pss(j) == 9999999.9_wp )  pss(j) = pss(j-1)
958       IF ( psn(j) == 9999999.9_wp )  psn(j) = psn(j-1)
959       IF ( psb(j) == 9999999.9_wp )  psb(j) = psb(j-1)
960       IF ( pst(j) == 9999999.9_wp )  pst(j) = pst(j-1)
961       IF ( pdx(j) == 9999999.9_wp  .OR.  pdx(j) == 0.0_wp )  pdx(j) = pdx(j-1)
962       IF ( pdy(j) == 9999999.9_wp  .OR.  pdy(j) == 0.0_wp )  pdy(j) = pdy(j-1)
963       IF ( pdz(j) == 9999999.9_wp  .OR.  pdz(j) == 0.0_wp )  pdz(j) = pdz(j-1)
964    ENDDO
965
966!
967!-- Allocate arrays required for calculating particle SGS velocities.
968!-- Initialize prefactor required for stoachastic Weil equation.
969    IF ( use_sgs_for_particles  .AND.  .NOT. cloud_droplets )  THEN
970       ALLOCATE( de_dx(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
971                 de_dy(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
972                 de_dz(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
973
974       de_dx = 0.0_wp
975       de_dy = 0.0_wp
976       de_dz = 0.0_wp             
977                 
978       sgs_wf_part = 1.0_wp / 3.0_wp
979    ENDIF
980
981!
982!-- Allocate array required for logarithmic vertical interpolation of
983!-- horizontal particle velocities between the surface and the first vertical
984!-- grid level. In order to avoid repeated CPU cost-intensive CALLS of
985!-- intrinsic FORTRAN procedure LOG(z/z0), LOG(z/z0) is precalculated for
986!-- several heights. Splitting into 20 sublayers turned out to be sufficient.
987!-- To obtain exact height levels of particles, linear interpolation is applied
988!-- (see lpm_advec.f90).
989    IF ( constant_flux_layer )  THEN
990
991       ALLOCATE ( log_z_z0(0:number_of_sublayers) )
992       z_p = zu(nzb+1) - zw(nzb)
993
994!
995!--    Calculate horizontal mean value of z0 used for logartihmic
996!--    interpolation. Note: this is not exact for heterogeneous z0.
997!--    However, sensitivity studies showed that the effect is
998!--    negligible.
999       z0_av_local  = SUM( surf_def_h(0)%z0 ) + SUM( surf_lsm_h%z0 ) +         &
1000                      SUM( surf_usm_h%z0 )
1001       z0_av_global = 0.0_wp
1002
1003#if defined( __parallel )
1004       CALL MPI_ALLREDUCE(z0_av_local, z0_av_global, 1, MPI_REAL, MPI_SUM, &
1005                          comm2d, ierr )
1006#else
1007       z0_av_global = z0_av_local
1008#endif
1009
1010       z0_av_global = z0_av_global  / ( ( ny + 1 ) * ( nx + 1 ) )
1011!
1012!--    Horizontal wind speed is zero below and at z0
1013       log_z_z0(0) = 0.0_wp
1014!
1015!--    Calculate vertical depth of the sublayers
1016       height_int  = ( z_p - z0_av_global ) / REAL( number_of_sublayers, KIND=wp )
1017!
1018!--    Precalculate LOG(z/z0)
1019       height_p    = z0_av_global
1020       DO  k = 1, number_of_sublayers
1021
1022          height_p    = height_p + height_int
1023          log_z_z0(k) = LOG( height_p / z0_av_global )
1024
1025       ENDDO
1026
1027    ENDIF
1028
1029!
1030!-- Check boundary condition and set internal variables
1031    SELECT CASE ( bc_par_b )
1032
1033       CASE ( 'absorb' )
1034          ibc_par_b = 1
1035
1036       CASE ( 'reflect' )
1037          ibc_par_b = 2
1038
1039       CASE DEFAULT
1040          WRITE( message_string, * )  'unknown boundary condition ',           &
1041                                       'bc_par_b = "', TRIM( bc_par_b ), '"'
1042          CALL message( 'lpm_init', 'PA0217', 1, 2, 0, 6, 0 )
1043
1044    END SELECT
1045    SELECT CASE ( bc_par_t )
1046
1047       CASE ( 'absorb' )
1048          ibc_par_t = 1
1049
1050       CASE ( 'reflect' )
1051          ibc_par_t = 2
1052         
1053       CASE ( 'nested' )
1054          ibc_par_t = 3
1055
1056       CASE DEFAULT
1057          WRITE( message_string, * ) 'unknown boundary condition ',            &
1058                                     'bc_par_t = "', TRIM( bc_par_t ), '"'
1059          CALL message( 'lpm_init', 'PA0218', 1, 2, 0, 6, 0 )
1060
1061    END SELECT
1062    SELECT CASE ( bc_par_lr )
1063
1064       CASE ( 'cyclic' )
1065          ibc_par_lr = 0
1066
1067       CASE ( 'absorb' )
1068          ibc_par_lr = 1
1069
1070       CASE ( 'reflect' )
1071          ibc_par_lr = 2
1072         
1073       CASE ( 'nested' )
1074          ibc_par_lr = 3
1075
1076       CASE DEFAULT
1077          WRITE( message_string, * ) 'unknown boundary condition ',   &
1078                                     'bc_par_lr = "', TRIM( bc_par_lr ), '"'
1079          CALL message( 'lpm_init', 'PA0219', 1, 2, 0, 6, 0 )
1080
1081    END SELECT
1082    SELECT CASE ( bc_par_ns )
1083
1084       CASE ( 'cyclic' )
1085          ibc_par_ns = 0
1086
1087       CASE ( 'absorb' )
1088          ibc_par_ns = 1
1089
1090       CASE ( 'reflect' )
1091          ibc_par_ns = 2
1092         
1093       CASE ( 'nested' )
1094          ibc_par_ns = 3
1095
1096       CASE DEFAULT
1097          WRITE( message_string, * ) 'unknown boundary condition ',   &
1098                                     'bc_par_ns = "', TRIM( bc_par_ns ), '"'
1099          CALL message( 'lpm_init', 'PA0220', 1, 2, 0, 6, 0 )
1100
1101    END SELECT
1102    SELECT CASE ( splitting_mode )
1103
1104       CASE ( 'const' )
1105          i_splitting_mode = 1
1106
1107       CASE ( 'cl_av' )
1108          i_splitting_mode = 2
1109
1110       CASE ( 'gb_av' )
1111          i_splitting_mode = 3
1112
1113       CASE DEFAULT
1114          WRITE( message_string, * )  'unknown splitting_mode = "',            &
1115                                      TRIM( splitting_mode ), '"'
1116          CALL message( 'lpm_init', 'PA0146', 1, 2, 0, 6, 0 )
1117
1118    END SELECT
1119    SELECT CASE ( splitting_function )
1120
1121       CASE ( 'gamma' )
1122          isf = 1
1123
1124       CASE ( 'log' )
1125          isf = 2
1126
1127       CASE ( 'exp' )
1128          isf = 3
1129
1130       CASE DEFAULT
1131          WRITE( message_string, * )  'unknown splitting function = "',        &
1132                                       TRIM( splitting_function ), '"'
1133          CALL message( 'lpm_init', 'PA0147', 1, 2, 0, 6, 0 )
1134
1135    END SELECT
1136!
1137!-- Initialize collision kernels
1138    IF ( collision_kernel /= 'none' )  CALL lpm_init_kernels
1139!
1140!-- For the first model run of a possible job chain initialize the
1141!-- particles, otherwise read the particle data from restart file.
1142    IF ( TRIM( initializing_actions ) == 'read_restart_data'  &
1143         .AND.  read_particles_from_restartfile )  THEN
1144       CALL lpm_rrd_local_particles
1145    ELSE
1146!
1147!--    Allocate particle arrays and set attributes of the initial set of
1148!--    particles, which can be also periodically released at later times.
1149       ALLOCATE( prt_count(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
1150                 grid_particles(nzb+1:nzt,nys:nyn,nxl:nxr) )
1151
1152       number_of_particles = 0
1153       prt_count           = 0
1154!
1155!--    initialize counter for particle IDs
1156       grid_particles%id_counter = 1
1157!
1158!--    Initialize all particles with dummy values (otherwise errors may
1159!--    occur within restart runs). The reason for this is still not clear
1160!--    and may be presumably caused by errors in the respective user-interface.
1161       zero_particle = particle_type( 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp,  &
1162                                      0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp,  &
1163                                      0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp,  &
1164                                      0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp,  &
1165                                      0, 0, 0_idp, .FALSE., -1 )
1166
1167       particle_groups = particle_groups_type( 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp )
1168!
1169!--    Set values for the density ratio and radius for all particle
1170!--    groups, if necessary
1171       IF ( density_ratio(1) == 9999999.9_wp )  density_ratio(1) = 0.0_wp
1172       IF ( radius(1)        == 9999999.9_wp )  radius(1) = 0.0_wp
1173       DO  i = 2, number_of_particle_groups
1174          IF ( density_ratio(i) == 9999999.9_wp )  THEN
1175             density_ratio(i) = density_ratio(i-1)
1176          ENDIF
1177          IF ( radius(i) == 9999999.9_wp )  radius(i) = radius(i-1)
1178       ENDDO
1179
1180       DO  i = 1, number_of_particle_groups
1181          IF ( density_ratio(i) /= 0.0_wp  .AND.  radius(i) == 0 )  THEN
1182             WRITE( message_string, * ) 'particle group #', i, ' has a',       &
1183                                        'density ratio /= 0 but radius = 0'
1184             CALL message( 'lpm_init', 'PA0215', 1, 2, 0, 6, 0 )
1185          ENDIF
1186          particle_groups(i)%density_ratio = density_ratio(i)
1187          particle_groups(i)%radius        = radius(i)
1188       ENDDO
1189!
1190!--    Set a seed value for the random number generator to be exclusively
1191!--    used for the particle code. The generated random numbers should be
1192!--    different on the different PEs.
1193       iran_part = iran_part + myid
1194!
1195!--    Create the particle set, and set the initial particles
1196       CALL lpm_create_particle( phase_init )
1197       last_particle_release_time = particle_advection_start
1198!
1199!--    User modification of initial particles
1200       CALL user_lpm_init
1201!
1202!--    Open file for statistical informations about particle conditions
1203       IF ( write_particle_statistics )  THEN
1204          CALL check_open( 80 )
1205          WRITE ( 80, 8000 )  current_timestep_number, simulated_time,         &
1206                              number_of_particles
1207          CALL close_file( 80 )
1208       ENDIF
1209
1210    ENDIF
1211
1212    IF ( nested_run )  CALL pmcp_g_init
1213!
1214!-- To avoid programm abort, assign particles array to the local version of
1215!-- first grid cell
1216    number_of_particles = prt_count(nzb+1,nys,nxl)
1217    particles => grid_particles(nzb+1,nys,nxl)%particles(1:number_of_particles)
1218!
1219!-- Formats
12208000 FORMAT (I6,1X,F7.2,4X,I10,71X,I10)
1221
1222 END SUBROUTINE lpm_init
1223 
1224!------------------------------------------------------------------------------!
1225! Description:
1226! ------------
1227!> Create Lagrangian particles
1228!------------------------------------------------------------------------------! 
1229 SUBROUTINE lpm_create_particle (phase)
1230
1231    INTEGER(iwp)               ::  alloc_size  !< relative increase of allocated memory for particles
1232    INTEGER(iwp)               ::  i           !< loop variable ( particle groups )
1233    INTEGER(iwp)               ::  ip          !< index variable along x
1234    INTEGER(iwp)               ::  j           !< loop variable ( particles per point )
1235    INTEGER(iwp)               ::  jp          !< index variable along y
1236    INTEGER(iwp)               ::  k           !< index variable along z
1237    INTEGER(iwp)               ::  k_surf      !< index of surface grid point
1238    INTEGER(iwp)               ::  kp          !< index variable along z
1239    INTEGER(iwp)               ::  loop_stride !< loop variable for initialization
1240    INTEGER(iwp)               ::  n           !< loop variable ( number of particles )
1241    INTEGER(iwp)               ::  new_size    !< new size of allocated memory for particles
1242
1243    INTEGER(iwp), INTENT(IN)   ::  phase       !< mode of inititialization
1244
1245    INTEGER(iwp), DIMENSION(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ::  local_count !< start address of new particle
1246    INTEGER(iwp), DIMENSION(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ::  local_start !< start address of new particle
1247
1248    LOGICAL                    ::  first_stride !< flag for initialization
1249
1250    REAL(wp)                   ::  pos_x      !< increment for particle position in x
1251    REAL(wp)                   ::  pos_y      !< increment for particle position in y
1252    REAL(wp)                   ::  pos_z      !< increment for particle position in z
1253    REAL(wp)                   ::  rand_contr !< dummy argument for random position
1254
1255    TYPE(particle_type),TARGET ::  tmp_particle !< temporary particle used for initialization
1256
1257!
1258!-- Calculate particle positions and store particle attributes, if
1259!-- particle is situated on this PE
1260    DO  loop_stride = 1, 2
1261       first_stride = (loop_stride == 1)
1262       IF ( first_stride )   THEN
1263          local_count = 0           ! count number of particles
1264       ELSE
1265          local_count = prt_count   ! Start address of new particles
1266       ENDIF
1267
1268!
1269!--    Calculate initial_weighting_factor diagnostically
1270       IF ( number_concentration /= -1.0_wp .AND. number_concentration > 0.0_wp ) THEN
1271          initial_weighting_factor =  number_concentration  *                        &
1272                                      pdx(1) * pdy(1) * pdz(1)
1273       END IF
1274
1275       n = 0
1276       DO  i = 1, number_of_particle_groups
1277          pos_z = psb(i)
1278          DO WHILE ( pos_z <= pst(i) )
1279             IF ( pos_z >= zw(0) .AND.  pos_z < zw(nzt) )  THEN
1280                pos_y = pss(i)
1281                DO WHILE ( pos_y <= psn(i) )
1282                   IF ( pos_y >= nys * dy  .AND.                  &
1283                        pos_y <  ( nyn + 1 ) * dy  ) THEN
1284                      pos_x = psl(i)
1285               xloop: DO WHILE ( pos_x <= psr(i) )
1286                         IF ( pos_x >= nxl * dx  .AND.            &
1287                              pos_x <  ( nxr + 1) * dx ) THEN
1288                            DO  j = 1, particles_per_point
1289                               n = n + 1
1290                               tmp_particle%x             = pos_x
1291                               tmp_particle%y             = pos_y
1292                               tmp_particle%z             = pos_z
1293                               tmp_particle%age           = 0.0_wp
1294                               tmp_particle%age_m         = 0.0_wp
1295                               tmp_particle%dt_sum        = 0.0_wp
1296                               tmp_particle%e_m           = 0.0_wp
1297                               tmp_particle%rvar1         = 0.0_wp
1298                               tmp_particle%rvar2         = 0.0_wp
1299                               tmp_particle%rvar3         = 0.0_wp
1300                               tmp_particle%speed_x       = 0.0_wp
1301                               tmp_particle%speed_y       = 0.0_wp
1302                               tmp_particle%speed_z       = 0.0_wp
1303                               tmp_particle%origin_x      = pos_x
1304                               tmp_particle%origin_y      = pos_y
1305                               tmp_particle%origin_z      = pos_z
1306                               IF ( curvature_solution_effects )  THEN
1307                                  tmp_particle%aux1      = 0.0_wp    ! dry aerosol radius
1308                                  tmp_particle%aux2      = dt_3d     ! last Rosenbrock timestep
1309                               ELSE
1310                                  tmp_particle%aux1      = 0.0_wp    ! free to use
1311                                  tmp_particle%aux2      = 0.0_wp    ! free to use
1312                               ENDIF
1313                               tmp_particle%radius        = particle_groups(i)%radius
1314                               tmp_particle%weight_factor = initial_weighting_factor
1315                               tmp_particle%class         = 1
1316                               tmp_particle%group         = i
1317                               tmp_particle%id            = 0_idp
1318                               tmp_particle%particle_mask = .TRUE.
1319                               tmp_particle%block_nr      = -1
1320!
1321!--                            Determine the grid indices of the particle position
1322                               ip = INT( tmp_particle%x * ddx )
1323                               jp = INT( tmp_particle%y * ddy )
1324                               kp = INT( tmp_particle%z / dz(1) + 1 + offset_ocean_nzt )
1325                               DO WHILE( zw(kp) < tmp_particle%z ) 
1326                                  kp = kp + 1
1327                               ENDDO
1328                               DO WHILE( zw(kp-1) > tmp_particle%z )
1329                                  kp = kp - 1
1330                               ENDDO 
1331!
1332!--                            Determine surface level. Therefore, check for
1333!--                            upward-facing wall on w-grid.
1334                               k_surf = get_topography_top_index_ji( jp, ip, 'w' )
1335                               IF ( seed_follows_topography )  THEN
1336!
1337!--                               Particle height is given relative to topography
1338                                  kp = kp + k_surf
1339                                  tmp_particle%z = tmp_particle%z + zw(k_surf)
1340!--                               Skip particle release if particle position is
1341!--                               above model top, or within topography in case
1342!--                               of overhanging structures.
1343                                  IF ( kp > nzt  .OR.                          &
1344                                 .NOT. BTEST( wall_flags_0(kp,jp,ip), 0 ) )  THEN
1345                                     pos_x = pos_x + pdx(i)
1346                                     CYCLE xloop
1347                                  ENDIF
1348!
1349!--                            Skip particle release if particle position is
1350!--                            below surface, or within topography in case
1351!--                            of overhanging structures.
1352                               ELSEIF ( .NOT. seed_follows_topography .AND.    &
1353                                         tmp_particle%z <= zw(k_surf)  .OR.    &
1354                                        .NOT. BTEST( wall_flags_0(kp,jp,ip), 0 ) )&
1355                               THEN
1356                                  pos_x = pos_x + pdx(i)
1357                                  CYCLE xloop
1358                               ENDIF
1359
1360                               local_count(kp,jp,ip) = local_count(kp,jp,ip) + 1
1361
1362                               IF ( .NOT. first_stride )  THEN
1363                                  IF ( ip < nxl  .OR.  jp < nys  .OR.  kp < nzb+1 )  THEN
1364                                     write(6,*) 'xl ',ip,jp,kp,nxl,nys,nzb+1
1365                                  ENDIF
1366                                  IF ( ip > nxr  .OR.  jp > nyn  .OR.  kp > nzt )  THEN
1367                                     write(6,*) 'xu ',ip,jp,kp,nxr,nyn,nzt
1368                                  ENDIF
1369                                  grid_particles(kp,jp,ip)%particles(local_count(kp,jp,ip)) = tmp_particle
1370                               ENDIF
1371                            ENDDO
1372                         ENDIF
1373                         pos_x = pos_x + pdx(i)
1374                      ENDDO xloop
1375                   ENDIF
1376                   pos_y = pos_y + pdy(i)
1377                ENDDO
1378             ENDIF
1379
1380             pos_z = pos_z + pdz(i)
1381          ENDDO
1382       ENDDO
1383
1384       IF ( first_stride )  THEN
1385          DO  ip = nxl, nxr
1386             DO  jp = nys, nyn
1387                DO  kp = nzb+1, nzt
1388                   IF ( phase == PHASE_INIT )  THEN
1389                      IF ( local_count(kp,jp,ip) > 0 )  THEN
1390                         alloc_size = MAX( INT( local_count(kp,jp,ip) *        &
1391                            ( 1.0_wp + alloc_factor / 100.0_wp ) ),            &
1392                            min_nr_particle )
1393                      ELSE
1394                         alloc_size = min_nr_particle
1395                      ENDIF
1396                      ALLOCATE(grid_particles(kp,jp,ip)%particles(1:alloc_size))
1397                      DO  n = 1, alloc_size
1398                         grid_particles(kp,jp,ip)%particles(n) = zero_particle
1399                      ENDDO
1400                   ELSEIF ( phase == PHASE_RELEASE )  THEN
1401                      IF ( local_count(kp,jp,ip) > 0 )  THEN
1402                         new_size   = local_count(kp,jp,ip) + prt_count(kp,jp,ip)
1403                         alloc_size = MAX( INT( new_size * ( 1.0_wp +          &
1404                            alloc_factor / 100.0_wp ) ), min_nr_particle )
1405                         IF( alloc_size > SIZE( grid_particles(kp,jp,ip)%particles) )  THEN
1406                            CALL realloc_particles_array(ip,jp,kp,alloc_size)
1407                         ENDIF
1408                      ENDIF
1409                   ENDIF
1410                ENDDO
1411             ENDDO
1412          ENDDO
1413       ENDIF
1414
1415    ENDDO
1416
1417    local_start = prt_count+1
1418    prt_count   = local_count
1419!
1420!-- Calculate particle IDs
1421    DO  ip = nxl, nxr
1422       DO  jp = nys, nyn
1423          DO  kp = nzb+1, nzt
1424             number_of_particles = prt_count(kp,jp,ip)
1425             IF ( number_of_particles <= 0 )  CYCLE
1426             particles => grid_particles(kp,jp,ip)%particles(1:number_of_particles)
1427
1428             DO  n = local_start(kp,jp,ip), number_of_particles  !only new particles
1429
1430                particles(n)%id = 10000_idp**3 * grid_particles(kp,jp,ip)%id_counter + &
1431                                  10000_idp**2 * kp + 10000_idp * jp + ip
1432!
1433!--             Count the number of particles that have been released before
1434                grid_particles(kp,jp,ip)%id_counter =                          &
1435                                         grid_particles(kp,jp,ip)%id_counter + 1
1436
1437             ENDDO
1438
1439          ENDDO
1440       ENDDO
1441    ENDDO
1442!
1443!-- Initialize aerosol background spectrum
1444    IF ( curvature_solution_effects )  THEN
1445       CALL lpm_init_aerosols(local_start)
1446    ENDIF
1447!
1448!-- Add random fluctuation to particle positions.
1449    IF ( random_start_position )  THEN
1450       DO  ip = nxl, nxr
1451          DO  jp = nys, nyn
1452             DO  kp = nzb+1, nzt
1453                number_of_particles = prt_count(kp,jp,ip)
1454                IF ( number_of_particles <= 0 )  CYCLE
1455                particles => grid_particles(kp,jp,ip)%particles(1:number_of_particles)
1456!
1457!--             Move only new particles. Moreover, limit random fluctuation
1458!--             in order to prevent that particles move more than one grid box,
1459!--             which would lead to problems concerning particle exchange
1460!--             between processors in case pdx/pdy are larger than dx/dy,
1461!--             respectively.
1462                DO  n = local_start(kp,jp,ip), number_of_particles
1463                   IF ( psl(particles(n)%group) /= psr(particles(n)%group) )  THEN
1464                      rand_contr = ( random_function( iran_part ) - 0.5_wp ) * &
1465                                     pdx(particles(n)%group)
1466                      particles(n)%x = particles(n)%x +                        &
1467                              MERGE( rand_contr, SIGN( dx, rand_contr ),       &
1468                                     ABS( rand_contr ) < dx                    &
1469                                   )
1470                   ENDIF
1471                   IF ( pss(particles(n)%group) /= psn(particles(n)%group) )  THEN
1472                      rand_contr = ( random_function( iran_part ) - 0.5_wp ) * &
1473                                     pdy(particles(n)%group)
1474                      particles(n)%y = particles(n)%y +                        &
1475                              MERGE( rand_contr, SIGN( dy, rand_contr ),       &
1476                                     ABS( rand_contr ) < dy                    &
1477                                   )
1478                   ENDIF
1479                   IF ( psb(particles(n)%group) /= pst(particles(n)%group) )  THEN
1480                      rand_contr = ( random_function( iran_part ) - 0.5_wp ) * &
1481                                     pdz(particles(n)%group)
1482                      particles(n)%z = particles(n)%z +                        &
1483                              MERGE( rand_contr, SIGN( dzw(kp), rand_contr ),  &
1484                                     ABS( rand_contr ) < dzw(kp)               &
1485                                   )
1486                   ENDIF
1487                ENDDO
1488!
1489!--             Identify particles located outside the model domain and reflect
1490!--             or absorb them if necessary.
1491                CALL lpm_boundary_conds( 'bottom/top', i, j, k )
1492!
1493!--             Furthermore, remove particles located in topography. Note, as
1494!--             the particle speed is still zero at this point, wall
1495!--             reflection boundary conditions will not work in this case.
1496                particles =>                                                   &
1497                       grid_particles(kp,jp,ip)%particles(1:number_of_particles)
1498                DO  n = local_start(kp,jp,ip), number_of_particles
1499                   i = particles(n)%x * ddx
1500                   j = particles(n)%y * ddy
1501                   k = particles(n)%z / dz(1) + 1 + offset_ocean_nzt
1502                   DO WHILE( zw(k) < particles(n)%z )
1503                      k = k + 1
1504                   ENDDO
1505                   DO WHILE( zw(k-1) > particles(n)%z )
1506                      k = k - 1
1507                   ENDDO
1508!
1509!--                Check if particle is within topography
1510                   IF ( .NOT. BTEST( wall_flags_0(k,j,i), 0 ) )  THEN
1511                      particles(n)%particle_mask = .FALSE.
1512                      deleted_particles = deleted_particles + 1
1513                   ENDIF
1514
1515                ENDDO
1516             ENDDO
1517          ENDDO
1518       ENDDO
1519!
1520!--    Exchange particles between grid cells and processors
1521       CALL lpm_move_particle
1522       CALL lpm_exchange_horiz
1523
1524    ENDIF
1525!
1526!-- In case of random_start_position, delete particles identified by
1527!-- lpm_exchange_horiz and lpm_boundary_conds. Then sort particles into blocks,
1528!-- which is needed for a fast interpolation of the LES fields on the particle
1529!-- position.
1530    CALL lpm_sort_in_subboxes
1531!
1532!-- Determine the current number of particles
1533    DO  ip = nxl, nxr
1534       DO  jp = nys, nyn
1535          DO  kp = nzb+1, nzt
1536             number_of_particles         = number_of_particles                 &
1537                                           + prt_count(kp,jp,ip)
1538          ENDDO
1539       ENDDO
1540    ENDDO
1541!
1542!-- Calculate the number of particles of the total domain
1543#if defined( __parallel )
1544    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1545    CALL MPI_ALLREDUCE( number_of_particles, total_number_of_particles, 1, &
1546    MPI_INTEGER, MPI_SUM, comm2d, ierr )
1547#else
1548    total_number_of_particles = number_of_particles
1549#endif
1550
1551    RETURN
1552
1553 END SUBROUTINE lpm_create_particle
1554 
1555 
1556!------------------------------------------------------------------------------!
1557! Description:
1558! ------------
1559!> This routine initialize the particles as aerosols with physio-chemical
1560!> properties.
1561!------------------------------------------------------------------------------!   
1562 SUBROUTINE lpm_init_aerosols(local_start)
1563
1564    REAL(wp)  :: afactor            !< curvature effects
1565    REAL(wp)  :: bfactor            !< solute effects
1566    REAL(wp)  :: dlogr              !< logarithmic width of radius bin
1567    REAL(wp)  :: e_a                !< vapor pressure
1568    REAL(wp)  :: e_s                !< saturation vapor pressure
1569    REAL(wp)  :: rmin = 0.005e-6_wp !< minimum aerosol radius
1570    REAL(wp)  :: rmax = 10.0e-6_wp  !< maximum aerosol radius
1571    REAL(wp)  :: r_mid              !< mean radius of bin
1572    REAL(wp)  :: r_l                !< left radius of bin
1573    REAL(wp)  :: r_r                !< right radius of bin
1574    REAL(wp)  :: sigma              !< surface tension
1575    REAL(wp)  :: t_int              !< temperature
1576
1577    INTEGER(iwp), DIMENSION(nzb:nzt+1,nysg:nyng,nxlg:nxrg), INTENT(IN) ::  local_start !<
1578
1579    INTEGER(iwp)  :: n              !<
1580    INTEGER(iwp)  :: ip             !<
1581    INTEGER(iwp)  :: jp             !<
1582    INTEGER(iwp)  :: kp             !<
1583
1584!
1585!-- Set constants for different aerosol species
1586    IF ( TRIM(aero_species) .EQ. 'nacl' ) THEN
1587       molecular_weight_of_solute = 0.05844_wp 
1588       rho_s                      = 2165.0_wp
1589       vanthoff                   = 2.0_wp
1590    ELSEIF ( TRIM(aero_species) .EQ. 'c3h4o4' ) THEN
1591       molecular_weight_of_solute = 0.10406_wp 
1592       rho_s                      = 1600.0_wp
1593       vanthoff                   = 1.37_wp
1594    ELSEIF ( TRIM(aero_species) .EQ. 'nh4o3' ) THEN
1595       molecular_weight_of_solute = 0.08004_wp 
1596       rho_s                      = 1720.0_wp
1597       vanthoff                   = 2.31_wp
1598    ELSE
1599       WRITE( message_string, * ) 'unknown aerosol species ',   &
1600                                'aero_species = "', TRIM( aero_species ), '"'
1601       CALL message( 'lpm_init', 'PA0470', 1, 2, 0, 6, 0 )
1602    ENDIF
1603!
1604!-- The following typical aerosol spectra are taken from Jaenicke (1993):
1605!-- Tropospheric aerosols. Published in Aerosol-Cloud-Climate Interactions.
1606    IF ( TRIM(aero_type) .EQ. 'polar' )  THEN
1607       na        = (/ 2.17e1, 1.86e-1, 3.04e-4 /) * 1.0E6
1608       rm        = (/ 0.0689, 0.375, 4.29 /) * 1.0E-6
1609       log_sigma = (/ 0.245, 0.300, 0.291 /)
1610    ELSEIF ( TRIM(aero_type) .EQ. 'background' )  THEN
1611       na        = (/ 1.29e2, 5.97e1, 6.35e1 /) * 1.0E6
1612       rm        = (/ 0.0036, 0.127, 0.259 /) * 1.0E-6
1613       log_sigma = (/ 0.645, 0.253, 0.425 /)
1614    ELSEIF ( TRIM(aero_type) .EQ. 'maritime' )  THEN
1615       na        = (/ 1.33e2, 6.66e1, 3.06e0 /) * 1.0E6
1616       rm        = (/ 0.0039, 0.133, 0.29 /) * 1.0E-6
1617       log_sigma = (/ 0.657, 0.210, 0.396 /)
1618    ELSEIF ( TRIM(aero_type) .EQ. 'continental' )  THEN
1619       na        = (/ 3.20e3, 2.90e3, 3.00e-1 /) * 1.0E6
1620       rm        = (/ 0.01, 0.058, 0.9 /) * 1.0E-6
1621       log_sigma = (/ 0.161, 0.217, 0.380 /)
1622    ELSEIF ( TRIM(aero_type) .EQ. 'desert' )  THEN
1623       na        = (/ 7.26e2, 1.14e3, 1.78e-1 /) * 1.0E6
1624       rm        = (/ 0.001, 0.0188, 10.8 /) * 1.0E-6
1625       log_sigma = (/ 0.247, 0.770, 0.438 /)
1626    ELSEIF ( TRIM(aero_type) .EQ. 'rural' )  THEN
1627       na        = (/ 6.65e3, 1.47e2, 1.99e3 /) * 1.0E6
1628       rm        = (/ 0.00739, 0.0269, 0.0419 /) * 1.0E-6
1629       log_sigma = (/ 0.225, 0.557, 0.266 /)
1630    ELSEIF ( TRIM(aero_type) .EQ. 'urban' )  THEN
1631       na        = (/ 9.93e4, 1.11e3, 3.64e4 /) * 1.0E6
1632       rm        = (/ 0.00651, 0.00714, 0.0248 /) * 1.0E-6
1633       log_sigma = (/ 0.245, 0.666, 0.337 /)
1634    ELSEIF ( TRIM(aero_type) .EQ. 'user' )  THEN
1635       CONTINUE
1636    ELSE
1637       WRITE( message_string, * ) 'unknown aerosol type ',   &
1638                                'aero_type = "', TRIM( aero_type ), '"'
1639       CALL message( 'lpm_init', 'PA0459', 1, 2, 0, 6, 0 )
1640    ENDIF
1641
1642    DO  ip = nxl, nxr
1643       DO  jp = nys, nyn
1644          DO  kp = nzb+1, nzt
1645
1646             number_of_particles = prt_count(kp,jp,ip)
1647             IF ( number_of_particles <= 0 )  CYCLE
1648             particles => grid_particles(kp,jp,ip)%particles(1:number_of_particles)
1649
1650             dlogr   = ( LOG10(rmax) - LOG10(rmin) ) / ( number_of_particles - local_start(kp,jp,ip) + 1 )
1651!
1652!--          Initialize the aerosols with a predefined spectral distribution
1653!--          of the dry radius (logarithmically increasing bins) and a varying
1654!--          weighting factor
1655             DO  n = local_start(kp,jp,ip), number_of_particles  !only new particles
1656
1657                r_l   = 10.0**( LOG10( rmin ) + (n-1) * dlogr )
1658                r_r   = 10.0**( LOG10( rmin ) + n * dlogr )
1659                r_mid = SQRT( r_l * r_r )
1660
1661                particles(n)%aux1          = r_mid
1662                particles(n)%weight_factor =                                           &
1663                   ( na(1) / ( SQRT( 2.0 * pi ) * log_sigma(1) ) *                     &
1664                     EXP( - LOG10( r_mid / rm(1) )**2 / ( 2.0 * log_sigma(1)**2 ) ) +  &
1665                     na(2) / ( SQRT( 2.0 * pi ) * log_sigma(2) ) *                     &
1666                     EXP( - LOG10( r_mid / rm(2) )**2 / ( 2.0 * log_sigma(2)**2 ) ) +  &
1667                     na(3) / ( SQRT( 2.0 * pi ) * log_sigma(3) ) *                     &
1668                     EXP( - LOG10( r_mid / rm(3) )**2 / ( 2.0 * log_sigma(3)**2 ) )    &
1669                   ) * ( LOG10(r_r) - LOG10(r_l) ) * ( dx * dy * dzw(kp) )
1670
1671!
1672!--             Multiply weight_factor with the namelist parameter aero_weight
1673!--             to increase or decrease the number of simulated aerosols
1674                particles(n)%weight_factor = particles(n)%weight_factor * aero_weight
1675
1676                IF ( particles(n)%weight_factor - FLOOR(particles(n)%weight_factor,KIND=wp) &
1677                     .GT. random_function( iran_part ) )  THEN
1678                   particles(n)%weight_factor = FLOOR(particles(n)%weight_factor,KIND=wp) + 1.0_wp
1679                ELSE
1680                   particles(n)%weight_factor = FLOOR(particles(n)%weight_factor,KIND=wp)
1681                ENDIF
1682!
1683!--             Unnecessary particles will be deleted
1684                IF ( particles(n)%weight_factor .LE. 0.0 )  particles(n)%particle_mask = .FALSE.
1685
1686             ENDDO
1687!
1688!--          Set particle radius to equilibrium radius based on the environmental
1689!--          supersaturation (Khvorostyanov and Curry, 2007, JGR). This avoids
1690!--          the sometimes lengthy growth toward their equilibrium radius within
1691!--          the simulation.
1692             t_int  = pt(kp,jp,ip) * exner(kp)
1693
1694             e_s = magnus( t_int )
1695             e_a = q(kp,jp,ip) * hyp(kp) / ( q(kp,jp,ip) + rd_d_rv )
1696
1697             sigma   = 0.0761_wp - 0.000155_wp * ( t_int - 273.15_wp )
1698             afactor = 2.0_wp * sigma / ( rho_l * r_v * t_int )
1699
1700             bfactor = vanthoff * molecular_weight_of_water *    &
1701                       rho_s / ( molecular_weight_of_solute * rho_l )
1702!
1703!--          The formula is only valid for subsaturated environments. For
1704!--          supersaturations higher than -5 %, the supersaturation is set to -5%.
1705             IF ( e_a / e_s >= 0.95_wp )  e_a = 0.95_wp * e_s
1706
1707             DO  n = local_start(kp,jp,ip), number_of_particles  !only new particles
1708!
1709!--             For details on this equation, see Eq. (14) of Khvorostyanov and
1710!--             Curry (2007, JGR)
1711                particles(n)%radius = bfactor**0.3333333_wp *                  &
1712                   particles(n)%aux1 / ( 1.0_wp - e_a / e_s )**0.3333333_wp / &
1713                   ( 1.0_wp + ( afactor / ( 3.0_wp * bfactor**0.3333333_wp *   &
1714                     particles(n)%aux1 ) ) /                                  &
1715                     ( 1.0_wp - e_a / e_s )**0.6666666_wp                      &
1716                   )
1717
1718             ENDDO
1719
1720          ENDDO
1721       ENDDO
1722    ENDDO
1723
1724 END SUBROUTINE lpm_init_aerosols
1725 
1726 
1727!------------------------------------------------------------------------------!
1728! Description:
1729! ------------
1730!> Calculates quantities required for considering the SGS velocity fluctuations
1731!> in the particle transport by a stochastic approach. The respective
1732!> quantities are: SGS-TKE gradients and horizontally averaged profiles of the
1733!> SGS TKE and the resolved-scale velocity variances.
1734!------------------------------------------------------------------------------!
1735 SUBROUTINE lpm_init_sgs_tke
1736   
1737    USE statistics,                                                            &
1738        ONLY:  flow_statistics_called, hom, sums, sums_l
1739
1740    INTEGER(iwp) ::  i      !< index variable along x
1741    INTEGER(iwp) ::  j      !< index variable along y
1742    INTEGER(iwp) ::  k      !< index variable along z
1743    INTEGER(iwp) ::  m      !< running index for the surface elements
1744
1745    REAL(wp) ::  flag1      !< flag to mask topography
1746
1747!
1748!-- TKE gradient along x and y
1749    DO  i = nxl, nxr
1750       DO  j = nys, nyn
1751          DO  k = nzb, nzt+1
1752
1753             IF ( .NOT. BTEST( wall_flags_0(k,j,i-1), 0 )  .AND.               &
1754                        BTEST( wall_flags_0(k,j,i), 0   )  .AND.               &
1755                        BTEST( wall_flags_0(k,j,i+1), 0 ) )                    &
1756             THEN
1757                de_dx(k,j,i) = 2.0_wp * sgs_wf_part *                          &
1758                               ( e(k,j,i+1) - e(k,j,i) ) * ddx
1759             ELSEIF ( BTEST( wall_flags_0(k,j,i-1), 0 )  .AND.                 &
1760                      BTEST( wall_flags_0(k,j,i), 0   )  .AND.                 &
1761                .NOT. BTEST( wall_flags_0(k,j,i+1), 0 ) )                      &
1762             THEN
1763                de_dx(k,j,i) = 2.0_wp * sgs_wf_part *                          &
1764                               ( e(k,j,i) - e(k,j,i-1) ) * ddx
1765             ELSEIF ( .NOT. BTEST( wall_flags_0(k,j,i), 22   )  .AND.          &
1766                      .NOT. BTEST( wall_flags_0(k,j,i+1), 22 ) )               &   
1767             THEN
1768                de_dx(k,j,i) = 0.0_wp
1769             ELSEIF ( .NOT. BTEST( wall_flags_0(k,j,i-1), 22 )  .AND.          &
1770                      .NOT. BTEST( wall_flags_0(k,j,i), 22   ) )               &
1771             THEN
1772                de_dx(k,j,i) = 0.0_wp
1773             ELSE
1774                de_dx(k,j,i) = sgs_wf_part * ( e(k,j,i+1) - e(k,j,i-1) ) * ddx
1775             ENDIF
1776
1777             IF ( .NOT. BTEST( wall_flags_0(k,j-1,i), 0 )  .AND.               &
1778                        BTEST( wall_flags_0(k,j,i), 0   )  .AND.               &
1779                        BTEST( wall_flags_0(k,j+1,i), 0 ) )                    &
1780             THEN
1781                de_dy(k,j,i) = 2.0_wp * sgs_wf_part *                          &
1782                               ( e(k,j+1,i) - e(k,j,i) ) * ddy
1783             ELSEIF ( BTEST( wall_flags_0(k,j-1,i), 0 )  .AND.                 &
1784                      BTEST( wall_flags_0(k,j,i), 0   )  .AND.                 &
1785                .NOT. BTEST( wall_flags_0(k,j+1,i), 0 ) )                      &
1786             THEN
1787                de_dy(k,j,i) = 2.0_wp * sgs_wf_part *                          &
1788                               ( e(k,j,i) - e(k,j-1,i) ) * ddy
1789             ELSEIF ( .NOT. BTEST( wall_flags_0(k,j,i), 22   )  .AND.          &
1790                      .NOT. BTEST( wall_flags_0(k,j+1,i), 22 ) )               &   
1791             THEN
1792                de_dy(k,j,i) = 0.0_wp
1793             ELSEIF ( .NOT. BTEST( wall_flags_0(k,j-1,i), 22 )  .AND.          &
1794                      .NOT. BTEST( wall_flags_0(k,j,i), 22   ) )               &
1795             THEN
1796                de_dy(k,j,i) = 0.0_wp
1797             ELSE
1798                de_dy(k,j,i) = sgs_wf_part * ( e(k,j+1,i) - e(k,j-1,i) ) * ddy
1799             ENDIF
1800
1801          ENDDO
1802       ENDDO
1803    ENDDO
1804
1805!
1806!-- TKE gradient along z at topograhy and  including bottom and top boundary conditions
1807    DO  i = nxl, nxr
1808       DO  j = nys, nyn
1809          DO  k = nzb+1, nzt-1
1810!
1811!--          Flag to mask topography
1812             flag1 = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_0(k,j,i), 0  ) )
1813
1814             de_dz(k,j,i) = 2.0_wp * sgs_wf_part *                             &
1815                           ( e(k+1,j,i) - e(k-1,j,i) ) / ( zu(k+1) - zu(k-1) ) &
1816                                                 * flag1 
1817          ENDDO
1818!
1819!--       upward-facing surfaces
1820          DO  m = bc_h(0)%start_index(j,i), bc_h(0)%end_index(j,i)
1821             k            = bc_h(0)%k(m)
1822             de_dz(k,j,i) = 2.0_wp * sgs_wf_part *                             &
1823                           ( e(k+1,j,i) - e(k,j,i)   ) / ( zu(k+1) - zu(k) )
1824          ENDDO
1825!
1826!--       downward-facing surfaces
1827          DO  m = bc_h(1)%start_index(j,i), bc_h(1)%end_index(j,i)
1828             k            = bc_h(1)%k(m)
1829             de_dz(k,j,i) = 2.0_wp * sgs_wf_part *                             &
1830                           ( e(k,j,i) - e(k-1,j,i)   ) / ( zu(k) - zu(k-1) )
1831          ENDDO
1832
1833          de_dz(nzb,j,i)   = 0.0_wp
1834          de_dz(nzt,j,i)   = 0.0_wp
1835          de_dz(nzt+1,j,i) = 0.0_wp
1836       ENDDO
1837    ENDDO
1838!
1839!-- Ghost point exchange
1840    CALL exchange_horiz( de_dx, nbgp )
1841    CALL exchange_horiz( de_dy, nbgp )
1842    CALL exchange_horiz( de_dz, nbgp )
1843    CALL exchange_horiz( diss, nbgp  )
1844!
1845!-- Set boundary conditions at non-periodic boundaries. Note, at non-period
1846!-- boundaries zero-gradient boundary conditions are set for the subgrid TKE.
1847!-- Thus, TKE gradients normal to the respective lateral boundaries are zero,
1848!-- while tangetial TKE gradients then must be the same as within the prognostic
1849!-- domain. 
1850    IF ( bc_dirichlet_l )  THEN
1851       de_dx(:,:,-1) = 0.0_wp
1852       de_dy(:,:,-1) = de_dy(:,:,0) 
1853       de_dz(:,:,-1) = de_dz(:,:,0)
1854    ENDIF
1855    IF ( bc_dirichlet_r )  THEN
1856       de_dx(:,:,nxr+1) = 0.0_wp
1857       de_dy(:,:,nxr+1) = de_dy(:,:,nxr) 
1858       de_dz(:,:,nxr+1) = de_dz(:,:,nxr)
1859    ENDIF
1860    IF ( bc_dirichlet_n )  THEN
1861       de_dx(:,nyn+1,:) = de_dx(:,nyn,:)
1862       de_dy(:,nyn+1,:) = 0.0_wp 
1863       de_dz(:,nyn+1,:) = de_dz(:,nyn,:)
1864    ENDIF
1865    IF ( bc_dirichlet_s )  THEN
1866       de_dx(:,nys-1,:) = de_dx(:,nys,:)
1867       de_dy(:,nys-1,:) = 0.0_wp 
1868       de_dz(:,nys-1,:) = de_dz(:,nys,:)
1869    ENDIF 
1870!
1871!-- Calculate the horizontally averaged profiles of SGS TKE and resolved
1872!-- velocity variances (they may have been already calculated in routine
1873!-- flow_statistics).
1874    IF ( .NOT. flow_statistics_called )  THEN
1875
1876!
1877!--    First calculate horizontally averaged profiles of the horizontal
1878!--    velocities.
1879       sums_l(:,1,0) = 0.0_wp
1880       sums_l(:,2,0) = 0.0_wp
1881
1882       DO  i = nxl, nxr
1883          DO  j =  nys, nyn
1884             DO  k = nzb, nzt+1
1885!
1886!--             Flag indicating vicinity of wall
1887                flag1 = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_0(k,j,i), 24 ) )
1888
1889                sums_l(k,1,0)  = sums_l(k,1,0)  + u(k,j,i) * flag1
1890                sums_l(k,2,0)  = sums_l(k,2,0)  + v(k,j,i) * flag1
1891             ENDDO
1892          ENDDO
1893       ENDDO
1894
1895#if defined( __parallel )
1896!
1897!--    Compute total sum from local sums
1898       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1899       CALL MPI_ALLREDUCE( sums_l(nzb,1,0), sums(nzb,1), nzt+2-nzb, &
1900                           MPI_REAL, MPI_SUM, comm2d, ierr )
1901       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1902       CALL MPI_ALLREDUCE( sums_l(nzb,2,0), sums(nzb,2), nzt+2-nzb, &
1903                              MPI_REAL, MPI_SUM, comm2d, ierr )
1904#else
1905       sums(:,1) = sums_l(:,1,0)
1906       sums(:,2) = sums_l(:,2,0)
1907#endif
1908
1909!
1910!--    Final values are obtained by division by the total number of grid
1911!--    points used for the summation.
1912       hom(:,1,1,0) = sums(:,1) / ngp_2dh_outer(:,0)   ! u
1913       hom(:,1,2,0) = sums(:,2) / ngp_2dh_outer(:,0)   ! v
1914
1915!
1916!--    Now calculate the profiles of SGS TKE and the resolved-scale
1917!--    velocity variances
1918       sums_l(:,8,0)  = 0.0_wp
1919       sums_l(:,30,0) = 0.0_wp
1920       sums_l(:,31,0) = 0.0_wp
1921       sums_l(:,32,0) = 0.0_wp
1922       DO  i = nxl, nxr
1923          DO  j = nys, nyn
1924             DO  k = nzb, nzt+1
1925!
1926!--             Flag indicating vicinity of wall
1927                flag1 = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_0(k,j,i), 24 ) )
1928
1929                sums_l(k,8,0)  = sums_l(k,8,0)  + e(k,j,i)                       * flag1
1930                sums_l(k,30,0) = sums_l(k,30,0) + ( u(k,j,i) - hom(k,1,1,0) )**2 * flag1
1931                sums_l(k,31,0) = sums_l(k,31,0) + ( v(k,j,i) - hom(k,1,2,0) )**2 * flag1
1932                sums_l(k,32,0) = sums_l(k,32,0) + w(k,j,i)**2                    * flag1
1933             ENDDO
1934          ENDDO
1935       ENDDO
1936
1937#if defined( __parallel )
1938!
1939!--    Compute total sum from local sums
1940       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1941       CALL MPI_ALLREDUCE( sums_l(nzb,8,0), sums(nzb,8), nzt+2-nzb, &
1942                           MPI_REAL, MPI_SUM, comm2d, ierr )
1943       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1944       CALL MPI_ALLREDUCE( sums_l(nzb,30,0), sums(nzb,30), nzt+2-nzb, &
1945                           MPI_REAL, MPI_SUM, comm2d, ierr )
1946       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1947       CALL MPI_ALLREDUCE( sums_l(nzb,31,0), sums(nzb,31), nzt+2-nzb, &
1948                           MPI_REAL, MPI_SUM, comm2d, ierr )
1949       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1950       CALL MPI_ALLREDUCE( sums_l(nzb,32,0), sums(nzb,32), nzt+2-nzb, &
1951                           MPI_REAL, MPI_SUM, comm2d, ierr )
1952
1953#else
1954       sums(:,8)  = sums_l(:,8,0)
1955       sums(:,30) = sums_l(:,30,0)
1956       sums(:,31) = sums_l(:,31,0)
1957       sums(:,32) = sums_l(:,32,0)
1958#endif
1959
1960!
1961!--    Final values are obtained by division by the total number of grid
1962!--    points used for the summation.
1963       hom(:,1,8,0)  = sums(:,8)  / ngp_2dh_outer(:,0)   ! e
1964       hom(:,1,30,0) = sums(:,30) / ngp_2dh_outer(:,0)   ! u*2
1965       hom(:,1,31,0) = sums(:,31) / ngp_2dh_outer(:,0)   ! v*2
1966       hom(:,1,32,0) = sums(:,32) / ngp_2dh_outer(:,0)   ! w*2
1967
1968    ENDIF
1969
1970 END SUBROUTINE lpm_init_sgs_tke
1971 
1972 
1973!------------------------------------------------------------------------------!
1974! Description:
1975! ------------
1976!> Sobroutine control lpm actions, i.e. all actions during one time step.
1977!------------------------------------------------------------------------------! 
1978 SUBROUTINE lpm_actions( location )
1979
1980    CHARACTER (LEN=*), INTENT(IN) ::  location !< call location string
1981
1982    INTEGER(iwp)       ::  i                  !<
1983    INTEGER(iwp)       ::  ie                 !<
1984    INTEGER(iwp)       ::  is                 !<
1985    INTEGER(iwp)       ::  j                  !<
1986    INTEGER(iwp)       ::  je                 !<
1987    INTEGER(iwp)       ::  js                 !<
1988    INTEGER(iwp), SAVE ::  lpm_count = 0      !<
1989    INTEGER(iwp)       ::  k                  !<
1990    INTEGER(iwp)       ::  ke                 !<
1991    INTEGER(iwp)       ::  ks                 !<
1992    INTEGER(iwp)       ::  m                  !<
1993    INTEGER(iwp), SAVE ::  steps = 0          !<
1994
1995    LOGICAL            ::  first_loop_stride  !<
1996
1997   
1998    SELECT CASE ( location )
1999
2000       CASE ( 'after_prognostic_equations' )
2001
2002          CALL cpu_log( log_point(25), 'lpm', 'start' )
2003!
2004!--       Write particle data at current time on file.
2005!--       This has to be done here, before particles are further processed,
2006!--       because they may be deleted within this timestep (in case that
2007!--       dt_write_particle_data = dt_prel = particle_maximum_age).
2008          time_write_particle_data = time_write_particle_data + dt_3d
2009          IF ( time_write_particle_data >= dt_write_particle_data )  THEN
2010
2011             CALL lpm_data_output_particles
2012!
2013!--       The MOD function allows for changes in the output interval with restart
2014!--       runs.
2015             time_write_particle_data = MOD( time_write_particle_data, &
2016                                        MAX( dt_write_particle_data, dt_3d ) )
2017          ENDIF
2018
2019!
2020!--       Initialize arrays for marking those particles to be deleted after the
2021!--       (sub-) timestep
2022          deleted_particles = 0
2023
2024!
2025!--       Initialize variables used for accumulating the number of particles
2026!--       xchanged between the subdomains during all sub-timesteps (if sgs
2027!--       velocities are included). These data are output further below on the
2028!--       particle statistics file.
2029          trlp_count_sum      = 0
2030          trlp_count_recv_sum = 0
2031          trrp_count_sum      = 0
2032          trrp_count_recv_sum = 0
2033          trsp_count_sum      = 0
2034          trsp_count_recv_sum = 0
2035          trnp_count_sum      = 0
2036          trnp_count_recv_sum = 0
2037!
2038!--       Calculate exponential term used in case of particle inertia for each
2039!--       of the particle groups
2040          DO  m = 1, number_of_particle_groups
2041             IF ( particle_groups(m)%density_ratio /= 0.0_wp )  THEN
2042                particle_groups(m)%exp_arg  =                                        &
2043                          4.5_wp * particle_groups(m)%density_ratio *                &
2044                          molecular_viscosity / ( particle_groups(m)%radius )**2
2045
2046                particle_groups(m)%exp_term = EXP( -particle_groups(m)%exp_arg *     &
2047                          dt_3d )
2048             ENDIF
2049          ENDDO
2050!
2051!--       If necessary, release new set of particles
2052          IF ( ( simulated_time - last_particle_release_time ) >= dt_prel  .AND. end_time_prel > simulated_time ) &
2053          THEN
2054             DO WHILE ( ( simulated_time - last_particle_release_time ) >= dt_prel )
2055                CALL lpm_create_particle( PHASE_RELEASE )
2056                last_particle_release_time = last_particle_release_time + dt_prel
2057             ENDDO
2058          ENDIF
2059!
2060!--       Reset summation arrays
2061          IF ( cloud_droplets )  THEN
2062             ql_c  = 0.0_wp
2063             ql_v  = 0.0_wp
2064             ql_vp = 0.0_wp
2065          ENDIF
2066
2067          first_loop_stride = .TRUE.
2068          grid_particles(:,:,:)%time_loop_done = .TRUE.
2069!
2070!--       Timestep loop for particle advection.
2071!--       This loop has to be repeated until the advection time of every particle
2072!--       (within the total domain!) has reached the LES timestep (dt_3d).
2073!--       In case of including the SGS velocities, the particle timestep may be
2074!--       smaller than the LES timestep (because of the Lagrangian timescale
2075!--       restriction) and particles may require to undergo several particle
2076!--       timesteps, before the LES timestep is reached. Because the number of these
2077!--       particle timesteps to be carried out is unknown at first, these steps are
2078!--       carried out in the following infinite loop with exit condition.
2079          DO
2080             CALL cpu_log( log_point_s(44), 'lpm_advec', 'start' )
2081             CALL cpu_log( log_point_s(44), 'lpm_advec', 'pause' )
2082             
2083!
2084!--          If particle advection includes SGS velocity components, calculate the
2085!--          required SGS quantities (i.e. gradients of the TKE, as well as
2086!--          horizontally averaged profiles of the SGS TKE and the resolved-scale
2087!--          velocity variances)
2088             IF ( use_sgs_for_particles  .AND.  .NOT. cloud_droplets )  THEN
2089                CALL lpm_init_sgs_tke
2090             ENDIF
2091!
2092!--          In case SGS-particle speed is considered, particles may carry out
2093!--          several particle timesteps. In order to prevent unnecessary
2094!--          treatment of particles that already reached the final time level,
2095!--          particles are sorted into contiguous blocks of finished and
2096!--          not-finished particles, in addition to their already sorting
2097!--          according to their sub-boxes.
2098             IF ( .NOT. first_loop_stride  .AND.  use_sgs_for_particles )            &
2099                CALL lpm_sort_timeloop_done
2100
2101             DO  i = nxl, nxr
2102                DO  j = nys, nyn
2103                   DO  k = nzb+1, nzt
2104
2105                      number_of_particles = prt_count(k,j,i)
2106!
2107!--                   If grid cell gets empty, flag must be true
2108                      IF ( number_of_particles <= 0 )  THEN
2109                         grid_particles(k,j,i)%time_loop_done = .TRUE.
2110                         CYCLE
2111                      ENDIF
2112
2113                      IF ( .NOT. first_loop_stride  .AND.  &
2114                           grid_particles(k,j,i)%time_loop_done ) CYCLE
2115
2116                      particles => grid_particles(k,j,i)%particles(1:number_of_particles)
2117
2118                      particles(1:number_of_particles)%particle_mask = .TRUE.
2119!
2120!--                   Initialize the variable storing the total time that a particle
2121!--                   has advanced within the timestep procedure
2122                      IF ( first_loop_stride )  THEN
2123                         particles(1:number_of_particles)%dt_sum = 0.0_wp
2124                      ENDIF
2125!
2126!--                   Particle (droplet) growth by condensation/evaporation and
2127!--                   collision
2128                      IF ( cloud_droplets  .AND.  first_loop_stride)  THEN
2129!
2130!--                      Droplet growth by condensation / evaporation
2131                         CALL lpm_droplet_condensation(i,j,k)
2132!
2133!--                      Particle growth by collision
2134                         IF ( collision_kernel /= 'none' )  THEN
2135                            CALL lpm_droplet_collision(i,j,k)
2136                         ENDIF
2137
2138                      ENDIF
2139!
2140!--                   Initialize the switch used for the loop exit condition checked
2141!--                   at the end of this loop. If at least one particle has failed to
2142!--                   reach the LES timestep, this switch will be set false in
2143!--                   lpm_advec.
2144                      dt_3d_reached_l = .TRUE.
2145
2146!
2147!--                   Particle advection
2148                      CALL lpm_advec(i,j,k)
2149!
2150!--                   Particle reflection from walls. Only applied if the particles
2151!--                   are in the vertical range of the topography. (Here, some
2152!--                   optimization is still possible.)
2153                      IF ( topography /= 'flat' .AND. k < nzb_max + 2 )  THEN
2154                         CALL  lpm_boundary_conds( 'walls', i, j, k )
2155                      ENDIF
2156!
2157!--                   User-defined actions after the calculation of the new particle
2158!--                   position
2159                      CALL user_lpm_advec(i,j,k)
2160!
2161!--                   Apply boundary conditions to those particles that have crossed
2162!--                   the top or bottom boundary and delete those particles, which are
2163!--                   older than allowed
2164                      CALL lpm_boundary_conds( 'bottom/top', i, j, k )
2165!
2166!---                  If not all particles of the actual grid cell have reached the
2167!--                   LES timestep, this cell has to do another loop iteration. Due to
2168!--                   the fact that particles can move into neighboring grid cells,
2169!--                   these neighbor cells also have to perform another loop iteration.
2170!--                   Please note, this realization does not work properly if
2171!--                   particles move into another subdomain.
2172                      IF ( .NOT. dt_3d_reached_l )  THEN
2173                         ks = MAX(nzb+1,k-1)
2174                         ke = MIN(nzt,k+1)
2175                         js = MAX(nys,j-1)
2176                         je = MIN(nyn,j+1)
2177                         is = MAX(nxl,i-1)
2178                         ie = MIN(nxr,i+1)
2179                         grid_particles(ks:ke,js:je,is:ie)%time_loop_done = .FALSE.
2180                      ELSE
2181                         grid_particles(k,j,i)%time_loop_done = .TRUE.
2182                      ENDIF
2183
2184                   ENDDO
2185                ENDDO
2186             ENDDO
2187
2188             steps = steps + 1
2189             dt_3d_reached_l = ALL(grid_particles(:,:,:)%time_loop_done)
2190!
2191!--          Find out, if all particles on every PE have completed the LES timestep
2192!--          and set the switch corespondingly
2193#if defined( __parallel )
2194             IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2195             CALL MPI_ALLREDUCE( dt_3d_reached_l, dt_3d_reached, 1, MPI_LOGICAL, &
2196                                 MPI_LAND, comm2d, ierr )
2197#else
2198             dt_3d_reached = dt_3d_reached_l
2199#endif
2200
2201             CALL cpu_log( log_point_s(44), 'lpm_advec', 'stop' )
2202
2203!
2204!--          Apply splitting and merging algorithm
2205             IF ( cloud_droplets )  THEN
2206                IF ( splitting ) THEN
2207                   CALL lpm_splitting
2208                ENDIF
2209                IF ( merging ) THEN
2210                   CALL lpm_merging
2211                ENDIF
2212             ENDIF
2213!
2214!--          Move Particles local to PE to a different grid cell
2215             CALL lpm_move_particle
2216!
2217!--          Horizontal boundary conditions including exchange between subdmains
2218             CALL lpm_exchange_horiz
2219
2220!
2221!--          IF .FALSE., lpm_sort_in_subboxes is done inside pcmp
2222             IF ( .NOT. dt_3d_reached .OR. .NOT. nested_run )   THEN   
2223!
2224!--             Pack particles (eliminate those marked for deletion),
2225!--             determine new number of particles
2226                CALL lpm_sort_in_subboxes
2227!
2228!--             Initialize variables for the next (sub-) timestep, i.e., for marking
2229!--             those particles to be deleted after the timestep
2230                deleted_particles = 0
2231             ENDIF
2232
2233             IF ( dt_3d_reached )  EXIT
2234
2235             first_loop_stride = .FALSE.
2236          ENDDO   ! timestep loop
2237!   
2238!--       in case of nested runs do the transfer of particles after every full model time step
2239          IF ( nested_run )   THEN
2240             CALL particles_from_parent_to_child
2241             CALL particles_from_child_to_parent
2242             CALL pmcp_p_delete_particles_in_fine_grid_area
2243
2244             CALL lpm_sort_in_subboxes
2245
2246             deleted_particles = 0
2247          ENDIF
2248
2249!
2250!--       Calculate the new liquid water content for each grid box
2251          IF ( cloud_droplets )  CALL lpm_calc_liquid_water_content
2252
2253!
2254!--       Deallocate unused memory
2255          IF ( deallocate_memory  .AND.  lpm_count == step_dealloc )  THEN
2256             CALL dealloc_particles_array
2257             lpm_count = 0
2258          ELSEIF ( deallocate_memory )  THEN
2259             lpm_count = lpm_count + 1
2260          ENDIF
2261
2262!
2263!--       Write particle statistics (in particular the number of particles
2264!--       exchanged between the subdomains) on file
2265          IF ( write_particle_statistics )  CALL lpm_write_exchange_statistics
2266
2267          CALL cpu_log( log_point(25), 'lpm', 'stop' )
2268         
2269! !
2270! !--       Output of particle time series
2271!           IF ( particle_advection )  THEN
2272!              IF ( time_dopts >= dt_dopts  .OR.                                                        &
2273!                   ( time_since_reference_point >= particle_advection_start  .AND.                     &
2274!                    first_call_lpm ) )  THEN
2275!                 CALL lpm_data_output_ptseries
2276!                 time_dopts = MOD( time_dopts, MAX( dt_dopts, dt_3d ) )
2277!              ENDIF
2278!           ENDIF
2279   
2280       CASE DEFAULT
2281          CONTINUE
2282
2283    END SELECT
2284
2285 END SUBROUTINE lpm_actions
2286 
2287 
2288!------------------------------------------------------------------------------!
2289! Description:
2290! ------------
2291!
2292!------------------------------------------------------------------------------!
2293 SUBROUTINE particles_from_parent_to_child
2294    IMPLICIT NONE
2295
2296    CALL pmcp_c_get_particle_from_parent                         ! Child actions
2297    CALL pmcp_p_fill_particle_win                                ! Parent actions
2298
2299    RETURN
2300 END SUBROUTINE particles_from_parent_to_child
2301
2302 
2303!------------------------------------------------------------------------------!
2304! Description:
2305! ------------
2306!
2307!------------------------------------------------------------------------------!
2308 SUBROUTINE particles_from_child_to_parent
2309    IMPLICIT NONE
2310
2311    CALL pmcp_c_send_particle_to_parent                         ! Child actions
2312    CALL pmcp_p_empty_particle_win                              ! Parent actions
2313
2314    RETURN
2315 END SUBROUTINE particles_from_child_to_parent
2316 
2317!------------------------------------------------------------------------------!
2318! Description:
2319! ------------
2320!> This routine write exchange statistics of the lpm in a ascii file.
2321!------------------------------------------------------------------------------!
2322 SUBROUTINE lpm_write_exchange_statistics
2323
2324    INTEGER(iwp) :: ip         !<
2325    INTEGER(iwp) :: jp         !<
2326    INTEGER(iwp) :: kp         !<
2327    INTEGER(iwp) :: tot_number_of_particles
2328
2329!
2330!-- Determine the current number of particles
2331    number_of_particles         = 0
2332    DO  ip = nxl, nxr
2333       DO  jp = nys, nyn
2334          DO  kp = nzb+1, nzt
2335             number_of_particles = number_of_particles                         &
2336                                     + prt_count(kp,jp,ip)
2337          ENDDO
2338       ENDDO
2339    ENDDO
2340
2341    CALL check_open( 80 )
2342#if defined( __parallel )
2343    WRITE ( 80, 8000 )  current_timestep_number+1, simulated_time+dt_3d, &
2344                        number_of_particles, pleft, trlp_count_sum,      &
2345                        trlp_count_recv_sum, pright, trrp_count_sum,     &
2346                        trrp_count_recv_sum, psouth, trsp_count_sum,     &
2347                        trsp_count_recv_sum, pnorth, trnp_count_sum,     &
2348                        trnp_count_recv_sum
2349#else
2350    WRITE ( 80, 8000 )  current_timestep_number+1, simulated_time+dt_3d, &
2351                        number_of_particles
2352#endif
2353    CALL close_file( 80 )
2354
2355    IF ( number_of_particles > 0 ) THEN
2356        WRITE(9,*) 'number_of_particles ', number_of_particles,                &
2357                    current_timestep_number + 1, simulated_time + dt_3d
2358    ENDIF
2359
2360#if defined( __parallel )
2361    CALL MPI_ALLREDUCE( number_of_particles, tot_number_of_particles, 1,       &
2362                        MPI_INTEGER, MPI_SUM, comm2d, ierr )
2363#else
2364    tot_number_of_particles = number_of_particles
2365#endif
2366
2367    IF ( nested_run )  THEN
2368       CALL pmcp_g_print_number_of_particles( simulated_time+dt_3d,            &
2369                                              tot_number_of_particles)
2370    ENDIF
2371
2372!
2373!-- Formats
23748000 FORMAT (I6,1X,F7.2,4X,I10,5X,4(I3,1X,I4,'/',I4,2X),6X,I10)
2375
2376
2377 END SUBROUTINE lpm_write_exchange_statistics
2378 
2379
2380!------------------------------------------------------------------------------!
2381! Description:
2382! ------------
2383!> Write particle data in FORTRAN binary and/or netCDF format
2384!------------------------------------------------------------------------------!
2385 SUBROUTINE lpm_data_output_particles
2386 
2387    INTEGER(iwp) ::  ip !<
2388    INTEGER(iwp) ::  jp !<
2389    INTEGER(iwp) ::  kp !<
2390
2391    CALL cpu_log( log_point_s(40), 'lpm_data_output', 'start' )
2392
2393!
2394!-- Attention: change version number for unit 85 (in routine check_open)
2395!--            whenever the output format for this unit is changed!
2396    CALL check_open( 85 )
2397
2398    WRITE ( 85 )  simulated_time
2399    WRITE ( 85 )  prt_count
2400         
2401    DO  ip = nxl, nxr
2402       DO  jp = nys, nyn
2403          DO  kp = nzb+1, nzt
2404             number_of_particles = prt_count(kp,jp,ip)
2405             particles => grid_particles(kp,jp,ip)%particles(1:number_of_particles)
2406             IF ( number_of_particles <= 0 )  CYCLE
2407             WRITE ( 85 )  particles
2408          ENDDO
2409       ENDDO
2410    ENDDO
2411
2412    CALL close_file( 85 )
2413
2414
2415#if defined( __netcdf )
2416! !
2417! !-- Output in netCDF format
2418!     CALL check_open( 108 )
2419!
2420! !
2421! !-- Update the NetCDF time axis
2422!     prt_time_count = prt_time_count + 1
2423!
2424!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_time_prt, &
2425!                             (/ simulated_time /),        &
2426!                             start = (/ prt_time_count /), count = (/ 1 /) )
2427!     CALL netcdf_handle_error( 'lpm_data_output_particles', 1 )
2428!
2429! !
2430! !-- Output the real number of particles used
2431!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_rnop_prt, &
2432!                             (/ number_of_particles /),   &
2433!                             start = (/ prt_time_count /), count = (/ 1 /) )
2434!     CALL netcdf_handle_error( 'lpm_data_output_particles', 2 )
2435!
2436! !
2437! !-- Output all particle attributes
2438!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(1), particles%age,      &
2439!                             start = (/ 1, prt_time_count /),               &
2440!                             count = (/ maximum_number_of_particles /) )
2441!     CALL netcdf_handle_error( 'lpm_data_output_particles', 3 )
2442!
2443!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(2), particles%user,     &
2444!                             start = (/ 1, prt_time_count /),               &
2445!                             count = (/ maximum_number_of_particles /) )
2446!     CALL netcdf_handle_error( 'lpm_data_output_particles', 4 )
2447!
2448!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(3), particles%origin_x, &
2449!                             start = (/ 1, prt_time_count /),               &
2450!                             count = (/ maximum_number_of_particles /) )
2451!     CALL netcdf_handle_error( 'lpm_data_output_particles', 5 )
2452!
2453!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(4), particles%origin_y, &
2454!                             start = (/ 1, prt_time_count /),               &
2455!                             count = (/ maximum_number_of_particles /) )
2456!     CALL netcdf_handle_error( 'lpm_data_output_particles', 6 )
2457!
2458!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(5), particles%origin_z, &
2459!                             start = (/ 1, prt_time_count /),               &
2460!                             count = (/ maximum_number_of_particles /) )
2461!     CALL netcdf_handle_error( 'lpm_data_output_particles', 7 )
2462!
2463!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(6), particles%radius,   &
2464!                             start = (/ 1, prt_time_count /),               &
2465!                             count = (/ maximum_number_of_particles /) )
2466!     CALL netcdf_handle_error( 'lpm_data_output_particles', 8 )
2467!
2468!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(7), particles%speed_x,  &
2469!                             start = (/ 1, prt_time_count /),               &
2470!                             count = (/ maximum_number_of_particles /) )
2471!     CALL netcdf_handle_error( 'lpm_data_output_particles', 9 )
2472!
2473!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(8), particles%speed_y,  &
2474!                             start = (/ 1, prt_time_count /),               &
2475!                             count = (/ maximum_number_of_particles /) )
2476!     CALL netcdf_handle_error( 'lpm_data_output_particles', 10 )
2477!
2478!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(9), particles%speed_z,  &
2479!                             start = (/ 1, prt_time_count /),               &
2480!                             count = (/ maximum_number_of_particles /) )
2481!     CALL netcdf_handle_error( 'lpm_data_output_particles', 11 )
2482!
2483!     nc_stat = NF90_PUT_VAR( id_set_prt,id_var_prt(10),                     &
2484!                             particles%weight_factor,                       &
2485!                             start = (/ 1, prt_time_count /),               &
2486!                             count = (/ maximum_number_of_particles /) )
2487!     CALL netcdf_handle_error( 'lpm_data_output_particles', 12 )
2488!
2489!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(11), particles%x,       &
2490!                             start = (/ 1, prt_time_count /),               &
2491!                             count = (/ maximum_number_of_particles /) )
2492!     CALL netcdf_handle_error( 'lpm_data_output_particles', 13 )
2493!
2494!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(12), particles%y,       &
2495!                             start = (/ 1, prt_time_count /),               &
2496!                             count = (/ maximum_number_of_particles /) )
2497!     CALL netcdf_handle_error( 'lpm_data_output_particles', 14 )
2498!
2499!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(13), particles%z,       &
2500!                             start = (/ 1, prt_time_count /),               &
2501!                             count = (/ maximum_number_of_particles /) )
2502!     CALL netcdf_handle_error( 'lpm_data_output_particles', 15 )
2503!
2504!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(14), particles%class,   &
2505!                             start = (/ 1, prt_time_count /),               &
2506!                             count = (/ maximum_number_of_particles /) )
2507!     CALL netcdf_handle_error( 'lpm_data_output_particles', 16 )
2508!
2509!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(15), particles%group,   &
2510!                             start = (/ 1, prt_time_count /),               &
2511!                             count = (/ maximum_number_of_particles /) )
2512!     CALL netcdf_handle_error( 'lpm_data_output_particles', 17 )
2513!
2514!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(16),                    &
2515!                             particles%id2,                                 &
2516!                             start = (/ 1, prt_time_count /),               &
2517!                             count = (/ maximum_number_of_particles /) )
2518!     CALL netcdf_handle_error( 'lpm_data_output_particles', 18 )
2519!
2520!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(17), particles%id1,     &
2521!                             start = (/ 1, prt_time_count /),               &
2522!                             count = (/ maximum_number_of_particles /) )
2523!     CALL netcdf_handle_error( 'lpm_data_output_particles', 19 )
2524!
2525#endif
2526
2527    CALL cpu_log( log_point_s(40), 'lpm_data_output', 'stop' )
2528
2529 END SUBROUTINE lpm_data_output_particles
2530 
2531!------------------------------------------------------------------------------!
2532! Description:
2533! ------------
2534!> This routine calculates and provide particle timeseries output.
2535!------------------------------------------------------------------------------!
2536 SUBROUTINE lpm_data_output_ptseries
2537 
2538    INTEGER(iwp) ::  i    !<
2539    INTEGER(iwp) ::  inum !<
2540    INTEGER(iwp) ::  j    !<
2541    INTEGER(iwp) ::  jg   !<
2542    INTEGER(iwp) ::  k    !<
2543    INTEGER(iwp) ::  n    !<
2544
2545    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  pts_value   !<
2546    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  pts_value_l !<
2547
2548
2549    CALL cpu_log( log_point(36), 'data_output_ptseries', 'start' )
2550
2551    IF ( myid == 0 )  THEN
2552!
2553!--    Open file for time series output in NetCDF format
2554       dopts_time_count = dopts_time_count + 1
2555       CALL check_open( 109 )
2556#if defined( __netcdf )
2557!
2558!--    Update the particle time series time axis
2559       nc_stat = NF90_PUT_VAR( id_set_pts, id_var_time_pts,      &
2560                               (/ time_since_reference_point /), &
2561                               start = (/ dopts_time_count /), count = (/ 1 /) )
2562       CALL netcdf_handle_error( 'data_output_ptseries', 391 )
2563#endif
2564
2565    ENDIF
2566
2567    ALLOCATE( pts_value(0:number_of_particle_groups,dopts_num), &
2568              pts_value_l(0:number_of_particle_groups,dopts_num) )
2569
2570    pts_value_l = 0.0_wp
2571    pts_value_l(:,16) = 9999999.9_wp    ! for calculation of minimum radius
2572
2573!
2574!-- Calculate or collect the particle time series quantities for all particles
2575!-- and seperately for each particle group (if there is more than one group)
2576    DO  i = nxl, nxr
2577       DO  j = nys, nyn
2578          DO  k = nzb, nzt
2579             number_of_particles = prt_count(k,j,i)
2580             IF (number_of_particles <= 0)  CYCLE
2581             particles => grid_particles(k,j,i)%particles(1:number_of_particles)
2582             DO  n = 1, number_of_particles
2583
2584                IF ( particles(n)%particle_mask )  THEN  ! Restrict analysis to active particles
2585
2586                   pts_value_l(0,1)  = pts_value_l(0,1) + 1.0_wp  ! total # of particles
2587                   pts_value_l(0,2)  = pts_value_l(0,2) +                      &
2588                          ( particles(n)%x - particles(n)%origin_x )  ! mean x
2589                   pts_value_l(0,3)  = pts_value_l(0,3) +                      &
2590                          ( particles(n)%y - particles(n)%origin_y )  ! mean y
2591                   pts_value_l(0,4)  = pts_value_l(0,4) +                      &
2592                          ( particles(n)%z - particles(n)%origin_z )  ! mean z
2593                   pts_value_l(0,5)  = pts_value_l(0,5) + particles(n)%z        ! mean z (absolute)
2594                   pts_value_l(0,6)  = pts_value_l(0,6) + particles(n)%speed_x  ! mean u
2595                   pts_value_l(0,7)  = pts_value_l(0,7) + particles(n)%speed_y  ! mean v
2596                   pts_value_l(0,8)  = pts_value_l(0,8) + particles(n)%speed_z  ! mean w
2597                   pts_value_l(0,9)  = pts_value_l(0,9)  + particles(n)%rvar1 ! mean sgsu
2598                   pts_value_l(0,10) = pts_value_l(0,10) + particles(n)%rvar2 ! mean sgsv
2599                   pts_value_l(0,11) = pts_value_l(0,11) + particles(n)%rvar3 ! mean sgsw
2600                   IF ( particles(n)%speed_z > 0.0_wp )  THEN
2601                      pts_value_l(0,12) = pts_value_l(0,12) + 1.0_wp  ! # of upward moving prts
2602                      pts_value_l(0,13) = pts_value_l(0,13) +                  &
2603                                              particles(n)%speed_z ! mean w upw.
2604                   ELSE
2605                      pts_value_l(0,14) = pts_value_l(0,14) +                  &
2606                                              particles(n)%speed_z ! mean w down
2607                   ENDIF
2608                   pts_value_l(0,15) = pts_value_l(0,15) + particles(n)%radius ! mean rad
2609                   pts_value_l(0,16) = MIN( pts_value_l(0,16), particles(n)%radius ) ! minrad
2610                   pts_value_l(0,17) = MAX( pts_value_l(0,17), particles(n)%radius ) ! maxrad
2611                   pts_value_l(0,18) = pts_value_l(0,18) + 1.0_wp
2612                   pts_value_l(0,19) = pts_value_l(0,18) + 1.0_wp
2613!
2614!--                Repeat the same for the respective particle group
2615                   IF ( number_of_particle_groups > 1 )  THEN
2616                      jg = particles(n)%group
2617
2618                      pts_value_l(jg,1)  = pts_value_l(jg,1) + 1.0_wp
2619                      pts_value_l(jg,2)  = pts_value_l(jg,2) +                   &
2620                           ( particles(n)%x - particles(n)%origin_x )
2621                      pts_value_l(jg,3)  = pts_value_l(jg,3) +                   &
2622                           ( particles(n)%y - particles(n)%origin_y )
2623                      pts_value_l(jg,4)  = pts_value_l(jg,4) +                   &
2624                           ( particles(n)%z - particles(n)%origin_z )
2625                      pts_value_l(jg,5)  = pts_value_l(jg,5) + particles(n)%z
2626                      pts_value_l(jg,6)  = pts_value_l(jg,6) + particles(n)%speed_x
2627                      pts_value_l(jg,7)  = pts_value_l(jg,7) + particles(n)%speed_y
2628                      pts_value_l(jg,8)  = pts_value_l(jg,8) + particles(n)%speed_z
2629                      pts_value_l(jg,9)  = pts_value_l(jg,9)  + particles(n)%rvar1
2630                      pts_value_l(jg,10) = pts_value_l(jg,10) + particles(n)%rvar2
2631                      pts_value_l(jg,11) = pts_value_l(jg,11) + particles(n)%rvar3
2632                      IF ( particles(n)%speed_z > 0.0_wp )  THEN
2633                         pts_value_l(jg,12) = pts_value_l(jg,12) + 1.0_wp
2634                         pts_value_l(jg,13) = pts_value_l(jg,13) + particles(n)%speed_z
2635                      ELSE
2636                         pts_value_l(jg,14) = pts_value_l(jg,14) + particles(n)%speed_z
2637                      ENDIF
2638                      pts_value_l(jg,15) = pts_value_l(jg,15) + particles(n)%radius
2639                      pts_value_l(jg,16) = MIN( pts_value(jg,16), particles(n)%radius )
2640                      pts_value_l(jg,17) = MAX( pts_value(jg,17), particles(n)%radius )
2641                      pts_value_l(jg,18) = pts_value_l(jg,18) + 1.0_wp
2642                      pts_value_l(jg,19) = pts_value_l(jg,19) + 1.0_wp
2643                   ENDIF
2644
2645                ENDIF
2646
2647             ENDDO
2648
2649          ENDDO
2650       ENDDO
2651    ENDDO
2652
2653
2654#if defined( __parallel )
2655!
2656!-- Sum values of the subdomains
2657    inum = number_of_particle_groups + 1
2658
2659    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2660    CALL MPI_ALLREDUCE( pts_value_l(0,1), pts_value(0,1), 15*inum, MPI_REAL, &
2661                        MPI_SUM, comm2d, ierr )
2662    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2663    CALL MPI_ALLREDUCE( pts_value_l(0,16), pts_value(0,16), inum, MPI_REAL, &
2664                        MPI_MIN, comm2d, ierr )
2665    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2666    CALL MPI_ALLREDUCE( pts_value_l(0,17), pts_value(0,17), inum, MPI_REAL, &
2667                        MPI_MAX, comm2d, ierr )
2668    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2669    CALL MPI_ALLREDUCE( pts_value_l(0,18), pts_value(0,18), inum, MPI_REAL, &
2670                        MPI_MAX, comm2d, ierr )
2671    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2672    CALL MPI_ALLREDUCE( pts_value_l(0,19), pts_value(0,19), inum, MPI_REAL, &
2673                        MPI_MIN, comm2d, ierr )
2674#else
2675    pts_value(:,1:19) = pts_value_l(:,1:19)
2676#endif
2677
2678!
2679!-- Normalize the above calculated quantities (except min/max values) with the
2680!-- total number of particles
2681    IF ( number_of_particle_groups > 1 )  THEN
2682       inum = number_of_particle_groups
2683    ELSE
2684       inum = 0
2685    ENDIF
2686
2687    DO  j = 0, inum
2688
2689       IF ( pts_value(j,1) > 0.0_wp )  THEN
2690
2691          pts_value(j,2:15) = pts_value(j,2:15) / pts_value(j,1)
2692          IF ( pts_value(j,12) > 0.0_wp  .AND.  pts_value(j,12) < 1.0_wp )  THEN
2693             pts_value(j,13) = pts_value(j,13) / pts_value(j,12)
2694             pts_value(j,14) = pts_value(j,14) / ( 1.0_wp - pts_value(j,12) )
2695          ELSEIF ( pts_value(j,12) == 0.0_wp )  THEN
2696             pts_value(j,13) = -1.0_wp
2697          ELSE
2698             pts_value(j,14) = -1.0_wp
2699          ENDIF
2700
2701       ENDIF
2702
2703    ENDDO
2704
2705!
2706!-- Calculate higher order moments of particle time series quantities,
2707!-- seperately for each particle group (if there is more than one group)
2708    DO  i = nxl, nxr
2709       DO  j = nys, nyn
2710          DO  k = nzb, nzt
2711             number_of_particles = prt_count(k,j,i)
2712             IF (number_of_particles <= 0)  CYCLE
2713             particles => grid_particles(k,j,i)%particles(1:number_of_particles)
2714             DO  n = 1, number_of_particles
2715
2716                pts_value_l(0,20) = pts_value_l(0,20) + ( particles(n)%x - &
2717                                    particles(n)%origin_x - pts_value(0,2) )**2 ! x*2
2718                pts_value_l(0,21) = pts_value_l(0,21) + ( particles(n)%y - &
2719                                    particles(n)%origin_y - pts_value(0,3) )**2 ! y*2
2720                pts_value_l(0,22) = pts_value_l(0,22) + ( particles(n)%z - &
2721                                    particles(n)%origin_z - pts_value(0,4) )**2 ! z*2
2722                pts_value_l(0,23) = pts_value_l(0,23) + ( particles(n)%speed_x - &
2723                                                         pts_value(0,6) )**2   ! u*2
2724                pts_value_l(0,24) = pts_value_l(0,24) + ( particles(n)%speed_y - &
2725                                                          pts_value(0,7) )**2   ! v*2
2726                pts_value_l(0,25) = pts_value_l(0,25) + ( particles(n)%speed_z - &
2727                                                          pts_value(0,8) )**2   ! w*2
2728                pts_value_l(0,26) = pts_value_l(0,26) + ( particles(n)%rvar1 - &
2729                                                          pts_value(0,9) )**2   ! u"2
2730                pts_value_l(0,27) = pts_value_l(0,27) + ( particles(n)%rvar2 - &
2731                                                          pts_value(0,10) )**2  ! v"2
2732                pts_value_l(0,28) = pts_value_l(0,28) + ( particles(n)%rvar3 - &
2733                                                          pts_value(0,11) )**2  ! w"2
2734!
2735!--             Repeat the same for the respective particle group
2736                IF ( number_of_particle_groups > 1 )  THEN
2737                   jg = particles(n)%group
2738
2739                   pts_value_l(jg,20) = pts_value_l(jg,20) + ( particles(n)%x - &
2740                                       particles(n)%origin_x - pts_value(jg,2) )**2
2741                   pts_value_l(jg,21) = pts_value_l(jg,21) + ( particles(n)%y - &
2742                                       particles(n)%origin_y - pts_value(jg,3) )**2
2743                   pts_value_l(jg,22) = pts_value_l(jg,22) + ( particles(n)%z - &
2744                                       particles(n)%origin_z - pts_value(jg,4) )**2
2745                   pts_value_l(jg,23) = pts_value_l(jg,23) + ( particles(n)%speed_x - &
2746                                                             pts_value(jg,6) )**2
2747                   pts_value_l(jg,24) = pts_value_l(jg,24) + ( particles(n)%speed_y - &
2748                                                             pts_value(jg,7) )**2
2749                   pts_value_l(jg,25) = pts_value_l(jg,25) + ( particles(n)%speed_z - &
2750                                                             pts_value(jg,8) )**2
2751                   pts_value_l(jg,26) = pts_value_l(jg,26) + ( particles(n)%rvar1 - &
2752                                                             pts_value(jg,9) )**2
2753                   pts_value_l(jg,27) = pts_value_l(jg,27) + ( particles(n)%rvar2 - &
2754                                                             pts_value(jg,10) )**2
2755                   pts_value_l(jg,28) = pts_value_l(jg,28) + ( particles(n)%rvar3 - &
2756                                                             pts_value(jg,11) )**2
2757                ENDIF
2758
2759             ENDDO
2760          ENDDO
2761       ENDDO
2762    ENDDO
2763
2764    pts_value_l(0,29) = ( number_of_particles - pts_value(0,1) / numprocs )**2
2765                                                 ! variance of particle numbers
2766    IF ( number_of_particle_groups > 1 )  THEN
2767       DO  j = 1, number_of_particle_groups
2768          pts_value_l(j,29) = ( pts_value_l(j,1) - &
2769                                pts_value(j,1) / numprocs )**2
2770       ENDDO
2771    ENDIF
2772
2773#if defined( __parallel )
2774!
2775!-- Sum values of the subdomains
2776    inum = number_of_particle_groups + 1
2777
2778    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2779    CALL MPI_ALLREDUCE( pts_value_l(0,20), pts_value(0,20), inum*10, MPI_REAL, &
2780                        MPI_SUM, comm2d, ierr )
2781#else
2782    pts_value(:,20:29) = pts_value_l(:,20:29)
2783#endif
2784
2785!
2786!-- Normalize the above calculated quantities with the total number of
2787!-- particles
2788    IF ( number_of_particle_groups > 1 )  THEN
2789       inum = number_of_particle_groups
2790    ELSE
2791       inum = 0
2792    ENDIF
2793
2794    DO  j = 0, inum
2795
2796       IF ( pts_value(j,1) > 0.0_wp )  THEN
2797          pts_value(j,20:28) = pts_value(j,20:28) / pts_value(j,1)
2798       ENDIF
2799       pts_value(j,29) = pts_value(j,29) / numprocs
2800
2801    ENDDO
2802
2803#if defined( __netcdf )
2804!
2805!-- Output particle time series quantities in NetCDF format
2806    IF ( myid == 0 )  THEN
2807       DO  j = 0, inum
2808          DO  i = 1, dopts_num
2809             nc_stat = NF90_PUT_VAR( id_set_pts, id_var_dopts(i,j),  &
2810                                     (/ pts_value(j,i) /),           &
2811                                     start = (/ dopts_time_count /), &
2812                                     count = (/ 1 /) )
2813             CALL netcdf_handle_error( 'data_output_ptseries', 392 )
2814          ENDDO
2815       ENDDO
2816    ENDIF
2817#endif
2818
2819    DEALLOCATE( pts_value, pts_value_l )
2820
2821    CALL cpu_log( log_point(36), 'data_output_ptseries', 'stop' )
2822
2823END SUBROUTINE lpm_data_output_ptseries
2824
2825 
2826!------------------------------------------------------------------------------!
2827! Description:
2828! ------------
2829!> This routine reads the respective restart data for the lpm.
2830!------------------------------------------------------------------------------!
2831 SUBROUTINE lpm_rrd_local_particles
2832
2833    CHARACTER (LEN=10) ::  particle_binary_version    !<
2834    CHARACTER (LEN=10) ::  version_on_file            !<
2835
2836    INTEGER(iwp) :: alloc_size !<
2837    INTEGER(iwp) :: ip         !<
2838    INTEGER(iwp) :: jp         !<
2839    INTEGER(iwp) :: kp         !<
2840
2841    TYPE(particle_type), DIMENSION(:), ALLOCATABLE :: tmp_particles !<
2842
2843!
2844!-- Read particle data from previous model run.
2845!-- First open the input unit.
2846    IF ( myid_char == '' )  THEN
2847       OPEN ( 90, FILE='PARTICLE_RESTART_DATA_IN'//myid_char,                  &
2848                  FORM='UNFORMATTED' )
2849    ELSE
2850       OPEN ( 90, FILE='PARTICLE_RESTART_DATA_IN/'//myid_char,                 &
2851                  FORM='UNFORMATTED' )
2852    ENDIF
2853
2854!
2855!-- First compare the version numbers
2856    READ ( 90 )  version_on_file
2857    particle_binary_version = '4.0'
2858    IF ( TRIM( version_on_file ) /= TRIM( particle_binary_version ) )  THEN
2859       message_string = 'version mismatch concerning data from prior ' //      &
2860                        'run &version on file = "' //                          &
2861                                      TRIM( version_on_file ) //               &
2862                        '&version in program = "' //                           &
2863                                      TRIM( particle_binary_version ) // '"'
2864       CALL message( 'lpm_read_restart_file', 'PA0214', 1, 2, 0, 6, 0 )
2865    ENDIF
2866
2867!
2868!-- If less particles are stored on the restart file than prescribed by
2869!-- min_nr_particle, the remainder is initialized by zero_particle to avoid
2870!-- errors.
2871    zero_particle = particle_type( 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp,     &
2872                                   0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp,     &
2873                                   0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp,     &
2874                                   0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp,     &
2875                                   0, 0, 0_idp, .FALSE., -1 )
2876!
2877!-- Read some particle parameters and the size of the particle arrays,
2878!-- allocate them and read their contents.
2879    READ ( 90 )  bc_par_b, bc_par_lr, bc_par_ns, bc_par_t,                     &
2880                 last_particle_release_time, number_of_particle_groups,        &
2881                 particle_groups, time_write_particle_data
2882
2883    ALLOCATE( prt_count(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                        &
2884              grid_particles(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
2885
2886    READ ( 90 )  prt_count
2887
2888    DO  ip = nxl, nxr
2889       DO  jp = nys, nyn
2890          DO  kp = nzb+1, nzt
2891
2892             number_of_particles = prt_count(kp,jp,ip)
2893             IF ( number_of_particles > 0 )  THEN
2894                alloc_size = MAX( INT( number_of_particles *                   &
2895                             ( 1.0_wp + alloc_factor / 100.0_wp ) ),           &
2896                             min_nr_particle )
2897             ELSE
2898                alloc_size = min_nr_particle
2899             ENDIF
2900
2901             ALLOCATE( grid_particles(kp,jp,ip)%particles(1:alloc_size) )
2902
2903             IF ( number_of_particles > 0 )  THEN
2904                ALLOCATE( tmp_particles(1:number_of_particles) )
2905                READ ( 90 )  tmp_particles
2906                grid_particles(kp,jp,ip)%particles(1:number_of_particles) = tmp_particles
2907                DEALLOCATE( tmp_particles )
2908                IF ( number_of_particles < alloc_size )  THEN
2909                   grid_particles(kp,jp,ip)%particles(number_of_particles+1:alloc_size) &
2910                      = zero_particle
2911                ENDIF
2912             ELSE
2913                grid_particles(kp,jp,ip)%particles(1:alloc_size) = zero_particle
2914             ENDIF
2915
2916          ENDDO
2917       ENDDO
2918    ENDDO
2919
2920    CLOSE ( 90 )
2921!
2922!-- Must be called to sort particles into blocks, which is needed for a fast
2923!-- interpolation of the LES fields on the particle position.
2924    CALL lpm_sort_in_subboxes
2925       
2926
2927 END SUBROUTINE lpm_rrd_local_particles
2928 
2929 
2930 SUBROUTINE lpm_rrd_local( k, nxlf, nxlc, nxl_on_file, nxrf, nxrc,          &
2931                              nxr_on_file, nynf, nync, nyn_on_file, nysf,  &
2932                              nysc, nys_on_file, tmp_3d, found )
2933                             
2934       
2935   USE control_parameters,                                                 &
2936       ONLY: length, restart_string           
2937
2938    INTEGER(iwp) ::  k               !<
2939    INTEGER(iwp) ::  nxlc            !<
2940    INTEGER(iwp) ::  nxlf            !<
2941    INTEGER(iwp) ::  nxl_on_file     !<
2942    INTEGER(iwp) ::  nxrc            !<
2943    INTEGER(iwp) ::  nxrf            !<
2944    INTEGER(iwp) ::  nxr_on_file     !<
2945    INTEGER(iwp) ::  nync            !<
2946    INTEGER(iwp) ::  nynf            !<
2947    INTEGER(iwp) ::  nyn_on_file     !<
2948    INTEGER(iwp) ::  nysc            !<
2949    INTEGER(iwp) ::  nysf            !<
2950    INTEGER(iwp) ::  nys_on_file     !<
2951
2952    LOGICAL, INTENT(OUT)  ::  found
2953
2954    REAL(wp), DIMENSION(nzb:nzt+1,nys_on_file-nbgp:nyn_on_file+nbgp,nxl_on_file-nbgp:nxr_on_file+nbgp) :: tmp_3d   !<
2955
2956   
2957    found = .TRUE.
2958
2959    SELECT CASE ( restart_string(1:length) )
2960   
2961       CASE ( 'iran' ) ! matching random numbers is still unresolved issue
2962          IF ( k == 1 )  READ ( 13 )  iran, iran_part
2963         
2964        CASE ( 'pc_av' )
2965           IF ( .NOT. ALLOCATED( pc_av ) )  THEN
2966              ALLOCATE( pc_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
2967           ENDIF
2968           IF ( k == 1 )  READ ( 13 )  tmp_3d
2969           pc_av(:,nysc-nbgp:nync+nbgp,nxlc-nbgp:nxrc+nbgp) =          &
2970              tmp_3d(:,nysf-nbgp:nynf+nbgp,nxlf-nbgp:nxrf+nbgp)
2971
2972        CASE ( 'pr_av' )
2973           IF ( .NOT. ALLOCATED( pr_av ) )  THEN
2974              ALLOCATE( pr_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
2975           ENDIF
2976           IF ( k == 1 )  READ ( 13 )  tmp_3d
2977           pr_av(:,nysc-nbgp:nync+nbgp,nxlc-nbgp:nxrc+nbgp) =          &
2978              tmp_3d(:,nysf-nbgp:nynf+nbgp,nxlf-nbgp:nxrf+nbgp)
2979 
2980         CASE ( 'ql_c_av' )
2981            IF ( .NOT. ALLOCATED( ql_c_av ) )  THEN
2982               ALLOCATE( ql_c_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
2983            ENDIF
2984            IF ( k == 1 )  READ ( 13 )  tmp_3d
2985            ql_c_av(:,nysc-nbgp:nync+nbgp,nxlc-nbgp:nxrc+nbgp) =        &
2986               tmp_3d(:,nysf-nbgp:nynf+nbgp,nxlf-nbgp:nxrf+nbgp)
2987
2988         CASE ( 'ql_v_av' )
2989            IF ( .NOT. ALLOCATED( ql_v_av ) )  THEN
2990               ALLOCATE( ql_v_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
2991            ENDIF
2992            IF ( k == 1 )  READ ( 13 )  tmp_3d
2993            ql_v_av(:,nysc-nbgp:nync+nbgp,nxlc-nbgp:nxrc+nbgp) =        &
2994               tmp_3d(:,nysf-nbgp:nynf+nbgp,nxlf-nbgp:nxrf+nbgp)
2995
2996         CASE ( 'ql_vp_av' )
2997            IF ( .NOT. ALLOCATED( ql_vp_av ) )  THEN
2998               ALLOCATE( ql_vp_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
2999            ENDIF
3000            IF ( k == 1 )  READ ( 13 )  tmp_3d
3001            ql_vp_av(:,nysc-nbgp:nync+nbgp,nxlc-nbgp:nxrc+nbgp) =       &
3002               tmp_3d(:,nysf-nbgp:nynf+nbgp,nxlf-nbgp:nxrf+nbgp)
3003
3004          CASE DEFAULT
3005
3006             found = .FALSE.
3007
3008       END SELECT
3009               
3010
3011 END SUBROUTINE lpm_rrd_local
3012 
3013!------------------------------------------------------------------------------!
3014! Description:
3015! ------------
3016!> This routine writes the respective restart data for the lpm.
3017!------------------------------------------------------------------------------!
3018 SUBROUTINE lpm_wrd_local
3019 
3020    CHARACTER (LEN=10) ::  particle_binary_version   !<
3021
3022    INTEGER(iwp) ::  ip                              !<
3023    INTEGER(iwp) ::  jp                              !<
3024    INTEGER(iwp) ::  kp                              !<
3025!
3026!-- First open the output unit.
3027    IF ( myid_char == '' )  THEN
3028       OPEN ( 90, FILE='PARTICLE_RESTART_DATA_OUT'//myid_char, &
3029                  FORM='UNFORMATTED')
3030    ELSE
3031       IF ( myid == 0 )  CALL local_system( 'mkdir PARTICLE_RESTART_DATA_OUT' )
3032#if defined( __parallel )
3033!
3034!--    Set a barrier in order to allow that thereafter all other processors
3035!--    in the directory created by PE0 can open their file
3036       CALL MPI_BARRIER( comm2d, ierr )
3037#endif
3038       OPEN ( 90, FILE='PARTICLE_RESTART_DATA_OUT/'//myid_char, &
3039                  FORM='UNFORMATTED' )
3040    ENDIF
3041
3042!
3043!-- Write the version number of the binary format.
3044!-- Attention: After changes to the following output commands the version
3045!-- ---------  number of the variable particle_binary_version must be
3046!--            changed! Also, the version number and the list of arrays
3047!--            to be read in lpm_read_restart_file must be adjusted
3048!--            accordingly.
3049    particle_binary_version = '4.0'
3050    WRITE ( 90 )  particle_binary_version
3051
3052!
3053!-- Write some particle parameters, the size of the particle arrays
3054    WRITE ( 90 )  bc_par_b, bc_par_lr, bc_par_ns, bc_par_t,                    &
3055                  last_particle_release_time, number_of_particle_groups,       &
3056                  particle_groups, time_write_particle_data
3057
3058    WRITE ( 90 )  prt_count
3059         
3060    DO  ip = nxl, nxr
3061       DO  jp = nys, nyn
3062          DO  kp = nzb+1, nzt
3063             number_of_particles = prt_count(kp,jp,ip)
3064             particles => grid_particles(kp,jp,ip)%particles(1:number_of_particles)
3065             IF ( number_of_particles <= 0 )  CYCLE
3066             WRITE ( 90 )  particles
3067          ENDDO
3068       ENDDO
3069    ENDDO
3070
3071    CLOSE ( 90 )
3072
3073#if defined( __parallel )
3074       CALL MPI_BARRIER( comm2d, ierr )
3075#endif
3076
3077    CALL wrd_write_string( 'iran' ) 
3078    WRITE ( 14 )  iran, iran_part 
3079
3080
3081 END SUBROUTINE lpm_wrd_local
3082 
3083 
3084!------------------------------------------------------------------------------!
3085! Description:
3086! ------------
3087!> This routine writes the respective restart data for the lpm.
3088!------------------------------------------------------------------------------!
3089 SUBROUTINE lpm_wrd_global
3090 
3091    CALL wrd_write_string( 'curvature_solution_effects' ) 
3092    WRITE ( 14 )  curvature_solution_effects
3093
3094 END SUBROUTINE lpm_wrd_global
3095 
3096
3097!------------------------------------------------------------------------------!
3098! Description:
3099! ------------
3100!> This routine writes the respective restart data for the lpm.
3101!------------------------------------------------------------------------------!
3102 SUBROUTINE lpm_rrd_global( found )
3103 
3104    USE control_parameters,                            &
3105        ONLY: length, restart_string
3106
3107    LOGICAL, INTENT(OUT)  ::  found
3108
3109    found = .TRUE.
3110
3111    SELECT CASE ( restart_string(1:length) )
3112
3113       CASE ( 'curvature_solution_effects' )
3114          READ ( 13 )  curvature_solution_effects
3115         
3116!          CASE ( 'global_paramter' )
3117!             READ ( 13 )  global_parameter
3118!          CASE ( 'global_array' )
3119!             IF ( .NOT. ALLOCATED( global_array ) )  ALLOCATE( global_array(1:10) )
3120!             READ ( 13 )  global_array
3121
3122       CASE DEFAULT
3123
3124          found = .FALSE.
3125
3126    END SELECT
3127   
3128 END SUBROUTINE lpm_rrd_global
3129
3130
3131!------------------------------------------------------------------------------!
3132! Description:
3133! ------------
3134!> This is a submodule of the lagrangian particle model. It contains all
3135!> dynamic processes of the lpm. This includes the advection (resolved and sub-
3136!> grid scale) as well as the boundary conditions of particles. As a next step
3137!> this submodule should be excluded as an own file.
3138!------------------------------------------------------------------------------!   
3139 SUBROUTINE lpm_advec (ip,jp,kp)
3140
3141    LOGICAL ::  subbox_at_wall !< flag to see if the current subgridbox is adjacent to a wall
3142
3143    INTEGER(iwp) ::  i                           !< index variable along x
3144    INTEGER(iwp) ::  ip                          !< index variable along x
3145    INTEGER(iwp) ::  j                           !< index variable along y
3146    INTEGER(iwp) ::  jp                          !< index variable along y
3147    INTEGER(iwp) ::  k                           !< index variable along z
3148    INTEGER(iwp) ::  k_wall                      !< vertical index of topography top
3149    INTEGER(iwp) ::  kp                          !< index variable along z
3150    INTEGER(iwp) ::  kw                          !< index variable along z
3151    INTEGER(iwp) ::  n                           !< loop variable over all particles in a grid box
3152    INTEGER(iwp) ::  nb                          !< block number particles are sorted in
3153    INTEGER(iwp) ::  surf_start                  !< Index on surface data-type for current grid box
3154
3155    INTEGER(iwp), DIMENSION(0:7) ::  start_index !< start particle index for current block
3156    INTEGER(iwp), DIMENSION(0:7) ::  end_index   !< start particle index for current block
3157
3158    REAL(wp) ::  aa                 !< dummy argument for horizontal particle interpolation
3159    REAL(wp) ::  bb                 !< dummy argument for horizontal particle interpolation
3160    REAL(wp) ::  cc                 !< dummy argument for horizontal particle interpolation
3161    REAL(wp) ::  d_z_p_z0           !< inverse of interpolation length for logarithmic interpolation
3162    REAL(wp) ::  dd                 !< dummy argument for horizontal particle interpolation
3163    REAL(wp) ::  de_dx_int_l        !< x/y-interpolated TKE gradient (x) at particle position at lower vertical level
3164    REAL(wp) ::  de_dx_int_u        !< x/y-interpolated TKE gradient (x) at particle position at upper vertical level
3165    REAL(wp) ::  de_dy_int_l        !< x/y-interpolated TKE gradient (y) at particle position at lower vertical level
3166    REAL(wp) ::  de_dy_int_u        !< x/y-interpolated TKE gradient (y) at particle position at upper vertical level
3167    REAL(wp) ::  de_dt              !< temporal derivative of TKE experienced by the particle
3168    REAL(wp) ::  de_dt_min          !< lower level for temporal TKE derivative
3169    REAL(wp) ::  de_dz_int_l        !< x/y-interpolated TKE gradient (z) at particle position at lower vertical level
3170    REAL(wp) ::  de_dz_int_u        !< x/y-interpolated TKE gradient (z) at particle position at upper vertical level
3171    REAL(wp) ::  diameter           !< diamter of droplet
3172    REAL(wp) ::  diss_int_l         !< x/y-interpolated dissipation at particle position at lower vertical level
3173    REAL(wp) ::  diss_int_u         !< x/y-interpolated dissipation at particle position at upper vertical level
3174    REAL(wp) ::  dt_particle_m      !< previous particle time step
3175    REAL(wp) ::  dz_temp            !< dummy for the vertical grid spacing
3176    REAL(wp) ::  e_int_l            !< x/y-interpolated TKE at particle position at lower vertical level
3177    REAL(wp) ::  e_int_u            !< x/y-interpolated TKE at particle position at upper vertical level
3178    REAL(wp) ::  e_mean_int         !< horizontal mean TKE at particle height
3179    REAL(wp) ::  exp_arg            !< argument in the exponent - particle radius
3180    REAL(wp) ::  exp_term           !< exponent term
3181    REAL(wp) ::  gg                 !< dummy argument for horizontal particle interpolation
3182    REAL(wp) ::  height_p           !< dummy argument for logarithmic interpolation
3183    REAL(wp) ::  log_z_z0_int       !< logarithmus used for surface_layer interpolation
3184    REAL(wp) ::  random_gauss       !< Gaussian-distributed random number used for SGS particle advection
3185    REAL(wp) ::  RL                 !< Lagrangian autocorrelation coefficient
3186    REAL(wp) ::  rg1                !< Gaussian distributed random number
3187    REAL(wp) ::  rg2                !< Gaussian distributed random number
3188    REAL(wp) ::  rg3                !< Gaussian distributed random number
3189    REAL(wp) ::  sigma              !< velocity standard deviation
3190    REAL(wp) ::  u_int_l            !< x/y-interpolated u-component at particle position at lower vertical level
3191    REAL(wp) ::  u_int_u            !< x/y-interpolated u-component at particle position at upper vertical level
3192    REAL(wp) ::  us_int             !< friction velocity at particle grid box
3193    REAL(wp) ::  usws_int           !< surface momentum flux (u component) at particle grid box
3194    REAL(wp) ::  v_int_l            !< x/y-interpolated v-component at particle position at lower vertical level
3195    REAL(wp) ::  v_int_u            !< x/y-interpolated v-component at particle position at upper vertical level
3196    REAL(wp) ::  vsws_int           !< surface momentum flux (u component) at particle grid box
3197    REAL(wp) ::  vv_int             !< dummy to compute interpolated mean SGS TKE, used to scale SGS advection
3198    REAL(wp) ::  w_int_l            !< x/y-interpolated w-component at particle position at lower vertical level
3199    REAL(wp) ::  w_int_u            !< x/y-interpolated w-component at particle position at upper vertical level
3200    REAL(wp) ::  w_s                !< terminal velocity of droplets
3201    REAL(wp) ::  x                  !< dummy argument for horizontal particle interpolation
3202    REAL(wp) ::  y                  !< dummy argument for horizontal particle interpolation
3203    REAL(wp) ::  z_p                !< surface layer height (0.5 dz)
3204
3205    REAL(wp), PARAMETER ::  a_rog = 9.65_wp      !< parameter for fall velocity
3206    REAL(wp), PARAMETER ::  b_rog = 10.43_wp     !< parameter for fall velocity
3207    REAL(wp), PARAMETER ::  c_rog = 0.6_wp       !< parameter for fall velocity
3208    REAL(wp), PARAMETER ::  k_cap_rog = 4.0_wp   !< parameter for fall velocity
3209    REAL(wp), PARAMETER ::  k_low_rog = 12.0_wp  !< parameter for fall velocity
3210    REAL(wp), PARAMETER ::  d0_rog = 0.745_wp    !< separation diameter
3211
3212    REAL(wp), DIMENSION(number_of_particles) ::  term_1_2       !< flag to communicate whether a particle is near topography or not
3213    REAL(wp), DIMENSION(number_of_particles) ::  dens_ratio     !< ratio between the density of the fluid and the density of the particles
3214    REAL(wp), DIMENSION(number_of_particles) ::  de_dx_int      !< horizontal TKE gradient along x at particle position
3215    REAL(wp), DIMENSION(number_of_particles) ::  de_dy_int      !< horizontal TKE gradient along y at particle position
3216    REAL(wp), DIMENSION(number_of_particles) ::  de_dz_int      !< horizontal TKE gradient along z at particle position
3217    REAL(wp), DIMENSION(number_of_particles) ::  diss_int       !< dissipation at particle position
3218    REAL(wp), DIMENSION(number_of_particles) ::  dt_gap         !< remaining time until particle time integration reaches LES time
3219    REAL(wp), DIMENSION(number_of_particles) ::  dt_particle    !< particle time step
3220    REAL(wp), DIMENSION(number_of_particles) ::  e_int          !< TKE at particle position
3221    REAL(wp), DIMENSION(number_of_particles) ::  fs_int         !< weighting factor for subgrid-scale particle speed
3222    REAL(wp), DIMENSION(number_of_particles) ::  lagr_timescale !< Lagrangian timescale
3223    REAL(wp), DIMENSION(number_of_particles) ::  rvar1_temp     !< SGS particle velocity - u-component
3224    REAL(wp), DIMENSION(number_of_particles) ::  rvar2_temp     !< SGS particle velocity - v-component
3225    REAL(wp), DIMENSION(number_of_particles) ::  rvar3_temp     !< SGS particle velocity - w-component
3226    REAL(wp), DIMENSION(number_of_particles) ::  u_int          !< u-component of particle speed
3227    REAL(wp), DIMENSION(number_of_particles) ::  v_int          !< v-component of particle speed
3228    REAL(wp), DIMENSION(number_of_particles) ::  w_int          !< w-component of particle speed
3229    REAL(wp), DIMENSION(number_of_particles) ::  xv             !< x-position
3230    REAL(wp), DIMENSION(number_of_particles) ::  yv             !< y-position
3231    REAL(wp), DIMENSION(number_of_particles) ::  zv             !< z-position
3232
3233    REAL(wp), DIMENSION(number_of_particles, 3) ::  rg !< vector of Gaussian distributed random numbers
3234
3235    CALL cpu_log( log_point_s(44), 'lpm_advec', 'continue' )
3236
3237!
3238!-- Determine height of Prandtl layer and distance between Prandtl-layer
3239!-- height and horizontal mean roughness height, which are required for
3240!-- vertical logarithmic interpolation of horizontal particle speeds
3241!-- (for particles below first vertical grid level).
3242    z_p      = zu(nzb+1) - zw(nzb)
3243    d_z_p_z0 = 1.0_wp / ( z_p - z0_av_global )
3244
3245    start_index = grid_particles(kp,jp,ip)%start_index
3246    end_index   = grid_particles(kp,jp,ip)%end_index
3247
3248    xv = particles(1:number_of_particles)%x
3249    yv = particles(1:number_of_particles)%y
3250    zv = particles(1:number_of_particles)%z
3251
3252    DO  nb = 0, 7
3253!
3254!--    Interpolate u velocity-component       
3255       i = ip
3256       j = jp + block_offset(nb)%j_off
3257       k = kp + block_offset(nb)%k_off
3258
3259       DO  n = start_index(nb), end_index(nb)
3260!
3261!--       Interpolation of the u velocity component onto particle position. 
3262!--       Particles are interpolation bi-linearly in the horizontal and a
3263!--       linearly in the vertical. An exception is made for particles below
3264!--       the first vertical grid level in case of a prandtl layer. In this
3265!--       case the horizontal particle velocity components are determined using
3266!--       Monin-Obukhov relations (if branch).
3267!--       First, check if particle is located below first vertical grid level
3268!--       above topography (Prandtl-layer height)
3269!--       Determine vertical index of topography top
3270          k_wall = get_topography_top_index_ji( jp, ip, 's' )
3271
3272          IF ( constant_flux_layer  .AND.  zv(n) - zw(k_wall) < z_p )  THEN
3273!
3274!--          Resolved-scale horizontal particle velocity is zero below z0.
3275             IF ( zv(n) - zw(k_wall) < z0_av_global )  THEN
3276                u_int(n) = 0.0_wp
3277             ELSE
3278!
3279!--             Determine the sublayer. Further used as index.
3280                height_p = ( zv(n) - zw(k_wall) - z0_av_global ) &
3281                                     * REAL( number_of_sublayers, KIND=wp )    &
3282                                     * d_z_p_z0 
3283!
3284!--             Calculate LOG(z/z0) for exact particle height. Therefore,   
3285!--             interpolate linearly between precalculated logarithm.
3286                log_z_z0_int = log_z_z0(INT(height_p))                         &
3287                                 + ( height_p - INT(height_p) )                &
3288                                 * ( log_z_z0(INT(height_p)+1)                 &
3289                                      - log_z_z0(INT(height_p))                &
3290                                   ) 
3291!
3292!--             Get friction velocity and momentum flux from new surface data
3293!--             types.
3294                IF ( surf_def_h(0)%start_index(jp,ip) <=                   &
3295                     surf_def_h(0)%end_index(jp,ip) )  THEN
3296                   surf_start = surf_def_h(0)%start_index(jp,ip)
3297!--                Limit friction velocity. In narrow canyons or holes the
3298!--                friction velocity can become very small, resulting in a too
3299!--                large particle speed.
3300                   us_int    = MAX( surf_def_h(0)%us(surf_start), 0.01_wp ) 
3301                   usws_int  = surf_def_h(0)%usws(surf_start)
3302                ELSEIF ( surf_lsm_h%start_index(jp,ip) <=                  &
3303                         surf_lsm_h%end_index(jp,ip) )  THEN
3304                   surf_start = surf_lsm_h%start_index(jp,ip)
3305                   us_int    = MAX( surf_lsm_h%us(surf_start), 0.01_wp ) 
3306                   usws_int  = surf_lsm_h%usws(surf_start)
3307                ELSEIF ( surf_usm_h%start_index(jp,ip) <=                  &
3308                         surf_usm_h%end_index(jp,ip) )  THEN
3309                   surf_start = surf_usm_h%start_index(jp,ip)
3310                   us_int    = MAX( surf_usm_h%us(surf_start), 0.01_wp ) 
3311                   usws_int  = surf_usm_h%usws(surf_start)
3312                ENDIF
3313
3314!
3315!--             Neutral solution is applied for all situations, e.g. also for
3316!--             unstable and stable situations. Even though this is not exact
3317!--             this saves a lot of CPU time since several calls of intrinsic
3318!--             FORTRAN procedures (LOG, ATAN) are avoided, This is justified
3319!--             as sensitivity studies revealed no significant effect of
3320!--             using the neutral solution also for un/stable situations.
3321                u_int(n) = -usws_int / ( us_int * kappa + 1E-10_wp )           & 
3322                            * log_z_z0_int - u_gtrans
3323               
3324             ENDIF
3325!
3326!--       Particle above the first grid level. Bi-linear interpolation in the
3327!--       horizontal and linear interpolation in the vertical direction.
3328          ELSE
3329
3330             x  = xv(n) - i * dx
3331             y  = yv(n) + ( 0.5_wp - j ) * dy
3332             aa = x**2          + y**2
3333             bb = ( dx - x )**2 + y**2
3334             cc = x**2          + ( dy - y )**2
3335             dd = ( dx - x )**2 + ( dy - y )**2
3336             gg = aa + bb + cc + dd
3337
3338             u_int_l = ( ( gg - aa ) * u(k,j,i)   + ( gg - bb ) * u(k,j,i+1)   &
3339                         + ( gg - cc ) * u(k,j+1,i) + ( gg - dd ) *            &
3340                         u(k,j+1,i+1) ) / ( 3.0_wp * gg ) - u_gtrans
3341
3342             IF ( k == nzt )  THEN
3343                u_int(n) = u_int_l
3344             ELSE
3345                u_int_u = ( ( gg-aa ) * u(k+1,j,i) + ( gg-bb ) * u(k+1,j,i+1)  &
3346                            + ( gg-cc ) * u(k+1,j+1,i) + ( gg-dd ) *           &
3347                            u(k+1,j+1,i+1) ) / ( 3.0_wp * gg ) - u_gtrans
3348                u_int(n) = u_int_l + ( zv(n) - zu(k) ) / dzw(k+1) *               &
3349                           ( u_int_u - u_int_l )
3350             ENDIF
3351
3352          ENDIF
3353
3354       ENDDO
3355!
3356!--    Same procedure for interpolation of the v velocity-component
3357       i = ip + block_offset(nb)%i_off
3358       j = jp
3359       k = kp + block_offset(nb)%k_off
3360
3361       DO  n = start_index(nb), end_index(nb)
3362
3363!
3364!--       Determine vertical index of topography top
3365          k_wall = get_topography_top_index_ji( jp,ip, 's' )
3366
3367          IF ( constant_flux_layer  .AND.  zv(n) - zw(k_wall) < z_p )  THEN
3368             IF ( zv(n) - zw(k_wall) < z0_av_global )  THEN
3369!
3370!--             Resolved-scale horizontal particle velocity is zero below z0.
3371                v_int(n) = 0.0_wp
3372             ELSE       
3373!
3374!--             Determine the sublayer. Further used as index. Please note,
3375!--             logarithmus can not be reused from above, as in in case of
3376!--             topography particle on u-grid can be above surface-layer height,
3377!--             whereas it can be below on v-grid.
3378                height_p = ( zv(n) - zw(k_wall) - z0_av_global ) &
3379                                  * REAL( number_of_sublayers, KIND=wp )       &
3380                                  * d_z_p_z0 
3381!
3382!--             Calculate LOG(z/z0) for exact particle height. Therefore,   
3383!--             interpolate linearly between precalculated logarithm.
3384                log_z_z0_int = log_z_z0(INT(height_p))                         &
3385                                 + ( height_p - INT(height_p) )                &
3386                                 * ( log_z_z0(INT(height_p)+1)                 &
3387                                      - log_z_z0(INT(height_p))                &
3388                                   ) 
3389!
3390!--             Get friction velocity and momentum flux from new surface data
3391!--             types.
3392                IF ( surf_def_h(0)%start_index(jp,ip) <=                   &
3393                     surf_def_h(0)%end_index(jp,ip) )  THEN
3394                   surf_start = surf_def_h(0)%start_index(jp,ip)
3395!--                Limit friction velocity. In narrow canyons or holes the
3396!--                friction velocity can become very small, resulting in a too
3397!--                large particle speed.
3398                   us_int    = MAX( surf_def_h(0)%us(surf_start), 0.01_wp ) 
3399                   vsws_int  = surf_def_h(0)%vsws(surf_start)
3400                ELSEIF ( surf_lsm_h%start_index(jp,ip) <=                  &
3401                         surf_lsm_h%end_index(jp,ip) )  THEN
3402                   surf_start = surf_lsm_h%start_index(jp,ip)
3403                   us_int    = MAX( surf_lsm_h%us(surf_start), 0.01_wp ) 
3404                   vsws_int  = surf_lsm_h%vsws(surf_start)
3405                ELSEIF ( surf_usm_h%start_index(jp,ip) <=                  &
3406                         surf_usm_h%end_index(jp,ip) )  THEN
3407                   surf_start = surf_usm_h%start_index(jp,ip)
3408                   us_int    = MAX( surf_usm_h%us(surf_start), 0.01_wp ) 
3409                   vsws_int  = surf_usm_h%vsws(surf_start)
3410                ENDIF 
3411!
3412!--             Neutral solution is applied for all situations, e.g. also for
3413!--             unstable and stable situations. Even though this is not exact
3414!--             this saves a lot of CPU time since several calls of intrinsic
3415!--             FORTRAN procedures (LOG, ATAN) are avoided, This is justified
3416!--             as sensitivity studies revealed no significant effect of
3417!--             using the neutral solution also for un/stable situations.
3418                v_int(n) = -vsws_int / ( us_int * kappa + 1E-10_wp )           &
3419                         * log_z_z0_int - v_gtrans
3420
3421             ENDIF
3422
3423          ELSE
3424             x  = xv(n) + ( 0.5_wp - i ) * dx
3425             y  = yv(n) - j * dy
3426             aa = x**2          + y**2
3427             bb = ( dx - x )**2 + y**2
3428             cc = x**2          + ( dy - y )**2
3429             dd = ( dx - x )**2 + ( dy - y )**2
3430             gg = aa + bb + cc + dd
3431
3432             v_int_l = ( ( gg - aa ) * v(k,j,i)   + ( gg - bb ) * v(k,j,i+1)   &
3433                       + ( gg - cc ) * v(k,j+1,i) + ( gg - dd ) * v(k,j+1,i+1) &
3434                       ) / ( 3.0_wp * gg ) - v_gtrans
3435
3436             IF ( k == nzt )  THEN
3437                v_int(n) = v_int_l
3438             ELSE
3439                v_int_u = ( ( gg-aa ) * v(k+1,j,i)   + ( gg-bb ) * v(k+1,j,i+1)   &
3440                          + ( gg-cc ) * v(k+1,j+1,i) + ( gg-dd ) * v(k+1,j+1,i+1) &
3441                          ) / ( 3.0_wp * gg ) - v_gtrans
3442                v_int(n) = v_int_l + ( zv(n) - zu(k) ) / dzw(k+1) *               &
3443                                  ( v_int_u - v_int_l )
3444             ENDIF
3445
3446          ENDIF
3447
3448       ENDDO
3449!
3450!--    Same procedure for interpolation of the w velocity-component
3451       i = ip + block_offset(nb)%i_off
3452       j = jp + block_offset(nb)%j_off
3453       k = kp - 1
3454
3455       DO  n = start_index(nb), end_index(nb)
3456
3457          IF ( vertical_particle_advection(particles(n)%group) )  THEN
3458
3459             x  = xv(n) + ( 0.5_wp - i ) * dx
3460             y  = yv(n) + ( 0.5_wp - j ) * dy
3461             aa = x**2          + y**2
3462             bb = ( dx - x )**2 + y**2
3463             cc = x**2          + ( dy - y )**2
3464             dd = ( dx - x )**2 + ( dy - y )**2
3465             gg = aa + bb + cc + dd
3466
3467             w_int_l = ( ( gg - aa ) * w(k,j,i)   + ( gg - bb ) * w(k,j,i+1)   &
3468                       + ( gg - cc ) * w(k,j+1,i) + ( gg - dd ) * w(k,j+1,i+1) &
3469                       ) / ( 3.0_wp * gg )
3470
3471             IF ( k == nzt )  THEN
3472                w_int(n) = w_int_l
3473             ELSE
3474                w_int_u = ( ( gg-aa ) * w(k+1,j,i)   + &
3475                            ( gg-bb ) * w(k+1,j,i+1) + &
3476                            ( gg-cc ) * w(k+1,j+1,i) + &
3477                            ( gg-dd ) * w(k+1,j+1,i+1) &
3478                          ) / ( 3.0_wp * gg )
3479                w_int(n) = w_int_l + ( zv(n) - zw(k) ) / dzw(k+1) *               &
3480                           ( w_int_u - w_int_l )
3481             ENDIF
3482
3483          ELSE
3484
3485             w_int(n) = 0.0_wp
3486
3487          ENDIF
3488
3489       ENDDO
3490
3491    ENDDO
3492
3493!-- Interpolate and calculate quantities needed for calculating the SGS
3494!-- velocities
3495    IF ( use_sgs_for_particles  .AND.  .NOT. cloud_droplets )  THEN
3496   
3497       DO  nb = 0,7
3498         
3499          subbox_at_wall = .FALSE.         
3500!
3501!--       In case of topography check if subbox is adjacent to a wall
3502          IF ( .NOT. topography == 'flat' ) THEN
3503             i = ip + MERGE( -1_iwp , 1_iwp, BTEST( nb, 2 ) )
3504             j = jp + MERGE( -1_iwp , 1_iwp, BTEST( nb, 1 ) )
3505             k = kp + MERGE( -1_iwp , 1_iwp, BTEST( nb, 0 ) )
3506             IF ( .NOT. BTEST(wall_flags_0(k,  jp, ip), 0) .OR.                &
3507                  .NOT. BTEST(wall_flags_0(kp, j,  ip), 0) .OR.                &
3508                  .NOT. BTEST(wall_flags_0(kp, jp, i ), 0) )                   &
3509             THEN
3510                subbox_at_wall = .TRUE.
3511             ENDIF
3512          ENDIF
3513          IF ( subbox_at_wall ) THEN
3514             e_int(start_index(nb):end_index(nb))     = e(kp,jp,ip) 
3515             diss_int(start_index(nb):end_index(nb))  = diss(kp,jp,ip)
3516             de_dx_int(start_index(nb):end_index(nb)) = de_dx(kp,jp,ip)
3517             de_dy_int(start_index(nb):end_index(nb)) = de_dy(kp,jp,ip)
3518             de_dz_int(start_index(nb):end_index(nb)) = de_dz(kp,jp,ip)
3519!
3520!--          Set flag for stochastic equation.
3521             term_1_2(start_index(nb):end_index(nb)) = 0.0_wp             
3522          ELSE
3523             i = ip + block_offset(nb)%i_off
3524             j = jp + block_offset(nb)%j_off
3525             k = kp + block_offset(nb)%k_off
3526
3527             DO  n = start_index(nb), end_index(nb)
3528!
3529!--             Interpolate TKE
3530                x  = xv(n) + ( 0.5_wp - i ) * dx
3531                y  = yv(n) + ( 0.5_wp - j ) * dy
3532                aa = x**2          + y**2
3533                bb = ( dx - x )**2 + y**2
3534                cc = x**2          + ( dy - y )**2
3535                dd = ( dx - x )**2 + ( dy - y )**2
3536                gg = aa + bb + cc + dd
3537
3538                e_int_l = ( ( gg-aa ) * e(k,j,i)   + ( gg-bb ) * e(k,j,i+1)   &
3539                          + ( gg-cc ) * e(k,j+1,i) + ( gg-dd ) * e(k,j+1,i+1) &
3540                          ) / ( 3.0_wp * gg )
3541
3542                IF ( k+1 == nzt+1 )  THEN
3543                   e_int(n) = e_int_l
3544                ELSE
3545                   e_int_u = ( ( gg - aa ) * e(k+1,j,i)   + &
3546                               ( gg - bb ) * e(k+1,j,i+1) + &
3547                               ( gg - cc ) * e(k+1,j+1,i) + &
3548                               ( gg - dd ) * e(k+1,j+1,i+1) &
3549                            ) / ( 3.0_wp * gg )
3550                   e_int(n) = e_int_l + ( zv(n) - zu(k) ) / dzw(k+1) *            &
3551                                     ( e_int_u - e_int_l )
3552                ENDIF
3553!
3554!--             Needed to avoid NaN particle velocities (this might not be
3555!--             required any more)
3556                IF ( e_int(n) <= 0.0_wp )  THEN
3557                   e_int(n) = 1.0E-20_wp
3558                ENDIF
3559!
3560!--             Interpolate the TKE gradient along x (adopt incides i,j,k and
3561!--             all position variables from above (TKE))
3562                de_dx_int_l = ( ( gg - aa ) * de_dx(k,j,i)   + &
3563                                ( gg - bb ) * de_dx(k,j,i+1) + &
3564                                ( gg - cc ) * de_dx(k,j+1,i) + &
3565                                ( gg - dd ) * de_dx(k,j+1,i+1) &
3566                               ) / ( 3.0_wp * gg )
3567
3568                IF ( ( k+1 == nzt+1 )  .OR.  ( k == nzb ) )  THEN
3569                   de_dx_int(n) = de_dx_int_l
3570                ELSE
3571                   de_dx_int_u = ( ( gg - aa ) * de_dx(k+1,j,i)   + &
3572                                   ( gg - bb ) * de_dx(k+1,j,i+1) + &
3573                                   ( gg - cc ) * de_dx(k+1,j+1,i) + &
3574                                   ( gg - dd ) * de_dx(k+1,j+1,i+1) &
3575                                  ) / ( 3.0_wp * gg )
3576                   de_dx_int(n) = de_dx_int_l + ( zv(n) - zu(k) ) / dzw(k+1) *    &
3577                                              ( de_dx_int_u - de_dx_int_l )
3578                ENDIF
3579!
3580!--             Interpolate the TKE gradient along y
3581                de_dy_int_l = ( ( gg - aa ) * de_dy(k,j,i)   + &
3582                                ( gg - bb ) * de_dy(k,j,i+1) + &
3583                                ( gg - cc ) * de_dy(k,j+1,i) + &
3584                                ( gg - dd ) * de_dy(k,j+1,i+1) &
3585                               ) / ( 3.0_wp * gg )
3586                IF ( ( k+1 == nzt+1 )  .OR.  ( k == nzb ) )  THEN
3587                   de_dy_int(n) = de_dy_int_l
3588                ELSE
3589                   de_dy_int_u = ( ( gg - aa ) * de_dy(k+1,j,i)   + &
3590                                   ( gg - bb ) * de_dy(k+1,j,i+1) + &
3591                                   ( gg - cc ) * de_dy(k+1,j+1,i) + &
3592                                   ( gg - dd ) * de_dy(k+1,j+1,i+1) &
3593                                  ) / ( 3.0_wp * gg )
3594                      de_dy_int(n) = de_dy_int_l + ( zv(n) - zu(k) ) / dzw(k+1) * &
3595                                                 ( de_dy_int_u - de_dy_int_l )
3596                ENDIF
3597
3598!
3599!--             Interpolate the TKE gradient along z
3600                IF ( zv(n) < 0.5_wp * dz(1) )  THEN
3601                   de_dz_int(n) = 0.0_wp
3602                ELSE
3603                   de_dz_int_l = ( ( gg - aa ) * de_dz(k,j,i)   + &
3604                                   ( gg - bb ) * de_dz(k,j,i+1) + &
3605                                   ( gg - cc ) * de_dz(k,j+1,i) + &
3606                                   ( gg - dd ) * de_dz(k,j+1,i+1) &
3607                                  ) / ( 3.0_wp * gg )
3608
3609                   IF ( ( k+1 == nzt+1 )  .OR.  ( k == nzb ) )  THEN
3610                      de_dz_int(n) = de_dz_int_l
3611                   ELSE
3612                      de_dz_int_u = ( ( gg - aa ) * de_dz(k+1,j,i)   + &
3613                                      ( gg - bb ) * de_dz(k+1,j,i+1) + &
3614                                      ( gg - cc ) * de_dz(k+1,j+1,i) + &
3615                                      ( gg - dd ) * de_dz(k+1,j+1,i+1) &
3616                                     ) / ( 3.0_wp * gg )
3617                      de_dz_int(n) = de_dz_int_l + ( zv(n) - zu(k) ) / dzw(k+1) * &
3618                                                 ( de_dz_int_u - de_dz_int_l )
3619                   ENDIF
3620                ENDIF
3621
3622!
3623!--             Interpolate the dissipation of TKE
3624                diss_int_l = ( ( gg - aa ) * diss(k,j,i)   + &
3625                               ( gg - bb ) * diss(k,j,i+1) + &
3626                               ( gg - cc ) * diss(k,j+1,i) + &
3627                               ( gg - dd ) * diss(k,j+1,i+1) &
3628                               ) / ( 3.0_wp * gg )
3629
3630                IF ( k == nzt )  THEN
3631                   diss_int(n) = diss_int_l
3632                ELSE
3633                   diss_int_u = ( ( gg - aa ) * diss(k+1,j,i)   + &
3634                                  ( gg - bb ) * diss(k+1,j,i+1) + &
3635                                  ( gg - cc ) * diss(k+1,j+1,i) + &
3636                                  ( gg - dd ) * diss(k+1,j+1,i+1) &
3637                                 ) / ( 3.0_wp * gg )
3638                   diss_int(n) = diss_int_l + ( zv(n) - zu(k) ) / dzw(k+1) *      &
3639                                            ( diss_int_u - diss_int_l )
3640                ENDIF
3641
3642!
3643!--             Set flag for stochastic equation.
3644                term_1_2(n) = 1.0_wp
3645             ENDDO
3646          ENDIF
3647       ENDDO
3648
3649       DO nb = 0,7
3650          i = ip + block_offset(nb)%i_off
3651          j = jp + block_offset(nb)%j_off
3652          k = kp + block_offset(nb)%k_off
3653
3654          DO  n = start_index(nb), end_index(nb)
3655!
3656!--          Vertical interpolation of the horizontally averaged SGS TKE and
3657!--          resolved-scale velocity variances and use the interpolated values
3658!--          to calculate the coefficient fs, which is a measure of the ratio
3659!--          of the subgrid-scale turbulent kinetic energy to the total amount
3660!--          of turbulent kinetic energy.
3661             IF ( k == 0 )  THEN
3662                e_mean_int = hom(0,1,8,0)
3663             ELSE
3664                e_mean_int = hom(k,1,8,0) +                                    &
3665                                           ( hom(k+1,1,8,0) - hom(k,1,8,0) ) / &
3666                                           ( zu(k+1) - zu(k) ) *               &
3667                                           ( zv(n) - zu(k) )
3668             ENDIF
3669
3670             kw = kp - 1
3671
3672             IF ( k == 0 )  THEN
3673                aa  = hom(k+1,1,30,0)  * ( zv(n) / &
3674                                         ( 0.5_wp * ( zu(k+1) - zu(k) ) ) )
3675                bb  = hom(k+1,1,31,0)  * ( zv(n) / &
3676                                         ( 0.5_wp * ( zu(k+1) - zu(k) ) ) )
3677                cc  = hom(kw+1,1,32,0) * ( zv(n) / &
3678                                         ( 1.0_wp * ( zw(kw+1) - zw(kw) ) ) )
3679             ELSE
3680                aa  = hom(k,1,30,0) + ( hom(k+1,1,30,0) - hom(k,1,30,0) ) *    &
3681                           ( ( zv(n) - zu(k) ) / ( zu(k+1) - zu(k) ) )
3682                bb  = hom(k,1,31,0) + ( hom(k+1,1,31,0) - hom(k,1,31,0) ) *    &
3683                           ( ( zv(n) - zu(k) ) / ( zu(k+1) - zu(k) ) )
3684                cc  = hom(kw,1,32,0) + ( hom(kw+1,1,32,0)-hom(kw,1,32,0) ) *   &
3685                           ( ( zv(n) - zw(kw) ) / ( zw(kw+1)-zw(kw) ) )
3686             ENDIF
3687
3688             vv_int = ( 1.0_wp / 3.0_wp ) * ( aa + bb + cc )
3689!
3690!--          Needed to avoid NaN particle velocities. The value of 1.0 is just
3691!--          an educated guess for the given case.
3692             IF ( vv_int + ( 2.0_wp / 3.0_wp ) * e_mean_int == 0.0_wp )  THEN
3693                fs_int(n) = 1.0_wp
3694             ELSE
3695                fs_int(n) = ( 2.0_wp / 3.0_wp ) * e_mean_int /                 &
3696                            ( vv_int + ( 2.0_wp / 3.0_wp ) * e_mean_int )
3697             ENDIF
3698
3699          ENDDO
3700       ENDDO
3701
3702       DO  nb = 0, 7
3703          DO  n = start_index(nb), end_index(nb)
3704             rg(n,1) = random_gauss( iran_part, 5.0_wp )
3705             rg(n,2) = random_gauss( iran_part, 5.0_wp )
3706             rg(n,3) = random_gauss( iran_part, 5.0_wp )
3707          ENDDO
3708       ENDDO
3709
3710       DO  nb = 0, 7
3711          DO  n = start_index(nb), end_index(nb)
3712
3713!
3714!--          Calculate the Lagrangian timescale according to Weil et al. (2004).
3715             lagr_timescale(n) = ( 4.0_wp * e_int(n) + 1E-20_wp ) / &
3716                              ( 3.0_wp * fs_int(n) * c_0 * diss_int(n) + 1E-20_wp )
3717
3718!
3719!--          Calculate the next particle timestep. dt_gap is the time needed to
3720!--          complete the current LES timestep.
3721             dt_gap(n) = dt_3d - particles(n)%dt_sum
3722             dt_particle(n) = MIN( dt_3d, 0.025_wp * lagr_timescale(n), dt_gap(n) )
3723             particles(n)%aux1 = lagr_timescale(n)
3724             particles(n)%aux2 = dt_gap(n)
3725!
3726!--          The particle timestep should not be too small in order to prevent
3727!--          the number of particle timesteps of getting too large
3728             IF ( dt_particle(n) < dt_min_part  .AND.  dt_min_part < dt_gap(n) )  THEN
3729                dt_particle(n) = dt_min_part
3730             ENDIF
3731             rvar1_temp(n) = particles(n)%rvar1
3732             rvar2_temp(n) = particles(n)%rvar2
3733             rvar3_temp(n) = particles(n)%rvar3
3734!
3735!--          Calculate the SGS velocity components
3736             IF ( particles(n)%age == 0.0_wp )  THEN
3737!
3738!--             For new particles the SGS components are derived from the SGS
3739!--             TKE. Limit the Gaussian random number to the interval
3740!--             [-5.0*sigma, 5.0*sigma] in order to prevent the SGS velocities
3741!--             from becoming unrealistically large.
3742                rvar1_temp(n) = SQRT( 2.0_wp * sgs_wf_part * e_int(n)          &
3743                                          + 1E-20_wp ) * ( rg(n,1) - 1.0_wp )
3744                rvar2_temp(n) = SQRT( 2.0_wp * sgs_wf_part * e_int(n)          &
3745                                          + 1E-20_wp ) * ( rg(n,2) - 1.0_wp )
3746                rvar3_temp(n) = SQRT( 2.0_wp * sgs_wf_part * e_int(n)          &
3747                                          + 1E-20_wp ) * ( rg(n,3) - 1.0_wp )
3748
3749             ELSE
3750!
3751!--             Restriction of the size of the new timestep: compared to the
3752!--             previous timestep the increase must not exceed 200%. First,
3753!--             check if age > age_m, in order to prevent that particles get zero
3754!--             timestep.
3755                dt_particle_m = MERGE( dt_particle(n),                         &
3756                                       particles(n)%age - particles(n)%age_m,  &
3757                                       particles(n)%age - particles(n)%age_m < &
3758                                       1E-8_wp )
3759                IF ( dt_particle(n) > 2.0_wp * dt_particle_m )  THEN
3760                   dt_particle(n) = 2.0_wp * dt_particle_m
3761                ENDIF
3762
3763!--             For old particles the SGS components are correlated with the
3764!--             values from the previous timestep. Random numbers have also to
3765!--             be limited (see above).
3766!--             As negative values for the subgrid TKE are not allowed, the
3767!--             change of the subgrid TKE with time cannot be smaller than
3768!--             -e_int(n)/dt_particle. This value is used as a lower boundary
3769!--             value for the change of TKE
3770                de_dt_min = - e_int(n) / dt_particle(n)
3771
3772                de_dt = ( e_int(n) - particles(n)%e_m ) / dt_particle_m
3773
3774                IF ( de_dt < de_dt_min )  THEN
3775                   de_dt = de_dt_min
3776                ENDIF
3777
3778                CALL weil_stochastic_eq(rvar1_temp(n), fs_int(n), e_int(n),& 
3779                                        de_dx_int(n), de_dt, diss_int(n),       &
3780                                        dt_particle(n), rg(n,1), term_1_2(n) )
3781
3782                CALL weil_stochastic_eq(rvar2_temp(n), fs_int(n), e_int(n),& 
3783                                        de_dy_int(n), de_dt, diss_int(n),       &
3784                                        dt_particle(n), rg(n,2), term_1_2(n) )
3785
3786                CALL weil_stochastic_eq(rvar3_temp(n), fs_int(n), e_int(n),& 
3787                                        de_dz_int(n), de_dt, diss_int(n),       &
3788                                        dt_particle(n), rg(n,3), term_1_2(n) )
3789
3790             ENDIF
3791
3792          ENDDO
3793       ENDDO
3794!
3795!--    Check if the added SGS velocities result in a violation of the CFL-
3796!--    criterion. If yes choose a smaller timestep based on the new velocities
3797!--    and calculate SGS velocities again
3798       dz_temp = zw(kp)-zw(kp-1)
3799       
3800       DO  nb = 0, 7
3801          DO  n = start_index(nb), end_index(nb)
3802             IF ( .NOT. particles(n)%age == 0.0_wp .AND.                       &
3803                (ABS( u_int(n) + rvar1_temp(n) ) > (dx/dt_particle(n))  .OR.   &
3804                 ABS( v_int(n) + rvar2_temp(n) ) > (dy/dt_particle(n))  .OR.   &
3805                 ABS( w_int(n) + rvar3_temp(n) ) > (dz_temp/dt_particle(n)))) THEN
3806               
3807                dt_particle(n) = 0.9_wp * MIN(                                 &
3808                                 ( dx / ABS( u_int(n) + rvar1_temp(n) ) ),     &
3809                                 ( dy / ABS( v_int(n) + rvar2_temp(n) ) ),     &
3810                                 ( dz_temp / ABS( w_int(n) + rvar3_temp(n) ) ) )
3811
3812!
3813!--             Reset temporary SGS velocites to "current" ones
3814                rvar1_temp(n) = particles(n)%rvar1
3815                rvar2_temp(n) = particles(n)%rvar2
3816                rvar3_temp(n) = particles(n)%rvar3
3817               
3818                de_dt_min = - e_int(n) / dt_particle(n)
3819
3820                de_dt = ( e_int(n) - particles(n)%e_m ) / dt_particle_m
3821
3822                IF ( de_dt < de_dt_min )  THEN
3823                   de_dt = de_dt_min
3824                ENDIF
3825
3826                CALL weil_stochastic_eq(rvar1_temp(n), fs_int(n), e_int(n),& 
3827                                        de_dx_int(n), de_dt, diss_int(n),       &
3828                                        dt_particle(n), rg(n,1), term_1_2(n) )
3829
3830                CALL weil_stochastic_eq(rvar2_temp(n), fs_int(n), e_int(n),& 
3831                                        de_dy_int(n), de_dt, diss_int(n),       &
3832                                        dt_particle(n), rg(n,2), term_1_2(n) )
3833
3834                CALL weil_stochastic_eq(rvar3_temp(n), fs_int(n), e_int(n),& 
3835                                        de_dz_int(n), de_dt, diss_int(n),       &
3836                                        dt_particle(n), rg(n,3), term_1_2(n) )
3837             ENDIF                           
3838             
3839!
3840!--          Update particle velocites
3841             particles(n)%rvar1 = rvar1_temp(n)
3842             particles(n)%rvar2 = rvar2_temp(n)
3843             particles(n)%rvar3 = rvar3_temp(n)
3844             u_int(n) = u_int(n) + particles(n)%rvar1
3845             v_int(n) = v_int(n) + particles(n)%rvar2
3846             w_int(n) = w_int(n) + particles(n)%rvar3
3847!
3848!--          Store the SGS TKE of the current timelevel which is needed for
3849!--          for calculating the SGS particle velocities at the next timestep
3850             particles(n)%e_m = e_int(n)
3851          ENDDO
3852       ENDDO
3853       
3854    ELSE
3855!
3856!--    If no SGS velocities are used, only the particle timestep has to
3857!--    be set
3858       dt_particle = dt_3d
3859
3860    ENDIF
3861
3862    dens_ratio = particle_groups(particles(1:number_of_particles)%group)%density_ratio
3863
3864    IF ( ANY( dens_ratio == 0.0_wp ) )  THEN
3865       DO  nb = 0, 7
3866          DO  n = start_index(nb), end_index(nb)
3867
3868!
3869!--          Particle advection
3870             IF ( dens_ratio(n) == 0.0_wp )  THEN
3871!
3872!--             Pure passive transport (without particle inertia)
3873                particles(n)%x = xv(n) + u_int(n) * dt_particle(n)
3874                particles(n)%y = yv(n) + v_int(n) * dt_particle(n)
3875                particles(n)%z = zv(n) + w_int(n) * dt_particle(n)
3876
3877                particles(n)%speed_x = u_int(n)
3878                particles(n)%speed_y = v_int(n)
3879                particles(n)%speed_z = w_int(n)
3880
3881             ELSE
3882!
3883!--             Transport of particles with inertia
3884                particles(n)%x = particles(n)%x + particles(n)%speed_x * &
3885                                                  dt_particle(n)
3886                particles(n)%y = particles(n)%y + particles(n)%speed_y * &
3887                                                  dt_particle(n)
3888                particles(n)%z = particles(n)%z + particles(n)%speed_z * &
3889                                                  dt_particle(n)
3890
3891!
3892!--             Update of the particle velocity
3893                IF ( cloud_droplets )  THEN
3894!
3895!--                Terminal velocity is computed for vertical direction (Rogers et
3896!--                al., 1993, J. Appl. Meteorol.)
3897                   diameter = particles(n)%radius * 2000.0_wp !diameter in mm
3898                   IF ( diameter <= d0_rog )  THEN
3899                      w_s = k_cap_rog * diameter * ( 1.0_wp - EXP( -k_low_rog * diameter ) )
3900                   ELSE
3901                      w_s = a_rog - b_rog * EXP( -c_rog * diameter )
3902                   ENDIF
3903
3904!
3905!--                If selected, add random velocities following Soelch and Kaercher
3906!--                (2010, Q. J. R. Meteorol. Soc.)
3907                   IF ( use_sgs_for_particles )  THEN
3908                      lagr_timescale(n) = km(kp,jp,ip) / MAX( e(kp,jp,ip), 1.0E-20_wp )
3909                      RL             = EXP( -1.0_wp * dt_3d / MAX( lagr_timescale(n), &
3910                                             1.0E-20_wp ) )
3911                      sigma          = SQRT( e(kp,jp,ip) )
3912
3913                      rg1 = random_gauss( iran_part, 5.0_wp ) - 1.0_wp
3914                      rg2 = random_gauss( iran_part, 5.0_wp ) - 1.0_wp
3915                      rg3 = random_gauss( iran_part, 5.0_wp ) - 1.0_wp
3916
3917                      particles(n)%rvar1 = RL * particles(n)%rvar1 +              &
3918                                           SQRT( 1.0_wp - RL**2 ) * sigma * rg1
3919                      particles(n)%rvar2 = RL * particles(n)%rvar2 +              &
3920                                           SQRT( 1.0_wp - RL**2 ) * sigma * rg2
3921                      particles(n)%rvar3 = RL * particles(n)%rvar3 +              &
3922                                           SQRT( 1.0_wp - RL**2 ) * sigma * rg3
3923
3924                      particles(n)%speed_x = u_int(n) + particles(n)%rvar1
3925                      particles(n)%speed_y = v_int(n) + particles(n)%rvar2
3926                      particles(n)%speed_z = w_int(n) + particles(n)%rvar3 - w_s
3927                   ELSE
3928                      particles(n)%speed_x = u_int(n)
3929                      particles(n)%speed_y = v_int(n)
3930                      particles(n)%speed_z = w_int(n) - w_s
3931                   ENDIF
3932
3933                ELSE
3934
3935                   IF ( use_sgs_for_particles )  THEN
3936                      exp_arg  = particle_groups(particles(n)%group)%exp_arg
3937                      exp_term = EXP( -exp_arg * dt_particle(n) )
3938                   ELSE
3939                      exp_arg  = particle_groups(particles(n)%group)%exp_arg
3940                      exp_term = particle_groups(particles(n)%group)%exp_term
3941                   ENDIF
3942                   particles(n)%speed_x = particles(n)%speed_x * exp_term +         &
3943                                          u_int(n) * ( 1.0_wp - exp_term )
3944                   particles(n)%speed_y = particles(n)%speed_y * exp_term +         &
3945                                          v_int(n) * ( 1.0_wp - exp_term )
3946                   particles(n)%speed_z = particles(n)%speed_z * exp_term +         &
3947                                          ( w_int(n) - ( 1.0_wp - dens_ratio(n) ) * &
3948                                          g / exp_arg ) * ( 1.0_wp - exp_term )
3949                ENDIF
3950
3951             ENDIF
3952          ENDDO
3953       ENDDO
3954   
3955    ELSE
3956
3957       DO  nb = 0, 7
3958          DO  n = start_index(nb), end_index(nb)
3959!
3960!--          Transport of particles with inertia
3961             particles(n)%x = xv(n) + particles(n)%speed_x * dt_particle(n)
3962             particles(n)%y = yv(n) + particles(n)%speed_y * dt_particle(n)
3963             particles(n)%z = zv(n) + particles(n)%speed_z * dt_particle(n)
3964!
3965!--          Update of the particle velocity
3966             IF ( cloud_droplets )  THEN
3967!
3968!--             Terminal velocity is computed for vertical direction (Rogers et al.,
3969!--             1993, J. Appl. Meteorol.)
3970                diameter = particles(n)%radius * 2000.0_wp !diameter in mm
3971                IF ( diameter <= d0_rog )  THEN
3972                   w_s = k_cap_rog * diameter * ( 1.0_wp - EXP( -k_low_rog * diameter ) )
3973                ELSE
3974                   w_s = a_rog - b_rog * EXP( -c_rog * diameter )
3975                ENDIF
3976
3977!
3978!--             If selected, add random velocities following Soelch and Kaercher
3979!--             (2010, Q. J. R. Meteorol. Soc.)
3980                IF ( use_sgs_for_particles )  THEN
3981                    lagr_timescale(n) = km(kp,jp,ip) / MAX( e(kp,jp,ip), 1.0E-20_wp )
3982                     RL             = EXP( -1.0_wp * dt_3d / MAX( lagr_timescale(n), &
3983                                             1.0E-20_wp ) )
3984                    sigma          = SQRT( e(kp,jp,ip) )
3985
3986                    rg1 = random_gauss( iran_part, 5.0_wp ) - 1.0_wp
3987                    rg2 = random_gauss( iran_part, 5.0_wp ) - 1.0_wp
3988                    rg3 = random_gauss( iran_part, 5.0_wp ) - 1.0_wp
3989
3990                    particles(n)%rvar1 = RL * particles(n)%rvar1 +                &
3991                                         SQRT( 1.0_wp - RL**2 ) * sigma * rg1
3992                    particles(n)%rvar2 = RL * particles(n)%rvar2 +                &
3993                                         SQRT( 1.0_wp - RL**2 ) * sigma * rg2
3994                    particles(n)%rvar3 = RL * particles(n)%rvar3 +                &
3995                                         SQRT( 1.0_wp - RL**2 ) * sigma * rg3
3996
3997                    particles(n)%speed_x = u_int(n) + particles(n)%rvar1
3998                    particles(n)%speed_y = v_int(n) + particles(n)%rvar2
3999                    particles(n)%speed_z = w_int(n) + particles(n)%rvar3 - w_s
4000                ELSE
4001                    particles(n)%speed_x = u_int(n)
4002                    particles(n)%speed_y = v_int(n)
4003                    particles(n)%speed_z = w_int(n) - w_s
4004                ENDIF
4005
4006             ELSE
4007
4008                IF ( use_sgs_for_particles )  THEN
4009                   exp_arg  = particle_groups(particles(n)%group)%exp_arg
4010                   exp_term = EXP( -exp_arg * dt_particle(n) )
4011                ELSE
4012                   exp_arg  = particle_groups(particles(n)%group)%exp_arg
4013                   exp_term = particle_groups(particles(n)%group)%exp_term
4014                ENDIF
4015                particles(n)%speed_x = particles(n)%speed_x * exp_term +             &
4016                                       u_int(n) * ( 1.0_wp - exp_term )
4017                particles(n)%speed_y = particles(n)%speed_y * exp_term +             &
4018                                       v_int(n) * ( 1.0_wp - exp_term )
4019                particles(n)%speed_z = particles(n)%speed_z * exp_term +             &
4020                                       ( w_int(n) - ( 1.0_wp - dens_ratio(n) ) * g / &
4021                                       exp_arg ) * ( 1.0_wp - exp_term )
4022             ENDIF
4023          ENDDO
4024       ENDDO
4025
4026    ENDIF
4027
4028!
4029!-- Store the old age of the particle ( needed to prevent that a
4030!-- particle crosses several PEs during one timestep, and for the
4031!-- evaluation of the subgrid particle velocity fluctuations )
4032    particles(1:number_of_particles)%age_m = particles(1:number_of_particles)%age
4033
4034    DO  nb = 0, 7
4035       DO  n = start_index(nb), end_index(nb)
4036!
4037!--       Increment the particle age and the total time that the particle
4038!--       has advanced within the particle timestep procedure
4039          particles(n)%age    = particles(n)%age    + dt_particle(n)
4040          particles(n)%dt_sum = particles(n)%dt_sum + dt_particle(n)
4041
4042!
4043!--       Check whether there is still a particle that has not yet completed
4044!--       the total LES timestep
4045          IF ( ( dt_3d - particles(n)%dt_sum ) > 1E-8_wp )  THEN
4046             dt_3d_reached_l = .FALSE.
4047          ENDIF
4048
4049       ENDDO
4050    ENDDO
4051
4052    CALL cpu_log( log_point_s(44), 'lpm_advec', 'pause' )
4053
4054
4055 END SUBROUTINE lpm_advec
4056
4057 
4058!------------------------------------------------------------------------------! 
4059! Description:
4060! ------------
4061!> Calculation of subgrid-scale particle speed using the stochastic model
4062!> of Weil et al. (2004, JAS, 61, 2877-2887).
4063!------------------------------------------------------------------------------!
4064 SUBROUTINE weil_stochastic_eq( v_sgs, fs_n, e_n, dedxi_n, dedt_n, diss_n,     &
4065                                dt_n, rg_n, fac )
4066                               
4067    REAL(wp) ::  a1      !< dummy argument
4068    REAL(wp) ::  dedt_n  !< time derivative of TKE at particle position
4069    REAL(wp) ::  dedxi_n !< horizontal derivative of TKE at particle position
4070    REAL(wp) ::  diss_n  !< dissipation at particle position
4071    REAL(wp) ::  dt_n    !< particle timestep
4072    REAL(wp) ::  e_n     !< TKE at particle position
4073    REAL(wp) ::  fac     !< flag to identify adjacent topography
4074    REAL(wp) ::  fs_n    !< weighting factor to prevent that subgrid-scale particle speed becomes too large
4075    REAL(wp) ::  rg_n    !< random number
4076    REAL(wp) ::  term1   !< memory term
4077    REAL(wp) ::  term2   !< drift correction term
4078    REAL(wp) ::  term3   !< random term
4079    REAL(wp) ::  v_sgs   !< subgrid-scale velocity component
4080
4081!-- At first, limit TKE to a small non-zero number, in order to prevent
4082!-- the occurrence of extremely large SGS-velocities in case TKE is zero,
4083!-- (could occur at the simulation begin).
4084    e_n = MAX( e_n, 1E-20_wp )
4085!
4086!-- Please note, terms 1 and 2 (drift and memory term, respectively) are
4087!-- multiplied by a flag to switch of both terms near topography.
4088!-- This is necessary, as both terms may cause a subgrid-scale velocity build up
4089!-- if particles are trapped in regions with very small TKE, e.g. in narrow street
4090!-- canyons resolved by only a few grid points. Hence, term 1 and term 2 are
4091!-- disabled if one of