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

Last change on this file since 2293 was 2274, checked in by Giersch, 7 years ago

Changed error messages

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