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

Last change on this file since 1614 was 1576, checked in by raasch, 9 years ago

last commit documented

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