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

Last change on this file since 1806 was 1784, checked in by raasch, 9 years ago

last commit documented

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