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

Last change on this file since 1327 was 1327, checked in by raasch, 10 years ago

Changed:


-s real64 removed (.mrun.config.hlrnIII)
-r8 removed (.mrun.config.imuk)
deleted: .mrun.config.imuk_ice2_netcdf4 .mrun.config.imuk_hlrn

REAL constants defined as wp-kind in modules

"baroclinicity" renamed "baroclinity", "ocean version" replaced by
"ocean mode"

code parts concerning old output formats "iso2d" and "avs" removed.
netCDF is the only remaining output format.

Errors:


bugfix: duplicate error message 56 removed

  • Property svn:keywords set to Id
File size: 26.0 KB
Line 
1 SUBROUTINE lpm_init
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! -netcdf_output
23!
24! Former revisions:
25! -----------------
26! $Id: lpm_init.f90 1327 2014-03-21 11:00:16Z raasch $
27!
28! 1322 2014-03-20 16:38:49Z raasch
29! REAL functions provided with KIND-attribute
30!
31! 1320 2014-03-20 08:40:49Z raasch
32! ONLY-attribute added to USE-statements,
33! kind-parameters added to all INTEGER and REAL declaration statements,
34! kinds are defined in new module kinds,
35! revision history before 2012 removed,
36! comment fields (!:) to be used for variable explanations added to
37! all variable declaration statements
38! bugfix: #if defined( __parallel ) added
39!
40! 1314 2014-03-14 18:25:17Z suehring
41! Vertical logarithmic interpolation of horizontal particle speed for particles
42! between roughness height and first vertical grid level.
43!
44! 1092 2013-02-02 11:24:22Z raasch
45! unused variables removed
46!
47! 1036 2012-10-22 13:43:42Z raasch
48! code put under GPL (PALM 3.9)
49!
50! 849 2012-03-15 10:35:09Z raasch
51! routine renamed: init_particles -> lpm_init
52! de_dx, de_dy, de_dz are allocated here (instead of automatic arrays in
53! advec_particles),
54! sort_particles renamed lpm_sort_arrays, user_init_particles renamed lpm_init
55!
56! 828 2012-02-21 12:00:36Z raasch
57! call of init_kernels, particle feature color renamed class
58!
59! 824 2012-02-17 09:09:57Z raasch
60! particle attributes speed_x|y|z_sgs renamed rvar1|2|3,
61! array particles implemented as pointer
62!
63! 667 2010-12-23 12:06:00Z suehring/gryschka
64! nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng for allocation
65! of arrays.
66!
67! Revision 1.1  1999/11/25 16:22:38  raasch
68! Initial revision
69!
70!
71! Description:
72! ------------
73! This routine initializes a set of particles and their attributes (position,
74! radius, ..) which are used by the Lagrangian particle model (see lpm).
75!------------------------------------------------------------------------------!
76
77    USE arrays_3d,                                                             &
78        ONLY:  de_dx, de_dy, de_dz, zu, zw, z0
79
80    USE cloud_parameters,                                                      &
81        ONLY:  curvature_solution_effects
82
83    USE control_parameters,                                                    &
84        ONLY:  cloud_droplets, current_timestep_number, initializing_actions,  &
85               message_string, netcdf_data_format, ocean,                      &
86               prandtl_layer, simulated_time
87
88    USE dvrp_variables,                                                        &
89        ONLY:  particle_color, particle_dvrpsize
90
91    USE grid_variables,                                                        &
92        ONLY:  dx, dy
93
94    USE indices,                                                               &
95        ONLY:  nx, nxl, nxlg, nxrg, nxr, ny, nyn, nys, nyng, nysg, nz, nzb, nzt
96
97    USE kinds
98
99    USE lpm_collision_kernels_mod,                                             &
100        ONLY:  init_kernels
101
102    USE particle_attributes,                                                   &
103        ONLY:   bc_par_b, bc_par_lr, bc_par_ns, bc_par_t, collision_kernel,    &
104                density_ratio, dvrp_psize, initial_weighting_factor, ibc_par_b,&
105                ibc_par_lr, ibc_par_ns, ibc_par_t, initial_particles,          &
106                iran_part, log_z_z0, max_number_of_particle_groups,            &
107                maximum_number_of_particles, maximum_number_of_tailpoints,     &
108                minimum_tailpoint_distance, maximum_number_of_tails,           &
109                mpi_particle_type, new_tail_id, number_of_initial_particles,   &
110                number_of_initial_tails, number_of_particles,                  &
111                number_of_particle_groups, number_of_sublayers,                &
112                number_of_tails, offset_ocean_nzt, offset_ocean_nzt_m1, part_1,&
113                part_2, particles, particle_advection_start, particle_groups,  &
114                particle_groups_type, particle_mask, particles_per_point,      &
115                particle_tail_coordinates,  particle_type, pdx, pdy, pdz,      &
116                prt_count, prt_start_index, psb, psl, psn, psr, pss, pst,      &
117                radius, random_start_position, read_particles_from_restartfile,&
118                skip_particles_for_tail, sort_count, tail_mask,                &
119                total_number_of_particles, total_number_of_tails,              &
120                use_particle_tails, use_sgs_for_particles,                     &
121                write_particle_statistics, uniform_particles, z0_av_global
122
123    USE pegrid
124
125    USE random_function_mod,                                                   &
126        ONLY:  random_function
127
128
129    IMPLICIT NONE
130
131    INTEGER(iwp) ::  i                           !:
132    INTEGER(iwp) ::  j                           !:
133    INTEGER(iwp) ::  k                           !:
134    INTEGER(iwp) ::  n                           !:
135    INTEGER(iwp) ::  nn                          !:
136
137#if defined( __parallel )
138    INTEGER(iwp), DIMENSION(3) ::  blocklengths  !:
139    INTEGER(iwp), DIMENSION(3) ::  displacements !:
140    INTEGER(iwp), DIMENSION(3) ::  types         !:
141#endif
142    LOGICAL ::  uniform_particles_l              !:
143
144    REAL(wp)    ::  height_int                   !:
145    REAL(wp)    ::  height_p                     !:
146    REAL(wp)    ::  pos_x                        !:
147    REAL(wp)    ::  pos_y                        !:
148    REAL(wp)    ::  pos_z                        !:
149    REAL(wp)    ::  z_p                          !:
150    REAL(wp)    ::  z0_av_local = 0.0            !:
151
152               
153
154
155#if defined( __parallel )
156!
157!-- Define MPI derived datatype for FORTRAN datatype particle_type (see module
158!-- particle_attributes). Integer length is 4 byte, Real is 8 byte
159    blocklengths(1)  = 19;  blocklengths(2)  =   4;  blocklengths(3)  =   1
160    displacements(1) =  0;  displacements(2) = 152;  displacements(3) = 168
161
162    types(1) = MPI_REAL
163    types(2) = MPI_INTEGER
164    types(3) = MPI_UB
165    CALL MPI_TYPE_STRUCT( 3, blocklengths, displacements, types, &
166                          mpi_particle_type, ierr )
167    CALL MPI_TYPE_COMMIT( mpi_particle_type, ierr )
168#endif
169
170!
171!-- In case of oceans runs, the vertical index calculations need an offset,
172!-- because otherwise the k indices will become negative
173    IF ( ocean )  THEN
174       offset_ocean_nzt    = nzt
175       offset_ocean_nzt_m1 = nzt - 1
176    ENDIF
177
178
179!
180!-- Check the number of particle groups.
181    IF ( number_of_particle_groups > max_number_of_particle_groups )  THEN
182       WRITE( message_string, * ) 'max_number_of_particle_groups =',      &
183                                  max_number_of_particle_groups ,         &
184                                  '&number_of_particle_groups reset to ', &
185                                  max_number_of_particle_groups
186       CALL message( 'lpm_init', 'PA0213', 0, 1, 0, 6, 0 )
187       number_of_particle_groups = max_number_of_particle_groups
188    ENDIF
189
190!
191!-- Set default start positions, if necessary
192    IF ( psl(1) == 9999999.9 )  psl(1) = -0.5 * dx
193    IF ( psr(1) == 9999999.9 )  psr(1) = ( nx + 0.5 ) * dx
194    IF ( pss(1) == 9999999.9 )  pss(1) = -0.5 * dy
195    IF ( psn(1) == 9999999.9 )  psn(1) = ( ny + 0.5 ) * dy
196    IF ( psb(1) == 9999999.9 )  psb(1) = zu(nz/2)
197    IF ( pst(1) == 9999999.9 )  pst(1) = psb(1)
198
199    IF ( pdx(1) == 9999999.9  .OR.  pdx(1) == 0.0 )  pdx(1) = dx
200    IF ( pdy(1) == 9999999.9  .OR.  pdy(1) == 0.0 )  pdy(1) = dy
201    IF ( pdz(1) == 9999999.9  .OR.  pdz(1) == 0.0 )  pdz(1) = zu(2) - zu(1)
202
203    DO  j = 2, number_of_particle_groups
204       IF ( psl(j) == 9999999.9 )  psl(j) = psl(j-1)
205       IF ( psr(j) == 9999999.9 )  psr(j) = psr(j-1)
206       IF ( pss(j) == 9999999.9 )  pss(j) = pss(j-1)
207       IF ( psn(j) == 9999999.9 )  psn(j) = psn(j-1)
208       IF ( psb(j) == 9999999.9 )  psb(j) = psb(j-1)
209       IF ( pst(j) == 9999999.9 )  pst(j) = pst(j-1)
210       IF ( pdx(j) == 9999999.9  .OR.  pdx(j) == 0.0 )  pdx(j) = pdx(j-1)
211       IF ( pdy(j) == 9999999.9  .OR.  pdy(j) == 0.0 )  pdy(j) = pdy(j-1)
212       IF ( pdz(j) == 9999999.9  .OR.  pdz(j) == 0.0 )  pdz(j) = pdz(j-1)
213    ENDDO
214
215!
216!-- Allocate arrays required for calculating particle SGS velocities
217    IF ( use_sgs_for_particles )  THEN
218       ALLOCATE( de_dx(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
219                 de_dy(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
220                 de_dz(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
221    ENDIF
222
223!
224!-- Allocate array required for logarithmic vertical interpolation of
225!-- horizontal particle velocities between the surface and the first vertical
226!-- grid level. In order to avoid repeated CPU cost-intensive CALLS of
227!-- intrinsic FORTRAN procedure LOG(z/z0), LOG(z/z0) is precalculated for
228!-- several heights. Splitting into 20 sublayers turned out to be sufficient.
229!-- To obtain exact height levels of particles, linear interpolation is applied
230!-- (see lpm_advec.f90).
231    IF ( prandtl_layer )  THEN
232       
233       ALLOCATE ( log_z_z0(0:number_of_sublayers) ) 
234       z_p         = zu(nzb+1) - zw(nzb)
235
236!
237!--    Calculate horizontal mean value of z0 used for logartihmic
238!--    interpolation. Note: this is not exact for heterogeneous z0.
239!--    However, sensitivity studies showed that the effect is
240!--    negligible.
241       z0_av_local  = SUM( z0(nys:nyn,nxl:nxr) )
242       z0_av_global = 0.0
243
244#if defined( __parallel )
245       CALL MPI_ALLREDUCE(z0_av_local, z0_av_global, 1, MPI_REAL, MPI_SUM, &
246                          comm2d, ierr )
247#else
248       z0_av_global = z0_av_local
249#endif
250
251       z0_av_global = z0_av_global  / ( ( ny + 1 ) * ( nx + 1 ) )
252!
253!--    Horizontal wind speed is zero below and at z0
254       log_z_z0(0) = 0.0   
255!
256!--    Calculate vertical depth of the sublayers
257       height_int  = ( z_p - z0_av_global ) / REAL( number_of_sublayers, KIND=wp )
258!
259!--    Precalculate LOG(z/z0)
260       height_p    = 0.0
261       DO  k = 1, number_of_sublayers
262
263          height_p    = height_p + height_int
264          log_z_z0(k) = LOG( height_p / z0_av_global )
265
266       ENDDO
267
268
269    ENDIF
270
271!
272!-- Initialize collision kernels
273    IF ( collision_kernel /= 'none' )  CALL init_kernels
274
275!
276!-- For the first model run of a possible job chain initialize the
277!-- particles, otherwise read the particle data from restart file.
278    IF ( TRIM( initializing_actions ) == 'read_restart_data'  &
279         .AND.  read_particles_from_restartfile )  THEN
280
281       CALL lpm_read_restart_file
282
283    ELSE
284
285!
286!--    Allocate particle arrays and set attributes of the initial set of
287!--    particles, which can be also periodically released at later times.
288!--    Also allocate array for particle tail coordinates, if needed.
289       ALLOCATE( prt_count(nzb:nzt+1,nysg:nyng,nxlg:nxrg),       &
290                 prt_start_index(nzb:nzt+1,nysg:nyng,nxlg:nxrg), &
291                 particle_mask(maximum_number_of_particles),     &
292                 part_1(maximum_number_of_particles),            &
293                 part_2(maximum_number_of_particles)  )
294
295       particles => part_1
296
297       sort_count = 0
298
299!
300!--    Initialize all particles with dummy values (otherwise errors may
301!--    occur within restart runs). The reason for this is still not clear
302!--    and may be presumably caused by errors in the respective user-interface.
303       particles = particle_type( 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, &
304                                  0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, &
305                                  0.0, 0, 0, 0, 0 )
306       particle_groups = particle_groups_type( 0.0, 0.0, 0.0, 0.0 )
307
308!
309!--    Set the default particle size used for dvrp plots
310       IF ( dvrp_psize == 9999999.9 )  dvrp_psize = 0.2 * dx
311
312!
313!--    Set values for the density ratio and radius for all particle
314!--    groups, if necessary
315       IF ( density_ratio(1) == 9999999.9 )  density_ratio(1) = 0.0
316       IF ( radius(1)        == 9999999.9 )  radius(1) = 0.0
317       DO  i = 2, number_of_particle_groups
318          IF ( density_ratio(i) == 9999999.9 )  THEN
319             density_ratio(i) = density_ratio(i-1)
320          ENDIF
321          IF ( radius(i) == 9999999.9 )  radius(i) = radius(i-1)
322       ENDDO
323
324       DO  i = 1, number_of_particle_groups
325          IF ( density_ratio(i) /= 0.0  .AND.  radius(i) == 0 )  THEN
326             WRITE( message_string, * ) 'particle group #', i, 'has a', &
327                                        'density ratio /= 0 but radius = 0'
328             CALL message( 'lpm_init', 'PA0215', 1, 2, 0, 6, 0 )
329          ENDIF
330          particle_groups(i)%density_ratio = density_ratio(i)
331          particle_groups(i)%radius        = radius(i)
332       ENDDO
333
334!
335!--    Calculate particle positions and store particle attributes, if
336!--    particle is situated on this PE
337       n = 0
338
339       DO  i = 1, number_of_particle_groups
340
341          pos_z = psb(i)
342
343          DO WHILE ( pos_z <= pst(i) )
344
345             pos_y = pss(i)
346
347             DO WHILE ( pos_y <= psn(i) )
348
349                IF ( pos_y >= ( nys - 0.5 ) * dy  .AND.  &
350                     pos_y <  ( nyn + 0.5 ) * dy )  THEN
351
352                   pos_x = psl(i)
353
354                   DO WHILE ( pos_x <= psr(i) )
355
356                      IF ( pos_x >= ( nxl - 0.5 ) * dx  .AND.  &
357                           pos_x <  ( nxr + 0.5 ) * dx )  THEN
358
359                         DO  j = 1, particles_per_point
360
361                            n = n + 1
362                            IF ( n > maximum_number_of_particles )  THEN
363                               WRITE( message_string, * ) 'number of initial', &
364                                      'particles (', n, ') exceeds',           &
365                                      '&maximum_number_of_particles (',        &
366                                      maximum_number_of_particles, ') on PE ', &
367                                             myid
368                               CALL message( 'lpm_init', 'PA0216', 2, 2, -1, 6,&
369                                             1 )
370                            ENDIF
371                            particles(n)%x             = pos_x
372                            particles(n)%y             = pos_y
373                            particles(n)%z             = pos_z
374                            particles(n)%age           = 0.0
375                            particles(n)%age_m         = 0.0
376                            particles(n)%dt_sum        = 0.0
377                            particles(n)%dvrp_psize    = dvrp_psize
378                            particles(n)%e_m           = 0.0
379                            IF ( curvature_solution_effects )  THEN
380!
381!--                            Initial values (internal timesteps, derivative)
382!--                            for Rosenbrock method
383                               particles(n)%rvar1      = 1.0E-12
384                               particles(n)%rvar2      = 1.0E-3
385                               particles(n)%rvar3      = -9999999.9
386                            ELSE
387!
388!--                            Initial values for SGS velocities
389                               particles(n)%rvar1      = 0.0
390                               particles(n)%rvar2      = 0.0
391                               particles(n)%rvar3      = 0.0
392                            ENDIF
393                            particles(n)%speed_x       = 0.0
394                            particles(n)%speed_y       = 0.0
395                            particles(n)%speed_z       = 0.0
396                            particles(n)%origin_x      = pos_x
397                            particles(n)%origin_y      = pos_y
398                            particles(n)%origin_z      = pos_z
399                            particles(n)%radius      = particle_groups(i)%radius
400                            particles(n)%weight_factor =initial_weighting_factor
401                            particles(n)%class         = 1
402                            particles(n)%group         = i
403                            particles(n)%tailpoints    = 0
404                            IF ( use_particle_tails  .AND. &
405                                 MOD( n, skip_particles_for_tail ) == 0 )  THEN
406                               number_of_tails         = number_of_tails + 1
407!
408!--                            This is a temporary provisional setting (see
409!--                            further below!)
410                               particles(n)%tail_id    = number_of_tails
411                            ELSE
412                               particles(n)%tail_id    = 0
413                            ENDIF
414
415                         ENDDO
416
417                      ENDIF
418
419                      pos_x = pos_x + pdx(i)
420
421                   ENDDO
422
423                ENDIF
424
425                pos_y = pos_y + pdy(i)
426
427             ENDDO
428
429             pos_z = pos_z + pdz(i)
430
431          ENDDO
432
433       ENDDO
434
435       number_of_initial_particles = n
436       number_of_particles         = n
437
438!
439!--    Calculate the number of particles and tails of the total domain
440#if defined( __parallel )
441       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
442       CALL MPI_ALLREDUCE( number_of_particles, total_number_of_particles, 1, &
443                           MPI_INTEGER, MPI_SUM, comm2d, ierr )
444       IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
445       CALL MPI_ALLREDUCE( number_of_tails, total_number_of_tails, 1, &
446                           MPI_INTEGER, MPI_SUM, comm2d, ierr )
447#else
448       total_number_of_particles = number_of_particles
449       total_number_of_tails     = number_of_tails
450#endif
451
452!
453!--    Set a seed value for the random number generator to be exclusively
454!--    used for the particle code. The generated random numbers should be
455!--    different on the different PEs.
456       iran_part = iran_part + myid
457
458!
459!--    User modification of initial particles
460       CALL user_lpm_init
461
462!
463!--    Store the initial set of particles for release at later times
464       IF ( number_of_initial_particles /= 0 )  THEN
465          ALLOCATE( initial_particles(1:number_of_initial_particles) )
466          initial_particles(1:number_of_initial_particles) = &
467                                        particles(1:number_of_initial_particles)
468       ENDIF
469
470!
471!--    Add random fluctuation to particle positions
472       IF ( random_start_position )  THEN
473
474          DO  n = 1, number_of_initial_particles
475             IF ( psl(particles(n)%group) /= psr(particles(n)%group) )  THEN
476                particles(n)%x = particles(n)%x + &
477                                 ( random_function( iran_part ) - 0.5 ) * &
478                                 pdx(particles(n)%group)
479                IF ( particles(n)%x  <=  ( nxl - 0.5 ) * dx )  THEN
480                   particles(n)%x = ( nxl - 0.4999999999 ) * dx
481                ELSEIF ( particles(n)%x  >=  ( nxr + 0.5 ) * dx )  THEN
482                   particles(n)%x = ( nxr + 0.4999999999 ) * dx
483                ENDIF
484             ENDIF
485             IF ( pss(particles(n)%group) /= psn(particles(n)%group) )  THEN
486                particles(n)%y = particles(n)%y + &
487                                 ( random_function( iran_part ) - 0.5 ) * &
488                                 pdy(particles(n)%group)
489                IF ( particles(n)%y  <=  ( nys - 0.5 ) * dy )  THEN
490                   particles(n)%y = ( nys - 0.4999999999 ) * dy
491                ELSEIF ( particles(n)%y  >=  ( nyn + 0.5 ) * dy )  THEN
492                   particles(n)%y = ( nyn + 0.4999999999 ) * dy
493                ENDIF
494             ENDIF
495             IF ( psb(particles(n)%group) /= pst(particles(n)%group) )  THEN
496                particles(n)%z = particles(n)%z + &
497                                 ( random_function( iran_part ) - 0.5 ) * &
498                                 pdz(particles(n)%group)
499             ENDIF
500          ENDDO
501       ENDIF
502
503!
504!--    Sort particles in the sequence the gridboxes are stored in the memory.
505!--    Only required if cloud droplets are used.
506       IF ( cloud_droplets )  CALL lpm_sort_arrays
507
508!
509!--    Open file for statistical informations about particle conditions
510       IF ( write_particle_statistics )  THEN
511          CALL check_open( 80 )
512          WRITE ( 80, 8000 )  current_timestep_number, simulated_time, &
513                              number_of_initial_particles,             &
514                              maximum_number_of_particles
515          CALL close_file( 80 )
516       ENDIF
517
518!
519!--    Check if particles are really uniform in color and radius (dvrp_size)
520!--    (uniform_particles is preset TRUE)
521       IF ( uniform_particles )  THEN
522          IF ( number_of_initial_particles == 0 )  THEN
523             uniform_particles_l = .TRUE.
524          ELSE
525             n = number_of_initial_particles
526             IF ( MINVAL( particles(1:n)%dvrp_psize  ) ==     &
527                  MAXVAL( particles(1:n)%dvrp_psize  )  .AND. &
528                  MINVAL( particles(1:n)%class ) ==     &
529                  MAXVAL( particles(1:n)%class ) )  THEN
530                uniform_particles_l = .TRUE.
531             ELSE
532                uniform_particles_l = .FALSE.
533             ENDIF
534          ENDIF
535
536#if defined( __parallel )
537          IF ( collective_wait )  CALL MPI_BARRIER( comm2d, ierr )
538          CALL MPI_ALLREDUCE( uniform_particles_l, uniform_particles, 1, &
539                              MPI_LOGICAL, MPI_LAND, comm2d, ierr )
540#else
541          uniform_particles = uniform_particles_l
542#endif
543
544       ENDIF
545
546!
547!--    Particles will probably become none-uniform, if their size and color
548!--    will be determined by flow variables
549       IF ( particle_color /= 'none'  .OR.  particle_dvrpsize /= 'none' )  THEN
550          uniform_particles = .FALSE.
551       ENDIF
552
553!
554!--    Set the beginning of the particle tails and their age
555       IF ( use_particle_tails )  THEN
556!
557!--       Choose the maximum number of tails with respect to the maximum number
558!--       of particles and skip_particles_for_tail
559          maximum_number_of_tails = maximum_number_of_particles / &
560                                    skip_particles_for_tail
561
562!
563!--       Create a minimum number of tails in case that there is no tail
564!--       initially (otherwise, index errors will occur when adressing the
565!--       arrays below)
566          IF ( maximum_number_of_tails == 0 )  maximum_number_of_tails = 10
567
568          ALLOCATE( particle_tail_coordinates(maximum_number_of_tailpoints,5, &
569                    maximum_number_of_tails),                                 &
570                    new_tail_id(maximum_number_of_tails),                     &
571                    tail_mask(maximum_number_of_tails) )
572
573          particle_tail_coordinates  = 0.0
574          minimum_tailpoint_distance = minimum_tailpoint_distance**2
575          number_of_initial_tails    = number_of_tails
576
577          nn = 0
578          DO  n = 1, number_of_particles
579!
580!--          Only for those particles marked above with a provisional tail_id
581!--          tails will be created. Particles now get their final tail_id.
582             IF ( particles(n)%tail_id /= 0 )  THEN
583
584                nn = nn + 1
585                particles(n)%tail_id = nn
586
587                particle_tail_coordinates(1,1,nn) = particles(n)%x
588                particle_tail_coordinates(1,2,nn) = particles(n)%y
589                particle_tail_coordinates(1,3,nn) = particles(n)%z
590                particle_tail_coordinates(1,4,nn) = particles(n)%class
591                particles(n)%tailpoints = 1
592                IF ( minimum_tailpoint_distance /= 0.0 )  THEN
593                   particle_tail_coordinates(2,1,nn) = particles(n)%x
594                   particle_tail_coordinates(2,2,nn) = particles(n)%y
595                   particle_tail_coordinates(2,3,nn) = particles(n)%z
596                   particle_tail_coordinates(2,4,nn) = particles(n)%class
597                   particle_tail_coordinates(1:2,5,nn) = 0.0
598                   particles(n)%tailpoints = 2
599                ENDIF
600
601             ENDIF
602          ENDDO
603       ENDIF
604
605!
606!--    Plot initial positions of particles (only if particle advection is
607!--    switched on from the beginning of the simulation (t=0))
608       IF ( particle_advection_start == 0.0 )  CALL data_output_dvrp
609
610    ENDIF
611
612!
613!-- Check boundary condition and set internal variables
614    SELECT CASE ( bc_par_b )
615   
616       CASE ( 'absorb' )
617          ibc_par_b = 1
618
619       CASE ( 'reflect' )
620          ibc_par_b = 2
621         
622       CASE DEFAULT
623          WRITE( message_string, * )  'unknown boundary condition ',   &
624                                       'bc_par_b = "', TRIM( bc_par_b ), '"'
625          CALL message( 'lpm_init', 'PA0217', 1, 2, 0, 6, 0 )
626         
627    END SELECT
628    SELECT CASE ( bc_par_t )
629   
630       CASE ( 'absorb' )
631          ibc_par_t = 1
632
633       CASE ( 'reflect' )
634          ibc_par_t = 2
635         
636       CASE DEFAULT
637          WRITE( message_string, * ) 'unknown boundary condition ',   &
638                                     'bc_par_t = "', TRIM( bc_par_t ), '"'
639          CALL message( 'lpm_init', 'PA0218', 1, 2, 0, 6, 0 )
640         
641    END SELECT
642    SELECT CASE ( bc_par_lr )
643
644       CASE ( 'cyclic' )
645          ibc_par_lr = 0
646
647       CASE ( 'absorb' )
648          ibc_par_lr = 1
649
650       CASE ( 'reflect' )
651          ibc_par_lr = 2
652         
653       CASE DEFAULT
654          WRITE( message_string, * ) 'unknown boundary condition ',   &
655                                     'bc_par_lr = "', TRIM( bc_par_lr ), '"'
656          CALL message( 'lpm_init', 'PA0219', 1, 2, 0, 6, 0 )
657         
658    END SELECT
659    SELECT CASE ( bc_par_ns )
660
661       CASE ( 'cyclic' )
662          ibc_par_ns = 0
663
664       CASE ( 'absorb' )
665          ibc_par_ns = 1
666
667       CASE ( 'reflect' )
668          ibc_par_ns = 2
669         
670       CASE DEFAULT
671          WRITE( message_string, * ) 'unknown boundary condition ',   &
672                                     'bc_par_ns = "', TRIM( bc_par_ns ), '"'
673          CALL message( 'lpm_init', 'PA0220', 1, 2, 0, 6, 0 )
674         
675    END SELECT
676!
677!-- Formats
6788000 FORMAT (I6,1X,F7.2,4X,I6,71X,I6)
679
680 END SUBROUTINE lpm_init
Note: See TracBrowser for help on using the repository browser.