PartMC  2.8.0
coagulation_dist.F90
Go to the documentation of this file.
1 ! Copyright (C) 2005-2012 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_coagulation_dist module.
7 
8 !> Parallel aerosol particle coagulation with MPI
10 
11  use pmc_bin_grid
12  use pmc_aero_data
13  use pmc_util
14  use pmc_env_state
15  use pmc_aero_state
17  use pmc_coagulation
18  use pmc_mpi
19 #ifdef PMC_USE_MPI
20  use mpi
21 #endif
22 
23  !> Size of the outgoing buffer for \c bsend (bytes).
24  !!
25  !! FIXME: check that this size is big enough. It must be large
26  !! enough to handle the required number of messages of the given
27  !! sizes, plus MPI_BSEND_OVERHEAD per message, plus some extra room
28  !! because it's only kind of a circular buffer --- the messages
29  !! themselves aren't allowed to wrap around then end, so we might
30  !! need extra space up to the size of the largest message type.
31  integer, parameter :: coag_dist_outgoing_buffer_size = 1000000
32  !> Size of send and receive buffer for each message (bytes).
33  !!
34  !! The biggest message type will be one of the particle-sending
35  !! types, for which we need pmc_mpi_pack_size_aero_particle(), plus
36  !! a couple of integers or something. At the moment this means
37  !! something like (10 + n_spec) reals, (3 + 2) integers, which for
38  !! n_spec = 20 gives a size of 260 bytes.
39  integer, parameter :: coag_dist_max_buffer_size = 10000
40  integer, parameter :: coag_dist_max_requests = 1
41  integer, parameter :: coag_dist_tag_request_particle = 5321
42  integer, parameter :: coag_dist_tag_return_req_particle = 5322
43  integer, parameter :: coag_dist_tag_return_unreq_particle = 5323
44  integer, parameter :: coag_dist_tag_return_no_particle = 5324
45  integer, parameter :: coag_dist_tag_done = 5325
46 
47  !> A single outstanding request for a remote particle.
48  type request_t
49  !> Local \c aero_particle to maybe coagulate with the received
50  !> particle.
51  type(aero_particle_t) :: local_aero_particle
52  !> Remote process number that we sent the request to
53  !> (-1 means this request is currently not used).
54  integer :: remote_proc
55  !> Local bin number from which we took \c local_aero_particle.
56  integer :: local_bin
57  !> Remote bin number from which we requested an \c aero_particle.
58  integer :: remote_bin
59  !> Whether this request is currently active
60  logical :: active
61  end type request_t
62 
63 contains
64 
65 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
66 
67  ! Allocate a request object and set its state to invalid.
68  subroutine request_allocate(request)
69 
70  !> Request object to allocate.
71  type(request_t), intent(out) :: request
72 
73  request%active = .false.
74 
75  end subroutine request_allocate
76 
77 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
78 
79  !> Deallocate a request object and set it to be invalid.
80  subroutine request_deallocate(request)
81 
82  !> Request object to deallocate
83  type(request_t), intent(inout) :: request
84 
85  request%active = .false.
86 
87  end subroutine request_deallocate
88 
89 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
90 
91  !> Whether the given reqest object is currectly active.
92  logical function request_is_active(request)
93 
94  !> Request object to test for activeness.
95  type(request_t), intent(in) :: request
96 
97  request_is_active = request%active
98 
99  end function request_is_active
100 
101 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
102 
103  !> Do coagulation for time del_t.
104  subroutine mc_coag_dist(coag_kernel_type, env_state, aero_data, &
105  aero_state, del_t, tot_n_samp, tot_n_coag)
106 
107  !> Coagulation kernel type.
108  integer, intent(in) :: coag_kernel_type
109  !> Environment state.
110  type(env_state_t), intent(in) :: env_state
111  !> Aerosol data.
112  type(aero_data_t), intent(in) :: aero_data
113  !> Aerosol state.
114  type(aero_state_t), intent(inout) :: aero_state
115  !> Timestep.
116  real(kind=dp), intent(in) :: del_t
117  !> Total number of samples tested.
118  integer, intent(out) :: tot_n_samp
119  !> Number of coagulation events.
120  integer, intent(out) :: tot_n_coag
121 
122  integer, parameter :: s1 = 1
123  integer, parameter :: s2 = 1
124  integer, parameter :: sc = 1
125 
126 #ifdef PMC_USE_MPI
127  logical :: samps_remaining, sent_dones
128  integer :: i_bin, j_bin, n_samp, i_samp, i_proc, n_proc
129  integer :: ierr, status(MPI_STATUS_SIZE), current_i, current_j, i_req
130  real(kind=dp) :: n_samp_real, f_max
131  integer, allocatable :: n_parts(:,:)
132  real(kind=dp), allocatable :: magnitudes(:,:)
133  type(request_t) :: requests(coag_dist_max_requests)
134  integer, allocatable :: n_samps(:,:)
135  real(kind=dp), allocatable :: accept_factors(:,:), k_max(:,:)
136  logical, allocatable :: procs_done(:)
137  integer, allocatable :: outgoing_buffer(:)
138  integer :: outgoing_buffer_size_check
139  type(aero_weight_array_t) :: aero_weight_total
140 
141  call assert_msg(667898741, &
142  aero_sorted_n_class(aero_state%aero_sorted) == 1, &
143  "FIXME: mc_coag_dist() can only handle one weight class")
144 
145  allocate(outgoing_buffer(coag_dist_outgoing_buffer_size))
146 
147  n_proc = pmc_mpi_size()
148 
149  call pmc_mpi_barrier()
150 
151  call aero_state_sort(aero_state, aero_data, all_procs_same=.true.)
152  if (.not. aero_state%aero_sorted%coag_kernel_bounds_valid) then
153  call est_k_minmax_binned_unweighted(aero_state%aero_sorted%bin_grid, &
154  coag_kernel_type, aero_data, env_state, &
155  aero_state%aero_sorted%coag_kernel_min, &
156  aero_state%aero_sorted%coag_kernel_max)
157  aero_state%aero_sorted%coag_kernel_bounds_valid = .true.
158  end if
159 
160  allocate(n_samps(bin_grid_size(aero_state%aero_sorted%bin_grid), &
161  bin_grid_size(aero_state%aero_sorted%bin_grid)))
162  allocate(accept_factors(bin_grid_size(aero_state%aero_sorted%bin_grid), &
163  bin_grid_size(aero_state%aero_sorted%bin_grid)))
164 
165  allocate(n_parts(bin_grid_size(aero_state%aero_sorted%bin_grid), n_proc))
166  call pmc_mpi_allgather_integer_array(integer_varray_n_entry( &
167  aero_state%aero_sorted%size_class%inverse(:, s1)), n_parts)
168 
169  allocate(magnitudes(size(aero_state%awa%weight), n_proc))
170  call pmc_mpi_allgather_real_array(aero_state%awa%weight(:, s1)%magnitude, &
171  magnitudes)
172 
173  aero_weight_total = aero_state%awa
174  aero_weight_total%weight(:, s1)%magnitude = 1d0 / sum(1d0 / magnitudes, 2)
175 
176  allocate(k_max(bin_grid_size(aero_state%aero_sorted%bin_grid), &
177  bin_grid_size(aero_state%aero_sorted%bin_grid)))
178  do i_bin = 1,bin_grid_size(aero_state%aero_sorted%bin_grid)
179  do j_bin = 1,bin_grid_size(aero_state%aero_sorted%bin_grid)
180  call max_coag_num_conc_factor(aero_weight_total, &
181  aero_data, aero_state%aero_sorted%bin_grid, &
182  i_bin, j_bin, s1, s2, sc, f_max)
183  k_max(i_bin, j_bin) &
184  = aero_state%aero_sorted%coag_kernel_max(i_bin, j_bin) * f_max
185  end do
186  end do
187 
188  call generate_n_samps(n_parts, del_t, aero_state%aero_sorted%bin_grid, &
189  aero_weight_total, k_max, n_samps, accept_factors)
190  tot_n_samp = sum(n_samps)
191  tot_n_coag = 0
192 
193  ! main loop
194  do i_req = 1,coag_dist_max_requests
195  call request_allocate(requests(i_req))
196  end do
197  samps_remaining = .true.
198  current_i = 1
199  current_j = 1
200  allocate(procs_done(n_proc))
201  procs_done = .false.
202  sent_dones = .false.
203  call mpi_buffer_attach(outgoing_buffer, &
205  call pmc_mpi_check_ierr(ierr)
206  do while (.not. all(procs_done))
207  ! add requests if we have any slots available call
208  call add_coagulation_requests(aero_state, requests, n_parts, &
209  current_i, current_j, n_samps, samps_remaining)
210 
211  ! if we have no outstanding requests, send done messages
212  if (.not. sent_dones) then
213  if (.not. any_requests_active(requests)) then
214  sent_dones = .true.
215  do i_proc = 0, (n_proc - 1)
216  call send_done(i_proc)
217  end do
218  end if
219  end if
220 
221  ! receive exactly one message
222  call coag_dist_recv(requests, env_state, aero_weight_total, aero_data, &
223  aero_state, accept_factors, coag_kernel_type, tot_n_coag, &
224  magnitudes, procs_done)
225  end do
226 
227  do i_req = 1,coag_dist_max_requests
228  call assert(502009333, .not. request_is_active(requests(i_req)))
229  call request_deallocate(requests(i_req))
230  end do
231  deallocate(procs_done)
232  deallocate(n_samps)
233  deallocate(accept_factors)
234  deallocate(n_parts)
235  deallocate(magnitudes)
236  call mpi_buffer_detach(outgoing_buffer, &
237  outgoing_buffer_size_check, ierr)
238  call pmc_mpi_check_ierr(ierr)
239  call assert(577822730, &
240  coag_dist_outgoing_buffer_size == outgoing_buffer_size_check)
241  deallocate(outgoing_buffer)
242 #endif
243 
244  end subroutine mc_coag_dist
245 
246 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
247 
248  subroutine coag_dist_recv(requests, env_state, aero_weight_total, &
249  aero_data, aero_state, accept_factors, coag_kernel_type, tot_n_coag, &
250  magnitudes, procs_done)
251 
252  !> Array of outstanding requests.
253  type(request_t), intent(inout) :: requests(COAG_DIST_MAX_REQUESTS)
254  !> Environment state.
255  type(env_state_t), intent(in) :: env_state
256  !> Total weighting functions.
257  type(aero_weight_array_t), intent(in) :: aero_weight_total
258  !> Aerosol data.
259  type(aero_data_t), intent(in) :: aero_data
260  !> Aerosol state.
261  type(aero_state_t), intent(inout) :: aero_state
262  !> Accept scale factors per bin pair (1).
263  real(kind=dp), intent(in) :: accept_factors(:,:)
264  !> Coagulation kernel type.
265  integer, intent(in) :: coag_kernel_type
266  !> Number of coagulation events.
267  integer, intent(inout) :: tot_n_coag
268  !> Computational volumes on all processes.
269  real(kind=dp), intent(in) :: magnitudes(:,:)
270  !> Which processes are finished with coagulation.
271  logical, intent(inout) :: procs_done(:)
272 
273 #ifdef PMC_USE_MPI
274  integer :: status(MPI_STATUS_SIZE), ierr
275 
276  call mpi_probe(mpi_any_source, mpi_any_tag, mpi_comm_world, &
277  status, ierr)
278  call pmc_mpi_check_ierr(ierr)
279  if (status(mpi_tag) == coag_dist_tag_request_particle) then
280  call recv_request_particle(aero_state)
281  elseif (status(mpi_tag) == coag_dist_tag_return_req_particle) then
282  call recv_return_req_particle(requests, env_state, aero_weight_total, &
283  aero_data, aero_state, accept_factors, coag_kernel_type, &
284  tot_n_coag, magnitudes)
285  elseif (status(mpi_tag) == coag_dist_tag_return_unreq_particle) then
286  call recv_return_unreq_particle(aero_state, aero_data)
287  elseif (status(mpi_tag) == coag_dist_tag_return_no_particle) then
288  call recv_return_no_particle(requests, aero_data, aero_state)
289  elseif (status(mpi_tag) == coag_dist_tag_done) then
290  call recv_done(procs_done)
291  else
292  call die_msg(856123972, &
293  'unknown tag: ' // trim(integer_to_string(status(mpi_tag))))
294  end if
295 #endif
296 
297  end subroutine coag_dist_recv
298 
299 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
300 
301  subroutine add_coagulation_requests(aero_state, requests, n_parts, &
302  local_bin, remote_bin, n_samps, samps_remaining)
303 
304  !> Aerosol state.
305  type(aero_state_t), intent(inout) :: aero_state
306  !> Array of outstanding requests.
307  type(request_t), intent(inout) :: requests(COAG_DIST_MAX_REQUESTS)
308  !> Number of particles per bin per process.
309  integer, intent(in) :: n_parts(:,:)
310  !> Bin index of first particle we need to coagulate.
311  integer, intent(inout) :: local_bin
312  !> Bin index of second particle we need to coagulate.
313  integer, intent(inout) :: remote_bin
314  !> Number of samples remaining per bin pair
315  integer, intent(inout) :: n_samps(:,:)
316  !> Whether there are still coagulation samples that need to be done.
317  logical, intent(inout) :: samps_remaining
318 
319  integer, parameter :: s1 = 1
320  integer, parameter :: s2 = 1
321  integer, parameter :: sc = 1
322 
323  integer :: i_req
324 
325  if (.not. samps_remaining) return
326 
327  outer: do i_req = 1,coag_dist_max_requests
328  if (.not. request_is_active(requests(i_req))) then
329  inner: do
330  call update_n_samps(n_samps, local_bin, remote_bin, &
331  samps_remaining)
332  if (.not. samps_remaining) exit outer
333  if (integer_varray_n_entry( &
334  aero_state%aero_sorted%size_class%inverse(local_bin, s2)) &
335  > 0) then
336  call find_rand_remote_proc(n_parts, remote_bin, &
337  requests(i_req)%remote_proc)
338  requests(i_req)%active = .true.
339  requests(i_req)%local_bin = local_bin
340  requests(i_req)%remote_bin = remote_bin
342  local_bin, s2, requests(i_req)%local_aero_particle)
343  call send_request_particle(requests(i_req)%remote_proc, &
344  requests(i_req)%remote_bin)
345  exit inner
346  end if
347  end do inner
348  end if
349  end do outer
350 
351  end subroutine add_coagulation_requests
352 
353 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
354 
355  !> Returns \c .true. if any of the requests are active, otherwise
356  !> returns \c .false.
357  logical function any_requests_active(requests)
358 
359  !> Array of outstanding requests.
360  type(request_t), intent(inout) :: requests(coag_dist_max_requests)
361 
362  integer :: i_req
363 
364  do i_req = 1,coag_dist_max_requests
365  if (request_is_active(requests(i_req))) then
366  any_requests_active = .true.
367  return
368  end if
369  end do
370  any_requests_active = .false.
371 
372  end function any_requests_active
373 
374 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
375 
376  subroutine find_rand_remote_proc(n_parts, remote_bin, remote_proc)
377 
378  !> Number of particles per bin per process.
379  integer, intent(in) :: n_parts(:,:)
380  !> Remote bin number.
381  integer, intent(in) :: remote_bin
382  !> Remote process number chosen at random.
383  integer, intent(out) :: remote_proc
384 
385 #ifdef PMC_USE_MPI
386  call assert(770964285, size(n_parts, 2) == pmc_mpi_size())
387  remote_proc = sample_disc_pdf(n_parts(remote_bin, :)) - 1
388 #endif
389 
390  end subroutine find_rand_remote_proc
391 
392 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
393 
394  subroutine update_n_samps(n_samps, local_bin, remote_bin, samps_remaining)
395 
396  !> Number of samples remaining per bin pair
397  integer, intent(inout) :: n_samps(:,:)
398  !> Bin index of first particle we need to coagulate.
399  integer, intent(inout) :: local_bin
400  !> Bin index of second particle we need to coagulate.
401  integer, intent(inout) :: remote_bin
402  !> Whether there are still coagulation samples that need to be done.
403  logical, intent(inout) :: samps_remaining
404 
405  integer :: n_bin
406 
407  if (.not. samps_remaining) return
408 
409  n_bin = size(n_samps, 1)
410  do
411  if (n_samps(local_bin, remote_bin) > 0) exit
412 
413  remote_bin = remote_bin + 1
414  if (remote_bin > n_bin) then
415  remote_bin = 1
416  local_bin = local_bin + 1
417  end if
418  if (local_bin > n_bin) exit
419  end do
420 
421  if (local_bin > n_bin) then
422  samps_remaining = .false.
423  else
424  n_samps(local_bin, remote_bin) = n_samps(local_bin, remote_bin) - 1
425  end if
426 
427  end subroutine update_n_samps
428 
429 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
430 
431  subroutine send_request_particle(remote_proc, remote_bin)
432 
433  !> Remote process number.
434  integer, intent(in) :: remote_proc
435  !> Remote bin number.
436  integer, intent(in) :: remote_bin
437 
438 #ifdef PMC_USE_MPI
439  character :: buffer(COAG_DIST_MAX_BUFFER_SIZE)
440  integer :: buffer_size, max_buffer_size, position, ierr
441 
442  max_buffer_size = 0
443  max_buffer_size = max_buffer_size + pmc_mpi_pack_size_integer(remote_bin)
444  call assert(893545122, max_buffer_size <= coag_dist_max_buffer_size)
445  position = 0
446  call pmc_mpi_pack_integer(buffer, position, remote_bin)
447  call assert(610314213, position <= max_buffer_size)
448  buffer_size = position
449  call mpi_bsend(buffer, buffer_size, mpi_character, remote_proc, &
450  coag_dist_tag_request_particle, mpi_comm_world, ierr)
451  call pmc_mpi_check_ierr(ierr)
452 #endif
453 
454  end subroutine send_request_particle
455 
456 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
457 
458  subroutine recv_request_particle(aero_state)
459 
460  !> Aero state.
461  type(aero_state_t), intent(inout) :: aero_state
462 
463  integer, parameter :: s1 = 1
464  integer, parameter :: s2 = 1
465  integer, parameter :: sc = 1
466 
467 #ifdef PMC_USE_MPI
468  integer :: buffer_size, position, request_bin, sent_proc
469  integer :: ierr, remote_proc, status(MPI_STATUS_SIZE)
470  character :: buffer(COAG_DIST_MAX_BUFFER_SIZE)
471  type(aero_particle_t) :: aero_particle
472 
473  ! get the message
474  call mpi_recv(buffer, coag_dist_max_buffer_size, mpi_character, &
475  mpi_any_source, coag_dist_tag_request_particle, mpi_comm_world, &
476  status, ierr)
477  call pmc_mpi_check_ierr(ierr)
478  call assert(920139874, status(mpi_tag) &
480  call mpi_get_count(status, mpi_character, buffer_size, ierr)
481  call pmc_mpi_check_ierr(ierr)
482  call assert(190658659, buffer_size <= coag_dist_max_buffer_size)
483  remote_proc = status(mpi_source)
484 
485  ! unpack it
486  position = 0
487  call pmc_mpi_unpack_integer(buffer, position, request_bin)
488  call assert(895128380, position == buffer_size)
489 
490  ! send the particle back if we have one
491  if (integer_varray_n_entry( &
492  aero_state%aero_sorted%size_class%inverse(request_bin, s1)) == 0) then
493  call send_return_no_particle(remote_proc, request_bin)
494  else
496  request_bin, s1, aero_particle)
497  call send_return_req_particle(aero_particle, request_bin, &
498  remote_proc)
499  end if
500 #endif
501 
502  end subroutine recv_request_particle
503 
504 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
505 
506  subroutine send_return_no_particle(dest_proc, i_bin)
507 
508  !> Process number to send message to.
509  integer, intent(in) :: dest_proc
510  !> Bin number where there was no particle.
511  integer, intent(in) :: i_bin
512 
513 #ifdef PMC_USE_MPI
514  character :: buffer(COAG_DIST_MAX_BUFFER_SIZE)
515  integer :: buffer_size, max_buffer_size, position, ierr
516 
517  max_buffer_size = 0
518  max_buffer_size = max_buffer_size + pmc_mpi_pack_size_integer(i_bin)
519  call assert(744787119, max_buffer_size <= coag_dist_max_buffer_size)
520  position = 0
521  call pmc_mpi_pack_integer(buffer, position, i_bin)
522  call assert(445960340, position <= max_buffer_size)
523  buffer_size = position
524  call mpi_bsend(buffer, buffer_size, mpi_character, dest_proc, &
525  coag_dist_tag_return_no_particle, mpi_comm_world, ierr)
526  call pmc_mpi_check_ierr(ierr)
527 #endif
528 
529  end subroutine send_return_no_particle
530 
531 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
532 
533  subroutine recv_return_no_particle(requests, aero_data, aero_state)
534 
535  !> Array of outstanding requests.
536  type(request_t), intent(inout) :: requests(COAG_DIST_MAX_REQUESTS)
537  !> Aerosol data.
538  type(aero_data_t), intent(in) :: aero_data
539  !> Aerosol state.
540  type(aero_state_t), intent(inout) :: aero_state
541 
542 #ifdef PMC_USE_MPI
543  logical :: found_request
544  integer :: buffer_size, position, sent_bin, sent_proc, i_req
545  integer :: ierr, status(MPI_STATUS_SIZE)
546  character :: buffer(COAG_DIST_MAX_BUFFER_SIZE)
547 
548  ! get the message
549  call mpi_recv(buffer, coag_dist_max_buffer_size, mpi_character, &
550  mpi_any_source, coag_dist_tag_return_no_particle, &
551  mpi_comm_world, status, ierr)
552  call pmc_mpi_check_ierr(ierr)
553  call assert(918153221, status(mpi_tag) &
555  call mpi_get_count(status, mpi_character, buffer_size, ierr)
556  call pmc_mpi_check_ierr(ierr)
557  call assert(461111487, buffer_size <= coag_dist_max_buffer_size)
558  sent_proc = status(mpi_source)
559 
560  ! unpack it
561  position = 0
562  call pmc_mpi_unpack_integer(buffer, position, sent_bin)
563  call assert(518172999, position == buffer_size)
564 
565  ! find the matching request
566  found_request = .false.
567  do i_req = 1,coag_dist_max_requests
568  if ((requests(i_req)%remote_proc == sent_proc) &
569  .and. (requests(i_req)%remote_bin == sent_bin)) then
570  found_request = .true.
571  exit
572  end if
573  end do
574  call assert(215612776, found_request)
575 
576  ! We can't do coagulation with the local particle, so store it
577  ! back. If we wanted to, we could use the knowledge that it should
578  ! go into bin requests(i_req)%local_bin
579  call aero_state_add_particle(aero_state, &
580  requests(i_req)%local_aero_particle, aero_data, &
581  allow_resort=.false.)
582  call request_deallocate(requests(i_req))
583  call request_allocate(requests(i_req))
584 #endif
585 
586  end subroutine recv_return_no_particle
587 
588 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
589 
590  subroutine send_return_req_particle(aero_particle, i_bin, dest_proc)
591 
592  !> Aero particle to send.
593  type(aero_particle_t), intent(in) :: aero_particle
594  !> Bin that the particle is in.
595  integer, intent(in) :: i_bin
596  !> Process number to send particle to.
597  integer, intent(in) :: dest_proc
598 
599 #ifdef PMC_USE_MPI
600  character :: buffer(COAG_DIST_MAX_BUFFER_SIZE)
601  integer :: buffer_size, max_buffer_size, position, ierr
602 
603  max_buffer_size = 0
604  max_buffer_size = max_buffer_size + pmc_mpi_pack_size_integer(i_bin)
605  max_buffer_size = max_buffer_size &
606  + pmc_mpi_pack_size_aero_particle(aero_particle)
607  call assert(496283814, max_buffer_size <= coag_dist_max_buffer_size)
608  position = 0
609  call pmc_mpi_pack_integer(buffer, position, i_bin)
610  call pmc_mpi_pack_aero_particle(buffer, position, aero_particle)
611  call assert(263666386, position <= max_buffer_size)
612  buffer_size = position
613  call mpi_bsend(buffer, buffer_size, mpi_character, dest_proc, &
614  coag_dist_tag_return_req_particle, mpi_comm_world, ierr)
615  call pmc_mpi_check_ierr(ierr)
616 #endif
617 
618  end subroutine send_return_req_particle
619 
620 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
621 
622  subroutine recv_return_req_particle(requests, env_state, aero_weight_total, &
623  aero_data, aero_state, accept_factors, coag_kernel_type, tot_n_coag, &
624  magnitudes)
625 
626  !> Array of outstanding requests.
627  type(request_t), intent(inout) :: requests(COAG_DIST_MAX_REQUESTS)
628  !> Environment state.
629  type(env_state_t), intent(in) :: env_state
630  !> Total weighting array.
631  type(aero_weight_array_t), intent(in) :: aero_weight_total
632  !> Aerosol data.
633  type(aero_data_t), intent(in) :: aero_data
634  !> Aerosol state.
635  type(aero_state_t), intent(inout) :: aero_state
636  !> Accept scale factors per bin pair (1).
637  real(kind=dp), intent(in) :: accept_factors(:,:)
638  !> Coagulation kernel type.
639  integer, intent(in) :: coag_kernel_type
640  !> Number of coagulation events.
641  integer, intent(inout) :: tot_n_coag
642  !> Computational volumes on all processes.
643  real(kind=dp), intent(in) :: magnitudes(:,:)
644 
645  integer, parameter :: s1 = 1
646  integer, parameter :: s2 = 1
647  integer, parameter :: sc = 1
648 
649 #ifdef PMC_USE_MPI
650  logical :: found_request, remove_1, remove_2
651  integer :: buffer_size, position, sent_bin, sent_proc, i_req
652  integer :: ierr, status(MPI_STATUS_SIZE)
653  character :: buffer(COAG_DIST_MAX_BUFFER_SIZE)
654  type(aero_particle_t) :: sent_aero_particle
655  real(kind=dp) :: k, p
656 
657  ! get the message
658  call mpi_recv(buffer, coag_dist_max_buffer_size, mpi_character, &
659  mpi_any_source, coag_dist_tag_return_req_particle, &
660  mpi_comm_world, status, ierr)
661  call pmc_mpi_check_ierr(ierr)
662  call assert(133285061, status(mpi_tag) &
664  call mpi_get_count(status, mpi_character, buffer_size, ierr)
665  call pmc_mpi_check_ierr(ierr)
666  call assert(563012836, buffer_size <= coag_dist_max_buffer_size)
667  sent_proc = status(mpi_source)
668 
669  ! unpack it
670  position = 0
671  call pmc_mpi_unpack_integer(buffer, position, sent_bin)
672  call pmc_mpi_unpack_aero_particle(buffer, position, sent_aero_particle)
673  call assert(753356021, position == buffer_size)
674 
675  ! find the matching request
676  found_request = .false.
677  do i_req = 1,coag_dist_max_requests
678  if ((requests(i_req)%remote_proc == sent_proc) &
679  .and. (requests(i_req)%remote_bin == sent_bin)) then
680  found_request = .true.
681  exit
682  end if
683  end do
684  call assert(579308475, found_request)
685 
686  ! maybe do coagulation
687  call num_conc_weighted_kernel(coag_kernel_type, &
688  requests(i_req)%local_aero_particle, sent_aero_particle, &
689  s1, s2, sc, aero_data, aero_weight_total, env_state, k)
690  p = k * accept_factors(requests(i_req)%local_bin, sent_bin)
691 
692  if (pmc_random() .lt. p) then
693  ! coagulation happened, do it
694  tot_n_coag = tot_n_coag + 1
695  call coagulate_dist(aero_data, aero_state, &
696  requests(i_req)%local_aero_particle, sent_aero_particle, &
697  sent_proc, aero_weight_total, magnitudes, remove_1, remove_2)
698  else
699  remove_1 = .false.
700  remove_2 = .false.
701  end if
702 
703  ! send the particles back
704  if (.not. remove_1) then
705  ! If we wanted to, we could use the knowledge that this will go
706  ! into bin requests(i_req)%local_bin
707  call aero_state_add_particle(aero_state, &
708  requests(i_req)%local_aero_particle, aero_data, &
709  allow_resort=.false.)
710  end if
711  if (.not. remove_2) then
712  call send_return_unreq_particle(sent_aero_particle, sent_proc)
713  end if
714 
715  call request_deallocate(requests(i_req))
716  call request_allocate(requests(i_req))
717 #endif
718 
719  end subroutine recv_return_req_particle
720 
721 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
722 
723  subroutine send_return_unreq_particle(aero_particle, dest_proc)
724 
725  !> Aero particle to send.
726  type(aero_particle_t), intent(in) :: aero_particle
727  !> Process to send the particle to.
728  integer, intent(in) :: dest_proc
729 
730 #ifdef PMC_USE_MPI
731  character :: buffer(coag_dist_max_buffer_size)
732  integer :: buffer_size, max_buffer_size, position, ierr
733 
734  max_buffer_size = 0
735  max_buffer_size = max_buffer_size &
736  + pmc_mpi_pack_size_aero_particle(aero_particle)
737  call assert(414990602, max_buffer_size <= coag_dist_max_buffer_size)
738  position = 0
739  call pmc_mpi_pack_aero_particle(buffer, position, aero_particle)
740  call assert(898537822, position <= max_buffer_size)
741  buffer_size = position
742  call mpi_bsend(buffer, buffer_size, mpi_character, dest_proc, &
743  coag_dist_tag_return_unreq_particle, mpi_comm_world, ierr)
744  call pmc_mpi_check_ierr(ierr)
745 #endif
746 
747  end subroutine send_return_unreq_particle
748 
749 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
750 
751  subroutine recv_return_unreq_particle(aero_state, aero_data)
752 
753  !> Aerosol state.
754  type(aero_state_t), intent(inout) :: aero_state
755  !> Aerosol data.
756  type(aero_data_t), intent(in) :: aero_data
757 
758 #ifdef PMC_USE_MPI
759  logical :: found_request
760  integer :: buffer_size, position, sent_proc, ierr
761  character :: buffer(COAG_DIST_MAX_BUFFER_SIZE)
762  type(aero_particle_t) :: aero_particle
763  integer :: status(MPI_STATUS_SIZE), send_proc
764 
765  ! get the message
766  call mpi_recv(buffer, coag_dist_max_buffer_size, mpi_character, &
767  mpi_any_source, coag_dist_tag_return_unreq_particle, &
768  mpi_comm_world, status, ierr)
769  call pmc_mpi_check_ierr(ierr)
770  call assert(496247788, status(mpi_tag) &
772  call mpi_get_count(status, mpi_character, buffer_size, ierr)
773  call pmc_mpi_check_ierr(ierr)
774  call assert(590644042, buffer_size <= coag_dist_max_buffer_size)
775  sent_proc = status(mpi_source)
776 
777  ! unpack it
778  position = 0
779  call pmc_mpi_unpack_aero_particle(buffer, position, aero_particle)
780  call assert(833588594, position == buffer_size)
781 
782  ! put it back
783  call aero_state_add_particle(aero_state, aero_particle, aero_data, &
784  allow_resort=.false.)
785 #endif
786 
787  end subroutine recv_return_unreq_particle
788 
789 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
790 
791  !> Send a message saying that this process is finished with its
792  !> coagulation.
793  subroutine send_done(dest_proc)
794 
795  !> Process to send the message to.
796  integer, intent(in) :: dest_proc
797 
798 #ifdef PMC_USE_MPI
799  character :: buffer(COAG_DIST_MAX_BUFFER_SIZE)
800  integer :: buffer_size, ierr
801 
802  buffer_size = 0
803  call mpi_bsend(buffer, buffer_size, mpi_character, dest_proc, &
804  coag_dist_tag_done, mpi_comm_world, ierr)
805  call pmc_mpi_check_ierr(ierr)
806 #endif
807 
808  end subroutine send_done
809 
810 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
811 
812  !> Receive a done message.
813  subroutine recv_done(procs_done)
814 
815  !> Which processes are finished with coagulation.
816  logical, intent(inout) :: procs_done(:)
817 
818 #ifdef PMC_USE_MPI
819  integer :: buffer_size, sent_proc, ierr
820  character :: buffer(COAG_DIST_MAX_BUFFER_SIZE)
821  integer :: status(MPI_STATUS_SIZE)
822 
823  ! get the message
824  call mpi_recv(buffer, coag_dist_max_buffer_size, mpi_character, &
825  mpi_any_source, coag_dist_tag_done, mpi_comm_world, &
826  status, ierr)
827  call pmc_mpi_check_ierr(ierr)
828  call assert(348737947, status(mpi_tag) &
830  call mpi_get_count(status, mpi_character, buffer_size, ierr)
831  call pmc_mpi_check_ierr(ierr)
832  call assert(214904056, buffer_size == 0)
833  sent_proc = status(mpi_source)
834 
835  ! process it
836  procs_done(sent_proc + 1) = .true.
837 #endif
838 
839  end subroutine recv_done
840 
841 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
842 
843  !> generate the number of samples to do per bin pair.
844  subroutine generate_n_samps(n_parts, del_t, bin_grid, aero_weight_array, &
845  k_max, n_samps, accept_factors)
846 
847  !> Number of particles per bin on all processes.
848  integer, intent(in) :: n_parts(:,:)
849  !> Timestep.
850  real(kind=dp), intent(in) :: del_t
851  !> Bin grid.
852  type(bin_grid_t), intent(in) :: bin_grid
853  !> Weighting function array.
854  type(aero_weight_array_t), intent(in) :: aero_weight_array
855  !> Maximum kernel.
856  real(kind=dp), intent(in) :: k_max(:,:)
857  !> Number of samples to do per bin pair.
858  integer, intent(out) :: n_samps(:,:)
859  !> Accept scale factors per bin pair (1).
860  real(kind=dp), intent(out) :: accept_factors(:,:)
861 
862  integer :: i_bin, j_bin, rank, n_bin
863  real(kind=dp) :: n_samp_mean
864 
865  n_bin = size(k_max, 1)
866  rank = pmc_mpi_rank()
867  n_samps = 0
868  do i_bin = 1,n_bin
869  if (n_parts(i_bin, rank + 1) == 0) &
870  cycle
871  do j_bin = i_bin,n_bin
872  call compute_n_samp(n_parts(i_bin, rank + 1), &
873  sum(n_parts(j_bin, :)), (i_bin == j_bin), &
874  k_max(i_bin, j_bin), del_t, n_samp_mean, &
875  n_samps(i_bin, j_bin), accept_factors(i_bin, j_bin))
876  end do
877  end do
878 
879  end subroutine generate_n_samps
880 
881 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
882 
883  subroutine coagulate_dist(aero_data, aero_state, aero_particle_1, &
884  aero_particle_2, remote_proc, aero_weight_total, magnitudes, &
885  remove_1, remove_2)
886 
887  !> Aerosol data.
888  type(aero_data_t), intent(in) :: aero_data
889  !> Aerosol state.
890  type(aero_state_t), intent(inout) :: aero_state
891  !> First particle to coagulate.
892  type(aero_particle_t), intent(in) :: aero_particle_1
893  !> Second particle to coagulate.
894  type(aero_particle_t), intent(in) :: aero_particle_2
895  !> Remote process that the particle came from.
896  integer, intent(in) :: remote_proc
897  !> Total weight across all processes.
898  type(aero_weight_array_t), intent(in) :: aero_weight_total
899  !> Computational volumes on all processes (m^3).
900  real(kind=dp), intent(in) :: magnitudes(:,:)
901  !> Whether to remove aero_particle_1 after the coagulation.
902  logical, intent(out) :: remove_1
903  !> Whether to remove aero_particle_2 after the coagulation.
904  logical, intent(out) :: remove_2
905 
906  integer, parameter :: s1 = 1
907  integer, parameter :: s2 = 1
908  integer, parameter :: sc = 1
909 
910  type(aero_particle_t) :: aero_particle_new
911  integer :: new_proc, new_group
912  type(aero_info_t) :: aero_info_1, aero_info_2
913  logical :: create_new, id_1_lost, id_2_lost
914 
915  call coagulate_weighting(aero_particle_1, aero_particle_2, &
916  aero_particle_new, s1, s2, sc, aero_data, aero_state%awa, &
917  remove_1, remove_2, create_new, id_1_lost, id_2_lost, &
918  aero_info_1, aero_info_2)
919 
920  if (id_1_lost) then
921  call aero_info_array_add_aero_info(aero_state%aero_info_array, &
922  aero_info_1)
923  end if
924  if (id_2_lost) then
925  call aero_info_array_add_aero_info(aero_state%aero_info_array, &
926  aero_info_2)
927  end if
928 
929  ! add new particle
930  if (create_new) then
931  new_group = aero_weight_array_rand_group(aero_weight_total, sc, &
932  aero_particle_radius(aero_particle_new, aero_data))
933  aero_particle_new%weight_group = new_group
934  new_proc = sample_cts_pdf(1d0 / magnitudes(new_group, :)) - 1
935  call send_return_unreq_particle(aero_particle_new, new_proc)
936  end if
937 
938  end subroutine coagulate_dist
939 
940 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
941 
942 end module pmc_coagulation_dist
pmc_coagulation_dist::send_return_no_particle
subroutine send_return_no_particle(dest_proc, i_bin)
Definition: coagulation_dist.F90:507
pmc_coagulation::coagulate_weighting
subroutine coagulate_weighting(pt1, pt2, ptc, c1, c2, cc, aero_data, aero_weight_array, remove_1, remove_2, create_new, id_1_lost, id_2_lost, aero_info_1, aero_info_2)
Actually coagulate pt1 and pt2 to form ptc and compute weighting effects, including which particles s...
Definition: coagulation.F90:813
pmc_mpi::pmc_mpi_size
integer function pmc_mpi_size()
Returns the total number of processes.
Definition: mpi.F90:134
pmc_mpi
Wrapper functions for MPI.
Definition: mpi.F90:13
pmc_coagulation_dist::send_return_req_particle
subroutine send_return_req_particle(aero_particle, i_bin, dest_proc)
Definition: coagulation_dist.F90:591
pmc_coagulation_dist::coag_dist_outgoing_buffer_size
integer, parameter coag_dist_outgoing_buffer_size
Size of the outgoing buffer for bsend (bytes).
Definition: coagulation_dist.F90:31
pmc_coagulation_dist::coagulate_dist
subroutine coagulate_dist(aero_data, aero_state, aero_particle_1, aero_particle_2, remote_proc, aero_weight_total, magnitudes, remove_1, remove_2)
Definition: coagulation_dist.F90:886
pmc_coagulation_dist::coag_dist_recv
subroutine coag_dist_recv(requests, env_state, aero_weight_total, aero_data, aero_state, accept_factors, coag_kernel_type, tot_n_coag, magnitudes, procs_done)
Definition: coagulation_dist.F90:251
pmc_coagulation_dist::coag_dist_max_buffer_size
integer, parameter coag_dist_max_buffer_size
Size of send and receive buffer for each message (bytes).
Definition: coagulation_dist.F90:39
pmc_coagulation_dist::recv_done
subroutine recv_done(procs_done)
Receive a done message.
Definition: coagulation_dist.F90:814
pmc_coagulation_dist::generate_n_samps
subroutine generate_n_samps(n_parts, del_t, bin_grid, aero_weight_array, k_max, n_samps, accept_factors)
generate the number of samples to do per bin pair.
Definition: coagulation_dist.F90:846
pmc_aero_state::aero_state_remove_rand_particle_from_bin
subroutine aero_state_remove_rand_particle_from_bin(aero_state, i_bin, i_class, aero_particle)
Remove a randomly chosen particle from the given bin and return it.
Definition: aero_state.F90:454
pmc_coagulation_dist::mc_coag_dist
subroutine mc_coag_dist(coag_kernel_type, env_state, aero_data, aero_state, del_t, tot_n_samp, tot_n_coag)
Do coagulation for time del_t.
Definition: coagulation_dist.F90:106
pmc_mpi::pmc_mpi_allgather_real_array
subroutine pmc_mpi_allgather_real_array(send, recv)
Does an allgather of real arrays (must be the same size on all processes).
Definition: mpi.F90:2010
pmc_util::die_msg
subroutine die_msg(code, error_msg)
Error immediately.
Definition: util.F90:135
pmc_mpi::pmc_mpi_rank
integer function pmc_mpi_rank()
Returns the rank of the current process.
Definition: mpi.F90:117
pmc_constants::dp
integer, parameter dp
Kind of a double precision real number.
Definition: constants.F90:12
pmc_aero_state::aero_state_add_particle
subroutine aero_state_add_particle(aero_state, aero_particle, aero_data, allow_resort)
Add the given particle.
Definition: aero_state.F90:367
pmc_aero_weight_array
The aero_weight_array_t structure and associated subroutines.
Definition: aero_weight_array.F90:10
pmc_coagulation_dist::send_return_unreq_particle
subroutine send_return_unreq_particle(aero_particle, dest_proc)
Definition: coagulation_dist.F90:724
pmc_bin_grid::bin_grid_size
elemental integer function bin_grid_size(bin_grid)
Return the number of bins in the grid, or -1 if the bin grid is not allocated.
Definition: bin_grid.F90:51
pmc_coagulation_dist::any_requests_active
logical function any_requests_active(requests)
Returns .true. if any of the requests are active, otherwise returns .false.
Definition: coagulation_dist.F90:358
pmc_env_state::env_state_t
Current environment state.
Definition: env_state.F90:29
pmc_util::assert
subroutine assert(code, condition_ok)
Errors unless condition_ok is true.
Definition: util.F90:104
pmc_coagulation_dist::update_n_samps
subroutine update_n_samps(n_samps, local_bin, remote_bin, samps_remaining)
Definition: coagulation_dist.F90:395
pmc_aero_state
The aero_state_t structure and assocated subroutines.
Definition: aero_state.F90:9
pmc_aero_weight_array::aero_weight_array_t
An array of aerosol size distribution weighting functions.
Definition: aero_weight_array.F90:35
pmc_coagulation_dist::find_rand_remote_proc
subroutine find_rand_remote_proc(n_parts, remote_bin, remote_proc)
Definition: coagulation_dist.F90:377
pmc_coagulation_dist::request_allocate
subroutine request_allocate(request)
Definition: coagulation_dist.F90:69
pmc_coagulation_dist::coag_dist_tag_request_particle
integer, parameter coag_dist_tag_request_particle
Definition: coagulation_dist.F90:41
pmc_coagulation_dist::recv_request_particle
subroutine recv_request_particle(aero_state)
Definition: coagulation_dist.F90:459
pmc_coagulation_dist::send_request_particle
subroutine send_request_particle(remote_proc, remote_bin)
Definition: coagulation_dist.F90:432
pmc_util::integer_to_string
character(len=pmc_util_convert_string_len) function integer_to_string(val)
Convert an integer to a string format.
Definition: util.F90:767
pmc_coagulation
Aerosol particle coagulation.
Definition: coagulation.F90:11
pmc_util::assert_msg
subroutine assert_msg(code, condition_ok, error_msg)
Errors unless condition_ok is true.
Definition: util.F90:78
pmc_coagulation_dist::coag_dist_tag_return_req_particle
integer, parameter coag_dist_tag_return_req_particle
Definition: coagulation_dist.F90:42
pmc_coagulation_dist::coag_dist_tag_done
integer, parameter coag_dist_tag_done
Definition: coagulation_dist.F90:45
pmc_coagulation_dist::coag_dist_tag_return_unreq_particle
integer, parameter coag_dist_tag_return_unreq_particle
Definition: coagulation_dist.F90:43
pmc_rand::sample_cts_pdf
integer function sample_cts_pdf(pdf)
Sample the given continuous probability density function.
Definition: rand.F90:551
pmc_coagulation_dist::send_done
subroutine send_done(dest_proc)
Send a message saying that this process is finished with its coagulation.
Definition: coagulation_dist.F90:794
pmc_coagulation::compute_n_samp
subroutine compute_n_samp(ni, nj, same_bin, k_max, del_t, n_samp_mean, n_samp, accept_factor)
Compute the number of samples required for the pair of bins.
Definition: coagulation.F90:650
pmc_coagulation_dist::recv_return_req_particle
subroutine recv_return_req_particle(requests, env_state, aero_weight_total, aero_data, aero_state, accept_factors, coag_kernel_type, tot_n_coag, magnitudes)
Definition: coagulation_dist.F90:625
pmc_aero_weight_array::aero_weight_array_rand_group
integer function aero_weight_array_rand_group(aero_weight_array, i_class, radius)
Choose a random group at the given radius, with probability inversely proportional to group weight at...
Definition: aero_weight_array.F90:376
pmc_env_state
The env_state_t structure and associated subroutines.
Definition: env_state.F90:9
pmc_coagulation_dist::request_deallocate
subroutine request_deallocate(request)
Deallocate a request object and set it to be invalid.
Definition: coagulation_dist.F90:81
pmc_coagulation_dist::coag_dist_tag_return_no_particle
integer, parameter coag_dist_tag_return_no_particle
Definition: coagulation_dist.F90:44
pmc_aero_data::aero_data_t
Aerosol material properties and associated data.
Definition: aero_data.F90:55
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:1139
pmc_coagulation_dist::add_coagulation_requests
subroutine add_coagulation_requests(aero_state, requests, n_parts, local_bin, remote_bin, n_samps, samps_remaining)
Definition: coagulation_dist.F90:303
pmc_mpi::pmc_mpi_check_ierr
subroutine pmc_mpi_check_ierr(ierr)
Dies if ierr is not ok.
Definition: mpi.F90:40
pmc_util
Common utility subroutines.
Definition: util.F90:9
pmc_mpi::pmc_mpi_allgather_integer_array
subroutine pmc_mpi_allgather_integer_array(send, recv)
Does an allgather of integer arrays (must be the same size on all processes).
Definition: mpi.F90:1975
pmc_mpi::pmc_mpi_barrier
subroutine pmc_mpi_barrier()
Synchronize all processes.
Definition: mpi.F90:103
pmc_coagulation_dist::request_is_active
logical function request_is_active(request)
Whether the given reqest object is currectly active.
Definition: coagulation_dist.F90:93
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_bin_grid
The bin_grid_t structure and associated subroutines.
Definition: bin_grid.F90:9
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:711
pmc_aero_data
The aero_data_t structure and associated subroutines.
Definition: aero_data.F90:9
pmc_bin_grid::bin_grid_t
1D grid, either logarithmic or linear.
Definition: bin_grid.F90:33
pmc_coagulation_dist::recv_return_no_particle
subroutine recv_return_no_particle(requests, aero_data, aero_state)
Definition: coagulation_dist.F90:534
pmc_aero_state::aero_state_t
The current collection of aerosol particles.
Definition: aero_state.F90:69
pmc_coagulation_dist::recv_return_unreq_particle
subroutine recv_return_unreq_particle(aero_state, aero_data)
Definition: coagulation_dist.F90:752
pmc_coagulation_dist
Parallel aerosol particle coagulation with MPI.
Definition: coagulation_dist.F90:9
pmc_coagulation_dist::coag_dist_max_requests
integer, parameter coag_dist_max_requests
Definition: coagulation_dist.F90:40
pmc_rand::sample_disc_pdf
integer function sample_disc_pdf(pdf)
Sample the given discrete probability density function.
Definition: rand.F90:585
pmc_aero_state::aero_state_sort
subroutine aero_state_sort(aero_state, aero_data, bin_grid, all_procs_same)
Sorts the particles if necessary.
Definition: aero_state.F90:3219
pmc_rand::pmc_random
real(kind=dp) function pmc_random()
Returns a random number between 0 and 1.
Definition: rand.F90:139
pmc_coagulation_dist::request_t
A single outstanding request for a remote particle.
Definition: coagulation_dist.F90:48