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

Last change on this file since 4232 was 4232, checked in by knoop, 6 years ago

Bugfix: wrong placement of INCLUDE "mpif.h" fixed.

  • Property svn:keywords set to Id
File size: 349.6 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 4232 2019-09-20 09:34:22Z knoop $
27! Removed INCLUDE "mpif.h", as it is not needed because of USE pegrid
28!
29! 4195 2019-08-28 13:44:27Z schwenkel
30! Bugfix for simple_corrector interpolation method in case of ocean runs and
31! output particle advection interpolation method into header
32!
33! 4182 2019-08-22 15:20:23Z scharf
34! Corrected "Former revisions" section
35!
36! 4168 2019-08-16 13:50:17Z suehring
37! Replace function get_topography_top_index by topo_top_ind
38!
39! 4145 2019-08-06 09:55:22Z schwenkel
40! Some reformatting
41!
42! 4144 2019-08-06 09:11:47Z raasch
43! relational operators .EQ., .NE., etc. replaced by ==, /=, etc.
44!
45! 4143 2019-08-05 15:14:53Z schwenkel
46! Rename variable and change select case to if statement
47!
48! 4122 2019-07-26 13:11:56Z schwenkel
49! Implement reset method as bottom boundary condition
50!
51! 4121 2019-07-26 10:01:22Z schwenkel
52! Implementation of an simple method for interpolating the velocities to
53! particle position
54!
55! 4114 2019-07-23 14:09:27Z schwenkel
56! Bugfix: Added working precision for if statement
57!
58! 4054 2019-06-27 07:42:18Z raasch
59! bugfix for calculating the minimum particle time step
60!
61! 4044 2019-06-19 12:28:27Z schwenkel
62! Bugfix in case of grid strecting: corrected calculation of k-Index
63!
64! 4043 2019-06-18 16:59:00Z schwenkel
65! Remove min_nr_particle, Add lpm_droplet_interactions_ptq into module
66!
67! 4028 2019-06-13 12:21:37Z schwenkel
68! Further modularization of particle code components
69!
70! 4020 2019-06-06 14:57:48Z schwenkel
71! Removing submodules
72!
73! 4018 2019-06-06 13:41:50Z eckhard
74! Bugfix for former revision
75!
76! 4017 2019-06-06 12:16:46Z schwenkel
77! Modularization of all lagrangian particle model code components
78!
79! 3655 2019-01-07 16:51:22Z knoop
80! bugfix to guarantee correct particle releases in case that the release
81! interval is smaller than the model timestep
82!
83! Revision 1.1  1999/11/25 16:16:06  raasch
84! Initial revision
85!
86!
87! Description:
88! ------------
89!> The embedded LPM allows for studying transport and dispersion processes within
90!> turbulent flows. This model including passive particles that do not show any
91!> feedback on the turbulent flow. Further also particles with inertia and
92!> cloud droplets ca be simulated explicitly.
93!>
94!> @todo test lcm
95!>       implement simple interpolation method for subgrid scale velocites
96!> @note <Enter notes on the module>
97!> @bug  <Enter bug on the module>
98!------------------------------------------------------------------------------!
99 MODULE lagrangian_particle_model_mod
100
101    USE, INTRINSIC ::  ISO_C_BINDING
102
103    USE arrays_3d,                                                             &
104        ONLY:  de_dx, de_dy, de_dz, dzw, zu, zw,  ql_c, ql_v, ql_vp, hyp,      &
105               pt, q, exner, ql, diss, e, u, v, w, km, ql_1, ql_2, pt_p, q_p,  &
106               d_exner, u_p, v_p, w_p
107 
108    USE averaging,                                                             &
109        ONLY:  ql_c_av, pr_av, pc_av, ql_vp_av, ql_v_av
110
111    USE basic_constants_and_equations_mod,                                     &
112        ONLY: molecular_weight_of_solute, molecular_weight_of_water, magnus,   &
113              pi, rd_d_rv, rho_l, r_v, rho_s, vanthoff, l_v, kappa, g, lv_d_cp
114
115    USE control_parameters,                                                    &
116        ONLY:  bc_dirichlet_l, bc_dirichlet_n, bc_dirichlet_r, bc_dirichlet_s, &
117               cloud_droplets, constant_flux_layer, current_timestep_number,   &
118               dt_3d, dt_3d_reached, humidity,                                 &
119               dt_3d_reached_l, dt_dopts, dz, initializing_actions,            &
120               message_string, molecular_viscosity, ocean_mode,                &
121               particle_maximum_age, iran,                                     & 
122               simulated_time, topography, dopts_time_count,                   &
123               time_since_reference_point, rho_surface, u_gtrans, v_gtrans,    &
124               dz_stretch_level, dz_stretch_level_start
125
126    USE cpulog,                                                                &
127        ONLY:  cpu_log, log_point, log_point_s
128
129    USE indices,                                                               &
130        ONLY:  nx, nxl, nxlg, nxrg, nxr, ny, nyn, nys, nyng, nysg, nz, nzb,    &
131               nzb_max, nzt,nbgp, ngp_2dh_outer,                               &
132               topo_top_ind,                                                   &
133               wall_flags_0
134
135    USE kinds
136
137    USE pegrid
138
139    USE particle_attributes
140
141    USE pmc_particle_interface,                                                &
142        ONLY: pmcp_c_get_particle_from_parent, pmcp_p_fill_particle_win,       &
143              pmcp_c_send_particle_to_parent, pmcp_p_empty_particle_win,       &
144              pmcp_p_delete_particles_in_fine_grid_area, pmcp_g_init,          &
145              pmcp_g_print_number_of_particles
146
147    USE pmc_interface,                                                         &
148        ONLY: nested_run
149
150    USE grid_variables,                                                        &
151        ONLY:  ddx, dx, ddy, dy
152
153    USE netcdf_interface,                                                      &
154        ONLY:  netcdf_data_format, netcdf_deflate, dopts_num, id_set_pts,      &
155               id_var_dopts, id_var_time_pts, nc_stat,                         &
156               netcdf_handle_error
157
158    USE random_function_mod,                                                   &
159        ONLY:  random_function
160
161    USE statistics,                                                            &
162        ONLY:  hom
163
164    USE surface_mod,                                                           &
165        ONLY:  bc_h,                                                           &
166               surf_def_h,                                                     &
167               surf_lsm_h,                                                     &
168               surf_usm_h
169
170#if defined( __parallel )  &&  !defined( __mpifh )
171    USE MPI
172#endif
173
174#if defined( __netcdf )
175    USE NETCDF
176#endif
177
178    IMPLICIT NONE
179
180    CHARACTER(LEN=15) ::  aero_species = 'nacl'                   !< aerosol species
181    CHARACTER(LEN=15) ::  aero_type    = 'maritime'               !< aerosol type
182    CHARACTER(LEN=15) ::  bc_par_lr    = 'cyclic'                 !< left/right boundary condition
183    CHARACTER(LEN=15) ::  bc_par_ns    = 'cyclic'                 !< north/south boundary condition
184    CHARACTER(LEN=15) ::  bc_par_b     = 'reflect'                !< bottom boundary condition
185    CHARACTER(LEN=15) ::  bc_par_t     = 'absorb'                 !< top boundary condition
186    CHARACTER(LEN=15) ::  collision_kernel   = 'none'             !< collision kernel
187
188    CHARACTER(LEN=5)  ::  splitting_function = 'gamma'            !< function for calculation critical weighting factor
189    CHARACTER(LEN=5)  ::  splitting_mode     = 'const'            !< splitting mode
190
191    CHARACTER(LEN=25) ::  particle_advection_interpolation = 'trilinear' !< interpolation method for calculatin the particle
192
193    INTEGER(iwp) ::  deleted_particles = 0                        !< number of deleted particles per time step   
194    INTEGER(iwp) ::  i_splitting_mode                             !< dummy for splitting mode
195    INTEGER(iwp) ::  iran_part = -1234567                         !< number for random generator   
196    INTEGER(iwp) ::  max_number_particles_per_gridbox = 100       !< namelist parameter (see documentation)
197    INTEGER(iwp) ::  isf                                          !< dummy for splitting function
198    INTEGER(iwp) ::  number_particles_per_gridbox = -1            !< namelist parameter (see documentation)
199    INTEGER(iwp) ::  number_of_sublayers = 20                     !< number of sublayers for particle velocities betwenn surface and first grid level
200    INTEGER(iwp) ::  offset_ocean_nzt = 0                         !< in case of oceans runs, the vertical index calculations need an offset
201    INTEGER(iwp) ::  offset_ocean_nzt_m1 = 0                      !< in case of oceans runs, the vertical index calculations need an offset
202    INTEGER(iwp) ::  particles_per_point = 1                      !< namelist parameter (see documentation)
203    INTEGER(iwp) ::  radius_classes = 20                          !< namelist parameter (see documentation)
204    INTEGER(iwp) ::  splitting_factor = 2                         !< namelist parameter (see documentation)
205    INTEGER(iwp) ::  splitting_factor_max = 5                     !< namelist parameter (see documentation)
206    INTEGER(iwp) ::  step_dealloc = 100                           !< namelist parameter (see documentation)
207    INTEGER(iwp) ::  total_number_of_particles                    !< total number of particles in the whole model domain
208    INTEGER(iwp) ::  trlp_count_sum                               !< parameter for particle exchange of PEs
209    INTEGER(iwp) ::  trlp_count_recv_sum                          !< parameter for particle exchange of PEs
210    INTEGER(iwp) ::  trrp_count_sum                               !< parameter for particle exchange of PEs
211    INTEGER(iwp) ::  trrp_count_recv_sum                          !< parameter for particle exchange of PEs
212    INTEGER(iwp) ::  trsp_count_sum                               !< parameter for particle exchange of PEs
213    INTEGER(iwp) ::  trsp_count_recv_sum                          !< parameter for particle exchange of PEs
214    INTEGER(iwp) ::  trnp_count_sum                               !< parameter for particle exchange of PEs
215    INTEGER(iwp) ::  trnp_count_recv_sum                          !< parameter for particle exchange of PEs
216
217    LOGICAL ::  lagrangian_particle_model = .FALSE.       !< namelist parameter (see documentation)
218    LOGICAL ::  curvature_solution_effects = .FALSE.      !< namelist parameter (see documentation)
219    LOGICAL ::  deallocate_memory = .TRUE.                !< namelist parameter (see documentation)
220    LOGICAL ::  hall_kernel = .FALSE.                     !< flag for collision kernel
221    LOGICAL ::  merging = .FALSE.                         !< namelist parameter (see documentation)
222    LOGICAL ::  random_start_position = .FALSE.           !< namelist parameter (see documentation)
223    LOGICAL ::  read_particles_from_restartfile = .TRUE.  !< namelist parameter (see documentation)
224    LOGICAL ::  seed_follows_topography = .FALSE.         !< namelist parameter (see documentation)
225    LOGICAL ::  splitting = .FALSE.                       !< namelist parameter (see documentation)
226    LOGICAL ::  use_kernel_tables = .FALSE.               !< parameter, which turns on the use of precalculated collision kernels
227    LOGICAL ::  write_particle_statistics = .FALSE.       !< namelist parameter (see documentation)
228    LOGICAL ::  interpolation_simple_predictor = .FALSE.  !< flag for simple particle advection interpolation with predictor step
229    LOGICAL ::  interpolation_simple_corrector = .FALSE.  !< flag for simple particle advection interpolation with corrector step
230    LOGICAL ::  interpolation_trilinear = .FALSE.         !< flag for trilinear particle advection interpolation
231
232    LOGICAL, DIMENSION(max_number_of_particle_groups) ::   vertical_particle_advection = .TRUE. !< Switch for vertical particle transport
233
234    REAL(wp) ::  aero_weight = 1.0_wp                      !< namelist parameter (see documentation)
235    REAL(wp) ::  dt_min_part = 0.0002_wp                   !< minimum particle time step when SGS velocities are used (s)
236    REAL(wp) ::  dt_prel = 9999999.9_wp                    !< namelist parameter (see documentation)
237    REAL(wp) ::  dt_write_particle_data = 9999999.9_wp     !< namelist parameter (see documentation)
238    REAL(wp) ::  end_time_prel = 9999999.9_wp              !< namelist parameter (see documentation)
239    REAL(wp) ::  initial_weighting_factor = 1.0_wp         !< namelist parameter (see documentation)
240    REAL(wp) ::  last_particle_release_time = 0.0_wp       !< last time of particle release
241    REAL(wp) ::  log_sigma(3) = 1.0_wp                     !< namelist parameter (see documentation)
242    REAL(wp) ::  na(3) = 0.0_wp                            !< namelist parameter (see documentation)
243    REAL(wp) ::  number_concentration = -1.0_wp            !< namelist parameter (see documentation)
244    REAL(wp) ::  radius_merge = 1.0E-7_wp                  !< namelist parameter (see documentation)
245    REAL(wp) ::  radius_split = 40.0E-6_wp                 !< namelist parameter (see documentation)
246    REAL(wp) ::  rm(3) = 1.0E-6_wp                         !< namelist parameter (see documentation)
247    REAL(wp) ::  sgs_wf_part                               !< parameter for sgs
248    REAL(wp) ::  time_write_particle_data = 0.0_wp         !< write particle data at current time on file
249    REAL(wp) ::  weight_factor_merge = -1.0_wp             !< namelist parameter (see documentation)
250    REAL(wp) ::  weight_factor_split = -1.0_wp             !< namelist parameter (see documentation)
251    REAL(wp) ::  z0_av_global                              !< horizontal mean value of z0
252
253    REAL(wp) ::  rclass_lbound !<
254    REAL(wp) ::  rclass_ubound !<
255
256    REAL(wp), PARAMETER ::  c_0 = 3.0_wp         !< parameter for lagrangian timescale
257
258    REAL(wp), DIMENSION(max_number_of_particle_groups) ::  density_ratio = 9999999.9_wp  !< namelist parameter (see documentation)
259    REAL(wp), DIMENSION(max_number_of_particle_groups) ::  pdx = 9999999.9_wp            !< namelist parameter (see documentation)
260    REAL(wp), DIMENSION(max_number_of_particle_groups) ::  pdy = 9999999.9_wp            !< namelist parameter (see documentation)
261    REAL(wp), DIMENSION(max_number_of_particle_groups) ::  pdz = 9999999.9_wp            !< namelist parameter (see documentation)
262    REAL(wp), DIMENSION(max_number_of_particle_groups) ::  psb = 9999999.9_wp            !< namelist parameter (see documentation)
263    REAL(wp), DIMENSION(max_number_of_particle_groups) ::  psl = 9999999.9_wp            !< namelist parameter (see documentation)
264    REAL(wp), DIMENSION(max_number_of_particle_groups) ::  psn = 9999999.9_wp            !< namelist parameter (see documentation)
265    REAL(wp), DIMENSION(max_number_of_particle_groups) ::  psr = 9999999.9_wp            !< namelist parameter (see documentation)
266    REAL(wp), DIMENSION(max_number_of_particle_groups) ::  pss = 9999999.9_wp            !< namelist parameter (see documentation)
267    REAL(wp), DIMENSION(max_number_of_particle_groups) ::  pst = 9999999.9_wp            !< namelist parameter (see documentation).
268    REAL(wp), DIMENSION(max_number_of_particle_groups) ::  radius = 9999999.9_wp         !< namelist parameter (see documentation)
269
270    REAL(wp), DIMENSION(:), ALLOCATABLE     ::  log_z_z0   !< Precalculate LOG(z/z0) 
271
272    INTEGER(iwp), PARAMETER ::  NR_2_direction_move = 10000 !<
273    INTEGER(iwp)            ::  nr_move_north               !<
274    INTEGER(iwp)            ::  nr_move_south               !<
275
276    TYPE(particle_type), DIMENSION(:), ALLOCATABLE ::  move_also_north
277    TYPE(particle_type), DIMENSION(:), ALLOCATABLE ::  move_also_south
278
279    REAL(wp) ::  epsilon_collision !<
280    REAL(wp) ::  urms              !<
281
282    REAL(wp), DIMENSION(:),   ALLOCATABLE ::  epsclass  !< dissipation rate class
283    REAL(wp), DIMENSION(:),   ALLOCATABLE ::  radclass  !< radius class
284    REAL(wp), DIMENSION(:),   ALLOCATABLE ::  winf      !<
285
286    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  ec        !<
287    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  ecf       !<
288    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  gck       !<
289    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  hkernel   !<
290    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  hwratio   !<
291
292    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  ckernel !<
293
294    INTEGER(iwp), PARAMETER         ::  PHASE_INIT    = 1  !<
295    INTEGER(iwp), PARAMETER, PUBLIC ::  PHASE_RELEASE = 2  !<
296
297    SAVE
298
299    PRIVATE
300
301    PUBLIC lpm_parin,     &
302           lpm_header,    &
303           lpm_init_arrays,&
304           lpm_init,      &
305           lpm_actions,   &
306           lpm_data_output_ptseries, &
307           lpm_interaction_droplets_ptq, &
308           lpm_rrd_local_particles, &
309           lpm_wrd_local, &
310           lpm_rrd_global, &
311           lpm_wrd_global, &
312           lpm_rrd_local, &
313           lpm_check_parameters
314
315    PUBLIC lagrangian_particle_model
316
317    INTERFACE lpm_check_parameters
318       MODULE PROCEDURE lpm_check_parameters
319    END INTERFACE lpm_check_parameters
320
321    INTERFACE lpm_parin
322       MODULE PROCEDURE lpm_parin
323    END INTERFACE lpm_parin
324
325    INTERFACE lpm_header
326       MODULE PROCEDURE lpm_header
327    END INTERFACE lpm_header
328
329    INTERFACE lpm_init_arrays
330       MODULE PROCEDURE lpm_init_arrays
331    END INTERFACE lpm_init_arrays
332 
333    INTERFACE lpm_init
334       MODULE PROCEDURE lpm_init
335    END INTERFACE lpm_init
336
337    INTERFACE lpm_actions
338       MODULE PROCEDURE lpm_actions
339    END INTERFACE lpm_actions
340
341    INTERFACE lpm_data_output_ptseries
342       MODULE PROCEDURE lpm_data_output_ptseries
343    END INTERFACE
344
345    INTERFACE lpm_rrd_local_particles
346       MODULE PROCEDURE lpm_rrd_local_particles
347    END INTERFACE lpm_rrd_local_particles
348
349    INTERFACE lpm_rrd_global
350       MODULE PROCEDURE lpm_rrd_global
351    END INTERFACE lpm_rrd_global
352
353    INTERFACE lpm_rrd_local
354       MODULE PROCEDURE lpm_rrd_local
355    END INTERFACE lpm_rrd_local
356
357    INTERFACE lpm_wrd_local
358       MODULE PROCEDURE lpm_wrd_local
359    END INTERFACE lpm_wrd_local
360
361    INTERFACE lpm_wrd_global
362       MODULE PROCEDURE lpm_wrd_global
363    END INTERFACE lpm_wrd_global
364
365    INTERFACE lpm_advec
366       MODULE PROCEDURE lpm_advec
367    END INTERFACE lpm_advec
368
369    INTERFACE lpm_calc_liquid_water_content
370       MODULE PROCEDURE lpm_calc_liquid_water_content
371    END INTERFACE
372
373    INTERFACE lpm_interaction_droplets_ptq
374       MODULE PROCEDURE lpm_interaction_droplets_ptq
375       MODULE PROCEDURE lpm_interaction_droplets_ptq_ij
376    END INTERFACE lpm_interaction_droplets_ptq
377
378    INTERFACE lpm_boundary_conds
379       MODULE PROCEDURE lpm_boundary_conds
380    END INTERFACE lpm_boundary_conds
381
382    INTERFACE lpm_droplet_condensation
383       MODULE PROCEDURE lpm_droplet_condensation
384    END INTERFACE
385
386    INTERFACE lpm_droplet_collision
387       MODULE PROCEDURE lpm_droplet_collision
388    END INTERFACE lpm_droplet_collision
389
390    INTERFACE lpm_init_kernels
391       MODULE PROCEDURE lpm_init_kernels
392    END INTERFACE lpm_init_kernels
393
394    INTERFACE lpm_splitting
395       MODULE PROCEDURE lpm_splitting
396    END INTERFACE lpm_splitting
397
398    INTERFACE lpm_merging
399       MODULE PROCEDURE lpm_merging
400    END INTERFACE lpm_merging
401
402    INTERFACE lpm_exchange_horiz
403       MODULE PROCEDURE lpm_exchange_horiz
404    END INTERFACE lpm_exchange_horiz
405
406    INTERFACE lpm_move_particle
407       MODULE PROCEDURE lpm_move_particle
408    END INTERFACE lpm_move_particle
409
410    INTERFACE realloc_particles_array
411       MODULE PROCEDURE realloc_particles_array
412    END INTERFACE realloc_particles_array
413
414    INTERFACE dealloc_particles_array
415       MODULE PROCEDURE dealloc_particles_array
416    END INTERFACE dealloc_particles_array
417
418    INTERFACE lpm_sort_and_delete
419       MODULE PROCEDURE lpm_sort_and_delete
420    END INTERFACE lpm_sort_and_delete
421
422    INTERFACE lpm_sort_timeloop_done
423       MODULE PROCEDURE lpm_sort_timeloop_done
424    END INTERFACE lpm_sort_timeloop_done
425
426    INTERFACE lpm_pack
427       MODULE PROCEDURE lpm_pack
428    END INTERFACE lpm_pack
429
430 CONTAINS
431 
432
433!------------------------------------------------------------------------------!
434! Description:
435! ------------
436!> Parin for &particle_parameters for the Lagrangian particle model
437!------------------------------------------------------------------------------!
438 SUBROUTINE lpm_parin
439 
440    CHARACTER (LEN=80) ::  line  !<
441
442    NAMELIST /particles_par/ &
443       aero_species, &
444       aero_type, &
445       aero_weight, &
446       alloc_factor, &
447       bc_par_b, &
448       bc_par_lr, &
449       bc_par_ns, &
450       bc_par_t, &
451       collision_kernel, &
452       curvature_solution_effects, &
453       deallocate_memory, &
454       density_ratio, &
455       dissipation_classes, &
456       dt_dopts, &
457       dt_min_part, &
458       dt_prel, &
459       dt_write_particle_data, &
460       end_time_prel, &
461       initial_weighting_factor, &
462       log_sigma, &
463       max_number_particles_per_gridbox, &
464       merging, &
465       na, &
466       number_concentration, &
467       number_of_particle_groups, &
468       number_particles_per_gridbox, &
469       particles_per_point, &
470       particle_advection_start, &
471       particle_advection_interpolation, &
472       particle_maximum_age, &
473       pdx, &
474       pdy, &
475       pdz, &
476       psb, &
477       psl, &
478       psn, &
479       psr, &
480       pss, &
481       pst, &
482       radius, &
483       radius_classes, &
484       radius_merge, &
485       radius_split, &
486       random_start_position, &
487       read_particles_from_restartfile, &
488       rm, &
489       seed_follows_topography, &
490       splitting, &
491       splitting_factor, &
492       splitting_factor_max, &
493       splitting_function, &
494       splitting_mode, &
495       step_dealloc, &
496       use_sgs_for_particles, &
497       vertical_particle_advection, &
498       weight_factor_merge, &
499       weight_factor_split, &
500       write_particle_statistics
501
502       NAMELIST /particle_parameters/ &
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       na, &
526       number_concentration, &
527       number_of_particle_groups, &
528       number_particles_per_gridbox, &
529       particles_per_point, &
530       particle_advection_start, &
531       particle_advection_interpolation, &
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!
563!-- Position the namelist-file at the beginning (it was already opened in
564!-- parin), search for the namelist-group of the package and position the
565!-- file at this line. Do the same for each optionally used package.
566    line = ' '
567   
568!
569!-- Try to find particles package
570    REWIND ( 11 )
571    line = ' '
572    DO   WHILE ( INDEX( line, '&particle_parameters' ) == 0 )
573       READ ( 11, '(A)', END=12 )  line
574    ENDDO
575    BACKSPACE ( 11 )
576!
577!-- Read user-defined namelist
578    READ ( 11, particle_parameters, ERR = 10 )
579!
580!-- Set flag that indicates that particles are switched on
581    particle_advection = .TRUE.
582   
583    GOTO 14
584
58510  BACKSPACE( 11 )
586    READ( 11 , '(A)') line
587    CALL parin_fail_message( 'particle_parameters', line )
588!
589!-- Try to find particles package (old namelist)
59012  REWIND ( 11 )
591    line = ' '
592    DO WHILE ( INDEX( line, '&particles_par' ) == 0 )
593       READ ( 11, '(A)', END=14 )  line
594    ENDDO
595    BACKSPACE ( 11 )
596!
597!-- Read user-defined namelist
598    READ ( 11, particles_par, ERR = 13, END = 14 )
599
600    message_string = 'namelist particles_par is deprecated and will be ' //    &
601                     'removed in near future. Please use namelist ' //         &
602                     'particle_parameters instead'
603    CALL message( 'package_parin', 'PA0487', 0, 1, 0, 6, 0 )
604
605!
606!-- Set flag that indicates that particles are switched on
607    particle_advection = .TRUE.
608
609    GOTO 14
610
61113    BACKSPACE( 11 )
612       READ( 11 , '(A)') line
613       CALL parin_fail_message( 'particles_par', line )
614
61514 CONTINUE
616
617 END SUBROUTINE lpm_parin
618 
619!------------------------------------------------------------------------------!
620! Description:
621! ------------
622!> Writes used particle attributes in header file.
623!------------------------------------------------------------------------------!
624 SUBROUTINE lpm_header ( io )
625
626    CHARACTER (LEN=40) ::  output_format       !< netcdf format
627
628    INTEGER(iwp) ::  i               !<
629    INTEGER(iwp), INTENT(IN) ::  io  !< Unit of the output file
630
631
632     IF ( humidity  .AND.  cloud_droplets )  THEN
633       WRITE ( io, 433 )
634       IF ( curvature_solution_effects )  WRITE ( io, 434 )
635       IF ( collision_kernel /= 'none' )  THEN
636          WRITE ( io, 435 )  TRIM( collision_kernel )
637          IF ( collision_kernel(6:9) == 'fast' )  THEN
638             WRITE ( io, 436 )  radius_classes, dissipation_classes
639          ENDIF
640       ELSE
641          WRITE ( io, 437 )
642       ENDIF
643    ENDIF
644 
645    IF ( particle_advection )  THEN
646!
647!--    Particle attributes
648       WRITE ( io, 480 )  particle_advection_start, TRIM(particle_advection_interpolation), &
649                          dt_prel, bc_par_lr, &
650                          bc_par_ns, bc_par_b, bc_par_t, particle_maximum_age, &
651                          end_time_prel
652       IF ( use_sgs_for_particles )  WRITE ( io, 488 )  dt_min_part
653       IF ( random_start_position )  WRITE ( io, 481 )
654       IF ( seed_follows_topography )  WRITE ( io, 496 )
655       IF ( particles_per_point > 1 )  WRITE ( io, 489 )  particles_per_point
656       WRITE ( io, 495 )  total_number_of_particles
657       IF ( dt_write_particle_data /= 9999999.9_wp )  THEN
658          WRITE ( io, 485 )  dt_write_particle_data
659          IF ( netcdf_data_format > 1 )  THEN
660             output_format = 'netcdf (64 bit offset) and binary'
661          ELSE
662             output_format = 'netcdf and binary'
663          ENDIF
664          IF ( netcdf_deflate == 0 )  THEN
665             WRITE ( io, 344 )  output_format
666          ELSE
667             WRITE ( io, 354 )  TRIM( output_format ), netcdf_deflate
668          ENDIF
669       ENDIF
670       IF ( dt_dopts /= 9999999.9_wp )  WRITE ( io, 494 )  dt_dopts
671       IF ( write_particle_statistics )  WRITE ( io, 486 )
672
673       WRITE ( io, 487 )  number_of_particle_groups
674
675       DO  i = 1, number_of_particle_groups
676          IF ( i == 1  .AND.  density_ratio(i) == 9999999.9_wp )  THEN
677             WRITE ( io, 490 )  i, 0.0_wp
678             WRITE ( io, 492 )
679          ELSE
680             WRITE ( io, 490 )  i, radius(i)
681             IF ( density_ratio(i) /= 0.0_wp )  THEN
682                WRITE ( io, 491 )  density_ratio(i)
683             ELSE
684                WRITE ( io, 492 )
685             ENDIF
686          ENDIF
687          WRITE ( io, 493 )  psl(i), psr(i), pss(i), psn(i), psb(i), pst(i), &
688                             pdx(i), pdy(i), pdz(i)
689          IF ( .NOT. vertical_particle_advection(i) )  WRITE ( io, 482 )
690       ENDDO
691
692    ENDIF
693   
694344 FORMAT ('       Output format: ',A/)
695354 FORMAT ('       Output format: ',A, '   compressed with level: ',I1/)
696
697433 FORMAT ('    Cloud droplets treated explicitly using the Lagrangian part', &
698                 'icle model')
699434 FORMAT ('    Curvature and solution effecs are considered for growth of', &
700                 ' droplets < 1.0E-6 m')
701435 FORMAT ('    Droplet collision is handled by ',A,'-kernel')
702436 FORMAT ('       Fast kernel with fixed radius- and dissipation classes ', &
703                    'are used'/ &
704            '          number of radius classes:       ',I3,'    interval ', &
705                       '[1.0E-6,2.0E-4] m'/ &
706            '          number of dissipation classes:   ',I2,'    interval ', &
707                       '[0,1000] cm**2/s**3')
708437 FORMAT ('    Droplet collision is switched off')
709
710480 FORMAT ('    Particles:'/ &
711            '    ---------'// &
712            '       Particle advection is active (switched on at t = ', F7.1, &
713                    ' s)'/ &
714            '       Interpolation of particle velocities is done by using ', A, &
715                    ' method'/ &
716            '       Start of new particle generations every  ',F6.1,' s'/ &
717            '       Boundary conditions: left/right: ', A, ' north/south: ', A/&
718            '                            bottom:     ', A, ' top:         ', A/&
719            '       Maximum particle age:                 ',F9.1,' s'/ &
720            '       Advection stopped at t = ',F9.1,' s'/)
721481 FORMAT ('       Particles have random start positions'/)
722482 FORMAT ('          Particles are advected only horizontally'/)
723485 FORMAT ('       Particle data are written on file every ', F9.1, ' s')
724486 FORMAT ('       Particle statistics are written on file'/)
725487 FORMAT ('       Number of particle groups: ',I2/)
726488 FORMAT ('       SGS velocity components are used for particle advection'/ &
727            '          minimum timestep for advection:', F8.5/)
728489 FORMAT ('       Number of particles simultaneously released at each ', &
729                    'point: ', I5/)
730490 FORMAT ('       Particle group ',I2,':'/ &
731            '          Particle radius: ',E10.3, 'm')
732491 FORMAT ('          Particle inertia is activated'/ &
733            '             density_ratio (rho_fluid/rho_particle) =',F6.3/)
734492 FORMAT ('          Particles are advected only passively (no inertia)'/)
735493 FORMAT ('          Boundaries of particle source: x:',F8.1,' - ',F8.1,' m'/&
736            '                                         y:',F8.1,' - ',F8.1,' m'/&
737            '                                         z:',F8.1,' - ',F8.1,' m'/&
738            '          Particle distances:  dx = ',F8.1,' m  dy = ',F8.1, &
739                       ' m  dz = ',F8.1,' m'/)
740494 FORMAT ('       Output of particle time series in NetCDF format every ', &
741                    F8.2,' s'/)
742495 FORMAT ('       Number of particles in total domain: ',I10/)
743496 FORMAT ('       Initial vertical particle positions are interpreted ', &
744                    'as relative to the given topography')
745   
746 END SUBROUTINE lpm_header
747 
748!------------------------------------------------------------------------------!
749! Description:
750! ------------
751!> Writes used particle attributes in header file.
752!------------------------------------------------------------------------------! 
753 SUBROUTINE lpm_check_parameters
754 
755!
756!-- Collision kernels:
757    SELECT CASE ( TRIM( collision_kernel ) )
758
759       CASE ( 'hall', 'hall_fast' )
760          hall_kernel = .TRUE.
761
762       CASE ( 'wang', 'wang_fast' )
763          wang_kernel = .TRUE.
764
765       CASE ( 'none' )
766
767
768       CASE DEFAULT
769          message_string = 'unknown collision kernel: collision_kernel = "' // &
770                           TRIM( collision_kernel ) // '"'
771          CALL message( 'check_parameters', 'PA0350', 1, 2, 0, 6, 0 )
772
773    END SELECT
774    IF ( collision_kernel(6:9) == 'fast' )  use_kernel_tables = .TRUE.
775
776!
777!-- Subgrid scale velocites with the simple interpolation method for resolved
778!-- velocites is not implemented for passive particles. However, for cloud
779!-- it can be combined as the sgs-velocites for active particles are
780!-- calculated differently, i.e. no subboxes are needed.
781    IF ( .NOT. TRIM( particle_advection_interpolation ) == 'trilinear'  .AND.  &
782       use_sgs_for_particles  .AND.  .NOT. cloud_droplets )  THEN
783          message_string = 'subrgrid scale velocities in combination with ' // &
784                           'simple interpolation method is not '            // &
785                           'implemented'
786          CALL message( 'check_parameters', 'PA0659', 1, 2, 0, 6, 0 )
787    ENDIF
788
789 END SUBROUTINE
790 
791!------------------------------------------------------------------------------!
792! Description:
793! ------------
794!> Initialize arrays for lpm
795!------------------------------------------------------------------------------!   
796 SUBROUTINE lpm_init_arrays
797 
798    IF ( cloud_droplets )  THEN
799!
800!--    Liquid water content, change in liquid water content
801       ALLOCATE ( ql_1(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                      &
802                  ql_2(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
803!
804!--    Real volume of particles (with weighting), volume of particles
805       ALLOCATE ( ql_v(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                      &
806                     ql_vp(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
807    ENDIF
808   
809!
810!--    Initial assignment of the pointers   
811    IF ( cloud_droplets )  THEN
812       ql   => ql_1
813       ql_c => ql_2
814    ENDIF
815   
816 END SUBROUTINE lpm_init_arrays
817 
818!------------------------------------------------------------------------------!
819! Description:
820! ------------
821!> Initialize Lagrangian particle model
822!------------------------------------------------------------------------------!
823 SUBROUTINE lpm_init
824
825    INTEGER(iwp) ::  i                           !<
826    INTEGER(iwp) ::  j                           !<
827    INTEGER(iwp) ::  k                           !<
828
829    REAL(wp) ::  div                             !<
830    REAL(wp) ::  height_int                      !<
831    REAL(wp) ::  height_p                        !<
832    REAL(wp) ::  z_p                             !<
833    REAL(wp) ::  z0_av_local                     !<
834
835!
836!-- In case of oceans runs, the vertical index calculations need an offset,
837!-- because otherwise the k indices will become negative
838    IF ( ocean_mode )  THEN
839       offset_ocean_nzt    = nzt
840       offset_ocean_nzt_m1 = nzt - 1
841    ENDIF
842
843!
844!-- Define block offsets for dividing a gridcell in 8 sub cells
845!-- See documentation for List of subgrid boxes
846!-- See pack_and_sort in lpm_pack_arrays.f90 for assignment of the subgrid boxes
847    block_offset(0) = block_offset_def ( 0, 0, 0)
848    block_offset(1) = block_offset_def ( 0, 0,-1)
849    block_offset(2) = block_offset_def ( 0,-1, 0)
850    block_offset(3) = block_offset_def ( 0,-1,-1)
851    block_offset(4) = block_offset_def (-1, 0, 0)
852    block_offset(5) = block_offset_def (-1, 0,-1)
853    block_offset(6) = block_offset_def (-1,-1, 0)
854    block_offset(7) = block_offset_def (-1,-1,-1)
855!
856!-- Check the number of particle groups.
857    IF ( number_of_particle_groups > max_number_of_particle_groups )  THEN
858       WRITE( message_string, * ) 'max_number_of_particle_groups =',           &
859                                  max_number_of_particle_groups ,              &
860                                  '&number_of_particle_groups reset to ',      &
861                                  max_number_of_particle_groups
862       CALL message( 'lpm_init', 'PA0213', 0, 1, 0, 6, 0 )
863       number_of_particle_groups = max_number_of_particle_groups
864    ENDIF
865!
866!-- Check if downward-facing walls exist. This case, reflection boundary
867!-- conditions (as well as subgrid-scale velocities) may do not work
868!-- propably (not realized so far).
869    IF ( surf_def_h(1)%ns >= 1 )  THEN
870       WRITE( message_string, * ) 'Overhanging topography do not work '//      &
871                                  'with particles'
872       CALL message( 'lpm_init', 'PA0212', 0, 1, 0, 6, 0 )
873
874    ENDIF
875
876!
877!-- Set default start positions, if necessary
878    IF ( psl(1) == 9999999.9_wp )  psl(1) = 0.0_wp
879    IF ( psr(1) == 9999999.9_wp )  psr(1) = ( nx +1 ) * dx
880    IF ( pss(1) == 9999999.9_wp )  pss(1) = 0.0_wp
881    IF ( psn(1) == 9999999.9_wp )  psn(1) = ( ny +1 ) * dy
882    IF ( psb(1) == 9999999.9_wp )  psb(1) = zu(nz/2)
883    IF ( pst(1) == 9999999.9_wp )  pst(1) = psb(1)
884
885    IF ( pdx(1) == 9999999.9_wp  .OR.  pdx(1) == 0.0_wp )  pdx(1) = dx
886    IF ( pdy(1) == 9999999.9_wp  .OR.  pdy(1) == 0.0_wp )  pdy(1) = dy
887    IF ( pdz(1) == 9999999.9_wp  .OR.  pdz(1) == 0.0_wp )  pdz(1) = zu(2) - zu(1)
888
889!
890!-- If number_particles_per_gridbox is set, the parametres pdx, pdy and pdz are
891!-- calculated diagnostically. Therfore an isotropic distribution is prescribed.
892    IF ( number_particles_per_gridbox /= -1 .AND.   &
893         number_particles_per_gridbox >= 1 )    THEN
894       pdx(1) = (( dx * dy * ( zu(2) - zu(1) ) ) /  &
895             REAL(number_particles_per_gridbox))**0.3333333_wp
896!
897!--    Ensure a smooth value (two significant digits) of distance between
898!--    particles (pdx, pdy, pdz).
899       div = 1000.0_wp
900       DO  WHILE ( pdx(1) < div )
901          div = div / 10.0_wp
902       ENDDO
903       pdx(1) = NINT( pdx(1) * 100.0_wp / div ) * div / 100.0_wp
904       pdy(1) = pdx(1)
905       pdz(1) = pdx(1)
906
907    ENDIF
908
909    DO  j = 2, number_of_particle_groups
910       IF ( psl(j) == 9999999.9_wp )  psl(j) = psl(j-1)
911       IF ( psr(j) == 9999999.9_wp )  psr(j) = psr(j-1)
912       IF ( pss(j) == 9999999.9_wp )  pss(j) = pss(j-1)
913       IF ( psn(j) == 9999999.9_wp )  psn(j) = psn(j-1)
914       IF ( psb(j) == 9999999.9_wp )  psb(j) = psb(j-1)
915       IF ( pst(j) == 9999999.9_wp )  pst(j) = pst(j-1)
916       IF ( pdx(j) == 9999999.9_wp  .OR.  pdx(j) == 0.0_wp )  pdx(j) = pdx(j-1)
917       IF ( pdy(j) == 9999999.9_wp  .OR.  pdy(j) == 0.0_wp )  pdy(j) = pdy(j-1)
918       IF ( pdz(j) == 9999999.9_wp  .OR.  pdz(j) == 0.0_wp )  pdz(j) = pdz(j-1)
919    ENDDO
920
921!
922!-- Allocate arrays required for calculating particle SGS velocities.
923!-- Initialize prefactor required for stoachastic Weil equation.
924    IF ( use_sgs_for_particles  .AND.  .NOT. cloud_droplets )  THEN
925       ALLOCATE( de_dx(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
926                 de_dy(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
927                 de_dz(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
928
929       de_dx = 0.0_wp
930       de_dy = 0.0_wp
931       de_dz = 0.0_wp
932
933       sgs_wf_part = 1.0_wp / 3.0_wp
934    ENDIF
935
936!
937!-- Allocate array required for logarithmic vertical interpolation of
938!-- horizontal particle velocities between the surface and the first vertical
939!-- grid level. In order to avoid repeated CPU cost-intensive CALLS of
940!-- intrinsic FORTRAN procedure LOG(z/z0), LOG(z/z0) is precalculated for
941!-- several heights. Splitting into 20 sublayers turned out to be sufficient.
942!-- To obtain exact height levels of particles, linear interpolation is applied
943!-- (see lpm_advec.f90).
944    IF ( constant_flux_layer )  THEN
945
946       ALLOCATE ( log_z_z0(0:number_of_sublayers) )
947       z_p = zu(nzb+1) - zw(nzb)
948
949!
950!--    Calculate horizontal mean value of z0 used for logartihmic
951!--    interpolation. Note: this is not exact for heterogeneous z0.
952!--    However, sensitivity studies showed that the effect is
953!--    negligible.
954       z0_av_local  = SUM( surf_def_h(0)%z0 ) + SUM( surf_lsm_h%z0 ) +         &
955                      SUM( surf_usm_h%z0 )
956       z0_av_global = 0.0_wp
957
958#if defined( __parallel )
959       CALL MPI_ALLREDUCE(z0_av_local, z0_av_global, 1, MPI_REAL, MPI_SUM, &
960                          comm2d, ierr )
961#else
962       z0_av_global = z0_av_local
963#endif
964
965       z0_av_global = z0_av_global  / ( ( ny + 1 ) * ( nx + 1 ) )
966!
967!--    Horizontal wind speed is zero below and at z0
968       log_z_z0(0) = 0.0_wp
969!
970!--    Calculate vertical depth of the sublayers
971       height_int  = ( z_p - z0_av_global ) / REAL( number_of_sublayers, KIND=wp )
972!
973!--    Precalculate LOG(z/z0)
974       height_p    = z0_av_global
975       DO  k = 1, number_of_sublayers
976
977          height_p    = height_p + height_int
978          log_z_z0(k) = LOG( height_p / z0_av_global )
979
980       ENDDO
981
982    ENDIF
983
984!
985!-- Check which particle interpolation method should be used
986    IF ( TRIM( particle_advection_interpolation )  ==  'trilinear' )  THEN
987       interpolation_simple_corrector = .FALSE.
988       interpolation_simple_predictor = .FALSE.
989       interpolation_trilinear        = .TRUE.
990    ELSEIF ( TRIM( particle_advection_interpolation )  ==  'simple_corrector' )  THEN
991       interpolation_simple_corrector = .TRUE.
992       interpolation_simple_predictor = .FALSE.
993       interpolation_trilinear        = .FALSE.
994    ELSEIF ( TRIM( particle_advection_interpolation )  ==  'simple_predictor' )  THEN
995       interpolation_simple_corrector = .FALSE.
996       interpolation_simple_predictor = .TRUE.
997       interpolation_trilinear        = .FALSE.
998    ENDIF
999
1000!
1001!-- Check boundary condition and set internal variables
1002    SELECT CASE ( bc_par_b )
1003
1004       CASE ( 'absorb' )
1005          ibc_par_b = 1
1006
1007       CASE ( 'reflect' )
1008          ibc_par_b = 2
1009
1010       CASE ( 'reset' )
1011          ibc_par_b = 3
1012
1013       CASE DEFAULT
1014          WRITE( message_string, * )  'unknown boundary condition ',           &
1015                                       'bc_par_b = "', TRIM( bc_par_b ), '"'
1016          CALL message( 'lpm_init', 'PA0217', 1, 2, 0, 6, 0 )
1017
1018    END SELECT
1019    SELECT CASE ( bc_par_t )
1020
1021       CASE ( 'absorb' )
1022          ibc_par_t = 1
1023
1024       CASE ( 'reflect' )
1025          ibc_par_t = 2
1026
1027       CASE ( 'nested' )
1028          ibc_par_t = 3
1029
1030       CASE DEFAULT
1031          WRITE( message_string, * ) 'unknown boundary condition ',            &
1032                                     'bc_par_t = "', TRIM( bc_par_t ), '"'
1033          CALL message( 'lpm_init', 'PA0218', 1, 2, 0, 6, 0 )
1034
1035    END SELECT
1036    SELECT CASE ( bc_par_lr )
1037
1038       CASE ( 'cyclic' )
1039          ibc_par_lr = 0
1040
1041       CASE ( 'absorb' )
1042          ibc_par_lr = 1
1043
1044       CASE ( 'reflect' )
1045          ibc_par_lr = 2
1046
1047       CASE ( 'nested' )
1048          ibc_par_lr = 3
1049
1050       CASE DEFAULT
1051          WRITE( message_string, * ) 'unknown boundary condition ',   &
1052                                     'bc_par_lr = "', TRIM( bc_par_lr ), '"'
1053          CALL message( 'lpm_init', 'PA0219', 1, 2, 0, 6, 0 )
1054
1055    END SELECT
1056    SELECT CASE ( bc_par_ns )
1057
1058       CASE ( 'cyclic' )
1059          ibc_par_ns = 0
1060
1061       CASE ( 'absorb' )
1062          ibc_par_ns = 1
1063
1064       CASE ( 'reflect' )
1065          ibc_par_ns = 2
1066
1067       CASE ( 'nested' )
1068          ibc_par_ns = 3
1069
1070       CASE DEFAULT
1071          WRITE( message_string, * ) 'unknown boundary condition ',   &
1072                                     'bc_par_ns = "', TRIM( bc_par_ns ), '"'
1073          CALL message( 'lpm_init', 'PA0220', 1, 2, 0, 6, 0 )
1074
1075    END SELECT
1076    SELECT CASE ( splitting_mode )
1077
1078       CASE ( 'const' )
1079          i_splitting_mode = 1
1080
1081       CASE ( 'cl_av' )
1082          i_splitting_mode = 2
1083
1084       CASE ( 'gb_av' )
1085          i_splitting_mode = 3
1086
1087       CASE DEFAULT
1088          WRITE( message_string, * )  'unknown splitting_mode = "',            &
1089                                      TRIM( splitting_mode ), '"'
1090          CALL message( 'lpm_init', 'PA0146', 1, 2, 0, 6, 0 )
1091
1092    END SELECT
1093    SELECT CASE ( splitting_function )
1094
1095       CASE ( 'gamma' )
1096          isf = 1
1097
1098       CASE ( 'log' )
1099          isf = 2
1100
1101       CASE ( 'exp' )
1102          isf = 3
1103
1104       CASE DEFAULT
1105          WRITE( message_string, * )  'unknown splitting function = "',        &
1106                                       TRIM( splitting_function ), '"'
1107          CALL message( 'lpm_init', 'PA0147', 1, 2, 0, 6, 0 )
1108
1109    END SELECT
1110!
1111!-- Initialize collision kernels
1112    IF ( collision_kernel /= 'none' )  CALL lpm_init_kernels
1113!
1114!-- For the first model run of a possible job chain initialize the
1115!-- particles, otherwise read the particle data from restart file.
1116    IF ( TRIM( initializing_actions ) == 'read_restart_data'  &
1117         .AND.  read_particles_from_restartfile )  THEN
1118       CALL lpm_rrd_local_particles
1119    ELSE
1120!
1121!--    Allocate particle arrays and set attributes of the initial set of
1122!--    particles, which can be also periodically released at later times.
1123       ALLOCATE( prt_count(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
1124                 grid_particles(nzb+1:nzt,nys:nyn,nxl:nxr) )
1125
1126       number_of_particles = 0
1127       prt_count           = 0
1128!
1129!--    initialize counter for particle IDs
1130       grid_particles%id_counter = 1
1131!
1132!--    Initialize all particles with dummy values (otherwise errors may
1133!--    occur within restart runs). The reason for this is still not clear
1134!--    and may be presumably caused by errors in the respective user-interface.
1135       zero_particle = particle_type( 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp,  &
1136                                      0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp,  &
1137                                      0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp,  &
1138                                      0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp,  &
1139                                      0, 0, 0_idp, .FALSE., -1 )
1140
1141       particle_groups = particle_groups_type( 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp )
1142!
1143!--    Set values for the density ratio and radius for all particle
1144!--    groups, if necessary
1145       IF ( density_ratio(1) == 9999999.9_wp )  density_ratio(1) = 0.0_wp
1146       IF ( radius(1)        == 9999999.9_wp )  radius(1) = 0.0_wp
1147       DO  i = 2, number_of_particle_groups
1148          IF ( density_ratio(i) == 9999999.9_wp )  THEN
1149             density_ratio(i) = density_ratio(i-1)
1150          ENDIF
1151          IF ( radius(i) == 9999999.9_wp )  radius(i) = radius(i-1)
1152       ENDDO
1153
1154       DO  i = 1, number_of_particle_groups
1155          IF ( density_ratio(i) /= 0.0_wp  .AND.  radius(i) == 0 )  THEN
1156             WRITE( message_string, * ) 'particle group #', i, ' has a',       &
1157                                        'density ratio /= 0 but radius = 0'
1158             CALL message( 'lpm_init', 'PA0215', 1, 2, 0, 6, 0 )
1159          ENDIF
1160          particle_groups(i)%density_ratio = density_ratio(i)
1161          particle_groups(i)%radius        = radius(i)
1162       ENDDO
1163!
1164!--    Set a seed value for the random number generator to be exclusively
1165!--    used for the particle code. The generated random numbers should be
1166!--    different on the different PEs.
1167       iran_part = iran_part + myid
1168!
1169!--    Create the particle set, and set the initial particles
1170       CALL lpm_create_particle( phase_init )
1171       last_particle_release_time = particle_advection_start
1172!
1173!--    User modification of initial particles
1174       CALL user_lpm_init
1175!
1176!--    Open file for statistical informations about particle conditions
1177       IF ( write_particle_statistics )  THEN
1178          CALL check_open( 80 )
1179          WRITE ( 80, 8000 )  current_timestep_number, simulated_time,         &
1180                              number_of_particles
1181          CALL close_file( 80 )
1182       ENDIF
1183
1184    ENDIF
1185
1186    IF ( nested_run )  CALL pmcp_g_init
1187!
1188!-- To avoid programm abort, assign particles array to the local version of
1189!-- first grid cell
1190    number_of_particles = prt_count(nzb+1,nys,nxl)
1191    particles => grid_particles(nzb+1,nys,nxl)%particles(1:number_of_particles)
1192!
1193!-- Formats
11948000 FORMAT (I6,1X,F7.2,4X,I10,71X,I10)
1195
1196 END SUBROUTINE lpm_init
1197 
1198!------------------------------------------------------------------------------!
1199! Description:
1200! ------------
1201!> Create Lagrangian particles
1202!------------------------------------------------------------------------------! 
1203 SUBROUTINE lpm_create_particle (phase)
1204
1205    INTEGER(iwp)               ::  alloc_size  !< relative increase of allocated memory for particles
1206    INTEGER(iwp)               ::  i           !< loop variable ( particle groups )
1207    INTEGER(iwp)               ::  ip          !< index variable along x
1208    INTEGER(iwp)               ::  j           !< loop variable ( particles per point )
1209    INTEGER(iwp)               ::  jp          !< index variable along y
1210    INTEGER(iwp)               ::  k           !< index variable along z
1211    INTEGER(iwp)               ::  k_surf      !< index of surface grid point
1212    INTEGER(iwp)               ::  kp          !< index variable along z
1213    INTEGER(iwp)               ::  loop_stride !< loop variable for initialization
1214    INTEGER(iwp)               ::  n           !< loop variable ( number of particles )
1215    INTEGER(iwp)               ::  new_size    !< new size of allocated memory for particles
1216
1217    INTEGER(iwp), INTENT(IN)   ::  phase       !< mode of inititialization
1218
1219    INTEGER(iwp), DIMENSION(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ::  local_count !< start address of new particle
1220    INTEGER(iwp), DIMENSION(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ::  local_start !< start address of new particle
1221
1222    LOGICAL                    ::  first_stride !< flag for initialization
1223
1224    REAL(wp)                   ::  pos_x      !< increment for particle position in x
1225    REAL(wp)                   ::  pos_y      !< increment for particle position in y
1226    REAL(wp)                   ::  pos_z      !< increment for particle position in z
1227    REAL(wp)                   ::  rand_contr !< dummy argument for random position
1228
1229    TYPE(particle_type),TARGET ::  tmp_particle !< temporary particle used for initialization
1230
1231
1232!
1233!-- Calculate particle positions and store particle attributes, if
1234!-- particle is situated on this PE
1235    DO  loop_stride = 1, 2
1236       first_stride = (loop_stride == 1)
1237       IF ( first_stride )   THEN
1238          local_count = 0           ! count number of particles
1239       ELSE
1240          local_count = prt_count   ! Start address of new particles
1241       ENDIF
1242
1243!
1244!--    Calculate initial_weighting_factor diagnostically
1245       IF ( number_concentration /= -1.0_wp  .AND.  number_concentration > 0.0_wp )  THEN
1246          initial_weighting_factor =  number_concentration  *                           &
1247                                      pdx(1) * pdy(1) * pdz(1)
1248       END IF
1249
1250       n = 0
1251       DO  i = 1, number_of_particle_groups
1252          pos_z = psb(i)
1253          DO WHILE ( pos_z <= pst(i) )
1254             IF ( pos_z >= zw(0) .AND.  pos_z < zw(nzt) )  THEN
1255                pos_y = pss(i)
1256                DO WHILE ( pos_y <= psn(i) )
1257                   IF ( pos_y >= nys * dy  .AND.                  &
1258                        pos_y <  ( nyn + 1 ) * dy  )  THEN
1259                      pos_x = psl(i)
1260               xloop: DO WHILE ( pos_x <= psr(i) )
1261                         IF ( pos_x >= nxl * dx  .AND.            &
1262                              pos_x <  ( nxr + 1) * dx )  THEN
1263                            DO  j = 1, particles_per_point
1264                               n = n + 1
1265                               tmp_particle%x             = pos_x
1266                               tmp_particle%y             = pos_y
1267                               tmp_particle%z             = pos_z
1268                               tmp_particle%age           = 0.0_wp
1269                               tmp_particle%age_m         = 0.0_wp
1270                               tmp_particle%dt_sum        = 0.0_wp
1271                               tmp_particle%e_m           = 0.0_wp
1272                               tmp_particle%rvar1         = 0.0_wp
1273                               tmp_particle%rvar2         = 0.0_wp
1274                               tmp_particle%rvar3         = 0.0_wp
1275                               tmp_particle%speed_x       = 0.0_wp
1276                               tmp_particle%speed_y       = 0.0_wp
1277                               tmp_particle%speed_z       = 0.0_wp
1278                               tmp_particle%origin_x      = pos_x
1279                               tmp_particle%origin_y      = pos_y
1280                               tmp_particle%origin_z      = pos_z
1281                               IF ( curvature_solution_effects )  THEN
1282                                  tmp_particle%aux1      = 0.0_wp    ! dry aerosol radius
1283                                  tmp_particle%aux2      = dt_3d     ! last Rosenbrock timestep
1284                               ELSE
1285                                  tmp_particle%aux1      = 0.0_wp    ! free to use
1286                                  tmp_particle%aux2      = 0.0_wp    ! free to use
1287                               ENDIF
1288                               tmp_particle%radius        = particle_groups(i)%radius
1289                               tmp_particle%weight_factor = initial_weighting_factor
1290                               tmp_particle%class         = 1
1291                               tmp_particle%group         = i
1292                               tmp_particle%id            = 0_idp
1293                               tmp_particle%particle_mask = .TRUE.
1294                               tmp_particle%block_nr      = -1
1295!
1296!--                            Determine the grid indices of the particle position
1297                               ip = INT( tmp_particle%x * ddx )
1298                               jp = INT( tmp_particle%y * ddy )
1299!
1300!--                            In case of stretching the actual k index is found iteratively
1301                               IF ( dz_stretch_level /= -9999999.9_wp  .OR.           &
1302                                    dz_stretch_level_start(1) /= -9999999.9_wp )  THEN
1303                                  kp = MINLOC( ABS( tmp_particle%z - zu ), DIM = 1 ) - 1
1304                               ELSE
1305                                  kp = INT( tmp_particle%z / dz(1) + 1 + offset_ocean_nzt )
1306                               ENDIF
1307!
1308!--                            Determine surface level. Therefore, check for
1309!--                            upward-facing wall on w-grid.
1310                               k_surf = topo_top_ind(jp,ip,3)
1311                               IF ( seed_follows_topography )  THEN
1312!
1313!--                               Particle height is given relative to topography
1314                                  kp = kp + k_surf
1315                                  tmp_particle%z = tmp_particle%z + zw(k_surf)
1316!--                               Skip particle release if particle position is
1317!--                               above model top, or within topography in case
1318!--                               of overhanging structures.
1319                                  IF ( kp > nzt  .OR.                          &
1320                                 .NOT. BTEST( wall_flags_0(kp,jp,ip), 0 ) )  THEN
1321                                     pos_x = pos_x + pdx(i)
1322                                     CYCLE xloop
1323                                  ENDIF
1324!
1325!--                            Skip particle release if particle position is
1326!--                            below surface, or within topography in case
1327!--                            of overhanging structures.
1328                               ELSEIF ( .NOT. seed_follows_topography .AND.    &
1329                                         tmp_particle%z <= zw(k_surf)  .OR.    &
1330                                        .NOT. BTEST( wall_flags_0(kp,jp,ip), 0 ) )&
1331                               THEN
1332                                  pos_x = pos_x + pdx(i)
1333                                  CYCLE xloop
1334                               ENDIF
1335
1336                               local_count(kp,jp,ip) = local_count(kp,jp,ip) + 1
1337
1338                               IF ( .NOT. first_stride )  THEN
1339                                  IF ( ip < nxl  .OR.  jp < nys  .OR.  kp < nzb+1 )  THEN
1340                                     write(6,*) 'xl ',ip,jp,kp,nxl,nys,nzb+1
1341                                  ENDIF
1342                                  IF ( ip > nxr  .OR.  jp > nyn  .OR.  kp > nzt )  THEN
1343                                     write(6,*) 'xu ',ip,jp,kp,nxr,nyn,nzt
1344                                  ENDIF
1345                                  grid_particles(kp,jp,ip)%particles(local_count(kp,jp,ip)) = tmp_particle
1346                               ENDIF
1347                            ENDDO
1348                         ENDIF
1349                         pos_x = pos_x + pdx(i)
1350                      ENDDO xloop
1351                   ENDIF
1352                   pos_y = pos_y + pdy(i)
1353                ENDDO
1354             ENDIF
1355
1356             pos_z = pos_z + pdz(i)
1357          ENDDO
1358       ENDDO
1359
1360       IF ( first_stride )  THEN
1361          DO  ip = nxl, nxr
1362             DO  jp = nys, nyn
1363                DO  kp = nzb+1, nzt
1364                   IF ( phase == PHASE_INIT )  THEN
1365                      IF ( local_count(kp,jp,ip) > 0 )  THEN
1366                         alloc_size = MAX( INT( local_count(kp,jp,ip) *        &
1367                            ( 1.0_wp + alloc_factor / 100.0_wp ) ),            &
1368                            1 )
1369                      ELSE
1370                         alloc_size = 1
1371                      ENDIF
1372                      ALLOCATE(grid_particles(kp,jp,ip)%particles(1:alloc_size))
1373                      DO  n = 1, alloc_size
1374                         grid_particles(kp,jp,ip)%particles(n) = zero_particle
1375                      ENDDO
1376                   ELSEIF ( phase == PHASE_RELEASE )  THEN
1377                      IF ( local_count(kp,jp,ip) > 0 )  THEN
1378                         new_size   = local_count(kp,jp,ip) + prt_count(kp,jp,ip)
1379                         alloc_size = MAX( INT( new_size * ( 1.0_wp +          &
1380                            alloc_factor / 100.0_wp ) ), 1 )
1381                         IF( alloc_size > SIZE( grid_particles(kp,jp,ip)%particles) )  THEN
1382                            CALL realloc_particles_array( ip, jp, kp, alloc_size )
1383                         ENDIF
1384                      ENDIF
1385                   ENDIF
1386                ENDDO
1387             ENDDO
1388          ENDDO
1389       ENDIF
1390
1391    ENDDO
1392
1393    local_start = prt_count+1
1394    prt_count   = local_count
1395!
1396!-- Calculate particle IDs
1397    DO  ip = nxl, nxr
1398       DO  jp = nys, nyn
1399          DO  kp = nzb+1, nzt
1400             number_of_particles = prt_count(kp,jp,ip)
1401             IF ( number_of_particles <= 0 )  CYCLE
1402             particles => grid_particles(kp,jp,ip)%particles(1:number_of_particles)
1403
1404             DO  n = local_start(kp,jp,ip), number_of_particles  !only new particles
1405
1406                particles(n)%id = 10000_idp**3 * grid_particles(kp,jp,ip)%id_counter + &
1407                                  10000_idp**2 * kp + 10000_idp * jp + ip
1408!
1409!--             Count the number of particles that have been released before
1410                grid_particles(kp,jp,ip)%id_counter =                          &
1411                                         grid_particles(kp,jp,ip)%id_counter + 1
1412
1413             ENDDO
1414
1415          ENDDO
1416       ENDDO
1417    ENDDO
1418!
1419!-- Initialize aerosol background spectrum
1420    IF ( curvature_solution_effects )  THEN
1421       CALL lpm_init_aerosols( local_start )
1422    ENDIF
1423!
1424!-- Add random fluctuation to particle positions.
1425    IF ( random_start_position )  THEN
1426       DO  ip = nxl, nxr
1427          DO  jp = nys, nyn
1428             DO  kp = nzb+1, nzt
1429                number_of_particles = prt_count(kp,jp,ip)
1430                IF ( number_of_particles <= 0 )  CYCLE
1431                particles => grid_particles(kp,jp,ip)%particles(1:number_of_particles)
1432!
1433!--             Move only new particles. Moreover, limit random fluctuation
1434!--             in order to prevent that particles move more than one grid box,
1435!--             which would lead to problems concerning particle exchange
1436!--             between processors in case pdx/pdy are larger than dx/dy,
1437!--             respectively.
1438                DO  n = local_start(kp,jp,ip), number_of_particles
1439                   IF ( psl(particles(n)%group) /= psr(particles(n)%group) )  THEN
1440                      rand_contr = ( random_function( iran_part ) - 0.5_wp ) * &
1441                                     pdx(particles(n)%group)
1442                      particles(n)%x = particles(n)%x +                        &
1443                              MERGE( rand_contr, SIGN( dx, rand_contr ),       &
1444                                     ABS( rand_contr ) < dx                    &
1445                                   )
1446                   ENDIF
1447                   IF ( pss(particles(n)%group) /= psn(particles(n)%group) )  THEN
1448                      rand_contr = ( random_function( iran_part ) - 0.5_wp ) * &
1449                                     pdy(particles(n)%group)
1450                      particles(n)%y = particles(n)%y +                        &
1451                              MERGE( rand_contr, SIGN( dy, rand_contr ),       &
1452                                     ABS( rand_contr ) < dy                    &
1453                                   )
1454                   ENDIF
1455                   IF ( psb(particles(n)%group) /= pst(particles(n)%group) )  THEN
1456                      rand_contr = ( random_function( iran_part ) - 0.5_wp ) * &
1457                                     pdz(particles(n)%group)
1458                      particles(n)%z = particles(n)%z +                        &
1459                              MERGE( rand_contr, SIGN( dzw(kp), rand_contr ),  &
1460                                     ABS( rand_contr ) < dzw(kp)               &
1461                                   )
1462                   ENDIF
1463                ENDDO
1464!
1465!--             Identify particles located outside the model domain and reflect
1466!--             or absorb them if necessary.
1467                CALL lpm_boundary_conds( 'bottom/top', i, j, k )
1468!
1469!--             Furthermore, remove particles located in topography. Note, as
1470!--             the particle speed is still zero at this point, wall
1471!--             reflection boundary conditions will not work in this case.
1472                particles =>                                                   &
1473                       grid_particles(kp,jp,ip)%particles(1:number_of_particles)
1474                DO  n = local_start(kp,jp,ip), number_of_particles
1475                   i = particles(n)%x * ddx
1476                   j = particles(n)%y * ddy
1477                   k = particles(n)%z / dz(1) + 1 + offset_ocean_nzt
1478                   DO WHILE( zw(k) < particles(n)%z )
1479                      k = k + 1
1480                   ENDDO
1481                   DO WHILE( zw(k-1) > particles(n)%z )
1482                      k = k - 1
1483                   ENDDO
1484!
1485!--                Check if particle is within topography
1486                   IF ( .NOT. BTEST( wall_flags_0(k,j,i), 0 ) )  THEN
1487                      particles(n)%particle_mask = .FALSE.
1488                      deleted_particles = deleted_particles + 1
1489                   ENDIF
1490
1491                ENDDO
1492             ENDDO
1493          ENDDO
1494       ENDDO
1495!
1496!--    Exchange particles between grid cells and processors
1497       CALL lpm_move_particle
1498       CALL lpm_exchange_horiz
1499
1500    ENDIF
1501!
1502!-- In case of random_start_position, delete particles identified by
1503!-- lpm_exchange_horiz and lpm_boundary_conds. Then sort particles into blocks,
1504!-- which is needed for a fast interpolation of the LES fields on the particle
1505!-- position.
1506    CALL lpm_sort_and_delete
1507!
1508!-- Determine the current number of particles
1509    DO  ip = nxl, nxr
1510       DO  jp = nys, nyn
1511          DO  kp = nzb+1, nzt
1512             number_of_particles         = number_of_particles                 &
1513                                           + prt_count(kp,jp,ip)
1514          ENDDO
1515       ENDDO
1516    ENDDO
1517!
1518!-- Calculate the number of particles of the total domain
1519#if defined( __parallel )
1520    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1521    CALL MPI_ALLREDUCE( number_of_particles, total_number_of_particles, 1, &
1522    MPI_INTEGER, MPI_SUM, comm2d, ierr )
1523#else
1524    total_number_of_particles = number_of_particles
1525#endif
1526
1527    RETURN
1528
1529 END SUBROUTINE lpm_create_particle
1530 
1531 
1532!------------------------------------------------------------------------------!
1533! Description:
1534! ------------
1535!> This routine initialize the particles as aerosols with physio-chemical
1536!> properties.
1537!------------------------------------------------------------------------------!   
1538 SUBROUTINE lpm_init_aerosols(local_start)
1539
1540    REAL(wp) ::  afactor            !< curvature effects
1541    REAL(wp) ::  bfactor            !< solute effects
1542    REAL(wp) ::  dlogr              !< logarithmic width of radius bin
1543    REAL(wp) ::  e_a                !< vapor pressure
1544    REAL(wp) ::  e_s                !< saturation vapor pressure
1545    REAL(wp) ::  rmin = 0.005e-6_wp !< minimum aerosol radius
1546    REAL(wp) ::  rmax = 10.0e-6_wp  !< maximum aerosol radius
1547    REAL(wp) ::  r_mid              !< mean radius of bin
1548    REAL(wp) ::  r_l                !< left radius of bin
1549    REAL(wp) ::  r_r                !< right radius of bin
1550    REAL(wp) ::  sigma              !< surface tension
1551    REAL(wp) ::  t_int              !< temperature
1552
1553    INTEGER(iwp), DIMENSION(nzb:nzt+1,nysg:nyng,nxlg:nxrg), INTENT(IN) ::  local_start !<
1554
1555    INTEGER(iwp) ::  n              !<
1556    INTEGER(iwp) ::  ip             !<
1557    INTEGER(iwp) ::  jp             !<
1558    INTEGER(iwp) ::  kp             !<
1559
1560!
1561!-- Set constants for different aerosol species
1562    IF ( TRIM( aero_species ) == 'nacl' )  THEN
1563       molecular_weight_of_solute = 0.05844_wp 
1564       rho_s                      = 2165.0_wp
1565       vanthoff                   = 2.0_wp
1566    ELSEIF ( TRIM( aero_species ) == 'c3h4o4' )  THEN
1567       molecular_weight_of_solute = 0.10406_wp 
1568       rho_s                      = 1600.0_wp
1569       vanthoff                   = 1.37_wp
1570    ELSEIF ( TRIM( aero_species ) == 'nh4o3' )  THEN
1571       molecular_weight_of_solute = 0.08004_wp 
1572       rho_s                      = 1720.0_wp
1573       vanthoff                   = 2.31_wp
1574    ELSE
1575       WRITE( message_string, * ) 'unknown aerosol species ',   &
1576                                'aero_species = "', TRIM( aero_species ), '"'
1577       CALL message( 'lpm_init', 'PA0470', 1, 2, 0, 6, 0 )
1578    ENDIF
1579!
1580!-- The following typical aerosol spectra are taken from Jaenicke (1993):
1581!-- Tropospheric aerosols. Published in Aerosol-Cloud-Climate Interactions.
1582    IF ( TRIM( aero_type ) == 'polar' )  THEN
1583       na        = (/ 2.17e1, 1.86e-1, 3.04e-4 /) * 1.0E6_wp
1584       rm        = (/ 0.0689, 0.375, 4.29 /) * 1.0E-6_wp
1585       log_sigma = (/ 0.245, 0.300, 0.291 /)
1586    ELSEIF ( TRIM( aero_type ) == 'background' )  THEN
1587       na        = (/ 1.29e2, 5.97e1, 6.35e1 /) * 1.0E6_wp
1588       rm        = (/ 0.0036, 0.127, 0.259 /) * 1.0E-6_wp
1589       log_sigma = (/ 0.645, 0.253, 0.425 /)
1590    ELSEIF ( TRIM( aero_type ) == 'maritime' )  THEN
1591       na        = (/ 1.33e2, 6.66e1, 3.06e0 /) * 1.0E6_wp
1592       rm        = (/ 0.0039, 0.133, 0.29 /) * 1.0E-6_wp
1593       log_sigma = (/ 0.657, 0.210, 0.396 /)
1594    ELSEIF ( TRIM( aero_type ) == 'continental' )  THEN
1595       na        = (/ 3.20e3, 2.90e3, 3.00e-1 /) * 1.0E6_wp
1596       rm        = (/ 0.01, 0.058, 0.9 /) * 1.0E-6_wp
1597       log_sigma = (/ 0.161, 0.217, 0.380 /)
1598    ELSEIF ( TRIM( aero_type ) == 'desert' )  THEN
1599       na        = (/ 7.26e2, 1.14e3, 1.78e-1 /) * 1.0E6_wp
1600       rm        = (/ 0.001, 0.0188, 10.8 /) * 1.0E-6_wp
1601       log_sigma = (/ 0.247, 0.770, 0.438 /)
1602    ELSEIF ( TRIM( aero_type ) == 'rural' )  THEN
1603       na        = (/ 6.65e3, 1.47e2, 1.99e3 /) * 1.0E6_wp
1604       rm        = (/ 0.00739, 0.0269, 0.0419 /) * 1.0E-6_wp
1605       log_sigma = (/ 0.225, 0.557, 0.266 /)
1606    ELSEIF ( TRIM( aero_type ) == 'urban' )  THEN
1607       na        = (/ 9.93e4, 1.11e3, 3.64e4 /) * 1.0E6_wp
1608       rm        = (/ 0.00651, 0.00714, 0.0248 /) * 1.0E-6_wp
1609       log_sigma = (/ 0.245, 0.666, 0.337 /)
1610    ELSEIF ( TRIM( aero_type ) == 'user' )  THEN
1611       CONTINUE
1612    ELSE
1613       WRITE( message_string, * ) 'unknown aerosol type ',   &
1614                                'aero_type = "', TRIM( aero_type ), '"'
1615       CALL message( 'lpm_init', 'PA0459', 1, 2, 0, 6, 0 )
1616    ENDIF
1617
1618    DO  ip = nxl, nxr
1619       DO  jp = nys, nyn
1620          DO  kp = nzb+1, nzt
1621
1622             number_of_particles = prt_count(kp,jp,ip)
1623             IF ( number_of_particles <= 0 )  CYCLE
1624             particles => grid_particles(kp,jp,ip)%particles(1:number_of_particles)
1625
1626             dlogr   = ( LOG10(rmax) - LOG10(rmin) ) / ( number_of_particles - local_start(kp,jp,ip) + 1 )
1627!
1628!--          Initialize the aerosols with a predefined spectral distribution
1629!--          of the dry radius (logarithmically increasing bins) and a varying
1630!--          weighting factor
1631             DO  n = local_start(kp,jp,ip), number_of_particles  !only new particles
1632
1633                r_l   = 10.0**( LOG10( rmin ) + (n-1) * dlogr )
1634                r_r   = 10.0**( LOG10( rmin ) + n * dlogr )
1635                r_mid = SQRT( r_l * r_r )
1636
1637                particles(n)%aux1          = r_mid
1638                particles(n)%weight_factor =                                           &
1639                   ( na(1) / ( SQRT( 2.0_wp * pi ) * log_sigma(1) ) *                     &
1640                     EXP( - LOG10( r_mid / rm(1) )**2 / ( 2.0_wp * log_sigma(1)**2 ) ) +  &
1641                     na(2) / ( SQRT( 2.0_wp * pi ) * log_sigma(2) ) *                     &
1642                     EXP( - LOG10( r_mid / rm(2) )**2 / ( 2.0_wp * log_sigma(2)**2 ) ) +  &
1643                     na(3) / ( SQRT( 2.0_wp * pi ) * log_sigma(3) ) *                     &
1644                     EXP( - LOG10( r_mid / rm(3) )**2 / ( 2.0_wp * log_sigma(3)**2 ) )    &
1645                   ) * ( LOG10(r_r) - LOG10(r_l) ) * ( dx * dy * dzw(kp) )
1646
1647!
1648!--             Multiply weight_factor with the namelist parameter aero_weight
1649!--             to increase or decrease the number of simulated aerosols
1650                particles(n)%weight_factor = particles(n)%weight_factor * aero_weight
1651
1652                IF ( particles(n)%weight_factor - FLOOR(particles(n)%weight_factor,KIND=wp) &
1653                     > random_function( iran_part ) )  THEN
1654                   particles(n)%weight_factor = FLOOR(particles(n)%weight_factor,KIND=wp) + 1.0_wp
1655                ELSE
1656                   particles(n)%weight_factor = FLOOR(particles(n)%weight_factor,KIND=wp)
1657                ENDIF
1658!
1659!--             Unnecessary particles will be deleted
1660                IF ( particles(n)%weight_factor <= 0.0_wp )  particles(n)%particle_mask = .FALSE.
1661
1662             ENDDO
1663!
1664!--          Set particle radius to equilibrium radius based on the environmental
1665!--          supersaturation (Khvorostyanov and Curry, 2007, JGR). This avoids
1666!--          the sometimes lengthy growth toward their equilibrium radius within
1667!--          the simulation.
1668             t_int  = pt(kp,jp,ip) * exner(kp)
1669
1670             e_s = magnus( t_int )
1671             e_a = q(kp,jp,ip) * hyp(kp) / ( q(kp,jp,ip) + rd_d_rv )
1672
1673             sigma   = 0.0761_wp - 0.000155_wp * ( t_int - 273.15_wp )
1674             afactor = 2.0_wp * sigma / ( rho_l * r_v * t_int )
1675
1676             bfactor = vanthoff * molecular_weight_of_water *    &
1677                       rho_s / ( molecular_weight_of_solute * rho_l )
1678!
1679!--          The formula is only valid for subsaturated environments. For
1680!--          supersaturations higher than -5 %, the supersaturation is set to -5%.
1681             IF ( e_a / e_s >= 0.95_wp )  e_a = 0.95_wp * e_s
1682
1683             DO  n = local_start(kp,jp,ip), number_of_particles  !only new particles
1684!
1685!--             For details on this equation, see Eq. (14) of Khvorostyanov and
1686!--             Curry (2007, JGR)
1687                particles(n)%radius = bfactor**0.3333333_wp *                  &
1688                   particles(n)%aux1 / ( 1.0_wp - e_a / e_s )**0.3333333_wp / &
1689                   ( 1.0_wp + ( afactor / ( 3.0_wp * bfactor**0.3333333_wp *   &
1690                     particles(n)%aux1 ) ) /                                  &
1691                     ( 1.0_wp - e_a / e_s )**0.6666666_wp                      &
1692                   )
1693
1694             ENDDO
1695
1696          ENDDO
1697       ENDDO
1698    ENDDO
1699
1700 END SUBROUTINE lpm_init_aerosols
1701
1702
1703!------------------------------------------------------------------------------!
1704! Description:
1705! ------------
1706!> Calculates quantities required for considering the SGS velocity fluctuations
1707!> in the particle transport by a stochastic approach. The respective
1708!> quantities are: SGS-TKE gradients and horizontally averaged profiles of the
1709!> SGS TKE and the resolved-scale velocity variances.
1710!------------------------------------------------------------------------------!
1711 SUBROUTINE lpm_init_sgs_tke
1712
1713    USE statistics,                                                            &
1714        ONLY:  flow_statistics_called, hom, sums, sums_l
1715
1716    INTEGER(iwp) ::  i      !< index variable along x
1717    INTEGER(iwp) ::  j      !< index variable along y
1718    INTEGER(iwp) ::  k      !< index variable along z
1719    INTEGER(iwp) ::  m      !< running index for the surface elements
1720
1721    REAL(wp) ::  flag1      !< flag to mask topography
1722
1723!
1724!-- TKE gradient along x and y
1725    DO  i = nxl, nxr
1726       DO  j = nys, nyn
1727          DO  k = nzb, nzt+1
1728
1729             IF ( .NOT. BTEST( wall_flags_0(k,j,i-1), 0 )  .AND.               &
1730                        BTEST( wall_flags_0(k,j,i), 0   )  .AND.               &
1731                        BTEST( wall_flags_0(k,j,i+1), 0 ) )                    &
1732             THEN
1733                de_dx(k,j,i) = 2.0_wp * sgs_wf_part *                          &
1734                               ( e(k,j,i+1) - e(k,j,i) ) * ddx
1735             ELSEIF ( BTEST( wall_flags_0(k,j,i-1), 0 )  .AND.                 &
1736                      BTEST( wall_flags_0(k,j,i), 0   )  .AND.                 &
1737                .NOT. BTEST( wall_flags_0(k,j,i+1), 0 ) )                      &
1738             THEN
1739                de_dx(k,j,i) = 2.0_wp * sgs_wf_part *                          &
1740                               ( e(k,j,i) - e(k,j,i-1) ) * ddx
1741             ELSEIF ( .NOT. BTEST( wall_flags_0(k,j,i), 22   )  .AND.          &
1742                      .NOT. BTEST( wall_flags_0(k,j,i+1), 22 ) )               &   
1743             THEN
1744                de_dx(k,j,i) = 0.0_wp
1745             ELSEIF ( .NOT. BTEST( wall_flags_0(k,j,i-1), 22 )  .AND.          &
1746                      .NOT. BTEST( wall_flags_0(k,j,i), 22   ) )               &
1747             THEN
1748                de_dx(k,j,i) = 0.0_wp
1749             ELSE
1750                de_dx(k,j,i) = sgs_wf_part * ( e(k,j,i+1) - e(k,j,i-1) ) * ddx
1751             ENDIF
1752
1753             IF ( .NOT. BTEST( wall_flags_0(k,j-1,i), 0 )  .AND.               &
1754                        BTEST( wall_flags_0(k,j,i), 0   )  .AND.               &
1755                        BTEST( wall_flags_0(k,j+1,i), 0 ) )                    &
1756             THEN
1757                de_dy(k,j,i) = 2.0_wp * sgs_wf_part *                          &
1758                               ( e(k,j+1,i) - e(k,j,i) ) * ddy
1759             ELSEIF ( BTEST( wall_flags_0(k,j-1,i), 0 )  .AND.                 &
1760                      BTEST( wall_flags_0(k,j,i), 0   )  .AND.                 &
1761                .NOT. BTEST( wall_flags_0(k,j+1,i), 0 ) )                      &
1762             THEN
1763                de_dy(k,j,i) = 2.0_wp * sgs_wf_part *                          &
1764                               ( e(k,j,i) - e(k,j-1,i) ) * ddy
1765             ELSEIF ( .NOT. BTEST( wall_flags_0(k,j,i), 22   )  .AND.          &
1766                      .NOT. BTEST( wall_flags_0(k,j+1,i), 22 ) )               &   
1767             THEN
1768                de_dy(k,j,i) = 0.0_wp
1769             ELSEIF ( .NOT. BTEST( wall_flags_0(k,j-1,i), 22 )  .AND.          &
1770                      .NOT. BTEST( wall_flags_0(k,j,i), 22   ) )               &
1771             THEN
1772                de_dy(k,j,i) = 0.0_wp
1773             ELSE
1774                de_dy(k,j,i) = sgs_wf_part * ( e(k,j+1,i) - e(k,j-1,i) ) * ddy
1775             ENDIF
1776
1777          ENDDO
1778       ENDDO
1779    ENDDO
1780
1781!
1782!-- TKE gradient along z at topograhy and  including bottom and top boundary conditions
1783    DO  i = nxl, nxr
1784       DO  j = nys, nyn
1785          DO  k = nzb+1, nzt-1
1786!
1787!--          Flag to mask topography
1788             flag1 = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_0(k,j,i), 0  ) )
1789
1790             de_dz(k,j,i) = 2.0_wp * sgs_wf_part *                             &
1791                           ( e(k+1,j,i) - e(k-1,j,i) ) / ( zu(k+1) - zu(k-1) ) &
1792                                                 * flag1 
1793          ENDDO
1794!
1795!--       upward-facing surfaces
1796          DO  m = bc_h(0)%start_index(j,i), bc_h(0)%end_index(j,i)
1797             k            = bc_h(0)%k(m)
1798             de_dz(k,j,i) = 2.0_wp * sgs_wf_part *                             &
1799                           ( e(k+1,j,i) - e(k,j,i)   ) / ( zu(k+1) - zu(k) )
1800          ENDDO
1801!
1802!--       downward-facing surfaces
1803          DO  m = bc_h(1)%start_index(j,i), bc_h(1)%end_index(j,i)
1804             k            = bc_h(1)%k(m)
1805             de_dz(k,j,i) = 2.0_wp * sgs_wf_part *                             &
1806                           ( e(k,j,i) - e(k-1,j,i)   ) / ( zu(k) - zu(k-1) )
1807          ENDDO
1808
1809          de_dz(nzb,j,i)   = 0.0_wp
1810          de_dz(nzt,j,i)   = 0.0_wp
1811          de_dz(nzt+1,j,i) = 0.0_wp
1812       ENDDO
1813    ENDDO
1814!
1815!-- Ghost point exchange
1816    CALL exchange_horiz( de_dx, nbgp )
1817    CALL exchange_horiz( de_dy, nbgp )
1818    CALL exchange_horiz( de_dz, nbgp )
1819    CALL exchange_horiz( diss, nbgp  )
1820!
1821!-- Set boundary conditions at non-periodic boundaries. Note, at non-period
1822!-- boundaries zero-gradient boundary conditions are set for the subgrid TKE.
1823!-- Thus, TKE gradients normal to the respective lateral boundaries are zero,
1824!-- while tangetial TKE gradients then must be the same as within the prognostic
1825!-- domain. 
1826    IF ( bc_dirichlet_l )  THEN
1827       de_dx(:,:,-1) = 0.0_wp
1828       de_dy(:,:,-1) = de_dy(:,:,0) 
1829       de_dz(:,:,-1) = de_dz(:,:,0)
1830    ENDIF
1831    IF ( bc_dirichlet_r )  THEN
1832       de_dx(:,:,nxr+1) = 0.0_wp
1833       de_dy(:,:,nxr+1) = de_dy(:,:,nxr) 
1834       de_dz(:,:,nxr+1) = de_dz(:,:,nxr)
1835    ENDIF
1836    IF ( bc_dirichlet_n )  THEN
1837       de_dx(:,nyn+1,:) = de_dx(:,nyn,:)
1838       de_dy(:,nyn+1,:) = 0.0_wp 
1839       de_dz(:,nyn+1,:) = de_dz(:,nyn,:)
1840    ENDIF
1841    IF ( bc_dirichlet_s )  THEN
1842       de_dx(:,nys-1,:) = de_dx(:,nys,:)
1843       de_dy(:,nys-1,:) = 0.0_wp 
1844       de_dz(:,nys-1,:) = de_dz(:,nys,:)
1845    ENDIF 
1846!
1847!-- Calculate the horizontally averaged profiles of SGS TKE and resolved
1848!-- velocity variances (they may have been already calculated in routine
1849!-- flow_statistics).
1850    IF ( .NOT. flow_statistics_called )  THEN
1851
1852!
1853!--    First calculate horizontally averaged profiles of the horizontal
1854!--    velocities.
1855       sums_l(:,1,0) = 0.0_wp
1856       sums_l(:,2,0) = 0.0_wp
1857
1858       DO  i = nxl, nxr
1859          DO  j =  nys, nyn
1860             DO  k = nzb, nzt+1
1861!
1862!--             Flag indicating vicinity of wall
1863                flag1 = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_0(k,j,i), 24 ) )
1864
1865                sums_l(k,1,0)  = sums_l(k,1,0)  + u(k,j,i) * flag1
1866                sums_l(k,2,0)  = sums_l(k,2,0)  + v(k,j,i) * flag1
1867             ENDDO
1868          ENDDO
1869       ENDDO
1870
1871#if defined( __parallel )
1872!
1873!--    Compute total sum from local sums
1874       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1875       CALL MPI_ALLREDUCE( sums_l(nzb,1,0), sums(nzb,1), nzt+2-nzb, &
1876                           MPI_REAL, MPI_SUM, comm2d, ierr )
1877       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1878       CALL MPI_ALLREDUCE( sums_l(nzb,2,0), sums(nzb,2), nzt+2-nzb, &
1879                              MPI_REAL, MPI_SUM, comm2d, ierr )
1880#else
1881       sums(:,1) = sums_l(:,1,0)
1882       sums(:,2) = sums_l(:,2,0)
1883#endif
1884
1885!
1886!--    Final values are obtained by division by the total number of grid
1887!--    points used for the summation.
1888       hom(:,1,1,0) = sums(:,1) / ngp_2dh_outer(:,0)   ! u
1889       hom(:,1,2,0) = sums(:,2) / ngp_2dh_outer(:,0)   ! v
1890
1891!
1892!--    Now calculate the profiles of SGS TKE and the resolved-scale
1893!--    velocity variances
1894       sums_l(:,8,0)  = 0.0_wp
1895       sums_l(:,30,0) = 0.0_wp
1896       sums_l(:,31,0) = 0.0_wp
1897       sums_l(:,32,0) = 0.0_wp
1898       DO  i = nxl, nxr
1899          DO  j = nys, nyn
1900             DO  k = nzb, nzt+1
1901!
1902!--             Flag indicating vicinity of wall
1903                flag1 = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_0(k,j,i), 24 ) )
1904
1905                sums_l(k,8,0)  = sums_l(k,8,0)  + e(k,j,i)                       * flag1
1906                sums_l(k,30,0) = sums_l(k,30,0) + ( u(k,j,i) - hom(k,1,1,0) )**2 * flag1
1907                sums_l(k,31,0) = sums_l(k,31,0) + ( v(k,j,i) - hom(k,1,2,0) )**2 * flag1
1908                sums_l(k,32,0) = sums_l(k,32,0) + w(k,j,i)**2                    * flag1
1909             ENDDO
1910          ENDDO
1911       ENDDO
1912
1913#if defined( __parallel )
1914!
1915!--    Compute total sum from local sums
1916       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1917       CALL MPI_ALLREDUCE( sums_l(nzb,8,0), sums(nzb,8), nzt+2-nzb, &
1918                           MPI_REAL, MPI_SUM, comm2d, ierr )
1919       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1920       CALL MPI_ALLREDUCE( sums_l(nzb,30,0), sums(nzb,30), nzt+2-nzb, &
1921                           MPI_REAL, MPI_SUM, comm2d, ierr )
1922       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1923       CALL MPI_ALLREDUCE( sums_l(nzb,31,0), sums(nzb,31), nzt+2-nzb, &
1924                           MPI_REAL, MPI_SUM, comm2d, ierr )
1925       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1926       CALL MPI_ALLREDUCE( sums_l(nzb,32,0), sums(nzb,32), nzt+2-nzb, &
1927                           MPI_REAL, MPI_SUM, comm2d, ierr )
1928
1929#else
1930       sums(:,8)  = sums_l(:,8,0)
1931       sums(:,30) = sums_l(:,30,0)
1932       sums(:,31) = sums_l(:,31,0)
1933       sums(:,32) = sums_l(:,32,0)
1934#endif
1935
1936!
1937!--    Final values are obtained by division by the total number of grid
1938!--    points used for the summation.
1939       hom(:,1,8,0)  = sums(:,8)  / ngp_2dh_outer(:,0)   ! e
1940       hom(:,1,30,0) = sums(:,30) / ngp_2dh_outer(:,0)   ! u*2
1941       hom(:,1,31,0) = sums(:,31) / ngp_2dh_outer(:,0)   ! v*2
1942       hom(:,1,32,0) = sums(:,32) / ngp_2dh_outer(:,0)   ! w*2
1943
1944    ENDIF
1945
1946 END SUBROUTINE lpm_init_sgs_tke
1947 
1948 
1949!------------------------------------------------------------------------------!
1950! Description:
1951! ------------
1952!> Sobroutine control lpm actions, i.e. all actions during one time step.
1953!------------------------------------------------------------------------------! 
1954 SUBROUTINE lpm_actions( location )
1955
1956    CHARACTER (LEN=*), INTENT(IN) ::  location !< call location string
1957
1958    INTEGER(iwp)       ::  i                  !<
1959    INTEGER(iwp)       ::  ie                 !<
1960    INTEGER(iwp)       ::  is                 !<
1961    INTEGER(iwp)       ::  j                  !<
1962    INTEGER(iwp)       ::  je                 !<
1963    INTEGER(iwp)       ::  js                 !<
1964    INTEGER(iwp), SAVE ::  lpm_count = 0      !<
1965    INTEGER(iwp)       ::  k                  !<
1966    INTEGER(iwp)       ::  ke                 !<
1967    INTEGER(iwp)       ::  ks                 !<
1968    INTEGER(iwp)       ::  m                  !<
1969    INTEGER(iwp), SAVE ::  steps = 0          !<
1970
1971    LOGICAL            ::  first_loop_stride  !<
1972
1973
1974    SELECT CASE ( location )
1975
1976       CASE ( 'after_prognostic_equations' )
1977
1978          CALL cpu_log( log_point(25), 'lpm', 'start' )
1979!
1980!--       Write particle data at current time on file.
1981!--       This has to be done here, before particles are further processed,
1982!--       because they may be deleted within this timestep (in case that
1983!--       dt_write_particle_data = dt_prel = particle_maximum_age).
1984          time_write_particle_data = time_write_particle_data + dt_3d
1985          IF ( time_write_particle_data >= dt_write_particle_data )  THEN
1986
1987             CALL lpm_data_output_particles
1988!
1989!--       The MOD function allows for changes in the output interval with restart
1990!--       runs.
1991             time_write_particle_data = MOD( time_write_particle_data, &
1992                                        MAX( dt_write_particle_data, dt_3d ) )
1993          ENDIF
1994
1995!
1996!--       Initialize arrays for marking those particles to be deleted after the
1997!--       (sub-) timestep
1998          deleted_particles = 0
1999
2000!
2001!--       Initialize variables used for accumulating the number of particles
2002!--       xchanged between the subdomains during all sub-timesteps (if sgs
2003!--       velocities are included). These data are output further below on the
2004!--       particle statistics file.
2005          trlp_count_sum      = 0
2006          trlp_count_recv_sum = 0
2007          trrp_count_sum      = 0
2008          trrp_count_recv_sum = 0
2009          trsp_count_sum      = 0
2010          trsp_count_recv_sum = 0
2011          trnp_count_sum      = 0
2012          trnp_count_recv_sum = 0
2013!
2014!--       Calculate exponential term used in case of particle inertia for each
2015!--       of the particle groups
2016          DO  m = 1, number_of_particle_groups
2017             IF ( particle_groups(m)%density_ratio /= 0.0_wp )  THEN
2018                particle_groups(m)%exp_arg  =                                        &
2019                          4.5_wp * particle_groups(m)%density_ratio *                &
2020                          molecular_viscosity / ( particle_groups(m)%radius )**2
2021
2022                particle_groups(m)%exp_term = EXP( -particle_groups(m)%exp_arg *     &
2023                          dt_3d )
2024             ENDIF
2025          ENDDO
2026!
2027!--       If necessary, release new set of particles
2028          IF ( ( simulated_time - last_particle_release_time ) >= dt_prel  .AND.     &
2029                 end_time_prel > simulated_time )  THEN
2030             DO WHILE ( ( simulated_time - last_particle_release_time ) >= dt_prel )
2031                CALL lpm_create_particle( PHASE_RELEASE )
2032                last_particle_release_time = last_particle_release_time + dt_prel
2033             ENDDO
2034          ENDIF
2035!
2036!--       Reset summation arrays
2037          IF ( cloud_droplets )  THEN
2038             ql_c  = 0.0_wp
2039             ql_v  = 0.0_wp
2040             ql_vp = 0.0_wp
2041          ENDIF
2042
2043          first_loop_stride = .TRUE.
2044          grid_particles(:,:,:)%time_loop_done = .TRUE.
2045!
2046!--       Timestep loop for particle advection.
2047!--       This loop has to be repeated until the advection time of every particle
2048!--       (within the total domain!) has reached the LES timestep (dt_3d).
2049!--       In case of including the SGS velocities, the particle timestep may be
2050!--       smaller than the LES timestep (because of the Lagrangian timescale
2051!--       restriction) and particles may require to undergo several particle
2052!--       timesteps, before the LES timestep is reached. Because the number of these
2053!--       particle timesteps to be carried out is unknown at first, these steps are
2054!--       carried out in the following infinite loop with exit condition.
2055          DO
2056             CALL cpu_log( log_point_s(44), 'lpm_advec', 'start' )
2057             CALL cpu_log( log_point_s(44), 'lpm_advec', 'pause' )
2058
2059!
2060!--          If particle advection includes SGS velocity components, calculate the
2061!--          required SGS quantities (i.e. gradients of the TKE, as well as
2062!--          horizontally averaged profiles of the SGS TKE and the resolved-scale
2063!--          velocity variances)
2064             IF ( use_sgs_for_particles  .AND.  .NOT. cloud_droplets )  THEN
2065                CALL lpm_init_sgs_tke
2066             ENDIF
2067!
2068!--          In case SGS-particle speed is considered, particles may carry out
2069!--          several particle timesteps. In order to prevent unnecessary
2070!--          treatment of particles that already reached the final time level,
2071!--          particles are sorted into contiguous blocks of finished and
2072!--          not-finished particles, in addition to their already sorting
2073!--          according to their sub-boxes.
2074             IF ( .NOT. first_loop_stride  .AND.  use_sgs_for_particles )            &
2075                CALL lpm_sort_timeloop_done
2076             DO  i = nxl, nxr
2077                DO  j = nys, nyn
2078                   DO  k = nzb+1, nzt
2079
2080                      number_of_particles = prt_count(k,j,i)
2081!
2082!--                   If grid cell gets empty, flag must be true
2083                      IF ( number_of_particles <= 0 )  THEN
2084                         grid_particles(k,j,i)%time_loop_done = .TRUE.
2085                         CYCLE
2086                      ENDIF
2087
2088                      IF ( .NOT. first_loop_stride  .AND.  &
2089                           grid_particles(k,j,i)%time_loop_done )  CYCLE
2090
2091                      particles => grid_particles(k,j,i)%particles(1:number_of_particles)
2092
2093                      particles(1:number_of_particles)%particle_mask = .TRUE.
2094!
2095!--                   Initialize the variable storing the total time that a particle
2096!--                   has advanced within the timestep procedure
2097                      IF ( first_loop_stride )  THEN
2098                         particles(1:number_of_particles)%dt_sum = 0.0_wp
2099                      ENDIF
2100!
2101!--                   Particle (droplet) growth by condensation/evaporation and
2102!--                   collision
2103                      IF ( cloud_droplets  .AND.  first_loop_stride)  THEN
2104!
2105!--                      Droplet growth by condensation / evaporation
2106                         CALL lpm_droplet_condensation(i,j,k)
2107!
2108!--                      Particle growth by collision
2109                         IF ( collision_kernel /= 'none' )  THEN
2110                            CALL lpm_droplet_collision(i,j,k)
2111                         ENDIF
2112
2113                      ENDIF
2114!
2115!--                   Initialize the switch used for the loop exit condition checked
2116!--                   at the end of this loop. If at least one particle has failed to
2117!--                   reach the LES timestep, this switch will be set false in
2118!--                   lpm_advec.
2119                      dt_3d_reached_l = .TRUE.
2120
2121!
2122!--                   Particle advection
2123                      CALL lpm_advec( i, j, k )
2124!
2125!--                   Particle reflection from walls. Only applied if the particles
2126!--                   are in the vertical range of the topography. (Here, some
2127!--                   optimization is still possible.)
2128                      IF ( topography /= 'flat'  .AND.  k < nzb_max + 2 )  THEN
2129                         CALL  lpm_boundary_conds( 'walls', i, j, k )
2130                      ENDIF
2131!
2132!--                   User-defined actions after the calculation of the new particle
2133!--                   position
2134                      CALL user_lpm_advec( i, j, k )
2135!
2136!--                   Apply boundary conditions to those particles that have crossed
2137!--                   the top or bottom boundary and delete those particles, which are
2138!--                   older than allowed
2139                      CALL lpm_boundary_conds( 'bottom/top', i, j, k )
2140!
2141!---                  If not all particles of the actual grid cell have reached the
2142!--                   LES timestep, this cell has to do another loop iteration. Due to
2143!--                   the fact that particles can move into neighboring grid cells,
2144!--                   these neighbor cells also have to perform another loop iteration.
2145!--                   Please note, this realization does not work properly if
2146!--                   particles move into another subdomain.
2147                      IF ( .NOT. dt_3d_reached_l )  THEN
2148                         ks = MAX(nzb+1,k-1)
2149                         ke = MIN(nzt,k+1)
2150                         js = MAX(nys,j-1)
2151                         je = MIN(nyn,j+1)
2152                         is = MAX(nxl,i-1)
2153                         ie = MIN(nxr,i+1)
2154                         grid_particles(ks:ke,js:je,is:ie)%time_loop_done = .FALSE.
2155                      ELSE
2156                         grid_particles(k,j,i)%time_loop_done = .TRUE.
2157                      ENDIF
2158
2159                   ENDDO
2160                ENDDO
2161             ENDDO
2162             steps = steps + 1
2163             dt_3d_reached_l = ALL(grid_particles(:,:,:)%time_loop_done)
2164!
2165!--          Find out, if all particles on every PE have completed the LES timestep
2166!--          and set the switch corespondingly
2167#if defined( __parallel )
2168             IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2169             CALL MPI_ALLREDUCE( dt_3d_reached_l, dt_3d_reached, 1, MPI_LOGICAL, &
2170                                 MPI_LAND, comm2d, ierr )
2171#else
2172             dt_3d_reached = dt_3d_reached_l
2173#endif
2174             CALL cpu_log( log_point_s(44), 'lpm_advec', 'stop' )
2175
2176!
2177!--          Apply splitting and merging algorithm
2178             IF ( cloud_droplets )  THEN
2179                IF ( splitting )  THEN
2180                   CALL lpm_splitting
2181                ENDIF
2182                IF ( merging )  THEN
2183                   CALL lpm_merging
2184                ENDIF
2185             ENDIF
2186!
2187!--          Move Particles local to PE to a different grid cell
2188             CALL lpm_move_particle
2189!
2190!--          Horizontal boundary conditions including exchange between subdmains
2191             CALL lpm_exchange_horiz
2192
2193!
2194!--          IF .FALSE., lpm_sort_and_delete is done inside pcmp
2195             IF ( .NOT. dt_3d_reached  .OR.  .NOT. nested_run )   THEN
2196!
2197!--             Pack particles (eliminate those marked for deletion),
2198!--             determine new number of particles
2199                CALL lpm_sort_and_delete
2200
2201!--             Initialize variables for the next (sub-) timestep, i.e., for marking
2202!--             those particles to be deleted after the timestep
2203                deleted_particles = 0
2204             ENDIF
2205
2206             IF ( dt_3d_reached )  EXIT
2207
2208             first_loop_stride = .FALSE.
2209          ENDDO   ! timestep loop
2210!
2211!--       in case of nested runs do the transfer of particles after every full model time step
2212          IF ( nested_run )   THEN
2213             CALL particles_from_parent_to_child
2214             CALL particles_from_child_to_parent
2215             CALL pmcp_p_delete_particles_in_fine_grid_area
2216
2217             CALL lpm_sort_and_delete
2218
2219             deleted_particles = 0
2220          ENDIF
2221
2222!
2223!--       Calculate the new liquid water content for each grid box
2224          IF ( cloud_droplets )  CALL lpm_calc_liquid_water_content
2225
2226!
2227!--       Deallocate unused memory
2228          IF ( deallocate_memory  .AND.  lpm_count == step_dealloc )  THEN
2229             CALL dealloc_particles_array
2230             lpm_count = 0
2231          ELSEIF ( deallocate_memory )  THEN
2232             lpm_count = lpm_count + 1
2233          ENDIF
2234
2235!
2236!--       Write particle statistics (in particular the number of particles
2237!--       exchanged between the subdomains) on file
2238          IF ( write_particle_statistics )  CALL lpm_write_exchange_statistics
2239
2240          CALL cpu_log( log_point(25), 'lpm', 'stop' )
2241
2242! !
2243! !--       Output of particle time series
2244!           IF ( particle_advection )  THEN
2245!              IF ( time_dopts >= dt_dopts  .OR.                                                        &
2246!                   ( time_since_reference_point >= particle_advection_start  .AND.                     &
2247!                    first_call_lpm ) )  THEN
2248!                 CALL lpm_data_output_ptseries
2249!                 time_dopts = MOD( time_dopts, MAX( dt_dopts, dt_3d ) )
2250!              ENDIF
2251!           ENDIF
2252
2253       CASE DEFAULT
2254          CONTINUE
2255
2256    END SELECT
2257
2258 END SUBROUTINE lpm_actions
2259 
2260 
2261!------------------------------------------------------------------------------!
2262! Description:
2263! ------------
2264!
2265!------------------------------------------------------------------------------!
2266 SUBROUTINE particles_from_parent_to_child
2267    IMPLICIT NONE
2268
2269    CALL pmcp_c_get_particle_from_parent                         ! Child actions
2270    CALL pmcp_p_fill_particle_win                                ! Parent actions
2271
2272    RETURN
2273 END SUBROUTINE particles_from_parent_to_child
2274
2275 
2276!------------------------------------------------------------------------------!
2277! Description:
2278! ------------
2279!
2280!------------------------------------------------------------------------------!
2281 SUBROUTINE particles_from_child_to_parent
2282    IMPLICIT NONE
2283
2284    CALL pmcp_c_send_particle_to_parent                         ! Child actions
2285    CALL pmcp_p_empty_particle_win                              ! Parent actions
2286
2287    RETURN
2288 END SUBROUTINE particles_from_child_to_parent
2289 
2290!------------------------------------------------------------------------------!
2291! Description:
2292! ------------
2293!> This routine write exchange statistics of the lpm in a ascii file.
2294!------------------------------------------------------------------------------!
2295 SUBROUTINE lpm_write_exchange_statistics
2296
2297    INTEGER(iwp) ::  ip         !<
2298    INTEGER(iwp) ::  jp         !<
2299    INTEGER(iwp) ::  kp         !<
2300    INTEGER(iwp) ::  tot_number_of_particles !<
2301
2302!
2303!-- Determine the current number of particles
2304    number_of_particles         = 0
2305    DO  ip = nxl, nxr
2306       DO  jp = nys, nyn
2307          DO  kp = nzb+1, nzt
2308             number_of_particles = number_of_particles                         &
2309                                     + prt_count(kp,jp,ip)
2310          ENDDO
2311       ENDDO
2312    ENDDO
2313
2314    CALL check_open( 80 )
2315#if defined( __parallel )
2316    WRITE ( 80, 8000 )  current_timestep_number+1, simulated_time+dt_3d, &
2317                        number_of_particles, pleft, trlp_count_sum,      &
2318                        trlp_count_recv_sum, pright, trrp_count_sum,     &
2319                        trrp_count_recv_sum, psouth, trsp_count_sum,     &
2320                        trsp_count_recv_sum, pnorth, trnp_count_sum,     &
2321                        trnp_count_recv_sum
2322#else
2323    WRITE ( 80, 8000 )  current_timestep_number+1, simulated_time+dt_3d, &
2324                        number_of_particles
2325#endif
2326    CALL close_file( 80 )
2327
2328    IF ( number_of_particles > 0 )  THEN
2329        WRITE(9,*) 'number_of_particles ', number_of_particles,                &
2330                    current_timestep_number + 1, simulated_time + dt_3d
2331    ENDIF
2332
2333#if defined( __parallel )
2334    CALL MPI_ALLREDUCE( number_of_particles, tot_number_of_particles, 1,       &
2335                        MPI_INTEGER, MPI_SUM, comm2d, ierr )
2336#else
2337    tot_number_of_particles = number_of_particles
2338#endif
2339
2340    IF ( nested_run )  THEN
2341       CALL pmcp_g_print_number_of_particles( simulated_time+dt_3d,            &
2342                                              tot_number_of_particles)
2343    ENDIF
2344
2345!
2346!-- Formats
23478000 FORMAT (I6,1X,F7.2,4X,I10,5X,4(I3,1X,I4,'/',I4,2X),6X,I10)
2348
2349
2350 END SUBROUTINE lpm_write_exchange_statistics
2351 
2352
2353!------------------------------------------------------------------------------!
2354! Description:
2355! ------------
2356!> Write particle data in FORTRAN binary and/or netCDF format
2357!------------------------------------------------------------------------------!
2358 SUBROUTINE lpm_data_output_particles
2359 
2360    INTEGER(iwp) ::  ip !<
2361    INTEGER(iwp) ::  jp !<
2362    INTEGER(iwp) ::  kp !<
2363
2364    CALL cpu_log( log_point_s(40), 'lpm_data_output', 'start' )
2365
2366!
2367!-- Attention: change version number for unit 85 (in routine check_open)
2368!--            whenever the output format for this unit is changed!
2369    CALL check_open( 85 )
2370
2371    WRITE ( 85 )  simulated_time
2372    WRITE ( 85 )  prt_count
2373
2374    DO  ip = nxl, nxr
2375       DO  jp = nys, nyn
2376          DO  kp = nzb+1, nzt
2377             number_of_particles = prt_count(kp,jp,ip)
2378             particles => grid_particles(kp,jp,ip)%particles(1:number_of_particles)
2379             IF ( number_of_particles <= 0 )  CYCLE
2380             WRITE ( 85 )  particles
2381          ENDDO
2382       ENDDO
2383    ENDDO
2384
2385    CALL close_file( 85 )
2386
2387
2388#if defined( __netcdf )
2389! !
2390! !-- Output in netCDF format
2391!     CALL check_open( 108 )
2392!
2393! !
2394! !-- Update the NetCDF time axis
2395!     prt_time_count = prt_time_count + 1
2396!
2397!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_time_prt, &
2398!                             (/ simulated_time /),        &
2399!                             start = (/ prt_time_count /), count = (/ 1 /) )
2400!     CALL netcdf_handle_error( 'lpm_data_output_particles', 1 )
2401!
2402! !
2403! !-- Output the real number of particles used
2404!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_rnop_prt, &
2405!                             (/ number_of_particles /),   &
2406!                             start = (/ prt_time_count /), count = (/ 1 /) )
2407!     CALL netcdf_handle_error( 'lpm_data_output_particles', 2 )
2408!
2409! !
2410! !-- Output all particle attributes
2411!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(1), particles%age,      &
2412!                             start = (/ 1, prt_time_count /),               &
2413!                             count = (/ maximum_number_of_particles /) )
2414!     CALL netcdf_handle_error( 'lpm_data_output_particles', 3 )
2415!
2416!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(2), particles%user,     &
2417!                             start = (/ 1, prt_time_count /),               &
2418!                             count = (/ maximum_number_of_particles /) )
2419!     CALL netcdf_handle_error( 'lpm_data_output_particles', 4 )
2420!
2421!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(3), particles%origin_x, &
2422!                             start = (/ 1, prt_time_count /),               &
2423!                             count = (/ maximum_number_of_particles /) )
2424!     CALL netcdf_handle_error( 'lpm_data_output_particles', 5 )
2425!
2426!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(4), particles%origin_y, &
2427!                             start = (/ 1, prt_time_count /),               &
2428!                             count = (/ maximum_number_of_particles /) )
2429!     CALL netcdf_handle_error( 'lpm_data_output_particles', 6 )
2430!
2431!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(5), particles%origin_z, &
2432!                             start = (/ 1, prt_time_count /),               &
2433!                             count = (/ maximum_number_of_particles /) )
2434!     CALL netcdf_handle_error( 'lpm_data_output_particles', 7 )
2435!
2436!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(6), particles%radius,   &
2437!                             start = (/ 1, prt_time_count /),               &
2438!                             count = (/ maximum_number_of_particles /) )
2439!     CALL netcdf_handle_error( 'lpm_data_output_particles', 8 )
2440!
2441!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(7), particles%speed_x,  &
2442!                             start = (/ 1, prt_time_count /),               &
2443!                             count = (/ maximum_number_of_particles /) )
2444!     CALL netcdf_handle_error( 'lpm_data_output_particles', 9 )
2445!
2446!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(8), particles%speed_y,  &
2447!                             start = (/ 1, prt_time_count /),               &
2448!                             count = (/ maximum_number_of_particles /) )
2449!     CALL netcdf_handle_error( 'lpm_data_output_particles', 10 )
2450!
2451!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(9), particles%speed_z,  &
2452!                             start = (/ 1, prt_time_count /),               &
2453!                             count = (/ maximum_number_of_particles /) )
2454!     CALL netcdf_handle_error( 'lpm_data_output_particles', 11 )
2455!
2456!     nc_stat = NF90_PUT_VAR( id_set_prt,id_var_prt(10),                     &
2457!                             particles%weight_factor,                       &
2458!                             start = (/ 1, prt_time_count /),               &
2459!                             count = (/ maximum_number_of_particles /) )
2460!     CALL netcdf_handle_error( 'lpm_data_output_particles', 12 )
2461!
2462!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(11), particles%x,       &
2463!                             start = (/ 1, prt_time_count /),               &
2464!                             count = (/ maximum_number_of_particles /) )
2465!     CALL netcdf_handle_error( 'lpm_data_output_particles', 13 )
2466!
2467!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(12), particles%y,       &
2468!                             start = (/ 1, prt_time_count /),               &
2469!                             count = (/ maximum_number_of_particles /) )
2470!     CALL netcdf_handle_error( 'lpm_data_output_particles', 14 )
2471!
2472!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(13), particles%z,       &
2473!                             start = (/ 1, prt_time_count /),               &
2474!                             count = (/ maximum_number_of_particles /) )
2475!     CALL netcdf_handle_error( 'lpm_data_output_particles', 15 )
2476!
2477!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(14), particles%class,   &
2478!                             start = (/ 1, prt_time_count /),               &
2479!                             count = (/ maximum_number_of_particles /) )
2480!     CALL netcdf_handle_error( 'lpm_data_output_particles', 16 )
2481!
2482!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(15), particles%group,   &
2483!                             start = (/ 1, prt_time_count /),               &
2484!                             count = (/ maximum_number_of_particles /) )
2485!     CALL netcdf_handle_error( 'lpm_data_output_particles', 17 )
2486!
2487!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(16),                    &
2488!                             particles%id2,                                 &
2489!                             start = (/ 1, prt_time_count /),               &
2490!                             count = (/ maximum_number_of_particles /) )
2491!     CALL netcdf_handle_error( 'lpm_data_output_particles', 18 )
2492!
2493!     nc_stat = NF90_PUT_VAR( id_set_prt, id_var_prt(17), particles%id1,     &
2494!                             start = (/ 1, prt_time_count /),               &
2495!                             count = (/ maximum_number_of_particles /) )
2496!     CALL netcdf_handle_error( 'lpm_data_output_particles', 19 )
2497!
2498#endif
2499
2500    CALL cpu_log( log_point_s(40), 'lpm_data_output', 'stop' )
2501
2502 END SUBROUTINE lpm_data_output_particles
2503 
2504!------------------------------------------------------------------------------!
2505! Description:
2506! ------------
2507!> This routine calculates and provide particle timeseries output.
2508!------------------------------------------------------------------------------!
2509 SUBROUTINE lpm_data_output_ptseries
2510 
2511    INTEGER(iwp) ::  i    !<
2512    INTEGER(iwp) ::  inum !<
2513    INTEGER(iwp) ::  j    !<
2514    INTEGER(iwp) ::  jg   !<
2515    INTEGER(iwp) ::  k    !<
2516    INTEGER(iwp) ::  n    !<
2517
2518    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  pts_value   !<
2519    REAL(wp), DIMENSION(:,:), ALLOCATABLE ::  pts_value_l !<
2520
2521
2522    CALL cpu_log( log_point(36), 'data_output_ptseries', 'start' )
2523
2524    IF ( myid == 0 )  THEN
2525!
2526!--    Open file for time series output in NetCDF format
2527       dopts_time_count = dopts_time_count + 1
2528       CALL check_open( 109 )
2529#if defined( __netcdf )
2530!
2531!--    Update the particle time series time axis
2532       nc_stat = NF90_PUT_VAR( id_set_pts, id_var_time_pts,      &
2533                               (/ time_since_reference_point /), &
2534                               start = (/ dopts_time_count /), count = (/ 1 /) )
2535       CALL netcdf_handle_error( 'data_output_ptseries', 391 )
2536#endif
2537
2538    ENDIF
2539
2540    ALLOCATE( pts_value(0:number_of_particle_groups,dopts_num), &
2541              pts_value_l(0:number_of_particle_groups,dopts_num) )
2542
2543    pts_value_l = 0.0_wp
2544    pts_value_l(:,16) = 9999999.9_wp    ! for calculation of minimum radius
2545
2546!
2547!-- Calculate or collect the particle time series quantities for all particles
2548!-- and seperately for each particle group (if there is more than one group)
2549    DO  i = nxl, nxr
2550       DO  j = nys, nyn
2551          DO  k = nzb, nzt
2552             number_of_particles = prt_count(k,j,i)
2553             IF (number_of_particles <= 0)  CYCLE
2554             particles => grid_particles(k,j,i)%particles(1:number_of_particles)
2555             DO  n = 1, number_of_particles
2556
2557                IF ( particles(n)%particle_mask )  THEN  ! Restrict analysis to active particles
2558
2559                   pts_value_l(0,1)  = pts_value_l(0,1) + 1.0_wp  ! total # of particles
2560                   pts_value_l(0,2)  = pts_value_l(0,2) +                      &
2561                          ( particles(n)%x - particles(n)%origin_x )  ! mean x
2562                   pts_value_l(0,3)  = pts_value_l(0,3) +                      &
2563                          ( particles(n)%y - particles(n)%origin_y )  ! mean y
2564                   pts_value_l(0,4)  = pts_value_l(0,4) +                      &
2565                          ( particles(n)%z - particles(n)%origin_z )  ! mean z
2566                   pts_value_l(0,5)  = pts_value_l(0,5) + particles(n)%z        ! mean z (absolute)
2567                   pts_value_l(0,6)  = pts_value_l(0,6) + particles(n)%speed_x  ! mean u
2568                   pts_value_l(0,7)  = pts_value_l(0,7) + particles(n)%speed_y  ! mean v
2569                   pts_value_l(0,8)  = pts_value_l(0,8) + particles(n)%speed_z  ! mean w
2570                   pts_value_l(0,9)  = pts_value_l(0,9)  + particles(n)%rvar1 ! mean sgsu
2571                   pts_value_l(0,10) = pts_value_l(0,10) + particles(n)%rvar2 ! mean sgsv
2572                   pts_value_l(0,11) = pts_value_l(0,11) + particles(n)%rvar3 ! mean sgsw
2573                   IF ( particles(n)%speed_z > 0.0_wp )  THEN
2574                      pts_value_l(0,12) = pts_value_l(0,12) + 1.0_wp  ! # of upward moving prts
2575                      pts_value_l(0,13) = pts_value_l(0,13) +                  &
2576                                              particles(n)%speed_z ! mean w upw.
2577                   ELSE
2578                      pts_value_l(0,14) = pts_value_l(0,14) +                  &
2579                                              particles(n)%speed_z ! mean w down
2580                   ENDIF
2581                   pts_value_l(0,15) = pts_value_l(0,15) + particles(n)%radius ! mean rad
2582                   pts_value_l(0,16) = MIN( pts_value_l(0,16), particles(n)%radius ) ! minrad
2583                   pts_value_l(0,17) = MAX( pts_value_l(0,17), particles(n)%radius ) ! maxrad
2584                   pts_value_l(0,18) = pts_value_l(0,18) + 1.0_wp
2585                   pts_value_l(0,19) = pts_value_l(0,18) + 1.0_wp
2586!
2587!--                Repeat the same for the respective particle group
2588                   IF ( number_of_particle_groups > 1 )  THEN
2589                      jg = particles(n)%group
2590
2591                      pts_value_l(jg,1)  = pts_value_l(jg,1) + 1.0_wp
2592                      pts_value_l(jg,2)  = pts_value_l(jg,2) +                   &
2593                           ( particles(n)%x - particles(n)%origin_x )
2594                      pts_value_l(jg,3)  = pts_value_l(jg,3) +                   &
2595                           ( particles(n)%y - particles(n)%origin_y )
2596                      pts_value_l(jg,4)  = pts_value_l(jg,4) +                   &
2597                           ( particles(n)%z - particles(n)%origin_z )
2598                      pts_value_l(jg,5)  = pts_value_l(jg,5) + particles(n)%z
2599                      pts_value_l(jg,6)  = pts_value_l(jg,6) + particles(n)%speed_x
2600                      pts_value_l(jg,7)  = pts_value_l(jg,7) + particles(n)%speed_y
2601                      pts_value_l(jg,8)  = pts_value_l(jg,8) + particles(n)%speed_z
2602                      pts_value_l(jg,9)  = pts_value_l(jg,9)  + particles(n)%rvar1
2603                      pts_value_l(jg,10) = pts_value_l(jg,10) + particles(n)%rvar2
2604                      pts_value_l(jg,11) = pts_value_l(jg,11) + particles(n)%rvar3
2605                      IF ( particles(n)%speed_z > 0.0_wp )  THEN
2606                         pts_value_l(jg,12) = pts_value_l(jg,12) + 1.0_wp
2607                         pts_value_l(jg,13) = pts_value_l(jg,13) + particles(n)%speed_z
2608                      ELSE
2609                         pts_value_l(jg,14) = pts_value_l(jg,14) + particles(n)%speed_z
2610                      ENDIF
2611                      pts_value_l(jg,15) = pts_value_l(jg,15) + particles(n)%radius
2612                      pts_value_l(jg,16) = MIN( pts_value(jg,16), particles(n)%radius )
2613                      pts_value_l(jg,17) = MAX( pts_value(jg,17), particles(n)%radius )
2614                      pts_value_l(jg,18) = pts_value_l(jg,18) + 1.0_wp
2615                      pts_value_l(jg,19) = pts_value_l(jg,19) + 1.0_wp
2616                   ENDIF
2617
2618                ENDIF
2619
2620             ENDDO
2621
2622          ENDDO
2623       ENDDO
2624    ENDDO
2625
2626
2627#if defined( __parallel )
2628!
2629!-- Sum values of the subdomains
2630    inum = number_of_particle_groups + 1
2631
2632    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2633    CALL MPI_ALLREDUCE( pts_value_l(0,1), pts_value(0,1), 15*inum, MPI_REAL, &
2634                        MPI_SUM, comm2d, ierr )
2635    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2636    CALL MPI_ALLREDUCE( pts_value_l(0,16), pts_value(0,16), inum, MPI_REAL, &
2637                        MPI_MIN, comm2d, ierr )
2638    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2639    CALL MPI_ALLREDUCE( pts_value_l(0,17), pts_value(0,17), inum, MPI_REAL, &
2640                        MPI_MAX, comm2d, ierr )
2641    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2642    CALL MPI_ALLREDUCE( pts_value_l(0,18), pts_value(0,18), inum, MPI_REAL, &
2643                        MPI_MAX, comm2d, ierr )
2644    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2645    CALL MPI_ALLREDUCE( pts_value_l(0,19), pts_value(0,19), inum, MPI_REAL, &
2646                        MPI_MIN, comm2d, ierr )
2647#else
2648    pts_value(:,1:19) = pts_value_l(:,1:19)
2649#endif
2650
2651!
2652!-- Normalize the above calculated quantities (except min/max values) with the
2653!-- total number of particles
2654    IF ( number_of_particle_groups > 1 )  THEN
2655       inum = number_of_particle_groups
2656    ELSE
2657       inum = 0
2658    ENDIF
2659
2660    DO  j = 0, inum
2661
2662       IF ( pts_value(j,1) > 0.0_wp )  THEN
2663
2664          pts_value(j,2:15) = pts_value(j,2:15) / pts_value(j,1)
2665          IF ( pts_value(j,12) > 0.0_wp  .AND.  pts_value(j,12) < 1.0_wp )  THEN
2666             pts_value(j,13) = pts_value(j,13) / pts_value(j,12)
2667             pts_value(j,14) = pts_value(j,14) / ( 1.0_wp - pts_value(j,12) )
2668          ELSEIF ( pts_value(j,12) == 0.0_wp )  THEN
2669             pts_value(j,13) = -1.0_wp
2670          ELSE
2671             pts_value(j,14) = -1.0_wp
2672          ENDIF
2673
2674       ENDIF
2675
2676    ENDDO
2677
2678!
2679!-- Calculate higher order moments of particle time series quantities,
2680!-- seperately for each particle group (if there is more than one group)
2681    DO  i = nxl, nxr
2682       DO  j = nys, nyn
2683          DO  k = nzb, nzt
2684             number_of_particles = prt_count(k,j,i)
2685             IF (number_of_particles <= 0)  CYCLE
2686             particles => grid_particles(k,j,i)%particles(1:number_of_particles)
2687             DO  n = 1, number_of_particles
2688
2689                pts_value_l(0,20) = pts_value_l(0,20) + ( particles(n)%x - &
2690                                    particles(n)%origin_x - pts_value(0,2) )**2 ! x*2
2691                pts_value_l(0,21) = pts_value_l(0,21) + ( particles(n)%y - &
2692                                    particles(n)%origin_y - pts_value(0,3) )**2 ! y*2
2693                pts_value_l(0,22) = pts_value_l(0,22) + ( particles(n)%z - &
2694                                    particles(n)%origin_z - pts_value(0,4) )**2 ! z*2
2695                pts_value_l(0,23) = pts_value_l(0,23) + ( particles(n)%speed_x - &
2696                                                         pts_value(0,6) )**2   ! u*2
2697                pts_value_l(0,24) = pts_value_l(0,24) + ( particles(n)%speed_y - &
2698                                                          pts_value(0,7) )**2   ! v*2
2699                pts_value_l(0,25) = pts_value_l(0,25) + ( particles(n)%speed_z - &
2700                                                          pts_value(0,8) )**2   ! w*2
2701                pts_value_l(0,26) = pts_value_l(0,26) + ( particles(n)%rvar1 - &
2702                                                          pts_value(0,9) )**2   ! u"2
2703                pts_value_l(0,27) = pts_value_l(0,27) + ( particles(n)%rvar2 - &
2704                                                          pts_value(0,10) )**2  ! v"2
2705                pts_value_l(0,28) = pts_value_l(0,28) + ( particles(n)%rvar3 - &
2706                                                          pts_value(0,11) )**2  ! w"2
2707!
2708!--             Repeat the same for the respective particle group
2709                IF ( number_of_particle_groups > 1 )  THEN
2710                   jg = particles(n)%group
2711
2712                   pts_value_l(jg,20) = pts_value_l(jg,20) + ( particles(n)%x - &
2713                                       particles(n)%origin_x - pts_value(jg,2) )**2
2714                   pts_value_l(jg,21) = pts_value_l(jg,21) + ( particles(n)%y - &
2715                                       particles(n)%origin_y - pts_value(jg,3) )**2
2716                   pts_value_l(jg,22) = pts_value_l(jg,22) + ( particles(n)%z - &
2717                                       particles(n)%origin_z - pts_value(jg,4) )**2
2718                   pts_value_l(jg,23) = pts_value_l(jg,23) + ( particles(n)%speed_x - &
2719                                                             pts_value(jg,6) )**2
2720                   pts_value_l(jg,24) = pts_value_l(jg,24) + ( particles(n)%speed_y - &
2721                                                             pts_value(jg,7) )**2
2722                   pts_value_l(jg,25) = pts_value_l(jg,25) + ( particles(n)%speed_z - &
2723                                                             pts_value(jg,8) )**2
2724                   pts_value_l(jg,26) = pts_value_l(jg,26) + ( particles(n)%rvar1 - &
2725                                                             pts_value(jg,9) )**2
2726                   pts_value_l(jg,27) = pts_value_l(jg,27) + ( particles(n)%rvar2 - &
2727                                                             pts_value(jg,10) )**2
2728                   pts_value_l(jg,28) = pts_value_l(jg,28) + ( particles(n)%rvar3 - &
2729                                                             pts_value(jg,11) )**2
2730                ENDIF
2731
2732             ENDDO
2733          ENDDO
2734       ENDDO
2735    ENDDO
2736
2737    pts_value_l(0,29) = ( number_of_particles - pts_value(0,1) / numprocs )**2
2738                                                 ! variance of particle numbers
2739    IF ( number_of_particle_groups > 1 )  THEN
2740       DO  j = 1, number_of_particle_groups
2741          pts_value_l(j,29) = ( pts_value_l(j,1) - &
2742                                pts_value(j,1) / numprocs )**2
2743       ENDDO
2744    ENDIF
2745
2746#if defined( __parallel )
2747!
2748!-- Sum values of the subdomains
2749    inum = number_of_particle_groups + 1
2750
2751    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
2752    CALL MPI_ALLREDUCE( pts_value_l(0,20), pts_value(0,20), inum*10, MPI_REAL, &
2753                        MPI_SUM, comm2d, ierr )
2754#else
2755    pts_value(:,20:29) = pts_value_l(:,20:29)
2756#endif
2757
2758!
2759!-- Normalize the above calculated quantities with the total number of
2760!-- particles
2761    IF ( number_of_particle_groups > 1 )  THEN
2762       inum = number_of_particle_groups
2763    ELSE
2764       inum = 0
2765    ENDIF
2766
2767    DO  j = 0, inum
2768
2769       IF ( pts_value(j,1) > 0.0_wp )  THEN
2770          pts_value(j,20:28) = pts_value(j,20:28) / pts_value(j,1)
2771       ENDIF
2772       pts_value(j,29) = pts_value(j,29) / numprocs
2773
2774    ENDDO
2775
2776#if defined( __netcdf )
2777!
2778!-- Output particle time series quantities in NetCDF format
2779    IF ( myid == 0 )  THEN
2780       DO  j = 0, inum
2781          DO  i = 1, dopts_num
2782             nc_stat = NF90_PUT_VAR( id_set_pts, id_var_dopts(i,j),  &
2783                                     (/ pts_value(j,i) /),           &
2784                                     start = (/ dopts_time_count /), &
2785                                     count = (/ 1 /) )
2786             CALL netcdf_handle_error( 'data_output_ptseries', 392 )
2787          ENDDO
2788       ENDDO
2789    ENDIF
2790#endif
2791
2792    DEALLOCATE( pts_value, pts_value_l )
2793
2794    CALL cpu_log( log_point(36), 'data_output_ptseries', 'stop' )
2795
2796END SUBROUTINE lpm_data_output_ptseries
2797
2798 
2799!------------------------------------------------------------------------------!
2800! Description:
2801! ------------
2802!> This routine reads the respective restart data for the lpm.
2803!------------------------------------------------------------------------------!
2804 SUBROUTINE lpm_rrd_local_particles
2805
2806    CHARACTER (LEN=10) ::  particle_binary_version    !<
2807    CHARACTER (LEN=10) ::  version_on_file            !<
2808
2809    INTEGER(iwp) ::  alloc_size !<
2810    INTEGER(iwp) ::  ip         !<
2811    INTEGER(iwp) ::  jp         !<
2812    INTEGER(iwp) ::  kp         !<
2813
2814    TYPE(particle_type), DIMENSION(:), ALLOCATABLE ::  tmp_particles !<
2815
2816!
2817!-- Read particle data from previous model run.
2818!-- First open the input unit.
2819    IF ( myid_char == '' )  THEN
2820       OPEN ( 90, FILE='PARTICLE_RESTART_DATA_IN'//myid_char,                  &
2821                  FORM='UNFORMATTED' )
2822    ELSE
2823       OPEN ( 90, FILE='PARTICLE_RESTART_DATA_IN/'//myid_char,                 &
2824                  FORM='UNFORMATTED' )
2825    ENDIF
2826
2827!
2828!-- First compare the version numbers
2829    READ ( 90 )  version_on_file
2830    particle_binary_version = '4.0'
2831    IF ( TRIM( version_on_file ) /= TRIM( particle_binary_version ) )  THEN
2832       message_string = 'version mismatch concerning data from prior ' //      &
2833                        'run &version on file = "' //                          &
2834                                      TRIM( version_on_file ) //               &
2835                        '&version in program = "' //                           &
2836                                      TRIM( particle_binary_version ) // '"'
2837       CALL message( 'lpm_read_restart_file', 'PA0214', 1, 2, 0, 6, 0 )
2838    ENDIF
2839
2840!
2841!-- If less particles are stored on the restart file than prescribed by
2842!-- 1, the remainder is initialized by zero_particle to avoid
2843!-- errors.
2844    zero_particle = particle_type( 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp,     &
2845                                   0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp,     &
2846                                   0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp,     &
2847                                   0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp,     &
2848                                   0, 0, 0_idp, .FALSE., -1 )
2849!
2850!-- Read some particle parameters and the size of the particle arrays,
2851!-- allocate them and read their contents.
2852    READ ( 90 )  bc_par_b, bc_par_lr, bc_par_ns, bc_par_t,                     &
2853                 last_particle_release_time, number_of_particle_groups,        &
2854                 particle_groups, time_write_particle_data
2855
2856    ALLOCATE( prt_count(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                        &
2857              grid_particles(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
2858
2859    READ ( 90 )  prt_count
2860
2861    DO  ip = nxl, nxr
2862       DO  jp = nys, nyn
2863          DO  kp = nzb+1, nzt
2864
2865             number_of_particles = prt_count(kp,jp,ip)
2866             IF ( number_of_particles > 0 )  THEN
2867                alloc_size = MAX( INT( number_of_particles *                   &
2868                             ( 1.0_wp + alloc_factor / 100.0_wp ) ),           &
2869                             1 )
2870             ELSE
2871                alloc_size = 1
2872             ENDIF
2873
2874             ALLOCATE( grid_particles(kp,jp,ip)%particles(1:alloc_size) )
2875
2876             IF ( number_of_particles > 0 )  THEN
2877                ALLOCATE( tmp_particles(1:number_of_particles) )
2878                READ ( 90 )  tmp_particles
2879                grid_particles(kp,jp,ip)%particles(1:number_of_particles) = tmp_particles
2880                DEALLOCATE( tmp_particles )
2881                IF ( number_of_particles < alloc_size )  THEN
2882                   grid_particles(kp,jp,ip)%particles(number_of_particles+1:alloc_size) &
2883                      = zero_particle
2884                ENDIF
2885             ELSE
2886                grid_particles(kp,jp,ip)%particles(1:alloc_size) = zero_particle
2887             ENDIF
2888
2889          ENDDO
2890       ENDDO
2891    ENDDO
2892
2893    CLOSE ( 90 )
2894!
2895!-- Must be called to sort particles into blocks, which is needed for a fast
2896!-- interpolation of the LES fields on the particle position.
2897    CALL lpm_sort_and_delete
2898
2899
2900 END SUBROUTINE lpm_rrd_local_particles
2901 
2902 
2903 SUBROUTINE lpm_rrd_local( k, nxlf, nxlc, nxl_on_file, nxrf, nxrc,          &
2904                              nxr_on_file, nynf, nync, nyn_on_file, nysf,  &
2905                              nysc, nys_on_file, tmp_3d, found )
2906
2907
2908   USE control_parameters,                                                 &
2909       ONLY: length, restart_string
2910
2911    INTEGER(iwp) ::  k               !<
2912    INTEGER(iwp) ::  nxlc            !<
2913    INTEGER(iwp) ::  nxlf            !<
2914    INTEGER(iwp) ::  nxl_on_file     !<
2915    INTEGER(iwp) ::  nxrc            !<
2916    INTEGER(iwp) ::  nxrf            !<
2917    INTEGER(iwp) ::  nxr_on_file     !<
2918    INTEGER(iwp) ::  nync            !<
2919    INTEGER(iwp) ::  nynf            !<
2920    INTEGER(iwp) ::  nyn_on_file     !<
2921    INTEGER(iwp) ::  nysc            !<
2922    INTEGER(iwp) ::  nysf            !<
2923    INTEGER(iwp) ::  nys_on_file     !<
2924
2925    LOGICAL, INTENT(OUT)  ::  found
2926
2927    REAL(wp), DIMENSION(nzb:nzt+1,nys_on_file-nbgp:nyn_on_file+nbgp,nxl_on_file-nbgp:nxr_on_file+nbgp) ::  tmp_3d   !<
2928
2929
2930    found = .TRUE.
2931
2932    SELECT CASE ( restart_string(1:length) )
2933
2934       CASE ( 'iran' ) ! matching random numbers is still unresolved issue
2935          IF ( k == 1 )  READ ( 13 )  iran, iran_part
2936
2937        CASE ( 'pc_av' )
2938           IF ( .NOT. ALLOCATED( pc_av ) )  THEN
2939              ALLOCATE( pc_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
2940           ENDIF
2941           IF ( k == 1 )  READ ( 13 )  tmp_3d
2942           pc_av(:,nysc-nbgp:nync+nbgp,nxlc-nbgp:nxrc+nbgp) =          &
2943              tmp_3d(:,nysf-nbgp:nynf+nbgp,nxlf-nbgp:nxrf+nbgp)
2944
2945        CASE ( 'pr_av' )
2946           IF ( .NOT. ALLOCATED( pr_av ) )  THEN
2947              ALLOCATE( pr_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
2948           ENDIF
2949           IF ( k == 1 )  READ ( 13 )  tmp_3d
2950           pr_av(:,nysc-nbgp:nync+nbgp,nxlc-nbgp:nxrc+nbgp) =          &
2951              tmp_3d(:,nysf-nbgp:nynf+nbgp,nxlf-nbgp:nxrf+nbgp)
2952 
2953         CASE ( 'ql_c_av' )
2954            IF ( .NOT. ALLOCATED( ql_c_av ) )  THEN
2955               ALLOCATE( ql_c_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
2956            ENDIF
2957            IF ( k == 1 )  READ ( 13 )  tmp_3d
2958            ql_c_av(:,nysc-nbgp:nync+nbgp,nxlc-nbgp:nxrc+nbgp) =        &
2959               tmp_3d(:,nysf-nbgp:nynf+nbgp,nxlf-nbgp:nxrf+nbgp)
2960
2961         CASE ( 'ql_v_av' )
2962            IF ( .NOT. ALLOCATED( ql_v_av ) )  THEN
2963               ALLOCATE( ql_v_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
2964            ENDIF
2965            IF ( k == 1 )  READ ( 13 )  tmp_3d
2966            ql_v_av(:,nysc-nbgp:nync+nbgp,nxlc-nbgp:nxrc+nbgp) =        &
2967               tmp_3d(:,nysf-nbgp:nynf+nbgp,nxlf-nbgp:nxrf+nbgp)
2968
2969         CASE ( 'ql_vp_av' )
2970            IF ( .NOT. ALLOCATED( ql_vp_av ) )  THEN
2971               ALLOCATE( ql_vp_av(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
2972            ENDIF
2973            IF ( k == 1 )  READ ( 13 )  tmp_3d
2974            ql_vp_av(:,nysc-nbgp:nync+nbgp,nxlc-nbgp:nxrc+nbgp) =       &
2975               tmp_3d(:,nysf-nbgp:nynf+nbgp,nxlf-nbgp:nxrf+nbgp)
2976
2977          CASE DEFAULT
2978
2979             found = .FALSE.
2980
2981       END SELECT
2982               
2983
2984 END SUBROUTINE lpm_rrd_local
2985 
2986!------------------------------------------------------------------------------!
2987! Description:
2988! ------------
2989!> This routine writes the respective restart data for the lpm.
2990!------------------------------------------------------------------------------!
2991 SUBROUTINE lpm_wrd_local
2992 
2993    CHARACTER (LEN=10) ::  particle_binary_version   !<
2994
2995    INTEGER(iwp) ::  ip                              !<
2996    INTEGER(iwp) ::  jp                              !<
2997    INTEGER(iwp) ::  kp                              !<
2998!
2999!-- First open the output unit.
3000    IF ( myid_char == '' )  THEN
3001       OPEN ( 90, FILE='PARTICLE_RESTART_DATA_OUT'//myid_char, &
3002                  FORM='UNFORMATTED')
3003    ELSE
3004       IF ( myid == 0 )  CALL local_system( 'mkdir PARTICLE_RESTART_DATA_OUT' )
3005#if defined( __parallel )
3006!
3007!--    Set a barrier in order to allow that thereafter all other processors
3008!--    in the directory created by PE0 can open their file
3009       CALL MPI_BARRIER( comm2d, ierr )
3010#endif
3011       OPEN ( 90, FILE='PARTICLE_RESTART_DATA_OUT/'//myid_char, &
3012                  FORM='UNFORMATTED' )
3013    ENDIF
3014
3015!
3016!-- Write the version number of the binary format.
3017!-- Attention: After changes to the following output commands the version
3018!-- ---------  number of the variable particle_binary_version must be
3019!--            changed! Also, the version number and the list of arrays
3020!--            to be read in lpm_read_restart_file must be adjusted
3021!--            accordingly.
3022    particle_binary_version = '4.0'
3023    WRITE ( 90 )  particle_binary_version
3024
3025!
3026!-- Write some particle parameters, the size of the particle arrays
3027    WRITE ( 90 )  bc_par_b, bc_par_lr, bc_par_ns, bc_par_t,                    &
3028                  last_particle_release_time, number_of_particle_groups,       &
3029                  particle_groups, time_write_particle_data
3030
3031    WRITE ( 90 )  prt_count
3032         
3033    DO  ip = nxl, nxr
3034       DO  jp = nys, nyn
3035          DO  kp = nzb+1, nzt
3036             number_of_particles = prt_count(kp,jp,ip)
3037             particles => grid_particles(kp,jp,ip)%particles(1:number_of_particles)
3038             IF ( number_of_particles <= 0 )  CYCLE
3039             WRITE ( 90 )  particles
3040          ENDDO
3041       ENDDO
3042    ENDDO
3043
3044    CLOSE ( 90 )
3045
3046#if defined( __parallel )
3047       CALL MPI_BARRIER( comm2d, ierr )
3048#endif
3049
3050    CALL wrd_write_string( 'iran' ) 
3051    WRITE ( 14 )  iran, iran_part 
3052
3053
3054 END SUBROUTINE lpm_wrd_local
3055
3056
3057!------------------------------------------------------------------------------!
3058! Description:
3059! ------------
3060!> This routine writes the respective restart data for the lpm.
3061!------------------------------------------------------------------------------!
3062 SUBROUTINE lpm_wrd_global
3063 
3064    CALL wrd_write_string( 'curvature_solution_effects' ) 
3065    WRITE ( 14 )  curvature_solution_effects
3066
3067    CALL wrd_write_string( 'interpolation_simple_corrector' )
3068    WRITE ( 14 )  interpolation_simple_corrector
3069
3070    CALL wrd_write_string( 'interpolation_simple_predictor' )
3071    WRITE ( 14 )  interpolation_simple_predictor
3072
3073    CALL wrd_write_string( 'interpolation_trilinear' )
3074    WRITE ( 14 )  interpolation_trilinear
3075
3076 END SUBROUTINE lpm_wrd_global
3077 
3078
3079!------------------------------------------------------------------------------!
3080! Description:
3081! ------------
3082!> This routine writes the respective restart data for the lpm.
3083!------------------------------------------------------------------------------!
3084 SUBROUTINE lpm_rrd_global( found )
3085 
3086    USE control_parameters,                            &
3087        ONLY: length, restart_string
3088
3089    LOGICAL, INTENT(OUT)  ::  found
3090
3091    found = .TRUE.
3092
3093    SELECT CASE ( restart_string(1:length) )
3094
3095       CASE ( 'curvature_solution_effects' )
3096          READ ( 13 )  curvature_solution_effects
3097
3098       CASE ( 'interpolation_simple_corrector' )
3099          READ ( 13 )  interpolation_simple_corrector
3100
3101       CASE ( 'interpolation_simple_predictor' )
3102          READ ( 13 )  interpolation_simple_predictor
3103
3104       CASE ( 'interpolation_trilinear' )
3105          READ ( 13 )  interpolation_trilinear
3106
3107!          CASE ( 'global_paramter' )
3108!             READ ( 13 )  global_parameter
3109!          CASE ( 'global_array' )
3110!             IF ( .NOT. ALLOCATED( global_array ) )  ALLOCATE( global_array(1:10) )
3111!             READ ( 13 )  global_array
3112
3113       CASE DEFAULT
3114
3115          found = .FALSE.
3116
3117    END SELECT
3118   
3119 END SUBROUTINE lpm_rrd_global
3120
3121
3122!------------------------------------------------------------------------------!
3123! Description:
3124! ------------
3125!> This is a submodule of the lagrangian particle model. It contains all
3126!> dynamic processes of the lpm. This includes the advection (resolved and sub-
3127!> grid scale) as well as the boundary conditions of particles. As a next step
3128!> this submodule should be excluded as an own file.
3129!------------------------------------------------------------------------------!
3130 SUBROUTINE lpm_advec (ip,jp,kp)
3131
3132    LOGICAL ::  subbox_at_wall !< flag to see if the current subgridbox is adjacent to a wall
3133
3134    INTEGER(iwp) ::  i                           !< index variable along x
3135    INTEGER(iwp) ::  i_next                      !< index variable along x
3136    INTEGER(iwp) ::  ip                          !< index variable along x
3137    INTEGER(iwp) ::  iteration_steps = 1         !< amount of iterations steps for corrector step
3138    INTEGER(iwp) ::  j                           !< index variable along y
3139    INTEGER(iwp) ::  j_next                      !< index variable along y
3140    INTEGER(iwp) ::  jp                          !< index variable along y
3141    INTEGER(iwp) ::  k                           !< index variable along z
3142    INTEGER(iwp) ::  k_wall                      !< vertical index of topography top
3143    INTEGER(iwp) ::  kp                          !< index variable along z
3144    INTEGER(iwp) ::  k_next                      !< index variable along z
3145    INTEGER(iwp) ::  kw                          !< index variable along z
3146    INTEGER(iwp) ::  kkw                         !< index variable along z
3147    INTEGER(iwp) ::  n                           !< loop variable over all particles in a grid box
3148    INTEGER(iwp) ::  nb                          !< block number particles are sorted in
3149    INTEGER(iwp) ::  particle_end                !< end index for partilce loop
3150    INTEGER(iwp) ::  particle_start              !< start index for particle loop
3151    INTEGER(iwp) ::  surf_start                  !< Index on surface data-type for current grid box
3152    INTEGER(iwp) ::  subbox_end                  !< end index for loop over subboxes in particle advection
3153    INTEGER(iwp) ::  subbox_start                !< start index for loop over subboxes in particle advection
3154    INTEGER(iwp) ::  nn                          !< loop variable over iterations steps
3155
3156    INTEGER(iwp), DIMENSION(0:7) ::  start_index !< start particle index for current block
3157    INTEGER(iwp), DIMENSION(0:7) ::  end_index   !< start particle index for current block
3158
3159    REAL(wp) ::  aa                 !< dummy argument for horizontal particle interpolation
3160    REAL(wp) ::  alpha              !< interpolation facor for x-direction
3161
3162    REAL(wp) ::  bb                 !< dummy argument for horizontal particle interpolation
3163    REAL(wp) ::  beta               !< interpolation facor for y-direction
3164    REAL(wp) ::  cc                 !< dummy argument for horizontal particle interpolation
3165    REAL(wp) ::  d_z_p_z0           !< inverse of interpolation length for logarithmic interpolation
3166    REAL(wp) ::  dd                 !< dummy argument for horizontal particle interpolation
3167    REAL(wp) ::  de_dx_int_l        !< x/y-interpolated TKE gradient (x) at particle position at lower vertical level
3168    REAL(wp) ::  de_dx_int_u        !< x/y-interpolated TKE gradient (x) at particle position at upper vertical level
3169    REAL(wp) ::  de_dy_int_l        !< x/y-interpolated TKE gradient (y) at particle position at lower vertical level
3170    REAL(wp) ::  de_dy_int_u        !< x/y-interpolated TKE gradient (y) at particle position at upper vertical level
3171    REAL(wp) ::  de_dt              !< temporal derivative of TKE experienced by the particle
3172    REAL(wp) ::  de_dt_min          !< lower level for temporal TKE derivative
3173    REAL(wp) ::  de_dz_int_l        !< x/y-interpolated TKE gradient (z) at particle position at lower vertical level
3174    REAL(wp) ::  de_dz_int_u        !< x/y-interpolated TKE gradient (z) at particle position at upper vertical level
3175    REAL(wp) ::  diameter           !< diamter of droplet
3176    REAL(wp) ::  diss_int_l         !< x/y-interpolated dissipation at particle position at lower vertical level
3177    REAL(wp) ::  diss_int_u         !< x/y-interpolated dissipation at particle position at upper vertical level
3178    REAL(wp) ::  dt_particle_m      !< previous particle time step
3179    REAL(wp) ::  dz_temp            !< dummy for the vertical grid spacing
3180    REAL(wp) ::  e_int_l            !< x/y-interpolated TKE at particle position at lower vertical level
3181    REAL(wp) ::  e_int_u            !< x/y-interpolated TKE at particle position at upper vertical level
3182    REAL(wp) ::  e_mean_int         !< horizontal mean TKE at particle height
3183    REAL(wp) ::  exp_arg            !< argument in the exponent - particle radius
3184    REAL(wp) ::  exp_term           !< exponent term
3185    REAL(wp) ::  gamma              !< interpolation facor for z-direction
3186    REAL(wp) ::  gg                 !< dummy argument for horizontal particle interpolation
3187    REAL(wp) ::  height_p           !< dummy argument for logarithmic interpolation
3188    REAL(wp) ::  log_z_z0_int       !< logarithmus used for surface_layer interpolation
3189    REAL(wp) ::  random_gauss       !< Gaussian-distributed random number used for SGS particle advection
3190    REAL(wp) ::  RL                 !< Lagrangian autocorrelation coefficient
3191    REAL(wp) ::  rg1                !< Gaussian distributed random number
3192    REAL(wp) ::  rg2                !< Gaussian distributed random number
3193    REAL(wp) ::  rg3                !< Gaussian distributed random number
3194    REAL(wp) ::  sigma              !< velocity standard deviation
3195    REAL(wp) ::  u_int_l            !< x/y-interpolated u-component at particle position at lower vertical level
3196    REAL(wp) ::  u_int_u            !< x/y-interpolated u-component at particle position at upper vertical level
3197    REAL(wp) ::  unext              !< calculated particle u-velocity of corrector step
3198    REAL(wp) ::  us_int             !< friction velocity at particle grid box
3199    REAL(wp) ::  usws_int           !< surface momentum flux (u component) at particle grid box
3200    REAL(wp) ::  v_int_l            !< x/y-interpolated v-component at particle position at lower vertical level
3201    REAL(wp) ::  v_int_u            !< x/y-interpolated v-component at particle position at upper vertical level
3202    REAL(wp) ::  vsws_int           !< surface momentum flux (u component) at particle grid box
3203    REAL(wp) ::  vnext              !< calculated particle v-velocity of corrector step
3204    REAL(wp) ::  vv_int             !< dummy to compute interpolated mean SGS TKE, used to scale SGS advection
3205    REAL(wp) ::  w_int_l            !< x/y-interpolated w-component at particle position at lower vertical level
3206    REAL(wp) ::  w_int_u            !< x/y-interpolated w-component at particle position at upper vertical level
3207    REAL(wp) ::  wnext              !< calculated particle w-velocity of corrector step
3208    REAL(wp) ::  w_s                !< terminal velocity of droplets
3209    REAL(wp) ::  x                  !< dummy argument for horizontal particle interpolation
3210    REAL(wp) ::  xp                 !< calculated particle position in x of predictor step
3211    REAL(wp) ::  y                  !< dummy argument for horizontal particle interpolation
3212    REAL(wp) ::  yp                 !< calculated particle position in y of predictor step
3213    REAL(wp) ::  z_p                !< surface layer height (0.5 dz)
3214    REAL(wp) ::  zp                 !< calculated particle position in z of predictor step
3215
3216    REAL(wp), PARAMETER ::  a_rog = 9.65_wp      !< parameter for fall velocity
3217    REAL(wp), PARAMETER ::  b_rog = 10.43_wp     !< parameter for fall velocity
3218    REAL(wp), PARAMETER ::  c_rog = 0.6_wp       !< parameter for fall velocity
3219    REAL(wp), PARAMETER ::  k_cap_rog = 4.0_wp   !< parameter for fall velocity
3220    REAL(wp), PARAMETER ::  k_low_rog = 12.0_wp  !< parameter for fall velocity
3221    REAL(wp), PARAMETER ::  d0_rog = 0.745_wp    !< separation diameter
3222
3223    REAL(wp), DIMENSION(number_of_particles) ::  term_1_2       !< flag to communicate whether a particle is near topography or not
3224    REAL(wp), DIMENSION(number_of_particles) ::  dens_ratio     !< ratio between the density of the fluid and the density of the particles
3225    REAL(wp), DIMENSION(number_of_particles) ::  de_dx_int      !< horizontal TKE gradient along x at particle position
3226    REAL(wp), DIMENSION(number_of_particles) ::  de_dy_int      !< horizontal TKE gradient along y at particle position
3227    REAL(wp), DIMENSION(number_of_particles) ::  de_dz_int      !< horizontal TKE gradient along z at particle position
3228    REAL(wp), DIMENSION(number_of_particles) ::  diss_int       !< dissipation at particle position
3229    REAL(wp), DIMENSION(number_of_particles) ::  dt_gap         !< remaining time until particle time integration reaches LES time
3230    REAL(wp), DIMENSION(number_of_particles) ::  dt_particle    !< particle time step
3231    REAL(wp), DIMENSION(number_of_particles) ::  e_int          !< TKE at particle position
3232    REAL(wp), DIMENSION(number_of_particles) ::  fs_int         !< weighting factor for subgrid-scale particle speed
3233    REAL(wp), DIMENSION(number_of_particles) ::  lagr_timescale !< Lagrangian timescale
3234    REAL(wp), DIMENSION(number_of_particles) ::  rvar1_temp     !< SGS particle velocity - u-component
3235    REAL(wp), DIMENSION(number_of_particles) ::  rvar2_temp     !< SGS particle velocity - v-component
3236    REAL(wp), DIMENSION(number_of_particles) ::  rvar3_temp     !< SGS particle velocity - w-component
3237    REAL(wp), DIMENSION(number_of_particles) ::  u_int          !< u-component of particle speed
3238    REAL(wp), DIMENSION(number_of_particles) ::  v_int          !< v-component of particle speed
3239    REAL(wp), DIMENSION(number_of_particles) ::  w_int          !< w-component of particle speed
3240    REAL(wp), DIMENSION(number_of_particles) ::  xv             !< x-position
3241    REAL(wp), DIMENSION(number_of_particles) ::  yv             !< y-position
3242    REAL(wp), DIMENSION(number_of_particles) ::  zv             !< z-position
3243
3244    REAL(wp), DIMENSION(number_of_particles, 3) ::  rg !< vector of Gaussian distributed random numbers
3245
3246    CALL cpu_log( log_point_s(44), 'lpm_advec', 'continue' )
3247!
3248!-- Determine height of Prandtl layer and distance between Prandtl-layer
3249!-- height and horizontal mean roughness height, which are required for
3250!-- vertical logarithmic interpolation of horizontal particle speeds
3251!-- (for particles below first vertical grid level).
3252    z_p      = zu(nzb+1) - zw(nzb)
3253    d_z_p_z0 = 1.0_wp / ( z_p - z0_av_global )
3254
3255    xv = particles(1:number_of_particles)%x
3256    yv = particles(1:number_of_particles)%y
3257    zv = particles(1:number_of_particles)%z
3258    dt_particle = dt_3d
3259
3260!
3261!-- This case uses a simple interpolation method for the particle velocites,
3262!-- and applying a predictor-corrector method. @attention: for the corrector
3263!-- step the velocities of t(n+1) are required. However, at this moment of
3264!-- the time integration they are not free of divergence. This interpolation
3265!-- method is described in more detail in Grabowski et al., 2018 (GMD).
3266    IF ( interpolation_simple_corrector )  THEN
3267!
3268!--    Predictor step
3269       kkw = kp - 1
3270       DO  n = 1, number_of_particles
3271
3272          alpha = MAX( MIN( ( particles(n)%x - ip * dx ) * ddx, 1.0_wp ), 0.0_wp )
3273          u_int(n) = u(kp,jp,ip) * ( 1.0_wp - alpha ) + u(kp,jp,ip+1) * alpha
3274
3275          beta  = MAX( MIN( ( particles(n)%y - jp * dy ) * ddy, 1.0_wp ), 0.0_wp )
3276          v_int(n) = v(kp,jp,ip) * ( 1.0_wp - beta ) + v(kp,jp+1,ip) * beta
3277
3278          gamma = MAX( MIN( ( particles(n)%z - zw(kkw) ) /                   &
3279                            ( zw(kkw+1) - zw(kkw) ), 1.0_wp ), 0.0_wp )
3280          w_int(n) = w(kkw,jp,ip) * ( 1.0_wp - gamma ) + w(kkw+1,jp,ip) * gamma
3281
3282       ENDDO
3283!
3284!--    Corrector step
3285       DO  n = 1, number_of_particles
3286
3287          IF ( .NOT. particles(n)%particle_mask )  CYCLE
3288
3289          DO  nn = 1, iteration_steps
3290
3291!
3292!--          Guess new position
3293             xp = particles(n)%x + u_int(n) * dt_particle(n)
3294             yp = particles(n)%y + v_int(n) * dt_particle(n)
3295             zp = particles(n)%z + w_int(n) * dt_particle(n)
3296!
3297!--          x direction
3298             i_next = FLOOR( xp * ddx , KIND=iwp)
3299             alpha  = MAX( MIN( ( xp - i_next * dx ) * ddx, 1.0_wp ), 0.0_wp )
3300!
3301!--          y direction
3302             j_next = FLOOR( yp * ddy )
3303             beta   = MAX( MIN( ( yp - j_next * dy ) * ddy, 1.0_wp ), 0.0_wp )
3304!
3305!--          z_direction
3306             k_next = MAX( MIN( FLOOR( zp / (zw(kkw+1)-zw(kkw)) + offset_ocean_nzt ), nzt ), 0)
3307             gamma = MAX( MIN( ( zp - zw(k_next) ) /                      &
3308                               ( zw(k_next+1) - zw(k_next) ), 1.0_wp ), 0.0_wp )
3309!
3310!--          Calculate part of the corrector step
3311             unext = u_p(k_next+1, j_next, i_next) * ( 1.0_wp - alpha ) +    &
3312                     u_p(k_next+1, j_next,   i_next+1) * alpha
3313
3314             vnext = v_p(k_next+1, j_next, i_next) * ( 1.0_wp - beta  ) +    &
3315                     v_p(k_next+1, j_next+1, i_next  ) * beta
3316
3317             wnext = w_p(k_next,   j_next, i_next) * ( 1.0_wp - gamma ) +    &
3318                     w_p(k_next+1, j_next, i_next  ) * gamma
3319
3320!
3321!--          Calculate interpolated particle velocity with predictor
3322!--          corrector step. u_int, v_int and w_int describes the part of
3323!--          the predictor step. unext, vnext and wnext is the part of the
3324!--          corrector step. The resulting new position is set below. The
3325!--          implementation is based on Grabowski et al., 2018 (GMD).
3326             u_int(n) = 0.5_wp * ( u_int(n) + unext )
3327             v_int(n) = 0.5_wp * ( v_int(n) + vnext )
3328             w_int(n) = 0.5_wp * ( w_int(n) + wnext )
3329
3330          ENDDO
3331       ENDDO
3332!
3333!-- This case uses a simple interpolation method for the particle velocites,
3334!-- and applying a predictor.
3335    ELSEIF ( interpolation_simple_predictor )  THEN
3336!
3337!--    The particle position for the w velociy is based on the value of kp and kp-1
3338       kkw = kp - 1
3339       DO  n = 1, number_of_particles
3340          IF ( .NOT. particles(n)%particle_mask )  CYCLE
3341
3342          alpha    = MAX( MIN( ( particles(n)%x - ip * dx ) * ddx, 1.0_wp ), 0.0_wp )
3343          u_int(n) = u(kp,jp,ip) * ( 1.0_wp - alpha ) + u(kp,jp,ip+1) * alpha
3344
3345          beta     = MAX( MIN( ( particles(n)%y - jp * dy ) * ddy, 1.0_wp ), 0.0_wp )
3346          v_int(n) = v(kp,jp,ip) * ( 1.0_wp - beta ) + v(kp,jp+1,ip) * beta
3347
3348          gamma    = MAX( MIN( ( particles(n)%z - zw(kkw) ) /                   &
3349                               ( zw(kkw+1) - zw(kkw) ), 1.0_wp ), 0.0_wp )
3350          w_int(n) = w(kkw,jp,ip) * ( 1.0_wp - gamma ) + w(kkw+1,jp,ip) * gamma
3351       ENDDO
3352!
3353!-- The trilinear interpolation.
3354    ELSEIF ( interpolation_trilinear )  THEN
3355
3356       start_index = grid_particles(kp,jp,ip)%start_index
3357       end_index   = grid_particles(kp,jp,ip)%end_index
3358
3359       DO  nb = 0, 7
3360!
3361!--       Interpolate u velocity-component
3362          i = ip
3363          j = jp + block_offset(nb)%j_off
3364          k = kp + block_offset(nb)%k_off
3365
3366          DO  n = start_index(nb), end_index(nb)
3367!
3368!--          Interpolation of the u velocity component onto particle position.
3369!--          Particles are interpolation bi-linearly in the horizontal and a
3370!--          linearly in the vertical. An exception is made for particles below
3371!--          the first vertical grid level in case of a prandtl layer. In this
3372!--          case the horizontal particle velocity components are determined using
3373!--          Monin-Obukhov relations (if branch).
3374!--          First, check if particle is located below first vertical grid level
3375!--          above topography (Prandtl-layer height)
3376!--          Determine vertical index of topography top
3377             k_wall = topo_top_ind(jp,ip,0)
3378
3379             IF ( constant_flux_layer  .AND.  zv(n) - zw(k_wall) < z_p )  THEN
3380!
3381!--             Resolved-scale horizontal particle velocity is zero below z0.
3382                IF ( zv(n) - zw(k_wall) < z0_av_global )  THEN
3383                   u_int(n) = 0.0_wp
3384                ELSE
3385!
3386!--                Determine the sublayer. Further used as index.
3387                   height_p = ( zv(n) - zw(k_wall) - z0_av_global ) &
3388                                        * REAL( number_of_sublayers, KIND=wp )    &
3389                                        * d_z_p_z0
3390!
3391!--                Calculate LOG(z/z0) for exact particle height. Therefore,
3392!--                interpolate linearly between precalculated logarithm.
3393                   log_z_z0_int = log_z_z0(INT(height_p))                         &
3394                                    + ( height_p - INT(height_p) )                &
3395                                    * ( log_z_z0(INT(height_p)+1)                 &
3396                                         - log_z_z0(INT(height_p))                &
3397                                      )
3398!
3399!--                Get friction velocity and momentum flux from new surface data
3400!--                types.
3401                   IF ( surf_def_h(0)%start_index(jp,ip) <=                   &
3402                        surf_def_h(0)%end_index(jp,ip) )  THEN
3403                      surf_start = surf_def_h(0)%start_index(jp,ip)
3404!--                   Limit friction velocity. In narrow canyons or holes the
3405!--                   friction velocity can become very small, resulting in a too
3406!--                   large particle speed.
3407                      us_int    = MAX( surf_def_h(0)%us(surf_start), 0.01_wp )
3408                      usws_int  = surf_def_h(0)%usws(surf_start)
3409                   ELSEIF ( surf_lsm_h%start_index(jp,ip) <=                  &
3410                            surf_lsm_h%end_index(jp,ip) )  THEN
3411                      surf_start = surf_lsm_h%start_index(jp,ip)
3412                      us_int    = MAX( surf_lsm_h%us(surf_start), 0.01_wp )
3413                      usws_int  = surf_lsm_h%usws(surf_start)
3414                   ELSEIF ( surf_usm_h%start_index(jp,ip) <=                  &
3415                            surf_usm_h%end_index(jp,ip) )  THEN
3416                      surf_start = surf_usm_h%start_index(jp,ip)
3417                      us_int    = MAX( surf_usm_h%us(surf_start), 0.01_wp )
3418                      usws_int  = surf_usm_h%usws(surf_start)
3419                   ENDIF
3420!
3421!--                Neutral solution is applied for all situations, e.g. also for
3422!--                unstable and stable situations. Even though this is not exact
3423!--                this saves a lot of CPU time since several calls of intrinsic
3424!--                FORTRAN procedures (LOG, ATAN) are avoided, This is justified
3425!--                as sensitivity studies revealed no significant effect of
3426!--                using the neutral solution also for un/stable situations.
3427                   u_int(n) = -usws_int / ( us_int * kappa + 1E-10_wp )           &
3428                               * log_z_z0_int - u_gtrans
3429                ENDIF
3430!
3431!--          Particle above the first grid level. Bi-linear interpolation in the
3432!--          horizontal and linear interpolation in the vertical direction.
3433             ELSE
3434                x  = xv(n) - i * dx
3435                y  = yv(n) + ( 0.5_wp - j ) * dy
3436                aa = x**2          + y**2
3437                bb = ( dx - x )**2 + y**2
3438                cc = x**2          + ( dy - y )**2
3439                dd = ( dx - x )**2 + ( dy - y )**2
3440                gg = aa + bb + cc + dd
3441
3442                u_int_l = ( ( gg - aa ) * u(k,j,i)   + ( gg - bb ) * u(k,j,i+1)   &
3443                            + ( gg - cc ) * u(k,j+1,i) + ( gg - dd ) *            &
3444                            u(k,j+1,i+1) ) / ( 3.0_wp * gg ) - u_gtrans
3445
3446                IF ( k == nzt )  THEN
3447                   u_int(n) = u_int_l
3448                ELSE
3449                   u_int_u = ( ( gg-aa ) * u(k+1,j,i) + ( gg-bb ) * u(k+1,j,i+1)  &
3450                               + ( gg-cc ) * u(k+1,j+1,i) + ( gg-dd ) *           &
3451                               u(k+1,j+1,i+1) ) / ( 3.0_wp * gg ) - u_gtrans
3452                   u_int(n) = u_int_l + ( zv(n) - zu(k) ) / dzw(k+1) *            &
3453                              ( u_int_u - u_int_l )
3454                ENDIF
3455             ENDIF
3456          ENDDO
3457!
3458!--       Same procedure for interpolation of the v velocity-component
3459          i = ip + block_offset(nb)%i_off
3460          j = jp
3461          k = kp + block_offset(nb)%k_off
3462
3463          DO  n = start_index(nb), end_index(nb)
3464!
3465!--          Determine vertical index of topography top
3466             k_wall = topo_top_ind(jp,ip,0)
3467
3468             IF ( constant_flux_layer  .AND.  zv(n) - zw(k_wall) < z_p )  THEN
3469                IF ( zv(n) - zw(k_wall) < z0_av_global )  THEN
3470!
3471!--                Resolved-scale horizontal particle velocity is zero below z0.
3472                   v_int(n) = 0.0_wp
3473                ELSE
3474!
3475!--                Determine the sublayer. Further used as index. Please note,
3476!--                logarithmus can not be reused from above, as in in case of
3477!--                topography particle on u-grid can be above surface-layer height,
3478!--                whereas it can be below on v-grid.
3479                   height_p = ( zv(n) - zw(k_wall) - z0_av_global ) &
3480                                     * REAL( number_of_sublayers, KIND=wp )       &
3481                                     * d_z_p_z0
3482!
3483!--                Calculate LOG(z/z0) for exact particle height. Therefore,
3484!--                interpolate linearly between precalculated logarithm.
3485                   log_z_z0_int = log_z_z0(INT(height_p))                         &
3486                                    + ( height_p - INT(height_p) )                &
3487                                    * ( log_z_z0(INT(height_p)+1)                 &
3488                                         - log_z_z0(INT(height_p))                &
3489                                      )
3490!
3491!--                Get friction velocity and momentum flux from new surface data
3492!--                types.
3493                   IF ( surf_def_h(0)%start_index(jp,ip) <=                   &
3494                        surf_def_h(0)%end_index(jp,ip) )  THEN
3495                      surf_start = surf_def_h(0)%start_index(jp,ip)
3496!--                   Limit friction velocity. In narrow canyons or holes the
3497!--                   friction velocity can become very small, resulting in a too
3498!--                   large particle speed.
3499                      us_int    = MAX( surf_def_h(0)%us(surf_start), 0.01_wp )
3500                      vsws_int  = surf_def_h(0)%vsws(surf_start)
3501                   ELSEIF ( surf_lsm_h%start_index(jp,ip) <=                  &
3502                            surf_lsm_h%end_index(jp,ip) )  THEN
3503                      surf_start = surf_lsm_h%start_index(jp,ip)
3504                      us_int    = MAX( surf_lsm_h%us(surf_start), 0.01_wp )
3505                      vsws_int  = surf_lsm_h%vsws(surf_start)
3506                   ELSEIF ( surf_usm_h%start_index(jp,ip) <=                  &
3507                            surf_usm_h%end_index(jp,ip) )  THEN
3508                      surf_start = surf_usm_h%start_index(jp,ip)
3509                      us_int    = MAX( surf_usm_h%us(surf_start), 0.01_wp )
3510                      vsws_int  = surf_usm_h%vsws(surf_start)
3511                   ENDIF
3512!
3513!--                Neutral solution is applied for all situations, e.g. also for
3514!--                unstable and stable situations. Even though this is not exact
3515!--                this saves a lot of CPU time since several calls of intrinsic
3516!--                FORTRAN procedures (LOG, ATAN) are avoided, This is justified
3517!--                as sensitivity studies revealed no significant effect of
3518!--                using the neutral solution also for un/stable situations.
3519                   v_int(n) = -vsws_int / ( us_int * kappa + 1E-10_wp )           &
3520                            * log_z_z0_int - v_gtrans
3521
3522                ENDIF
3523             ELSE
3524                x  = xv(n) + ( 0.5_wp - i ) * dx
3525                y  = yv(n) - j * dy
3526                aa = x**2          + y**2
3527                bb = ( dx - x )**2 + y**2
3528                cc = x**2          + ( dy - y )**2
3529                dd = ( dx - x )**2 + ( dy - y )**2
3530                gg = aa + bb + cc + dd
3531
3532                v_int_l = ( ( gg - aa ) * v(k,j,i)   + ( gg - bb ) * v(k,j,i+1)   &
3533                          + ( gg - cc ) * v(k,j+1,i) + ( gg - dd ) * v(k,j+1,i+1) &
3534                          ) / ( 3.0_wp * gg ) - v_gtrans
3535
3536                IF ( k == nzt )  THEN
3537                   v_int(n) = v_int_l
3538                ELSE
3539                   v_int_u = ( ( gg-aa ) * v(k+1,j,i)   + ( gg-bb ) * v(k+1,j,i+1)   &
3540                             + ( gg-cc ) * v(k+1,j+1,i) + ( gg-dd ) * v(k+1,j+1,i+1) &
3541                             ) / ( 3.0_wp * gg ) - v_gtrans
3542                   v_int(n) = v_int_l + ( zv(n) - zu(k) ) / dzw(k+1) *               &
3543                                     ( v_int_u - v_int_l )
3544                ENDIF
3545             ENDIF
3546          ENDDO
3547!
3548!--       Same procedure for interpolation of the w velocity-component
3549          i = ip + block_offset(nb)%i_off
3550          j = jp + block_offset(nb)%j_off
3551          k = kp - 1
3552
3553          DO  n = start_index(nb), end_index(nb)
3554             IF ( vertical_particle_advection(particles(n)%group) )  THEN
3555                x  = xv(n) + ( 0.5_wp - i ) * dx
3556                y  = yv(n) + ( 0.5_wp - j ) * dy
3557                aa = x**2          + y**2
3558                bb = ( dx - x )**2 + y**2
3559                cc = x**2          + ( dy - y )**2
3560                dd = ( dx - x )**2 + ( dy - y )**2
3561                gg = aa + bb + cc + dd
3562
3563                w_int_l = ( ( gg - aa ) * w(k,j,i)   + ( gg - bb ) * w(k,j,i+1)   &
3564                          + ( gg - cc ) * w(k,j+1,i) + ( gg - dd ) * w(k,j+1,i+1) &
3565                          ) / ( 3.0_wp * gg )
3566
3567                IF ( k == nzt )  THEN
3568                   w_int(n) = w_int_l
3569                ELSE
3570                   w_int_u = ( ( gg-aa ) * w(k+1,j,i)   + &
3571                               ( gg-bb ) * w(k+1,j,i+1) + &
3572                               ( gg-cc ) * w(k+1,j+1,i) + &
3573                               ( gg-dd ) * w(k+1,j+1,i+1) &
3574                             ) / ( 3.0_wp * gg )
3575                   w_int(n) = w_int_l + ( zv(n) - zw(k) ) / dzw(k+1) *               &
3576                              ( w_int_u - w_int_l )
3577                ENDIF
3578             ELSE
3579                w_int(n) = 0.0_wp
3580             ENDIF
3581          ENDDO
3582       ENDDO
3583    ENDIF
3584
3585!-- Interpolate and calculate quantities needed for calculating the SGS
3586!-- velocities
3587    IF ( use_sgs_for_particles  .AND.  .NOT. cloud_droplets )  THEN
3588
3589       DO  nb = 0,7
3590
3591          subbox_at_wall = .FALSE.
3592!
3593!--       In case of topography check if subbox is adjacent to a wall
3594          IF ( .NOT. topography == 'flat' )  THEN
3595             i = ip + MERGE( -1_iwp , 1_iwp, BTEST( nb, 2 ) )
3596             j = jp + MERGE( -1_iwp , 1_iwp, BTEST( nb, 1 ) )
3597             k = kp + MERGE( -1_iwp , 1_iwp, BTEST( nb, 0 ) )
3598             IF ( .NOT. BTEST(wall_flags_0(k,  jp, ip), 0) .OR.                &
3599                  .NOT. BTEST(wall_flags_0(kp, j,  ip), 0) .OR.                &
3600                  .NOT. BTEST(wall_flags_0(kp, jp, i ), 0) )                   &
3601             THEN
3602                subbox_at_wall = .TRUE.
3603             ENDIF
3604          ENDIF
3605          IF ( subbox_at_wall )  THEN
3606             e_int(start_index(nb):end_index(nb))     = e(kp,jp,ip) 
3607             diss_int(start_index(nb):end_index(nb))  = diss(kp,jp,ip)
3608             de_dx_int(start_index(nb):end_index(nb)) = de_dx(kp,jp,ip)
3609             de_dy_int(start_index(nb):end_index(nb)) = de_dy(kp,jp,ip)
3610             de_dz_int(start_index(nb):end_index(nb)) = de_dz(kp,jp,ip)
3611!
3612!--          Set flag for stochastic equation.
3613             term_1_2(start_index(nb):end_index(nb)) = 0.0_wp
3614          ELSE
3615             i = ip + block_offset(nb)%i_off
3616             j = jp + block_offset(nb)%j_off
3617             k = kp + block_offset(nb)%k_off
3618
3619             DO  n = start_index(nb), end_index(nb)
3620!
3621!--             Interpolate TKE
3622                x  = xv(n) + ( 0.5_wp - i ) * dx
3623                y  = yv(n) + ( 0.5_wp - j ) * dy
3624                aa = x**2          + y**2
3625                bb = ( dx - x )**2 + y**2
3626                cc = x**2          + ( dy - y )**2
3627                dd = ( dx - x )**2 + ( dy - y )**2
3628                gg = aa + bb + cc + dd
3629
3630                e_int_l = ( ( gg-aa ) * e(k,j,i)   + ( gg-bb ) * e(k,j,i+1)   &
3631                          + ( gg-cc ) * e(k,j+1,i) + ( gg-dd ) * e(k,j+1,i+1) &
3632                          ) / ( 3.0_wp * gg )
3633
3634                IF ( k+1 == nzt+1 )  THEN
3635                   e_int(n) = e_int_l
3636                ELSE
3637                   e_int_u = ( ( gg - aa ) * e(k+1,j,i)   + &
3638                               ( gg - bb ) * e(k+1,j,i+1) + &
3639                               ( gg - cc ) * e(k+1,j+1,i) + &
3640                               ( gg - dd ) * e(k+1,j+1,i+1) &
3641                            ) / ( 3.0_wp * gg )
3642                   e_int(n) = e_int_l + ( zv(n) - zu(k) ) / dzw(k+1) *            &
3643                                     ( e_int_u - e_int_l )
3644                ENDIF
3645!
3646!--             Needed to avoid NaN particle velocities (this might not be
3647!--             required any more)
3648                IF ( e_int(n) <= 0.0_wp )  THEN
3649                   e_int(n) = 1.0E-20_wp
3650                ENDIF
3651!
3652!--             Interpolate the TKE gradient along x (adopt incides i,j,k and
3653!--             all position variables from above (TKE))
3654                de_dx_int_l = ( ( gg - aa ) * de_dx(k,j,i)   + &
3655                                ( gg - bb ) * de_dx(k,j,i+1) + &
3656                                ( gg - cc ) * de_dx(k,j+1,i) + &
3657                                ( gg - dd ) * de_dx(k,j+1,i+1) &
3658                               ) / ( 3.0_wp * gg )
3659
3660                IF ( ( k+1 == nzt+1 )  .OR.  ( k == nzb ) )  THEN
3661                   de_dx_int(n) = de_dx_int_l
3662                ELSE
3663                   de_dx_int_u = ( ( gg - aa ) * de_dx(k+1,j,i)   + &
3664                                   ( gg - bb ) * de_dx(k+1,j,i+1) + &
3665                                   ( gg - cc ) * de_dx(k+1,j+1,i) + &
3666                                   ( gg - dd ) * de_dx(k+1,j+1,i+1) &
3667                                  ) / ( 3.0_wp * gg )
3668                   de_dx_int(n) = de_dx_int_l + ( zv(n) - zu(k) ) / dzw(k+1) *    &
3669                                              ( de_dx_int_u - de_dx_int_l )
3670                ENDIF
3671!
3672!--             Interpolate the TKE gradient along y
3673                de_dy_int_l = ( ( gg - aa ) * de_dy(k,j,i)   + &
3674                                ( gg - bb ) * de_dy(k,j,i+1) + &
3675                                ( gg - cc ) * de_dy(k,j+1,i) + &
3676                                ( gg - dd ) * de_dy(k,j+1,i+1) &
3677                               ) / ( 3.0_wp * gg )
3678                IF ( ( k+1 == nzt+1 )  .OR.  ( k == nzb ) )  THEN
3679                   de_dy_int(n) = de_dy_int_l
3680                ELSE
3681                   de_dy_int_u = ( ( gg - aa ) * de_dy(k+1,j,i)   + &
3682                                   ( gg - bb ) * de_dy(k+1,j,i+1) + &
3683                                   ( gg - cc ) * de_dy(k+1,j+1,i) + &
3684                                   ( gg - dd ) * de_dy(k+1,j+1,i+1) &
3685                                  ) / ( 3.0_wp * gg )
3686                      de_dy_int(n) = de_dy_int_l + ( zv(n) - zu(k) ) / dzw(k+1) * &
3687                                                 ( de_dy_int_u - de_dy_int_l )
3688                ENDIF
3689
3690!
3691!--             Interpolate the TKE gradient along z
3692                IF ( zv(n) < 0.5_wp * dz(1) )  THEN
3693                   de_dz_int(n) = 0.0_wp
3694                ELSE
3695                   de_dz_int_l = ( ( gg - aa ) * de_dz(k,j,i)   + &
3696                                   ( gg - bb ) * de_dz(k,j,i+1) + &
3697                                   ( gg - cc ) * de_dz(k,j+1,i) + &
3698                                   ( gg - dd ) * de_dz(k,j+1,i+1) &
3699                                  ) / ( 3.0_wp * gg )
3700
3701                   IF ( ( k+1 == nzt+1 )  .OR.  ( k == nzb ) )  THEN
3702                      de_dz_int(n) = de_dz_int_l
3703                   ELSE
3704                      de_dz_int_u = ( ( gg - aa ) * de_dz(k+1,j,i)   + &
3705                                      ( gg - bb ) * de_dz(k+1,j,i+1) + &
3706                                      ( gg - cc ) * de_dz(k+1,j+1,i) + &
3707                                      ( gg - dd ) * de_dz(k+1,j+1,i+1) &
3708                                     ) / ( 3.0_wp * gg )
3709                      de_dz_int(n) = de_dz_int_l + ( zv(n) - zu(k) ) / dzw(k+1) * &
3710                                                 ( de_dz_int_u - de_dz_int_l )
3711                   ENDIF
3712                ENDIF
3713
3714!
3715!--             Interpolate the dissipation of TKE
3716                diss_int_l = ( ( gg - aa ) * diss(k,j,i)   + &
3717                               ( gg - bb ) * diss(k,j,i+1) + &
3718                               ( gg - cc ) * diss(k,j+1,i) + &
3719                               ( gg - dd ) * diss(k,j+1,i+1) &
3720                               ) / ( 3.0_wp * gg )
3721
3722                IF ( k == nzt )  THEN
3723                   diss_int(n) = diss_int_l
3724                ELSE
3725                   diss_int_u = ( ( gg - aa ) * diss(k+1,j,i)   + &
3726                                  ( gg - bb ) * diss(k+1,j,i+1) + &
3727                                  ( gg - cc ) * diss(k+1,j+1,i) + &
3728                                  ( gg - dd ) * diss(k+1,j+1,i+1) &
3729                                 ) / ( 3.0_wp * gg )
3730                   diss_int(n) = diss_int_l + ( zv(n) - zu(k) ) / dzw(k+1) *      &
3731                                            ( diss_int_u - diss_int_l )
3732                ENDIF
3733
3734!
3735!--             Set flag for stochastic equation.
3736                term_1_2(n) = 1.0_wp
3737             ENDDO
3738          ENDIF
3739       ENDDO
3740
3741       DO  nb = 0,7
3742          i = ip + block_offset(nb)%i_off
3743          j = jp + block_offset(nb)%j_off
3744          k = kp + block_offset(nb)%k_off
3745
3746          DO  n = start_index(nb), end_index(nb)
3747!
3748!--          Vertical interpolation of the horizontally averaged SGS TKE and
3749!--          resolved-scale velocity variances and use the interpolated values
3750!--          to calculate the coefficient fs, which is a measure of the ratio
3751!--          of the subgrid-scale turbulent kinetic energy to the total amount
3752!--          of turbulent kinetic energy.
3753             IF ( k == 0 )  THEN
3754                e_mean_int = hom(0,1,8,0)
3755             ELSE
3756                e_mean_int = hom(k,1,8,0) +                                    &
3757                                           ( hom(k+1,1,8,0) - hom(k,1,8,0) ) / &
3758                                           ( zu(k+1) - zu(k) ) *               &
3759                                           ( zv(n) - zu(k) )
3760             ENDIF
3761
3762             kw = kp - 1
3763
3764             IF ( k == 0 )  THEN
3765                aa  = hom(k+1,1,30,0)  * ( zv(n) / &
3766                                         ( 0.5_wp * ( zu(k+1) - zu(k) ) ) )
3767                bb  = hom(k+1,1,31,0)  * ( zv(n) / &
3768                                         ( 0.5_wp * ( zu(k+1) - zu(k) ) ) )
3769                cc  = hom(kw+1,1,32,0) * ( zv(n) / &
3770                                         ( 1.0_wp * ( zw(kw+1) - zw(kw) ) ) )
3771             ELSE
3772                aa  = hom(k,1,30,0) + ( hom(k+1,1,30,0) - hom(k,1,30,0) ) *    &
3773                           ( ( zv(n) - zu(k) ) / ( zu(k+1) - zu(k) ) )
3774                bb  = hom(k,1,31,0) + ( hom(k+1,1,31,0) - hom(k,1,31,0) ) *    &
3775                           ( ( zv(n) - zu(k) ) / ( zu(k+1) - zu(k) ) )
3776                cc  = hom(kw,1,32,0) + ( hom(kw+1,1,32,0)-hom(kw,1,32,0) ) *   &
3777                           ( ( zv(n) - zw(kw) ) / ( zw(kw+1)-zw(kw) ) )
3778             ENDIF
3779
3780             vv_int = ( 1.0_wp / 3.0_wp ) * ( aa + bb + cc )
3781!
3782!--          Needed to avoid NaN particle velocities. The value of 1.0 is just
3783!--          an educated guess for the given case.
3784             IF ( vv_int + ( 2.0_wp / 3.0_wp ) * e_mean_int == 0.0_wp )  THEN
3785                fs_int(n) = 1.0_wp
3786             ELSE
3787                fs_int(n) = ( 2.0_wp / 3.0_wp ) * e_mean_int /                 &
3788                            ( vv_int + ( 2.0_wp / 3.0_wp ) * e_mean_int )
3789             ENDIF
3790
3791          ENDDO
3792       ENDDO
3793
3794       DO  nb = 0, 7
3795          DO  n = start_index(nb), end_index(nb)
3796             rg(n,1) = random_gauss( iran_part, 5.0_wp )
3797             rg(n,2) = random_gauss( iran_part, 5.0_wp )
3798             rg(n,3) = random_gauss( iran_part, 5.0_wp )
3799          ENDDO
3800       ENDDO
3801
3802       DO  nb = 0, 7
3803          DO  n = start_index(nb), end_index(nb)
3804
3805!
3806!--          Calculate the Lagrangian timescale according to Weil et al. (2004).
3807             lagr_timescale(n) = ( 4.0_wp * e_int(n) + 1E-20_wp ) / &
3808                              ( 3.0_wp * fs_int(n) * c_0 * diss_int(n) + 1E-20_wp )
3809
3810!
3811!--          Calculate the next particle timestep. dt_gap is the time needed to
3812!--          complete the current LES timestep.
3813             dt_gap(n) = dt_3d - particles(n)%dt_sum
3814             dt_particle(n) = MIN( dt_3d, 0.025_wp * lagr_timescale(n), dt_gap(n) )
3815             particles(n)%aux1 = lagr_timescale(n)
3816             particles(n)%aux2 = dt_gap(n)
3817!
3818!--          The particle timestep should not be too small in order to prevent
3819!--          the number of particle timesteps of getting too large
3820             IF ( dt_particle(n) < dt_min_part )  THEN
3821                IF ( dt_min_part < dt_gap(n) )  THEN
3822                   dt_particle(n) = dt_min_part
3823                ELSE
3824                   dt_particle(n) = dt_gap(n)
3825                ENDIF
3826             ENDIF
3827             rvar1_temp(n) = particles(n)%rvar1
3828             rvar2_temp(n) = particles(n)%rvar2
3829             rvar3_temp(n) = particles(n)%rvar3
3830!
3831!--          Calculate the SGS velocity components
3832             IF ( particles(n)%age == 0.0_wp )  THEN
3833!
3834!--             For new particles the SGS components are derived from the SGS
3835!--             TKE. Limit the Gaussian random number to the interval
3836!--             [-5.0*sigma, 5.0*sigma] in order to prevent the SGS velocities
3837!--             from becoming unrealistically large.
3838                rvar1_temp(n) = SQRT( 2.0_wp * sgs_wf_part * e_int(n)          &
3839                                          + 1E-20_wp ) * ( rg(n,1) - 1.0_wp )
3840                rvar2_temp(n) = SQRT( 2.0_wp * sgs_wf_part * e_int(n)          &
3841                                          + 1E-20_wp ) * ( rg(n,2) - 1.0_wp )
3842                rvar3_temp(n) = SQRT( 2.0_wp * sgs_wf_part * e_int(n)          &
3843                                          + 1E-20_wp ) * ( rg(n,3) - 1.0_wp )
3844
3845             ELSE
3846!
3847!--             Restriction of the size of the new timestep: compared to the
3848!--             previous timestep the increase must not exceed 200%. First,
3849!--             check if age > age_m, in order to prevent that particles get zero
3850!--             timestep.
3851                dt_particle_m = MERGE( dt_particle(n),                         &
3852                                       particles(n)%age - particles(n)%age_m,  &
3853                                       particles(n)%age - particles(n)%age_m < &
3854                                       1E-8_wp )
3855                IF ( dt_particle(n) > 2.0_wp * dt_particle_m )  THEN
3856                   dt_particle(n) = 2.0_wp * dt_particle_m
3857                ENDIF
3858
3859!--             For old particles the SGS components are correlated with the
3860!--             values from the previous timestep. Random numbers have also to
3861!--             be limited (see above).
3862!--             As negative values for the subgrid TKE are not allowed, the
3863!--             change of the subgrid TKE with time cannot be smaller than
3864!--             -e_int(n)/dt_particle. This value is used as a lower boundary
3865!--             value for the change of TKE
3866                de_dt_min = - e_int(n) / dt_particle(n)
3867
3868                de_dt = ( e_int(n) - particles(n)%e_m ) / dt_particle_m
3869
3870                IF ( de_dt < de_dt_min )  THEN
3871                   de_dt = de_dt_min
3872                ENDIF
3873
3874                CALL weil_stochastic_eq( rvar1_temp(n), fs_int(n), e_int(n),    &
3875                                        de_dx_int(n), de_dt, diss_int(n),       &
3876                                        dt_particle(n), rg(n,1), term_1_2(n) )
3877
3878                CALL weil_stochastic_eq( rvar2_temp(n), fs_int(n), e_int(n),    &
3879                                        de_dy_int(n), de_dt, diss_int(n),       &
3880                                        dt_particle(n), rg(n,2), term_1_2(n) )
3881
3882                CALL weil_stochastic_eq( rvar3_temp(n), fs_int(n), e_int(n),    &
3883                                        de_dz_int(n), de_dt, diss_int(n),       &
3884                                        dt_particle(n), rg(n,3), term_1_2(n) )
3885
3886             ENDIF
3887
3888          ENDDO
3889       ENDDO
3890!
3891!--    Check if the added SGS velocities result in a violation of the CFL-
3892!--    criterion. If yes choose a smaller timestep based on the new velocities
3893!--    and calculate SGS velocities again
3894       dz_temp = zw(kp)-zw(kp-1)
3895
3896       DO  nb = 0, 7
3897          DO  n = start_index(nb), end_index(nb)
3898             IF ( .NOT. particles(n)%age == 0.0_wp .AND.                       &
3899                (ABS( u_int(n) + rvar1_temp(n) ) > (dx/dt_particle(n))  .OR.   &
3900                 ABS( v_int(n) + rvar2_temp(n) ) > (dy/dt_particle(n))  .OR.   &
3901                 ABS( w_int(n) + rvar3_temp(n) ) > (dz_temp/dt_particle(n))))  THEN
3902
3903                dt_particle(n) = 0.9_wp * MIN(                                 &
3904                                 ( dx / ABS( u_int(n) + rvar1_temp(n) ) ),     &
3905                                 ( dy / ABS( v_int(n) + rvar2_temp(n) ) ),     &
3906                                 ( dz_temp / ABS( w_int(n) + rvar3_temp(n) ) ) )
3907
3908!
3909!--             Reset temporary SGS velocites to "current" ones
3910                rvar1_temp(n) = particles(n)%rvar1
3911                rvar2_temp(n) = particles(n)%rvar2
3912                rvar3_temp(n) = particles(n)%rvar3
3913
3914                de_dt_min = - e_int(n) / dt_particle(n)
3915
3916                de_dt = ( e_int(n) - particles(n)%e_m ) / dt_particle_m
3917
3918                IF ( de_dt < de_dt_min )  THEN
3919                   de_dt = de_dt_min
3920                ENDIF
3921
3922                CALL weil_stochastic_eq( rvar1_temp(n), fs_int(n), e_int(n),    &
3923                                        de_dx_int(n), de_dt, diss_int(n),       &
3924                                        dt_particle(n), rg(n,1), term_1_2(n) )
3925
3926                CALL weil_stochastic_eq( rvar2_temp(n), fs_int(n), e_int(n),    &
3927                                        de_dy_int(n), de_dt, diss_int(n),       &
3928                                        dt_particle(n), rg(n,2), term_1_2(n) )
3929
3930                CALL weil_stochastic_eq( rvar3_temp(n), fs_int(n), e_int(n),    &
3931                                        de_dz_int(n), de_dt, diss_int(n),       &
3932                                        dt_particle(n), rg(n,3), term_1_2(n) )
3933             ENDIF
3934
3935!
3936!--          Update particle velocites
3937             particles(n)%rvar1 = rvar1_temp(n)
3938             particles(n)%rvar2 = rvar2_temp(n)
3939             particles(n)%rvar3 = rvar3_temp(n)
3940             u_int(n) = u_int(n) + particles(n)%rvar1
3941             v_int(n) = v_int(n) + particles(n)%rvar2
3942             w_int(n) = w_int(n) + particles(n)%rvar3
3943!
3944!--          Store the SGS TKE of the current timelevel which is needed for
3945!--          for calculating the SGS particle velocities at the next timestep
3946             particles(n)%e_m = e_int(n)
3947          ENDDO
3948       ENDDO
3949
3950    ELSE
3951!
3952!--    If no SGS velocities are used, only the particle timestep has to
3953!--    be set
3954       dt_particle = dt_3d
3955
3956    ENDIF
3957
3958    dens_ratio = particle_groups(particles(1:number_of_particles)%group)%density_ratio
3959    IF ( ANY( dens_ratio == 0.0_wp ) )  THEN
3960!
3961!--    Decide whether the particle loop runs over the subboxes or only over 1,
3962!--    number_of_particles. This depends on the selected interpolation method.
3963!--    If particle interpolation method is not trilinear, then the sorting within
3964!--    subboxes is not required. However, therefore the index start_index(nb) and
3965!--    end_index(nb) are not defined and the loops are still over
3966!--    number_of_particles. @todo find a more generic way to write this loop or
3967!--    delete trilinear interpolation
3968       IF ( interpolation_trilinear )  THEN
3969          subbox_start = 0
3970          subbox_end   = 7
3971       ELSE
3972          subbox_start = 1
3973          subbox_end   = 1
3974       ENDIF
3975!
3976!--    loop over subboxes. In case of simple interpolation scheme no subboxes
3977!--    are introduced, as they are not required. Accordingly, this loops goes
3978!--    from 1 to 1.
3979       DO  nb = subbox_start, subbox_end
3980          IF ( interpolation_trilinear )  THEN
3981             particle_start = start_index(nb)
3982             particle_end   = end_index(nb)
3983          ELSE
3984             particle_start = 1
3985             particle_end   = number_of_particles
3986          ENDIF
3987!
3988!--         Loop from particle start to particle end
3989            DO  n = particle_start, particle_end
3990
3991!
3992!--          Particle advection
3993             IF ( dens_ratio(n) == 0.0_wp )  THEN
3994!
3995!--             Pure passive transport (without particle inertia)
3996                particles(n)%x = xv(n) + u_int(n) * dt_particle(n)
3997                particles(n)%y = yv(n) + v_int(n) * dt_particle(n)
3998                particles(n)%z = zv(n) + w_int(n) * dt_particle(n)
3999
4000                particles(n)%speed_x = u_int(n)
4001                particles(n)%speed_y = v_int(n)
4002                particles(n)%speed_z = w_int(n)
4003
4004             ELSE
4005!
4006!--             Transport of particles with inertia
4007                particles(n)%x = particles(n)%x + particles(n)%speed_x * &
4008                                                  dt_particle(n)
4009                particles(n)%y = particles(n)%y + particles(n)%speed_y * &
4010                                                  dt_particle(n)
4011                particles(n)%z = particles(n)%z + particles(n)%speed_z * &
4012                                                  dt_particle(n)
4013
4014!
4015!--             Update of the particle velocity
4016                IF ( cloud_droplets )  THEN
4017!
4018!--                Terminal velocity is computed for vertical direction (Rogers et
4019!--                al., 1993, J. Appl. Meteorol.)
4020                   diameter = particles(n)%radius * 2000.0_wp !diameter in mm
4021                   IF ( diameter <= d0_rog )  THEN
4022                      w_s = k_cap_rog * diameter * ( 1.0_wp - EXP( -k_low_rog * diameter ) )
4023                   ELSE
4024                      w_s = a_rog - b_rog * EXP( -c_rog * diameter )
4025                   ENDIF
4026
4027!
4028!--                If selected, add random velocities following Soelch and Kaercher
4029!--                (2010, Q. J. R. Meteorol. Soc.)
4030                   IF ( use_sgs_for_particles )  THEN
4031                      lagr_timescale(n) = km(kp,jp,ip) / MAX( e(kp,jp,ip), 1.0E-20_wp )
4032                      RL             = EXP( -1.0_wp * dt_3d / MAX( lagr_timescale(n), &
4033                                             1.0E-20_wp ) )
4034                      sigma          = SQRT( e(kp,jp,ip) )
4035
4036                      rg1 = random_gauss( iran_part, 5.0_wp ) - 1.0_wp
4037                      rg2 = random_gauss( iran_part, 5.0_wp ) - 1.0_wp
4038                      rg3 = random_gauss( iran_part, 5.0_wp ) - 1.0_wp
4039
4040                      particles(n)%rvar1 = RL * particles(n)%rvar1 +              &
4041                                           SQRT( 1.0_wp - RL**2 ) * sigma * rg1
4042                      particles(n)%rvar2 = RL * particles(n)%rvar2 +              &
4043                                           SQRT( 1.0_wp - RL**2 ) * sigma * rg2
4044                      particles(n)%rvar3 = RL * particles(n)%rvar3 +              &
4045                                           SQRT( 1.0_wp - RL**2 ) * sigma * rg3
4046
4047                      particles(n)%speed_x = u_int(n) + particles(n)%rvar1
4048                      particles(n)%speed_y = v_int(n) + particles(n)%rvar2
4049                      particles(n)%speed_z = w_int(n) + particles(n)%rvar3 - w_s
4050                   ELSE
4051                      particles(n)%speed_x = u_int(n)
4052                      particles(n)%speed_y = v_int(n)
4053                      particles(n)%speed_z = w_int(n) - w_s
4054                   ENDIF
4055
4056                ELSE
4057
4058                   IF ( use_sgs_for_particles )  THEN
4059                      exp_arg  = particle_groups(particles(n)%group)%exp_arg
4060                      exp_term = EXP( -exp_arg * dt_particle(n) )
4061                   ELSE
4062                      exp_arg  = particle_groups(particles(n)%group)%exp_arg
4063                      exp_term = particle_groups(particles(n)%group)%exp_term
4064                   ENDIF
4065                   particles(n)%speed_x = particles(n)%speed_x * exp_term +         &
4066                                          u_int(n) * ( 1.0_wp - exp_term )
4067                   particles(n)%speed_y = particles(n)%speed_y * exp_term +         &
4068                                          v_int(n) * ( 1.0_wp - exp_term )
4069                   particles(n)%speed_z = particles(n)%speed_z * exp_term +         &
4070                                          ( w_int(n) - ( 1.0_wp - dens_ratio(n) ) * &
4071                                          g / exp_arg ) * ( 1.0_wp - exp_term )
4072                ENDIF
4073
4074             ENDIF
4075          ENDDO
4076       ENDDO
4077
4078    ELSE
4079!
4080!--    Decide whether the particle loop runs over the subboxes or only over 1,
4081!--    number_of_particles. This depends on the selected interpolation method.
4082       IF ( interpolation_trilinear )  THEN
4083          subbox_start = 0
4084          subbox_end   = 7
4085       ELSE
4086          subbox_start = 1
4087          subbox_end   = 1
4088       ENDIF
4089!--    loop over subboxes. In case of simple interpolation scheme no subboxes
4090!--    are introduced, as they are not required. Accordingly, this loops goes
4091!--    from 1 to 1.
4092       DO  nb = subbox_start, subbox_end
4093          IF ( interpolation_trilinear )  THEN
4094             particle_start = start_index(nb)
4095             particle_end   = end_index(nb)
4096          ELSE
4097             particle_start = 1
4098             particle_end   = number_of_particles
4099          ENDIF
4100!
4101!--         Loop from particle start to particle end
4102            DO  n = particle_start, particle_end
4103
4104!
4105!--          Transport of particles with inertia
4106             particles(n)%x = xv(n) + particles(n)%speed_x * dt_particle(n)
4107             particles(n)%y = yv(n) + particles(n)%speed_y * dt_particle(n)
4108             particles(n)%z = zv(n) + particles(n)%speed_z * dt_particle(n)
4109!
4110!--          Update of the particle velocity
4111             IF ( cloud_droplets )  THEN
4112!
4113!--             Terminal velocity is computed for vertical direction (Rogers et al.,
4114!--             1993, J. Appl. Meteorol.)
4115                diameter = particles(n)%radius * 2000.0_wp !diameter in mm
4116                IF ( diameter <= d0_rog )  THEN
4117                   w_s = k_cap_rog * diameter * ( 1.0_wp - EXP( -k_low_rog * diameter ) )
4118                ELSE
4119                   w_s = a_rog - b_rog * EXP( -c_rog * diameter )
4120                ENDIF
4121
4122!
4123!--             If selected, add random velocities following Soelch and Kaercher
4124!--             (2010, Q. J. R. Meteorol. Soc.)
4125                IF ( use_sgs_for_particles )  THEN
4126                    lagr_timescale(n) = km(kp,jp,ip) / MAX( e(kp,jp,ip), 1.0E-20_wp )
4127                     RL             = EXP( -1.0_wp * dt_3d / MAX( lagr_timescale(n), &
4128                                             1.0E-20_wp ) )
4129                    sigma          = SQRT( e(kp,jp,ip) )
4130
4131                    rg1 = random_gauss( iran_part, 5.0_wp ) - 1.0_wp
4132                    rg2 = random_gauss( iran_part, 5.0_wp ) - 1.0_wp
4133                    rg3 = random_gauss( iran_part, 5.0_wp ) - 1.0_wp
4134
4135                    particles(n)%rvar1 = RL * particles(n)%rvar1 +                &
4136                                         SQRT( 1.0_wp - RL**2 ) * sigma * rg1
4137                    particles(n)%rvar2 = RL * particles(n)%rvar2 +                &
4138                                         SQRT( 1.0_wp - RL**2 ) * sigma * rg2
4139                    particles(n)%rvar3 = RL * particles(n)%rvar3 +                &
4140                                         SQRT( 1.0_wp - RL**2 ) * sigma * rg3
4141
4142                    particles(n)%speed_x = u_int(n) + particles(n)%rvar1
4143                    particles(n)%speed_y = v_int(n) + particles(n)%rvar2
4144                    particles(n)%speed_z = w_int(n) + particles(n)%rvar3 - w_s
4145                ELSE
4146                    particles(n)%speed_x = u_int(n)
4147                    particles(n)%speed_y = v_int(n)
4148                    particles(n)%speed_z = w_int(n) - w_s
4149                ENDIF
4150
4151             ELSE
4152
4153                IF ( use_sgs_for_particles )  THEN
4154                   exp_arg  = particle_groups(particles(n)%group)%exp_arg
4155                   exp_term = EXP( -exp_arg * dt_particle(n) )
4156                ELSE
4157                   exp_arg  = particle_groups(particles(n)%group)%exp_arg
4158                   exp_term = particle_groups(particles(n)%group)%exp_term
4159                ENDIF
4160                particles(n)%speed_x = particles(n)%speed_x * exp_term +             &
4161                                       u_int(n) * ( 1.0_wp - exp_term )
4162                particles(n)%speed_y = particles(n)%speed_y * exp_term +             &
4163                                       v_int(n) * ( 1.0_wp - exp_term )
4164                particles(n)%speed_z = particles(n)%speed_z * exp_term +             &
4165                                       ( w_int(n) - ( 1.0_wp - dens_ratio(n) ) * g / &
4166                                       exp_arg ) * ( 1.0_wp - exp_term )
4167             ENDIF
4168          ENDDO
4169       ENDDO
4170
4171    ENDIF
4172
4173!
4174!-- Store the old age of the particle ( needed to prevent that a
4175!-- particle crosses several PEs during one timestep, and for the
4176!-- evaluation of the subgrid particle velocity fluctuations )
4177    particles(1:number_of_particles)%age_m = particles(1:number_of_particles)%age
4178
4179!
4180!--    loop over subboxes. In case of simple interpolation scheme no subboxes
4181!--    are introduced, as they are not required. Accordingly, this loops goes
4182!--    from 1 to 1.
4183!
4184!-- Decide whether the particle loop runs over the subboxes or only over 1,
4185!-- number_of_particles. This depends on the selected interpolation method.
4186    IF ( interpolation_trilinear )  THEN
4187       subbox_start = 0
4188       subbox_end   = 7
4189    ELSE
4190       subbox_start = 1
4191       subbox_end   = 1
4192    ENDIF
4193    DO  nb = subbox_start, subbox_end
4194       IF ( interpolation_trilinear )  THEN
4195          particle_start = start_index(nb)
4196          particle_end   = end_index(nb)
4197       ELSE
4198          particle_start = 1
4199          particle_end   = number_of_particles
4200       ENDIF
4201!
4202!--    Loop from particle start to particle end
4203       DO  n = particle_start, particle_end
4204!
4205!--       Increment the particle age and the total time that the particle
4206!--       has advanced within the particle timestep procedure
4207          particles(n)%age    = particles(n)%age    + dt_particle(n)
4208          particles(n)%dt_sum = particles(n)%dt_sum + dt_particle(n)
4209
4210!
4211!--       Check whether there is still a particle that has not yet completed
4212!--       the total LES timestep
4213          IF ( ( dt_3d - particles(n)%dt_sum ) > 1E-8_wp )  THEN
4214             dt_3d_reached_l = .FALSE.
4215          ENDIF
4216
4217       ENDDO
4218    ENDDO
4219
4220    CALL cpu_log( log_point_s(44), 'lpm_advec', 'pause' )
4221
4222
4223 END SUBROUTINE lpm_advec
4224
4225 
4226!------------------------------------------------------------------------------! 
4227! Description:
4228! ------------
4229!> Calculation of subgrid-scale particle speed using the stochastic model
4230!> of Weil et al. (2004, JAS, 61, 2877-2887).
4231!------------------------------------------------------------------------------!
4232 SUBROUTINE weil_stochastic_eq( v_sgs, fs_n, e_n, dedxi_n, dedt_n, diss_n,     &
4233                                dt_n, rg_n, fac )
4234
4235    REAL(wp) ::  a1      !< dummy argument
4236    REAL(wp) ::  dedt_n  !< time derivative of TKE at particle position
4237    REAL(wp) ::  dedxi_n !< horizontal derivative of TKE at particle position
4238    REAL(wp) ::  diss_n  !< dissipation at particle position
4239    REAL(wp) ::  dt_n    !< particle timestep
4240    REAL(wp) ::  e_n     !< TKE at particle position
4241    REAL(wp) ::  fac     !< flag to identify adjacent topography
4242    REAL(wp) ::  fs_n    !< weighting factor to prevent that subgrid-scale particle speed becomes too large
4243    REAL(wp) ::  rg_n    !< random number
4244    REAL(wp) ::  term1   !< memory term
4245    REAL(wp) ::  term2   !< drift correction term
4246    REAL(wp) ::  term3   !< random term
4247    REAL(wp) ::  v_sgs   !< subgrid-scale velocity component
4248
4249!-- At first, limit TKE to a small non-zero number, in order to prevent
4250!-- the occurrence of extremely large SGS-velocities in case TKE is zero,
4251!-- (could occur at the simulation begin).
4252    e_n = MAX( e_n, 1E-20_wp )
4253!
4254!-- Please note, terms 1 and 2 (drift and memory term, respectively) are
4255!-- multiplied by a flag to switch of both terms near topography.
4256!-- This is necessary, as both terms may cause a subgrid-scale velocity build up
4257!-- if particles are trapped in regions with very small TKE, e.g. in narrow street
4258!-- canyons resolved by only a few grid points. Hence, term 1 and term 2 are
4259!-- disabled if one of the adjacent grid points belongs to topography.
4260!-- Moreover, in this case, the  previous subgrid-scale component is also set
4261!-- to zero.
4262
4263    a1 = fs_n * c_0 * diss_n
4264!
4265!-- Memory term
4266    term1 = - a1 * v_sgs * dt_n / ( 4.0_wp * sgs_wf_part * e_n + 1E-20_wp )    &
4267                 * fac
4268!
4269!-- Drift correction term
4270    term2 = ( ( dedt_n * v_sgs / e_n ) + dedxi_n ) * 0.5_wp * dt_n              &
4271                 * fac
4272!
4273!-- Random term
4274    term3 = SQRT( MAX( a1, 1E-20_wp ) ) * ( rg_n - 1.0_wp ) * SQRT( dt_n )
4275!
4276!-- In cese one of the adjacent grid-boxes belongs to topograhy, the previous
4277!-- subgrid-scale velocity component is set to zero, in order to prevent a
4278!-- velocity build-up.
4279!-- This case, set also previous subgrid-scale component to zero.
4280    v_sgs = v_sgs * fac + term1 + term2 + term3
4281
4282 END SUBROUTINE weil_stochastic_eq 
4283 
4284 
4285!------------------------------------------------------------------------------! 
4286! Description:
4287! ------------
4288!> Boundary conditions for the Lagrangian particles.
4289!> The routine consists of two different parts. One handles the bottom (flat)
4290!> and top boundary. In this part, also particles which exceeded their lifetime
4291!> are deleted.
4292!> The other part handles the reflection of particles from vertical walls.
4293!> This part was developed by Jin Zhang during 2006-2007.
4294!>
4295!> To do: Code structure for finding the t_index values and for checking the
4296!> -----  reflection conditions is basically the same for all four cases, so it
4297!>        should be possible to further simplify/shorten it.
4298!>
4299!> THE WALLS PART OF THIS ROUTINE HAS NOT BEEN TESTED FOR OCEAN RUNS SO FAR!!!!
4300!> (see offset_ocean_*)
4301!------------------------------------------------------------------------------!
4302 SUBROUTINE lpm_boundary_conds( location_bc , i, j, k )
4303
4304    CHARACTER (LEN=*), INTENT(IN) ::  location_bc !< general mode: boundary conditions at bottom/top of the model domain
4305                                   !< or at vertical surfaces (buildings, terrain steps)   
4306    INTEGER(iwp), INTENT(IN) ::  i !< grid index of particle box along x
4307    INTEGER(iwp), INTENT(IN) ::  j !< grid index of particle box along y
4308    INTEGER(iwp), INTENT(IN) ::  k !< grid index of particle box along z
4309
4310    INTEGER(iwp) ::  inc            !< dummy for sorting algorithmus
4311    INTEGER(iwp) ::  ir             !< dummy for sorting algorithmus
4312    INTEGER(iwp) ::  i1             !< grid index (x) of old particle position
4313    INTEGER(iwp) ::  i2             !< grid index (x) of current particle position
4314    INTEGER(iwp) ::  i3             !< grid index (x) of intermediate particle position
4315    INTEGER(iwp) ::  index_reset    !< index reset height
4316    INTEGER(iwp) ::  jr             !< dummy for sorting algorithmus
4317    INTEGER(iwp) ::  j1             !< grid index (y) of old particle position
4318    INTEGER(iwp) ::  j2             !< grid index (y) of current particle position
4319    INTEGER(iwp) ::  j3             !< grid index (y) of intermediate particle position
4320    INTEGER(iwp) ::  k1             !< grid index (z) of old particle position
4321    INTEGER(iwp) ::  k2             !< grid index (z) of current particle position
4322    INTEGER(iwp) ::  k3             !< grid index (z) of intermediate particle position
4323    INTEGER(iwp) ::  n              !< particle number
4324    INTEGER(iwp) ::  particles_top  !< maximum reset height
4325    INTEGER(iwp) ::  t_index        !< running index for intermediate particle timesteps in reflection algorithmus
4326    INTEGER(iwp) ::  t_index_number !< number of intermediate particle timesteps in reflection algorithmus
4327    INTEGER(iwp) ::  tmp_x          !< dummy for sorting algorithm
4328    INTEGER(iwp) ::  tmp_y          !< dummy for sorting algorithm
4329    INTEGER(iwp) ::  tmp_z          !< dummy for sorting algorithm
4330
4331    INTEGER(iwp), DIMENSION(0:10) ::  x_ind(0:10) = 0 !< index array (x) of intermediate particle positions
4332    INTEGER(iwp), DIMENSION(0:10) ::  y_ind(0:10) = 0 !< index array (y) of intermediate particle positions
4333    INTEGER(iwp), DIMENSION(0:10) ::  z_ind(0:10) = 0 !< index array (z) of intermediate particle positions
4334
4335    LOGICAL  ::  cross_wall_x    !< flag to check if particle reflection along x is necessary
4336    LOGICAL  ::  cross_wall_y    !< flag to check if particle reflection along y is necessary
4337    LOGICAL  ::  cross_wall_z    !< flag to check if particle reflection along z is necessary
4338    LOGICAL  ::  reflect_x       !< flag to check if particle is already reflected along x
4339    LOGICAL  ::  reflect_y       !< flag to check if particle is already reflected along y
4340    LOGICAL  ::  reflect_z       !< flag to check if particle is already reflected along z
4341    LOGICAL  ::  tmp_reach_x     !< dummy for sorting algorithmus
4342    LOGICAL  ::  tmp_reach_y     !< dummy for sorting algorithmus
4343    LOGICAL  ::  tmp_reach_z     !< dummy for sorting algorithmus
4344    LOGICAL  ::  x_wall_reached  !< flag to check if particle has already reached wall
4345    LOGICAL  ::  y_wall_reached  !< flag to check if particle has already reached wall
4346    LOGICAL  ::  z_wall_reached  !< flag to check if particle has already reached wall
4347
4348    LOGICAL, DIMENSION(0:10) ::  reach_x  !< flag to check if particle is at a yz-wall
4349    LOGICAL, DIMENSION(0:10) ::  reach_y  !< flag to check if particle is at a xz-wall
4350    LOGICAL, DIMENSION(0:10) ::  reach_z  !< flag to check if particle is at a xy-wall
4351
4352    REAL(wp) ::  dt_particle    !< particle timestep
4353    REAL(wp) ::  eps = 1E-10_wp !< security number to check if particle has reached a wall
4354    REAL(wp) ::  pos_x          !< intermediate particle position (x)
4355    REAL(wp) ::  pos_x_old      !< particle position (x) at previous particle timestep
4356    REAL(wp) ::  pos_y          !< intermediate particle position (y)
4357    REAL(wp) ::  pos_y_old      !< particle position (y) at previous particle timestep
4358    REAL(wp) ::  pos_z          !< intermediate particle position (z)
4359    REAL(wp) ::  pos_z_old      !< particle position (z) at previous particle timestep
4360    REAL(wp) ::  prt_x          !< current particle position (x)
4361    REAL(wp) ::  prt_y          !< current particle position (y)
4362    REAL(wp) ::  prt_z          !< current particle position (z)
4363    REAL(wp) ::  ran_val        !< location of wall in z
4364    REAL(wp) ::  reset_top      !< location of wall in z
4365    REAL(wp) ::  t_old          !< previous reflection time
4366    REAL(wp) ::  tmp_t          !< dummy for sorting algorithmus
4367    REAL(wp) ::  xwall          !< location of wall in x
4368    REAL(wp) ::  ywall          !< location of wall in y
4369    REAL(wp) ::  zwall          !< location of wall in z
4370
4371    REAL(wp), DIMENSION(0:10) ::  t  !< reflection time
4372
4373    SELECT CASE ( location_bc )
4374
4375       CASE ( 'bottom/top' )
4376
4377!
4378!--    Apply boundary conditions to those particles that have crossed the top or
4379!--    bottom boundary and delete those particles, which are older than allowed
4380       DO  n = 1, number_of_particles
4381
4382!
4383!--       Stop if particles have moved further than the length of one
4384!--       PE subdomain (newly released particles have age = age_m!)
4385          IF ( particles(n)%age /= particles(n)%age_m )  THEN
4386             IF ( ABS(particles(n)%speed_x) >                                  &
4387                  ((nxr-nxl+2)*dx)/(particles(n)%age-particles(n)%age_m)  .OR. &
4388                  ABS(particles(n)%speed_y) >                                  &
4389                  ((nyn-nys+2)*dy)/(particles(n)%age-particles(n)%age_m) )  THEN
4390
4391                  WRITE( message_string, * )  'particle too fast.  n = ',  n 
4392                  CALL message( 'lpm_boundary_conds', 'PA0148', 2, 2, -1, 6, 1 )
4393             ENDIF
4394          ENDIF
4395
4396          IF ( particles(n)%age > particle_maximum_age  .AND.  &
4397               particles(n)%particle_mask )                              &
4398          THEN
4399             particles(n)%particle_mask  = .FALSE.
4400             deleted_particles = deleted_particles + 1
4401          ENDIF
4402
4403          IF ( particles(n)%z >= zw(nz)  .AND.  particles(n)%particle_mask )  THEN
4404             IF ( ibc_par_t == 1 )  THEN
4405!
4406!--             Particle absorption
4407                particles(n)%particle_mask  = .FALSE.
4408                deleted_particles = deleted_particles + 1
4409             ELSEIF ( ibc_par_t == 2 )  THEN
4410!
4411!--             Particle reflection
4412                particles(n)%z       = 2.0_wp * zw(nz) - particles(n)%z
4413                particles(n)%speed_z = -particles(n)%speed_z
4414                IF ( use_sgs_for_particles  .AND. &
4415                     particles(n)%rvar3 > 0.0_wp )  THEN
4416                   particles(n)%rvar3 = -particles(n)%rvar3
4417                ENDIF
4418             ENDIF
4419          ENDIF
4420
4421          IF ( particles(n)%z < zw(0)  .AND.  particles(n)%particle_mask )  THEN
4422             IF ( ibc_par_b == 1 )  THEN
4423!
4424!--             Particle absorption
4425                particles(n)%particle_mask  = .FALSE.
4426                deleted_particles = deleted_particles + 1
4427             ELSEIF ( ibc_par_b == 2 )  THEN
4428!
4429!--             Particle reflection
4430                particles(n)%z       = 2.0_wp * zw(0) - particles(n)%z
4431                particles(n)%speed_z = -particles(n)%speed_z
4432                IF ( use_sgs_for_particles  .AND. &
4433                     particles(n)%rvar3 < 0.0_wp )  THEN
4434                   particles(n)%rvar3 = -particles(n)%rvar3
4435                ENDIF
4436             ELSEIF ( ibc_par_b == 3 )  THEN
4437!
4438!--             Find reset height. @note this works only in non-strechted cases
4439                particles_top = INT( pst(1) / dz(1) )
4440                index_reset = MINLOC( prt_count(nzb+1:particles_top,j,i), DIM = 1 )
4441                reset_top = zu(index_reset)
4442                iran_part = iran_part + myid
4443                ran_val = random_function( iran_part )
4444                particles(n)%z       = reset_top *  ( 1.0  + ( ran_val / 10.0_wp) )
4445                particles(n)%speed_z = 0.0_wp
4446                IF ( curvature_solution_effects )  THEN
4447                   particles(n)%radius = particles(n)%aux1
4448                ELSE
4449                   particles(n)%radius = 1.0E-8
4450                ENDIF
4451             ENDIF
4452          ENDIF
4453       ENDDO
4454
4455      CASE ( 'walls' )
4456
4457       CALL cpu_log( log_point_s(48), 'lpm_wall_reflect', 'start' )
4458
4459       DO  n = 1, number_of_particles
4460!
4461!--       Recalculate particle timestep
4462          dt_particle = particles(n)%age - particles(n)%age_m
4463!
4464!--       Obtain x/y indices for current particle position
4465          i2 = particles(n)%x * ddx
4466          j2 = particles(n)%y * ddy
4467          IF ( zw(k)   < particles(n)%z ) k2 = k + 1
4468          IF ( zw(k)   > particles(n)%z  .AND.  zw(k-1) < particles(n)%z ) k2 = k
4469          IF ( zw(k-1) > particles(n)%z ) k2 = k - 1
4470!
4471!--       Save current particle positions
4472          prt_x = particles(n)%x
4473          prt_y = particles(n)%y
4474          prt_z = particles(n)%z
4475!
4476!--       Recalculate old particle positions
4477          pos_x_old = particles(n)%x - particles(n)%speed_x * dt_particle
4478          pos_y_old = particles(n)%y - particles(n)%speed_y * dt_particle
4479          pos_z_old = particles(n)%z - particles(n)%speed_z * dt_particle
4480!
4481!--       Obtain x/y indices for old particle positions
4482          i1 = i
4483          j1 = j
4484          k1 = k
4485!
4486!--       Determine horizontal as well as vertical walls at which particle can
4487!--       be potentially reflected.
4488!--       Start with walls aligned in yz layer.
4489!--       Wall to the right
4490          IF ( prt_x > pos_x_old )  THEN
4491             xwall = ( i1 + 1 ) * dx
4492!
4493!--       Wall to the left
4494          ELSE
4495             xwall = i1 * dx
4496          ENDIF
4497!
4498!--       Walls aligned in xz layer
4499!--       Wall to the north
4500          IF ( prt_y > pos_y_old )  THEN
4501             ywall = ( j1 + 1 ) * dy
4502!--       Wall to the south
4503          ELSE
4504             ywall = j1 * dy
4505          ENDIF
4506
4507          IF ( prt_z > pos_z_old )  THEN
4508             zwall = zw(k)
4509          ELSE
4510             zwall = zw(k-1)
4511          ENDIF
4512!
4513!--       Initialize flags to check if particle reflection is necessary
4514          cross_wall_x = .FALSE.
4515          cross_wall_y = .FALSE.
4516          cross_wall_z = .FALSE.
4517!
4518!--       Initialize flags to check if a wall is reached
4519          reach_x      = .FALSE.
4520          reach_y      = .FALSE.
4521          reach_z      = .FALSE.
4522!
4523!--       Initialize flags to check if a particle was already reflected
4524          reflect_x    = .FALSE.
4525          reflect_y    = .FALSE.
4526          reflect_z    = .FALSE.
4527!
4528!--       Initialize flags to check if a wall is already crossed.
4529!--       ( Required to obtain correct indices. )
4530          x_wall_reached = .FALSE.
4531          y_wall_reached = .FALSE.
4532          z_wall_reached = .FALSE.
4533!
4534!--       Initialize time array
4535          t     = 0.0_wp
4536!
4537!--       Check if particle can reach any wall. This case, calculate the
4538!--       fractional time needed to reach this wall. Store this fractional
4539!--       timestep in array t. Moreover, store indices for these grid
4540!--       boxes where the respective wall belongs to. 
4541!--       Start with x-direction.
4542          t_index    = 1
4543          t(t_index) = ( xwall - pos_x_old )                                   &
4544                     / MERGE( MAX( prt_x - pos_x_old,  1E-30_wp ),             &
4545                              MIN( prt_x - pos_x_old, -1E-30_wp ),             &
4546                              prt_x > pos_x_old )
4547          x_ind(t_index)   = i2
4548          y_ind(t_index)   = j1
4549          z_ind(t_index)   = k1
4550          reach_x(t_index) = .TRUE.
4551          reach_y(t_index) = .FALSE.
4552          reach_z(t_index) = .FALSE.
4553!
4554!--       Store these values only if particle really reaches any wall. t must
4555!--       be in a interval between [0:1].
4556          IF ( t(t_index) <= 1.0_wp  .AND.  t(t_index) >= 0.0_wp )  THEN
4557             t_index      = t_index + 1
4558             cross_wall_x = .TRUE.
4559          ENDIF
4560!
4561!--       y-direction
4562          t(t_index) = ( ywall - pos_y_old )                                   &
4563                     / MERGE( MAX( prt_y - pos_y_old,  1E-30_wp ),             &
4564                              MIN( prt_y - pos_y_old, -1E-30_wp ),             &
4565                              prt_y > pos_y_old )
4566          x_ind(t_index)   = i1
4567          y_ind(t_index)   = j2
4568          z_ind(t_index)   = k1
4569          reach_x(t_index) = .FALSE.
4570          reach_y(t_index) = .TRUE.
4571          reach_z(t_index) = .FALSE.
4572          IF ( t(t_index) <= 1.0_wp  .AND.  t(t_index) >= 0.0_wp )  THEN
4573             t_index      = t_index + 1
4574             cross_wall_y = .TRUE.
4575          ENDIF
4576!
4577!--       z-direction
4578          t(t_index) = (zwall - pos_z_old )                                    &
4579                     / MERGE( MAX( prt_z - pos_z_old,  1E-30_wp ),             &
4580                              MIN( prt_z - pos_z_old, -1E-30_wp ),             &
4581                              prt_z > pos_z_old )
4582
4583          x_ind(t_index)   = i1
4584          y_ind(t_index)   = j1
4585          z_ind(t_index)   = k2
4586          reach_x(t_index) = .FALSE.
4587          reach_y(t_index) = .FALSE.
4588          reach_z(t_index) = .TRUE.
4589          IF( t(t_index) <= 1.0_wp  .AND.  t(t_index) >= 0.0_wp)  THEN
4590             t_index      = t_index + 1
4591             cross_wall_z = .TRUE.
4592          ENDIF
4593
4594          t_index_number = t_index - 1
4595!
4596!--       Carry out reflection only if particle reaches any wall
4597          IF ( cross_wall_x  .OR.  cross_wall_y  .OR.  cross_wall_z )  THEN
4598!
4599!--          Sort fractional timesteps in ascending order. Also sort the
4600!--          corresponding indices and flag according to the time interval a 
4601!--          particle reaches the respective wall.
4602             inc = 1
4603             jr  = 1
4604             DO WHILE ( inc <= t_index_number )
4605                inc = 3 * inc + 1
4606             ENDDO
4607
4608             DO WHILE ( inc > 1 )
4609                inc = inc / 3
4610                DO  ir = inc+1, t_index_number
4611                   tmp_t       = t(ir)
4612                   tmp_x       = x_ind(ir)
4613                   tmp_y       = y_ind(ir)
4614                   tmp_z       = z_ind(ir)
4615                   tmp_reach_x = reach_x(ir)
4616                   tmp_reach_y = reach_y(ir)
4617                   tmp_reach_z = reach_z(ir)
4618                   jr    = ir
4619                   DO WHILE ( t(jr-inc) > tmp_t )
4620                      t(jr)       = t(jr-inc)
4621                      x_ind(jr)   = x_ind(jr-inc)
4622                      y_ind(jr)   = y_ind(jr-inc)
4623                      z_ind(jr)   = z_ind(jr-inc)
4624                      reach_x(jr) = reach_x(jr-inc)
4625                      reach_y(jr) = reach_y(jr-inc)
4626                      reach_z(jr) = reach_z(jr-inc)
4627                      jr    = jr - inc
4628                      IF ( jr <= inc )  EXIT
4629                   ENDDO
4630                   t(jr)       = tmp_t
4631                   x_ind(jr)   = tmp_x
4632                   y_ind(jr)   = tmp_y
4633                   z_ind(jr)   = tmp_z
4634                   reach_x(jr) = tmp_reach_x
4635                   reach_y(jr) = tmp_reach_y
4636                   reach_z(jr) = tmp_reach_z
4637                ENDDO
4638             ENDDO
4639!
4640!--          Initialize temporary particle positions
4641             pos_x = pos_x_old
4642             pos_y = pos_y_old
4643             pos_z = pos_z_old
4644!
4645!--          Loop over all times a particle possibly moves into a new grid box
4646             t_old = 0.0_wp
4647             DO t_index = 1, t_index_number 
4648!
4649!--             Calculate intermediate particle position according to the
4650!--             timesteps a particle reaches any wall.
4651                pos_x = pos_x + ( t(t_index) - t_old ) * dt_particle           &
4652                                                       * particles(n)%speed_x
4653                pos_y = pos_y + ( t(t_index) - t_old ) * dt_particle           &
4654                                                       * particles(n)%speed_y
4655                pos_z = pos_z + ( t(t_index) - t_old ) * dt_particle           &
4656                                                       * particles(n)%speed_z
4657!
4658!--             Obtain x/y grid indices for intermediate particle position from
4659!--             sorted index array
4660                i3 = x_ind(t_index)
4661                j3 = y_ind(t_index)
4662                k3 = z_ind(t_index)
4663!
4664!--             Check which wall is already reached
4665                IF ( .NOT. x_wall_reached )  x_wall_reached = reach_x(t_index) 
4666                IF ( .NOT. y_wall_reached )  y_wall_reached = reach_y(t_index)
4667                IF ( .NOT. z_wall_reached )  z_wall_reached = reach_z(t_index)
4668!
4669!--             Check if a particle needs to be reflected at any yz-wall. If
4670!--             necessary, carry out reflection. Please note, a security
4671!--             constant is required, as the particle position does not
4672!--             necessarily exactly match the wall location due to rounding
4673!--             errors.
4674                IF ( reach_x(t_index)                      .AND.               & 
4675                     ABS( pos_x - xwall ) < eps            .AND.               &
4676                     .NOT. BTEST(wall_flags_0(k3,j3,i3),0) .AND.               &
4677                     .NOT. reflect_x )  THEN
4678!
4679!
4680!--                Reflection in x-direction.
4681!--                Ensure correct reflection by MIN/MAX functions, depending on
4682!--                direction of particle transport.
4683!--                Due to rounding errors pos_x does not exactly match the wall
4684!--                location, leading to erroneous reflection.             
4685                   pos_x = MERGE( MIN( 2.0_wp * xwall - pos_x, xwall ),        &
4686                                  MAX( 2.0_wp * xwall - pos_x, xwall ),        &
4687                                  particles(n)%x > xwall )
4688!
4689!--                Change sign of particle speed                     
4690                   particles(n)%speed_x = - particles(n)%speed_x
4691!
4692!--                Also change sign of subgrid-scale particle speed
4693                   particles(n)%rvar1 = - particles(n)%rvar1
4694!
4695!--                Set flag that reflection along x is already done
4696                   reflect_x          = .TRUE.
4697!
4698!--                As the particle does not cross any further yz-wall during
4699!--                this timestep, set further x-indices to the current one.
4700                   x_ind(t_index:t_index_number) = i1
4701!
4702!--             If particle already reached the wall but was not reflected,
4703!--             set further x-indices to the new one.
4704                ELSEIF ( x_wall_reached .AND. .NOT. reflect_x )  THEN
4705                    x_ind(t_index:t_index_number) = i2
4706                ENDIF !particle reflection in x direction done
4707
4708!
4709!--             Check if a particle needs to be reflected at any xz-wall. If
4710!--             necessary, carry out reflection. Please note, a security
4711!--             constant is required, as the particle position does not
4712!--             necessarily exactly match the wall location due to rounding
4713!--             errors.
4714                IF ( reach_y(t_index)                      .AND.               & 
4715                     ABS( pos_y - ywall ) < eps            .AND.               &
4716                     .NOT. BTEST(wall_flags_0(k3,j3,i3),0) .AND.               &
4717                     .NOT. reflect_y )  THEN
4718!
4719!
4720!--                Reflection in y-direction.
4721!--                Ensure correct reflection by MIN/MAX functions, depending on
4722!--                direction of particle transport.
4723!--                Due to rounding errors pos_y does not exactly match the wall
4724!--                location, leading to erroneous reflection.             
4725                   pos_y = MERGE( MIN( 2.0_wp * ywall - pos_y, ywall ),        &
4726                                  MAX( 2.0_wp * ywall - pos_y, ywall ),        &
4727                                  particles(n)%y > ywall )
4728!
4729!--                Change sign of particle speed                     
4730                   particles(n)%speed_y = - particles(n)%speed_y
4731!
4732!--                Also change sign of subgrid-scale particle speed
4733                   particles(n)%rvar2 = - particles(n)%rvar2
4734!
4735!--                Set flag that reflection along y is already done
4736                   reflect_y          = .TRUE.
4737!
4738!--                As the particle does not cross any further xz-wall during
4739!--                this timestep, set further y-indices to the current one.
4740                   y_ind(t_index:t_index_number) = j1
4741!
4742!--             If particle already reached the wall but was not reflected,
4743!--             set further y-indices to the new one.
4744                ELSEIF ( y_wall_reached .AND. .NOT. reflect_y )  THEN
4745                    y_ind(t_index:t_index_number) = j2
4746                ENDIF !particle reflection in y direction done
4747
4748!
4749!--             Check if a particle needs to be reflected at any xy-wall. If
4750!--             necessary, carry out reflection. Please note, a security
4751!--             constant is required, as the particle position does not
4752!--             necessarily exactly match the wall location due to rounding
4753!--             errors.
4754                IF ( reach_z(t_index)                      .AND.               & 
4755                     ABS( pos_z - zwall ) < eps            .AND.               &
4756                     .NOT. BTEST(wall_flags_0(k3,j3,i3),0) .AND.               &
4757                     .NOT. reflect_z )  THEN
4758!
4759!
4760!--                Reflection in z-direction.
4761!--                Ensure correct reflection by MIN/MAX functions, depending on
4762!--                direction of particle transport.
4763!--                Due to rounding errors pos_z does not exactly match the wall
4764!--                location, leading to erroneous reflection.             
4765                   pos_z = MERGE( MIN( 2.0_wp * zwall - pos_z, zwall ),        &
4766                                  MAX( 2.0_wp * zwall - pos_z, zwall ),        &
4767                                  particles(n)%z > zwall )
4768!
4769!--                Change sign of particle speed                     
4770                   particles(n)%speed_z = - particles(n)%speed_z
4771!
4772!--                Also change sign of subgrid-scale particle speed
4773                   particles(n)%rvar3 = - particles(n)%rvar3
4774!
4775!--                Set flag that reflection along z is already done
4776                   reflect_z          = .TRUE.
4777!
4778!--                As the particle does not cross any further xy-wall during
4779!--                this timestep, set further z-indices to the current one.
4780                   z_ind(t_index:t_index_number) = k1
4781!
4782!--             If particle already reached the wall but was not reflected,
4783!--             set further z-indices to the new one.
4784                ELSEIF ( z_wall_reached .AND. .NOT. reflect_z )  THEN
4785                    z_ind(t_index:t_index_number) = k2
4786                ENDIF !particle reflection in z direction done               
4787
4788!
4789!--             Swap time
4790                t_old = t(t_index)
4791
4792             ENDDO
4793!
4794!--          If a particle was reflected, calculate final position from last
4795!--          intermediate position.
4796             IF ( reflect_x  .OR.  reflect_y  .OR.  reflect_z )  THEN
4797
4798                particles(n)%x = pos_x + ( 1.0_wp - t_old ) * dt_particle      &
4799                                                         * particles(n)%speed_x
4800                particles(n)%y = pos_y + ( 1.0_wp - t_old ) * dt_particle      &
4801                                                         * particles(n)%speed_y
4802                particles(n)%z = pos_z + ( 1.0_wp - t_old ) * dt_particle      &
4803                                                         * particles(n)%speed_z
4804
4805             ENDIF
4806
4807          ENDIF
4808
4809       ENDDO
4810
4811       CALL cpu_log( log_point_s(48), 'lpm_wall_reflect', 'stop' )
4812
4813       CASE DEFAULT
4814          CONTINUE
4815
4816    END SELECT
4817
4818 END SUBROUTINE lpm_boundary_conds 
4819
4820
4821!------------------------------------------------------------------------------!
4822! Description:
4823! ------------
4824!> Calculates change in droplet radius by condensation/evaporation, using
4825!> either an analytic formula or by numerically integrating the radius growth
4826!> equation including curvature and solution effects using Rosenbrocks method
4827!> (see Numerical recipes in FORTRAN, 2nd edition, p. 731).
4828!> The analytical formula and growth equation follow those given in
4829!> Rogers and Yau (A short course in cloud physics, 3rd edition, p. 102/103).
4830!------------------------------------------------------------------------------!
4831 SUBROUTINE lpm_droplet_condensation (i,j,k)
4832
4833    INTEGER(iwp), INTENT(IN) ::  i              !<
4834    INTEGER(iwp), INTENT(IN) ::  j              !<
4835    INTEGER(iwp), INTENT(IN) ::  k              !<
4836    INTEGER(iwp) ::  n                          !<
4837
4838    REAL(wp) ::  afactor                       !< curvature effects
4839    REAL(wp) ::  arg                           !<
4840    REAL(wp) ::  bfactor                       !< solute effects
4841    REAL(wp) ::  ddenom                        !<
4842    REAL(wp) ::  delta_r                       !<
4843    REAL(wp) ::  diameter                      !< diameter of cloud droplets
4844    REAL(wp) ::  diff_coeff                    !< diffusivity for water vapor
4845    REAL(wp) ::  drdt                          !<
4846    REAL(wp) ::  dt_ros                        !<
4847    REAL(wp) ::  dt_ros_sum                    !<
4848    REAL(wp) ::  d2rdtdr                       !<
4849    REAL(wp) ::  e_a                           !< current vapor pressure
4850    REAL(wp) ::  e_s                           !< current saturation vapor pressure
4851    REAL(wp) ::  error                         !< local truncation error in Rosenbrock
4852    REAL(wp) ::  k1                            !<
4853    REAL(wp) ::  k2                            !<
4854    REAL(wp) ::  r_err                         !< First order estimate of Rosenbrock radius
4855    REAL(wp) ::  r_ros                         !< Rosenbrock radius
4856    REAL(wp) ::  r_ros_ini                     !< initial Rosenbrock radius
4857    REAL(wp) ::  r0                            !< gas-kinetic lengthscale
4858    REAL(wp) ::  sigma                         !< surface tension of water
4859    REAL(wp) ::  thermal_conductivity          !< thermal conductivity for water
4860    REAL(wp) ::  t_int                         !< temperature
4861    REAL(wp) ::  w_s                           !< terminal velocity of droplets
4862    REAL(wp) ::  re_p                          !< particle Reynolds number
4863!
4864!-- Parameters for Rosenbrock method (see Verwer et al., 1999)
4865    REAL(wp), PARAMETER ::  prec = 1.0E-3_wp     !< precision of Rosenbrock solution
4866    REAL(wp), PARAMETER ::  q_increase = 1.5_wp  !< increase factor in timestep
4867    REAL(wp), PARAMETER ::  q_decrease = 0.9_wp  !< decrease factor in timestep
4868    REAL(wp), PARAMETER ::  gamma = 0.292893218814_wp !< = 1.0 - 1.0 / SQRT(2.0)
4869!
4870!-- Parameters for terminal velocity
4871    REAL(wp), PARAMETER ::  a_rog = 9.65_wp      !< parameter for fall velocity
4872    REAL(wp), PARAMETER ::  b_rog = 10.43_wp     !< parameter for fall velocity
4873    REAL(wp), PARAMETER ::  c_rog = 0.6_wp       !< parameter for fall velocity
4874    REAL(wp), PARAMETER ::  k_cap_rog = 4.0_wp   !< parameter for fall velocity
4875    REAL(wp), PARAMETER ::  k_low_rog = 12.0_wp  !< parameter for fall velocity
4876    REAL(wp), PARAMETER ::  d0_rog = 0.745_wp    !< separation diameter
4877
4878    REAL(wp), DIMENSION(number_of_particles) ::  ventilation_effect     !<
4879    REAL(wp), DIMENSION(number_of_particles) ::  new_r                  !<
4880
4881    CALL cpu_log( log_point_s(42), 'lpm_droplet_condens', 'start' )
4882
4883!
4884!-- Absolute temperature
4885    t_int = pt(k,j,i) * exner(k)
4886!
4887!-- Saturation vapor pressure (Eq. 10 in Bolton, 1980)
4888    e_s = magnus( t_int )
4889!
4890!-- Current vapor pressure
4891    e_a = q(k,j,i) * hyp(k) / ( q(k,j,i) + rd_d_rv )
4892!
4893!-- Thermal conductivity for water (from Rogers and Yau, Table 7.1)
4894    thermal_conductivity = 7.94048E-05_wp * t_int + 0.00227011_wp
4895!
4896!-- Moldecular diffusivity of water vapor in air (Hall und Pruppacher, 1976)
4897    diff_coeff           = 0.211E-4_wp * ( t_int / 273.15_wp )**1.94_wp * &
4898                           ( 101325.0_wp / hyp(k) )
4899!
4900!-- Lengthscale for gas-kinetic effects (from Mordy, 1959, p. 23):
4901    r0 = diff_coeff / 0.036_wp * SQRT( 2.0_wp * pi / ( r_v * t_int ) )
4902!
4903!-- Calculate effects of heat conductivity and diffusion of water vapor on the
4904!-- diffusional growth process (usually known as 1.0 / (F_k + F_d) )
4905    ddenom  = 1.0_wp / ( rho_l * r_v * t_int / ( e_s * diff_coeff ) +          &
4906                         ( l_v / ( r_v * t_int ) - 1.0_wp ) * rho_l *          &
4907                         l_v / ( thermal_conductivity * t_int )                &
4908                       )
4909    new_r = 0.0_wp
4910!
4911!-- Determine ventilation effect on evaporation of large drops
4912    DO  n = 1, number_of_particles
4913
4914       IF ( particles(n)%radius >= 4.0E-5_wp  .AND.  e_a / e_s < 1.0_wp )  THEN
4915!
4916!--       Terminal velocity is computed for vertical direction (Rogers et al.,
4917!--       1993, J. Appl. Meteorol.)
4918          diameter = particles(n)%radius * 2000.0_wp !diameter in mm
4919          IF ( diameter <= d0_rog )  THEN
4920             w_s = k_cap_rog * diameter * ( 1.0_wp - EXP( -k_low_rog * diameter ) )
4921          ELSE
4922             w_s = a_rog - b_rog * EXP( -c_rog * diameter )
4923          ENDIF
4924!
4925!--       Calculate droplet's Reynolds number
4926          re_p = 2.0_wp * particles(n)%radius * w_s / molecular_viscosity
4927!
4928!--       Ventilation coefficient (Rogers and Yau, 1989):
4929          IF ( re_p > 2.5_wp )  THEN
4930             ventilation_effect(n) = 0.78_wp + 0.28_wp * SQRT( re_p )
4931          ELSE
4932             ventilation_effect(n) = 1.0_wp + 0.09_wp * re_p
4933          ENDIF
4934       ELSE
4935!
4936!--       For small droplets or in supersaturated environments, the ventilation
4937!--       effect does not play a role
4938          ventilation_effect(n) = 1.0_wp
4939       ENDIF
4940    ENDDO
4941
4942    IF( .NOT. curvature_solution_effects )  THEN
4943!
4944!--    Use analytic model for diffusional growth including gas-kinetic
4945!--    effects (Mordy, 1959) but without the impact of aerosols.
4946       DO  n = 1, number_of_particles
4947          arg      = ( particles(n)%radius + r0 )**2 + 2.0_wp * dt_3d * ddenom * &
4948                                                       ventilation_effect(n) *   &
4949                                                       ( e_a / e_s - 1.0_wp )
4950          arg      = MAX( arg, ( 0.01E-6 + r0 )**2 )
4951          new_r(n) = SQRT( arg ) - r0
4952       ENDDO
4953
4954    ELSE
4955!
4956!--    Integrate the diffusional growth including gas-kinetic (Mordy, 1959),
4957!--    as well as curvature and solute effects (e.g., Köhler, 1936).
4958!
4959!--    Curvature effect (afactor) with surface tension (sigma) by Straka (2009)
4960       sigma = 0.0761_wp - 0.000155_wp * ( t_int - 273.15_wp )
4961!
4962!--    Solute effect (afactor)
4963       afactor = 2.0_wp * sigma / ( rho_l * r_v * t_int )
4964
4965       DO  n = 1, number_of_particles
4966!
4967!--       Solute effect (bfactor)
4968          bfactor = vanthoff * rho_s * particles(n)%aux1**3 *                    &
4969                    molecular_weight_of_water / ( rho_l * molecular_weight_of_solute )
4970
4971          dt_ros     = particles(n)%aux2  ! use previously stored Rosenbrock timestep
4972          dt_ros_sum = 0.0_wp
4973
4974          r_ros     = particles(n)%radius  ! initialize Rosenbrock particle radius
4975          r_ros_ini = r_ros
4976!
4977!--       Integrate growth equation using a 2nd-order Rosenbrock method
4978!--       (see Verwer et al., 1999, Eq. (3.2)). The Rosenbrock method adjusts
4979!--       its with internal timestep to minimize the local truncation error.
4980          DO WHILE ( dt_ros_sum < dt_3d )
4981
4982             dt_ros = MIN( dt_ros, dt_3d - dt_ros_sum )
4983
4984             DO
4985
4986                drdt = ddenom * ventilation_effect(n) * ( e_a / e_s - 1.0_wp - &
4987                                                          afactor / r_ros +    &
4988                                                          bfactor / r_ros**3   &
4989                                                        ) / ( r_ros + r0 )
4990
4991                d2rdtdr = -ddenom * ventilation_effect(n) * (                  &
4992                                            (e_a / e_s - 1.0_wp ) * r_ros**4 - &
4993                                            afactor * r0 * r_ros**2 -          &
4994                                            2.0_wp * afactor * r_ros**3 +      &
4995                                            3.0_wp * bfactor * r0 +            &
4996                                            4.0_wp * bfactor * r_ros           &
4997                                                            )                  &
4998                          / ( r_ros**4 * ( r_ros + r0 )**2 )
4999
5000                k1    = drdt / ( 1.0_wp - gamma * dt_ros * d2rdtdr )
5001
5002                r_ros = MAX(r_ros_ini + k1 * dt_ros, particles(n)%aux1)
5003                r_err = r_ros
5004
5005                drdt  = ddenom * ventilation_effect(n) * ( e_a / e_s - 1.0_wp - &
5006                                                           afactor / r_ros +    &
5007                                                           bfactor / r_ros**3   &
5008                                                         ) / ( r_ros + r0 )
5009
5010                k2 = ( drdt - dt_ros * 2.0 * gamma * d2rdtdr * k1 ) / &
5011                     ( 1.0_wp - dt_ros * gamma * d2rdtdr )
5012
5013                r_ros = MAX(r_ros_ini + dt_ros * ( 1.5_wp * k1 + 0.5_wp * k2), particles(n)%aux1)
5014   !
5015   !--          Check error of the solution, and reduce dt_ros if necessary.
5016                error = ABS(r_err - r_ros) / r_ros
5017                IF ( error > prec )  THEN
5018                   dt_ros = SQRT( q_decrease * prec / error ) * dt_ros
5019                   r_ros  = r_ros_ini
5020                ELSE
5021                   dt_ros_sum = dt_ros_sum + dt_ros
5022                   dt_ros     = q_increase * dt_ros
5023                   r_ros_ini  = r_ros
5024                   EXIT
5025                ENDIF
5026
5027             END DO
5028
5029          END DO !Rosenbrock loop
5030!
5031!--       Store new particle radius
5032          new_r(n) = r_ros
5033!
5034!--       Store internal time step value for next PALM step
5035          particles(n)%aux2 = dt_ros
5036
5037       ENDDO !Particle loop
5038
5039    ENDIF
5040
5041    DO  n = 1, number_of_particles
5042!
5043!--    Sum up the change in liquid water for the respective grid
5044!--    box for the computation of the release/depletion of water vapor
5045!--    and heat.
5046       ql_c(k,j,i) = ql_c(k,j,i) + particles(n)%weight_factor *          &
5047                                   rho_l * 1.33333333_wp * pi *                &
5048                                   ( new_r(n)**3 - particles(n)%radius**3 ) /  &
5049                                   ( rho_surface * dx * dy * dzw(k) )
5050!
5051!--    Check if the increase in liqid water is not too big. If this is the case,
5052!--    the model timestep might be too long.
5053       IF ( ql_c(k,j,i) > 100.0_wp )  THEN
5054          WRITE( message_string, * ) 'k=',k,' j=',j,' i=',i,                &
5055                       ' ql_c=',ql_c(k,j,i), '&part(',n,')%wf=',            &
5056                       particles(n)%weight_factor,' delta_r=',delta_r
5057          CALL message( 'lpm_droplet_condensation', 'PA0143', 2, 2, -1, 6, 1 )
5058       ENDIF
5059!
5060!--    Check if the change in the droplet radius is not too big. If this is the
5061!--    case, the model timestep might be too long.
5062       delta_r = new_r(n) - particles(n)%radius
5063       IF ( delta_r < 0.0_wp  .AND.  new_r(n) < 0.0_wp )  THEN
5064          WRITE( message_string, * ) '#1 k=',k,' j=',j,' i=',i,                &
5065                       ' e_s=',e_s, ' e_a=',e_a,' t_int=',t_int,               &
5066                       '&delta_r=',delta_r,                                    &
5067                       ' particle_radius=',particles(n)%radius
5068          CALL message( 'lpm_droplet_condensation', 'PA0144', 2, 2, -1, 6, 1 )
5069       ENDIF
5070!
5071!--    Sum up the total volume of liquid water (needed below for
5072!--    re-calculating the weighting factors)
5073       ql_v(k,j,i) = ql_v(k,j,i) + particles(n)%weight_factor * new_r(n)**3
5074!
5075!--    Determine radius class of the particle needed for collision
5076       IF ( use_kernel_tables )  THEN
5077          particles(n)%class = ( LOG( new_r(n) ) - rclass_lbound ) /           &
5078                               ( rclass_ubound - rclass_lbound ) *             &
5079                               radius_classes
5080          particles(n)%class = MIN( particles(n)%class, radius_classes )
5081          particles(n)%class = MAX( particles(n)%class, 1 )
5082       ENDIF
5083!
5084!--    Store new radius to particle features
5085       particles(n)%radius = new_r(n)
5086
5087    ENDDO
5088
5089    CALL cpu_log( log_point_s(42), 'lpm_droplet_condens', 'stop' )
5090
5091
5092 END SUBROUTINE lpm_droplet_condensation
5093
5094
5095!------------------------------------------------------------------------------!
5096! Description:
5097! ------------
5098!> Release of latent heat and change of mixing ratio due to condensation /
5099!> evaporation of droplets.
5100!------------------------------------------------------------------------------!
5101 SUBROUTINE lpm_interaction_droplets_ptq
5102
5103    INTEGER(iwp) ::  i    !< running index x direction
5104    INTEGER(iwp) ::  j    !< running index y direction
5105    INTEGER(iwp) ::  k    !< running index z direction
5106
5107    REAL(wp) ::  flag     !< flag to mask topography grid points
5108
5109    DO  i = nxl, nxr
5110       DO  j = nys, nyn
5111          DO  k = nzb+1, nzt
5112!
5113!--          Predetermine flag to mask topography
5114             flag = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_0(k,j,i), 0 ) )
5115
5116             q_p(k,j,i)  = q_p(k,j,i)  - ql_c(k,j,i) * flag
5117             pt_p(k,j,i) = pt_p(k,j,i) + lv_d_cp * ql_c(k,j,i) * d_exner(k) &
5118                                                     * flag
5119          ENDDO
5120       ENDDO
5121    ENDDO
5122
5123 END SUBROUTINE lpm_interaction_droplets_ptq
5124
5125
5126!------------------------------------------------------------------------------!
5127! Description:
5128! ------------
5129!> Release of latent heat and change of mixing ratio due to condensation /
5130!> evaporation of droplets. Call for grid point i,j
5131!------------------------------------------------------------------------------!
5132 SUBROUTINE lpm_interaction_droplets_ptq_ij( i, j )
5133
5134    INTEGER(iwp) ::  i    !< running index x direction
5135    INTEGER(iwp) ::  j    !< running index y direction
5136    INTEGER(iwp) ::  k    !< running index z direction
5137
5138    REAL(wp) ::  flag     !< flag to mask topography grid points
5139
5140
5141    DO  k = nzb+1, nzt
5142!
5143!--    Predetermine flag to mask topography
5144       flag = MERGE( 1.0_wp, 0.0_wp, BTEST( wall_flags_0(k,j,i), 0 ) )
5145
5146       q_p(k,j,i)  = q_p(k,j,i)  - ql_c(k,j,i) * flag
5147       pt_p(k,j,i) = pt_p(k,j,i) + lv_d_cp * ql_c(k,j,i) * d_exner(k) * flag
5148    ENDDO
5149
5150 END SUBROUTINE lpm_interaction_droplets_ptq_ij
5151
5152
5153!------------------------------------------------------------------------------!
5154! Description:
5155! ------------
5156!> Calculate the liquid water content for each grid box.
5157!------------------------------------------------------------------------------!
5158 SUBROUTINE lpm_calc_liquid_water_content
5159
5160
5161    INTEGER(iwp) ::  i   !<
5162    INTEGER(iwp) ::  j   !<
5163    INTEGER(iwp) ::  k   !<
5164    INTEGER(iwp) ::  n   !<
5165
5166    CALL cpu_log( log_point_s(45), 'lpm_calc_ql', 'start' )
5167
5168!
5169!-- Set water content initially to zero
5170    ql = 0.0_wp;  ql_v = 0.0_wp;  ql_vp = 0.0_wp
5171
5172!
5173!-- Calculate for each grid box
5174    DO  i = nxl, nxr
5175       DO  j = nys, nyn
5176          DO  k = nzb+1, nzt
5177             number_of_particles = prt_count(k,j,i)
5178             IF ( number_of_particles <= 0 )  CYCLE
5179             particles => grid_particles(k,j,i)%particles(1:number_of_particles)
5180!
5181!--          Calculate the total volume in the boxes (ql_v, weighting factor
5182!--          has to beincluded)
5183             DO  n = 1, prt_count(k,j,i)
5184                ql_v(k,j,i)  = ql_v(k,j,i)  + particles(n)%weight_factor *     &
5185                                              particles(n)%radius**3
5186             ENDDO
5187!
5188!--          Calculate the liquid water content
5189             IF ( ql_v(k,j,i) /= 0.0_wp )  THEN
5190                ql(k,j,i) = ql(k,j,i) + rho_l * 1.33333333_wp * pi *           &
5191                                        ql_v(k,j,i) /                          &
5192                                        ( rho_surface * dx * dy * dzw(k) )
5193                IF ( ql(k,j,i) < 0.0_wp )  THEN
5194                   WRITE( message_string, * )  'LWC out of range: ' , &
5195                                               ql(k,j,i),i,j,k
5196                   CALL message( 'lpm_calc_liquid_water_content', '', 2, 2,    &
5197                                 -1, 6, 1 )
5198                ENDIF
5199             ELSE
5200                ql(k,j,i) = 0.0_wp
5201             ENDIF
5202          ENDDO
5203       ENDDO
5204    ENDDO
5205
5206    CALL cpu_log( log_point_s(45), 'lpm_calc_ql', 'stop' )
5207
5208 END SUBROUTINE lpm_calc_liquid_water_content
5209
5210
5211!------------------------------------------------------------------------------!
5212! Description:
5213! ------------
5214!> Calculates change in droplet radius by collision. Droplet collision is
5215!> calculated for each grid box seperately. Collision is parameterized by
5216!> using collision kernels. Two different kernels are available:
5217!> Hall kernel: Kernel from Hall (1980, J. Atmos. Sci., 2486-2507), which
5218!>              considers collision due to pure gravitational effects.
5219!> Wang kernel: Beside gravitational effects (treated with the Hall-kernel) also
5220!>              the effects of turbulence on the collision are considered using
5221!>              parameterizations of Ayala et al. (2008, New J. Phys., 10,
5222!>              075015) and Wang and Grabowski (2009, Atmos. Sci. Lett., 10,
5223!>              1-8). This kernel includes three possible effects of turbulence:
5224!>              the modification of the relative velocity between the droplets,
5225!>              the effect of preferential concentration, and the enhancement of
5226!>              collision efficiencies.
5227!------------------------------------------------------------------------------!
5228 SUBROUTINE lpm_droplet_collision (i,j,k)
5229
5230    INTEGER(iwp), INTENT(IN) ::  i        !<
5231    INTEGER(iwp), INTENT(IN) ::  j        !<
5232    INTEGER(iwp), INTENT(IN) ::  k        !<
5233
5234    INTEGER(iwp) ::  eclass   !<
5235    INTEGER(iwp) ::  n        !<
5236    INTEGER(iwp) ::  m        !<
5237    INTEGER(iwp) ::  rclass_l !<
5238    INTEGER(iwp) ::  rclass_s !<
5239
5240    REAL(wp) ::  collection_probability  !< probability for collection
5241    REAL(wp) ::  ddV                     !< inverse grid box volume
5242    REAL(wp) ::  epsilon_collision       !< dissipation rate
5243    REAL(wp) ::  factor_volume_to_mass   !< 4.0 / 3.0 * pi * rho_l
5244    REAL(wp) ::  xm                      !< droplet mass of super-droplet m
5245    REAL(wp) ::  xn                      !< droplet mass of super-droplet n
5246    REAL(wp) ::  xsm                     !< aerosol mass of super-droplet m
5247    REAL(wp) ::  xsn                     !< aerosol mass of super-droplet n
5248
5249    REAL(wp), DIMENSION(:), ALLOCATABLE ::  weight    !< weighting factor
5250    REAL(wp), DIMENSION(:), ALLOCATABLE ::  mass      !< total mass of super droplet
5251    REAL(wp), DIMENSION(:), ALLOCATABLE ::  aero_mass !< total aerosol mass of super droplet
5252
5253    CALL cpu_log( log_point_s(43), 'lpm_droplet_coll', 'start' )
5254
5255    number_of_particles   = prt_count(k,j,i)
5256    factor_volume_to_mass = 4.0_wp / 3.0_wp * pi * rho_l
5257    ddV                   = 1.0_wp / ( dx * dy * dzw(k) )
5258!
5259!-- Collision requires at least one super droplet inside the box
5260    IF ( number_of_particles > 0 )  THEN
5261
5262       IF ( use_kernel_tables )  THEN
5263!
5264!--       Fast method with pre-calculated collection kernels for
5265!--       discrete radius- and dissipation-classes.
5266          IF ( wang_kernel )  THEN
5267             eclass = INT( diss(k,j,i) * 1.0E4_wp / 600.0_wp * &
5268                           dissipation_classes ) + 1
5269             epsilon_collision = diss(k,j,i)
5270          ELSE
5271             epsilon_collision = 0.0_wp
5272          ENDIF
5273
5274          IF ( hall_kernel  .OR.  epsilon_collision * 1.0E4_wp < 0.001_wp )  THEN
5275             eclass = 0   ! Hall kernel is used
5276          ELSE
5277             eclass = MIN( dissipation_classes, eclass )
5278          ENDIF
5279
5280       ELSE
5281!
5282!--       Collection kernels are re-calculated for every new
5283!--       grid box. First, allocate memory for kernel table.
5284!--       Third dimension is 1, because table is re-calculated for
5285!--       every new dissipation value.
5286          ALLOCATE( ckernel(1:number_of_particles,1:number_of_particles,1:1) )
5287!
5288!--       Now calculate collection kernel for this box. Note that
5289!--       the kernel is based on the previous time step
5290          CALL recalculate_kernel( i, j, k )
5291
5292       ENDIF
5293!
5294!--    Temporary fields for total mass of super-droplet, aerosol mass, and
5295!--    weighting factor are allocated.
5296       ALLOCATE(mass(1:number_of_particles), weight(1:number_of_particles))
5297       IF ( curvature_solution_effects )  ALLOCATE(aero_mass(1:number_of_particles))
5298
5299       mass(1:number_of_particles)   = particles(1:number_of_particles)%weight_factor * &
5300                                       particles(1:number_of_particles)%radius**3     * &
5301                                       factor_volume_to_mass
5302
5303       weight(1:number_of_particles) = particles(1:number_of_particles)%weight_factor
5304
5305       IF ( curvature_solution_effects )  THEN
5306          aero_mass(1:number_of_particles) = particles(1:number_of_particles)%weight_factor * &
5307                                             particles(1:number_of_particles)%aux1**3       * &
5308                                             4.0_wp / 3.0_wp * pi * rho_s
5309       ENDIF
5310!
5311!--    Calculate collision/coalescence
5312       DO  n = 1, number_of_particles
5313
5314          DO  m = n, number_of_particles
5315!
5316!--          For collisions, the weighting factor of at least one super-droplet
5317!--          needs to be larger or equal to one.
5318             IF ( MIN( weight(n), weight(m) ) < 1.0_wp )  CYCLE
5319!
5320!--          Get mass of individual droplets (aerosols)
5321             xn = mass(n) / weight(n)
5322             xm = mass(m) / weight(m)
5323             IF ( curvature_solution_effects )  THEN
5324                xsn = aero_mass(n) / weight(n)
5325                xsm = aero_mass(m) / weight(m)
5326             ENDIF
5327!
5328!--          Probability that the necessary collisions take place
5329             IF ( use_kernel_tables )  THEN
5330                rclass_l = particles(n)%class
5331                rclass_s = particles(m)%class
5332
5333                collection_probability  = MAX( weight(n), weight(m) ) *     &
5334                                          ckernel(rclass_l,rclass_s,eclass) * ddV * dt_3d
5335             ELSE
5336                collection_probability  = MAX( weight(n), weight(m) ) *     &
5337                                          ckernel(n,m,1) * ddV * dt_3d
5338             ENDIF
5339!
5340!--          Calculate the number of collections and consider multiple collections.
5341!--          (Accordingly, p_crit will be 0.0, 1.0, 2.0, ...)
5342             IF ( collection_probability - FLOOR(collection_probability)    &
5343                  > random_function( iran_part ) )  THEN
5344                collection_probability = FLOOR(collection_probability) + 1.0_wp
5345             ELSE
5346                collection_probability = FLOOR(collection_probability)
5347             ENDIF
5348
5349             IF ( collection_probability > 0.0_wp )  THEN
5350!
5351!--             Super-droplet n collects droplets of super-droplet m
5352                IF ( weight(n) < weight(m) )  THEN
5353
5354                   mass(n)   = mass(n)   + weight(n) * xm * collection_probability
5355                   weight(m) = weight(m) - weight(n)      * collection_probability
5356                   mass(m)   = mass(m)   - weight(n) * xm * collection_probability
5357                   IF ( curvature_solution_effects )  THEN
5358                      aero_mass(n) = aero_mass(n) + weight(n) * xsm * collection_probability
5359                      aero_mass(m) = aero_mass(m) - weight(n) * xsm * collection_probability
5360                   ENDIF
5361
5362                ELSEIF ( weight(m) < weight(n) )  THEN
5363
5364                   mass(m)   = mass(m)   + weight(m) * xn * collection_probability
5365                   weight(n) = weight(n) - weight(m)      * collection_probability
5366                   mass(n)   = mass(n)   - weight(m) * xn * collection_probability
5367                   IF ( curvature_solution_effects )  THEN
5368                      aero_mass(m) = aero_mass(m) + weight(m) * xsn * collection_probability
5369                      aero_mass(n) = aero_mass(n) - weight(m) * xsn * collection_probability
5370                   ENDIF
5371
5372                ELSE
5373!
5374!--                Collisions of particles of the same weighting factor.
5375!--                Particle n collects 1/2 weight(n) droplets of particle m,
5376!--                particle m collects 1/2 weight(m) droplets of particle n.
5377!--                The total mass mass changes accordingly.
5378!--                If n = m, the first half of the droplets coalesces with the
5379!--                second half of the droplets; mass is unchanged because
5380!--                xm = xn for n = m.
5381!--
5382!--                Note: For m = n this equation is an approximation only
5383!--                valid for weight >> 1 (which is usually the case). The
5384!--                approximation is weight(n)-1 = weight(n).
5385                   mass(n)   = mass(n)   + 0.5_wp * weight(n) * ( xm - xn )
5386                   mass(m)   = mass(m)   + 0.5_wp * weight(m) * ( xn - xm )
5387                   IF ( curvature_solution_effects )  THEN
5388                      aero_mass(n) = aero_mass(n) + 0.5_wp * weight(n) * ( xsm - xsn )
5389                      aero_mass(m) = aero_mass(m) + 0.5_wp * weight(m) * ( xsn - xsm )
5390                   ENDIF
5391                   weight(n) = weight(n) - 0.5_wp * weight(m)
5392                   weight(m) = weight(n)
5393
5394                ENDIF
5395
5396             ENDIF
5397
5398          ENDDO
5399
5400          ql_vp(k,j,i) = ql_vp(k,j,i) + mass(n) / factor_volume_to_mass
5401
5402       ENDDO
5403
5404       IF ( ANY(weight < 0.0_wp) )  THEN
5405             WRITE( message_string, * ) 'negative weighting factor'
5406             CALL message( 'lpm_droplet_collision', 'PA0028',      &
5407                            2, 2, -1, 6, 1 )
5408       ENDIF
5409
5410       particles(1:number_of_particles)%radius = ( mass(1:number_of_particles) /   &
5411                                                   ( weight(1:number_of_particles) &
5412                                                     * factor_volume_to_mass       &
5413                                                   )                               &
5414                                                 )**0.33333333333333_wp
5415
5416       IF ( curvature_solution_effects )  THEN
5417          particles(1:number_of_particles)%aux1 = ( aero_mass(1:number_of_particles) / &
5418                                                    ( weight(1:number_of_particles)    &
5419                                                      * 4.0_wp / 3.0_wp * pi * rho_s   &
5420                                                    )                                  &
5421                                                  )**0.33333333333333_wp
5422       ENDIF
5423
5424       particles(1:number_of_particles)%weight_factor = weight(1:number_of_particles)
5425
5426       DEALLOCATE( weight, mass )
5427       IF ( curvature_solution_effects )  DEALLOCATE( aero_mass )
5428       IF ( .NOT. use_kernel_tables )  DEALLOCATE( ckernel )
5429
5430!
5431!--    Check if LWC is conserved during collision process
5432       IF ( ql_v(k,j,i) /= 0.0_wp )  THEN
5433          IF ( ql_vp(k,j,i) / ql_v(k,j,i) >= 1.0001_wp  .OR.                      &
5434               ql_vp(k,j,i) / ql_v(k,j,i) <= 0.9999_wp )  THEN
5435             WRITE( message_string, * ) ' LWC is not conserved during',           &
5436                                        ' collision! ',                           &
5437                                        ' LWC after condensation: ', ql_v(k,j,i), &
5438                                        ' LWC after collision: ', ql_vp(k,j,i)
5439             CALL message( 'lpm_droplet_collision', 'PA0040', 2, 2, -1, 6, 1 )
5440          ENDIF
5441       ENDIF
5442
5443    ENDIF
5444
5445    CALL cpu_log( log_point_s(43), 'lpm_droplet_coll', 'stop' )
5446
5447 END SUBROUTINE lpm_droplet_collision
5448 
5449!------------------------------------------------------------------------------!
5450! Description:
5451! ------------
5452!> Initialization of the collision efficiency matrix with fixed radius and
5453!> dissipation classes, calculated at simulation start only.
5454!------------------------------------------------------------------------------!
5455 SUBROUTINE lpm_init_kernels
5456
5457    INTEGER(iwp) ::  i !<
5458    INTEGER(iwp) ::  j !<
5459    INTEGER(iwp) ::  k !<
5460   
5461!
5462!-- Calculate collision efficiencies for fixed radius- and dissipation
5463!-- classes
5464    IF ( collision_kernel(6:9) == 'fast' )  THEN
5465
5466       ALLOCATE( ckernel(1:radius_classes,1:radius_classes,                 &
5467                 0:dissipation_classes), epsclass(1:dissipation_classes),   &
5468                 radclass(1:radius_classes) )
5469
5470!
5471!--    Calculate the radius class bounds with logarithmic distances
5472!--    in the interval [1.0E-6, 1000.0E-6] m
5473       rclass_lbound = LOG( 1.0E-6_wp )
5474       rclass_ubound = LOG( 1000.0E-6_wp )
5475       radclass(1)   = EXP( rclass_lbound )
5476       DO  i = 2, radius_classes
5477          radclass(i) = EXP( rclass_lbound +                                &
5478                             ( rclass_ubound - rclass_lbound ) *            &
5479                             ( i - 1.0_wp ) / ( radius_classes - 1.0_wp ) )
5480       ENDDO
5481
5482!
5483!--    Set the class bounds for dissipation in interval [0.0, 600.0] cm**2/s**3
5484       DO  i = 1, dissipation_classes
5485          epsclass(i) = 0.06_wp * REAL( i, KIND=wp ) / dissipation_classes
5486       ENDDO
5487!
5488!--    Calculate collision efficiencies of the Wang/ayala kernel
5489       ALLOCATE( ec(1:radius_classes,1:radius_classes),  &
5490                 ecf(1:radius_classes,1:radius_classes), &
5491                 gck(1:radius_classes,1:radius_classes), &
5492                 winf(1:radius_classes) )
5493
5494       DO  k = 1, dissipation_classes
5495
5496          epsilon_collision = epsclass(k)
5497          urms    = 2.02_wp * ( epsilon_collision / 0.04_wp )**( 1.0_wp / 3.0_wp )
5498
5499          CALL turbsd
5500          CALL turb_enhance_eff
5501          CALL effic
5502
5503          DO  j = 1, radius_classes
5504             DO  i = 1, radius_classes
5505                ckernel(i,j,k) = ec(i,j) * gck(i,j) * ecf(i,j)
5506             ENDDO
5507          ENDDO
5508
5509       ENDDO
5510
5511!
5512!--    Calculate collision efficiencies of the Hall kernel
5513       ALLOCATE( hkernel(1:radius_classes,1:radius_classes), &
5514                 hwratio(1:radius_classes,1:radius_classes) )
5515
5516       CALL fallg
5517       CALL effic
5518
5519       DO  j = 1, radius_classes
5520          DO  i =  1, radius_classes
5521             hkernel(i,j) = pi * ( radclass(j) + radclass(i) )**2 &
5522                               * ec(i,j) * ABS( winf(j) - winf(i) )
5523             ckernel(i,j,0) = hkernel(i,j)  ! hall kernel stored on index 0
5524           ENDDO
5525       ENDDO
5526
5527!
5528!--    Test output of efficiencies
5529       IF ( j == -1 )  THEN
5530          PRINT*, '*** Hall kernel'
5531          WRITE ( *,'(5X,20(F4.0,1X))' ) ( radclass(i)*1.0E6_wp, &
5532                                           i = 1,radius_classes )
5533          DO  j = 1, radius_classes
5534             WRITE ( *,'(F4.0,1X,20(F8.4,1X))' ) radclass(j),  &
5535                                       ( hkernel(i,j), i = 1,radius_classes )
5536          ENDDO
5537
5538          DO  k = 1, dissipation_classes
5539             DO  i = 1, radius_classes
5540                DO  j = 1, radius_classes
5541                   IF ( hkernel(i,j) == 0.0_wp )  THEN
5542                      hwratio(i,j) = 9999999.9_wp
5543                   ELSE
5544                      hwratio(i,j) = ckernel(i,j,k) / hkernel(i,j)
5545                   ENDIF
5546                ENDDO
5547             ENDDO
5548
5549             PRINT*, '*** epsilon = ', epsclass(k)
5550             WRITE ( *,'(5X,20(F4.0,1X))' ) ( radclass(i) * 1.0E6_wp, &
5551                                              i = 1,radius_classes )
5552             DO  j = 1, radius_classes
5553                WRITE ( *,'(F4.0,1X,20(F8.4,1X))' ) radclass(j) * 1.0E6_wp, &
5554                                       ( hwratio(i,j), i = 1,radius_classes )
5555             ENDDO
5556          ENDDO
5557       ENDIF
5558
5559       DEALLOCATE( ec, ecf, epsclass, gck, hkernel, winf )
5560
5561    ENDIF
5562
5563 END SUBROUTINE lpm_init_kernels
5564 
5565!------------------------------------------------------------------------------!
5566! Description:
5567! ------------
5568!> Calculation of collision kernels during each timestep and for each grid box
5569!------------------------------------------------------------------------------!
5570 SUBROUTINE recalculate_kernel( i1, j1, k1 )
5571
5572
5573    INTEGER(iwp) ::  i      !<
5574    INTEGER(iwp) ::  i1     !<
5575    INTEGER(iwp) ::  j      !<
5576    INTEGER(iwp) ::  j1     !<
5577    INTEGER(iwp) ::  k1     !<
5578
5579
5580    number_of_particles = prt_count(k1,j1,i1)
5581    radius_classes      = number_of_particles   ! necessary to use the same
5582                                                ! subroutines as for
5583                                                ! precalculated kernels
5584
5585    ALLOCATE( ec(1:number_of_particles,1:number_of_particles), &
5586              radclass(1:number_of_particles), winf(1:number_of_particles) )
5587
5588!
5589!-- Store particle radii on the radclass array
5590    radclass(1:number_of_particles) = particles(1:number_of_particles)%radius
5591
5592    IF ( wang_kernel )  THEN
5593       epsilon_collision = diss(k1,j1,i1)   ! dissipation rate in m**2/s**3
5594    ELSE
5595       epsilon_collision = 0.0_wp
5596    ENDIF
5597    urms    = 2.02_wp * ( epsilon_collision / 0.04_wp )**( 0.33333333333_wp )
5598
5599    IF ( wang_kernel  .AND.  epsilon_collision > 1.0E-7_wp )  THEN
5600!
5601!--    Call routines to calculate efficiencies for the Wang kernel
5602       ALLOCATE( gck(1:number_of_particles,1:number_of_particles), &
5603                 ecf(1:number_of_particles,1:number_of_particles) )
5604
5605       CALL turbsd
5606       CALL turb_enhance_eff
5607       CALL effic
5608
5609       DO  j = 1, number_of_particles
5610          DO  i =  1, number_of_particles
5611             ckernel(1+i-1,1+j-1,1) = ec(i,j) * gck(i,j) * ecf(i,j)
5612          ENDDO
5613       ENDDO
5614
5615       DEALLOCATE( gck, ecf )
5616    ELSE
5617!
5618!--    Call routines to calculate efficiencies for the Hall kernel
5619       CALL fallg
5620       CALL effic
5621
5622       DO  j = 1, number_of_particles
5623          DO  i =  1, number_of_particles
5624             ckernel(i,j,1) = pi * ( radclass(j) + radclass(i) )**2         &
5625                                 * ec(i,j) * ABS( winf(j) - winf(i) )
5626          ENDDO
5627       ENDDO
5628    ENDIF
5629
5630    DEALLOCATE( ec, radclass, winf )
5631
5632 END SUBROUTINE recalculate_kernel
5633
5634!------------------------------------------------------------------------------!
5635! Description:
5636! ------------
5637!> Calculation of effects of turbulence on the geometric collision kernel
5638!> (by including the droplets' average radial relative velocities and their
5639!> radial distribution function) following the analytic model by Aayala et al.
5640!> (2008, New J. Phys.). For details check the second part 2 of the publication,
5641!> page 37ff.
5642!>
5643!> Input parameters, which need to be replaced by PALM parameters:
5644!>    water density, air density
5645!------------------------------------------------------------------------------!
5646 SUBROUTINE turbsd
5647
5648    INTEGER(iwp) ::  i     !<
5649    INTEGER(iwp) ::  j     !<
5650
5651    REAL(wp) ::  ao        !<
5652    REAL(wp) ::  ao_gr     !<
5653    REAL(wp) ::  bbb       !<
5654    REAL(wp) ::  be        !<
5655    REAL(wp) ::  b1        !<
5656    REAL(wp) ::  b2        !<
5657    REAL(wp) ::  ccc       !<
5658    REAL(wp) ::  c1        !<
5659    REAL(wp) ::  c1_gr     !<
5660    REAL(wp) ::  c2        !<
5661    REAL(wp) ::  d1        !<
5662    REAL(wp) ::  d2        !<
5663    REAL(wp) ::  eta       !<
5664    REAL(wp) ::  e1        !<
5665    REAL(wp) ::  e2        !<
5666    REAL(wp) ::  fao_gr    !<
5667    REAL(wp) ::  fr        !<
5668    REAL(wp) ::  grfin     !<
5669    REAL(wp) ::  lambda    !<
5670    REAL(wp) ::  lambda_re !<
5671    REAL(wp) ::  lf        !<
5672    REAL(wp) ::  rc        !<
5673    REAL(wp) ::  rrp       !<
5674    REAL(wp) ::  sst       !<
5675    REAL(wp) ::  tauk      !<
5676    REAL(wp) ::  tl        !<
5677    REAL(wp) ::  t2        !<
5678    REAL(wp) ::  tt        !<
5679    REAL(wp) ::  t1        !<
5680    REAL(wp) ::  vk        !<
5681    REAL(wp) ::  vrms1xy   !<
5682    REAL(wp) ::  vrms2xy   !<
5683    REAL(wp) ::  v1        !<
5684    REAL(wp) ::  v1v2xy    !<
5685    REAL(wp) ::  v1xysq    !<
5686    REAL(wp) ::  v2        !<
5687    REAL(wp) ::  v2xysq    !<
5688    REAL(wp) ::  wrfin     !<
5689    REAL(wp) ::  wrgrav2   !<
5690    REAL(wp) ::  wrtur2xy  !<
5691    REAL(wp) ::  xx        !<
5692    REAL(wp) ::  yy        !<
5693    REAL(wp) ::  z         !<
5694
5695    REAL(wp), DIMENSION(1:radius_classes) ::  st  !< Stokes number
5696    REAL(wp), DIMENSION(1:radius_classes) ::  tau !< inertial time scale
5697
5698    lambda    = urms * SQRT( 15.0_wp * molecular_viscosity / epsilon_collision )
5699    lambda_re = urms**2 * SQRT( 15.0_wp / epsilon_collision / molecular_viscosity )
5700    tl        = urms**2 / epsilon_collision
5701    lf        = 0.5_wp * urms**3 / epsilon_collision
5702    tauk      = SQRT( molecular_viscosity / epsilon_collision )
5703    eta       = ( molecular_viscosity**3 / epsilon_collision )**0.25_wp
5704    vk        = eta / tauk
5705
5706    ao = ( 11.0_wp + 7.0_wp * lambda_re ) / ( 205.0_wp + lambda_re )
5707    tt = SQRT( 2.0_wp * lambda_re / ( SQRT( 15.0_wp ) * ao ) ) * tauk
5708
5709!
5710!-- Get terminal velocity of droplets
5711    CALL fallg
5712
5713    DO  i = 1, radius_classes
5714       tau(i) = winf(i) / g    ! inertial time scale
5715       st(i)  = tau(i) / tauk  ! Stokes number
5716    ENDDO
5717
5718!
5719!-- Calculate average radial relative velocity at contact (wrfin)
5720    z   = tt / tl
5721    be  = SQRT( 2.0_wp ) * lambda / lf
5722    bbb = SQRT( 1.0_wp - 2.0_wp * be**2 )
5723    d1  = ( 1.0_wp + bbb ) / ( 2.0_wp * bbb )
5724    e1  = lf * ( 1.0_wp + bbb ) * 0.5_wp
5725    d2  = ( 1.0_wp - bbb ) * 0.5_wp / bbb
5726    e2  = lf * ( 1.0_wp - bbb ) * 0.5_wp
5727    ccc = SQRT( 1.0_wp - 2.0_wp * z**2 )
5728    b1  = ( 1.0_wp + ccc ) * 0.5_wp / ccc
5729    c1  = tl * ( 1.0_wp + ccc ) * 0.5_wp
5730    b2  = ( 1.0_wp - ccc ) * 0.5_wp / ccc
5731    c2  = tl * ( 1.0_wp - ccc ) * 0.5_wp
5732
5733    DO  i = 1, radius_classes
5734
5735       v1 = winf(i)
5736       t1 = tau(i)
5737
5738       DO  j = 1, i
5739          rrp = radclass(i) + radclass(j)
5740          v2  = winf(j)
5741          t2  = tau(j)
5742
5743          v1xysq  = b1 * d1 * phi_w(c1,e1,v1,t1) - b1 * d2 * phi_w(c1,e2,v1,t1) &
5744                  - b2 * d1 * phi_w(c2,e1,v1,t1) + b2 * d2 * phi_w(c2,e2,v1,t1)
5745          v1xysq  = v1xysq * urms**2 / t1
5746          vrms1xy = SQRT( v1xysq )
5747
5748          v2xysq  = b1 * d1 * phi_w(c1,e1,v2,t2) - b1 * d2 * phi_w(c1,e2,v2,t2) &
5749                  - b2 * d1 * phi_w(c2,e1,v2,t2) + b2 * d2 * phi_w(c2,e2,v2,t2)
5750          v2xysq  = v2xysq * urms**2 / t2
5751          vrms2xy = SQRT( v2xysq )
5752
5753          IF ( winf(i) >= winf(j) )  THEN
5754             v1 = winf(i)
5755             t1 = tau(i)
5756             v2 = winf(j)
5757             t2 = tau(j)
5758          ELSE
5759             v1 = winf(j)
5760             t1 = tau(j)
5761             v2 = winf(i)
5762             t2 = tau(i)
5763          ENDIF
5764
5765          v1v2xy   =  b1 * d1 * zhi(c1,e1,v1,t1,v2,t2) - &
5766                      b1 * d2 * zhi(c1,e2,v1,t1,v2,t2) - &
5767                      b2 * d1 * zhi(c2,e1,v1,t1,v2,t2) + &
5768                      b2 * d2* zhi(c2,e2,v1,t1,v2,t2)
5769          fr       = d1 * EXP( -rrp / e1 ) - d2 * EXP( -rrp / e2 )
5770          v1v2xy   = v1v2xy * fr * urms**2 / tau(i) / tau(j)
5771          wrtur2xy = vrms1xy**2 + vrms2xy**2 - 2.0_wp * v1v2xy
5772          IF ( wrtur2xy < 0.0_wp )  wrtur2xy = 0.0_wp
5773          wrgrav2  = pi / 8.0_wp * ( winf(j) - winf(i) )**2
5774          wrfin    = SQRT( ( 2.0_wp / pi ) * ( wrtur2xy + wrgrav2) )
5775
5776!
5777!--       Calculate radial distribution function (grfin)
5778          IF ( st(j) > st(i) )  THEN
5779             sst = st(j)
5780          ELSE
5781             sst = st(i)
5782          ENDIF
5783
5784          xx = -0.1988_wp * sst**4 + 1.5275_wp * sst**3 - 4.2942_wp *       &
5785                sst**2 + 5.3406_wp * sst
5786          IF ( xx < 0.0_wp )  xx = 0.0_wp
5787          yy = 0.1886_wp * EXP( 20.306_wp / lambda_re )
5788
5789          c1_gr  =  xx / ( g / vk * tauk )**yy
5790
5791          ao_gr  = ao + ( pi / 8.0_wp) * ( g / vk * tauk )**2
5792          fao_gr = 20.115_wp * SQRT( ao_gr / lambda_re )
5793          rc     = SQRT( fao_gr * ABS( st(j) - st(i) ) ) * eta
5794
5795          grfin  = ( ( eta**2 + rc**2 ) / ( rrp**2 + rc**2) )**( c1_gr*0.5_wp )
5796          IF ( grfin < 1.0_wp )  grfin = 1.0_wp
5797
5798!
5799!--       Calculate general collection kernel (without the consideration of
5800!--       collection efficiencies)
5801          gck(i,j) = 2.0_wp * pi * rrp**2 * wrfin * grfin
5802          gck(j,i) = gck(i,j)
5803
5804       ENDDO
5805    ENDDO
5806
5807 END SUBROUTINE turbsd
5808
5809 REAL(wp) FUNCTION phi_w( a, b, vsett, tau0 )
5810!
5811!-- Function used in the Ayala et al. (2008) analytical model for turbulent
5812!-- effects on the collision kernel
5813   
5814
5815    REAL(wp) ::  a     !<
5816    REAL(wp) ::  aa1   !<
5817    REAL(wp) ::  b     !<
5818    REAL(wp) ::  tau0  !<
5819    REAL(wp) ::  vsett !<
5820
5821    aa1 = 1.0_wp / tau0 + 1.0_wp / a + vsett / b
5822    phi_w = 1.0_wp / aa1  - 0.5_wp * vsett / b / aa1**2
5823
5824 END FUNCTION phi_w
5825
5826 REAL(wp) FUNCTION zhi( a, b, vsett1, tau1, vsett2, tau2 )
5827!
5828!-- Function used in the Ayala et al. (2008) analytical model for turbulent
5829!-- effects on the collision kernel
5830
5831    REAL(wp) ::  a      !<
5832    REAL(wp) ::  aa1    !<
5833    REAL(wp) ::  aa2    !<
5834    REAL(wp) ::  aa3    !<
5835    REAL(wp) ::  aa4    !<
5836    REAL(wp) ::  aa5    !<
5837    REAL(wp) ::  aa6    !<
5838    REAL(wp) ::  b      !<
5839    REAL(wp) ::  tau1   !<
5840    REAL(wp) ::  tau2   !<
5841    REAL(wp) ::  vsett1 !<
5842    REAL(wp) ::  vsett2 !<
5843
5844    aa1 = vsett2 / b - 1.0_wp / tau2 - 1.0_wp / a
5845    aa2 = vsett1 / b + 1.0_wp / tau1 + 1.0_wp / a
5846    aa3 = ( vsett1 - vsett2 ) / b + 1.0_wp / tau1 + 1.0_wp / tau2
5847    aa4 = ( vsett2 / b )**2 - ( 1.0_wp / tau2 + 1.0_wp / a )**2
5848    aa5 = vsett2 / b + 1.0_wp / tau2 + 1.0_wp / a
5849    aa6 = 1.0_wp / tau1 - 1.0_wp / a + ( 1.0_wp / tau2 + 1.0_wp / a) *      &
5850          vsett1 / vsett2
5851    zhi = (1.0_wp / aa1 - 1.0_wp / aa2 ) * ( vsett1 - vsett2 ) * 0.5_wp /   &
5852          b / aa3**2 + ( 4.0_wp / aa4 - 1.0_wp / aa5**2 - 1.0_wp / aa1**2 ) &
5853          * vsett2 * 0.5_wp / b /aa6 + ( 2.0_wp * ( b / aa2 - b / aa1 ) -   &
5854          vsett1 / aa2**2 + vsett2 / aa1**2 ) * 0.5_wp / b / aa3
5855
5856 END FUNCTION zhi
5857
5858
5859!------------------------------------------------------------------------------!
5860! Description:
5861! ------------
5862!> Parameterization of terminal velocity following Rogers et al. (1993, J. Appl.
5863!> Meteorol.)
5864!------------------------------------------------------------------------------!
5865 SUBROUTINE fallg
5866
5867    INTEGER(iwp) ::  j                            !<
5868
5869    REAL(wp), PARAMETER ::  k_cap_rog = 4.0_wp    !< parameter
5870    REAL(wp), PARAMETER ::  k_low_rog = 12.0_wp   !< parameter
5871    REAL(wp), PARAMETER ::  a_rog     = 9.65_wp   !< parameter
5872    REAL(wp), PARAMETER ::  b_rog     = 10.43_wp  !< parameter
5873    REAL(wp), PARAMETER ::  c_rog     = 0.6_wp    !< parameter
5874    REAL(wp), PARAMETER ::  d0_rog    = 0.745_wp  !< seperation diameter
5875
5876    REAL(wp)            ::  diameter              !< droplet diameter in mm
5877
5878
5879    DO  j = 1, radius_classes
5880
5881       diameter = radclass(j) * 2000.0_wp
5882
5883       IF ( diameter <= d0_rog )  THEN
5884          winf(j) = k_cap_rog * diameter * ( 1.0_wp -                       &
5885                                             EXP( -k_low_rog * diameter ) )
5886       ELSE
5887          winf(j) = a_rog - b_rog * EXP( -c_rog * diameter )
5888       ENDIF
5889
5890    ENDDO
5891
5892 END SUBROUTINE fallg
5893
5894
5895!------------------------------------------------------------------------------!
5896! Description:
5897! ------------
5898!> Interpolation of collision efficiencies (Hall, 1980, J. Atmos. Sci.)
5899!------------------------------------------------------------------------------!
5900 SUBROUTINE effic
5901 
5902    INTEGER(iwp) ::  i  !<
5903    INTEGER(iwp) ::  iq !<
5904    INTEGER(iwp) ::  ir !<
5905    INTEGER(iwp) ::  j  !<
5906    INTEGER(iwp) ::  k  !<
5907
5908    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  ira !<
5909
5910    LOGICAL, SAVE ::  first = .TRUE. !<
5911
5912    REAL(wp) ::  ek              !<
5913    REAL(wp) ::  particle_radius !<
5914    REAL(wp) ::  pp              !<
5915    REAL(wp) ::  qq              !<
5916    REAL(wp) ::  rq              !<
5917
5918    REAL(wp), DIMENSION(1:21), SAVE ::  rat        !<
5919   
5920    REAL(wp), DIMENSION(1:15), SAVE ::  r0         !<
5921   
5922    REAL(wp), DIMENSION(1:15,1:21), SAVE ::  ecoll !<
5923
5924!
5925!-- Initial assignment of constants
5926    IF ( first )  THEN
5927
5928      first = .FALSE.
5929      r0  = (/   6.0_wp,   8.0_wp,  10.0_wp, 15.0_wp,  20.0_wp,  25.0_wp,   &
5930                30.0_wp,  40.0_wp,  50.0_wp, 60.0_wp,  70.0_wp, 100.0_wp,   &
5931               150.0_wp, 200.0_wp, 300.0_wp /)
5932
5933      rat = (/ 0.00_wp, 0.05_wp, 0.10_wp, 0.15_wp, 0.20_wp, 0.25_wp,        &
5934               0.30_wp, 0.35_wp, 0.40_wp, 0.45_wp, 0.50_wp, 0.55_wp,        &
5935               0.60_wp, 0.65_wp, 0.70_wp, 0.75_wp, 0.80_wp, 0.85_wp,        &
5936               0.90_wp, 0.95_wp, 1.00_wp /)
5937
5938      ecoll(:,1)  = (/ 0.001_wp, 0.001_wp, 0.001_wp, 0.001_wp, 0.001_wp,    &
5939                       0.001_wp, 0.001_wp, 0.001_wp, 0.001_wp, 0.001_wp,    &
5940                       0.001_wp, 0.001_wp, 0.001_wp, 0.001_wp, 0.001_wp /)
5941      ecoll(:,2)  = (/ 0.003_wp, 0.003_wp, 0.003_wp, 0.004_wp, 0.005_wp,    &
5942                       0.005_wp, 0.005_wp, 0.010_wp, 0.100_wp, 0.050_wp,    &
5943                       0.200_wp, 0.500_wp, 0.770_wp, 0.870_wp, 0.970_wp /)
5944      ecoll(:,3)  = (/ 0.007_wp, 0.007_wp, 0.007_wp, 0.008_wp, 0.009_wp,    &
5945                       0.010_wp, 0.010_wp, 0.070_wp, 0.400_wp, 0.430_wp,    &
5946                       0.580_wp, 0.790_wp, 0.930_wp, 0.960_wp, 1.000_wp /)
5947      ecoll(:,4)  = (/ 0.009_wp, 0.009_wp, 0.009_wp, 0.012_wp, 0.015_wp,    &
5948                       0.010_wp, 0.020_wp, 0.280_wp, 0.600_wp, 0.640_wp,    &
5949                       0.750_wp, 0.910_wp, 0.970_wp, 0.980_wp, 1.000_wp /)
5950      ecoll(:,5)  = (/ 0.014_wp, 0.014_wp, 0.014_wp, 0.015_wp, 0.016_wp,    &
5951                       0.030_wp, 0.060_wp, 0.500_wp, 0.700_wp, 0.770_wp,    &
5952                       0.840_wp, 0.950_wp, 0.970_wp, 1.000_wp, 1.000_wp /)
5953      ecoll(:,6)  = (/ 0.017_wp, 0.017_wp, 0.017_wp, 0.020_wp, 0.022_wp,    &
5954                       0.060_wp, 0.100_wp, 0.620_wp, 0.780_wp, 0.840_wp,    &
5955                       0.880_wp, 0.950_wp, 1.000_wp, 1.000_wp, 1.000_wp /)
5956      ecoll(:,7)  = (/ 0.030_wp, 0.030_wp, 0.024_wp, 0.022_wp, 0.032_wp,    &
5957                       0.062_wp, 0.200_wp, 0.680_wp, 0.830_wp, 0.870_wp,    &
5958                       0.900_wp, 0.950_wp, 1.000_wp, 1.000_wp, 1.000_wp /)
5959      ecoll(:,8)  = (/ 0.025_wp, 0.025_wp, 0.025_wp, 0.036_wp, 0.043_wp,    &
5960                       0.130_wp, 0.270_wp, 0.740_wp, 0.860_wp, 0.890_wp,    &
5961                       0.920_wp, 1.000_wp, 1.000_wp, 1.000_wp, 1.000_wp /)
5962      ecoll(:,9)  = (/ 0.027_wp, 0.027_wp, 0.027_wp, 0.040_wp, 0.052_wp,    &
5963                       0.200_wp, 0.400_wp, 0.780_wp, 0.880_wp, 0.900_wp,    &
5964                       0.940_wp, 1.000_wp, 1.000_wp, 1.000_wp, 1.000_wp /)
5965      ecoll(:,10) = (/ 0.030_wp, 0.030_wp, 0.030_wp, 0.047_wp, 0.064_wp,    &
5966                       0.250_wp, 0.500_wp, 0.800_wp, 0.900_wp, 0.910_wp,    &
5967                       0.950_wp, 1.000_wp, 1.000_wp, 1.000_wp, 1.000_wp /)
5968      ecoll(:,11) = (/ 0.040_wp, 0.040_wp, 0.033_wp, 0.037_wp, 0.068_wp,    &
5969                       0.240_wp, 0.550_wp, 0.800_wp, 0.900_wp, 0.910_wp,    &
5970                       0.950_wp, 1.000_wp, 1.000_wp, 1.000_wp, 1.000_wp /)
5971      ecoll(:,12) = (/ 0.035_wp, 0.035_wp, 0.035_wp, 0.055_wp, 0.079_wp,    &
5972                       0.290_wp, 0.580_wp, 0.800_wp, 0.900_wp, 0.910_wp,    &
5973                       0.950_wp, 1.000_wp, 1.000_wp, 1.000_wp, 1.000_wp /)
5974      ecoll(:,13) = (/ 0.037_wp, 0.037_wp, 0.037_wp, 0.062_wp, 0.082_wp,    &
5975                       0.290_wp, 0.590_wp, 0.780_wp, 0.900_wp, 0.910_wp,    &
5976                       0.950_wp, 1.000_wp, 1.000_wp, 1.000_wp, 1.000_wp /)
5977      ecoll(:,14) = (/ 0.037_wp, 0.037_wp, 0.037_wp, 0.060_wp, 0.080_wp,    &
5978                       0.290_wp, 0.580_wp, 0.770_wp, 0.890_wp, 0.910_wp,    &
5979                       0.950_wp, 1.000_wp, 1.000_wp, 1.000_wp, 1.000_wp /)
5980      ecoll(:,15) = (/ 0.037_wp, 0.037_wp, 0.037_wp, 0.041_wp, 0.075_wp,    &
5981                       0.250_wp, 0.540_wp, 0.760_wp, 0.880_wp, 0.920_wp,    &
5982                       0.950_wp, 1.000_wp, 1.000_wp, 1.000_wp, 1.000_wp /)
5983      ecoll(:,16) = (/ 0.037_wp, 0.037_wp, 0.037_wp, 0.052_wp, 0.067_wp,    &
5984                       0.250_wp, 0.510_wp, 0.770_wp, 0.880_wp, 0.930_wp,    &
5985                       0.970_wp, 1.000_wp, 1.000_wp, 1.000_wp, 1.000_wp /)
5986      ecoll(:,17) = (/ 0.037_wp, 0.037_wp, 0.037_wp, 0.047_wp, 0.057_wp,    &
5987                       0.250_wp, 0.490_wp, 0.770_wp, 0.890_wp, 0.950_wp,    &
5988                       1.000_wp, 1.000_wp, 1.000_wp, 1.000_wp, 1.000_wp /)
5989      ecoll(:,18) = (/ 0.036_wp, 0.036_wp, 0.036_wp, 0.042_wp, 0.048_wp,    &
5990                       0.230_wp, 0.470_wp, 0.780_wp, 0.920_wp, 1.000_wp,    &
5991                       1.020_wp, 1.020_wp, 1.020_wp, 1.020_wp, 1.020_wp /)
5992      ecoll(:,19) = (/ 0.040_wp, 0.040_wp, 0.035_wp, 0.033_wp, 0.040_wp,    &
5993                       0.112_wp, 0.450_wp, 0.790_wp, 1.010_wp, 1.030_wp,    &
5994                       1.040_wp, 1.040_wp, 1.040_wp, 1.040_wp, 1.040_wp /)
5995      ecoll(:,20) = (/ 0.033_wp, 0.033_wp, 0.033_wp, 0.033_wp, 0.033_wp,    &
5996                       0.119_wp, 0.470_wp, 0.950_wp, 1.300_wp, 1.700_wp,    &
5997                       2.300_wp, 2.300_wp, 2.300_wp, 2.300_wp, 2.300_wp /)
5998      ecoll(:,21) = (/ 0.027_wp, 0.027_wp, 0.027_wp, 0.027_wp, 0.027_wp,    &
5999                       0.125_wp, 0.520_wp, 1.400_wp, 2.300_wp, 3.000_wp,    &
6000                       4.000_wp, 4.000_wp, 4.000_wp, 4.000_wp, 4.000_wp /)
6001    ENDIF
6002
6003!
6004!-- Calculate the radius class index of particles with respect to array r
6005!-- Radius has to be in microns
6006    ALLOCATE( ira(1:radius_classes) )
6007    DO  j = 1, radius_classes
6008       particle_radius = radclass(j) * 1.0E6_wp
6009       DO  k = 1, 15
6010          IF ( particle_radius < r0(k) )  THEN
6011             ira(j) = k
6012             EXIT
6013          ENDIF
6014       ENDDO
6015       IF ( particle_radius >= r0(15) )  ira(j) = 16
6016    ENDDO
6017
6018!
6019!-- Two-dimensional linear interpolation of the collision efficiency.
6020!-- Radius has to be in microns
6021    DO  j = 1, radius_classes
6022       DO  i = 1, j
6023
6024          ir = MAX( ira(i), ira(j) )
6025          rq = MIN( radclass(i) / radclass(j), radclass(j) / radclass(i) )
6026          iq = INT( rq * 20 ) + 1
6027          iq = MAX( iq , 2)
6028
6029          IF ( ir < 16 )  THEN
6030             IF ( ir >= 2 )  THEN
6031                pp = ( ( MAX( radclass(j), radclass(i) ) * 1.0E6_wp ) -     &
6032                       r0(ir-1) ) / ( r0(ir) - r0(ir-1) )
6033                qq = ( rq - rat(iq-1) ) / ( rat(iq) - rat(iq-1) )
6034                ec(j,i) = ( 1.0_wp - pp ) * ( 1.0_wp - qq )                 &
6035                          * ecoll(ir-1,iq-1)                                &
6036                          + pp * ( 1.0_wp - qq ) * ecoll(ir,iq-1)           &
6037                          + qq * ( 1.0_wp - pp ) * ecoll(ir-1,iq)           &
6038                          + pp * qq * ecoll(ir,iq)
6039             ELSE
6040                qq = ( rq - rat(iq-1) ) / ( rat(iq) - rat(iq-1) )
6041                ec(j,i) = ( 1.0_wp - qq ) * ecoll(1,iq-1) + qq * ecoll(1,iq)
6042             ENDIF
6043          ELSE
6044             qq = ( rq - rat(iq-1) ) / ( rat(iq) - rat(iq-1) )
6045             ek = ( 1.0_wp - qq ) * ecoll(15,iq-1) + qq * ecoll(15,iq)
6046             ec(j,i) = MIN( ek, 1.0_wp )
6047          ENDIF
6048
6049          IF ( ec(j,i) < 1.0E-20_wp )  ec(j,i) = 0.0_wp
6050
6051          ec(i,j) = ec(j,i)
6052
6053       ENDDO
6054    ENDDO
6055
6056    DEALLOCATE( ira )
6057
6058 END SUBROUTINE effic
6059
6060
6061!------------------------------------------------------------------------------!
6062! Description:
6063! ------------
6064!> Interpolation of turbulent enhancement factor for collision efficencies
6065!> following Wang and Grabowski (2009, Atmos. Sci. Let.)
6066!------------------------------------------------------------------------------!
6067 SUBROUTINE turb_enhance_eff
6068
6069    INTEGER(iwp) ::  i  !<
6070    INTEGER(iwp) ::  iq !<
6071    INTEGER(iwp) ::  ir !<
6072    INTEGER(iwp) ::  j  !<
6073    INTEGER(iwp) ::  k  !<
6074    INTEGER(iwp) ::  kk !<
6075
6076    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  ira !<
6077
6078    LOGICAL, SAVE ::  first = .TRUE. !<
6079
6080    REAL(wp) ::  particle_radius !<
6081    REAL(wp) ::  pp              !<
6082    REAL(wp) ::  qq              !<
6083    REAL(wp) ::  rq              !<
6084    REAL(wp) ::  y1              !<
6085    REAL(wp) ::  y2              !<
6086    REAL(wp) ::  y3              !<
6087
6088    REAL(wp), DIMENSION(1:11), SAVE ::  rat           !<
6089    REAL(wp), DIMENSION(1:7), SAVE  ::  r0            !<
6090
6091    REAL(wp), DIMENSION(1:7,1:11), SAVE ::  ecoll_100 !<
6092    REAL(wp), DIMENSION(1:7,1:11), SAVE ::  ecoll_400 !<
6093
6094!
6095!-- Initial assignment of constants
6096    IF ( first )  THEN
6097
6098       first = .FALSE.
6099
6100       r0  = (/  10.0_wp, 20.0_wp, 30.0_wp, 40.0_wp, 50.0_wp, 60.0_wp,  &
6101                100.0_wp /)
6102
6103       rat = (/ 0.0_wp, 0.1_wp, 0.2_wp, 0.3_wp, 0.4_wp, 0.5_wp, 0.6_wp, &
6104                0.7_wp, 0.8_wp, 0.9_wp, 1.0_wp /)
6105!
6106!--    Tabulated turbulent enhancement factor at 100 cm**2/s**3
6107       ecoll_100(:,1)  = (/  1.74_wp,   1.74_wp,   1.773_wp, 1.49_wp,  &
6108                             1.207_wp,  1.207_wp,  1.0_wp /)
6109       ecoll_100(:,2)  = (/  1.46_wp,   1.46_wp,   1.421_wp, 1.245_wp, &
6110                             1.069_wp,  1.069_wp,  1.0_wp /)
6111       ecoll_100(:,3)  = (/  1.32_wp,   1.32_wp,   1.245_wp, 1.123_wp, &
6112                             1.000_wp,  1.000_wp,  1.0_wp /)
6113       ecoll_100(:,4)  = (/  1.250_wp,  1.250_wp,  1.148_wp, 1.087_wp, &
6114                             1.025_wp,  1.025_wp,  1.0_wp /)
6115       ecoll_100(:,5)  = (/  1.186_wp,  1.186_wp,  1.066_wp, 1.060_wp, &
6116                             1.056_wp,  1.056_wp,  1.0_wp /)
6117       ecoll_100(:,6)  = (/  1.045_wp,  1.045_wp,  1.000_wp, 1.014_wp, &
6118                             1.028_wp,  1.028_wp,  1.0_wp /)
6119       ecoll_100(:,7)  = (/  1.070_wp,  1.070_wp,  1.030_wp, 1.038_wp, &
6120                             1.046_wp,  1.046_wp,  1.0_wp /)
6121       ecoll_100(:,8)  = (/  1.000_wp,  1.000_wp,  1.054_wp, 1.042_wp, &
6122                             1.029_wp,  1.029_wp,  1.0_wp /)
6123       ecoll_100(:,9)  = (/  1.223_wp,  1.223_wp,  1.117_wp, 1.069_wp, &
6124                             1.021_wp,  1.021_wp,  1.0_wp /)
6125       ecoll_100(:,10) = (/  1.570_wp,  1.570_wp,  1.244_wp, 1.166_wp, &
6126                             1.088_wp,  1.088_wp,  1.0_wp /)
6127       ecoll_100(:,11) = (/ 20.3_wp,   20.3_wp,   14.6_wp,   8.61_wp,  &
6128                             2.60_wp,   2.60_wp,   1.0_wp /)
6129!
6130!--    Tabulated turbulent enhancement factor at 400 cm**2/s**3
6131       ecoll_400(:,1)  = (/  4.976_wp,  4.976_wp,  3.593_wp,  2.519_wp, &
6132                             1.445_wp,  1.445_wp,  1.0_wp /)
6133       ecoll_400(:,2)  = (/  2.984_wp,  2.984_wp,  2.181_wp,  1.691_wp, &
6134                             1.201_wp,  1.201_wp,  1.0_wp /)
6135       ecoll_400(:,3)  = (/  1.988_wp,  1.988_wp,  1.475_wp,  1.313_wp, &
6136                             1.150_wp,  1.150_wp,  1.0_wp /)
6137       ecoll_400(:,4)  = (/  1.490_wp,  1.490_wp,  1.187_wp,  1.156_wp, &
6138                             1.126_wp,  1.126_wp,  1.0_wp /)
6139       ecoll_400(:,5)  = (/  1.249_wp,  1.249_wp,  1.088_wp,  1.090_wp, &
6140                             1.092_wp,  1.092_wp,  1.0_wp /)
6141       ecoll_400(:,6)  = (/  1.139_wp,  1.139_wp,  1.130_wp,  1.091_wp, &
6142                             1.051_wp,  1.051_wp,  1.0_wp /)
6143       ecoll_400(:,7)  = (/  1.220_wp,  1.220_wp,  1.190_wp,  1.138_wp, &
6144                             1.086_wp,  1.086_wp,  1.0_wp /)
6145       ecoll_400(:,8)  = (/  1.325_wp,  1.325_wp,  1.267_wp,  1.165_wp, &
6146                             1.063_wp,  1.063_wp,  1.0_wp /)
6147       ecoll_400(:,9)  = (/  1.716_wp,  1.716_wp,  1.345_wp,  1.223_wp, &
6148                             1.100_wp,  1.100_wp,  1.0_wp /)
6149       ecoll_400(:,10) = (/  3.788_wp,  3.788_wp,  1.501_wp,  1.311_wp, &
6150                             1.120_wp,  1.120_wp,  1.0_wp /)
6151       ecoll_400(:,11) = (/ 36.52_wp,  36.52_wp,  19.16_wp,  22.80_wp,  &
6152                            26.0_wp,   26.0_wp,    1.0_wp /)
6153
6154    ENDIF
6155
6156!
6157!-- Calculate the radius class index of particles with respect to array r0
6158!-- The droplet radius has to be given in microns.
6159    ALLOCATE( ira(1:radius_classes) )
6160
6161    DO  j = 1, radius_classes
6162       particle_radius = radclass(j) * 1.0E6_wp
6163       DO  k = 1, 7
6164          IF ( particle_radius < r0(k) )  THEN
6165             ira(j) = k
6166             EXIT
6167          ENDIF
6168       ENDDO
6169       IF ( particle_radius >= r0(7) )  ira(j) = 8
6170    ENDDO
6171
6172!
6173!-- Two-dimensional linear interpolation of the turbulent enhancement factor.
6174!-- The droplet radius has to be given in microns.
6175    DO  j =  1, radius_classes
6176       DO  i = 1, j
6177
6178          ir = MAX( ira(i), ira(j) )
6179          rq = MIN( radclass(i) / radclass(j), radclass(j) / radclass(i) )
6180
6181          DO  kk = 2, 11
6182             IF ( rq <= rat(kk) )  THEN
6183                iq = kk
6184                EXIT
6185             ENDIF
6186          ENDDO
6187
6188          y1 = 1.0_wp  ! turbulent enhancement factor at 0 m**2/s**3
6189
6190          IF ( ir < 8 )  THEN
6191             IF ( ir >= 2 )  THEN
6192                pp = ( MAX( radclass(j), radclass(i) ) * 1.0E6_wp -  &
6193                       r0(ir-1) ) / ( r0(ir) - r0(ir-1) )
6194                qq = ( rq - rat(iq-1) ) / ( rat(iq) - rat(iq-1) )
6195                y2 = ( 1.0_wp - pp ) * ( 1.0_wp - qq ) * ecoll_100(ir-1,iq-1) + &
6196                             pp * ( 1.0_wp - qq ) * ecoll_100(ir,iq-1)        + &
6197                             qq * ( 1.0_wp - pp ) * ecoll_100(ir-1,iq)        + &
6198                             pp * qq              * ecoll_100(ir,iq)
6199                y3 = ( 1.0-pp ) * ( 1.0_wp - qq ) * ecoll_400(ir-1,iq-1)      + &
6200                             pp * ( 1.0_wp - qq ) * ecoll_400(ir,iq-1)        + &
6201                             qq * ( 1.0_wp - pp ) * ecoll_400(ir-1,iq)        + &
6202                             pp * qq              * ecoll_400(ir,iq)
6203             ELSE
6204                qq = ( rq - rat(iq-1) ) / ( rat(iq) - rat(iq-1) )
6205                y2 = ( 1.0_wp - qq ) * ecoll_100(1,iq-1) + qq * ecoll_100(1,iq)
6206                y3 = ( 1.0_wp - qq ) * ecoll_400(1,iq-1) + qq * ecoll_400(1,iq)
6207             ENDIF
6208          ELSE
6209             qq = ( rq - rat(iq-1) ) / ( rat(iq) - rat(iq-1) )
6210             y2 = ( 1.0_wp - qq ) * ecoll_100(7,iq-1) + qq * ecoll_100(7,iq)
6211             y3 = ( 1.0_wp - qq ) * ecoll_400(7,iq-1) + qq * ecoll_400(7,iq)
6212          ENDIF
6213!
6214!--       Linear interpolation of turbulent enhancement factor
6215          IF ( epsilon_collision <= 0.01_wp )  THEN
6216             ecf(j,i) = ( epsilon_collision - 0.01_wp ) / ( 0.0_wp  - 0.01_wp ) * y1 &
6217                      + ( epsilon_collision - 0.0_wp  ) / ( 0.01_wp - 0.0_wp  ) * y2
6218          ELSEIF ( epsilon_collision <= 0.06_wp )  THEN
6219             ecf(j,i) = ( epsilon_collision - 0.04_wp ) / ( 0.01_wp - 0.04_wp ) * y2 &
6220                      + ( epsilon_collision - 0.01_wp ) / ( 0.04_wp - 0.01_wp ) * y3
6221          ELSE
6222             ecf(j,i) = ( 0.06_wp - 0.04_wp ) / ( 0.01_wp - 0.04_wp ) * y2 &
6223                      + ( 0.06_wp - 0.01_wp ) / ( 0.04_wp - 0.01_wp ) * y3
6224          ENDIF
6225
6226          IF ( ecf(j,i) < 1.0_wp )  ecf(j,i) = 1.0_wp
6227
6228          ecf(i,j) = ecf(j,i)
6229
6230       ENDDO
6231    ENDDO
6232
6233 END SUBROUTINE turb_enhance_eff
6234 
6235 
6236 !------------------------------------------------------------------------------!
6237! Description:
6238! ------------
6239! This routine is a part of the Lagrangian particle model. Super droplets which
6240! fulfill certain criterion's (e.g. a big weighting factor and a large radius)
6241! can be split into several super droplets with a reduced number of
6242! represented particles of every super droplet. This mechanism ensures an
6243! improved representation of the right tail of the drop size distribution with
6244! a feasible amount of computational costs. The limits of particle creation
6245! should be chosen carefully! The idea of this algorithm is based on
6246! Unterstrasser and Soelch, 2014.
6247!------------------------------------------------------------------------------!
6248 SUBROUTINE lpm_splitting
6249
6250    INTEGER(iwp) ::  i                !<
6251    INTEGER(iwp) ::  j                !<
6252    INTEGER(iwp) ::  jpp              !<
6253    INTEGER(iwp) ::  k                !<
6254    INTEGER(iwp) ::  n                !<
6255    INTEGER(iwp) ::  new_particles_gb !< counter of created particles within one grid box
6256    INTEGER(iwp) ::  new_size         !< new particle array size
6257    INTEGER(iwp) ::  np               !<
6258    INTEGER(iwp) ::  old_size         !< old particle array size
6259
6260    INTEGER(iwp), PARAMETER ::  n_max = 100 !< number of radii bin for splitting functions   
6261
6262    LOGICAL ::  first_loop_stride_sp = .TRUE. !< flag to calculate constants only once
6263
6264    REAL(wp) ::  diameter                 !< diameter of droplet
6265    REAL(wp) ::  dlog                     !< factor for DSD calculation
6266    REAL(wp) ::  factor_volume_to_mass    !< pre calculate factor volume to mass
6267    REAL(wp) ::  lambda                   !< slope parameter of gamma-distribution
6268    REAL(wp) ::  lwc                      !< liquid water content of grid box
6269    REAL(wp) ::  lwc_total                !< average liquid water content of cloud
6270    REAL(wp) ::  m1                       !< first moment of DSD
6271    REAL(wp) ::  m1_total                 !< average over all PEs of first moment of DSD
6272    REAL(wp) ::  m2                       !< second moment of DSD
6273    REAL(wp) ::  m2_total                 !< average average over all PEs second moment of DSD
6274    REAL(wp) ::  m3                       !< third moment of DSD
6275    REAL(wp) ::  m3_total                 !< average average over all PEs third moment of DSD
6276    REAL(wp) ::  mu                       !< spectral shape parameter of gamma distribution
6277    REAL(wp) ::  nrclgb                   !< number of cloudy grid boxes (ql >= 1.0E-5 kg/kg)
6278    REAL(wp) ::  nrclgb_total             !< average over all PEs of number of cloudy grid boxes
6279    REAL(wp) ::  nr                       !< number concentration of cloud droplets
6280    REAL(wp) ::  nr_total                 !< average over all PEs of number of cloudy grid boxes
6281    REAL(wp) ::  nr0                      !< intercept parameter of gamma distribution
6282    REAL(wp) ::  pirho_l                  !< pi * rho_l / 6.0
6283    REAL(wp) ::  ql_crit = 1.0E-5_wp      !< threshold lwc for cloudy grid cells
6284                                          !< (Siebesma et al 2003, JAS, 60)
6285    REAL(wp) ::  rm                       !< volume averaged mean radius
6286    REAL(wp) ::  rm_total                 !< average over all PEs of volume averaged mean radius
6287    REAL(wp) ::  r_min = 1.0E-6_wp        !< minimum radius of approximated spectra
6288    REAL(wp) ::  r_max = 1.0E-3_wp        !< maximum radius of approximated spectra
6289    REAL(wp) ::  sigma_log = 1.5_wp       !< standard deviation of the LOG-distribution
6290    REAL(wp) ::  zeta                     !< Parameter for DSD calculation of Seifert
6291
6292    REAL(wp), DIMENSION(0:n_max-1) ::  an_spl     !< size dependent critical weight factor
6293    REAL(wp), DIMENSION(0:n_max-1) ::  r_bin_mid  !< mass weighted mean radius of a bin
6294    REAL(wp), DIMENSION(0:n_max)   ::  r_bin      !< boundaries of a radius bin
6295
6296    TYPE(particle_type) ::  tmp_particle   !< temporary particle TYPE
6297
6298    CALL cpu_log( log_point_s(80), 'lpm_splitting', 'start' )
6299
6300    IF ( first_loop_stride_sp )  THEN
6301       IF ( i_splitting_mode == 2  .OR.  i_splitting_mode == 3 )  THEN
6302          dlog   = ( LOG10(r_max) - LOG10(r_min) ) / ( n_max - 1 )
6303          DO  i = 0, n_max-1
6304             r_bin(i) = 10.0_wp**( LOG10(r_min) + i * dlog - 0.5_wp * dlog )
6305             r_bin_mid(i) = 10.0_wp**( LOG10(r_min) + i * dlog )
6306          ENDDO
6307          r_bin(n_max) = 10.0_wp**( LOG10(r_min) + n_max * dlog - 0.5_wp * dlog )
6308       ENDIF   
6309       factor_volume_to_mass =  4.0_wp / 3.0_wp * pi * rho_l
6310       pirho_l  = pi * rho_l / 6.0_wp
6311       IF ( weight_factor_split == -1.0_wp )  THEN
6312          weight_factor_split = 0.1_wp * initial_weighting_factor 
6313       ENDIF
6314    ENDIF
6315
6316
6317    IF ( i_splitting_mode == 1 )  THEN
6318
6319       DO  i = nxl, nxr
6320          DO  j = nys, nyn
6321             DO  k = nzb+1, nzt
6322
6323                new_particles_gb = 0
6324                number_of_particles = prt_count(k,j,i)
6325                IF ( number_of_particles <= 0  .OR.                            & 
6326                     ql(k,j,i) < ql_crit )  CYCLE
6327                particles => grid_particles(k,j,i)%particles(1:number_of_particles)
6328!
6329!--             Start splitting operations. Each particle is checked if it
6330!--             fulfilled the splitting criterion's. In splitting mode 'const'   
6331!--             a critical radius  (radius_split) a critical weighting factor
6332!--             (weight_factor_split) and a splitting factor (splitting_factor)
6333!--             must  be prescribed (see particle_parameters). Super droplets
6334!--             which have a larger radius and larger weighting factor are split
6335!--             into 'splitting_factor' super droplets. Therefore, the weighting
6336!--             factor of  the super droplet and all created clones is reduced
6337!--             by the factor of 'splitting_factor'.
6338                DO  n = 1, number_of_particles
6339                   IF ( particles(n)%particle_mask  .AND.                      &
6340                        particles(n)%radius >= radius_split  .AND.             & 
6341                        particles(n)%weight_factor >= weight_factor_split )    &
6342                   THEN
6343!
6344!--                   Calculate the new number of particles.
6345                      new_size = prt_count(k,j,i) + splitting_factor - 1
6346!
6347!--                   Cycle if maximum number of particles per grid box
6348!--                   is greater than the allowed maximum number.
6349                      IF ( new_size >= max_number_particles_per_gridbox )  CYCLE
6350!
6351!--                   Reallocate particle array if necessary.
6352                      IF ( new_size > SIZE(particles) )  THEN
6353                         CALL realloc_particles_array( i, j, k, new_size )
6354                      ENDIF
6355                      old_size = prt_count(k,j,i)
6356!
6357!--                   Calculate new weighting factor.
6358                      particles(n)%weight_factor =  & 
6359                         particles(n)%weight_factor / splitting_factor
6360                      tmp_particle = particles(n)
6361!
6362!--                   Create splitting_factor-1 new particles.
6363                      DO  jpp = 1, splitting_factor-1
6364                         grid_particles(k,j,i)%particles(jpp+old_size) =       & 
6365                            tmp_particle
6366                      ENDDO 
6367                      new_particles_gb = new_particles_gb + splitting_factor - 1
6368!   
6369!--                   Save the new number of super droplets for every grid box.
6370                      prt_count(k,j,i) = prt_count(k,j,i) +                    &
6371                                         splitting_factor - 1
6372                   ENDIF
6373                ENDDO
6374
6375             ENDDO
6376          ENDDO
6377       ENDDO
6378
6379    ELSEIF ( i_splitting_mode == 2 )  THEN 
6380!
6381!--    Initialize summing variables.
6382       lwc          = 0.0_wp
6383       lwc_total    = 0.0_wp 
6384       m1           = 0.0_wp
6385       m1_total     = 0.0_wp
6386       m2           = 0.0_wp
6387       m2_total     = 0.0_wp
6388       m3           = 0.0_wp
6389       m3_total     = 0.0_wp
6390       nr           = 0.0_wp
6391       nrclgb       = 0.0_wp
6392       nrclgb_total = 0.0_wp
6393       nr_total     = 0.0_wp
6394       rm           = 0.0_wp
6395       rm_total     = 0.0_wp
6396
6397       DO  i = nxl, nxr
6398          DO  j = nys, nyn
6399             DO  k = nzb+1, nzt
6400                number_of_particles = prt_count(k,j,i)
6401                IF ( number_of_particles <= 0  .OR.                            & 
6402                     ql(k,j,i) < ql_crit )  CYCLE
6403                particles => grid_particles(k,j,i)%particles(1:number_of_particles)
6404                nrclgb = nrclgb + 1.0_wp
6405!
6406!--             Calculate moments of DSD.
6407                DO  n = 1, number_of_particles
6408                   IF ( particles(n)%particle_mask  .AND.                      &
6409                        particles(n)%radius >= r_min )                         &
6410                   THEN
6411                      nr  = nr  + particles(n)%weight_factor
6412                      rm  = rm  + factor_volume_to_mass  *                     &
6413                                 particles(n)%radius**3  *                     &
6414                                 particles(n)%weight_factor
6415                      IF ( isf == 1 )  THEN           
6416                         diameter   = particles(n)%radius * 2.0_wp
6417                         lwc = lwc + factor_volume_to_mass *                   &
6418                                     particles(n)%radius**3 *                  & 
6419                                     particles(n)%weight_factor 
6420                         m1  = m1  + particles(n)%weight_factor * diameter
6421                         m2  = m2  + particles(n)%weight_factor * diameter**2
6422                         m3  = m3  + particles(n)%weight_factor * diameter**3
6423                      ENDIF
6424                   ENDIF
6425                ENDDO 
6426             ENDDO
6427          ENDDO
6428       ENDDO
6429
6430#if defined( __parallel )
6431       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
6432       CALL MPI_ALLREDUCE( nr, nr_total, 1 , &
6433       MPI_REAL, MPI_SUM, comm2d, ierr )
6434       CALL MPI_ALLREDUCE( rm, rm_total, 1 , &
6435       MPI_REAL, MPI_SUM, comm2d, ierr )
6436       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
6437       CALL MPI_ALLREDUCE( nrclgb, nrclgb_total, 1 , &
6438       MPI_REAL, MPI_SUM, comm2d, ierr )
6439       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
6440       CALL MPI_ALLREDUCE( lwc, lwc_total, 1 , &
6441       MPI_REAL, MPI_SUM, comm2d, ierr )
6442       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
6443       CALL MPI_ALLREDUCE( m1, m1_total, 1 , &
6444       MPI_REAL, MPI_SUM, comm2d, ierr )
6445       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
6446       CALL MPI_ALLREDUCE( m2, m2_total, 1 , &
6447       MPI_REAL, MPI_SUM, comm2d, ierr )
6448       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
6449       CALL MPI_ALLREDUCE( m3, m3_total, 1 , &
6450       MPI_REAL, MPI_SUM, comm2d, ierr )
6451#endif 
6452
6453!
6454!--    Calculate number concentration and mean volume averaged radius.
6455       nr_total = MERGE( nr_total / nrclgb_total,                              &
6456                         0.0_wp, nrclgb_total > 0.0_wp                         &
6457                       )
6458       rm_total = MERGE( ( rm_total /                                          &
6459                            ( nr_total * factor_volume_to_mass )               &
6460                          )**0.3333333_wp, 0.0_wp, nrclgb_total > 0.0_wp       &
6461                       )
6462!
6463!--    Check which function should be used to approximate the DSD.
6464       IF ( isf == 1 )  THEN
6465          lwc_total = MERGE( lwc_total / nrclgb_total,                         &
6466                             0.0_wp, nrclgb_total > 0.0_wp                     &
6467                           )
6468          m1_total  = MERGE( m1_total / nrclgb_total,                          &
6469                             0.0_wp, nrclgb_total > 0.0_wp                     &
6470                           )
6471          m2_total  = MERGE( m2_total / nrclgb_total,                          &
6472                             0.0_wp, nrclgb_total > 0.0_wp                     &
6473                           )
6474          m3_total  = MERGE( m3_total / nrclgb_total,                          &
6475                             0.0_wp, nrclgb_total > 0.0_wp                     &
6476                           )
6477          zeta = m1_total * m3_total / m2_total**2
6478          mu   = MAX( ( ( 1.0_wp - zeta ) * 2.0_wp + 1.0_wp ) /                &
6479                        ( zeta - 1.0_wp ), 0.0_wp                              &
6480                    )
6481
6482          lambda = ( pirho_l * nr_total / lwc_total *                          &
6483                     ( mu + 3.0_wp ) * ( mu + 2.0_wp ) * ( mu + 1.0_wp )       &
6484                   )**0.3333333_wp
6485          nr0 = nr_total / gamma( mu + 1.0_wp ) * lambda**( mu + 1.0_wp ) 
6486
6487          DO  n = 0, n_max-1
6488             diameter  = r_bin_mid(n) * 2.0_wp
6489             an_spl(n) = nr0 * diameter**mu * EXP( -lambda * diameter ) *      & 
6490                         ( r_bin(n+1) - r_bin(n) ) * 2.0_wp 
6491          ENDDO
6492       ELSEIF ( isf == 2 )  THEN
6493          DO  n = 0, n_max-1
6494             an_spl(n) = nr_total / ( SQRT( 2.0_wp * pi ) *                    &
6495                                     LOG(sigma_log) * r_bin_mid(n)             &
6496                                     ) *                                       &
6497                         EXP( -( LOG( r_bin_mid(n) / rm_total )**2 ) /         &
6498                               ( 2.0_wp * LOG(sigma_log)**2 )                  & 
6499                             ) *                                               & 
6500                         ( r_bin(n+1) - r_bin(n) )
6501          ENDDO
6502       ELSEIF( isf == 3 )  THEN
6503          DO  n = 0, n_max-1 
6504             an_spl(n) = 3.0_wp * nr_total * r_bin_mid(n)**2 / rm_total**3  *  &
6505                         EXP( - ( r_bin_mid(n)**3 / rm_total**3 ) )         *  &
6506                         ( r_bin(n+1) - r_bin(n) )
6507          ENDDO
6508       ENDIF
6509!
6510!--    Criterion to avoid super droplets with a weighting factor < 1.0.
6511       an_spl = MAX(an_spl, 1.0_wp)
6512
6513       DO  i = nxl, nxr
6514          DO  j = nys, nyn
6515             DO  k = nzb+1, nzt
6516                number_of_particles = prt_count(k,j,i)
6517                IF ( number_of_particles <= 0  .OR.                            &
6518                     ql(k,j,i) < ql_crit )  CYCLE
6519                particles => grid_particles(k,j,i)%particles(1:number_of_particles)
6520                new_particles_gb = 0
6521!
6522!--             Start splitting operations. Each particle is checked if it
6523!--             fulfilled the splitting criterion's. In splitting mode 'cl_av'
6524!--             a critical radius (radius_split) and a splitting function must
6525!--             be prescribed (see particles_par). The critical weighting factor
6526!--             is calculated while approximating a 'gamma', 'log' or 'exp'-
6527!--             drop size distribution. In this mode the DSD is calculated as
6528!--             an average over all cloudy grid boxes. Super droplets which
6529!--             have a larger radius and larger weighting factor are split into
6530!--             'splitting_factor' super droplets. In this case the splitting
6531!--             factor is calculated of weighting factor of the super droplet
6532!--             and the approximated number concentration for droplet of such
6533!--             a size. Due to the splitting, the weighting factor of the
6534!--             super droplet and all created clones is reduced by the factor
6535!--             of 'splitting_facor'.
6536                DO  n = 1, number_of_particles
6537                   DO  np = 0, n_max-1
6538                      IF ( r_bin(np) >= radius_split  .AND.                    &
6539                           particles(n)%particle_mask  .AND.                   &
6540                           particles(n)%radius >= r_bin(np)  .AND.             &
6541                           particles(n)%radius < r_bin(np+1)  .AND.            &
6542                           particles(n)%weight_factor >= an_spl(np)  )         &
6543                      THEN
6544!
6545!--                      Calculate splitting factor
6546                         splitting_factor =                                    & 
6547                             MIN( INT( particles(n)%weight_factor /            &
6548                                        an_spl(np)                             &
6549                                     ), splitting_factor_max                   &
6550                                )
6551                         IF ( splitting_factor < 2 )  CYCLE
6552!
6553!--                      Calculate the new number of particles.
6554                         new_size = prt_count(k,j,i) + splitting_factor - 1
6555!
6556!--                      Cycle if maximum number of particles per grid box
6557!--                      is greater than the allowed maximum number.
6558                         IF ( new_size >= max_number_particles_per_gridbox )   & 
6559                         CYCLE
6560!
6561!--                      Reallocate particle array if necessary.
6562                         IF ( new_size > SIZE(particles) )  THEN
6563                            CALL realloc_particles_array( i, j, k, new_size )
6564                         ENDIF
6565                         old_size  = prt_count(k,j,i)
6566                         new_particles_gb = new_particles_gb +                 &
6567                                            splitting_factor - 1
6568!
6569!--                      Calculate new weighting factor.
6570                         particles(n)%weight_factor =                          & 
6571                            particles(n)%weight_factor / splitting_factor
6572                         tmp_particle = particles(n)
6573!
6574!--                      Create splitting_factor-1 new particles.
6575                         DO  jpp = 1, splitting_factor-1
6576                            grid_particles(k,j,i)%particles(jpp+old_size) =    &
6577                                                                    tmp_particle
6578                         ENDDO
6579!
6580!--                      Save the new number of super droplets.
6581                         prt_count(k,j,i) = prt_count(k,j,i) +                 &
6582                                            splitting_factor - 1
6583                      ENDIF
6584                   ENDDO
6585                ENDDO 
6586
6587             ENDDO
6588          ENDDO
6589       ENDDO
6590
6591    ELSEIF ( i_splitting_mode == 3 )  THEN
6592
6593       DO  i = nxl, nxr
6594          DO  j = nys, nyn
6595             DO  k = nzb+1, nzt
6596
6597!
6598!--             Initialize summing variables.
6599                lwc = 0.0_wp
6600                m1  = 0.0_wp
6601                m2  = 0.0_wp
6602                m3  = 0.0_wp
6603                nr  = 0.0_wp
6604                rm  = 0.0_wp 
6605
6606                new_particles_gb = 0
6607                number_of_particles = prt_count(k,j,i)
6608                IF ( number_of_particles <= 0  .OR.                            & 
6609                     ql(k,j,i) < ql_crit )  CYCLE
6610                particles => grid_particles(k,j,i)%particles
6611!
6612!--             Calculate moments of DSD.
6613                DO  n = 1, number_of_particles
6614                   IF ( particles(n)%particle_mask  .AND.                      &
6615                        particles(n)%radius >= r_min )                         &
6616                   THEN
6617                      nr  = nr + particles(n)%weight_factor
6618                      rm  = rm + factor_volume_to_mass  *                      &
6619                                 particles(n)%radius**3  *                     &
6620                                 particles(n)%weight_factor
6621                      IF ( isf == 1 )  THEN
6622                         diameter   = particles(n)%radius * 2.0_wp
6623                         lwc = lwc + factor_volume_to_mass *                   &
6624                                     particles(n)%radius**3 *                  &
6625                                     particles(n)%weight_factor 
6626                         m1  = m1 + particles(n)%weight_factor * diameter
6627                         m2  = m2 + particles(n)%weight_factor * diameter**2
6628                         m3  = m3 + particles(n)%weight_factor * diameter**3
6629                      ENDIF
6630                   ENDIF
6631                ENDDO
6632
6633                IF ( nr <= 0.0_wp  .OR.  rm <= 0.0_wp )  CYCLE
6634!
6635!--             Calculate mean volume averaged radius.
6636                rm = ( rm / ( nr * factor_volume_to_mass ) )**0.3333333_wp
6637!
6638!--             Check which function should be used to approximate the DSD.
6639                IF ( isf == 1 )  THEN
6640!
6641!--                Gamma size distribution to calculate 
6642!--                critical weight_factor (e.g. Marshall + Palmer, 1948).
6643                   zeta = m1 * m3 / m2**2
6644                   mu   = MAX( ( ( 1.0_wp - zeta ) * 2.0_wp + 1.0_wp ) /       &
6645                                ( zeta - 1.0_wp ), 0.0_wp                      &
6646                             )   
6647                   lambda = ( pirho_l * nr / lwc *                             &
6648                              ( mu + 3.0_wp ) * ( mu + 2.0_wp ) *              &
6649                              ( mu + 1.0_wp )                                  &
6650                            )**0.3333333_wp
6651                   nr0 =  ( nr / (gamma( mu + 1.0_wp ) ) ) *                   &
6652                          lambda**( mu + 1.0_wp ) 
6653
6654                   DO  n = 0, n_max-1
6655                      diameter         = r_bin_mid(n) * 2.0_wp
6656                      an_spl(n) = nr0 * diameter**mu *                         &
6657                                  EXP( -lambda * diameter ) *                  & 
6658                                  ( r_bin(n+1) - r_bin(n) ) * 2.0_wp 
6659                   ENDDO
6660                ELSEIF ( isf == 2 )  THEN
6661!
6662!--                Lognormal size distribution to calculate critical
6663!--                weight_factor (e.g. Levin, 1971, Bradley + Stow, 1974).
6664                   DO  n = 0, n_max-1
6665                      an_spl(n) = nr / ( SQRT( 2.0_wp * pi ) *                 &
6666                                              LOG(sigma_log) * r_bin_mid(n)    &
6667                                        ) *                                    &
6668                                  EXP( -( LOG( r_bin_mid(n) / rm )**2 ) /      &
6669                                        ( 2.0_wp * LOG(sigma_log)**2 )         &
6670                                      ) *                                      &
6671                                  ( r_bin(n+1) - r_bin(n) )
6672                   ENDDO
6673                ELSEIF ( isf == 3 )  THEN
6674!
6675!--                Exponential size distribution to calculate critical
6676!--                weight_factor (e.g. Berry + Reinhardt, 1974). 
6677                   DO  n = 0, n_max-1
6678                      an_spl(n) = 3.0_wp * nr * r_bin_mid(n)**2 / rm**3 *     &
6679                                  EXP( - ( r_bin_mid(n)**3 / rm**3 ) ) *      &
6680                                  ( r_bin(n+1) - r_bin(n) )
6681                   ENDDO
6682                ENDIF
6683
6684!
6685!--             Criterion to avoid super droplets with a weighting factor < 1.0.
6686                an_spl = MAX(an_spl, 1.0_wp)
6687!
6688!--             Start splitting operations. Each particle is checked if it
6689!--             fulfilled the splitting criterion's. In splitting mode 'gb_av'
6690!--             a critical radius (radius_split) and a splitting function must
6691!--             be prescribed (see particles_par). The critical weighting factor
6692!--             is calculated while appoximating a 'gamma', 'log' or 'exp'-
6693!--             drop size distribution. In this mode a DSD is calculated for
6694!--             every cloudy grid box. Super droplets which have a larger
6695!--             radius and larger weighting factor are split into
6696!--             'splitting_factor' super droplets. In this case the splitting 
6697!--             factor is calculated of weighting factor of the super droplet 
6698!--             and theapproximated number concentration for droplet of such
6699!--             a size. Due to the splitting, the weighting factor of the 
6700!--             super droplet and all created clones is reduced by the factor 
6701!--             of 'splitting_facor'.
6702                DO  n = 1, number_of_particles
6703                   DO  np = 0, n_max-1
6704                      IF ( r_bin(np) >= radius_split  .AND.                    &
6705                           particles(n)%particle_mask  .AND.                   &
6706                           particles(n)%radius >= r_bin(np)    .AND.           &
6707                           particles(n)%radius < r_bin(np+1)   .AND.           &
6708                           particles(n)%weight_factor >= an_spl(np) )          &
6709                      THEN
6710!
6711!--                      Calculate splitting factor.
6712                         splitting_factor =                                    & 
6713                             MIN( INT( particles(n)%weight_factor /            &
6714                                        an_spl(np)                             &
6715                                     ), splitting_factor_max                   &
6716                                )
6717                         IF ( splitting_factor < 2 )  CYCLE
6718
6719!
6720!--                      Calculate the new number of particles.
6721                         new_size = prt_count(k,j,i) + splitting_factor - 1
6722!
6723!--                      Cycle if maximum number of particles per grid box
6724!--                      is greater than the allowed maximum number.
6725                         IF ( new_size >= max_number_particles_per_gridbox )   &
6726                         CYCLE
6727!
6728!--                      Reallocate particle array if necessary.
6729                         IF ( new_size > SIZE(particles) )  THEN
6730                            CALL realloc_particles_array( i, j, k, new_size )
6731                         ENDIF
6732!
6733!--                      Calculate new weighting factor.
6734                         particles(n)%weight_factor = & 
6735                            particles(n)%weight_factor / splitting_factor
6736                         tmp_particle               = particles(n)
6737                         old_size                   = prt_count(k,j,i)
6738!
6739!--                      Create splitting_factor-1 new particles.
6740                         DO  jpp = 1, splitting_factor-1
6741                            grid_particles(k,j,i)%particles( jpp + old_size ) = &
6742                               tmp_particle
6743                         ENDDO
6744!
6745!--                      Save the new number of droplets for every grid box.
6746                         prt_count(k,j,i)    = prt_count(k,j,i) +              &
6747                                               splitting_factor - 1
6748                         new_particles_gb    = new_particles_gb +              &
6749                                               splitting_factor - 1
6750                      ENDIF
6751                   ENDDO
6752                ENDDO
6753             ENDDO
6754          ENDDO
6755       ENDDO
6756    ENDIF
6757
6758    CALL cpu_log( log_point_s(80), 'lpm_splitting', 'stop' )
6759
6760 END SUBROUTINE lpm_splitting
6761 
6762
6763!------------------------------------------------------------------------------!
6764! Description:
6765! ------------
6766! This routine is a part of the Lagrangian particle model. Two Super droplets
6767! which fulfill certain criterion's (e.g. a big weighting factor and a small
6768! radius) can be merged into one super droplet with a increased number of
6769! represented particles of the super droplet. This mechanism ensures an
6770! improved a feasible amount of computational costs. The limits of particle
6771! creation should be chosen carefully! The idea of this algorithm is based on
6772! Unterstrasser and Soelch, 2014.
6773!------------------------------------------------------------------------------!
6774 SUBROUTINE lpm_merging
6775
6776    INTEGER(iwp) ::  i         !<
6777    INTEGER(iwp) ::  j         !<
6778    INTEGER(iwp) ::  k         !<
6779    INTEGER(iwp) ::  n         !<
6780    INTEGER(iwp) ::  merge_drp = 0     !< number of merged droplets
6781
6782
6783    REAL(wp) ::  ql_crit = 1.0E-5_wp  !< threshold lwc for cloudy grid cells
6784                                      !< (e.g. Siebesma et al 2003, JAS, 60)
6785
6786    CALL cpu_log( log_point_s(81), 'lpm_merging', 'start' )
6787
6788    merge_drp  = 0
6789
6790    IF ( weight_factor_merge == -1.0_wp )  THEN
6791       weight_factor_merge = 0.5_wp * initial_weighting_factor 
6792    ENDIF
6793
6794    DO  i = nxl, nxr
6795       DO  j = nys, nyn
6796          DO  k = nzb+1, nzt
6797
6798             number_of_particles = prt_count(k,j,i)
6799             IF ( number_of_particles <= 0  .OR.                               &
6800                   ql(k,j,i) >= ql_crit )  CYCLE
6801             particles => grid_particles(k,j,i)%particles(1:number_of_particles)
6802!
6803!--          Start merging operations: This routine delete super droplets with
6804!--          a small radius (radius <= radius_merge) and a low weighting
6805!--          factor (weight_factor  <= weight_factor_merge). The number of
6806!--          represented particles will be added to the next particle of the
6807!--          particle array. Tests showed that this simplified method can be
6808!--          used because it will only take place outside of cloudy grid
6809!--          boxes where ql <= 1.0E-5 kg/kg. Therefore, especially former cloned
6810!--          and subsequent evaporated super droplets will be merged.
6811             DO  n = 1, number_of_particles-1
6812                IF ( particles(n)%particle_mask                    .AND.       &
6813                     particles(n+1)%particle_mask                  .AND.       &
6814                     particles(n)%radius        <= radius_merge    .AND.       &
6815                     particles(n)%weight_factor <= weight_factor_merge )       &
6816                THEN
6817                   particles(n+1)%weight_factor  =                             &
6818                                       particles(n+1)%weight_factor +          &
6819                                       ( particles(n)%radius**3     /          &
6820                                         particles(n+1)%radius**3   *          &
6821                                         particles(n)%weight_factor            &
6822                                       )
6823                   particles(n)%particle_mask = .FALSE.
6824                   deleted_particles          = deleted_particles + 1 
6825                   merge_drp                  = merge_drp + 1
6826
6827                ENDIF
6828             ENDDO
6829          ENDDO
6830       ENDDO
6831    ENDDO
6832
6833
6834    CALL cpu_log( log_point_s(81), 'lpm_merging', 'stop' )
6835
6836 END SUBROUTINE lpm_merging
6837
6838 
6839
6840 
6841!------------------------------------------------------------------------------!
6842! Description:
6843! ------------
6844!> Exchange between subdomains.
6845!> As soon as one particle has moved beyond the boundary of the domain, it
6846!> is included in the relevant transfer arrays and marked for subsequent
6847!> deletion on this PE.
6848!> First sweep for crossings in x direction. Find out first the number of
6849!> particles to be transferred and allocate temporary arrays needed to store
6850!> them.
6851!> For a one-dimensional decomposition along y, no transfer is necessary,
6852!> because the particle remains on the PE, but the particle coordinate has to
6853!> be adjusted.
6854!------------------------------------------------------------------------------!
6855 SUBROUTINE lpm_exchange_horiz
6856
6857    INTEGER(iwp) ::  i                 !< grid index (x) of particle positition
6858    INTEGER(iwp) ::  ip                !< index variable along x
6859    INTEGER(iwp) ::  j                 !< grid index (y) of particle positition
6860    INTEGER(iwp) ::  jp                !< index variable along y
6861    INTEGER(iwp) ::  kp                !< index variable along z
6862    INTEGER(iwp) ::  n                 !< particle index variable
6863    INTEGER(iwp) ::  par_size          !< Particle size in bytes
6864    INTEGER(iwp) ::  trlp_count        !< number of particles send to left PE
6865    INTEGER(iwp) ::  trlp_count_recv   !< number of particles receive from right PE
6866    INTEGER(iwp) ::  trnp_count        !< number of particles send to north PE
6867    INTEGER(iwp) ::  trnp_count_recv   !< number of particles receive from south PE
6868    INTEGER(iwp) ::  trrp_count        !< number of particles send to right PE
6869    INTEGER(iwp) ::  trrp_count_recv   !< number of particles receive from left PE
6870    INTEGER(iwp) ::  trsp_count        !< number of particles send to south PE
6871    INTEGER(iwp) ::  trsp_count_recv   !< number of particles receive from north PE
6872
6873    TYPE(particle_type), DIMENSION(:), ALLOCATABLE ::  rvlp  !< particles received from right PE
6874    TYPE(particle_type), DIMENSION(:), ALLOCATABLE ::  rvnp  !< particles received from south PE
6875    TYPE(particle_type), DIMENSION(:), ALLOCATABLE ::  rvrp  !< particles received from left PE
6876    TYPE(particle_type), DIMENSION(:), ALLOCATABLE ::  rvsp  !< particles received from north PE
6877    TYPE(particle_type), DIMENSION(:), ALLOCATABLE ::  trlp  !< particles send to left PE
6878    TYPE(particle_type), DIMENSION(:), ALLOCATABLE ::  trnp  !< particles send to north PE
6879    TYPE(particle_type), DIMENSION(:), ALLOCATABLE ::  trrp  !< particles send to right PE
6880    TYPE(particle_type), DIMENSION(:), ALLOCATABLE ::  trsp  !< particles send to south PE
6881
6882    CALL cpu_log( log_point_s(23), 'lpm_exchange_horiz', 'start' )
6883
6884#if defined( __parallel )
6885
6886!
6887!-- Exchange between subdomains.
6888!-- As soon as one particle has moved beyond the boundary of the domain, it
6889!-- is included in the relevant transfer arrays and marked for subsequent
6890!-- deletion on this PE.
6891!-- First sweep for crossings in x direction. Find out first the number of
6892!-- particles to be transferred and allocate temporary arrays needed to store
6893!-- them.
6894!-- For a one-dimensional decomposition along y, no transfer is necessary,
6895!-- because the particle remains on the PE, but the particle coordinate has to
6896!-- be adjusted.
6897    trlp_count  = 0
6898    trrp_count  = 0
6899
6900    trlp_count_recv   = 0
6901    trrp_count_recv   = 0
6902
6903    IF ( pdims(1) /= 1 )  THEN
6904!
6905!--    First calculate the storage necessary for sending and receiving the data.
6906!--    Compute only first (nxl) and last (nxr) loop iterration.
6907       DO  ip = nxl, nxr, nxr - nxl
6908          DO  jp = nys, nyn
6909             DO  kp = nzb+1, nzt
6910
6911                number_of_particles = prt_count(kp,jp,ip)
6912                IF ( number_of_particles <= 0 )  CYCLE
6913                particles => grid_particles(kp,jp,ip)%particles(1:number_of_particles)
6914                DO  n = 1, number_of_particles
6915                   IF ( particles(n)%particle_mask )  THEN
6916                      i = particles(n)%x * ddx
6917!
6918!--                   Above calculation does not work for indices less than zero
6919                      IF ( particles(n)%x < 0.0_wp)  i = -1
6920
6921                      IF ( i < nxl )  THEN
6922                         trlp_count = trlp_count + 1
6923                      ELSEIF ( i > nxr )  THEN
6924                         trrp_count = trrp_count + 1
6925                      ENDIF
6926                   ENDIF
6927                ENDDO
6928
6929             ENDDO
6930          ENDDO
6931       ENDDO
6932
6933       IF ( trlp_count  == 0 )  trlp_count  = 1
6934       IF ( trrp_count  == 0 )  trrp_count  = 1
6935
6936       ALLOCATE( trlp(trlp_count), trrp(trrp_count) )
6937
6938       trlp = zero_particle
6939       trrp = zero_particle
6940
6941       trlp_count  = 0
6942       trrp_count  = 0
6943
6944    ENDIF
6945!
6946!-- Compute only first (nxl) and last (nxr) loop iterration
6947    DO  ip = nxl, nxr, nxr-nxl
6948       DO  jp = nys, nyn
6949          DO  kp = nzb+1, nzt
6950             number_of_particles = prt_count(kp,jp,ip)
6951             IF ( number_of_particles <= 0 )  CYCLE
6952             particles => grid_particles(kp,jp,ip)%particles(1:number_of_particles)
6953             DO  n = 1, number_of_particles
6954!
6955!--             Only those particles that have not been marked as 'deleted' may
6956!--             be moved.
6957                IF ( particles(n)%particle_mask )  THEN
6958
6959                   i = particles(n)%x * ddx
6960!
6961!--                Above calculation does not work for indices less than zero
6962                   IF ( particles(n)%x < 0.0_wp )  i = -1
6963
6964                   IF ( i <  nxl )  THEN
6965                      IF ( i < 0 )  THEN
6966!
6967!--                   Apply boundary condition along x
6968                         IF ( ibc_par_lr == 0 )  THEN
6969!
6970!--                         Cyclic condition
6971                            IF ( pdims(1) == 1 )  THEN
6972                               particles(n)%x        = ( nx + 1 ) * dx + particles(n)%x
6973                               particles(n)%origin_x = ( nx + 1 ) * dx + &
6974                               particles(n)%origin_x
6975                            ELSE
6976                               trlp_count = trlp_count + 1
6977                               trlp(trlp_count)   = particles(n)
6978                               trlp(trlp_count)%x = ( nx + 1 ) * dx + trlp(trlp_count)%x
6979                               trlp(trlp_count)%origin_x = trlp(trlp_count)%origin_x + &
6980                               ( nx + 1 ) * dx
6981                               particles(n)%particle_mask  = .FALSE.
6982                               deleted_particles = deleted_particles + 1
6983
6984                               IF ( trlp(trlp_count)%x >= (nx + 1)* dx - 1.0E-12_wp )  THEN
6985                                  trlp(trlp_count)%x = trlp(trlp_count)%x - 1.0E-10_wp
6986                                  !++ why is 1 subtracted in next statement???
6987                                  trlp(trlp_count)%origin_x = trlp(trlp_count)%origin_x - 1
6988                               ENDIF
6989
6990                            ENDIF
6991
6992                         ELSEIF ( ibc_par_lr == 1 )  THEN
6993!
6994!--                         Particle absorption
6995                            particles(n)%particle_mask = .FALSE.
6996                            deleted_particles = deleted_particles + 1
6997
6998                         ELSEIF ( ibc_par_lr == 2 )  THEN
6999!
7000!--                         Particle reflection
7001                            particles(n)%x       = -particles(n)%x
7002                            particles(n)%speed_x = -particles(n)%speed_x
7003
7004                         ENDIF
7005                      ELSE
7006!
7007!--                      Store particle data in the transfer array, which will be
7008!--                      send to the neighbouring PE
7009                         trlp_count = trlp_count + 1
7010                         trlp(trlp_count) = particles(n)
7011                         particles(n)%particle_mask = .FALSE.
7012                         deleted_particles = deleted_particles + 1
7013
7014                      ENDIF
7015
7016                   ELSEIF ( i > nxr )  THEN
7017                      IF ( i > nx )  THEN
7018!
7019!--                      Apply boundary condition along x
7020                         IF ( ibc_par_lr == 0 )  THEN
7021!
7022!--                         Cyclic condition
7023                            IF ( pdims(1) == 1 )  THEN
7024                               particles(n)%x = particles(n)%x - ( nx + 1 ) * dx
7025                               particles(n)%origin_x = particles(n)%origin_x - &
7026                               ( nx + 1 ) * dx
7027                            ELSE
7028                               trrp_count = trrp_count + 1
7029                               trrp(trrp_count) = particles(n)
7030                               trrp(trrp_count)%x = trrp(trrp_count)%x - ( nx + 1 ) * dx
7031                               trrp(trrp_count)%origin_x = trrp(trrp_count)%origin_x - &
7032                               ( nx + 1 ) * dx
7033                               particles(n)%particle_mask = .FALSE.
7034                               deleted_particles = deleted_particles + 1
7035
7036                            ENDIF
7037
7038                         ELSEIF ( ibc_par_lr == 1 )  THEN
7039!
7040!--                         Particle absorption
7041                            particles(n)%particle_mask = .FALSE.
7042                            deleted_particles = deleted_particles + 1
7043
7044                         ELSEIF ( ibc_par_lr == 2 )  THEN
7045!
7046!--                         Particle reflection
7047                            particles(n)%x       = 2 * ( nx * dx ) - particles(n)%x
7048                            particles(n)%speed_x = -particles(n)%speed_x
7049
7050                         ENDIF
7051                      ELSE
7052!
7053!--                      Store particle data in the transfer array, which will be send
7054!--                      to the neighbouring PE
7055                         trrp_count = trrp_count + 1
7056                         trrp(trrp_count) = particles(n)
7057                         particles(n)%particle_mask = .FALSE.
7058                         deleted_particles = deleted_particles + 1
7059
7060                      ENDIF
7061
7062                   ENDIF
7063                ENDIF
7064
7065             ENDDO
7066          ENDDO
7067       ENDDO
7068    ENDDO
7069
7070!
7071!-- STORAGE_SIZE returns the storage size of argument A in bits. However , it
7072!-- is needed in bytes. The function C_SIZEOF which produces this value directly
7073!-- causes problems with gfortran. For this reason the use of C_SIZEOF is avoided
7074    par_size = STORAGE_SIZE(trlp(1))/8
7075
7076
7077!
7078!-- Allocate arrays required for north-south exchange, as these
7079!-- are used directly after particles are exchange along x-direction.
7080    ALLOCATE( move_also_north(1:NR_2_direction_move) )
7081    ALLOCATE( move_also_south(1:NR_2_direction_move) )
7082
7083    nr_move_north = 0
7084    nr_move_south = 0
7085!
7086!-- Send left boundary, receive right boundary (but first exchange how many
7087!-- and check, if particle storage must be extended)
7088    IF ( pdims(1) /= 1 )  THEN
7089
7090       CALL MPI_SENDRECV( trlp_count,      1, MPI_INTEGER, pleft,  0, &
7091                          trrp_count_recv, 1, MPI_INTEGER, pright, 0, &
7092                          comm2d, status, ierr )
7093
7094       ALLOCATE(rvrp(MAX(1,trrp_count_recv)))
7095
7096       CALL MPI_SENDRECV( trlp, max(1,trlp_count)*par_size, MPI_BYTE,&
7097                          pleft, 1, rvrp,                            &
7098                          max(1,trrp_count_recv)*par_size, MPI_BYTE, pright, 1,&
7099                          comm2d, status, ierr )
7100
7101       IF ( trrp_count_recv > 0 )  CALL lpm_add_particles_to_gridcell(rvrp(1:trrp_count_recv))
7102
7103       DEALLOCATE(rvrp)
7104
7105!
7106!--    Send right boundary, receive left boundary
7107       CALL MPI_SENDRECV( trrp_count,      1, MPI_INTEGER, pright, 0, &
7108                          trlp_count_recv, 1, MPI_INTEGER, pleft,  0, &
7109                          comm2d, status, ierr )
7110
7111       ALLOCATE(rvlp(MAX(1,trlp_count_recv)))
7112!
7113!--    This MPI_SENDRECV should work even with odd mixture on 32 and 64 Bit
7114!--    variables in structure particle_type (due to the calculation of par_size)
7115       CALL MPI_SENDRECV( trrp, max(1,trrp_count)*par_size, MPI_BYTE,&
7116                          pright, 1, rvlp,                           &
7117                          max(1,trlp_count_recv)*par_size, MPI_BYTE, pleft, 1, &
7118                          comm2d, status, ierr )
7119
7120       IF ( trlp_count_recv > 0 )  CALL lpm_add_particles_to_gridcell(rvlp(1:trlp_count_recv))
7121
7122       DEALLOCATE( rvlp )
7123       DEALLOCATE( trlp, trrp )
7124
7125    ENDIF
7126
7127!
7128!-- Check whether particles have crossed the boundaries in y direction. Note
7129!-- that this case can also apply to particles that have just been received
7130!-- from the adjacent right or left PE.
7131!-- Find out first the number of particles to be transferred and allocate
7132!-- temporary arrays needed to store them.
7133!-- For a one-dimensional decomposition along y, no transfer is necessary,
7134!-- because the particle remains on the PE.
7135    trsp_count  = nr_move_south
7136    trnp_count  = nr_move_north
7137
7138    trsp_count_recv   = 0
7139    trnp_count_recv   = 0
7140
7141    IF ( pdims(2) /= 1 )  THEN
7142!
7143!--    First calculate the storage necessary for sending and receiving the
7144!--    data
7145       DO  ip = nxl, nxr
7146          DO  jp = nys, nyn, nyn-nys    !compute only first (nys) and last (nyn) loop iterration
7147             DO  kp = nzb+1, nzt
7148                number_of_particles = prt_count(kp,jp,ip)
7149                IF ( number_of_particles <= 0 )  CYCLE
7150                particles => grid_particles(kp,jp,ip)%particles(1:number_of_particles)
7151                DO  n = 1, number_of_particles
7152                   IF ( particles(n)%particle_mask )  THEN
7153                      j = particles(n)%y * ddy
7154!
7155!--                   Above calculation does not work for indices less than zero
7156                      IF ( particles(n)%y < 0.0_wp)  j = -1
7157
7158                      IF ( j < nys )  THEN
7159                         trsp_count = trsp_count + 1
7160                      ELSEIF ( j > nyn )  THEN
7161                         trnp_count = trnp_count + 1
7162                      ENDIF
7163                   ENDIF
7164                ENDDO
7165             ENDDO
7166          ENDDO
7167       ENDDO
7168
7169       IF ( trsp_count  == 0 )  trsp_count  = 1
7170       IF ( trnp_count  == 0 )  trnp_count  = 1
7171
7172       ALLOCATE( trsp(trsp_count), trnp(trnp_count) )
7173
7174       trsp = zero_particle
7175       trnp = zero_particle
7176
7177       trsp_count  = nr_move_south
7178       trnp_count  = nr_move_north
7179
7180       trsp(1:nr_move_south) = move_also_south(1:nr_move_south)
7181       trnp(1:nr_move_north) = move_also_north(1:nr_move_north)
7182
7183    ENDIF
7184
7185    DO  ip = nxl, nxr
7186       DO  jp = nys, nyn, nyn-nys ! compute only first (nys) and last (nyn) loop iterration
7187          DO  kp = nzb+1, nzt
7188             number_of_particles = prt_count(kp,jp,ip)
7189             IF ( number_of_particles <= 0 )  CYCLE
7190             particles => grid_particles(kp,jp,ip)%particles(1:number_of_particles)
7191             DO  n = 1, number_of_particles
7192!
7193!--             Only those particles that have not been marked as 'deleted' may
7194!--             be moved.
7195                IF ( particles(n)%particle_mask )  THEN
7196
7197                   j = particles(n)%y * ddy
7198!
7199!--                Above calculation does not work for indices less than zero
7200                   IF ( particles(n)%y < 0.0_wp )  j = -1
7201
7202                   IF ( j < nys )  THEN
7203                      IF ( j < 0 )  THEN
7204!
7205!--                      Apply boundary condition along y
7206                         IF ( ibc_par_ns == 0 )  THEN
7207!
7208!--                         Cyclic condition
7209                            IF ( pdims(2) == 1 )  THEN
7210                               particles(n)%y = ( ny + 1 ) * dy + particles(n)%y
7211                               particles(n)%origin_y = ( ny + 1 ) * dy + &
7212                                                     particles(n)%origin_y
7213                            ELSE
7214                               trsp_count         = trsp_count + 1
7215                               trsp(trsp_count)   = particles(n)
7216                               trsp(trsp_count)%y = ( ny + 1 ) * dy + &
7217                                                 trsp(trsp_count)%y
7218                               trsp(trsp_count)%origin_y = trsp(trsp_count)%origin_y &
7219                                                + ( ny + 1 ) * dy
7220                               particles(n)%particle_mask = .FALSE.
7221                               deleted_particles = deleted_particles + 1
7222
7223                               IF ( trsp(trsp_count)%y >= (ny+1)* dy - 1.0E-12_wp )  THEN
7224                                  trsp(trsp_count)%y = trsp(trsp_count)%y - 1.0E-10_wp
7225                                  !++ why is 1 subtracted in next statement???
7226                                  trsp(trsp_count)%origin_y =                        &
7227                                                  trsp(trsp_count)%origin_y - 1
7228                               ENDIF
7229
7230                            ENDIF
7231
7232                         ELSEIF ( ibc_par_ns == 1 )  THEN
7233!
7234!--                         Particle absorption
7235                            particles(n)%particle_mask = .FALSE.
7236                            deleted_particles          = deleted_particles + 1
7237
7238                         ELSEIF ( ibc_par_ns == 2 )  THEN
7239!
7240!--                         Particle reflection
7241                            particles(n)%y       = -particles(n)%y
7242                            particles(n)%speed_y = -particles(n)%speed_y
7243
7244                         ENDIF
7245                      ELSE
7246!
7247!--                      Store particle data in the transfer array, which will
7248!--                      be send to the neighbouring PE
7249                         trsp_count = trsp_count + 1
7250                         trsp(trsp_count) = particles(n)
7251                         particles(n)%particle_mask = .FALSE.
7252                         deleted_particles = deleted_particles + 1
7253
7254                      ENDIF
7255
7256                   ELSEIF ( j > nyn )  THEN
7257                      IF ( j > ny )  THEN
7258!
7259!--                       Apply boundary condition along y
7260                         IF ( ibc_par_ns == 0 )  THEN
7261!
7262!--                         Cyclic condition
7263                            IF ( pdims(2) == 1 )  THEN
7264                               particles(n)%y        = particles(n)%y - ( ny + 1 ) * dy
7265                               particles(n)%origin_y =                         &
7266                                          particles(n)%origin_y - ( ny + 1 ) * dy
7267                            ELSE
7268                               trnp_count         = trnp_count + 1
7269                               trnp(trnp_count)   = particles(n)
7270                               trnp(trnp_count)%y =                            &
7271                                          trnp(trnp_count)%y - ( ny + 1 ) * dy
7272                               trnp(trnp_count)%origin_y =                     &
7273                                         trnp(trnp_count)%origin_y - ( ny + 1 ) * dy
7274                               particles(n)%particle_mask = .FALSE.
7275                               deleted_particles          = deleted_particles + 1
7276                            ENDIF
7277
7278                         ELSEIF ( ibc_par_ns == 1 )  THEN
7279!
7280!--                         Particle absorption
7281                            particles(n)%particle_mask = .FALSE.
7282                            deleted_particles = deleted_particles + 1
7283
7284                         ELSEIF ( ibc_par_ns == 2 )  THEN
7285!
7286!--                         Particle reflection
7287                            particles(n)%y       = 2 * ( ny * dy ) - particles(n)%y
7288                            particles(n)%speed_y = -particles(n)%speed_y
7289
7290                         ENDIF
7291                      ELSE
7292!
7293!--                      Store particle data in the transfer array, which will
7294!--                      be send to the neighbouring PE
7295                         trnp_count = trnp_count + 1
7296                         trnp(trnp_count) = particles(n)
7297                         particles(n)%particle_mask = .FALSE.
7298                         deleted_particles = deleted_particles + 1
7299
7300                      ENDIF
7301
7302                   ENDIF
7303                ENDIF
7304             ENDDO
7305          ENDDO
7306       ENDDO
7307    ENDDO
7308
7309!
7310!-- Send front boundary, receive back boundary (but first exchange how many
7311!-- and check, if particle storage must be extended)
7312    IF ( pdims(2) /= 1 )  THEN
7313
7314       CALL MPI_SENDRECV( trsp_count,      1, MPI_INTEGER, psouth, 0, &
7315                          trnp_count_recv, 1, MPI_INTEGER, pnorth, 0, &
7316                          comm2d, status, ierr )
7317
7318       ALLOCATE(rvnp(MAX(1,trnp_count_recv)))
7319!
7320!--    This MPI_SENDRECV should work even with odd mixture on 32 and 64 Bit
7321!--    variables in structure particle_type (due to the calculation of par_size)
7322       CALL MPI_SENDRECV( trsp, trsp_count*par_size, MPI_BYTE,      &
7323                          psouth, 1, rvnp,                             &
7324                          trnp_count_recv*par_size, MPI_BYTE, pnorth, 1,   &
7325                          comm2d, status, ierr )
7326
7327       IF ( trnp_count_recv  > 0 )  CALL lpm_add_particles_to_gridcell(rvnp(1:trnp_count_recv))
7328
7329       DEALLOCATE(rvnp)
7330
7331!
7332!--    Send back boundary, receive front boundary
7333       CALL MPI_SENDRECV( trnp_count,      1, MPI_INTEGER, pnorth, 0, &
7334                          trsp_count_recv, 1, MPI_INTEGER, psouth, 0, &
7335                          comm2d, status, ierr )
7336
7337       ALLOCATE(rvsp(MAX(1,trsp_count_recv)))
7338!
7339!--    This MPI_SENDRECV should work even with odd mixture on 32 and 64 Bit
7340!--    variables in structure particle_type (due to the calculation of par_size)
7341       CALL MPI_SENDRECV( trnp, trnp_count*par_size, MPI_BYTE,      &
7342                          pnorth, 1, rvsp,                          &
7343                          trsp_count_recv*par_size, MPI_BYTE, psouth, 1,   &
7344                          comm2d, status, ierr )
7345
7346       IF ( trsp_count_recv > 0 )  CALL lpm_add_particles_to_gridcell(rvsp(1:trsp_count_recv))
7347
7348       DEALLOCATE(rvsp)
7349
7350       number_of_particles = number_of_particles + trsp_count_recv
7351
7352       DEALLOCATE( trsp, trnp )
7353
7354    ENDIF
7355
7356    DEALLOCATE( move_also_north )
7357    DEALLOCATE( move_also_south )
7358
7359#else
7360
7361    DO  ip = nxl, nxr, nxr-nxl
7362       DO  jp = nys, nyn
7363          DO  kp = nzb+1, nzt
7364             number_of_particles = prt_count(kp,jp,ip)
7365             IF ( number_of_particles <= 0 )  CYCLE
7366             particles => grid_particles(kp,jp,ip)%particles(1:number_of_particles)
7367             DO  n = 1, number_of_particles
7368!
7369!--             Apply boundary conditions
7370
7371                IF ( particles(n)%x < 0.0_wp )  THEN
7372
7373                   IF ( ibc_par_lr == 0 )  THEN
7374!
7375!--                   Cyclic boundary. Relevant coordinate has to be changed.
7376                      particles(n)%x = ( nx + 1 ) * dx + particles(n)%x
7377                      particles(n)%origin_x = ( nx + 1 ) * dx + &
7378                               particles(n)%origin_x
7379                   ELSEIF ( ibc_par_lr == 1 )  THEN
7380!
7381!--                   Particle absorption
7382                      particles(n)%particle_mask = .FALSE.
7383                      deleted_particles = deleted_particles + 1
7384
7385                   ELSEIF ( ibc_par_lr == 2 )  THEN
7386!
7387!--                   Particle reflection
7388                      particles(n)%x       = -dx - particles(n)%x
7389                      particles(n)%speed_x = -particles(n)%speed_x
7390                   ENDIF
7391
7392                ELSEIF ( particles(n)%x >= ( nx + 1) * dx )  THEN
7393
7394                   IF ( ibc_par_lr == 0 )  THEN
7395!
7396!--                   Cyclic boundary. Relevant coordinate has to be changed.
7397                      particles(n)%x = particles(n)%x - ( nx + 1 ) * dx
7398                      particles(n)%origin_x = particles(n)%origin_x - &
7399                               ( nx + 1 ) * dx
7400
7401                   ELSEIF ( ibc_par_lr == 1 )  THEN
7402!
7403!--                   Particle absorption
7404                      particles(n)%particle_mask = .FALSE.
7405                      deleted_particles = deleted_particles + 1
7406
7407                   ELSEIF ( ibc_par_lr == 2 )  THEN
7408!
7409!--                   Particle reflection
7410                      particles(n)%x       = ( nx + 1 ) * dx - particles(n)%x
7411                      particles(n)%speed_x = -particles(n)%speed_x
7412                   ENDIF
7413
7414                ENDIF
7415             ENDDO
7416          ENDDO
7417       ENDDO
7418    ENDDO
7419
7420    DO  ip = nxl, nxr
7421       DO  jp = nys, nyn, nyn-nys
7422          DO  kp = nzb+1, nzt
7423             number_of_particles = prt_count(kp,jp,ip)
7424             IF ( number_of_particles <= 0 )  CYCLE
7425             particles => grid_particles(kp,jp,ip)%particles(1:number_of_particles)
7426             DO  n = 1, number_of_particles
7427
7428                IF ( particles(n)%y < 0.0_wp)  THEN
7429
7430                   IF ( ibc_par_ns == 0 )  THEN
7431!
7432!--                   Cyclic boundary. Relevant coordinate has to be changed.
7433                      particles(n)%y = ( ny + 1 ) * dy + particles(n)%y
7434                      particles(n)%origin_y = ( ny + 1 ) * dy + &
7435                           particles(n)%origin_y
7436
7437                   ELSEIF ( ibc_par_ns == 1 )  THEN
7438!
7439!--                   Particle absorption
7440                      particles(n)%particle_mask = .FALSE.
7441                      deleted_particles = deleted_particles + 1
7442
7443                   ELSEIF ( ibc_par_ns == 2 )  THEN
7444!
7445!--                   Particle reflection
7446                      particles(n)%y       = -dy - particles(n)%y
7447                      particles(n)%speed_y = -particles(n)%speed_y
7448                   ENDIF
7449
7450                ELSEIF ( particles(n)%y >= ( ny + 1) * dy )  THEN
7451
7452                   IF ( ibc_par_ns == 0 )  THEN
7453!
7454!--                   Cyclic boundary. Relevant coordinate has to be changed.
7455                      particles(n)%y = particles(n)%y - ( ny + 1 ) * dy
7456                      particles(n)%origin_y = particles(n)%origin_y - &
7457                                ( ny + 1 ) * dy
7458
7459                   ELSEIF ( ibc_par_ns == 1 )  THEN
7460!
7461!--                   Particle absorption
7462                      particles(n)%particle_mask = .FALSE.
7463                      deleted_particles = deleted_particles + 1
7464
7465                   ELSEIF ( ibc_par_ns == 2 )  THEN
7466!
7467!--                   Particle reflection
7468                      particles(n)%y       = ( ny + 1 ) * dy - particles(n)%y
7469                      particles(n)%speed_y = -particles(n)%speed_y
7470                   ENDIF
7471
7472                ENDIF
7473
7474             ENDDO
7475          ENDDO
7476       ENDDO
7477    ENDDO
7478#endif
7479
7480!
7481!-- Accumulate the number of particles transferred between the subdomains
7482#if defined( __parallel )
7483    trlp_count_sum      = trlp_count_sum      + trlp_count
7484    trlp_count_recv_sum = trlp_count_recv_sum + trlp_count_recv
7485    trrp_count_sum      = trrp_count_sum      + trrp_count
7486    trrp_count_recv_sum = trrp_count_recv_sum + trrp_count_recv
7487    trsp_count_sum      = trsp_count_sum      + trsp_count
7488    trsp_count_recv_sum = trsp_count_recv_sum + trsp_count_recv
7489    trnp_count_sum      = trnp_count_sum      + trnp_count
7490    trnp_count_recv_sum = trnp_count_recv_sum + trnp_count_recv
7491#endif
7492
7493    CALL cpu_log( log_point_s(23), 'lpm_exchange_horiz', 'stop' )
7494
7495 END SUBROUTINE lpm_exchange_horiz
7496
7497!------------------------------------------------------------------------------!
7498! Description:
7499! ------------
7500!> If a particle moves from one processor to another, this subroutine moves
7501!> the corresponding elements from the particle arrays of the old grid cells
7502!> to the particle arrays of the new grid cells.
7503!------------------------------------------------------------------------------!
7504 SUBROUTINE lpm_add_particles_to_gridcell (particle_array)
7505
7506    IMPLICIT NONE
7507
7508    INTEGER(iwp)        ::  ip        !< grid index (x) of particle
7509    INTEGER(iwp)        ::  jp        !< grid index (x) of particle
7510    INTEGER(iwp)        ::  kp        !< grid index (x) of particle
7511    INTEGER(iwp)        ::  n         !< index variable of particle
7512    INTEGER(iwp)        ::  pindex    !< dummy argument for new number of particles per grid box
7513
7514    LOGICAL             ::  pack_done !<
7515
7516    TYPE(particle_type), DIMENSION(:), INTENT(IN)  ::  particle_array !< new particles in a grid box
7517    TYPE(particle_type), DIMENSION(:), ALLOCATABLE ::  temp_ns        !< temporary particle array for reallocation
7518
7519    pack_done     = .FALSE.
7520
7521    DO  n = 1, SIZE(particle_array)
7522
7523       IF ( .NOT. particle_array(n)%particle_mask )  CYCLE
7524
7525       ip = particle_array(n)%x * ddx
7526       jp = particle_array(n)%y * ddy
7527!
7528!--    In case of stretching the actual k index must be found
7529       IF ( dz_stretch_level /= -9999999.9_wp  .OR.         &
7530            dz_stretch_level_start(1) /= -9999999.9_wp )  THEN
7531          kp = MINLOC( ABS( particle_array(n)%z - zu ), DIM = 1 ) - 1
7532       ELSE
7533          kp = INT( particle_array(n)%z / dz(1) + 1 + offset_ocean_nzt )
7534       ENDIF
7535
7536       IF ( ip >= nxl  .AND.  ip <= nxr  .AND.  jp >= nys  .AND.  jp <= nyn    &
7537            .AND.  kp >= nzb+1  .AND.  kp <= nzt)  THEN ! particle stays on processor
7538          number_of_particles = prt_count(kp,jp,ip)
7539          particles => grid_particles(kp,jp,ip)%particles(1:number_of_particles)
7540
7541          pindex = prt_count(kp,jp,ip)+1
7542          IF( pindex > SIZE(grid_particles(kp,jp,ip)%particles) )  THEN
7543             IF ( pack_done )  THEN
7544                CALL realloc_particles_array ( ip, jp, kp )
7545             ELSE
7546                CALL lpm_pack
7547                prt_count(kp,jp,ip) = number_of_particles
7548                pindex = prt_count(kp,jp,ip)+1
7549                IF ( pindex > SIZE(grid_particles(kp,jp,ip)%particles) )  THEN
7550                   CALL realloc_particles_array ( ip, jp, kp )
7551                ENDIF
7552                pack_done = .TRUE.
7553             ENDIF
7554          ENDIF
7555          grid_particles(kp,jp,ip)%particles(pindex) = particle_array(n)
7556          prt_count(kp,jp,ip) = pindex
7557       ELSE
7558          IF ( jp <= nys - 1 )  THEN
7559             nr_move_south = nr_move_south+1
7560!
7561!--          Before particle information is swapped to exchange-array, check
7562!--          if enough memory is allocated. If required, reallocate exchange
7563!--          array.
7564             IF ( nr_move_south > SIZE(move_also_south) )  THEN
7565!
7566!--             At first, allocate further temporary array to swap particle
7567!--             information.
7568                ALLOCATE( temp_ns(SIZE(move_also_south)+NR_2_direction_move) )
7569                temp_ns(1:nr_move_south-1) = move_also_south(1:nr_move_south-1)
7570                DEALLOCATE( move_also_south )
7571                ALLOCATE( move_also_south(SIZE(temp_ns)) )
7572                move_also_south(1:nr_move_south-1) = temp_ns(1:nr_move_south-1)
7573                DEALLOCATE( temp_ns )
7574
7575             ENDIF
7576
7577             move_also_south(nr_move_south) = particle_array(n)
7578
7579             IF ( jp == -1 )  THEN
7580!
7581!--             Apply boundary condition along y
7582                IF ( ibc_par_ns == 0 )  THEN
7583                   move_also_south(nr_move_south)%y =                          &
7584                      move_also_south(nr_move_south)%y + ( ny + 1 ) * dy
7585                   move_also_south(nr_move_south)%origin_y =                   &
7586                      move_also_south(nr_move_south)%origin_y + ( ny + 1 ) * dy
7587                ELSEIF ( ibc_par_ns == 1 )  THEN
7588!
7589!--                Particle absorption
7590                   move_also_south(nr_move_south)%particle_mask = .FALSE.
7591                   deleted_particles = deleted_particles + 1
7592
7593                ELSEIF ( ibc_par_ns == 2 )  THEN
7594!
7595!--                Particle reflection
7596                   move_also_south(nr_move_south)%y       =                    &
7597                      -move_also_south(nr_move_south)%y
7598                   move_also_south(nr_move_south)%speed_y =                    &
7599                      -move_also_south(nr_move_south)%speed_y
7600
7601                ENDIF
7602             ENDIF
7603          ELSEIF ( jp >= nyn+1 )  THEN
7604             nr_move_north = nr_move_north+1
7605!
7606!--          Before particle information is swapped to exchange-array, check
7607!--          if enough memory is allocated. If required, reallocate exchange
7608!--          array.
7609             IF ( nr_move_north > SIZE(move_also_north) )  THEN
7610!
7611!--             At first, allocate further temporary array to swap particle
7612!--             information.
7613                ALLOCATE( temp_ns(SIZE(move_also_north)+NR_2_direction_move) )
7614                temp_ns(1:nr_move_north-1) = move_also_south(1:nr_move_north-1)
7615                DEALLOCATE( move_also_north )
7616                ALLOCATE( move_also_north(SIZE(temp_ns)) )
7617                move_also_north(1:nr_move_north-1) = temp_ns(1:nr_move_north-1)
7618                DEALLOCATE( temp_ns )
7619
7620             ENDIF
7621
7622             move_also_north(nr_move_north) = particle_array(n)
7623             IF ( jp == ny+1 )  THEN
7624!
7625!--             Apply boundary condition along y
7626                IF ( ibc_par_ns == 0 )  THEN
7627
7628                   move_also_north(nr_move_north)%y =                          &
7629                      move_also_north(nr_move_north)%y - ( ny + 1 ) * dy
7630                   move_also_north(nr_move_north)%origin_y =                   &
7631                      move_also_north(nr_move_north)%origin_y - ( ny + 1 ) * dy
7632                ELSEIF ( ibc_par_ns == 1 )  THEN
7633!
7634!--                Particle absorption
7635                   move_also_north(nr_move_north)%particle_mask = .FALSE.
7636                   deleted_particles = deleted_particles + 1
7637
7638                ELSEIF ( ibc_par_ns == 2 )  THEN
7639!
7640!--                Particle reflection
7641                   move_also_north(nr_move_north)%y       =                    &
7642                      -move_also_north(nr_move_north)%y
7643                   move_also_north(nr_move_north)%speed_y =                    &
7644                      -move_also_north(nr_move_north)%speed_y
7645
7646                ENDIF
7647             ENDIF
7648          ELSE
7649             WRITE(0,'(a,8i7)') 'particle out of range ',myid,ip,jp,kp,nxl,nxr,nys,nyn
7650          ENDIF
7651       ENDIF
7652    ENDDO
7653
7654    RETURN
7655
7656 END SUBROUTINE lpm_add_particles_to_gridcell
7657 
7658 
7659!------------------------------------------------------------------------------!
7660! Description:
7661! ------------
7662!> If a particle moves from one grid cell to another (on the current
7663!> processor!), this subroutine moves the corresponding element from the
7664!> particle array of the old grid cell to the particle array of the new grid
7665!> cell.
7666!------------------------------------------------------------------------------!
7667 SUBROUTINE lpm_move_particle
7668 
7669    INTEGER(iwp)        ::  i           !< grid index (x) of particle position
7670    INTEGER(iwp)        ::  ip          !< index variable along x
7671    INTEGER(iwp)        ::  j           !< grid index (y) of particle position
7672    INTEGER(iwp)        ::  jp          !< index variable along y
7673    INTEGER(iwp)        ::  k           !< grid index (z) of particle position
7674    INTEGER(iwp)        ::  kp          !< index variable along z
7675    INTEGER(iwp)        ::  n           !< index variable for particle array
7676    INTEGER(iwp)        ::  np_before_move !< number of particles per grid box before moving
7677    INTEGER(iwp)        ::  pindex      !< dummy argument for number of new particle per grid box
7678
7679    TYPE(particle_type), DIMENSION(:), POINTER  ::  particles_before_move !< particles before moving
7680
7681    CALL cpu_log( log_point_s(41), 'lpm_move_particle', 'start' )
7682    CALL lpm_check_cfl
7683    DO  ip = nxl, nxr
7684       DO  jp = nys, nyn
7685          DO  kp = nzb+1, nzt
7686
7687             np_before_move = prt_count(kp,jp,ip)
7688             IF ( np_before_move <= 0 )  CYCLE
7689             particles_before_move => grid_particles(kp,jp,ip)%particles(1:np_before_move)
7690
7691             DO  n = 1, np_before_move
7692                i = particles_before_move(n)%x * ddx
7693                j = particles_before_move(n)%y * ddy
7694                k = kp
7695!
7696!--             Find correct vertical particle grid box (necessary in case of grid stretching)
7697!--             Due to the CFL limitations only the neighbouring grid boxes are considered.
7698                IF( zw(k)   < particles_before_move(n)%z ) k = k + 1
7699                IF( zw(k-1) > particles_before_move(n)%z ) k = k - 1 
7700
7701!--             For lpm_exchange_horiz to work properly particles need to be moved to the outermost gridboxes
7702!--             of the respective processor. If the particle index is inside the processor the following lines
7703!--             will not change the index
7704                i = MIN ( i , nxr )
7705                i = MAX ( i , nxl )
7706                j = MIN ( j , nyn )
7707                j = MAX ( j , nys )
7708
7709                k = MIN ( k , nzt )
7710                k = MAX ( k , nzb+1 )
7711
7712!
7713!--             Check, if particle has moved to another grid cell.
7714                IF ( i /= ip  .OR.  j /= jp  .OR.  k /= kp )  THEN
7715!!
7716!--                If the particle stays on the same processor, the particle
7717!--                will be added to the particle array of the new processor.
7718                   number_of_particles = prt_count(k,j,i)
7719                   particles => grid_particles(k,j,i)%particles(1:number_of_particles)
7720
7721                   pindex = prt_count(k,j,i)+1
7722                   IF (  pindex > SIZE(grid_particles(k,j,i)%particles)  )     &
7723                   THEN
7724                      CALL realloc_particles_array( i, j, k )
7725                   ENDIF
7726
7727                   grid_particles(k,j,i)%particles(pindex) = particles_before_move(n)
7728                   prt_count(k,j,i) = pindex
7729
7730                   particles_before_move(n)%particle_mask = .FALSE.
7731                ENDIF
7732             ENDDO
7733
7734          ENDDO
7735       ENDDO
7736    ENDDO
7737
7738    CALL cpu_log( log_point_s(41), 'lpm_move_particle', 'stop' )
7739
7740    RETURN
7741
7742 END SUBROUTINE lpm_move_particle
7743 
7744
7745!------------------------------------------------------------------------------!
7746! Description:
7747! ------------
7748!> Check CFL-criterion for each particle. If one particle violated the
7749!> criterion the particle will be deleted and a warning message is given.
7750!------------------------------------------------------------------------------!
7751 SUBROUTINE lpm_check_cfl 
7752
7753    IMPLICIT NONE
7754
7755    INTEGER(iwp)  ::  i !< running index, x-direction
7756    INTEGER(iwp)  ::  j !< running index, y-direction
7757    INTEGER(iwp)  ::  k !< running index, z-direction
7758    INTEGER(iwp)  ::  n !< running index, number of particles
7759
7760    DO  i = nxl, nxr
7761       DO  j = nys, nyn
7762          DO  k = nzb+1, nzt
7763             number_of_particles = prt_count(k,j,i)
7764             IF ( number_of_particles <= 0 )  CYCLE
7765             particles => grid_particles(k,j,i)%particles(1:number_of_particles)         
7766             DO  n = 1, number_of_particles
7767!
7768!--             Note, check for CFL does not work at first particle timestep
7769!--             when both, age and age_m are zero.
7770                IF ( particles(n)%age - particles(n)%age_m > 0.0_wp )  THEN 
7771                   IF(ABS(particles(n)%speed_x) >                              &
7772                      (dx/(particles(n)%age-particles(n)%age_m))  .OR.         &
7773                      ABS(particles(n)%speed_y) >                              & 
7774                      (dy/(particles(n)%age-particles(n)%age_m))  .OR.         &
7775                      ABS(particles(n)%speed_z) >                              &
7776                      ((zw(k)-zw(k-1))/(particles(n)%age-particles(n)%age_m))) &
7777                   THEN
7778                      WRITE( message_string, * )                               &
7779                      'Particle violated CFL-criterion: &particle with id ',   &
7780                      particles(n)%id, ' will be deleted!'   
7781                      CALL message( 'lpm_check_cfl', 'PA0475', 0, 1, -1, 6, 0 )
7782                      particles(n)%particle_mask= .FALSE.
7783                   ENDIF
7784                ENDIF
7785             ENDDO
7786          ENDDO
7787       ENDDO
7788    ENDDO   
7789
7790 END SUBROUTINE lpm_check_cfl
7791 
7792 
7793!------------------------------------------------------------------------------!
7794! Description:
7795! ------------
7796!> If the allocated memory for the particle array do not suffice to add arriving
7797!> particles from neighbour grid cells, this subrouting reallocates the
7798!> particle array to assure enough memory is available.
7799!------------------------------------------------------------------------------!
7800 SUBROUTINE realloc_particles_array ( i, j, k, size_in )
7801
7802    INTEGER(iwp), INTENT(IN)                       ::  i              !<
7803    INTEGER(iwp), INTENT(IN)                       ::  j              !<
7804    INTEGER(iwp), INTENT(IN)                       ::  k              !<
7805    INTEGER(iwp), INTENT(IN), OPTIONAL             ::  size_in        !<
7806
7807    INTEGER(iwp)                                   ::  old_size        !<
7808    INTEGER(iwp)                                   ::  new_size        !<
7809    TYPE(particle_type), DIMENSION(:), ALLOCATABLE ::  tmp_particles_d !<
7810    TYPE(particle_type), DIMENSION(500)            ::  tmp_particles_s !<
7811
7812    old_size = SIZE(grid_particles(k,j,i)%particles)
7813
7814    IF ( PRESENT(size_in) )   THEN
7815       new_size = size_in
7816    ELSE
7817       new_size = old_size * ( 1.0_wp + alloc_factor / 100.0_wp )
7818    ENDIF
7819
7820    new_size = MAX( new_size, 1, old_size + 1 )
7821
7822    IF ( old_size <= 500 )  THEN
7823
7824       tmp_particles_s(1:old_size) = grid_particles(k,j,i)%particles(1:old_size)
7825
7826       DEALLOCATE(grid_particles(k,j,i)%particles)
7827       ALLOCATE(grid_particles(k,j,i)%particles(new_size))
7828
7829       grid_particles(k,j,i)%particles(1:old_size)          = tmp_particles_s(1:old_size)
7830       grid_particles(k,j,i)%particles(old_size+1:new_size) = zero_particle
7831
7832    ELSE
7833
7834       ALLOCATE(tmp_particles_d(new_size))
7835       tmp_particles_d(1:old_size) = grid_particles(k,j,i)%particles
7836
7837       DEALLOCATE(grid_particles(k,j,i)%particles)
7838       ALLOCATE(grid_particles(k,j,i)%particles(new_size))
7839
7840       grid_particles(k,j,i)%particles(1:old_size)          = tmp_particles_d(1:old_size)
7841       grid_particles(k,j,i)%particles(old_size+1:new_size) = zero_particle
7842
7843       DEALLOCATE(tmp_particles_d)
7844
7845    ENDIF
7846    particles => grid_particles(k,j,i)%particles(1:new_size)
7847
7848    RETURN
7849   
7850 END SUBROUTINE realloc_particles_array
7851 
7852 
7853!------------------------------------------------------------------------------!
7854! Description:
7855! ------------
7856!> Not needed but allocated space for particles is dealloced.
7857!------------------------------------------------------------------------------!
7858 SUBROUTINE dealloc_particles_array
7859
7860 
7861    INTEGER(iwp) ::  i               !<
7862    INTEGER(iwp) ::  j               !<
7863    INTEGER(iwp) ::  k               !<
7864    INTEGER(iwp) ::  old_size        !<
7865    INTEGER(iwp) ::  new_size        !<
7866
7867    LOGICAL ::  dealloc
7868
7869    TYPE(particle_type), DIMENSION(:), ALLOCATABLE ::  tmp_particles_d !<
7870    TYPE(particle_type), DIMENSION(500)            ::  tmp_particles_s !<
7871
7872    DO  i = nxl, nxr
7873       DO  j = nys, nyn
7874          DO  k = nzb+1, nzt
7875!
7876!--          Determine number of active particles
7877             number_of_particles = prt_count(k,j,i)
7878!
7879!--          Determine allocated memory size
7880             old_size = SIZE( grid_particles(k,j,i)%particles )
7881!
7882!--          Check for large unused memory
7883             dealloc = ( ( number_of_particles < 1 .AND.         &
7884                           old_size            > 1 )  .OR.       &
7885                         ( number_of_particles > 1 .AND.         &
7886                           old_size - number_of_particles *                    &
7887                              ( 1.0_wp + 0.01_wp * alloc_factor ) > 0.0_wp ) )
7888
7889             IF ( dealloc )  THEN
7890                IF ( number_of_particles < 1 )  THEN
7891                   new_size = 1
7892                ELSE
7893                   new_size = INT( number_of_particles * ( 1.0_wp + 0.01_wp * alloc_factor ) )
7894                ENDIF
7895
7896                IF ( number_of_particles <= 500 )  THEN
7897
7898                   tmp_particles_s(1:number_of_particles) = grid_particles(k,j,i)%particles(1:number_of_particles)
7899
7900                   DEALLOCATE(grid_particles(k,j,i)%particles)
7901                   ALLOCATE(grid_particles(k,j,i)%particles(new_size))
7902
7903                   grid_particles(k,j,i)%particles(1:number_of_particles)          = tmp_particles_s(1:number_of_particles)
7904                   grid_particles(k,j,i)%particles(number_of_particles+1:new_size) = zero_particle
7905
7906                ELSE
7907
7908                   ALLOCATE(tmp_particles_d(number_of_particles))
7909                   tmp_particles_d(1:number_of_particles) = grid_particles(k,j,i)%particles(1:number_of_particles)
7910
7911                   DEALLOCATE(grid_particles(k,j,i)%particles)
7912                   ALLOCATE(grid_particles(k,j,i)%particles(new_size))
7913
7914                   grid_particles(k,j,i)%particles(1:number_of_particles)          = tmp_particles_d(1:number_of_particles)
7915                   grid_particles(k,j,i)%particles(number_of_particles+1:new_size) = zero_particle
7916
7917                   DEALLOCATE(tmp_particles_d)
7918
7919                ENDIF
7920
7921             ENDIF
7922          ENDDO
7923       ENDDO
7924    ENDDO
7925
7926 END SUBROUTINE dealloc_particles_array 
7927 
7928 
7929!------------------------------------------------------------------------------!
7930! Description:
7931! -----------
7932!> Routine for the whole processor
7933!> Sort all particles into the 8 respective subgrid boxes (in case of trilinear
7934!> interpolation method) and free space of particles which has been marked for
7935!> deletion.
7936!------------------------------------------------------------------------------!
7937   SUBROUTINE lpm_sort_and_delete
7938
7939       INTEGER(iwp) ::  i  !<
7940       INTEGER(iwp) ::  ip !<
7941       INTEGER(iwp) ::  is !<
7942       INTEGER(iwp) ::  j  !<
7943       INTEGER(iwp) ::  jp !<
7944       INTEGER(iwp) ::  kp !<
7945       INTEGER(iwp) ::  m  !<
7946       INTEGER(iwp) ::  n  !<
7947       INTEGER(iwp) ::  nn !<
7948       INTEGER(iwp) ::  sort_index  !<
7949
7950       INTEGER(iwp), DIMENSION(0:7) ::  sort_count  !<
7951
7952       TYPE(particle_type), DIMENSION(:,:), ALLOCATABLE ::  sort_particles    !<
7953
7954       CALL cpu_log( log_point_s(51), 'lpm_sort_and_delete', 'start' )
7955       IF ( interpolation_trilinear )  THEN
7956          DO  ip = nxl, nxr
7957             DO  jp = nys, nyn
7958                DO  kp = nzb+1, nzt
7959                   number_of_particles = prt_count(kp,jp,ip)
7960                   IF ( number_of_particles <= 0 )  CYCLE
7961                   particles => grid_particles(kp,jp,ip)%particles(1:number_of_particles)
7962                   nn = 0
7963                   sort_count = 0
7964                   ALLOCATE( sort_particles(number_of_particles, 0:7) )
7965
7966                   DO  n = 1, number_of_particles
7967                      sort_index = 0
7968
7969                      IF ( particles(n)%particle_mask )  THEN
7970                         nn = nn + 1
7971!
7972!--                      Sorting particles with a binary scheme
7973!--                      sort_index=111_2=7_10 -> particle at the left,south,bottom subgridbox
7974!--                      sort_index=000_2=0_10 -> particle at the right,north,top subgridbox
7975!--                      For this the center of the gridbox is calculated
7976                         i = (particles(n)%x + 0.5_wp * dx) * ddx
7977                         j = (particles(n)%y + 0.5_wp * dy) * ddy
7978
7979                         IF ( i == ip )  sort_index = sort_index + 4
7980                         IF ( j == jp )  sort_index = sort_index + 2
7981                         IF ( zu(kp) > particles(n)%z ) sort_index = sort_index + 1
7982
7983                         sort_count(sort_index) = sort_count(sort_index) + 1
7984                         m = sort_count(sort_index)
7985                         sort_particles(m,sort_index) = particles(n)
7986                         sort_particles(m,sort_index)%block_nr = sort_index
7987                      ENDIF
7988                   ENDDO
7989!
7990!--                Delete and resort particles by overwritting and set
7991!--                the number_of_particles to the actual value.
7992                   nn = 0
7993                   DO  is = 0,7
7994                      grid_particles(kp,jp,ip)%start_index(is) = nn + 1
7995                      DO  n = 1,sort_count(is)
7996                         nn = nn + 1
7997                         particles(nn) = sort_particles(n,is)
7998                      ENDDO
7999                      grid_particles(kp,jp,ip)%end_index(is) = nn
8000                   ENDDO
8001
8002                   number_of_particles = nn
8003                   prt_count(kp,jp,ip) = number_of_particles
8004                   DEALLOCATE(sort_particles)
8005                ENDDO
8006             ENDDO
8007          ENDDO
8008
8009!--    In case of the simple interpolation method the particles must not
8010!--    be sorted in subboxes. Particles marked for deletion however, must be
8011!--    deleted and number of particles must be recalculated as it is also
8012!--    done for the trilinear particle advection interpolation method.
8013       ELSE
8014
8015          DO  ip = nxl, nxr
8016             DO  jp = nys, nyn
8017                DO  kp = nzb+1, nzt
8018
8019                   number_of_particles = prt_count(kp,jp,ip)
8020                   IF ( number_of_particles <= 0 )  CYCLE
8021                   particles => grid_particles(kp,jp,ip)%particles(1:number_of_particles)
8022!
8023!--                Repack particles array, i.e. delete particles and recalculate
8024!--                number of particles
8025                   CALL lpm_pack
8026                   prt_count(kp,jp,ip) = number_of_particles
8027                ENDDO
8028             ENDDO
8029          ENDDO
8030       ENDIF
8031       CALL cpu_log( log_point_s(51), 'lpm_sort_and_delete', 'stop' )
8032
8033    END SUBROUTINE lpm_sort_and_delete
8034
8035 
8036!------------------------------------------------------------------------------!
8037! Description:
8038! ------------
8039!> Move all particles not marked for deletion to lowest indices (packing)
8040!------------------------------------------------------------------------------!
8041    SUBROUTINE lpm_pack
8042
8043       INTEGER(iwp) ::  n       !<
8044       INTEGER(iwp) ::  nn      !<
8045!
8046!--    Find out elements marked for deletion and move data from highest index
8047!--    values to these free indices
8048       nn = number_of_particles
8049
8050       DO WHILE ( .NOT. particles(nn)%particle_mask )
8051          nn = nn-1
8052          IF ( nn == 0 )  EXIT
8053       ENDDO
8054
8055       IF ( nn > 0 )  THEN
8056          DO  n = 1, number_of_particles
8057             IF ( .NOT. particles(n)%particle_mask )  THEN
8058                particles(n) = particles(nn)
8059                nn = nn - 1
8060                DO WHILE ( .NOT. particles(nn)%particle_mask )
8061                   nn = nn-1
8062                   IF ( n == nn )  EXIT
8063                ENDDO
8064             ENDIF
8065             IF ( n == nn )  EXIT
8066          ENDDO
8067       ENDIF
8068
8069!
8070!--    The number of deleted particles has been determined in routines
8071!--    lpm_boundary_conds, lpm_droplet_collision, and lpm_exchange_horiz
8072       number_of_particles = nn
8073
8074    END SUBROUTINE lpm_pack 
8075
8076
8077!------------------------------------------------------------------------------!
8078! Description:
8079! ------------
8080!> Sort particles in each sub-grid box into two groups: particles that already
8081!> completed the LES timestep, and particles that need further timestepping to
8082!> complete the LES timestep.
8083!------------------------------------------------------------------------------!
8084    SUBROUTINE lpm_sort_timeloop_done
8085
8086       INTEGER(iwp) ::  end_index     !< particle end index for each sub-box
8087       INTEGER(iwp) ::  i             !< index of particle grid box in x-direction
8088       INTEGER(iwp) ::  j             !< index of particle grid box in y-direction
8089       INTEGER(iwp) ::  k             !< index of particle grid box in z-direction
8090       INTEGER(iwp) ::  n             !< running index for number of particles
8091       INTEGER(iwp) ::  nb            !< index of subgrid boux
8092       INTEGER(iwp) ::  nf            !< indices for particles in each sub-box that already finalized their substeps
8093       INTEGER(iwp) ::  nnf           !< indices for particles in each sub-box that need further treatment
8094       INTEGER(iwp) ::  num_finalized !< number of particles in each sub-box that already finalized their substeps
8095       INTEGER(iwp) ::  start_index   !< particle start index for each sub-box
8096
8097       TYPE(particle_type), DIMENSION(:), ALLOCATABLE ::  sort_particles  !< temporary particle array
8098
8099       DO  i = nxl, nxr
8100          DO  j = nys, nyn
8101             DO  k = nzb+1, nzt
8102
8103                number_of_particles = prt_count(k,j,i)
8104                IF ( number_of_particles <= 0 )  CYCLE
8105                particles => grid_particles(k,j,i)%particles(1:number_of_particles)
8106
8107                DO  nb = 0, 7
8108!
8109!--                Obtain start and end index for each subgrid box
8110                   start_index = grid_particles(k,j,i)%start_index(nb)
8111                   end_index   = grid_particles(k,j,i)%end_index(nb)
8112!
8113!--                Allocate temporary array used for sorting.
8114                   ALLOCATE( sort_particles(start_index:end_index) )
8115!
8116!--                Determine number of particles already completed the LES
8117!--                timestep, and write them into a temporary array.
8118                   nf = start_index
8119                   num_finalized = 0
8120                   DO  n = start_index, end_index
8121                      IF ( dt_3d - particles(n)%dt_sum < 1E-8_wp )  THEN
8122                         sort_particles(nf) = particles(n)
8123                         nf                 = nf + 1
8124                         num_finalized      = num_finalized + 1
8125                      ENDIF
8126                   ENDDO
8127!
8128!--                Determine number of particles that not completed the LES
8129!--                timestep, and write them into a temporary array.
8130                   nnf = nf
8131                   DO  n = start_index, end_index
8132                      IF ( dt_3d - particles(n)%dt_sum > 1E-8_wp )  THEN
8133                         sort_particles(nnf) = particles(n)
8134                         nnf                 = nnf + 1
8135                      ENDIF
8136                   ENDDO
8137!
8138!--                Write back sorted particles
8139                   particles(start_index:end_index) =                          &
8140                                           sort_particles(start_index:end_index)
8141!
8142!--                Determine updated start_index, used to masked already
8143!--                completed particles.
8144                   grid_particles(k,j,i)%start_index(nb) =                     &
8145                                      grid_particles(k,j,i)%start_index(nb)    &
8146                                    + num_finalized
8147!
8148!--                Deallocate dummy array
8149                   DEALLOCATE ( sort_particles )
8150!
8151!--                Finally, if number of non-completed particles is non zero
8152!--                in any of the sub-boxes, set control flag appropriately.
8153                   IF ( nnf > nf )                                             &
8154                      grid_particles(k,j,i)%time_loop_done = .FALSE.
8155
8156                ENDDO
8157             ENDDO
8158          ENDDO
8159       ENDDO
8160
8161    END SUBROUTINE lpm_sort_timeloop_done 
8162
8163END MODULE lagrangian_particle_model_mod
Note: See TracBrowser for help on using the repository browser.