PartMC  2.3.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  call aero_particle_allocate(request%local_aero_particle)
74  request%active = .false.
75 
76  end subroutine request_allocate
77 
78 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
79 
80  !> Deallocate a request object and set it to be invalid.
81  subroutine request_deallocate(request)
82 
83  !> Request object to deallocate
84  type(request_t), intent(inout) :: request
85 
86  call aero_particle_deallocate(request%local_aero_particle)
87  request%active = .false.
88 
89  end subroutine request_deallocate
90 
91 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
92 
93  !> Whether the given reqest object is currectly active.
94  logical function request_is_active(request)
95 
96  !> Request object to test for activeness.
97  type(request_t), intent(in) :: request
98 
99  request_is_active = request%active
100 
101  end function request_is_active
102 
103 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
104 
105  !> Do coagulation for time del_t.
106  subroutine mc_coag_dist(coag_kernel_type, env_state, aero_data, &
107  aero_state, del_t, tot_n_samp, tot_n_coag)
108 
109  !> Coagulation kernel type.
110  integer, intent(in) :: coag_kernel_type
111  !> Environment state.
112  type(env_state_t), intent(in) :: env_state
113  !> Aerosol data.
114  type(aero_data_t), intent(in) :: aero_data
115  !> Aerosol state.
116  type(aero_state_t), intent(inout) :: aero_state
117  !> Timestep.
118  real(kind=dp), intent(in) :: del_t
119  !> Total number of samples tested.
120  integer, intent(out) :: tot_n_samp
121  !> Number of coagulation events.
122  integer, intent(out) :: tot_n_coag
123 
124  integer, parameter :: s1 = 1
125  integer, parameter :: s2 = 1
126  integer, parameter :: sc = 1
127 
128 #ifdef PMC_USE_MPI
129  logical :: samps_remaining, sent_dones
130  integer :: i_bin, j_bin, n_samp, i_samp, i_proc, n_proc
131  integer :: ierr, status(mpi_status_size), current_i, current_j, i_req
132  real(kind=dp) :: n_samp_real, f_max
133  integer, allocatable :: n_parts(:,:)
134  real(kind=dp), allocatable :: magnitudes(:,:)
135  type(request_t) :: requests(coag_dist_max_requests)
136  integer, allocatable :: n_samps(:,:)
137  real(kind=dp), allocatable :: accept_factors(:,:), k_max(:,:)
138  logical, allocatable :: procs_done(:)
139  integer :: outgoing_buffer(coag_dist_outgoing_buffer_size)
140  integer :: outgoing_buffer_size_check
141  type(aero_weight_array_t) :: aero_weight_total
142 
143  call assert_msg(667898741, &
144  aero_sorted_n_class(aero_state%aero_sorted) == 1, &
145  "FIXME: mc_coag_dist() can only handle one weight class")
146 
147  n_proc = pmc_mpi_size()
148 
149  call pmc_mpi_barrier()
150 
151  call aero_state_sort(aero_state, 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(aero_state%aero_sorted%bin_grid%n_bin, &
161  aero_state%aero_sorted%bin_grid%n_bin))
162  allocate(accept_factors(aero_state%aero_sorted%bin_grid%n_bin, &
163  aero_state%aero_sorted%bin_grid%n_bin))
164 
165  allocate(n_parts(aero_state%aero_sorted%bin_grid%n_bin, n_proc))
167  aero_state%aero_sorted%size_class%inverse(:, s1)%n_entry, 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  call aero_weight_array_allocate(aero_weight_total)
174  call aero_weight_array_copy(aero_state%awa, aero_weight_total)
175  aero_weight_total%weight(:, s1)%magnitude = 1d0 / sum(1d0 / magnitudes, 2)
176 
177  allocate(k_max(aero_state%aero_sorted%bin_grid%n_bin, &
178  aero_state%aero_sorted%bin_grid%n_bin))
179  do i_bin = 1,aero_state%aero_sorted%bin_grid%n_bin
180  do j_bin = 1,aero_state%aero_sorted%bin_grid%n_bin
181  call max_coag_num_conc_factor(aero_weight_total, &
182  aero_state%aero_sorted%bin_grid, i_bin, j_bin, s1, s2, sc, &
183  f_max)
184  k_max(i_bin, j_bin) &
185  = aero_state%aero_sorted%coag_kernel_max(i_bin, j_bin) * f_max
186  end do
187  end do
188 
189  call generate_n_samps(n_parts, del_t, aero_state%aero_sorted%bin_grid, &
190  aero_weight_total, k_max, n_samps, accept_factors)
191  tot_n_samp = sum(n_samps)
192  tot_n_coag = 0
193 
194  ! main loop
195  do i_req = 1,coag_dist_max_requests
196  call request_allocate(requests(i_req))
197  end do
198  samps_remaining = .true.
199  current_i = 1
200  current_j = 1
201  allocate(procs_done(n_proc))
202  procs_done = .false.
203  sent_dones = .false.
204  call mpi_buffer_attach(outgoing_buffer, &
205  coag_dist_outgoing_buffer_size, ierr)
206  call pmc_mpi_check_ierr(ierr)
207  do while (.not. all(procs_done))
208  ! add requests if we have any slots available call
209  call add_coagulation_requests(aero_state, requests, n_parts, &
210  current_i, current_j, n_samps, samps_remaining)
211 
212  ! if we have no outstanding requests, send done messages
213  if (.not. sent_dones) then
214  if (.not. any_requests_active(requests)) then
215  sent_dones = .true.
216  do i_proc = 0, (n_proc - 1)
217  call send_done(i_proc)
218  end do
219  end if
220  end if
221 
222  ! receive exactly one message
223  call coag_dist_recv(requests, env_state, aero_weight_total, aero_data, &
224  aero_state, accept_factors, coag_kernel_type, tot_n_coag, &
225  magnitudes, procs_done)
226  end do
227 
228  do i_req = 1,coag_dist_max_requests
229  call assert(502009333, .not. request_is_active(requests(i_req)))
230  call request_deallocate(requests(i_req))
231  end do
232  deallocate(procs_done)
233  deallocate(n_samps)
234  deallocate(accept_factors)
235  deallocate(n_parts)
236  deallocate(magnitudes)
237  call mpi_buffer_detach(outgoing_buffer, &
238  outgoing_buffer_size_check, ierr)
239  call pmc_mpi_check_ierr(ierr)
240  call assert(577822730, &
241  coag_dist_outgoing_buffer_size == outgoing_buffer_size_check)
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)
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 (aero_state%aero_sorted%size_class%inverse(local_bin, &
334  s2)%n_entry > 0) then
335  call find_rand_remote_proc(n_parts, remote_bin, &
336  requests(i_req)%remote_proc)
337  requests(i_req)%active = .true.
338  requests(i_req)%local_bin = local_bin
339  requests(i_req)%remote_bin = remote_bin
341  local_bin, s2, requests(i_req)%local_aero_particle)
342  call send_request_particle(requests(i_req)%remote_proc, &
343  requests(i_req)%remote_bin)
344  exit inner
345  end if
346  end do inner
347  end if
348  end do outer
349 
350  end subroutine add_coagulation_requests
351 
352 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
353 
354  !> Returns \c .true. if any of the requests are active, otherwise
355  !> returns \c .false.
356  logical function any_requests_active(requests)
357 
358  !> Array of outstanding requests.
359  type(request_t), intent(inout) :: requests(coag_dist_max_requests)
360 
361  integer :: i_req
362 
363  do i_req = 1,coag_dist_max_requests
364  if (request_is_active(requests(i_req))) then
365  any_requests_active = .true.
366  return
367  end if
368  end do
369  any_requests_active = .false.
370 
371  end function any_requests_active
372 
373 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
374 
375  subroutine find_rand_remote_proc(n_parts, remote_bin, remote_proc)
376 
377  !> Number of particles per bin per process.
378  integer, intent(in) :: n_parts(:,:)
379  !> Remote bin number.
380  integer, intent(in) :: remote_bin
381  !> Remote process number chosen at random.
382  integer, intent(out) :: remote_proc
383 
384 #ifdef PMC_USE_MPI
385  call assert(770964285, size(n_parts, 2) == pmc_mpi_size())
386  remote_proc = sample_disc_pdf(n_parts(remote_bin, :)) - 1
387 #endif
388 
389  end subroutine find_rand_remote_proc
390 
391 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
392 
393  subroutine update_n_samps(n_samps, local_bin, remote_bin, samps_remaining)
394 
395  !> Number of samples remaining per bin pair
396  integer, intent(inout) :: n_samps(:,:)
397  !> Bin index of first particle we need to coagulate.
398  integer, intent(inout) :: local_bin
399  !> Bin index of second particle we need to coagulate.
400  integer, intent(inout) :: remote_bin
401  !> Whether there are still coagulation samples that need to be done.
402  logical, intent(inout) :: samps_remaining
403 
404  integer :: n_bin
405 
406  if (.not. samps_remaining) return
407 
408  n_bin = size(n_samps, 1)
409  do
410  if (n_samps(local_bin, remote_bin) > 0) exit
411 
412  remote_bin = remote_bin + 1
413  if (remote_bin > n_bin) then
414  remote_bin = 1
415  local_bin = local_bin + 1
416  end if
417  if (local_bin > n_bin) exit
418  end do
419 
420  if (local_bin > n_bin) then
421  samps_remaining = .false.
422  else
423  n_samps(local_bin, remote_bin) = n_samps(local_bin, remote_bin) - 1
424  end if
425 
426  end subroutine update_n_samps
427 
428 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
429 
430  subroutine send_request_particle(remote_proc, remote_bin)
431 
432  !> Remote process number.
433  integer, intent(in) :: remote_proc
434  !> Remote bin number.
435  integer, intent(in) :: remote_bin
436 
437 #ifdef PMC_USE_MPI
438  character :: buffer(coag_dist_max_buffer_size)
439  integer :: buffer_size, max_buffer_size, position, ierr
440 
441  max_buffer_size = 0
442  max_buffer_size = max_buffer_size + pmc_mpi_pack_size_integer(remote_bin)
443  call assert(893545122, max_buffer_size <= coag_dist_max_buffer_size)
444  position = 0
445  call pmc_mpi_pack_integer(buffer, position, remote_bin)
446  call assert(610314213, position <= max_buffer_size)
447  buffer_size = position
448  call mpi_bsend(buffer, buffer_size, mpi_character, remote_proc, &
449  coag_dist_tag_request_particle, mpi_comm_world, ierr)
450  call pmc_mpi_check_ierr(ierr)
451 #endif
452 
453  end subroutine send_request_particle
454 
455 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
456 
457  subroutine recv_request_particle(aero_state)
458 
459  !> Aero state.
460  type(aero_state_t), intent(inout) :: aero_state
461 
462  integer, parameter :: s1 = 1
463  integer, parameter :: s2 = 1
464  integer, parameter :: sc = 1
465 
466 #ifdef PMC_USE_MPI
467  integer :: buffer_size, position, request_bin, sent_proc
468  integer :: ierr, remote_proc, status(mpi_status_size)
469  character :: buffer(coag_dist_max_buffer_size)
470  type(aero_particle_t) :: aero_particle
471 
472  ! get the message
473  call mpi_recv(buffer, coag_dist_max_buffer_size, mpi_character, &
474  mpi_any_source, coag_dist_tag_request_particle, mpi_comm_world, &
475  status, ierr)
476  call pmc_mpi_check_ierr(ierr)
477  call assert(920139874, status(mpi_tag) &
478  == coag_dist_tag_request_particle)
479  call mpi_get_count(status, mpi_character, buffer_size, ierr)
480  call pmc_mpi_check_ierr(ierr)
481  call assert(190658659, buffer_size <= coag_dist_max_buffer_size)
482  remote_proc = status(mpi_source)
483 
484  ! unpack it
485  position = 0
486  call pmc_mpi_unpack_integer(buffer, position, request_bin)
487  call assert(895128380, position == buffer_size)
488 
489  ! send the particle back if we have one
490  if (aero_state%aero_sorted%size_class%inverse(request_bin, &
491  s1)%n_entry == 0) then
492  call send_return_no_particle(remote_proc, request_bin)
493  else
494  call aero_particle_allocate(aero_particle)
496  request_bin, s1, aero_particle)
497  call send_return_req_particle(aero_particle, request_bin, &
498  remote_proc)
499  call aero_particle_deallocate(aero_particle)
500  end if
501 #endif
502 
503  end subroutine recv_request_particle
504 
505 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
506 
507  subroutine send_return_no_particle(dest_proc, i_bin)
508 
509  !> Process number to send message to.
510  integer, intent(in) :: dest_proc
511  !> Bin number where there was no particle.
512  integer, intent(in) :: i_bin
513 
514 #ifdef PMC_USE_MPI
515  character :: buffer(coag_dist_max_buffer_size)
516  integer :: buffer_size, max_buffer_size, position, ierr
517 
518  max_buffer_size = 0
519  max_buffer_size = max_buffer_size + pmc_mpi_pack_size_integer(i_bin)
520  call assert(744787119, max_buffer_size <= coag_dist_max_buffer_size)
521  position = 0
522  call pmc_mpi_pack_integer(buffer, position, i_bin)
523  call assert(445960340, position <= max_buffer_size)
524  buffer_size = position
525  call mpi_bsend(buffer, buffer_size, mpi_character, dest_proc, &
526  coag_dist_tag_return_no_particle, mpi_comm_world, ierr)
527  call pmc_mpi_check_ierr(ierr)
528 #endif
529 
530  end subroutine send_return_no_particle
531 
532 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
533 
534  subroutine recv_return_no_particle(requests, aero_data, aero_state)
535 
536  !> Array of outstanding requests.
537  type(request_t), intent(inout) :: requests(coag_dist_max_requests)
538  !> Aerosol data.
539  type(aero_data_t), intent(in) :: aero_data
540  !> Aerosol state.
541  type(aero_state_t), intent(inout) :: aero_state
542 
543 #ifdef PMC_USE_MPI
544  logical :: found_request
545  integer :: buffer_size, position, sent_bin, sent_proc, i_req
546  integer :: ierr, status(mpi_status_size)
547  character :: buffer(coag_dist_max_buffer_size)
548 
549  ! get the message
550  call mpi_recv(buffer, coag_dist_max_buffer_size, mpi_character, &
551  mpi_any_source, coag_dist_tag_return_no_particle, &
552  mpi_comm_world, status, ierr)
553  call pmc_mpi_check_ierr(ierr)
554  call assert(918153221, status(mpi_tag) &
555  == coag_dist_tag_return_no_particle)
556  call mpi_get_count(status, mpi_character, buffer_size, ierr)
557  call pmc_mpi_check_ierr(ierr)
558  call assert(461111487, buffer_size <= coag_dist_max_buffer_size)
559  sent_proc = status(mpi_source)
560 
561  ! unpack it
562  position = 0
563  call pmc_mpi_unpack_integer(buffer, position, sent_bin)
564  call assert(518172999, position == buffer_size)
565 
566  ! find the matching request
567  found_request = .false.
568  do i_req = 1,coag_dist_max_requests
569  if ((requests(i_req)%remote_proc == sent_proc) &
570  .and. (requests(i_req)%remote_bin == sent_bin)) then
571  found_request = .true.
572  exit
573  end if
574  end do
575  call assert(215612776, found_request)
576 
577  ! We can't do coagulation with the local particle, so store it
578  ! back. If we wanted to, we could use the knowledge that it should
579  ! go into bin requests(i_req)%local_bin
580  call aero_state_add_particle(aero_state, &
581  requests(i_req)%local_aero_particle, 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) &
663  == coag_dist_tag_return_req_particle)
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 aero_particle_allocate(sent_aero_particle)
673  call pmc_mpi_unpack_aero_particle(buffer, position, sent_aero_particle)
674  call assert(753356021, position == buffer_size)
675 
676  ! find the matching request
677  found_request = .false.
678  do i_req = 1,coag_dist_max_requests
679  if ((requests(i_req)%remote_proc == sent_proc) &
680  .and. (requests(i_req)%remote_bin == sent_bin)) then
681  found_request = .true.
682  exit
683  end if
684  end do
685  call assert(579308475, found_request)
686 
687  ! maybe do coagulation
688  call num_conc_weighted_kernel(coag_kernel_type, &
689  requests(i_req)%local_aero_particle, sent_aero_particle, &
690  s1, s2, sc, aero_data, aero_weight_total, env_state, k)
691  p = k * accept_factors(requests(i_req)%local_bin, sent_bin)
692 
693  if (pmc_random() .lt. p) then
694  ! coagulation happened, do it
695  tot_n_coag = tot_n_coag + 1
696  call coagulate_dist(aero_data, aero_state, &
697  requests(i_req)%local_aero_particle, sent_aero_particle, &
698  sent_proc, aero_weight_total, magnitudes, remove_1, remove_2)
699  else
700  remove_1 = .false.
701  remove_2 = .false.
702  end if
703 
704  ! send the particles back
705  if (.not. remove_1) then
706  ! If we wanted to, we could use the knowledge that this will go
707  ! into bin requests(i_req)%local_bin
708  call aero_state_add_particle(aero_state, &
709  requests(i_req)%local_aero_particle, 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  call aero_particle_deallocate(sent_aero_particle)
718 #endif
719 
720  end subroutine recv_return_req_particle
721 
722 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
723 
724  subroutine send_return_unreq_particle(aero_particle, dest_proc)
725 
726  !> Aero particle to send.
727  type(aero_particle_t), intent(in) :: aero_particle
728  !> Process to send the particle to.
729  integer, intent(in) :: dest_proc
730 
731 #ifdef PMC_USE_MPI
732  character :: buffer(coag_dist_max_buffer_size)
733  integer :: buffer_size, max_buffer_size, position, ierr
734 
735  max_buffer_size = 0
736  max_buffer_size = max_buffer_size &
737  + pmc_mpi_pack_size_aero_particle(aero_particle)
738  call assert(414990602, max_buffer_size <= coag_dist_max_buffer_size)
739  position = 0
740  call pmc_mpi_pack_aero_particle(buffer, position, aero_particle)
741  call assert(898537822, position <= max_buffer_size)
742  buffer_size = position
743  call mpi_bsend(buffer, buffer_size, mpi_character, dest_proc, &
744  coag_dist_tag_return_unreq_particle, mpi_comm_world, ierr)
745  call pmc_mpi_check_ierr(ierr)
746 #endif
747 
748  end subroutine send_return_unreq_particle
749 
750 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
751 
752  subroutine recv_return_unreq_particle(aero_state)
753 
754  !> Aerosol state.
755  type(aero_state_t), intent(inout) :: aero_state
756 
757 #ifdef PMC_USE_MPI
758  logical :: found_request
759  integer :: buffer_size, position, sent_proc, ierr
760  character :: buffer(coag_dist_max_buffer_size)
761  type(aero_particle_t) :: aero_particle
762  integer :: status(mpi_status_size), send_proc
763 
764  ! get the message
765  call mpi_recv(buffer, coag_dist_max_buffer_size, mpi_character, &
766  mpi_any_source, coag_dist_tag_return_unreq_particle, &
767  mpi_comm_world, status, ierr)
768  call pmc_mpi_check_ierr(ierr)
769  call assert(496247788, status(mpi_tag) &
770  == coag_dist_tag_return_unreq_particle)
771  call mpi_get_count(status, mpi_character, buffer_size, ierr)
772  call pmc_mpi_check_ierr(ierr)
773  call assert(590644042, buffer_size <= coag_dist_max_buffer_size)
774  sent_proc = status(mpi_source)
775 
776  ! unpack it
777  position = 0
778  call aero_particle_allocate(aero_particle)
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, &
784  allow_resort=.false.)
785  call aero_particle_deallocate(aero_particle)
786 #endif
787 
788  end subroutine recv_return_unreq_particle
789 
790 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
791 
792  !> Send a message saying that this process is finished with its
793  !> coagulation.
794  subroutine send_done(dest_proc)
795 
796  !> Process to send the message to.
797  integer, intent(in) :: dest_proc
798 
799 #ifdef PMC_USE_MPI
800  character :: buffer(coag_dist_max_buffer_size)
801  integer :: buffer_size, ierr
802 
803  buffer_size = 0
804  call mpi_bsend(buffer, buffer_size, mpi_character, dest_proc, &
805  coag_dist_tag_done, mpi_comm_world, ierr)
806  call pmc_mpi_check_ierr(ierr)
807 #endif
808 
809  end subroutine send_done
810 
811 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
812 
813  !> Receive a done message.
814  subroutine recv_done(procs_done)
815 
816  !> Which processes are finished with coagulation.
817  logical, intent(inout) :: procs_done(:)
818 
819 #ifdef PMC_USE_MPI
820  integer :: buffer_size, sent_proc, ierr
821  character :: buffer(coag_dist_max_buffer_size)
822  integer :: status(mpi_status_size)
823 
824  ! get the message
825  call mpi_recv(buffer, coag_dist_max_buffer_size, mpi_character, &
826  mpi_any_source, coag_dist_tag_done, mpi_comm_world, &
827  status, ierr)
828  call pmc_mpi_check_ierr(ierr)
829  call assert(348737947, status(mpi_tag) &
830  == coag_dist_tag_done)
831  call mpi_get_count(status, mpi_character, buffer_size, ierr)
832  call pmc_mpi_check_ierr(ierr)
833  call assert(214904056, buffer_size == 0)
834  sent_proc = status(mpi_source)
835 
836  ! process it
837  procs_done(sent_proc + 1) = .true.
838 #endif
839 
840  end subroutine recv_done
841 
842 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
843 
844  !> generate the number of samples to do per bin pair.
845  subroutine generate_n_samps(n_parts, del_t, bin_grid, aero_weight_array, &
846  k_max, n_samps, accept_factors)
847 
848  !> Number of particles per bin on all processes.
849  integer, intent(in) :: n_parts(:,:)
850  !> Timestep.
851  real(kind=dp), intent(in) :: del_t
852  !> Bin grid.
853  type(bin_grid_t), intent(in) :: bin_grid
854  !> Weighting function array.
855  type(aero_weight_array_t), intent(in) :: aero_weight_array
856  !> Maximum kernel.
857  real(kind=dp), intent(in) :: k_max(:,:)
858  !> Number of samples to do per bin pair.
859  integer, intent(out) :: n_samps(:,:)
860  !> Accept scale factors per bin pair (1).
861  real(kind=dp), intent(out) :: accept_factors(:,:)
862 
863  integer :: i_bin, j_bin, rank, n_bin
864  real(kind=dp) :: n_samp_mean
865 
866  n_bin = size(k_max, 1)
867  rank = pmc_mpi_rank()
868  n_samps = 0
869  do i_bin = 1,n_bin
870  if (n_parts(i_bin, rank + 1) == 0) &
871  cycle
872  do j_bin = i_bin,n_bin
873  call compute_n_samp(n_parts(i_bin, rank + 1), &
874  sum(n_parts(j_bin, :)), (i_bin == j_bin), &
875  k_max(i_bin, j_bin), del_t, n_samp_mean, &
876  n_samps(i_bin, j_bin), accept_factors(i_bin, j_bin))
877  end do
878  end do
879 
880  end subroutine generate_n_samps
881 
882 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
883 
884  subroutine coagulate_dist(aero_data, aero_state, aero_particle_1, &
885  aero_particle_2, remote_proc, aero_weight_total, magnitudes, &
886  remove_1, remove_2)
887 
888  !> Aerosol data.
889  type(aero_data_t), intent(in) :: aero_data
890  !> Aerosol state.
891  type(aero_state_t), intent(inout) :: aero_state
892  !> First particle to coagulate.
893  type(aero_particle_t), intent(in) :: aero_particle_1
894  !> Second particle to coagulate.
895  type(aero_particle_t), intent(in) :: aero_particle_2
896  !> Remote process that the particle came from.
897  integer, intent(in) :: remote_proc
898  !> Total weight across all processes.
899  type(aero_weight_array_t), intent(in) :: aero_weight_total
900  !> Computational volumes on all processes (m^3).
901  real(kind=dp), intent(in) :: magnitudes(:,:)
902  !> Whether to remove aero_particle_1 after the coagulation.
903  logical, intent(out) :: remove_1
904  !> Whether to remove aero_particle_2 after the coagulation.
905  logical, intent(out) :: remove_2
906 
907  integer, parameter :: s1 = 1
908  integer, parameter :: s2 = 1
909  integer, parameter :: sc = 1
910 
911  type(aero_particle_t) :: aero_particle_new
912  integer :: new_proc, new_group
913  type(aero_info_t) :: aero_info_1, aero_info_2
914  logical :: create_new, id_1_lost, id_2_lost
915 
916  call aero_particle_allocate(aero_particle_new)
917  call aero_info_allocate(aero_info_1)
918  call aero_info_allocate(aero_info_2)
919 
920  call coagulate_weighting(aero_particle_1, aero_particle_2, &
921  aero_particle_new, s1, s2, sc, aero_data, aero_state%awa, &
922  remove_1, remove_2, create_new, id_1_lost, id_2_lost, &
923  aero_info_1, aero_info_2)
924 
925  if (id_1_lost) then
926  call aero_info_array_add_aero_info(aero_state%aero_info_array, &
927  aero_info_1)
928  end if
929  if (id_2_lost) then
930  call aero_info_array_add_aero_info(aero_state%aero_info_array, &
931  aero_info_2)
932  end if
933 
934  ! add new particle
935  if (create_new) then
936  new_group = aero_weight_array_rand_group(aero_weight_total, sc, &
937  aero_particle_radius(aero_particle_new))
938  aero_particle_new%weight_group = new_group
939  new_proc = sample_cts_pdf(1d0 / magnitudes(new_group, :)) - 1
940  call send_return_unreq_particle(aero_particle_new, new_proc)
941  end if
942 
943  call aero_particle_deallocate(aero_particle_new)
944  call aero_info_deallocate(aero_info_1)
945  call aero_info_deallocate(aero_info_2)
946 
947  end subroutine coagulate_dist
948 
949 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
950 
951 end module pmc_coagulation_dist
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)
subroutine send_request_particle(remote_proc, remote_bin)
subroutine pmc_mpi_pack_aero_particle(buffer, position, val)
Packs the given value into the buffer, advancing position.
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.
The aero_data_t structure and associated subroutines.
Definition: aero_data.F90:9
integer function pmc_mpi_pack_size_integer(val)
Determines the number of bytes required to pack the given value.
Definition: mpi.F90:347
subroutine die_msg(code, error_msg)
Error immediately.
Definition: util.F90:133
subroutine aero_info_deallocate(aero_info)
Deallocates.
Definition: aero_info.F90:75
subroutine pmc_mpi_unpack_aero_particle(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
subroutine add_coagulation_requests(aero_state, requests, n_parts, local_bin, remote_bin, n_samps, samps_remaining)
logical function any_requests_active(requests)
Returns .true. if any of the requests are active, otherwise returns .false.
subroutine aero_info_allocate(aero_info)
Allocates and initializes.
Definition: aero_info.F90:63
subroutine num_conc_weighted_kernel(coag_kernel_type, aero_particle_1, aero_particle_2, i_class, j_class, ij_class, aero_data, aero_weight_array, env_state, k)
Compute the kernel value with the given number concentration weighting.
subroutine max_coag_num_conc_factor(aero_weight_array, bin_grid, i_bin, j_bin, i_class, j_class, ij_class, f_max)
Determine the minimum and maximum number concentration factors for coagulation.
subroutine aero_particle_deallocate(aero_particle)
Deallocates memory associated with an aero_particle_t.
The env_state_t structure and associated subroutines.
Definition: env_state.F90:9
Aerosol particle coagulation.
Definition: coagulation.F90:11
subroutine aero_state_add_particle(aero_state, aero_particle, allow_resort)
Add the given particle.
Definition: aero_state.F90:330
integer function sample_cts_pdf(pdf)
Sample the given continuous probability density function.
Definition: rand.F90:489
subroutine assert_msg(code, condition_ok, error_msg)
Errors unless condition_ok is true.
Definition: util.F90:76
subroutine request_deallocate(request)
Deallocate a request object and set it to be invalid.
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...
subroutine aero_info_array_add_aero_info(aero_info_array, aero_info)
Adds the given aero_info to the end of the array.
elemental real(kind=dp) function aero_particle_radius(aero_particle)
Total radius of the particle (m).
subroutine est_k_minmax_binned_unweighted(bin_grid, coag_kernel_type, aero_data, env_state, k_min, k_max)
Estimate an array of minimum and maximum kernel values. Given particles v1 in bin b1 and v2 in bin b2...
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...
Parallel aerosol particle coagulation with MPI.
subroutine send_return_req_particle(aero_particle, i_bin, dest_proc)
subroutine pmc_mpi_pack_integer(buffer, position, val)
Packs the given value into the buffer, advancing position.
Definition: mpi.F90:534
integer function pmc_mpi_size()
Returns the total number of processes.
Definition: mpi.F90:133
subroutine coagulate_dist(aero_data, aero_state, aero_particle_1, aero_particle_2, remote_proc, aero_weight_total, magnitudes, remove_1, remove_2)
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.
Common utility subroutines.
Definition: util.F90:9
subroutine request_allocate(request)
Current environment state.
Definition: env_state.F90:26
subroutine aero_state_sort(aero_state, bin_grid, all_procs_same)
Sorts the particles if necessary.
subroutine aero_weight_array_allocate(aero_weight_array)
Allocates an aero_weight_array.
subroutine send_return_unreq_particle(aero_particle, dest_proc)
subroutine pmc_mpi_unpack_integer(buffer, position, val)
Unpacks the given value from the buffer, advancing position.
Definition: mpi.F90:771
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:414
An array of aerosol size distribution weighting functions.
The aero_state_t structure and assocated subroutines.
Definition: aero_state.F90:9
integer function pmc_mpi_rank()
Returns the rank of the current process.
Definition: mpi.F90:116
subroutine recv_return_no_particle(requests, aero_data, aero_state)
subroutine send_done(dest_proc)
Send a message saying that this process is finished with its coagulation.
Wrapper functions for MPI.
Definition: mpi.F90:13
integer function aero_sorted_n_class(aero_sorted)
Returns the number of weight classes.
The aero_weight_array_t structure and associated subroutines.
subroutine update_n_samps(n_samps, local_bin, remote_bin, samps_remaining)
subroutine aero_particle_allocate(aero_particle)
Allocates memory in an aero_particle_t.
character(len=pmc_util_convert_string_len) function integer_to_string(val)
Convert an integer to a string format.
Definition: util.F90:743
The current collection of aerosol particles.
Definition: aero_state.F90:63
Single aerosol particle data structure.
1D grid, either logarithmic or linear.
Definition: bin_grid.F90:33
subroutine pmc_mpi_check_ierr(ierr)
Dies if ierr is not ok.
Definition: mpi.F90:39
The bin_grid_t structure and associated subroutines.
Definition: bin_grid.F90:9
subroutine find_rand_remote_proc(n_parts, remote_bin, remote_proc)
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.
subroutine pmc_mpi_barrier()
Synchronize all processes.
Definition: mpi.F90:102
A single outstanding request for a remote particle.
subroutine recv_return_unreq_particle(aero_state)
integer function sample_disc_pdf(pdf)
Sample the given discrete probability density function.
Definition: rand.F90:523
subroutine recv_request_particle(aero_state)
real(kind=dp) function pmc_random()
Returns a random number between 0 and 1.
Definition: rand.F90:138
integer function pmc_mpi_pack_size_aero_particle(val)
Determines the number of bytes required to pack the given value.
Aerosol material properties and associated data.
Definition: aero_data.F90:40
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:1464
subroutine assert(code, condition_ok)
Errors unless condition_ok is true.
Definition: util.F90:102
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:1429
subroutine recv_done(procs_done)
Receive a done message.
subroutine recv_return_req_particle(requests, env_state, aero_weight_total, aero_data, aero_state, accept_factors, coag_kernel_type, tot_n_coag, magnitudes)
logical function request_is_active(request)
Whether the given reqest object is currectly active.
subroutine send_return_no_particle(dest_proc, i_bin)
subroutine aero_weight_array_copy(aero_weight_array_from, aero_weight_array_to)
Copy an aero_weight_array.