PartMC  2.6.1
aero_info_array.F90
Go to the documentation of this file.
1 ! Copyright (C) 2007-2012, 2017 Nicole Riemer and Matthew West
2 ! Licensed under the GNU General Public License version 2 or (at your
3 ! option) any later version. See the file COPYING for details.
4 
5 !> \file
6 !> The pmc_aero_info_array module.
7 
8 !> The aero_info_array_t structure and assoicated subroutines.
10 
11  use pmc_aero_info
12  use pmc_util
13  use pmc_spec_file
14  use pmc_mpi
15 #ifdef PMC_USE_MPI
16  use mpi
17 #endif
18 
19  !> 1-D arrays of aero_info_t structure.
20  !!
21  !! This type implements a variable-length array of aero_info_t
22  !! structures. To give a reasonable tradeoff between frequent
23  !! re-allocs and memory usage, the length of an aero_info_array is
24  !! generally a bit longer than the number of particles stored in
25  !! it. When the array is full then a larger array is allocated with
26  !! new extra space. As a balance between memory usage and frequency
27  !! of re-allocs the length of the array is currently doubled when
28  !! necessary and halved when possible.
29  !!
30  !! The true allocated length of the aero_info_array can be obtained
31  !! by size(aero_info_array%%aero_info), while the number of used
32  !! particle slots in it is given by aero_info_array_n_item().
34  !> Number of items in the array (not the same as the length of
35  !> the allocated memory).
36  integer :: n_item
37  !> Array of aero_info_t structures.
38  type(aero_info_t), allocatable :: aero_info(:)
39  end type aero_info_array_t
40 
41 contains
42 
43 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
44 
45  !> Return the current number of items.
46  elemental integer function aero_info_array_n_item(aero_info_array)
47 
48  !> Aero info array.
49  type(aero_info_array_t), intent(in) :: aero_info_array
50 
51  if (allocated(aero_info_array%aero_info)) then
52  aero_info_array_n_item = aero_info_array%n_item
53  else
55  end if
56 
57  end function aero_info_array_n_item
58 
59 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
60 
61  !> Sets an aero_info_array to contain zero data.
62  subroutine aero_info_array_zero(aero_info_array)
63 
64  !> Structure to reset.
65  type(aero_info_array_t), intent(inout) :: aero_info_array
66 
67  aero_info_array%n_item = 0
68  if (allocated(aero_info_array%aero_info)) then
69  deallocate(aero_info_array%aero_info)
70  end if
71  allocate(aero_info_array%aero_info(0))
72 
73  end subroutine aero_info_array_zero
74 
75 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
76 
77  !> Changes the given aero_info_array to exactly the given
78  !> new_length.
79  !!
80  !! This function should not be called directly, but rather use
81  !! aero_info_array_enlarge(), aero_info_array_enlarge_to()
82  !! or aero_info_array_shrink().
83  subroutine aero_info_array_realloc(aero_info_array, new_length)
84 
85  !> Array to reallocate (must already be allocated on entry).
86  type(aero_info_array_t), intent(inout) :: aero_info_array
87  !> New length of the array.
88  integer, intent(in) :: new_length
89 
90  integer :: i
91  type(aero_info_t), allocatable :: new_items(:)
92 
93  if (.not. allocated(aero_info_array%aero_info)) then
94  allocate(aero_info_array%aero_info(new_length))
95  aero_info_array%n_item = 0
96  return
97  end if
98 
99  call assert(955874877, new_length >= aero_info_array%n_item)
100  allocate(new_items(new_length))
101  do i = 1,aero_info_array%n_item
102  new_items(i) = aero_info_array%aero_info(i)
103  end do
104  call move_alloc(new_items, aero_info_array%aero_info)
105 
106  end subroutine aero_info_array_realloc
107 
108 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
109 
110  !> Possibly enlarges the given array, ensuring that it is at least of size n.
111  subroutine aero_info_array_enlarge_to(aero_info_array, n)
112 
113  !> Array to enlarge.
114  type(aero_info_array_t), intent(inout) :: aero_info_array
115  !> Minimum new size of array.
116  integer, intent(in) :: n
117 
118  if (.not. allocated(aero_info_array%aero_info)) then
119  call aero_info_array_realloc(aero_info_array, pow2_above(n))
120  return
121  end if
122 
123  if (n <= size(aero_info_array%aero_info)) return
124 
125  call aero_info_array_realloc(aero_info_array, pow2_above(n))
126 
127  end subroutine aero_info_array_enlarge_to
128 
129 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
130 
131  !> Possibly shrinks the storage of the given array, ensuring that
132  !> it can still store the allocated particles.
133  subroutine aero_info_array_shrink(aero_info_array)
134 
135  !> Array to shrink.
136  type(aero_info_array_t), intent(inout) :: aero_info_array
137 
138  integer :: new_length
139 
140  if (.not. allocated(aero_info_array%aero_info)) return
141 
142  new_length = pow2_above(aero_info_array%n_item)
143  if (new_length < size(aero_info_array%aero_info)) then
144  call aero_info_array_realloc(aero_info_array, new_length)
145  end if
146 
147  end subroutine aero_info_array_shrink
148 
149 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
150 
151  !> Adds the given aero_info to the end of the array.
152  subroutine aero_info_array_add_aero_info(aero_info_array, &
153  aero_info)
154 
155  !> Array to add to.
156  type(aero_info_array_t), intent(inout) :: aero_info_array
157  !> Aero_info to add.
158  type(aero_info_t), intent(in) :: aero_info
159 
160  integer :: n
161 
162  n = aero_info_array_n_item(aero_info_array) + 1
163  call aero_info_array_enlarge_to(aero_info_array, n)
164  aero_info_array%aero_info(n) = aero_info
165  aero_info_array%n_item = n
166 
167  end subroutine aero_info_array_add_aero_info
168 
169 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
170 
171  !> Removes the aero_info at the given index.
172  subroutine aero_info_array_remove_aero_info(aero_info_array, &
173  index)
174 
175  !> Array to remove from.
176  type(aero_info_array_t), intent(inout) :: aero_info_array
177  !> Index of aero_info to remove.
178  integer, intent(in) :: index
179 
180  call assert(578870706, allocated(aero_info_array%aero_info))
181  call assert(213892348, index >= 1)
182  call assert(953927392, index <= aero_info_array%n_item)
183  if (index < aero_info_array%n_item) then
184  ! shift last aero_info into empty slot to preserve dense packing
185  aero_info_array%aero_info(index) &
186  = aero_info_array%aero_info(aero_info_array%n_item)
187  end if
188  aero_info_array%n_item = aero_info_array%n_item - 1
189  call aero_info_array_shrink(aero_info_array)
190 
191  end subroutine aero_info_array_remove_aero_info
192 
193 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
194 
195  !> Adds \c aero_info_array_delta to the end of \c aero_info_array.
196  subroutine aero_info_array_add(aero_info_array, &
197  aero_info_array_delta)
198 
199  !> Array to add to.
200  type(aero_info_array_t), intent(inout) :: aero_info_array
201  !> Aero_info to add.
202  type(aero_info_array_t), intent(in) :: aero_info_array_delta
203 
204  integer :: i, n, n_delta, n_new
205 
206  n = aero_info_array%n_item
207  n_delta = aero_info_array_delta%n_item
208  n_new = n + n_delta
209  call aero_info_array_enlarge_to(aero_info_array, n_new)
210  do i = 1,n_delta
211  aero_info_array%aero_info(n + i) = aero_info_array_delta%aero_info(i)
212  end do
213  aero_info_array%n_item = n_new
214 
215  end subroutine aero_info_array_add
216 
217 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
218 
219  !> Determines the number of bytes required to pack the given value.
220  integer function pmc_mpi_pack_size_aia(val)
221 
222  !> Value to pack.
223  type(aero_info_array_t), intent(in) :: val
224 
225  integer :: i, total_size
226 
227  total_size = 0
228  total_size = total_size &
230  do i = 1,aero_info_array_n_item(val)
231  total_size = total_size &
232  + pmc_mpi_pack_size_aero_info(val%aero_info(i))
233  end do
234  pmc_mpi_pack_size_aia = total_size
235 
236  end function pmc_mpi_pack_size_aia
237 
238 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
239 
240  !> Packs the given value into the buffer, advancing position.
241  subroutine pmc_mpi_pack_aero_info_array(buffer, position, val)
242 
243  !> Memory buffer.
244  character, intent(inout) :: buffer(:)
245  !> Current buffer position.
246  integer, intent(inout) :: position
247  !> Value to pack.
248  type(aero_info_array_t), intent(in) :: val
249 
250 #ifdef PMC_USE_MPI
251  integer :: prev_position, i
252 
253  prev_position = position
254  call pmc_mpi_pack_integer(buffer, position, aero_info_array_n_item(val))
255  do i = 1,aero_info_array_n_item(val)
256  call pmc_mpi_pack_aero_info(buffer, position, val%aero_info(i))
257  end do
258  call assert(732927292, &
259  position - prev_position <= pmc_mpi_pack_size_aia(val))
260 #endif
261 
262  end subroutine pmc_mpi_pack_aero_info_array
263 
264 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
265 
266  !> Unpacks the given value from the buffer, advancing position.
267  subroutine pmc_mpi_unpack_aero_info_array(buffer, position, val)
268 
269  !> Memory buffer.
270  character, intent(inout) :: buffer(:)
271  !> Current buffer position.
272  integer, intent(inout) :: position
273  !> Value to pack.
274  type(aero_info_array_t), intent(inout) :: val
275 
276 #ifdef PMC_USE_MPI
277  integer :: prev_position, i, n
278 
279  prev_position = position
280  call pmc_mpi_unpack_integer(buffer, position, n)
281  call aero_info_array_realloc(val, n)
282  val%n_item = n
283  do i = 1,n
284  call pmc_mpi_unpack_aero_info(buffer, position, val%aero_info(i))
285  end do
286  call assert(262838429, &
287  position - prev_position <= pmc_mpi_pack_size_aia(val))
288 #endif
289 
290  end subroutine pmc_mpi_unpack_aero_info_array
291 
292 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
293 
294 end module pmc_aero_info_array
pmc_aero_info_array::pmc_mpi_pack_aero_info_array
subroutine pmc_mpi_pack_aero_info_array(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: aero_info_array.F90:242
pmc_mpi
Wrapper functions for MPI.
Definition: mpi.F90:13
pmc_aero_info_array::aero_info_array_enlarge_to
subroutine aero_info_array_enlarge_to(aero_info_array, n)
Possibly enlarges the given array, ensuring that it is at least of size n.
Definition: aero_info_array.F90:112
pmc_aero_info_array::aero_info_array_add_aero_info
subroutine aero_info_array_add_aero_info(aero_info_array, aero_info)
Adds the given aero_info to the end of the array.
Definition: aero_info_array.F90:154
pmc_aero_info_array::aero_info_array_realloc
subroutine aero_info_array_realloc(aero_info_array, new_length)
Changes the given aero_info_array to exactly the given new_length.
Definition: aero_info_array.F90:84
pmc_aero_info::pmc_mpi_pack_size_aero_info
integer function pmc_mpi_pack_size_aero_info(val)
Determines the number of bytes required to pack the given value.
Definition: aero_info.F90:64
pmc_aero_info::pmc_mpi_unpack_aero_info
subroutine pmc_mpi_unpack_aero_info(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: aero_info.F90:107
pmc_aero_info_array::aero_info_array_shrink
subroutine aero_info_array_shrink(aero_info_array)
Possibly shrinks the storage of the given array, ensuring that it can still store the allocated parti...
Definition: aero_info_array.F90:134
pmc_spec_file
Reading formatted text input.
Definition: spec_file.F90:43
pmc_util::pow2_above
integer function pow2_above(n)
Return the least power-of-2 that is at least equal to n.
Definition: util.F90:1642
pmc_util::assert
subroutine assert(code, condition_ok)
Errors unless condition_ok is true.
Definition: util.F90:103
pmc_aero_info::pmc_mpi_pack_aero_info
subroutine pmc_mpi_pack_aero_info(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: aero_info.F90:82
pmc_aero_info_array::pmc_mpi_pack_size_aia
integer function pmc_mpi_pack_size_aia(val)
Determines the number of bytes required to pack the given value.
Definition: aero_info_array.F90:221
pmc_aero_info_array::aero_info_array_add
subroutine aero_info_array_add(aero_info_array, aero_info_array_delta)
Adds aero_info_array_delta to the end of aero_info_array.
Definition: aero_info_array.F90:198
pmc_aero_info_array::aero_info_array_zero
subroutine aero_info_array_zero(aero_info_array)
Sets an aero_info_array to contain zero data.
Definition: aero_info_array.F90:63
pmc_aero_info_array
The aero_info_array_t structure and assoicated subroutines.
Definition: aero_info_array.F90:9
pmc_aero_info_array::aero_info_array_n_item
elemental integer function aero_info_array_n_item(aero_info_array)
Return the current number of items.
Definition: aero_info_array.F90:47
pmc_mpi::pmc_mpi_unpack_integer
subroutine pmc_mpi_unpack_integer(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:818
pmc_aero_info::aero_info_t
Information about removed particles describing the sink.
Definition: aero_info.F90:48
pmc_util
Common utility subroutines.
Definition: util.F90:9
pmc_aero_info
The aero_info_t structure and associated subroutines.
Definition: aero_info.F90:9
pmc_aero_info_array::aero_info_array_t
1-D arrays of aero_info_t structure.
Definition: aero_info_array.F90:33
pmc_aero_info_array::aero_info_array_remove_aero_info
subroutine aero_info_array_remove_aero_info(aero_info_array, index)
Removes the aero_info at the given index.
Definition: aero_info_array.F90:174
pmc_mpi::pmc_mpi_pack_size_integer
integer function pmc_mpi_pack_size_integer(val)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:345
pmc_mpi::pmc_mpi_pack_integer
subroutine pmc_mpi_pack_integer(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:561
pmc_aero_info_array::pmc_mpi_unpack_aero_info_array
subroutine pmc_mpi_unpack_aero_info_array(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: aero_info_array.F90:268