source: palm/trunk/SOURCE/lpm_init.f90 @ 3573

Last change on this file since 3573 was 3560, checked in by raasch, 6 years ago

bugfix to guarantee correct particle releases in case that the release interval is smaller than the model timestep

  • Property svn:keywords set to Id
File size: 47.8 KB
RevLine 
[1873]1!> @file lpm_init.f90
[2000]2!------------------------------------------------------------------------------!
[2696]3! This file is part of the PALM model system.
[1036]4!
[2000]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.
[1036]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!
[2718]17! Copyright 1997-2018 Leibniz Universitaet Hannover
[2000]18!------------------------------------------------------------------------------!
[1036]19!
[254]20! Current revisions:
[1]21! -----------------
[2318]22!
[3049]23!
[1930]24! Former revisions:
25! -----------------
26! $Id: lpm_init.f90 3560 2018-11-23 09:20:21Z raasch $
[3560]27! set the first particle release time
28!
29! 3524 2018-11-14 13:36:44Z raasch
[3524]30! added missing working precision
31!
32! 3361 2018-10-16 20:39:37Z knoop
[3294]33! ocean renamed ocean_mode
34!
35! 3274 2018-09-24 15:42:55Z knoop
[3274]36! Modularization of all bulk cloud physics code components
37!
38! 3241 2018-09-12 15:02:00Z raasch
[3241]39! unused variables removed
40!
41! 3065 2018-06-12 07:03:02Z Giersch
[3065]42! dz was replaced by dzw or dz(1) to allow for right vertical stretching
43!
44! 3049 2018-05-29 13:52:36Z Giersch
[3049]45! Error messages revised
46!
47! 3045 2018-05-28 07:55:41Z Giersch
[3045]48! Error message revised
49!
50! 3039 2018-05-24 13:13:11Z schwenkel
[3039]51! bugfix for lcm with grid stretching
52!
53! 2967 2018-04-13 11:22:08Z raasch
[2967]54! nesting routine is only called if nesting is switched on
55!
56! 2954 2018-04-09 14:35:46Z schwenkel
[2954]57! Bugfix for particle initialization in case of ocean
58!
59! 2801 2018-02-14 16:01:55Z thiele
[2801]60! Introduce particle transfer in nested models.
61!
62! 2718 2018-01-02 08:49:38Z maronga
[2716]63! Corrected "Former revisions" section
[2701]64!
[2716]65! 2701 2017-12-15 15:40:50Z suehring
66! Changes from last commit documented
67!
[2701]68! 2698 2017-12-14 18:46:24Z suehring
[2716]69! Grid indices passed to lpm_boundary_conds. (responsible Philipp Thiele)
70!
71! 2696 2017-12-14 17:12:51Z kanani
72! Change in file header (GPL part)
73!
74! 2628 2017-11-20 12:40:38Z schwenkel
[2628]75! Enabled particle advection with grid stretching.
76!
77! 2608 2017-11-13 14:04:26Z schwenkel
[2608]78! Calculation of magnus equation in external module (diagnostic_quantities_mod).
79!
80! 2606 2017-11-10 10:36:31Z schwenkel
[2606]81! Changed particle box locations: center of particle box now coincides
82! with scalar grid point of same index.
83! Renamed module and subroutines: lpm_pack_arrays_mod -> lpm_pack_and_sort_mod
84! lpm_pack_all_arrays -> lpm_sort_in_subboxes, lpm_pack_arrays -> lpm_pack
85! lpm_sort -> lpm_sort_timeloop_done
86!
87! 2375 2017-08-29 14:10:28Z schwenkel
[2375]88! Initialization of chemical aerosol composition
89!
90! 2346 2017-08-09 16:39:17Z suehring
[2346]91! Bugfix, correct determination of topography top index
92!
93! 2318 2017-07-20 17:27:44Z suehring
[2318]94! Get topography top index via Function call
95!
96! 2317 2017-07-20 17:27:19Z suehring
[2312]97! Extended particle data type. Aerosol initialization improved.
98!
99! 2305 2017-07-06 11:18:47Z hoffmann
[2305]100! Improved calculation of particle IDs.
[2312]101!
[2305]102! 2274 2017-06-09 13:27:48Z Giersch
[2274]103!  Changed error messages
[2312]104!
[2274]105! 2265 2017-06-08 16:58:28Z schwenkel
[2265]106! Unused variables removed.
[2312]107!
[2265]108! 2263 2017-06-08 14:59:01Z schwenkel
[2263]109! Implemented splitting and merging algorithm
[2312]110!
[2263]111! 2233 2017-05-30 18:08:54Z suehring
[1930]112!
[2233]113! 2232 2017-05-30 17:47:52Z suehring
114! Adjustments according to new topography realization
[2312]115!
116!
[2224]117! 2223 2017-05-15 16:38:09Z suehring
118! Add check for particle release at model top
[2312]119!
[2183]120! 2182 2017-03-17 14:27:40Z schwenkel
121! Added parameters for simplified particle initialization.
[2305]122!
[2123]123! 2122 2017-01-18 12:22:54Z hoffmann
[2312]124! Improved initialization of equilibrium aerosol radii
[2123]125! Calculation of particle ID
126!
[2001]127! 2000 2016-08-20 18:09:15Z knoop
128! Forced header and separation lines into 80 columns
[2312]129!
[1930]130! 2016-06-09 16:25:25Z suehring
[2312]131! Bugfix in determining initial particle height and grid index in case of
[1929]132! seed_follows_topography.
[2312]133! Bugfix concerning random positions, ensure that particles do not move more
[1929]134! than one grid length.
135! Bugfix logarithmic interpolation.
136! Initial setting of sgs_wf_part.
[1321]137!
[1891]138! 1890 2016-04-22 08:52:11Z hoffmann
[2312]139! Initialization of aerosol equilibrium radius not possible in supersaturated
[1891]140! environments. Therefore, a maximum supersaturation of -1 % is assumed during
141! initialization.
142!
[1874]143! 1873 2016-04-18 14:50:06Z maronga
[1929]144! Module renamed (removed _mod
[2312]145!
[1872]146! 1871 2016-04-15 11:46:09Z hoffmann
147! Initialization of aerosols added.
148!
[1851]149! 1850 2016-04-08 13:29:27Z maronga
150! Module renamed
151!
[1832]152! 1831 2016-04-07 13:15:51Z hoffmann
153! curvature_solution_effects moved to particle_attributes
154!
[1823]155! 1822 2016-04-07 07:49:42Z hoffmann
156! Unused variables removed.
157!
[1784]158! 1783 2016-03-06 18:36:17Z raasch
159! netcdf module added
160!
[2312]161! 1725 2015-11-17 13:01:51Z hoffmann
162! Bugfix: Processor-dependent seed for random function is generated before it is
[1726]163! used.
164!
[1692]165! 1691 2015-10-26 16:17:44Z maronga
166! Renamed prandtl_layer to constant_flux_layer.
167!
[1686]168! 1685 2015-10-08 07:32:13Z raasch
169! bugfix concerning vertical index offset in case of ocean
170!
[1683]171! 1682 2015-10-07 23:56:08Z knoop
[2312]172! Code annotations made doxygen readable
[1683]173!
[1576]174! 1575 2015-03-27 09:56:27Z raasch
175! initial vertical particle position is allowed to follow the topography
176!
[1360]177! 1359 2014-04-11 17:15:14Z hoffmann
[2312]178! New particle structure integrated.
[1360]179! Kind definition added to all floating point numbers.
180! lpm_init changed form a subroutine to a module.
[2312]181!
[1329]182! 1327 2014-03-21 11:00:16Z raasch
183! -netcdf_output
184!
[1323]185! 1322 2014-03-20 16:38:49Z raasch
186! REAL functions provided with KIND-attribute
187!
[1321]188! 1320 2014-03-20 08:40:49Z raasch
[1320]189! ONLY-attribute added to USE-statements,
190! kind-parameters added to all INTEGER and REAL declaration statements,
191! kinds are defined in new module kinds,
192! revision history before 2012 removed,
193! comment fields (!:) to be used for variable explanations added to
194! all variable declaration statements
195! bugfix: #if defined( __parallel ) added
[850]196!
[1315]197! 1314 2014-03-14 18:25:17Z suehring
[2312]198! Vertical logarithmic interpolation of horizontal particle speed for particles
[1315]199! between roughness height and first vertical grid level.
200!
[1093]201! 1092 2013-02-02 11:24:22Z raasch
202! unused variables removed
203!
[1037]204! 1036 2012-10-22 13:43:42Z raasch
205! code put under GPL (PALM 3.9)
206!
[850]207! 849 2012-03-15 10:35:09Z raasch
[849]208! routine renamed: init_particles -> lpm_init
209! de_dx, de_dy, de_dz are allocated here (instead of automatic arrays in
210! advec_particles),
211! sort_particles renamed lpm_sort_arrays, user_init_particles renamed lpm_init
[392]212!
[829]213! 828 2012-02-21 12:00:36Z raasch
214! call of init_kernels, particle feature color renamed class
215!
[826]216! 824 2012-02-17 09:09:57Z raasch
217! particle attributes speed_x|y|z_sgs renamed rvar1|2|3,
218! array particles implemented as pointer
219!
[668]220! 667 2010-12-23 12:06:00Z suehring/gryschka
[2312]221! nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng for allocation
[668]222! of arrays.
223!
[1]224! Revision 1.1  1999/11/25 16:22:38  raasch
225! Initial revision
226!
227!
228! Description:
229! ------------
[1682]230!> This routine initializes a set of particles and their attributes (position,
231!> radius, ..) which are used by the Lagrangian particle model (see lpm).
[1]232!------------------------------------------------------------------------------!
[1682]233 MODULE lpm_init_mod
[2312]234
[2305]235    USE, INTRINSIC ::  ISO_C_BINDING
[1]236
[1320]237    USE arrays_3d,                                                             &
[3039]238        ONLY:  de_dx, de_dy, de_dz, dzw, zu, zw
[1320]239
240    USE control_parameters,                                                    &
[1691]241        ONLY:  cloud_droplets, constant_flux_layer, current_timestep_number,   &
[3294]242               dt_3d, dz, initializing_actions, message_string, ocean_mode,    &
[2312]243               simulated_time
[1320]244
245    USE grid_variables,                                                        &
[1359]246        ONLY:  ddx, dx, ddy, dy
[1320]247
248    USE indices,                                                               &
[1575]249        ONLY:  nx, nxl, nxlg, nxrg, nxr, ny, nyn, nys, nyng, nysg, nz, nzb,    &
[2317]250               nzt, wall_flags_0
[1320]251
252    USE kinds
253
254    USE lpm_collision_kernels_mod,                                             &
255        ONLY:  init_kernels
256
[1783]257    USE netcdf_interface,                                                      &
258        ONLY:  netcdf_data_format
259
[1320]260    USE particle_attributes,                                                   &
[1359]261        ONLY:   alloc_factor, bc_par_b, bc_par_lr, bc_par_ns, bc_par_t,        &
262                block_offset, block_offset_def, collision_kernel,              &
[2265]263                curvature_solution_effects, density_ratio, grid_particles,     &
264                isf,i_splitting_mode, initial_weighting_factor, ibc_par_b,     &
[3560]265                ibc_par_lr, ibc_par_ns, ibc_par_t, iran_part,                  &
266                last_particle_release_time, log_z_z0,                          &
[2265]267                max_number_of_particle_groups, min_nr_particle,                &
[2305]268                number_concentration,                                          &
[2265]269                number_particles_per_gridbox,  number_of_particles,            &
[1320]270                number_of_particle_groups, number_of_sublayers,                &
[1822]271                offset_ocean_nzt, offset_ocean_nzt_m1,                         &
[1359]272                particles, particle_advection_start, particle_groups,          &
273                particle_groups_type, particles_per_point,                     &
[2312]274                particle_type, pdx, pdy, pdz,  prt_count, psb, psl, psn, psr,  &
275                pss, pst, radius, random_start_position,                       &
[2265]276                read_particles_from_restartfile, seed_follows_topography,      &
277                sgs_wf_part, sort_count, splitting_function, splitting_mode,   &
278                total_number_of_particles, use_sgs_for_particles,              &
279                write_particle_statistics, zero_particle, z0_av_global
[1320]280
[1]281    USE pegrid
282
[1320]283    USE random_function_mod,                                                   &
284        ONLY:  random_function
[1]285
[2232]286    USE surface_mod,                                                           &
[2698]287        ONLY:  get_topography_top_index_ji, surf_def_h, surf_lsm_h, surf_usm_h
[2232]288
[2801]289    USE pmc_particle_interface,                                                &
290        ONLY:  pmcp_g_init
291
[1359]292    IMPLICIT NONE
[1320]293
[1359]294    PRIVATE
295
[1682]296    INTEGER(iwp), PARAMETER         :: PHASE_INIT    = 1  !<
297    INTEGER(iwp), PARAMETER, PUBLIC :: PHASE_RELEASE = 2  !<
[1359]298
299    INTERFACE lpm_init
300       MODULE PROCEDURE lpm_init
301    END INTERFACE lpm_init
302
303    INTERFACE lpm_create_particle
304       MODULE PROCEDURE lpm_create_particle
305    END INTERFACE lpm_create_particle
306
307    PUBLIC lpm_init, lpm_create_particle
308
[1929]309 CONTAINS
[1359]310
[1682]311!------------------------------------------------------------------------------!
312! Description:
313! ------------
314!> @todo Missing subroutine description.
315!------------------------------------------------------------------------------!
[1359]316 SUBROUTINE lpm_init
317
318    USE lpm_collision_kernels_mod,                                             &
319        ONLY:  init_kernels
320
[2967]321    USE pmc_interface,                                                         &
322        ONLY: nested_run
323
[1]324    IMPLICIT NONE
325
[1682]326    INTEGER(iwp) ::  i                           !<
327    INTEGER(iwp) ::  j                           !<
328    INTEGER(iwp) ::  k                           !<
[1320]329
[2312]330    REAL(wp) ::  div                             !<
[1682]331    REAL(wp) ::  height_int                      !<
332    REAL(wp) ::  height_p                        !<
333    REAL(wp) ::  z_p                             !<
334    REAL(wp) ::  z0_av_local                     !<
[1]335
[1359]336
[1]337!
[150]338!-- In case of oceans runs, the vertical index calculations need an offset,
339!-- because otherwise the k indices will become negative
[3294]340    IF ( ocean_mode )  THEN
[150]341       offset_ocean_nzt    = nzt
342       offset_ocean_nzt_m1 = nzt - 1
343    ENDIF
344
[1359]345!
346!-- Define block offsets for dividing a gridcell in 8 sub cells
[2606]347!-- See documentation for List of subgrid boxes
348!-- See pack_and_sort in lpm_pack_arrays.f90 for assignment of the subgrid boxes
349    block_offset(0) = block_offset_def ( 0, 0, 0)
350    block_offset(1) = block_offset_def ( 0, 0,-1)
351    block_offset(2) = block_offset_def ( 0,-1, 0)
352    block_offset(3) = block_offset_def ( 0,-1,-1)
353    block_offset(4) = block_offset_def (-1, 0, 0)
354    block_offset(5) = block_offset_def (-1, 0,-1)
355    block_offset(6) = block_offset_def (-1,-1, 0)
356    block_offset(7) = block_offset_def (-1,-1,-1)
[150]357!
[1]358!-- Check the number of particle groups.
359    IF ( number_of_particle_groups > max_number_of_particle_groups )  THEN
[3045]360       WRITE( message_string, * ) 'max_number_of_particle_groups =',           &
361                                  max_number_of_particle_groups ,              &
[3046]362                                  '&number_of_particle_groups reset to ',      &
[254]363                                  max_number_of_particle_groups
[849]364       CALL message( 'lpm_init', 'PA0213', 0, 1, 0, 6, 0 )
[1]365       number_of_particle_groups = max_number_of_particle_groups
366    ENDIF
[2232]367!
368!-- Check if downward-facing walls exist. This case, reflection boundary
369!-- conditions (as well as subgrid-scale velocities) may do not work
[2312]370!-- propably (not realized so far).
[2232]371    IF ( surf_def_h(1)%ns >= 1 )  THEN
[3045]372       WRITE( message_string, * ) 'Overhanging topography do not work '//      &
[2312]373                                  'with particles'
[2232]374       CALL message( 'lpm_init', 'PA0212', 0, 1, 0, 6, 0 )
[1]375
[2312]376    ENDIF
[2232]377
[1]378!
379!-- Set default start positions, if necessary
[2606]380    IF ( psl(1) == 9999999.9_wp )  psl(1) = 0.0_wp
381    IF ( psr(1) == 9999999.9_wp )  psr(1) = ( nx +1 ) * dx
382    IF ( pss(1) == 9999999.9_wp )  pss(1) = 0.0_wp
383    IF ( psn(1) == 9999999.9_wp )  psn(1) = ( ny +1 ) * dy
[1359]384    IF ( psb(1) == 9999999.9_wp )  psb(1) = zu(nz/2)
385    IF ( pst(1) == 9999999.9_wp )  pst(1) = psb(1)
[1]386
[1359]387    IF ( pdx(1) == 9999999.9_wp  .OR.  pdx(1) == 0.0_wp )  pdx(1) = dx
388    IF ( pdy(1) == 9999999.9_wp  .OR.  pdy(1) == 0.0_wp )  pdy(1) = dy
389    IF ( pdz(1) == 9999999.9_wp  .OR.  pdz(1) == 0.0_wp )  pdz(1) = zu(2) - zu(1)
[1]390
[2182]391!
392!-- If number_particles_per_gridbox is set, the parametres pdx, pdy and pdz are
393!-- calculated diagnostically. Therfore an isotropic distribution is prescribed.
[2312]394    IF ( number_particles_per_gridbox /= -1 .AND.   &
[2182]395         number_particles_per_gridbox >= 1 )    THEN
[2312]396       pdx(1) = (( dx * dy * ( zu(2) - zu(1) ) ) /  &
[2182]397             REAL(number_particles_per_gridbox))**0.3333333_wp
398!
[2312]399!--    Ensure a smooth value (two significant digits) of distance between
[2182]400!--    particles (pdx, pdy, pdz).
401       div = 1000.0_wp
402       DO  WHILE ( pdx(1) < div )
403          div = div / 10.0_wp
404       ENDDO
405       pdx(1) = NINT( pdx(1) * 100.0_wp / div ) * div / 100.0_wp
406       pdy(1) = pdx(1)
407       pdz(1) = pdx(1)
408
409    ENDIF
410
[1]411    DO  j = 2, number_of_particle_groups
[1359]412       IF ( psl(j) == 9999999.9_wp )  psl(j) = psl(j-1)
413       IF ( psr(j) == 9999999.9_wp )  psr(j) = psr(j-1)
414       IF ( pss(j) == 9999999.9_wp )  pss(j) = pss(j-1)
415       IF ( psn(j) == 9999999.9_wp )  psn(j) = psn(j-1)
416       IF ( psb(j) == 9999999.9_wp )  psb(j) = psb(j-1)
417       IF ( pst(j) == 9999999.9_wp )  pst(j) = pst(j-1)
418       IF ( pdx(j) == 9999999.9_wp  .OR.  pdx(j) == 0.0_wp )  pdx(j) = pdx(j-1)
419       IF ( pdy(j) == 9999999.9_wp  .OR.  pdy(j) == 0.0_wp )  pdy(j) = pdy(j-1)
420       IF ( pdz(j) == 9999999.9_wp  .OR.  pdz(j) == 0.0_wp )  pdz(j) = pdz(j-1)
[1]421    ENDDO
422
423!
[2312]424!-- Allocate arrays required for calculating particle SGS velocities.
[1929]425!-- Initialize prefactor required for stoachastic Weil equation.
[1822]426    IF ( use_sgs_for_particles  .AND.  .NOT. cloud_droplets )  THEN
[849]427       ALLOCATE( de_dx(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
428                 de_dy(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
429                 de_dz(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
[1929]430
[2312]431       sgs_wf_part = 1.0_wp / 3.0_wp
[849]432    ENDIF
433
434!
[2312]435!-- Allocate array required for logarithmic vertical interpolation of
[1314]436!-- horizontal particle velocities between the surface and the first vertical
[2312]437!-- grid level. In order to avoid repeated CPU cost-intensive CALLS of
438!-- intrinsic FORTRAN procedure LOG(z/z0), LOG(z/z0) is precalculated for
439!-- several heights. Splitting into 20 sublayers turned out to be sufficient.
[1314]440!-- To obtain exact height levels of particles, linear interpolation is applied
441!-- (see lpm_advec.f90).
[1691]442    IF ( constant_flux_layer )  THEN
[2312]443
444       ALLOCATE ( log_z_z0(0:number_of_sublayers) )
[2232]445       z_p = zu(nzb+1) - zw(nzb)
[1314]446
447!
448!--    Calculate horizontal mean value of z0 used for logartihmic
[2312]449!--    interpolation. Note: this is not exact for heterogeneous z0.
450!--    However, sensitivity studies showed that the effect is
451!--    negligible.
[2232]452       z0_av_local  = SUM( surf_def_h(0)%z0 ) + SUM( surf_lsm_h%z0 ) +         &
453                      SUM( surf_usm_h%z0 )
[1359]454       z0_av_global = 0.0_wp
[1314]455
[1320]456#if defined( __parallel )
[1314]457       CALL MPI_ALLREDUCE(z0_av_local, z0_av_global, 1, MPI_REAL, MPI_SUM, &
458                          comm2d, ierr )
[1320]459#else
460       z0_av_global = z0_av_local
461#endif
[1314]462
463       z0_av_global = z0_av_global  / ( ( ny + 1 ) * ( nx + 1 ) )
464!
465!--    Horizontal wind speed is zero below and at z0
[1359]466       log_z_z0(0) = 0.0_wp
[1314]467!
468!--    Calculate vertical depth of the sublayers
[1322]469       height_int  = ( z_p - z0_av_global ) / REAL( number_of_sublayers, KIND=wp )
[1314]470!
471!--    Precalculate LOG(z/z0)
[1929]472       height_p    = z0_av_global
[1314]473       DO  k = 1, number_of_sublayers
474
475          height_p    = height_p + height_int
476          log_z_z0(k) = LOG( height_p / z0_av_global )
477
478       ENDDO
479
480    ENDIF
481
482!
[1359]483!-- Check boundary condition and set internal variables
484    SELECT CASE ( bc_par_b )
[2312]485
[1359]486       CASE ( 'absorb' )
487          ibc_par_b = 1
488
489       CASE ( 'reflect' )
490          ibc_par_b = 2
[2312]491
[1359]492       CASE DEFAULT
[3045]493          WRITE( message_string, * )  'unknown boundary condition ',           &
[1359]494                                       'bc_par_b = "', TRIM( bc_par_b ), '"'
495          CALL message( 'lpm_init', 'PA0217', 1, 2, 0, 6, 0 )
[2312]496
[1359]497    END SELECT
498    SELECT CASE ( bc_par_t )
[2312]499
[1359]500       CASE ( 'absorb' )
501          ibc_par_t = 1
502
503       CASE ( 'reflect' )
504          ibc_par_t = 2
[2801]505         
506       CASE ( 'nested' )
507          ibc_par_t = 3
[2312]508
[1359]509       CASE DEFAULT
[3045]510          WRITE( message_string, * ) 'unknown boundary condition ',            &
[1359]511                                     'bc_par_t = "', TRIM( bc_par_t ), '"'
512          CALL message( 'lpm_init', 'PA0218', 1, 2, 0, 6, 0 )
[2312]513
[1359]514    END SELECT
515    SELECT CASE ( bc_par_lr )
516
517       CASE ( 'cyclic' )
518          ibc_par_lr = 0
519
520       CASE ( 'absorb' )
521          ibc_par_lr = 1
522
523       CASE ( 'reflect' )
524          ibc_par_lr = 2
[2801]525         
526       CASE ( 'nested' )
527          ibc_par_lr = 3
[2312]528
[1359]529       CASE DEFAULT
530          WRITE( message_string, * ) 'unknown boundary condition ',   &
531                                     'bc_par_lr = "', TRIM( bc_par_lr ), '"'
532          CALL message( 'lpm_init', 'PA0219', 1, 2, 0, 6, 0 )
[2312]533
[1359]534    END SELECT
535    SELECT CASE ( bc_par_ns )
536
537       CASE ( 'cyclic' )
538          ibc_par_ns = 0
539
540       CASE ( 'absorb' )
541          ibc_par_ns = 1
542
543       CASE ( 'reflect' )
544          ibc_par_ns = 2
[2801]545         
546       CASE ( 'nested' )
547          ibc_par_ns = 3
[2312]548
[1359]549       CASE DEFAULT
550          WRITE( message_string, * ) 'unknown boundary condition ',   &
551                                     'bc_par_ns = "', TRIM( bc_par_ns ), '"'
552          CALL message( 'lpm_init', 'PA0220', 1, 2, 0, 6, 0 )
[2312]553
[1359]554    END SELECT
[2263]555    SELECT CASE ( splitting_mode )
[2312]556
[2263]557       CASE ( 'const' )
558          i_splitting_mode = 1
[1359]559
[2263]560       CASE ( 'cl_av' )
561          i_splitting_mode = 2
562
563       CASE ( 'gb_av' )
564          i_splitting_mode = 3
[2312]565
[2263]566       CASE DEFAULT
[3045]567          WRITE( message_string, * )  'unknown splitting_mode = "',            &
568                                      TRIM( splitting_mode ), '"'
[2263]569          CALL message( 'lpm_init', 'PA0146', 1, 2, 0, 6, 0 )
[2312]570
[2263]571    END SELECT
572    SELECT CASE ( splitting_function )
[2312]573
[2263]574       CASE ( 'gamma' )
575          isf = 1
576
577       CASE ( 'log' )
578          isf = 2
579
580       CASE ( 'exp' )
581          isf = 3
[2312]582
[2263]583       CASE DEFAULT
[3045]584          WRITE( message_string, * )  'unknown splitting function = "',        &
585                                       TRIM( splitting_function ), '"'
[2263]586          CALL message( 'lpm_init', 'PA0147', 1, 2, 0, 6, 0 )
[2312]587
[2263]588    END SELECT
589
[2312]590
[1359]591!
[828]592!-- Initialize collision kernels
593    IF ( collision_kernel /= 'none' )  CALL init_kernels
594
595!
[1]596!-- For the first model run of a possible job chain initialize the
[849]597!-- particles, otherwise read the particle data from restart file.
[1]598    IF ( TRIM( initializing_actions ) == 'read_restart_data'  &
599         .AND.  read_particles_from_restartfile )  THEN
600
[849]601       CALL lpm_read_restart_file
[1]602
603    ELSE
604
605!
606!--    Allocate particle arrays and set attributes of the initial set of
607!--    particles, which can be also periodically released at later times.
[1359]608       ALLOCATE( prt_count(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
609                 grid_particles(nzb+1:nzt,nys:nyn,nxl:nxr) )
[1]610
[1359]611       number_of_particles         = 0
[792]612
613       sort_count = 0
[1359]614       prt_count  = 0
[792]615
[1]616!
[2312]617!--    initialize counter for particle IDs
[2305]618       grid_particles%id_counter = 1
[2122]619
620!
[1]621!--    Initialize all particles with dummy values (otherwise errors may
622!--    occur within restart runs). The reason for this is still not clear
623!--    and may be presumably caused by errors in the respective user-interface.
[1359]624       zero_particle = particle_type( 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp,  &
625                                      0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp,  &
626                                      0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp,  &
[2312]627                                      0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp,  &
[2305]628                                      0, 0, 0_idp, .FALSE., -1 )
[1822]629
[1359]630       particle_groups = particle_groups_type( 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp )
[1]631
632!
633!--    Set values for the density ratio and radius for all particle
634!--    groups, if necessary
[1359]635       IF ( density_ratio(1) == 9999999.9_wp )  density_ratio(1) = 0.0_wp
636       IF ( radius(1)        == 9999999.9_wp )  radius(1) = 0.0_wp
[1]637       DO  i = 2, number_of_particle_groups
[1359]638          IF ( density_ratio(i) == 9999999.9_wp )  THEN
[1]639             density_ratio(i) = density_ratio(i-1)
640          ENDIF
[1359]641          IF ( radius(i) == 9999999.9_wp )  radius(i) = radius(i-1)
[1]642       ENDDO
643
644       DO  i = 1, number_of_particle_groups
[1359]645          IF ( density_ratio(i) /= 0.0_wp  .AND.  radius(i) == 0 )  THEN
[3045]646             WRITE( message_string, * ) 'particle group #', i, ' has a',       &
[254]647                                        'density ratio /= 0 but radius = 0'
[849]648             CALL message( 'lpm_init', 'PA0215', 1, 2, 0, 6, 0 )
[1]649          ENDIF
650          particle_groups(i)%density_ratio = density_ratio(i)
651          particle_groups(i)%radius        = radius(i)
652       ENDDO
653
654!
[1359]655!--    Set a seed value for the random number generator to be exclusively
656!--    used for the particle code. The generated random numbers should be
657!--    different on the different PEs.
658       iran_part = iran_part + myid
659
660!
[3560]661!--    Create the particle set, and set the initial particles
662       CALL lpm_create_particle( phase_init )
663       last_particle_release_time = particle_advection_start
664
665!
[1359]666!--    User modification of initial particles
667       CALL user_lpm_init
668
669!
670!--    Open file for statistical informations about particle conditions
671       IF ( write_particle_statistics )  THEN
672          CALL check_open( 80 )
673          WRITE ( 80, 8000 )  current_timestep_number, simulated_time,         &
[2312]674                              number_of_particles
[1359]675          CALL close_file( 80 )
676       ENDIF
677
678    ENDIF
679
[2967]680    IF ( nested_run )  CALL pmcp_g_init
[2801]681
[1359]682!
[2312]683!-- To avoid programm abort, assign particles array to the local version of
[1359]684!-- first grid cell
685    number_of_particles = prt_count(nzb+1,nys,nxl)
686    particles => grid_particles(nzb+1,nys,nxl)%particles(1:number_of_particles)
687!
688!-- Formats
6898000 FORMAT (I6,1X,F7.2,4X,I10,71X,I10)
690
691 END SUBROUTINE lpm_init
692
[1682]693!------------------------------------------------------------------------------!
694! Description:
695! ------------
696!> @todo Missing subroutine description.
697!------------------------------------------------------------------------------!
[1359]698 SUBROUTINE lpm_create_particle (phase)
[2628]699   
700    USE arrays_3d,                                                             &
701       ONLY:  zw
[1359]702    USE lpm_exchange_horiz_mod,                                                &
703        ONLY: lpm_exchange_horiz, lpm_move_particle, realloc_particles_array
704
[2628]705    USE lpm_pack_and_sort_mod,                                                 &
[2606]706        ONLY: lpm_sort_in_subboxes
[1359]707
[1871]708    USE particle_attributes,                                                   &
[2312]709        ONLY: deleted_particles
[1871]710
[1359]711    IMPLICIT  NONE
712
[1929]713    INTEGER(iwp)               ::  alloc_size  !< relative increase of allocated memory for particles
714    INTEGER(iwp)               ::  i           !< loop variable ( particle groups )
715    INTEGER(iwp)               ::  ip          !< index variable along x
716    INTEGER(iwp)               ::  j           !< loop variable ( particles per point )
717    INTEGER(iwp)               ::  jp          !< index variable along y
[2232]718    INTEGER(iwp)               ::  k           !< index variable along z
719    INTEGER(iwp)               ::  k_surf      !< index of surface grid point
[1929]720    INTEGER(iwp)               ::  kp          !< index variable along z
721    INTEGER(iwp)               ::  loop_stride !< loop variable for initialization
722    INTEGER(iwp)               ::  n           !< loop variable ( number of particles )
723    INTEGER(iwp)               ::  new_size    !< new size of allocated memory for particles
[1359]724
[1929]725    INTEGER(iwp), INTENT(IN)   ::  phase       !< mode of inititialization
[1359]726
[1929]727    INTEGER(iwp), DIMENSION(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ::  local_count !< start address of new particle
728    INTEGER(iwp), DIMENSION(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ::  local_start !< start address of new particle
[1359]729
[1929]730    LOGICAL                    ::  first_stride !< flag for initialization
[1359]731
[2312]732    REAL(wp)                   ::  pos_x      !< increment for particle position in x
733    REAL(wp)                   ::  pos_y      !< increment for particle position in y
734    REAL(wp)                   ::  pos_z      !< increment for particle position in z
[1929]735    REAL(wp)                   ::  rand_contr !< dummy argument for random position
[1359]736
[1929]737    TYPE(particle_type),TARGET ::  tmp_particle !< temporary particle used for initialization
[1359]738
739!
740!-- Calculate particle positions and store particle attributes, if
741!-- particle is situated on this PE
742    DO  loop_stride = 1, 2
743       first_stride = (loop_stride == 1)
744       IF ( first_stride )   THEN
745          local_count = 0           ! count number of particles
746       ELSE
747          local_count = prt_count   ! Start address of new particles
748       ENDIF
749
[2182]750!
751!--    Calculate initial_weighting_factor diagnostically
752       IF ( number_concentration /= -1.0_wp .AND. number_concentration > 0.0_wp ) THEN
[2312]753          initial_weighting_factor =  number_concentration * 1.0E6_wp *             &
754                                      pdx(1) * pdy(1) * pdz(1)
[2182]755       END IF
756
[1]757       n = 0
758       DO  i = 1, number_of_particle_groups
759
760          pos_z = psb(i)
761
762          DO WHILE ( pos_z <= pst(i) )
763
[2954]764             IF ( pos_z >= zw(0) .AND.  pos_z < zw(nzt) )  THEN
[1]765
766
[2223]767                pos_y = pss(i)
[1]768
[2223]769                DO WHILE ( pos_y <= psn(i) )
[1]770
[2606]771                   IF ( pos_y >= nys * dy  .AND.                  &
772                        pos_y <  ( nyn + 1 ) * dy  ) THEN
[1]773
[2223]774                      pos_x = psl(i)
[1]775
[2223]776               xloop: DO WHILE ( pos_x <= psr(i) )
[1]777
[2606]778                         IF ( pos_x >= nxl * dx  .AND.            &
779                              pos_x <  ( nxr + 1) * dx ) THEN
[2223]780
781                            DO  j = 1, particles_per_point
782
[2305]783
[2223]784                               n = n + 1
785                               tmp_particle%x             = pos_x
786                               tmp_particle%y             = pos_y
787                               tmp_particle%z             = pos_z
788                               tmp_particle%age           = 0.0_wp
789                               tmp_particle%age_m         = 0.0_wp
790                               tmp_particle%dt_sum        = 0.0_wp
791                               tmp_particle%e_m           = 0.0_wp
[2312]792                               tmp_particle%rvar1         = 0.0_wp
793                               tmp_particle%rvar2         = 0.0_wp
794                               tmp_particle%rvar3         = 0.0_wp
[2223]795                               tmp_particle%speed_x       = 0.0_wp
796                               tmp_particle%speed_y       = 0.0_wp
797                               tmp_particle%speed_z       = 0.0_wp
798                               tmp_particle%origin_x      = pos_x
799                               tmp_particle%origin_y      = pos_y
800                               tmp_particle%origin_z      = pos_z
[2312]801                               IF ( curvature_solution_effects )  THEN
802                                  tmp_particle%aux1      = 0.0_wp    ! dry aerosol radius
803                                  tmp_particle%aux2      = dt_3d     ! last Rosenbrock timestep
804                               ELSE
805                                  tmp_particle%aux1      = 0.0_wp    ! free to use
806                                  tmp_particle%aux2      = 0.0_wp    ! free to use
807                               ENDIF
[2223]808                               tmp_particle%radius        = particle_groups(i)%radius
809                               tmp_particle%weight_factor = initial_weighting_factor
810                               tmp_particle%class         = 1
811                               tmp_particle%group         = i
[2305]812                               tmp_particle%id            = 0_idp
[2223]813                               tmp_particle%particle_mask = .TRUE.
814                               tmp_particle%block_nr      = -1
[1]815!
[2223]816!--                            Determine the grid indices of the particle position
[3524]817                               ip = INT( tmp_particle%x * ddx )
818                               jp = INT( tmp_particle%y * ddy )
819                               kp = INT( tmp_particle%z / dz(1) + 1 + offset_ocean_nzt )
[2628]820                               DO WHILE( zw(kp) < tmp_particle%z ) 
821                                  kp = kp + 1
822                               ENDDO
823                               DO WHILE( zw(kp-1) > tmp_particle%z )
824                                  kp = kp - 1
825                               ENDDO 
[2232]826!
827!--                            Determine surface level. Therefore, check for
[2317]828!--                            upward-facing wall on w-grid.
[2698]829                               k_surf = get_topography_top_index_ji( jp, ip, 'w' )
[1]830
[2223]831                               IF ( seed_follows_topography )  THEN
[1575]832!
[2223]833!--                               Particle height is given relative to topography
[2232]834                                  kp = kp + k_surf
835                                  tmp_particle%z = tmp_particle%z + zw(k_surf)
836!--                               Skip particle release if particle position is
837!--                               above model top, or within topography in case
838!--                               of overhanging structures.
839                                  IF ( kp > nzt  .OR.                          &
840                                 .NOT. BTEST( wall_flags_0(kp,jp,ip), 0 ) )  THEN
[2223]841                                     pos_x = pos_x + pdx(i)
842                                     CYCLE xloop
843                                  ENDIF
[2232]844!
[2312]845!--                            Skip particle release if particle position is
[2232]846!--                            below surface, or within topography in case
847!--                            of overhanging structures.
[2223]848                               ELSEIF ( .NOT. seed_follows_topography .AND.    &
[2232]849                                         tmp_particle%z <= zw(k_surf)  .OR.    &
850                                        .NOT. BTEST( wall_flags_0(kp,jp,ip), 0 ) )&
851                               THEN
[1575]852                                  pos_x = pos_x + pdx(i)
[2312]853                                  CYCLE xloop
[1575]854                               ENDIF
855
[2223]856                               local_count(kp,jp,ip) = local_count(kp,jp,ip) + 1
[2182]857
[2223]858                               IF ( .NOT. first_stride )  THEN
859                                  IF ( ip < nxl  .OR.  jp < nys  .OR.  kp < nzb+1 )  THEN
860                                     write(6,*) 'xl ',ip,jp,kp,nxl,nys,nzb+1
861                                  ENDIF
862                                  IF ( ip > nxr  .OR.  jp > nyn  .OR.  kp > nzt )  THEN
863                                     write(6,*) 'xu ',ip,jp,kp,nxr,nyn,nzt
864                                  ENDIF
865                                  grid_particles(kp,jp,ip)%particles(local_count(kp,jp,ip)) = tmp_particle
866
[1359]867                               ENDIF
[2223]868                            ENDDO
[1929]869
[2223]870                         ENDIF
[1]871
[2223]872                         pos_x = pos_x + pdx(i)
[1]873
[2223]874                      ENDDO xloop
[1]875
[2223]876                   ENDIF
[1]877
[2223]878                   pos_y = pos_y + pdy(i)
[1]879
[2223]880                ENDDO
[1]881
[2223]882             ENDIF
[1]883
884             pos_z = pos_z + pdz(i)
885
886          ENDDO
887
888       ENDDO
889
[1359]890       IF ( first_stride )  THEN
891          DO  ip = nxl, nxr
892             DO  jp = nys, nyn
893                DO  kp = nzb+1, nzt
894                   IF ( phase == PHASE_INIT )  THEN
895                      IF ( local_count(kp,jp,ip) > 0 )  THEN
896                         alloc_size = MAX( INT( local_count(kp,jp,ip) *        &
897                            ( 1.0_wp + alloc_factor / 100.0_wp ) ),            &
898                            min_nr_particle )
899                      ELSE
900                         alloc_size = min_nr_particle
901                      ENDIF
902                      ALLOCATE(grid_particles(kp,jp,ip)%particles(1:alloc_size))
903                      DO  n = 1, alloc_size
904                         grid_particles(kp,jp,ip)%particles(n) = zero_particle
905                      ENDDO
906                   ELSEIF ( phase == PHASE_RELEASE )  THEN
907                      IF ( local_count(kp,jp,ip) > 0 )  THEN
908                         new_size   = local_count(kp,jp,ip) + prt_count(kp,jp,ip)
909                         alloc_size = MAX( INT( new_size * ( 1.0_wp +          &
910                            alloc_factor / 100.0_wp ) ), min_nr_particle )
911                         IF( alloc_size > SIZE( grid_particles(kp,jp,ip)%particles) )  THEN
[2182]912                            CALL realloc_particles_array(ip,jp,kp,alloc_size)
[1359]913                         ENDIF
914                      ENDIF
915                   ENDIF
916                ENDDO
917             ENDDO
918          ENDDO
919       ENDIF
[1929]920
[1359]921    ENDDO
[1]922
[2182]923
924
[1359]925    local_start = prt_count+1
926    prt_count   = local_count
[1871]927
[1]928!
[2122]929!-- Calculate particle IDs
930    DO  ip = nxl, nxr
931       DO  jp = nys, nyn
932          DO  kp = nzb+1, nzt
933             number_of_particles = prt_count(kp,jp,ip)
934             IF ( number_of_particles <= 0 )  CYCLE
935             particles => grid_particles(kp,jp,ip)%particles(1:number_of_particles)
936
[2312]937             DO  n = local_start(kp,jp,ip), number_of_particles  !only new particles
[2122]938
[2305]939                particles(n)%id = 10000_idp**3 * grid_particles(kp,jp,ip)%id_counter + &
940                                  10000_idp**2 * kp + 10000_idp * jp + ip
941!
942!--             Count the number of particles that have been released before
[2122]943                grid_particles(kp,jp,ip)%id_counter =                          &
944                                         grid_particles(kp,jp,ip)%id_counter + 1
945
946             ENDDO
947
948          ENDDO
949       ENDDO
950    ENDDO
951
952!
[1871]953!-- Initialize aerosol background spectrum
[2312]954    IF ( curvature_solution_effects )  THEN
[1871]955       CALL lpm_init_aerosols(local_start)
956    ENDIF
957
958!
[1929]959!-- Add random fluctuation to particle positions.
[1359]960    IF ( random_start_position )  THEN
961       DO  ip = nxl, nxr
962          DO  jp = nys, nyn
963             DO  kp = nzb+1, nzt
964                number_of_particles = prt_count(kp,jp,ip)
965                IF ( number_of_particles <= 0 )  CYCLE
966                particles => grid_particles(kp,jp,ip)%particles(1:number_of_particles)
[1929]967!
[2312]968!--             Move only new particles. Moreover, limit random fluctuation
969!--             in order to prevent that particles move more than one grid box,
970!--             which would lead to problems concerning particle exchange
971!--             between processors in case pdx/pdy are larger than dx/dy,
972!--             respectively.
[1929]973                DO  n = local_start(kp,jp,ip), number_of_particles
[1359]974                   IF ( psl(particles(n)%group) /= psr(particles(n)%group) )  THEN
[1929]975                      rand_contr = ( random_function( iran_part ) - 0.5_wp ) * &
976                                     pdx(particles(n)%group)
[1359]977                      particles(n)%x = particles(n)%x +                        &
[2232]978                              MERGE( rand_contr, SIGN( dx, rand_contr ),       &
[1929]979                                     ABS( rand_contr ) < dx                    &
[2312]980                                   )
[1359]981                   ENDIF
982                   IF ( pss(particles(n)%group) /= psn(particles(n)%group) )  THEN
[1929]983                      rand_contr = ( random_function( iran_part ) - 0.5_wp ) * &
984                                     pdy(particles(n)%group)
[1359]985                      particles(n)%y = particles(n)%y +                        &
[2232]986                              MERGE( rand_contr, SIGN( dy, rand_contr ),       &
[1929]987                                     ABS( rand_contr ) < dy                    &
[2312]988                                   )
[1359]989                   ENDIF
990                   IF ( psb(particles(n)%group) /= pst(particles(n)%group) )  THEN
[1929]991                      rand_contr = ( random_function( iran_part ) - 0.5_wp ) * &
992                                     pdz(particles(n)%group)
[1359]993                      particles(n)%z = particles(n)%z +                        &
[3065]994                              MERGE( rand_contr, SIGN( dzw(kp), rand_contr ),  &
995                                     ABS( rand_contr ) < dzw(kp)               &
[2312]996                                   )
[1359]997                   ENDIF
998                ENDDO
[1]999!
[1929]1000!--             Identify particles located outside the model domain and reflect
1001!--             or absorb them if necessary.
[2698]1002                CALL lpm_boundary_conds( 'bottom/top', i, j, k )
[1929]1003!
1004!--             Furthermore, remove particles located in topography. Note, as
[2312]1005!--             the particle speed is still zero at this point, wall
[1929]1006!--             reflection boundary conditions will not work in this case.
1007                particles =>                                                   &
1008                       grid_particles(kp,jp,ip)%particles(1:number_of_particles)
1009                DO  n = local_start(kp,jp,ip), number_of_particles
[2606]1010                   i = particles(n)%x * ddx
1011                   j = particles(n)%y * ddy
[3065]1012                   k = particles(n)%z / dz(1) + 1 + offset_ocean_nzt
[2628]1013                   DO WHILE( zw(k) < particles(n)%z )
1014                      k = k + 1
1015                   ENDDO
1016                   DO WHILE( zw(k-1) > particles(n)%z )
1017                      k = k - 1
1018                   ENDDO
[2232]1019!
1020!--                Check if particle is within topography
1021                   IF ( .NOT. BTEST( wall_flags_0(k,j,i), 0 ) )  THEN
[1929]1022                      particles(n)%particle_mask = .FALSE.
1023                      deleted_particles = deleted_particles + 1
1024                   ENDIF
[2232]1025
[1929]1026                ENDDO
[1359]1027             ENDDO
1028          ENDDO
1029       ENDDO
[1]1030!
[1359]1031!--    Exchange particles between grid cells and processors
1032       CALL lpm_move_particle
1033       CALL lpm_exchange_horiz
[1]1034
[1359]1035    ENDIF
[1]1036!
[2312]1037!-- In case of random_start_position, delete particles identified by
1038!-- lpm_exchange_horiz and lpm_boundary_conds. Then sort particles into blocks,
1039!-- which is needed for a fast interpolation of the LES fields on the particle
[1359]1040!-- position.
[2606]1041    CALL lpm_sort_in_subboxes
[1]1042
1043!
[2265]1044!-- Determine the current number of particles
[1359]1045    DO  ip = nxl, nxr
1046       DO  jp = nys, nyn
1047          DO  kp = nzb+1, nzt
1048             number_of_particles         = number_of_particles                 &
1049                                           + prt_count(kp,jp,ip)
[1]1050          ENDDO
[1359]1051       ENDDO
1052    ENDDO
[1]1053!
[1822]1054!-- Calculate the number of particles of the total domain
[1]1055#if defined( __parallel )
[1359]1056    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
1057    CALL MPI_ALLREDUCE( number_of_particles, total_number_of_particles, 1, &
1058    MPI_INTEGER, MPI_SUM, comm2d, ierr )
[1]1059#else
[1359]1060    total_number_of_particles = number_of_particles
[1]1061#endif
1062
[1359]1063    RETURN
[1]1064
[1359]1065 END SUBROUTINE lpm_create_particle
[336]1066
[1871]1067 SUBROUTINE lpm_init_aerosols(local_start)
1068
1069    USE arrays_3d,                                                             &
[3274]1070        ONLY: hyp, pt, q, exner
[1871]1071
[3274]1072    USE basic_constants_and_equations_mod,                                     &
1073        ONLY: molecular_weight_of_solute, molecular_weight_of_water, magnus,   &
[3361]1074              pi, rd_d_rv, rho_l, r_v, rho_s, vanthoff
[1871]1075
1076    USE kinds
1077
1078    USE particle_attributes,                                                   &
[2375]1079        ONLY: aero_species, aero_type, aero_weight, log_sigma, na, rm
[1871]1080
1081    IMPLICIT NONE
1082
[2122]1083    REAL(wp)  :: afactor            !< curvature effects
[1871]1084    REAL(wp)  :: bfactor            !< solute effects
[2312]1085    REAL(wp)  :: dlogr              !< logarithmic width of radius bin
[1871]1086    REAL(wp)  :: e_a                !< vapor pressure
1087    REAL(wp)  :: e_s                !< saturation vapor pressure
[2312]1088    REAL(wp)  :: rmin = 0.005e-6_wp !< minimum aerosol radius
1089    REAL(wp)  :: rmax = 10.0e-6_wp  !< maximum aerosol radius
1090    REAL(wp)  :: r_mid              !< mean radius of bin
1091    REAL(wp)  :: r_l                !< left radius of bin
1092    REAL(wp)  :: r_r                !< right radius of bin
[2122]1093    REAL(wp)  :: sigma              !< surface tension
[1871]1094    REAL(wp)  :: t_int              !< temperature
1095
[1890]1096    INTEGER(iwp), DIMENSION(nzb:nzt+1,nysg:nyng,nxlg:nxrg), INTENT(IN) ::  local_start !<
[1871]1097
1098    INTEGER(iwp)  :: n              !<
1099    INTEGER(iwp)  :: ip             !<
1100    INTEGER(iwp)  :: jp             !<
1101    INTEGER(iwp)  :: kp             !<
1102
1103!
[2375]1104!-- Set constants for different aerosol species
1105    IF ( TRIM(aero_species) .EQ. 'nacl' ) THEN
1106       molecular_weight_of_solute = 0.05844_wp 
1107       rho_s                      = 2165.0_wp
1108       vanthoff                   = 2.0_wp
1109    ELSEIF ( TRIM(aero_species) .EQ. 'c3h4o4' ) THEN
1110       molecular_weight_of_solute = 0.10406_wp 
1111       rho_s                      = 1600.0_wp
1112       vanthoff                   = 1.37_wp
1113    ELSEIF ( TRIM(aero_species) .EQ. 'nh4o3' ) THEN
1114       molecular_weight_of_solute = 0.08004_wp 
1115       rho_s                      = 1720.0_wp
1116       vanthoff                   = 2.31_wp
1117    ELSE
1118       WRITE( message_string, * ) 'unknown aerosol species ',   &
1119                                'aero_species = "', TRIM( aero_species ), '"'
1120       CALL message( 'lpm_init', 'PA0470', 1, 2, 0, 6, 0 )
1121    ENDIF
1122!
[2312]1123!-- The following typical aerosol spectra are taken from Jaenicke (1993):
1124!-- Tropospheric aerosols. Published in Aerosol-Cloud-Climate Interactions.
1125    IF ( TRIM(aero_type) .EQ. 'polar' )  THEN
1126       na        = (/ 2.17e1, 1.86e-1, 3.04e-4 /) * 1.0E6
1127       rm        = (/ 0.0689, 0.375, 4.29 /) * 1.0E-6
1128       log_sigma = (/ 0.245, 0.300, 0.291 /)
1129    ELSEIF ( TRIM(aero_type) .EQ. 'background' )  THEN
1130       na        = (/ 1.29e2, 5.97e1, 6.35e1 /) * 1.0E6
1131       rm        = (/ 0.0036, 0.127, 0.259 /) * 1.0E-6
1132       log_sigma = (/ 0.645, 0.253, 0.425 /)
1133    ELSEIF ( TRIM(aero_type) .EQ. 'maritime' )  THEN
1134       na        = (/ 1.33e2, 6.66e1, 3.06e0 /) * 1.0E6
1135       rm        = (/ 0.0039, 0.133, 0.29 /) * 1.0E-6
1136       log_sigma = (/ 0.657, 0.210, 0.396 /)
1137    ELSEIF ( TRIM(aero_type) .EQ. 'continental' )  THEN
1138       na        = (/ 3.20e3, 2.90e3, 3.00e-1 /) * 1.0E6
1139       rm        = (/ 0.01, 0.058, 0.9 /) * 1.0E-6
1140       log_sigma = (/ 0.161, 0.217, 0.380 /)
1141    ELSEIF ( TRIM(aero_type) .EQ. 'desert' )  THEN
1142       na        = (/ 7.26e2, 1.14e3, 1.78e-1 /) * 1.0E6
1143       rm        = (/ 0.001, 0.0188, 10.8 /) * 1.0E-6
1144       log_sigma = (/ 0.247, 0.770, 0.438 /)
1145    ELSEIF ( TRIM(aero_type) .EQ. 'rural' )  THEN
1146       na        = (/ 6.65e3, 1.47e2, 1.99e3 /) * 1.0E6
1147       rm        = (/ 0.00739, 0.0269, 0.0419 /) * 1.0E-6
1148       log_sigma = (/ 0.225, 0.557, 0.266 /)
1149    ELSEIF ( TRIM(aero_type) .EQ. 'urban' )  THEN
1150       na        = (/ 9.93e4, 1.11e3, 3.64e4 /) * 1.0E6
1151       rm        = (/ 0.00651, 0.00714, 0.0248 /) * 1.0E-6
1152       log_sigma = (/ 0.245, 0.666, 0.337 /)
1153    ELSEIF ( TRIM(aero_type) .EQ. 'user' )  THEN
1154       CONTINUE
1155    ELSE
1156       WRITE( message_string, * ) 'unknown aerosol type ',   &
1157                                'aero_type = "', TRIM( aero_type ), '"'
1158       CALL message( 'lpm_init', 'PA0459', 1, 2, 0, 6, 0 )
[1871]1159    ENDIF
1160
1161    DO  ip = nxl, nxr
1162       DO  jp = nys, nyn
1163          DO  kp = nzb+1, nzt
1164
1165             number_of_particles = prt_count(kp,jp,ip)
1166             IF ( number_of_particles <= 0 )  CYCLE
1167             particles => grid_particles(kp,jp,ip)%particles(1:number_of_particles)
[2312]1168
1169             dlogr   = ( LOG10(rmax) - LOG10(rmin) ) / ( number_of_particles - local_start(kp,jp,ip) + 1 )
[1871]1170!
1171!--          Initialize the aerosols with a predefined spectral distribution
[2312]1172!--          of the dry radius (logarithmically increasing bins) and a varying
[1871]1173!--          weighting factor
[2312]1174             DO  n = local_start(kp,jp,ip), number_of_particles  !only new particles
[1871]1175
[2312]1176                r_l   = 10.0**( LOG10( rmin ) + (n-1) * dlogr )
1177                r_r   = 10.0**( LOG10( rmin ) + n * dlogr )
1178                r_mid = SQRT( r_l * r_r )
[1871]1179
[2312]1180                particles(n)%aux1          = r_mid
1181                particles(n)%weight_factor =                                           &
1182                   ( na(1) / ( SQRT( 2.0 * pi ) * log_sigma(1) ) *                     &
1183                     EXP( - LOG10( r_mid / rm(1) )**2 / ( 2.0 * log_sigma(1)**2 ) ) +  &
1184                     na(2) / ( SQRT( 2.0 * pi ) * log_sigma(2) ) *                     &
1185                     EXP( - LOG10( r_mid / rm(2) )**2 / ( 2.0 * log_sigma(2)**2 ) ) +  &
1186                     na(3) / ( SQRT( 2.0 * pi ) * log_sigma(3) ) *                     &
1187                     EXP( - LOG10( r_mid / rm(3) )**2 / ( 2.0 * log_sigma(3)**2 ) )    &
[3039]1188                   ) * ( LOG10(r_r) - LOG10(r_l) ) * ( dx * dy * dzw(kp) )
[1871]1189
[2312]1190!
1191!--             Multiply weight_factor with the namelist parameter aero_weight
1192!--             to increase or decrease the number of simulated aerosols
1193                particles(n)%weight_factor = particles(n)%weight_factor * aero_weight
[1871]1194
[2312]1195                IF ( particles(n)%weight_factor - FLOOR(particles(n)%weight_factor,KIND=wp) &
1196                     .GT. random_function( iran_part ) )  THEN
[3524]1197                   particles(n)%weight_factor = FLOOR(particles(n)%weight_factor,KIND=wp) + 1.0_wp
[2312]1198                ELSE
1199                   particles(n)%weight_factor = FLOOR(particles(n)%weight_factor,KIND=wp)
[1871]1200                ENDIF
1201!
[2312]1202!--             Unnecessary particles will be deleted
1203                IF ( particles(n)%weight_factor .LE. 0.0 )  particles(n)%particle_mask = .FALSE.
[1871]1204
[2312]1205             ENDDO
[1871]1206!
1207!--          Set particle radius to equilibrium radius based on the environmental
1208!--          supersaturation (Khvorostyanov and Curry, 2007, JGR). This avoids
[2312]1209!--          the sometimes lengthy growth toward their equilibrium radius within
[1871]1210!--          the simulation.
[3274]1211             t_int  = pt(kp,jp,ip) * exner(kp)
[1871]1212
[2608]1213             e_s = magnus( t_int )
[3361]1214             e_a = q(kp,jp,ip) * hyp(kp) / ( q(kp,jp,ip) + rd_d_rv )
[1871]1215
[2122]1216             sigma   = 0.0761_wp - 0.000155_wp * ( t_int - 273.15_wp )
1217             afactor = 2.0_wp * sigma / ( rho_l * r_v * t_int )
1218
1219             bfactor = vanthoff * molecular_weight_of_water *    &
1220                       rho_s / ( molecular_weight_of_solute * rho_l )
[1871]1221!
[2312]1222!--          The formula is only valid for subsaturated environments. For
[2122]1223!--          supersaturations higher than -5 %, the supersaturation is set to -5%.
1224             IF ( e_a / e_s >= 0.95_wp )  e_a = 0.95_wp * e_s
[1871]1225
[2312]1226             DO  n = local_start(kp,jp,ip), number_of_particles  !only new particles
[2122]1227!
[2312]1228!--             For details on this equation, see Eq. (14) of Khvorostyanov and
1229!--             Curry (2007, JGR)
[2122]1230                particles(n)%radius = bfactor**0.3333333_wp *                  &
[2312]1231                   particles(n)%aux1 / ( 1.0_wp - e_a / e_s )**0.3333333_wp / &
[2122]1232                   ( 1.0_wp + ( afactor / ( 3.0_wp * bfactor**0.3333333_wp *   &
[2312]1233                     particles(n)%aux1 ) ) /                                  &
[2122]1234                     ( 1.0_wp - e_a / e_s )**0.6666666_wp                      &
1235                   )
[1871]1236
[2122]1237             ENDDO
[1871]1238
1239          ENDDO
1240       ENDDO
1241    ENDDO
1242
1243 END SUBROUTINE lpm_init_aerosols
1244
[1359]1245END MODULE lpm_init_mod
Note: See TracBrowser for help on using the repository browser.