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

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

bugfix for lcm with grid stretching

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