PartMC  2.2.1
coagulation_dist.F90
Go to the documentation of this file.
00001 ! Copyright (C) 2005-2011 Nicole Riemer and Matthew West
00002 ! Licensed under the GNU General Public License version 2 or (at your
00003 ! option) any later version. See the file COPYING for details.
00004 
00005 !> \file
00006 !> The pmc_coagulation_dist module.
00007 
00008 !> Parallel aerosol particle coagulation with MPI
00009 module pmc_coagulation_dist
00010 
00011   use pmc_bin_grid
00012   use pmc_aero_data
00013   use pmc_util
00014   use pmc_env_state
00015   use pmc_aero_state
00016   use pmc_coagulation
00017   use pmc_mpi
00018 #ifdef PMC_USE_MPI
00019   use mpi
00020 #endif
00021 
00022   !> Size of the outgoing buffer for \c bsend (bytes).
00023   !!
00024   !! FIXME: check that this size is big enough. It must be large
00025   !! enough to handle the required number of messages of the given
00026   !! sizes, plus MPI_BSEND_OVERHEAD per message, plus some extra room
00027   !! because it's only kind of a circular buffer --- the messages
00028   !! themselves aren't allowed to wrap around then end, so we might
00029   !! need extra space up to the size of the largest message type.
00030   integer, parameter :: COAG_DIST_OUTGOING_BUFFER_SIZE      = 1000000
00031   !> Size of send and receive buffer for each message (bytes).
00032   !!
00033   !! The biggest message type will be one of the particle-sending
00034   !! types, for which we need pmc_mpi_pack_size_aero_particle(), plus
00035   !! a couple of integers or something. At the moment this means
00036   !! something like (10 + n_spec) reals, (3 + 2) integers, which for
00037   !! n_spec = 20 gives a size of 260 bytes.
00038   integer, parameter :: COAG_DIST_MAX_BUFFER_SIZE           = 10000
00039   integer, parameter :: COAG_DIST_MAX_REQUESTS              = 1
00040   integer, parameter :: COAG_DIST_TAG_REQUEST_PARTICLE      = 5321
00041   integer, parameter :: COAG_DIST_TAG_RETURN_REQ_PARTICLE   = 5322
00042   integer, parameter :: COAG_DIST_TAG_RETURN_UNREQ_PARTICLE = 5323
00043   integer, parameter :: COAG_DIST_TAG_RETURN_NO_PARTICLE    = 5324
00044   integer, parameter :: COAG_DIST_TAG_DONE                  = 5325
00045 
00046   !> A single outstanding request for a remote particle.
00047   type request_t
00048      !> Local \c aero_particle to maybe coagulate with the received
00049      !> particle.
00050      type(aero_particle_t) :: local_aero_particle
00051      !> Remote process number that we sent the request to
00052      !> (-1 means this request is currently not used).
00053      integer :: remote_proc
00054      !> Local bin number from which we took \c local_aero_particle.
00055      integer :: local_bin
00056      !> Remote bin number from which we requested an \c aero_particle.
00057      integer :: remote_bin
00058      !> Whether this request is currently active
00059      logical :: active
00060   end type request_t
00061   
00062 contains
00063 
00064 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00065 
00066   ! Allocate a request object and set its state to invalid.
00067   subroutine request_allocate(request)
00068 
00069     !> Request object to allocate.
00070     type(request_t), intent(out) :: request
00071 
00072     call aero_particle_allocate(request%local_aero_particle)
00073     request%active = .false.
00074 
00075   end subroutine request_allocate
00076 
00077 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00078 
00079   !> Deallocate a request object and set it to be invalid.
00080   subroutine request_deallocate(request)
00081    
00082     !> Request object to deallocate
00083     type(request_t), intent(inout) :: request
00084 
00085     call aero_particle_deallocate(request%local_aero_particle)
00086     request%active = .false.
00087 
00088   end subroutine request_deallocate
00089 
00090 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00091 
00092   !> Whether the given reqest object is currectly active.
00093   logical function request_is_active(request)
00094 
00095     !> Request object to test for activeness.
00096     type(request_t), intent(in) :: request
00097 
00098     request_is_active = request%active
00099 
00100   end function request_is_active
00101 
00102 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00103 
00104   !> Do coagulation for time del_t.
00105   subroutine mc_coag_dist(coag_kernel_type, env_state, aero_data, &
00106        aero_state, del_t, tot_n_samp, tot_n_coag)
00107 
00108     !> Coagulation kernel type.
00109     integer, intent(in) :: coag_kernel_type
00110     !> Environment state.
00111     type(env_state_t), intent(in) :: env_state
00112     !> Aerosol data.
00113     type(aero_data_t), intent(in) :: aero_data
00114     !> Aerosol state.
00115     type(aero_state_t), intent(inout) :: aero_state
00116     !> Timestep.
00117     real(kind=dp), intent(in) :: del_t
00118     !> Total number of samples tested.
00119     integer, intent(out) :: tot_n_samp
00120     !> Number of coagulation events.
00121     integer, intent(out) :: tot_n_coag
00122 
00123 #ifdef PMC_USE_MPI
00124     logical :: samps_remaining, sent_dones
00125     integer :: i_bin, j_bin, n_samp, i_samp, i_proc, n_proc
00126     integer :: ierr, status(MPI_STATUS_SIZE), current_i, current_j, i_req
00127     real(kind=dp) :: n_samp_real, f_max
00128     integer, allocatable :: n_parts(:,:)
00129     real(kind=dp), allocatable :: comp_vols(:,:)
00130     type(request_t) :: requests(COAG_DIST_MAX_REQUESTS)
00131     integer, allocatable :: n_samps(:,:)
00132     real(kind=dp), allocatable :: accept_factors(:,:), k_max(:,:)
00133     logical, allocatable :: procs_done(:)
00134     integer :: outgoing_buffer(COAG_DIST_OUTGOING_BUFFER_SIZE)
00135     integer :: outgoing_buffer_size_check
00136     type(aero_weight_t), allocatable :: aero_weight_total(:)
00137     
00138     n_proc = pmc_mpi_size()
00139 
00140     call pmc_mpi_barrier()
00141 
00142     call aero_state_sort(aero_state, all_procs_same=.true.)
00143     if (.not. aero_state%aero_sorted%coag_kernel_bounds_valid) then
00144        call est_k_minmax_binned_unweighted(aero_state%aero_sorted%bin_grid, &
00145             coag_kernel_type, aero_data, env_state, &
00146             aero_state%aero_sorted%coag_kernel_min, &
00147             aero_state%aero_sorted%coag_kernel_max)
00148        aero_state%aero_sorted%coag_kernel_bounds_valid = .true.
00149     end if
00150 
00151     allocate(n_samps(aero_state%aero_sorted%bin_grid%n_bin, &
00152          aero_state%aero_sorted%bin_grid%n_bin))
00153     allocate(accept_factors(aero_state%aero_sorted%bin_grid%n_bin, &
00154          aero_state%aero_sorted%bin_grid%n_bin))
00155 
00156     allocate(n_parts(aero_state%aero_sorted%bin_grid%n_bin, n_proc))
00157     call pmc_mpi_allgather_integer_array( &
00158          aero_state%aero_sorted%size%inverse(:)%n_entry, n_parts)
00159 
00160     allocate(comp_vols(size(aero_state%aero_weight), n_proc))
00161     call pmc_mpi_allgather_real_array(aero_state%aero_weight%comp_vol, &
00162          comp_vols)
00163 
00164     allocate(aero_weight_total(size(aero_state%aero_weight)))
00165     call aero_weight_allocate(aero_weight_total)
00166     call aero_weight_array_copy(aero_state%aero_weight, aero_weight_total)
00167     aero_weight_total%comp_vol = sum(comp_vols, 2)
00168 
00169     allocate(k_max(aero_state%aero_sorted%bin_grid%n_bin, &
00170          aero_state%aero_sorted%bin_grid%n_bin))
00171     do i_bin = 1,aero_state%aero_sorted%bin_grid%n_bin
00172        do j_bin = 1,aero_state%aero_sorted%bin_grid%n_bin
00173           call max_coag_num_conc_factor_better(aero_weight_total, &
00174                aero_state%aero_sorted%bin_grid, i_bin, j_bin, f_max)
00175           k_max(i_bin, j_bin) &
00176                = aero_state%aero_sorted%coag_kernel_max(i_bin, j_bin) * f_max
00177        end do
00178     end do
00179 
00180     call generate_n_samps(n_parts, del_t, aero_state%aero_sorted%bin_grid, &
00181          aero_weight_total, k_max, n_samps, accept_factors)
00182     tot_n_samp = sum(n_samps)
00183     tot_n_coag = 0
00184 
00185     ! main loop
00186     do i_req = 1,COAG_DIST_MAX_REQUESTS
00187        call request_allocate(requests(i_req))
00188     end do
00189     samps_remaining = .true.
00190     current_i = 1
00191     current_j = 1
00192     allocate(procs_done(n_proc))
00193     procs_done = .false.
00194     sent_dones = .false.
00195     call mpi_buffer_attach(outgoing_buffer, &
00196          COAG_DIST_OUTGOING_BUFFER_SIZE, ierr)
00197     call pmc_mpi_check_ierr(ierr)
00198     do while (.not. all(procs_done))
00199        ! add requests if we have any slots available call
00200        call add_coagulation_requests(aero_state, requests, n_parts, &
00201             current_i, current_j, n_samps, samps_remaining)
00202 
00203        ! if we have no outstanding requests, send done messages
00204        if (.not. sent_dones) then
00205           if (.not. any_requests_active(requests)) then
00206              sent_dones = .true.
00207              do i_proc = 0, (n_proc - 1)
00208                 call send_done(i_proc)
00209              end do
00210           end if
00211        end if
00212 
00213        ! receive exactly one message
00214        call coag_dist_recv(requests, env_state, aero_weight_total, aero_data, &
00215             aero_state, accept_factors, coag_kernel_type, tot_n_coag, &
00216             comp_vols, procs_done)
00217     end do
00218 
00219     do i_req = 1,COAG_DIST_MAX_REQUESTS
00220        call assert(502009333, .not. request_is_active(requests(i_req)))
00221        call request_deallocate(requests(i_req))
00222     end do
00223     deallocate(procs_done)
00224     deallocate(n_samps)
00225     deallocate(accept_factors)
00226     deallocate(n_parts)
00227     deallocate(comp_vols)
00228     call mpi_buffer_detach(outgoing_buffer, &
00229          outgoing_buffer_size_check, ierr)
00230     call pmc_mpi_check_ierr(ierr)
00231     call assert(577822730, &
00232          COAG_DIST_OUTGOING_BUFFER_SIZE == outgoing_buffer_size_check)
00233 #endif
00234 
00235   end subroutine mc_coag_dist
00236 
00237 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00238 
00239   subroutine coag_dist_recv(requests, env_state, aero_weight_total, &
00240        aero_data, aero_state, accept_factors, coag_kernel_type, tot_n_coag, &
00241        comp_vols, procs_done)
00242 
00243     !> Array of outstanding requests.
00244     type(request_t), intent(inout) :: requests(COAG_DIST_MAX_REQUESTS)
00245     !> Environment state.
00246     type(env_state_t), intent(in) :: env_state
00247     !> Total weighting functions.
00248     type(aero_weight_t), intent(in) :: aero_weight_total(:)
00249     !> Aerosol data.
00250     type(aero_data_t), intent(in) :: aero_data
00251     !> Aerosol state.
00252     type(aero_state_t), intent(inout) :: aero_state
00253     !> Accept scale factors per bin pair (1).
00254     real(kind=dp), intent(in) :: accept_factors(:,:)
00255     !> Coagulation kernel type.
00256     integer, intent(in) :: coag_kernel_type
00257     !> Number of coagulation events.
00258     integer, intent(inout) :: tot_n_coag
00259     !> Computational volumes on all processes.
00260     real(kind=dp), intent(in) :: comp_vols(:,:)
00261     !> Which processes are finished with coagulation.
00262     logical, intent(inout) :: procs_done(:)
00263 
00264 #ifdef PMC_USE_MPI
00265     integer :: status(MPI_STATUS_SIZE), ierr
00266 
00267     call mpi_probe(MPI_ANY_SOURCE, MPI_ANY_TAG, MPI_COMM_WORLD, &
00268          status, ierr)
00269     call pmc_mpi_check_ierr(ierr)
00270     if (status(MPI_TAG) == COAG_DIST_TAG_REQUEST_PARTICLE) then
00271        call recv_request_particle(aero_state)
00272     elseif (status(MPI_TAG) == COAG_DIST_TAG_RETURN_REQ_PARTICLE) then
00273        call recv_return_req_particle(requests, env_state, aero_weight_total, &
00274             aero_data, aero_state, accept_factors, coag_kernel_type, &
00275             tot_n_coag, comp_vols)
00276     elseif (status(MPI_TAG) == COAG_DIST_TAG_RETURN_UNREQ_PARTICLE) then
00277        call recv_return_unreq_particle(aero_state)
00278     elseif (status(MPI_TAG) == COAG_DIST_TAG_RETURN_NO_PARTICLE) then
00279        call recv_return_no_particle(requests, aero_data, aero_state)
00280     elseif (status(MPI_TAG) == COAG_DIST_TAG_DONE) then
00281        call recv_done(procs_done)
00282     else
00283        call die_msg(856123972, &
00284             'unknown tag: ' // trim(integer_to_string(status(MPI_TAG))))
00285     end if
00286 #endif
00287     
00288   end subroutine coag_dist_recv
00289 
00290 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00291 
00292   subroutine add_coagulation_requests(aero_state, requests, n_parts, &
00293        local_bin, remote_bin, n_samps, samps_remaining)
00294 
00295     !> Aerosol state.
00296     type(aero_state_t), intent(inout) :: aero_state
00297     !> Array of outstanding requests.
00298     type(request_t), intent(inout) :: requests(COAG_DIST_MAX_REQUESTS)
00299     !> Number of particles per bin per process.
00300     integer, intent(in) :: n_parts(:,:)
00301     !> Bin index of first particle we need to coagulate.
00302     integer, intent(inout) :: local_bin
00303     !> Bin index of second particle we need to coagulate.
00304     integer, intent(inout) :: remote_bin
00305     !> Number of samples remaining per bin pair
00306     integer, intent(inout) :: n_samps(:,:)
00307     !> Whether there are still coagulation samples that need to be done.
00308     logical, intent(inout) :: samps_remaining
00309 
00310     integer :: i_req
00311 
00312     if (.not. samps_remaining) return
00313 
00314     outer: do i_req = 1,COAG_DIST_MAX_REQUESTS
00315        if (.not. request_is_active(requests(i_req))) then
00316           inner: do
00317              call update_n_samps(n_samps, local_bin, remote_bin, &
00318                   samps_remaining)
00319              if (.not. samps_remaining) exit outer
00320              if (aero_state%aero_sorted%size%inverse(local_bin)%n_entry > 0) then
00321                 call find_rand_remote_proc(n_parts, remote_bin, &
00322                      requests(i_req)%remote_proc)
00323                 requests(i_req)%active = .true.
00324                 requests(i_req)%local_bin = local_bin
00325                 requests(i_req)%remote_bin = remote_bin
00326                 call aero_state_remove_rand_particle_from_bin(aero_state, &
00327                      local_bin, requests(i_req)%local_aero_particle)
00328                 call send_request_particle(requests(i_req)%remote_proc, &
00329                      requests(i_req)%remote_bin)
00330                 exit inner
00331              end if
00332           end do inner
00333        end if
00334     end do outer
00335 
00336   end subroutine add_coagulation_requests
00337 
00338 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00339 
00340   !> Returns \c .true. if any of the requests are active, otherwise
00341   !> returns \c .false.
00342   logical function any_requests_active(requests)
00343 
00344     !> Array of outstanding requests.
00345     type(request_t), intent(inout) :: requests(COAG_DIST_MAX_REQUESTS)
00346     
00347     integer :: i_req
00348 
00349     do i_req = 1,COAG_DIST_MAX_REQUESTS
00350        if (request_is_active(requests(i_req))) then
00351           any_requests_active = .true.
00352           return
00353        end if
00354     end do
00355     any_requests_active = .false.
00356     
00357   end function any_requests_active
00358 
00359 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00360 
00361   subroutine find_rand_remote_proc(n_parts, remote_bin, remote_proc)
00362 
00363     !> Number of particles per bin per process.
00364     integer, intent(in) :: n_parts(:,:)
00365     !> Remote bin number.
00366     integer, intent(in) :: remote_bin
00367     !> Remote process number chosen at random.
00368     integer, intent(out) :: remote_proc
00369 
00370 #ifdef PMC_USE_MPI
00371     call assert(770964285, size(n_parts, 2) == pmc_mpi_size())
00372     remote_proc = sample_disc_pdf(n_parts(remote_bin, :)) - 1
00373 #endif
00374 
00375   end subroutine find_rand_remote_proc
00376 
00377 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00378 
00379   subroutine update_n_samps(n_samps, local_bin, remote_bin, samps_remaining)
00380 
00381     !> Number of samples remaining per bin pair
00382     integer, intent(inout) :: n_samps(:,:)
00383     !> Bin index of first particle we need to coagulate.
00384     integer, intent(inout) :: local_bin
00385     !> Bin index of second particle we need to coagulate.
00386     integer, intent(inout) :: remote_bin
00387     !> Whether there are still coagulation samples that need to be done.
00388     logical, intent(inout) :: samps_remaining
00389 
00390     integer :: n_bin
00391 
00392     if (.not. samps_remaining) return
00393 
00394     n_bin = size(n_samps, 1)
00395     do
00396        if (n_samps(local_bin, remote_bin) > 0) exit
00397 
00398        remote_bin = remote_bin + 1
00399        if (remote_bin > n_bin) then
00400           remote_bin = 1
00401           local_bin = local_bin + 1
00402        end if
00403        if (local_bin > n_bin) exit
00404     end do
00405     
00406     if (local_bin > n_bin) then
00407        samps_remaining = .false.
00408     else
00409        n_samps(local_bin, remote_bin) = n_samps(local_bin, remote_bin) - 1
00410     end if
00411 
00412   end subroutine update_n_samps
00413 
00414 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00415 
00416   subroutine send_request_particle(remote_proc, remote_bin)
00417 
00418     !> Remote process number.
00419     integer, intent(in) :: remote_proc
00420     !> Remote bin number.
00421     integer, intent(in) :: remote_bin
00422 
00423 #ifdef PMC_USE_MPI
00424     character :: buffer(COAG_DIST_MAX_BUFFER_SIZE)
00425     integer :: buffer_size, max_buffer_size, position, ierr
00426 
00427     max_buffer_size = 0
00428     max_buffer_size = max_buffer_size + pmc_mpi_pack_size_integer(remote_bin)
00429     call assert(893545122, max_buffer_size <= COAG_DIST_MAX_BUFFER_SIZE)
00430     position = 0
00431     call pmc_mpi_pack_integer(buffer, position, remote_bin)
00432     call assert(610314213, position <= max_buffer_size)
00433     buffer_size = position
00434     call mpi_bsend(buffer, buffer_size, MPI_CHARACTER, remote_proc, &
00435          COAG_DIST_TAG_REQUEST_PARTICLE, MPI_COMM_WORLD, ierr)
00436     call pmc_mpi_check_ierr(ierr)
00437 #endif
00438 
00439   end subroutine send_request_particle
00440 
00441 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00442 
00443   subroutine recv_request_particle(aero_state)
00444 
00445     !> Aero state.
00446     type(aero_state_t), intent(inout) :: aero_state
00447 
00448 #ifdef PMC_USE_MPI
00449     integer :: buffer_size, position, request_bin, sent_proc
00450     integer :: ierr, remote_proc, status(MPI_STATUS_SIZE)
00451     character :: buffer(COAG_DIST_MAX_BUFFER_SIZE)
00452     type(aero_particle_t) :: aero_particle
00453 
00454     ! get the message
00455     call mpi_recv(buffer, COAG_DIST_MAX_BUFFER_SIZE, MPI_CHARACTER, &
00456          MPI_ANY_SOURCE, COAG_DIST_TAG_REQUEST_PARTICLE, MPI_COMM_WORLD, &
00457          status, ierr)
00458     call pmc_mpi_check_ierr(ierr)
00459     call assert(920139874, status(MPI_TAG) &
00460          == COAG_DIST_TAG_REQUEST_PARTICLE)
00461     call mpi_get_count(status, MPI_CHARACTER, buffer_size, ierr)
00462     call pmc_mpi_check_ierr(ierr)
00463     call assert(190658659, buffer_size <= COAG_DIST_MAX_BUFFER_SIZE)
00464     remote_proc = status(MPI_SOURCE)
00465 
00466     ! unpack it
00467     position = 0
00468     call pmc_mpi_unpack_integer(buffer, position, request_bin)
00469     call assert(895128380, position == buffer_size)
00470 
00471     ! send the particle back if we have one
00472     if (aero_state%aero_sorted%size%inverse(request_bin)%n_entry == 0) &
00473          then
00474        call send_return_no_particle(remote_proc, request_bin)
00475     else
00476        call aero_particle_allocate(aero_particle)
00477        call aero_state_remove_rand_particle_from_bin(aero_state, &
00478             request_bin, aero_particle)
00479        call send_return_req_particle(aero_particle, request_bin, &
00480             remote_proc)
00481        call aero_particle_deallocate(aero_particle)
00482     end if
00483 #endif
00484 
00485   end subroutine recv_request_particle
00486 
00487 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00488 
00489   subroutine send_return_no_particle(dest_proc, i_bin)
00490 
00491     !> Process number to send message to.
00492     integer, intent(in) :: dest_proc
00493     !> Bin number where there was no particle.
00494     integer, intent(in) :: i_bin
00495 
00496 #ifdef PMC_USE_MPI
00497     character :: buffer(COAG_DIST_MAX_BUFFER_SIZE)
00498     integer :: buffer_size, max_buffer_size, position, ierr
00499 
00500     max_buffer_size = 0
00501     max_buffer_size = max_buffer_size + pmc_mpi_pack_size_integer(i_bin)
00502     call assert(744787119, max_buffer_size <= COAG_DIST_MAX_BUFFER_SIZE)
00503     position = 0
00504     call pmc_mpi_pack_integer(buffer, position, i_bin)
00505     call assert(445960340, position <= max_buffer_size)
00506     buffer_size = position
00507     call mpi_bsend(buffer, buffer_size, MPI_CHARACTER, dest_proc, &
00508          COAG_DIST_TAG_RETURN_NO_PARTICLE, MPI_COMM_WORLD, ierr)
00509     call pmc_mpi_check_ierr(ierr)
00510 #endif
00511 
00512   end subroutine send_return_no_particle
00513 
00514 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00515 
00516   subroutine recv_return_no_particle(requests, aero_data, aero_state)
00517 
00518     !> Array of outstanding requests.
00519     type(request_t), intent(inout) :: requests(COAG_DIST_MAX_REQUESTS)
00520     !> Aerosol data.
00521     type(aero_data_t), intent(in) :: aero_data
00522     !> Aerosol state.
00523     type(aero_state_t), intent(inout) :: aero_state
00524 
00525 #ifdef PMC_USE_MPI
00526     logical :: found_request
00527     integer :: buffer_size, position, sent_bin, sent_proc, i_req
00528     integer :: ierr, status(MPI_STATUS_SIZE)
00529     character :: buffer(COAG_DIST_MAX_BUFFER_SIZE)
00530 
00531     ! get the message
00532     call mpi_recv(buffer, COAG_DIST_MAX_BUFFER_SIZE, MPI_CHARACTER, &
00533          MPI_ANY_SOURCE, COAG_DIST_TAG_RETURN_NO_PARTICLE, &
00534          MPI_COMM_WORLD, status, ierr)
00535     call pmc_mpi_check_ierr(ierr)
00536     call assert(918153221, status(MPI_TAG) &
00537          == COAG_DIST_TAG_RETURN_NO_PARTICLE)
00538     call mpi_get_count(status, MPI_CHARACTER, buffer_size, ierr)
00539     call pmc_mpi_check_ierr(ierr)
00540     call assert(461111487, buffer_size <= COAG_DIST_MAX_BUFFER_SIZE)
00541     sent_proc = status(MPI_SOURCE)
00542 
00543     ! unpack it
00544     position = 0
00545     call pmc_mpi_unpack_integer(buffer, position, sent_bin)
00546     call assert(518172999, position == buffer_size)
00547 
00548     ! find the matching request
00549     found_request = .false.
00550     do i_req = 1,COAG_DIST_MAX_REQUESTS
00551        if ((requests(i_req)%remote_proc == sent_proc) &
00552             .and. (requests(i_req)%remote_bin == sent_bin)) then
00553           found_request = .true.
00554           exit
00555        end if
00556     end do
00557     call assert(215612776, found_request)
00558 
00559     ! We can't do coagulation with the local particle, so store it
00560     ! back. If we wanted to, we could use the knowledge that it should
00561     ! go into bin requests(i_req)%local_bin
00562     call aero_state_add_particle(aero_state, &
00563          requests(i_req)%local_aero_particle, allow_resort=.false.)
00564     call request_deallocate(requests(i_req))
00565     call request_allocate(requests(i_req))
00566 #endif
00567 
00568   end subroutine recv_return_no_particle
00569 
00570 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00571 
00572   subroutine send_return_req_particle(aero_particle, i_bin, dest_proc)
00573 
00574     !> Aero particle to send.
00575     type(aero_particle_t), intent(in) :: aero_particle
00576     !> Bin that the particle is in.
00577     integer, intent(in) :: i_bin
00578     !> Process number to send particle to.
00579     integer, intent(in) :: dest_proc
00580 
00581 #ifdef PMC_USE_MPI
00582     character :: buffer(COAG_DIST_MAX_BUFFER_SIZE)
00583     integer :: buffer_size, max_buffer_size, position, ierr
00584 
00585     max_buffer_size = 0
00586     max_buffer_size = max_buffer_size + pmc_mpi_pack_size_integer(i_bin)
00587     max_buffer_size = max_buffer_size &
00588          + pmc_mpi_pack_size_aero_particle(aero_particle)
00589     call assert(496283814, max_buffer_size <= COAG_DIST_MAX_BUFFER_SIZE)
00590     position = 0
00591     call pmc_mpi_pack_integer(buffer, position, i_bin)
00592     call pmc_mpi_pack_aero_particle(buffer, position, aero_particle)
00593     call assert(263666386, position <= max_buffer_size)
00594     buffer_size = position
00595     call mpi_bsend(buffer, buffer_size, MPI_CHARACTER, dest_proc, &
00596          COAG_DIST_TAG_RETURN_REQ_PARTICLE, MPI_COMM_WORLD, ierr)
00597     call pmc_mpi_check_ierr(ierr)
00598 #endif
00599 
00600   end subroutine send_return_req_particle
00601 
00602 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00603 
00604   subroutine recv_return_req_particle(requests, env_state, aero_weight_total, &
00605        aero_data, aero_state, accept_factors, coag_kernel_type, tot_n_coag, &
00606        comp_vols)
00607 
00608     !> Array of outstanding requests.
00609     type(request_t), intent(inout) :: requests(COAG_DIST_MAX_REQUESTS)
00610     !> Environment state.
00611     type(env_state_t), intent(in) :: env_state
00612     !> Total weighting array.
00613     type(aero_weight_t), intent(in) :: aero_weight_total(:)
00614     !> Aerosol data.
00615     type(aero_data_t), intent(in) :: aero_data
00616     !> Aerosol state.
00617     type(aero_state_t), intent(inout) :: aero_state
00618     !> Accept scale factors per bin pair (1).
00619     real(kind=dp), intent(in) :: accept_factors(:,:)
00620     !> Coagulation kernel type.
00621     integer, intent(in) :: coag_kernel_type
00622     !> Number of coagulation events.
00623     integer, intent(inout) :: tot_n_coag
00624     !> Computational volumes on all processes.
00625     real(kind=dp), intent(in) :: comp_vols(:,:)
00626 
00627 #ifdef PMC_USE_MPI
00628     logical :: found_request, remove_1, remove_2
00629     integer :: buffer_size, position, sent_bin, sent_proc, i_req
00630     integer :: ierr, status(MPI_STATUS_SIZE)
00631     character :: buffer(COAG_DIST_MAX_BUFFER_SIZE)
00632     type(aero_particle_t) :: sent_aero_particle
00633     real(kind=dp) :: k, p
00634 
00635     ! get the message
00636     call mpi_recv(buffer, COAG_DIST_MAX_BUFFER_SIZE, MPI_CHARACTER, &
00637          MPI_ANY_SOURCE, COAG_DIST_TAG_RETURN_REQ_PARTICLE, &
00638          MPI_COMM_WORLD, status, ierr)
00639     call pmc_mpi_check_ierr(ierr)
00640     call assert(133285061, status(MPI_TAG) &
00641          == COAG_DIST_TAG_RETURN_REQ_PARTICLE)
00642     call mpi_get_count(status, MPI_CHARACTER, buffer_size, ierr)
00643     call pmc_mpi_check_ierr(ierr)
00644     call assert(563012836, buffer_size <= COAG_DIST_MAX_BUFFER_SIZE)
00645     sent_proc = status(MPI_SOURCE)
00646 
00647     ! unpack it
00648     position = 0
00649     call pmc_mpi_unpack_integer(buffer, position, sent_bin)
00650     call aero_particle_allocate(sent_aero_particle)
00651     call pmc_mpi_unpack_aero_particle(buffer, position, sent_aero_particle)
00652     call assert(753356021, position == buffer_size)
00653 
00654     ! find the matching request
00655     found_request = .false.
00656     do i_req = 1,COAG_DIST_MAX_REQUESTS
00657        if ((requests(i_req)%remote_proc == sent_proc) &
00658             .and. (requests(i_req)%remote_bin == sent_bin)) then
00659           found_request = .true.
00660           exit
00661        end if
00662     end do
00663     call assert(579308475, found_request)
00664     
00665     ! maybe do coagulation
00666     call num_conc_weighted_kernel(coag_kernel_type, &
00667          requests(i_req)%local_aero_particle, sent_aero_particle, &
00668          aero_data, aero_weight_total, env_state, k)
00669     p = k * accept_factors(requests(i_req)%local_bin, sent_bin)
00670 
00671     if (pmc_random() .lt. p) then
00672        ! coagulation happened, do it
00673        tot_n_coag = tot_n_coag + 1
00674        call coagulate_dist(aero_data, aero_state, &
00675             requests(i_req)%local_aero_particle, sent_aero_particle, &
00676             sent_proc, aero_weight_total, comp_vols, remove_1, remove_2)
00677     else
00678        remove_1 = .false.
00679        remove_2 = .false.
00680     end if
00681 
00682     ! send the particles back
00683     if (.not. remove_1) then
00684        ! If we wanted to, we could use the knowledge that this will go
00685        ! into bin requests(i_req)%local_bin
00686        call aero_state_add_particle(aero_state, &
00687             requests(i_req)%local_aero_particle, allow_resort=.false.)
00688     end if
00689     if (.not. remove_2) then
00690        call send_return_unreq_particle(sent_aero_particle, sent_proc)
00691     end if
00692 
00693     call request_deallocate(requests(i_req))
00694     call request_allocate(requests(i_req))
00695     call aero_particle_deallocate(sent_aero_particle)
00696 #endif
00697 
00698   end subroutine recv_return_req_particle
00699 
00700 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00701 
00702   subroutine send_return_unreq_particle(aero_particle, dest_proc)
00703 
00704     !> Aero particle to send.
00705     type(aero_particle_t), intent(in) :: aero_particle
00706     !> Process to send the particle to.
00707     integer, intent(in) :: dest_proc
00708 
00709 #ifdef PMC_USE_MPI
00710     character :: buffer(COAG_DIST_MAX_BUFFER_SIZE)
00711     integer :: buffer_size, max_buffer_size, position, ierr
00712 
00713     max_buffer_size = 0
00714     max_buffer_size = max_buffer_size &
00715          + pmc_mpi_pack_size_aero_particle(aero_particle)
00716     call assert(414990602, max_buffer_size <= COAG_DIST_MAX_BUFFER_SIZE)
00717     position = 0
00718     call pmc_mpi_pack_aero_particle(buffer, position, aero_particle)
00719     call assert(898537822, position <= max_buffer_size)
00720     buffer_size = position
00721     call mpi_bsend(buffer, buffer_size, MPI_CHARACTER, dest_proc, &
00722          COAG_DIST_TAG_RETURN_UNREQ_PARTICLE, MPI_COMM_WORLD, ierr)
00723     call pmc_mpi_check_ierr(ierr)
00724 #endif
00725 
00726   end subroutine send_return_unreq_particle
00727 
00728 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00729 
00730   subroutine recv_return_unreq_particle(aero_state)
00731 
00732     !> Aerosol state.
00733     type(aero_state_t), intent(inout) :: aero_state
00734 
00735 #ifdef PMC_USE_MPI
00736     logical :: found_request
00737     integer :: buffer_size, position, sent_proc, ierr
00738     character :: buffer(COAG_DIST_MAX_BUFFER_SIZE)
00739     type(aero_particle_t) :: aero_particle
00740     integer :: status(MPI_STATUS_SIZE), send_proc
00741 
00742     ! get the message
00743     call mpi_recv(buffer, COAG_DIST_MAX_BUFFER_SIZE, MPI_CHARACTER, &
00744          MPI_ANY_SOURCE, COAG_DIST_TAG_RETURN_UNREQ_PARTICLE, &
00745          MPI_COMM_WORLD, status, ierr)
00746     call pmc_mpi_check_ierr(ierr)
00747     call assert(496247788, status(MPI_TAG) &
00748          == COAG_DIST_TAG_RETURN_UNREQ_PARTICLE)
00749     call mpi_get_count(status, MPI_CHARACTER, buffer_size, ierr)
00750     call pmc_mpi_check_ierr(ierr)
00751     call assert(590644042, buffer_size <= COAG_DIST_MAX_BUFFER_SIZE)
00752     sent_proc = status(MPI_SOURCE)
00753 
00754     ! unpack it
00755     position = 0
00756     call aero_particle_allocate(aero_particle)
00757     call pmc_mpi_unpack_aero_particle(buffer, position, aero_particle)
00758     call assert(833588594, position == buffer_size)
00759 
00760     ! put it back
00761     call aero_state_add_particle(aero_state, aero_particle, &
00762          allow_resort=.false.)
00763     call aero_particle_deallocate(aero_particle)
00764 #endif
00765 
00766   end subroutine recv_return_unreq_particle
00767 
00768 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00769 
00770   !> Send a message saying that this process is finished with its
00771   !> coagulation.
00772   subroutine send_done(dest_proc)
00773 
00774     !> Process to send the message to.
00775     integer, intent(in) :: dest_proc
00776 
00777 #ifdef PMC_USE_MPI
00778     character :: buffer(COAG_DIST_MAX_BUFFER_SIZE)
00779     integer :: buffer_size, ierr
00780 
00781     buffer_size = 0
00782     call mpi_bsend(buffer, buffer_size, MPI_CHARACTER, dest_proc, &
00783          COAG_DIST_TAG_DONE, MPI_COMM_WORLD, ierr)
00784     call pmc_mpi_check_ierr(ierr)
00785 #endif
00786 
00787   end subroutine send_done
00788 
00789 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00790 
00791   !> Receive a done message.
00792   subroutine recv_done(procs_done)
00793 
00794     !> Which processes are finished with coagulation.
00795     logical, intent(inout) :: procs_done(:)
00796     
00797 #ifdef PMC_USE_MPI
00798     integer :: buffer_size, sent_proc, ierr
00799     character :: buffer(COAG_DIST_MAX_BUFFER_SIZE)
00800     integer :: status(MPI_STATUS_SIZE)
00801 
00802     ! get the message
00803     call mpi_recv(buffer, COAG_DIST_MAX_BUFFER_SIZE, MPI_CHARACTER, &
00804          MPI_ANY_SOURCE, COAG_DIST_TAG_DONE, MPI_COMM_WORLD, &
00805          status, ierr)
00806     call pmc_mpi_check_ierr(ierr)
00807     call assert(348737947, status(MPI_TAG) &
00808          == COAG_DIST_TAG_DONE)
00809     call mpi_get_count(status, MPI_CHARACTER, buffer_size, ierr)
00810     call pmc_mpi_check_ierr(ierr)
00811     call assert(214904056, buffer_size == 0)
00812     sent_proc = status(MPI_SOURCE)
00813 
00814     ! process it
00815     procs_done(sent_proc + 1) = .true.
00816 #endif
00817 
00818   end subroutine recv_done
00819 
00820 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00821 
00822   !> generate the number of samples to do per bin pair.
00823   subroutine generate_n_samps(n_parts, del_t, bin_grid, aero_weight_array, &
00824        k_max, n_samps, accept_factors)
00825     
00826     !> Number of particles per bin on all processes.
00827     integer, intent(in) :: n_parts(:,:)
00828     !> Timestep.
00829     real(kind=dp), intent(in) :: del_t
00830     !> Bin grid.
00831     type(bin_grid_t), intent(in) :: bin_grid
00832     !> Weighting function array.
00833     type(aero_weight_t), intent(in) :: aero_weight_array(:)
00834     !> Maximum kernel.
00835     real(kind=dp), intent(in) :: k_max(:,:)
00836     !> Number of samples to do per bin pair.
00837     integer, intent(out) :: n_samps(:,:)
00838     !> Accept scale factors per bin pair (1).
00839     real(kind=dp), intent(out) :: accept_factors(:,:)
00840 
00841     integer :: i_bin, j_bin, rank, n_bin
00842     real(kind=dp) :: n_samp_mean
00843 
00844     n_bin = size(k_max, 1)
00845     rank = pmc_mpi_rank()
00846     n_samps = 0
00847     do i_bin = 1,n_bin
00848        if (n_parts(i_bin, rank + 1) == 0) &
00849             cycle
00850        do j_bin = i_bin,n_bin
00851           call compute_n_samp(n_parts(i_bin, rank + 1), &
00852                sum(n_parts(j_bin, :)), (i_bin == j_bin), &
00853                k_max(i_bin, j_bin), del_t, n_samp_mean, &
00854                n_samps(i_bin, j_bin), accept_factors(i_bin, j_bin))
00855        end do
00856     end do
00857 
00858   end subroutine generate_n_samps
00859 
00860 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00861 
00862   subroutine coagulate_dist(aero_data, aero_state, aero_particle_1, &
00863        aero_particle_2, remote_proc, aero_weight_total, comp_vols, &
00864        remove_1, remove_2)
00865 
00866     !> Aerosol data.
00867     type(aero_data_t), intent(in) :: aero_data
00868     !> Aerosol state.
00869     type(aero_state_t), intent(inout) :: aero_state
00870     !> First particle to coagulate.
00871     type(aero_particle_t), intent(in) :: aero_particle_1
00872     !> Second particle to coagulate.
00873     type(aero_particle_t), intent(in) :: aero_particle_2
00874     !> Remote process that the particle came from.
00875     integer, intent(in) :: remote_proc
00876     !> Total weight across all processes.
00877     type(aero_weight_t), intent(in) :: aero_weight_total(:)
00878     !> Computational volumes on all processes (m^3).
00879     real(kind=dp), intent(in) :: comp_vols(:,:)
00880     !> Whether to remove aero_particle_1 after the coagulation.
00881     logical, intent(out) :: remove_1
00882     !> Whether to remove aero_particle_2 after the coagulation.
00883     logical, intent(out) :: remove_2
00884     
00885     type(aero_particle_t) :: aero_particle_new
00886     integer :: new_proc, new_group
00887     type(aero_info_t) :: aero_info_1, aero_info_2
00888     logical :: create_new, id_1_lost, id_2_lost
00889 
00890     call aero_particle_allocate(aero_particle_new)
00891     call aero_info_allocate(aero_info_1)
00892     call aero_info_allocate(aero_info_2)
00893 
00894     call coagulate_weighting(aero_particle_1, aero_particle_2, &
00895          aero_particle_new, aero_data, aero_state%aero_weight, &
00896          remove_1, remove_2, create_new, id_1_lost, id_2_lost, &
00897          aero_info_1, aero_info_2)
00898 
00899     if (id_1_lost) then
00900        call aero_info_array_add_aero_info(aero_state%aero_info_array, &
00901             aero_info_1)
00902     end if
00903     if (id_2_lost) then
00904        call aero_info_array_add_aero_info(aero_state%aero_info_array, &
00905             aero_info_2)
00906     end if
00907 
00908     ! add new particle
00909     if (create_new) then
00910        new_group = aero_weight_array_rand_group(aero_weight_total, &
00911             aero_particle_radius(aero_particle_new))
00912        aero_particle_new%weight_group = new_group
00913        new_proc = sample_cts_pdf(comp_vols(new_group, :)) - 1
00914        call send_return_unreq_particle(aero_particle_new, new_proc)
00915     end if
00916 
00917     call aero_particle_deallocate(aero_particle_new)
00918     call aero_info_deallocate(aero_info_1)
00919     call aero_info_deallocate(aero_info_2)
00920 
00921   end subroutine coagulate_dist
00922 
00923 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
00924   
00925 end module pmc_coagulation_dist