source: palm/trunk/SOURCE/lpm_read_restart_file.f90 @ 1682

Last change on this file since 1682 was 1682, checked in by knoop, 9 years ago

Code annotations made doxygen readable

  • Property svn:keywords set to Id
File size: 7.6 KB
RevLine 
[1682]1!> @file lpm_read_restart_file.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!
[1310]16! Copyright 1997-2014 Leibniz Universitaet Hannover
[1036]17!--------------------------------------------------------------------------------!
18!
[849]19! Current revisions:
20! ------------------
[1682]21! Code annotations made doxygen readable
[1360]22!
[849]23! Former revisions:
24! -----------------
25! $Id: lpm_read_restart_file.f90 1682 2015-10-07 23:56:08Z knoop $
26!
[1360]27! 1359 2014-04-11 17:15:14Z hoffmann
28! New particle structure integrated.
29!
[1321]30! 1320 2014-03-20 08:40:49Z raasch
31! ONLY-attribute added to USE-statements,
32! comment fields (!:) to be used for variable explanations added to
33! all variable declaration statements
34!
[1037]35! 1036 2012-10-22 13:43:42Z raasch
36! code put under GPL (PALM 3.9)
37!
[850]38! 849 2012-03-15 10:35:09Z raasch
39! initial revision (former part of init_particles)
[849]40!
[850]41!
[849]42! Description:
43! ------------
[1682]44!> Read particle data from the restart file.
[849]45!------------------------------------------------------------------------------!
[1682]46 SUBROUTINE lpm_read_restart_file
47 
[849]48
[1320]49    USE control_parameters,                                                    &
50        ONLY:  message_string
51
52    USE indices,                                                               &
53        ONLY:  nx, nxl, nxlg, nxr, nxrg, ny, nyn, nyng, nys, nysg, nz, nzb, nzt
54
[1359]55    USE kinds
56
57    USE lpm_pack_arrays_mod,                                                   &
58        ONLY:  lpm_pack_all_arrays
59
[1320]60    USE particle_attributes,                                                   &
[1359]61        ONLY:  alloc_factor, bc_par_b, bc_par_lr, bc_par_ns, bc_par_t,         &
62               grid_particles, maximum_number_of_particles,                    &
[1320]63               maximum_number_of_tailpoints, maximum_number_of_tails,          &
[1359]64               min_nr_particle, new_tail_id, number_of_particles,              &
65               number_of_particle_groups, number_of_tails, particles,          &
66               particle_groups, particle_tail_coordinates, particle_type,      &
67               prt_count, sort_count, tail_mask, time_prel,                    &
68               time_write_particle_data, uniform_particles,                    &
69               use_particle_tails, zero_particle
[1320]70
[849]71    USE pegrid
72
73    IMPLICIT NONE
74
[1682]75    CHARACTER (LEN=10) ::  particle_binary_version    !<
76    CHARACTER (LEN=10) ::  version_on_file            !<
[849]77
[1682]78    INTEGER(iwp) :: alloc_size !<
79    INTEGER(iwp) :: ip         !<
80    INTEGER(iwp) :: jp         !<
81    INTEGER(iwp) :: kp         !<
[1359]82
[1682]83    TYPE(particle_type), DIMENSION(:), ALLOCATABLE :: tmp_particles !<
[1359]84
[849]85!
86!-- Read particle data from previous model run.
87!-- First open the input unit.
88    IF ( myid_char == '' )  THEN
89       OPEN ( 90, FILE='PARTICLE_RESTART_DATA_IN'//myid_char, &
90                  FORM='UNFORMATTED' )
91    ELSE
92       OPEN ( 90, FILE='PARTICLE_RESTART_DATA_IN/'//myid_char, &
93                  FORM='UNFORMATTED' )
94    ENDIF
95
96!
97!-- First compare the version numbers
98    READ ( 90 )  version_on_file
[1359]99    particle_binary_version = '3.2'
[849]100    IF ( TRIM( version_on_file ) /= TRIM( particle_binary_version ) )  THEN
101       message_string = 'version mismatch concerning data from prior ' // &
102                        'run &version on file    = "' //                  &
103                                      TRIM( version_on_file ) //          &
104                        '&version in program = "' //                      &
105                                      TRIM( particle_binary_version ) // '"'
106       CALL message( 'lpm_read_restart_file', 'PA0214', 1, 2, 0, 6, 0 )
107    ENDIF
108
109!
[1359]110!-- If less particles are stored on the restart file than prescribed by
111!-- min_nr_particle, the remainder is initialized by zero_particle to avoid
112!-- errors.
113#if defined( __twocachelines )
114    zero_particle = particle_type( 0.0_wp, 0.0_sp, 0.0_sp, 0.0_sp, 0.0_sp,  &
115                                   0.0_sp, 0.0_sp, 0.0_wp, 0.0_wp, 0.0_wp,  &
116                                   0, .FALSE., 0.0_wp, 0.0_wp, 0.0_wp,      &
117                                   0.0_sp, 0.0_sp, 0.0_sp, 0.0_sp, 0.0_sp,  &
118                                   0.0_sp, 0, 0, 0, -1)
119#else
120    zero_particle = particle_type( 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp,  &
121                                   0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp,  &
122                                   0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp,  &
123                                   0.0_wp, 0.0_wp, 0.0_wp, 0.0_wp, 0, 0, 0, &
124                                   0, .FALSE., -1)
125#endif
126
127!
[849]128!-- Read some particle parameters and the size of the particle arrays,
129!-- allocate them and read their contents.
130    READ ( 90 )  bc_par_b, bc_par_lr, bc_par_ns, bc_par_t,                  &
[1359]131                 maximum_number_of_tailpoints, maximum_number_of_tails,     &
132                 number_of_particle_groups, number_of_tails,                &
133                 particle_groups, time_prel, time_write_particle_data,      &
134                 uniform_particles
[849]135
[1359]136    ALLOCATE( prt_count(nzb:nzt+1,nysg:nyng,nxlg:nxrg),                     &
137              grid_particles(nzb:nzt+1,nysg:nyng,nxlg:nxrg) )
[849]138
[1359]139    READ ( 90 )  prt_count
[849]140
[1359]141    maximum_number_of_particles = 0
142    DO  ip = nxl, nxr
143       DO  jp = nys, nyn
144          DO  kp = nzb+1, nzt
[849]145
[1359]146             number_of_particles = prt_count(kp,jp,ip)
147             IF ( number_of_particles > 0 )  THEN
148                alloc_size = MAX( INT( number_of_particles *                   &
149                             ( 1.0_wp + alloc_factor / 100.0_wp ) ),           &
150                             min_nr_particle )
151             ELSE
152                alloc_size = min_nr_particle
153             ENDIF
[849]154
[1359]155             ALLOCATE( grid_particles(kp,jp,ip)%particles(1:alloc_size) )
[849]156
[1359]157             IF ( number_of_particles > 0 )  THEN
158                ALLOCATE( tmp_particles(1:number_of_particles) )
159                READ ( 90 )  tmp_particles
160                grid_particles(kp,jp,ip)%particles(1:number_of_particles) = tmp_particles
161                DEALLOCATE( tmp_particles )
162                IF ( number_of_particles < alloc_size )  THEN
163                   grid_particles(kp,jp,ip)%particles(number_of_particles+1:alloc_size) &
164                      = zero_particle
165                ENDIF
166             ELSE
167                grid_particles(kp,jp,ip)%particles(1:alloc_size) = zero_particle
168             ENDIF
[849]169
[1359]170             maximum_number_of_particles = maximum_number_of_particles + alloc_size
[849]171
[1359]172          ENDDO
173       ENDDO
174    ENDDO
[849]175
[1359]176!
177!-- particle tails currently not available
178!    IF ( use_particle_tails )  THEN
179!       ALLOCATE( particle_tail_coordinates(maximum_number_of_tailpoints,5, &
180!                 maximum_number_of_tails),                                 &
181!                 new_tail_id(maximum_number_of_tails),                     &
182!                 tail_mask(maximum_number_of_tails) )
183!       READ ( 90 )  particle_tail_coordinates
184!    ENDIF
185
[849]186    CLOSE ( 90 )
[1359]187!
188!-- Must be called to sort particles into blocks, which is needed for a fast
189!-- interpolation of the LES fields on the particle position.
190    CALL lpm_pack_all_arrays
[849]191
192
193 END SUBROUTINE lpm_read_restart_file
Note: See TracBrowser for help on using the repository browser.