source: palm/trunk/SOURCE/lpm_sort_arrays.f90 @ 1282

Last change on this file since 1282 was 1037, checked in by raasch, 12 years ago

last commit documented

  • Property svn:keywords set to Id
File size: 4.1 KB
Line 
1 SUBROUTINE lpm_sort_arrays
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-2012  Leibniz University Hannover
18!--------------------------------------------------------------------------------!
19!
20! Current revisions:
21! ------------------
22!
23!
24! Former revisions:
25! -----------------
26! $Id: lpm_sort_arrays.f90 1037 2012-10-22 14:10:22Z raasch $
27!
28! 1036 2012-10-22 13:43:42Z raasch
29! code put under GPL (PALM 3.9)
30!
31! 849 2012-03-15 10:35:09Z raasch
32! initial revision (former part of advec_particles)
33!
34!
35! Description:
36! ------------
37! Sort particles in the sequence the grid boxes are stored in memory.
38!------------------------------------------------------------------------------!
39
40    USE arrays_3d
41    USE control_parameters
42    USE cpulog
43    USE grid_variables
44    USE indices
45    USE interfaces
46    USE particle_attributes
47
48    IMPLICIT NONE
49
50    INTEGER ::  i, ilow, j, k, n
51
52    TYPE(particle_type), DIMENSION(:), POINTER ::  particles_temp
53
54
55    CALL cpu_log( log_point_s(47), 'lpm_sort_arrays', 'start' )
56
57!
58!-- Initialize counters and set pointer of the temporary array into which
59!-- particles are sorted to free memory
60    prt_count  = 0
61    sort_count = sort_count +1
62
63    SELECT CASE ( MOD( sort_count, 2 ) )
64
65       CASE ( 0 )
66
67          particles_temp => part_1
68
69       CASE ( 1 )
70
71          particles_temp => part_2
72
73    END SELECT
74
75!
76!-- Count the particles per gridbox
77    DO  n = 1, number_of_particles
78
79       i = ( particles(n)%x + 0.5 * dx ) * ddx
80       j = ( particles(n)%y + 0.5 * dy ) * ddy
81       k = particles(n)%z / dz + 1 + offset_ocean_nzt
82           ! only exact if equidistant
83
84       prt_count(k,j,i) = prt_count(k,j,i) + 1
85
86       IF ( i < nxl .OR. i > nxr .OR. j < nys .OR. j > nyn .OR. k < nzb+1 .OR. &
87            k > nzt )  THEN
88          WRITE( message_string, * ) ' particle out of range: i=', i, ' j=', &
89                          j, ' k=', k,                                       &
90                          ' nxl=', nxl, ' nxr=', nxr,                        &
91                          ' nys=', nys, ' nyn=', nyn,                        &
92                          ' nzb=', nzb, ' nzt=', nzt
93         CALL message( 'lpm_sort_arrays', 'PA0149', 1, 2, 0, 6, 0 ) 
94       ENDIF
95
96    ENDDO
97
98!
99!-- Calculate the lower indices of those ranges of the particles-array
100!-- containing particles which belong to the same gridpox i,j,k
101    ilow = 1
102    DO  i = nxl, nxr
103       DO  j = nys, nyn
104          DO  k = nzb+1, nzt
105             prt_start_index(k,j,i) = ilow
106             ilow = ilow + prt_count(k,j,i)
107          ENDDO
108       ENDDO
109    ENDDO
110
111!
112!-- Sorting the particles
113    DO  n = 1, number_of_particles
114
115       i = ( particles(n)%x + 0.5 * dx ) * ddx
116       j = ( particles(n)%y + 0.5 * dy ) * ddy
117       k = particles(n)%z / dz + 1 + offset_ocean_nzt
118           ! only exact if equidistant
119
120       particles_temp(prt_start_index(k,j,i)) = particles(n)
121
122       prt_start_index(k,j,i) = prt_start_index(k,j,i) + 1
123
124    ENDDO
125
126!
127!-- Redirect the particles pointer to the sorted array
128    SELECT CASE ( MOD( sort_count, 2 ) )
129
130       CASE ( 0 )
131
132          particles => part_1
133
134       CASE ( 1 )
135
136          particles => part_2
137
138    END SELECT
139
140!
141!-- Reset the index array to the actual start position
142    DO  i = nxl, nxr
143       DO  j = nys, nyn
144          DO  k = nzb+1, nzt
145             prt_start_index(k,j,i) = prt_start_index(k,j,i) - prt_count(k,j,i)
146          ENDDO
147       ENDDO
148    ENDDO
149
150    CALL cpu_log( log_point_s(47), 'lpm_sort_arrays', 'stop' )
151
152
153 END SUBROUTINE lpm_sort_arrays
Note: See TracBrowser for help on using the repository browser.