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

Last change on this file since 1890 was 1890, checked in by hoffmann, 8 years ago

improvements for the consideration of aerosols in the LCM

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