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

Last change on this file since 2271 was 2265, checked in by schwenkel, 7 years ago

unused variables removed

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