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

Last change on this file since 1704 was 1692, checked in by maronga, 9 years ago

last commit documented

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