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

Last change on this file since 1755 was 1726, checked in by hoffmann, 9 years ago

last commit documented

  • Property svn:keywords set to Id
File size: 36.4 KB
RevLine 
[1682]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!
[1691]16! Copyright 1997-2015 Leibniz Universitaet Hannover
[1036]17!--------------------------------------------------------------------------------!
18!
[254]19! Current revisions:
[1]20! -----------------
[1726]21!
[1692]22!
[1321]23! Former revisions:
24! -----------------
25! $Id: lpm_init.f90 1726 2015-11-17 13:03:53Z maronga $
26!
[1726]27! 1725 2015-11-17 13:01:51Z hoffmann
28! Bugfix: Processor-dependent seed for random function is generated before it is
29! used.
30!
[1692]31! 1691 2015-10-26 16:17:44Z maronga
32! Renamed prandtl_layer to constant_flux_layer.
33!
[1686]34! 1685 2015-10-08 07:32:13Z raasch
35! bugfix concerning vertical index offset in case of ocean
36!
[1683]37! 1682 2015-10-07 23:56:08Z knoop
38! Code annotations made doxygen readable
39!
[1576]40! 1575 2015-03-27 09:56:27Z raasch
41! initial vertical particle position is allowed to follow the topography
42!
[1360]43! 1359 2014-04-11 17:15:14Z hoffmann
44! New particle structure integrated.
45! Kind definition added to all floating point numbers.
46! lpm_init changed form a subroutine to a module.
47!
[1329]48! 1327 2014-03-21 11:00:16Z raasch
49! -netcdf_output
50!
[1323]51! 1322 2014-03-20 16:38:49Z raasch
52! REAL functions provided with KIND-attribute
53!
[1321]54! 1320 2014-03-20 08:40:49Z raasch
[1320]55! ONLY-attribute added to USE-statements,
56! kind-parameters added to all INTEGER and REAL declaration statements,
57! kinds are defined in new module kinds,
58! revision history before 2012 removed,
59! comment fields (!:) to be used for variable explanations added to
60! all variable declaration statements
61! bugfix: #if defined( __parallel ) added
[850]62!
[1315]63! 1314 2014-03-14 18:25:17Z suehring
64! Vertical logarithmic interpolation of horizontal particle speed for particles
65! between roughness height and first vertical grid level.
66!
[1093]67! 1092 2013-02-02 11:24:22Z raasch
68! unused variables removed
69!
[1037]70! 1036 2012-10-22 13:43:42Z raasch
71! code put under GPL (PALM 3.9)
72!
[850]73! 849 2012-03-15 10:35:09Z raasch
[849]74! routine renamed: init_particles -> lpm_init
75! de_dx, de_dy, de_dz are allocated here (instead of automatic arrays in
76! advec_particles),
77! sort_particles renamed lpm_sort_arrays, user_init_particles renamed lpm_init
[392]78!
[829]79! 828 2012-02-21 12:00:36Z raasch
80! call of init_kernels, particle feature color renamed class
81!
[826]82! 824 2012-02-17 09:09:57Z raasch
83! particle attributes speed_x|y|z_sgs renamed rvar1|2|3,
84! array particles implemented as pointer
85!
[668]86! 667 2010-12-23 12:06:00Z suehring/gryschka
87! nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng for allocation
88! of arrays.
89!
[1]90! Revision 1.1  1999/11/25 16:22:38  raasch
91! Initial revision
92!
93!
94! Description:
95! ------------
[1682]96!> This routine initializes a set of particles and their attributes (position,
97!> radius, ..) which are used by the Lagrangian particle model (see lpm).
[1]98!------------------------------------------------------------------------------!
[1682]99 MODULE lpm_init_mod
100 
[1]101
[1320]102    USE arrays_3d,                                                             &
103        ONLY:  de_dx, de_dy, de_dz, zu, zw, z0
104
105    USE cloud_parameters,                                                      &
106        ONLY:  curvature_solution_effects
107
108    USE control_parameters,                                                    &
[1691]109        ONLY:  cloud_droplets, constant_flux_layer, current_timestep_number,   &
110               dz, initializing_actions, message_string, netcdf_data_format,   &
111               ocean, simulated_time
[1320]112
113    USE dvrp_variables,                                                        &
114        ONLY:  particle_color, particle_dvrpsize
115
116    USE grid_variables,                                                        &
[1359]117        ONLY:  ddx, dx, ddy, dy
[1320]118
119    USE indices,                                                               &
[1575]120        ONLY:  nx, nxl, nxlg, nxrg, nxr, ny, nyn, nys, nyng, nysg, nz, nzb,    &
121               nzb_w_inner, nzt
[1320]122
123    USE kinds
124
125    USE lpm_collision_kernels_mod,                                             &
126        ONLY:  init_kernels
127
128    USE particle_attributes,                                                   &
[1359]129        ONLY:   alloc_factor, bc_par_b, bc_par_lr, bc_par_ns, bc_par_t,        &
130                block_offset, block_offset_def, collision_kernel,              &
131                density_ratio, dvrp_psize, grid_particles,                     &
132                initial_weighting_factor, ibc_par_b, ibc_par_lr, ibc_par_ns,   &
133                ibc_par_t, iran_part, log_z_z0,                                &
134                max_number_of_particle_groups, maximum_number_of_particles,    &
135                maximum_number_of_tailpoints, maximum_number_of_tails,         &
136                minimum_tailpoint_distance, min_nr_particle,                   &
137                mpi_particle_type, new_tail_id,                                &
[1320]138                number_of_initial_tails, number_of_particles,                  &
139                number_of_particle_groups, number_of_sublayers,                &
[1685]140                number_of_tails, offset_ocean_nzt, offset_ocean_nzt_m1,        &
[1359]141                particles, particle_advection_start, particle_groups,          &
142                particle_groups_type, particles_per_point,                     &
[1320]143                particle_tail_coordinates,  particle_type, pdx, pdy, pdz,      &
[1359]144                prt_count, psb, psl, psn, psr, pss, pst,                       &
[1320]145                radius, random_start_position, read_particles_from_restartfile,&
[1575]146                seed_follows_topography, skip_particles_for_tail, sort_count,  &
[1359]147                tail_mask, total_number_of_particles, total_number_of_tails,   &
[1320]148                use_particle_tails, use_sgs_for_particles,                     &
[1359]149                write_particle_statistics, uniform_particles, zero_particle,   &
150                z0_av_global
[1320]151
[1]152    USE pegrid
153
[1320]154    USE random_function_mod,                                                   &
155        ONLY:  random_function
[1]156
[1359]157    IMPLICIT NONE
[1320]158
[1359]159    PRIVATE
160
[1682]161    INTEGER(iwp), PARAMETER         :: PHASE_INIT    = 1  !<
162    INTEGER(iwp), PARAMETER, PUBLIC :: PHASE_RELEASE = 2  !<
[1359]163
164    INTERFACE lpm_init
165       MODULE PROCEDURE lpm_init
166    END INTERFACE lpm_init
167
168    INTERFACE lpm_create_particle
169       MODULE PROCEDURE lpm_create_particle
170    END INTERFACE lpm_create_particle
171
172    PUBLIC lpm_init, lpm_create_particle
173
174CONTAINS
175
[1682]176!------------------------------------------------------------------------------!
177! Description:
178! ------------
179!> @todo Missing subroutine description.
180!------------------------------------------------------------------------------!
[1359]181 SUBROUTINE lpm_init
182
183    USE lpm_collision_kernels_mod,                                             &
184        ONLY:  init_kernels
185
[1]186    IMPLICIT NONE
187
[1682]188    INTEGER(iwp) ::  i                           !<
189    INTEGER(iwp) ::  ip                          !<
190    INTEGER(iwp) ::  j                           !<
191    INTEGER(iwp) ::  jp                          !<
192    INTEGER(iwp) ::  k                           !<
193    INTEGER(iwp) ::  kp                          !<
194    INTEGER(iwp) ::  n                           !<
195    INTEGER(iwp) ::  nn                          !<
[1320]196
[1]197#if defined( __parallel )
[1682]198    INTEGER(iwp), DIMENSION(3) ::  blocklengths  !<
199    INTEGER(iwp), DIMENSION(3) ::  displacements !<
200    INTEGER(iwp), DIMENSION(3) ::  types         !<
[1]201#endif
[1682]202    LOGICAL ::  uniform_particles_l              !<
[1]203
[1682]204    REAL(wp) ::  height_int                      !<
205    REAL(wp) ::  height_p                        !<
206    REAL(wp) ::  z_p                             !<
207    REAL(wp) ::  z0_av_local                     !<
[1]208
209#if defined( __parallel )
210!
211!-- Define MPI derived datatype for FORTRAN datatype particle_type (see module
[82]212!-- particle_attributes). Integer length is 4 byte, Real is 8 byte
[1359]213#if defined( __twocachelines )
214    blocklengths(1)  =  7;  blocklengths(2)  =  18;  blocklengths(3)  =   1
215    displacements(1) =  0;  displacements(2) =  64;  displacements(3) = 128
[82]216
[1359]217    types(1) = MPI_REAL                               ! 64 bit words
218    types(2) = MPI_INTEGER                            ! 32 Bit words
219    types(3) = MPI_UB
220#else
221    blocklengths(1)  = 19;  blocklengths(2)  =   6;  blocklengths(3)  =   1
222    displacements(1) =  0;  displacements(2) = 152;  displacements(3) = 176
223
[1]224    types(1) = MPI_REAL
225    types(2) = MPI_INTEGER
226    types(3) = MPI_UB
[1359]227#endif
[1]228    CALL MPI_TYPE_STRUCT( 3, blocklengths, displacements, types, &
229                          mpi_particle_type, ierr )
230    CALL MPI_TYPE_COMMIT( mpi_particle_type, ierr )
231#endif
232
233!
[150]234!-- In case of oceans runs, the vertical index calculations need an offset,
235!-- because otherwise the k indices will become negative
236    IF ( ocean )  THEN
237       offset_ocean_nzt    = nzt
238       offset_ocean_nzt_m1 = nzt - 1
239    ENDIF
240
[1359]241!
242!-- Define block offsets for dividing a gridcell in 8 sub cells
[150]243
[1359]244    block_offset(0) = block_offset_def (-1,-1,-1)
245    block_offset(1) = block_offset_def (-1,-1, 0)
246    block_offset(2) = block_offset_def (-1, 0,-1)
247    block_offset(3) = block_offset_def (-1, 0, 0)
248    block_offset(4) = block_offset_def ( 0,-1,-1)
249    block_offset(5) = block_offset_def ( 0,-1, 0)
250    block_offset(6) = block_offset_def ( 0, 0,-1)
251    block_offset(7) = block_offset_def ( 0, 0, 0)
[150]252!
[1]253!-- Check the number of particle groups.
254    IF ( number_of_particle_groups > max_number_of_particle_groups )  THEN
[274]255       WRITE( message_string, * ) 'max_number_of_particle_groups =',      &
256                                  max_number_of_particle_groups ,         &
[254]257                                  '&number_of_particle_groups reset to ', &
258                                  max_number_of_particle_groups
[849]259       CALL message( 'lpm_init', 'PA0213', 0, 1, 0, 6, 0 )
[1]260       number_of_particle_groups = max_number_of_particle_groups
261    ENDIF
262
263!
264!-- Set default start positions, if necessary
[1359]265    IF ( psl(1) == 9999999.9_wp )  psl(1) = -0.5_wp * dx
266    IF ( psr(1) == 9999999.9_wp )  psr(1) = ( nx + 0.5_wp ) * dx
267    IF ( pss(1) == 9999999.9_wp )  pss(1) = -0.5_wp * dy
268    IF ( psn(1) == 9999999.9_wp )  psn(1) = ( ny + 0.5_wp ) * dy
269    IF ( psb(1) == 9999999.9_wp )  psb(1) = zu(nz/2)
270    IF ( pst(1) == 9999999.9_wp )  pst(1) = psb(1)
[1]271
[1359]272    IF ( pdx(1) == 9999999.9_wp  .OR.  pdx(1) == 0.0_wp )  pdx(1) = dx
273    IF ( pdy(1) == 9999999.9_wp  .OR.  pdy(1) == 0.0_wp )  pdy(1) = dy
274    IF ( pdz(1) == 9999999.9_wp  .OR.  pdz(1) == 0.0_wp )  pdz(1) = zu(2) - zu(1)
[1]275
276    DO  j = 2, number_of_particle_groups
[1359]277       IF ( psl(j) == 9999999.9_wp )  psl(j) = psl(j-1)
278       IF ( psr(j) == 9999999.9_wp )  psr(j) = psr(j-1)
279       IF ( pss(j) == 9999999.9_wp )  pss(j) = pss(j-1)
280       IF ( psn(j) == 9999999.9_wp )  psn(j) = psn(j-1)
281       IF ( psb(j) == 9999999.9_wp )  psb(j) = psb(j-1)
282       IF ( pst(j) == 9999999.9_wp )  pst(j) = pst(j-1)
283       IF ( pdx(j) == 9999999.9_wp  .OR.  pdx(j) == 0.0_wp )  pdx(j) = pdx(j-1)
284       IF ( pdy(j) == 9999999.9_wp  .OR.  pdy(j) == 0.0_wp )  pdy(j) = pdy(j-1)
285       IF ( pdz(j) == 9999999.9_wp  .OR.  pdz(j) == 0.0_wp )  pdz(j) = pdz(j-1)
[1]286    ENDDO
287
288!
[849]289!-- Allocate arrays required for calculating particle SGS velocities
290    IF ( use_sgs_for_particles )  THEN
291       ALLOCATE( de_dx(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
292                 de_dy(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
293                 de_dz(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
294    ENDIF
295
296!
[1314]297!-- Allocate array required for logarithmic vertical interpolation of
298!-- horizontal particle velocities between the surface and the first vertical
299!-- grid level. In order to avoid repeated CPU cost-intensive CALLS of
300!-- intrinsic FORTRAN procedure LOG(z/z0), LOG(z/z0) is precalculated for
301!-- several heights. Splitting into 20 sublayers turned out to be sufficient.
302!-- To obtain exact height levels of particles, linear interpolation is applied
303!-- (see lpm_advec.f90).
[1691]304    IF ( constant_flux_layer )  THEN
[1314]305       
306       ALLOCATE ( log_z_z0(0:number_of_sublayers) ) 
307       z_p         = zu(nzb+1) - zw(nzb)
308
309!
310!--    Calculate horizontal mean value of z0 used for logartihmic
311!--    interpolation. Note: this is not exact for heterogeneous z0.
312!--    However, sensitivity studies showed that the effect is
313!--    negligible.
314       z0_av_local  = SUM( z0(nys:nyn,nxl:nxr) )
[1359]315       z0_av_global = 0.0_wp
[1314]316
[1320]317#if defined( __parallel )
[1314]318       CALL MPI_ALLREDUCE(z0_av_local, z0_av_global, 1, MPI_REAL, MPI_SUM, &
319                          comm2d, ierr )
[1320]320#else
321       z0_av_global = z0_av_local
322#endif
[1314]323
324       z0_av_global = z0_av_global  / ( ( ny + 1 ) * ( nx + 1 ) )
325!
326!--    Horizontal wind speed is zero below and at z0
[1359]327       log_z_z0(0) = 0.0_wp
[1314]328!
329!--    Calculate vertical depth of the sublayers
[1322]330       height_int  = ( z_p - z0_av_global ) / REAL( number_of_sublayers, KIND=wp )
[1314]331!
332!--    Precalculate LOG(z/z0)
[1359]333       height_p    = 0.0_wp
[1314]334       DO  k = 1, number_of_sublayers
335
336          height_p    = height_p + height_int
337          log_z_z0(k) = LOG( height_p / z0_av_global )
338
339       ENDDO
340
341
342    ENDIF
343
344!
[1359]345!-- Check boundary condition and set internal variables
346    SELECT CASE ( bc_par_b )
347   
348       CASE ( 'absorb' )
349          ibc_par_b = 1
350
351       CASE ( 'reflect' )
352          ibc_par_b = 2
353         
354       CASE DEFAULT
355          WRITE( message_string, * )  'unknown boundary condition ',   &
356                                       'bc_par_b = "', TRIM( bc_par_b ), '"'
357          CALL message( 'lpm_init', 'PA0217', 1, 2, 0, 6, 0 )
358         
359    END SELECT
360    SELECT CASE ( bc_par_t )
361   
362       CASE ( 'absorb' )
363          ibc_par_t = 1
364
365       CASE ( 'reflect' )
366          ibc_par_t = 2
367         
368       CASE DEFAULT
369          WRITE( message_string, * ) 'unknown boundary condition ',   &
370                                     'bc_par_t = "', TRIM( bc_par_t ), '"'
371          CALL message( 'lpm_init', 'PA0218', 1, 2, 0, 6, 0 )
372         
373    END SELECT
374    SELECT CASE ( bc_par_lr )
375
376       CASE ( 'cyclic' )
377          ibc_par_lr = 0
378
379       CASE ( 'absorb' )
380          ibc_par_lr = 1
381
382       CASE ( 'reflect' )
383          ibc_par_lr = 2
384         
385       CASE DEFAULT
386          WRITE( message_string, * ) 'unknown boundary condition ',   &
387                                     'bc_par_lr = "', TRIM( bc_par_lr ), '"'
388          CALL message( 'lpm_init', 'PA0219', 1, 2, 0, 6, 0 )
389         
390    END SELECT
391    SELECT CASE ( bc_par_ns )
392
393       CASE ( 'cyclic' )
394          ibc_par_ns = 0
395
396       CASE ( 'absorb' )
397          ibc_par_ns = 1
398
399       CASE ( 'reflect' )
400          ibc_par_ns = 2
401         
402       CASE DEFAULT
403          WRITE( message_string, * ) 'unknown boundary condition ',   &
404                                     'bc_par_ns = "', TRIM( bc_par_ns ), '"'
405          CALL message( 'lpm_init', 'PA0220', 1, 2, 0, 6, 0 )
406         
407    END SELECT
408
409!
[828]410!-- Initialize collision kernels
411    IF ( collision_kernel /= 'none' )  CALL init_kernels
412
413!
[1]414!-- For the first model run of a possible job chain initialize the
[849]415!-- particles, otherwise read the particle data from restart file.
[1]416    IF ( TRIM( initializing_actions ) == 'read_restart_data'  &
417         .AND.  read_particles_from_restartfile )  THEN
418
[849]419       CALL lpm_read_restart_file
[1]420
421    ELSE
422
423!
424!--    Allocate particle arrays and set attributes of the initial set of
425!--    particles, which can be also periodically released at later times.
426!--    Also allocate array for particle tail coordinates, if needed.
[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#if defined( __twocachelines )
441       zero_particle = particle_type( 0.0_wp, 0.0_sp, 0.0_sp, 0.0_sp, 0.0_sp,  &
442                                      0.0_sp, 0.0_sp, 0.0_wp, 0.0_wp, 0.0_wp,  &
443                                      0, .FALSE., 0.0_wp, 0.0_wp, 0.0_wp,      &
444                                      0.0_sp, 0.0_sp, 0.0_sp, 0.0_sp, 0.0_sp,  &
445                                      0.0_sp, 0, 0, 0, -1)
446#else
447       zero_particle = particle_type( 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp,  &
448                                      0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp,  &
449                                      0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp,  &
450                                      0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0, 0, 0, &
451                                      0, .FALSE., -1)
452#endif
453       particle_groups = particle_groups_type( 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp )
[1]454
455!
456!--    Set the default particle size used for dvrp plots
[1359]457       IF ( dvrp_psize == 9999999.9_wp )  dvrp_psize = 0.2_wp * dx
[1]458
459!
460!--    Set values for the density ratio and radius for all particle
461!--    groups, if necessary
[1359]462       IF ( density_ratio(1) == 9999999.9_wp )  density_ratio(1) = 0.0_wp
463       IF ( radius(1)        == 9999999.9_wp )  radius(1) = 0.0_wp
[1]464       DO  i = 2, number_of_particle_groups
[1359]465          IF ( density_ratio(i) == 9999999.9_wp )  THEN
[1]466             density_ratio(i) = density_ratio(i-1)
467          ENDIF
[1359]468          IF ( radius(i) == 9999999.9_wp )  radius(i) = radius(i-1)
[1]469       ENDDO
470
471       DO  i = 1, number_of_particle_groups
[1359]472          IF ( density_ratio(i) /= 0.0_wp  .AND.  radius(i) == 0 )  THEN
[254]473             WRITE( message_string, * ) 'particle group #', i, 'has a', &
474                                        'density ratio /= 0 but radius = 0'
[849]475             CALL message( 'lpm_init', 'PA0215', 1, 2, 0, 6, 0 )
[1]476          ENDIF
477          particle_groups(i)%density_ratio = density_ratio(i)
478          particle_groups(i)%radius        = radius(i)
479       ENDDO
480
481!
[1359]482!--    Set a seed value for the random number generator to be exclusively
483!--    used for the particle code. The generated random numbers should be
484!--    different on the different PEs.
485       iran_part = iran_part + myid
486
[1725]487       CALL lpm_create_particle (PHASE_INIT)
488
[1359]489!
490!--    User modification of initial particles
491       CALL user_lpm_init
492
493!
494!--    Open file for statistical informations about particle conditions
495       IF ( write_particle_statistics )  THEN
496          CALL check_open( 80 )
497          WRITE ( 80, 8000 )  current_timestep_number, simulated_time,         &
498                              number_of_particles,                             &
499                              maximum_number_of_particles
500          CALL close_file( 80 )
501       ENDIF
502
503!
504!--    Check if particles are really uniform in color and radius (dvrp_size)
505!--    (uniform_particles is preset TRUE)
506       IF ( uniform_particles )  THEN
507          DO  ip = nxl, nxr
508             DO  jp = nys, nyn
509                DO  kp = nzb+1, nzt
510
511                   n = prt_count(kp,jp,ip)
512                   IF ( MINVAL( grid_particles(kp,jp,ip)%particles(1:n)%dvrp_psize  ) ==     &
513                        MAXVAL( grid_particles(kp,jp,ip)%particles(1:n)%dvrp_psize  )  .AND. &
514                        MINVAL( grid_particles(kp,jp,ip)%particles(1:n)%class ) ==           &
515                        MAXVAL( grid_particles(kp,jp,ip)%particles(1:n)%class ) )  THEN
516                      uniform_particles_l = .TRUE.
517                   ELSE
518                      uniform_particles_l = .FALSE.
519                   ENDIF
520
521                ENDDO
522             ENDDO
523          ENDDO
524
525#if defined( __parallel )
526          IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
527          CALL MPI_ALLREDUCE( uniform_particles_l, uniform_particles, 1,       &
528                              MPI_LOGICAL, MPI_LAND, comm2d, ierr )
529#else
530          uniform_particles = uniform_particles_l
531#endif
532
533       ENDIF
534
535!
536!--    Particles will probably become none-uniform, if their size and color
537!--    will be determined by flow variables
538       IF ( particle_color /= 'none'  .OR.  particle_dvrpsize /= 'none' )  THEN
539          uniform_particles = .FALSE.
540       ENDIF
541
542! !kk    Not implemented aft individual particle array fort every gridcell
543! !
544! !--    Set the beginning of the particle tails and their age
545!        IF ( use_particle_tails )  THEN
546! !
547! !--       Choose the maximum number of tails with respect to the maximum number
548! !--       of particles and skip_particles_for_tail
549!           maximum_number_of_tails = maximum_number_of_particles / &
550!                                     skip_particles_for_tail
551!
552! !
553! !--       Create a minimum number of tails in case that there is no tail
554! !--       initially (otherwise, index errors will occur when adressing the
555! !--       arrays below)
556!           IF ( maximum_number_of_tails == 0 )  maximum_number_of_tails = 10
557!
558!           ALLOCATE( particle_tail_coordinates(maximum_number_of_tailpoints,5, &
559!                     maximum_number_of_tails),                                 &
560!                     new_tail_id(maximum_number_of_tails),                     &
561!                     tail_mask(maximum_number_of_tails) )
562!
563!           particle_tail_coordinates  = 0.0_wp
564!           minimum_tailpoint_distance = minimum_tailpoint_distance**2
565!           number_of_initial_tails    = number_of_tails
566!
567!           nn = 0
568!           DO  n = 1, number_of_particles
569! !
570! !--          Only for those particles marked above with a provisional tail_id
571! !--          tails will be created. Particles now get their final tail_id.
572!              IF ( particles(n)%tail_id /= 0 )  THEN
573!
574!                 nn = nn + 1
575!                 particles(n)%tail_id = nn
576!
577!                 particle_tail_coordinates(1,1,nn) = particles(n)%x
578!                 particle_tail_coordinates(1,2,nn) = particles(n)%y
579!                 particle_tail_coordinates(1,3,nn) = particles(n)%z
580!                 particle_tail_coordinates(1,4,nn) = particles(n)%class
581!                 particles(n)%tailpoints = 1
582!                 IF ( minimum_tailpoint_distance /= 0.0_wp )  THEN
583!                    particle_tail_coordinates(2,1,nn) = particles(n)%x
584!                    particle_tail_coordinates(2,2,nn) = particles(n)%y
585!                    particle_tail_coordinates(2,3,nn) = particles(n)%z
586!                    particle_tail_coordinates(2,4,nn) = particles(n)%class
587!                    particle_tail_coordinates(1:2,5,nn) = 0.0_wp
588!                    particles(n)%tailpoints = 2
589!                 ENDIF
590!
591!              ENDIF
592!           ENDDO
593!        ENDIF
594!
595! !
596! !--    Plot initial positions of particles (only if particle advection is
597! !--    switched on from the beginning of the simulation (t=0))
598!        IF ( particle_advection_start == 0.0_wp )  CALL data_output_dvrp
599
600    ENDIF
601
602!
603!-- To avoid programm abort, assign particles array to the local version of
604!-- first grid cell
605    number_of_particles = prt_count(nzb+1,nys,nxl)
606    particles => grid_particles(nzb+1,nys,nxl)%particles(1:number_of_particles)
607
608!
609!-- Formats
6108000 FORMAT (I6,1X,F7.2,4X,I10,71X,I10)
611
612 END SUBROUTINE lpm_init
613
[1682]614!------------------------------------------------------------------------------!
615! Description:
616! ------------
617!> @todo Missing subroutine description.
618!------------------------------------------------------------------------------!
[1359]619 SUBROUTINE lpm_create_particle (phase)
620
621    USE lpm_exchange_horiz_mod,                                                &
622        ONLY: lpm_exchange_horiz, lpm_move_particle, realloc_particles_array
623
624    USE lpm_pack_arrays_mod,                                                   &
625        ONLY: lpm_pack_all_arrays
626
627    IMPLICIT  NONE
628
[1682]629    INTEGER(iwp)               ::  alloc_size  !<
630    INTEGER(iwp)               ::  i           !<
631    INTEGER(iwp)               ::  ip          !<
632    INTEGER(iwp)               ::  j           !<
633    INTEGER(iwp)               ::  jp          !<
634    INTEGER(iwp)               ::  kp          !<
635    INTEGER(iwp)               ::  loop_stride !<
636    INTEGER(iwp)               ::  n           !<
637    INTEGER(iwp)               ::  new_size    !<
638    INTEGER(iwp)               ::  nn          !<
[1359]639
[1682]640    INTEGER(iwp), INTENT(IN)   ::  phase       !<
[1359]641
[1682]642    INTEGER(iwp), DIMENSION(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ::  local_count !<
643    INTEGER(iwp), DIMENSION(nzb:nzt+1,nysg:nyng,nxlg:nxrg) ::  local_start !<
[1359]644
[1682]645    LOGICAL                    ::  first_stride !<
[1359]646
[1682]647    REAL(wp)                   ::  pos_x !<
648    REAL(wp)                   ::  pos_y !<
649    REAL(wp)                   ::  pos_z !<
[1359]650
[1682]651    TYPE(particle_type),TARGET ::  tmp_particle !<
[1359]652
653!
654!-- Calculate particle positions and store particle attributes, if
655!-- particle is situated on this PE
656    DO  loop_stride = 1, 2
657       first_stride = (loop_stride == 1)
658       IF ( first_stride )   THEN
659          local_count = 0           ! count number of particles
660       ELSE
661          local_count = prt_count   ! Start address of new particles
662       ENDIF
663
[1]664       n = 0
665       DO  i = 1, number_of_particle_groups
666
667          pos_z = psb(i)
668
669          DO WHILE ( pos_z <= pst(i) )
670
671             pos_y = pss(i)
672
673             DO WHILE ( pos_y <= psn(i) )
674
[1359]675                IF ( pos_y >= ( nys - 0.5_wp ) * dy  .AND.  &
676                     pos_y <  ( nyn + 0.5_wp ) * dy )  THEN
[1]677
678                   pos_x = psl(i)
679
[1575]680            xloop: DO WHILE ( pos_x <= psr(i) )
[1]681
[1359]682                      IF ( pos_x >= ( nxl - 0.5_wp ) * dx  .AND.  &
683                           pos_x <  ( nxr + 0.5_wp ) * dx )  THEN
[1]684
685                         DO  j = 1, particles_per_point
686
687                            n = n + 1
[1359]688#if defined( __twocachelines )
689                            tmp_particle%x             = pos_x
690                            tmp_particle%y             = pos_y
691                            tmp_particle%z             = pos_z
692                            tmp_particle%age           = 0.0_sp
693                            tmp_particle%age_m         = 0.0_sp
694                            tmp_particle%dt_sum        = 0.0_wp
695                            tmp_particle%dvrp_psize    = dvrp_psize
696                            tmp_particle%e_m           = 0.0_sp
697                            IF ( curvature_solution_effects )  THEN
698!
699!--                            Initial values (internal timesteps, derivative)
700!--                            for Rosenbrock method
701                               tmp_particle%rvar1      = 1.0E-12_wp
702                               tmp_particle%rvar2      = 1.0E-3_wp
703                               tmp_particle%rvar3      = -9999999.9_wp
704                            ELSE
705!
706!--                            Initial values for SGS velocities
707                               tmp_particle%rvar1      = 0.0_wp
708                               tmp_particle%rvar2      = 0.0_wp
709                               tmp_particle%rvar3      = 0.0_wp
[1]710                            ENDIF
[1359]711                            tmp_particle%speed_x       = 0.0_sp
712                            tmp_particle%speed_y       = 0.0_sp
713                            tmp_particle%speed_z       = 0.0_sp
714                            tmp_particle%origin_x      = pos_x
715                            tmp_particle%origin_y      = pos_y
716                            tmp_particle%origin_z      = pos_z
717                            tmp_particle%radius        = particle_groups(i)%radius
718                            tmp_particle%weight_factor = initial_weighting_factor
719                            tmp_particle%class         = 1
720                            tmp_particle%group         = i
721                            tmp_particle%tailpoints    = 0
722                            tmp_particle%particle_mask = .TRUE.
723#else
724                            tmp_particle%x             = pos_x
725                            tmp_particle%y             = pos_y
726                            tmp_particle%z             = pos_z
727                            tmp_particle%age           = 0.0_wp
728                            tmp_particle%age_m         = 0.0_wp
729                            tmp_particle%dt_sum        = 0.0_wp
730                            tmp_particle%dvrp_psize    = dvrp_psize
731                            tmp_particle%e_m           = 0.0_wp
[824]732                            IF ( curvature_solution_effects )  THEN
733!
734!--                            Initial values (internal timesteps, derivative)
735!--                            for Rosenbrock method
[1359]736                               tmp_particle%rvar1      = 1.0E-12_wp
737                               tmp_particle%rvar2      = 1.0E-3_wp
738                               tmp_particle%rvar3      = -9999999.9_wp
[824]739                            ELSE
740!
741!--                            Initial values for SGS velocities
[1359]742                               tmp_particle%rvar1      = 0.0_wp
743                               tmp_particle%rvar2      = 0.0_wp
744                               tmp_particle%rvar3      = 0.0_wp
[824]745                            ENDIF
[1359]746                            tmp_particle%speed_x       = 0.0_wp
747                            tmp_particle%speed_y       = 0.0_wp
748                            tmp_particle%speed_z       = 0.0_wp
749                            tmp_particle%origin_x      = pos_x
750                            tmp_particle%origin_y      = pos_y
751                            tmp_particle%origin_z      = pos_z
752                            tmp_particle%radius        = particle_groups(i)%radius
753                            tmp_particle%weight_factor = initial_weighting_factor
754                            tmp_particle%class         = 1
755                            tmp_particle%group         = i
756                            tmp_particle%tailpoints    = 0
757                            tmp_particle%particle_mask = .TRUE.
758#endif
[1]759                            IF ( use_particle_tails  .AND. &
760                                 MOD( n, skip_particles_for_tail ) == 0 )  THEN
761                               number_of_tails         = number_of_tails + 1
762!
763!--                            This is a temporary provisional setting (see
764!--                            further below!)
[1359]765                               tmp_particle%tail_id    = number_of_tails
[1]766                            ELSE
[1359]767                               tmp_particle%tail_id    = 0
[1]768                            ENDIF
[1575]769!
770!--                         Determine the grid indices of the particle position
[1359]771                            ip = ( tmp_particle%x + 0.5_wp * dx ) * ddx
772                            jp = ( tmp_particle%y + 0.5_wp * dy ) * ddy
[1685]773                            kp = tmp_particle%z / dz + 1 + offset_ocean_nzt
[1]774
[1575]775                            IF ( seed_follows_topography )  THEN
776!
777!--                            Particle height is given relative to topography
778                               kp = kp + nzb_w_inner(jp,ip)
779                               tmp_particle%z = tmp_particle%z + zw(kp)
780                               IF ( kp > nzt )  THEN
781                                  pos_x = pos_x + pdx(i)
782                                  CYCLE xloop
783                               ENDIF
784                            ENDIF
785
[1359]786                            local_count(kp,jp,ip) = local_count(kp,jp,ip) + 1
787                            IF ( .NOT. first_stride )  THEN
788                               IF ( ip < nxl  .OR.  jp < nys  .OR.  kp < nzb+1 )  THEN
789                                  write(6,*) 'xl ',ip,jp,kp,nxl,nys,nzb+1
790                               ENDIF
791                               IF ( ip > nxr  .OR.  jp > nyn  .OR.  kp > nzt )  THEN
792                                  write(6,*) 'xu ',ip,jp,kp,nxr,nyn,nzt
793                               ENDIF
794                               grid_particles(kp,jp,ip)%particles(local_count(kp,jp,ip)) = tmp_particle
795                            ENDIF
[1]796                         ENDDO
797
798                      ENDIF
799
800                      pos_x = pos_x + pdx(i)
801
[1575]802                   ENDDO xloop
[1]803
804                ENDIF
805
806                pos_y = pos_y + pdy(i)
807
808             ENDDO
809
810             pos_z = pos_z + pdz(i)
811
812          ENDDO
813
814       ENDDO
815
[1359]816       IF ( first_stride )  THEN
817          DO  ip = nxl, nxr
818             DO  jp = nys, nyn
819                DO  kp = nzb+1, nzt
820                   IF ( phase == PHASE_INIT )  THEN
821                      IF ( local_count(kp,jp,ip) > 0 )  THEN
822                         alloc_size = MAX( INT( local_count(kp,jp,ip) *        &
823                            ( 1.0_wp + alloc_factor / 100.0_wp ) ),            &
824                            min_nr_particle )
825                      ELSE
826                         alloc_size = min_nr_particle
827                      ENDIF
828                      ALLOCATE(grid_particles(kp,jp,ip)%particles(1:alloc_size))
829                      DO  n = 1, alloc_size
830                         grid_particles(kp,jp,ip)%particles(n) = zero_particle
831                      ENDDO
832                   ELSEIF ( phase == PHASE_RELEASE )  THEN
833                      IF ( local_count(kp,jp,ip) > 0 )  THEN
834                         new_size   = local_count(kp,jp,ip) + prt_count(kp,jp,ip)
835                         alloc_size = MAX( INT( new_size * ( 1.0_wp +          &
836                            alloc_factor / 100.0_wp ) ), min_nr_particle )
837                         IF( alloc_size > SIZE( grid_particles(kp,jp,ip)%particles) )  THEN
838                           CALL realloc_particles_array(ip,jp,kp,alloc_size)
839                         ENDIF
840                      ENDIF
841                   ENDIF
842                ENDDO
843             ENDDO
844          ENDDO
845       ENDIF
846    ENDDO
[1]847
[1359]848    local_start = prt_count+1
849    prt_count   = local_count
[1]850!
[1359]851!-- Add random fluctuation to particle positions
852    IF ( random_start_position )  THEN
853       DO  ip = nxl, nxr
854          DO  jp = nys, nyn
855             DO  kp = nzb+1, nzt
856                number_of_particles = prt_count(kp,jp,ip)
857                IF ( number_of_particles <= 0 )  CYCLE
858                particles => grid_particles(kp,jp,ip)%particles(1:number_of_particles)
[1]859
[1359]860                DO  n = local_start(kp,jp,ip), number_of_particles              !Move only new particles
861                   IF ( psl(particles(n)%group) /= psr(particles(n)%group) )  THEN
862                      particles(n)%x = particles(n)%x +                        &
863                                   ( random_function( iran_part ) - 0.5_wp ) * &
864                                   pdx(particles(n)%group)
865                   ENDIF
866                   IF ( pss(particles(n)%group) /= psn(particles(n)%group) )  THEN
867                      particles(n)%y = particles(n)%y +                        &
868                                   ( random_function( iran_part ) - 0.5_wp ) * &
869                                   pdy(particles(n)%group)
870                   ENDIF
871                   IF ( psb(particles(n)%group) /= pst(particles(n)%group) )  THEN
872                      particles(n)%z = particles(n)%z +                        &
873                                   ( random_function( iran_part ) - 0.5_wp ) * &
874                                   pdz(particles(n)%group)
875                   ENDIF
876                ENDDO
[1]877!
[1359]878!--             Identify particles located outside the model domain
879                CALL lpm_boundary_conds( 'bottom/top' )
880             ENDDO
881          ENDDO
882       ENDDO
[1]883!
[1359]884!--    Exchange particles between grid cells and processors
885       CALL lpm_move_particle
886       CALL lpm_exchange_horiz
[1]887
[1359]888    ENDIF
[1]889!
[1359]890!-- In case of random_start_position, delete particles identified by
891!-- lpm_exchange_horiz and lpm_boundary_conds. Then sort particles into blocks,
892!-- which is needed for a fast interpolation of the LES fields on the particle
893!-- position.
894    CALL lpm_pack_all_arrays
[1]895
896!
[1359]897!-- Determine maximum number of particles (i.e., all possible particles that
898!-- have been allocated) and the current number of particles
899    DO  ip = nxl, nxr
900       DO  jp = nys, nyn
901          DO  kp = nzb+1, nzt
902             maximum_number_of_particles = maximum_number_of_particles         &
903                                           + SIZE(grid_particles(kp,jp,ip)%particles)
904             number_of_particles         = number_of_particles                 &
905                                           + prt_count(kp,jp,ip)
[1]906          ENDDO
[1359]907       ENDDO
908    ENDDO
[1]909!
[1359]910!-- Calculate the number of particles and tails of the total domain
[1]911#if defined( __parallel )
[1359]912    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
913    CALL MPI_ALLREDUCE( number_of_particles, total_number_of_particles, 1, &
914    MPI_INTEGER, MPI_SUM, comm2d, ierr )
915    IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
916    CALL MPI_ALLREDUCE( number_of_tails, total_number_of_tails, 1, &
917    MPI_INTEGER, MPI_SUM, comm2d, ierr )
[1]918#else
[1359]919    total_number_of_particles = number_of_particles
920    total_number_of_tails     = number_of_tails
[1]921#endif
922
[1359]923    RETURN
[1]924
[1359]925 END SUBROUTINE lpm_create_particle
[336]926
[1359]927END MODULE lpm_init_mod
Note: See TracBrowser for help on using the repository browser.