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

Last change on this file since 3657 was 3655, checked in by knoop, 6 years ago

Bugfix: made "unit" and "found" intend INOUT in module interface subroutines + automatic copyright update

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