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

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

enable particle advection with grid stretching and some formatation changes in lpm

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