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

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

Introduced global constant rd_d_rv=0.622

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