Changeset 1002


Ignore:
Timestamp:
01/10/20 12:20:52 (4 years ago)
Author:
adurocher
Message:

transfert_mpi : Moved buffer copies in separate subroutines

new routines to transfer data between Halos to Buffers are 'copy_HtoB', 'copy_HtoH', 'copy_BtoH'

File:
1 edited

Legend:

Unmodified
Added
Removed
  • codes/icosagcm/trunk/src/parallel/transfert_mpi.f90

    r1001 r1002  
    414414  end subroutine 
    415415 
    416   subroutine send_message(field, message) 
    417     use mpi_mod 
     416  ! Halo to Buffer : copy outbound message to MPI buffers 
     417  subroutine copy_HtoB(message) 
    418418    use domain_mod, only : assigned_domain 
    419419    use omp_para, only : distrib_level 
    420     type(t_field),pointer :: field(:) 
    421     type(t_message), target :: message 
    422     integer :: i, k, d3, d4, local_displ 
    423     integer :: ierr, d3_begin, d3_end, dim4, dim3 
    424  
    425     call enter_profile(profile_mpi) 
    426  
    427     ! Needed because rval4d is set in init_message 
    428     if( .not. associated( message%field, field ) ) & 
    429       call dynamico_abort("send_message must be called with the same field used in init_message") 
    430  
    431     !Prepare 'message' for on-device copies if field is on device 
    432     !$omp master 
    433     if( field(1)%ondevice .and. .not. message%ondevice ) call message_create_ondevice(message) 
    434     if( field(1)%ondevice .neqv. message%ondevice ) call dynamico_abort("send_message : internal device/host memory synchronization error") 
    435     ! Check if previous message has been waited 
    436     if(message%send_seq /= message%wait_seq) & 
    437       call dynamico_abort("No matching wait_message before new send_message") 
    438     message%send_seq = message%send_seq + 1 
    439     !$omp end master 
    440  
    441     call enter_profile(profile_mpi_barrier) 
    442     !$omp barrier 
    443     call exit_profile(profile_mpi_barrier) 
     420    type(t_message), intent(inout) :: message 
     421    integer :: dim3, dim4, d3_begin, d3_end 
     422    integer :: k, d3, d4, i 
     423    integer :: local_displ 
    444424 
    445425    dim4 = size(message%field(1)%rval4d, 3) 
    446426    dim3 = size(message%field(1)%rval4d, 2) 
    447427    CALL distrib_level( 1, dim3, d3_begin, d3_end ) 
    448  
    449     call enter_profile(profile_mpi_copies) 
    450     !$acc data present(message) async if(message%ondevice) 
    451  
    452     ! Halo to Buffer : copy outbound message to MPI buffers 
     428     
    453429    !$acc parallel loop default(none) present(message, assigned_domain) async if(message%ondevice) 
    454430    do i = 1, size( message%message_out ) 
     
    460436            !$acc loop 
    461437             do k = 1, message%message_out(i)%npoints 
    462               !print *, "buff", message%message_out(i)%remote_rank, k+local_displ, "<- field", message%message_out(i)%ind_loc, message%message_out(i)%displs(k), d3, d4, message%field(message%message_out(i)%ind_loc)%rval4d(message%message_out(i)%displs(k),d3, d4) 
    463438              message%mpi_buffer_out(message%message_out(i)%remote_rank)%buff(k+local_displ) = message%field(message%message_out(i)%ind_loc)%rval4d(message%message_out(i)%displs(k),d3, d4) 
    464439            end do 
     
    467442      endif 
    468443    end do 
    469  
    470     ! Halo to Halo : copy local messages from source field to destination field 
     444  end subroutine 
     445 
     446  ! Halo to Halo : copy local messages from source field to destination field 
     447  subroutine copy_HtoH(message) 
     448    use domain_mod, only : assigned_domain 
     449    use omp_para, only : distrib_level 
     450    type(t_message), intent(inout) :: message 
     451    integer :: dim3, dim4, d3_begin, d3_end 
     452    integer :: k, d3, d4, i 
     453 
     454    dim4 = size(message%field(1)%rval4d, 3) 
     455    dim3 = size(message%field(1)%rval4d, 2) 
     456    CALL distrib_level( 1, dim3, d3_begin, d3_end ) 
     457 
    471458    !$acc parallel loop default(none) present(message, assigned_domain) async if(message%ondevice) 
    472459    do i = 1, size( message%message_local ) 
     
    485472      endif 
    486473    end do 
    487  
    488     !$acc end data 
    489     call exit_profile(profile_mpi_copies) 
    490  
    491  
    492     call enter_profile(profile_mpi_barrier) 
    493     !$acc wait 
    494     !$omp barrier 
    495     call exit_profile(profile_mpi_barrier) 
    496  
    497     !$omp master 
    498     call MPI_Startall( size(message%mpi_requests_out), message%mpi_requests_out, ierr ) 
    499     call MPI_Startall( size(message%mpi_requests_in), message%mpi_requests_in, ierr ) 
    500     !$omp end master 
    501  
    502     call exit_profile(profile_mpi) 
    503   end subroutine 
    504  
    505   subroutine test_message(message) 
    506     use mpi_mod 
    507     type(t_message) :: message 
    508     integer :: ierr 
    509     logical :: completed 
    510  
    511     !$omp master 
    512     call MPI_Testall( size(message%mpi_requests_out), message%mpi_requests_out, completed, MPI_STATUSES_IGNORE, ierr ) 
    513     call MPI_Testall( size(message%mpi_requests_in), message%mpi_requests_in, completed, MPI_STATUSES_IGNORE, ierr ) 
    514     !$omp end master 
    515   end subroutine 
    516  
    517   subroutine wait_message(message) 
    518     use mpi_mod 
     474  end subroutine 
     475 
     476  ! Buffer to Halo : copy inbound message to field 
     477  subroutine copy_BtoH(message) 
    519478    use domain_mod, only : assigned_domain 
    520479    use omp_para, only : distrib_level 
    521     type(t_message), target :: message 
    522     integer :: d3_begin, d3_end, dim4, dim3, dim12 
    523     integer :: i, ierr, k, d3, d4, local_displ 
    524  
    525     ! Check if message has been sent and not recieved yet 
    526     ! note : barrier needed between this and send_seq increment, and this and wait_seq increment 
    527     ! note : watch out for integer overflow a = b+1 doesn't imply a>b 
    528     if(message%send_seq /= message%wait_seq+1) then 
    529       print*, "WARNING : wait_message called multiple times for one send_message, skipping" 
    530       return ! Don't recieve message if already recieved 
    531     end if 
    532  
    533     call enter_profile(profile_mpi) 
     480    type(t_message), intent(inout) :: message 
     481    integer :: dim3, dim4, d3_begin, d3_end 
     482    integer :: k, d3, d4, i 
     483    integer :: local_displ 
    534484 
    535485    dim4 = size(message%field(1)%rval4d, 3) 
    536486    dim3 = size(message%field(1)%rval4d, 2) 
    537     dim12 = size(message%field(1)%rval4d, 1) 
    538487    CALL distrib_level( 1, dim3, d3_begin, d3_end ) 
    539  
    540     call enter_profile(profile_mpi_waitall) 
    541     !$omp master 
    542     call MPI_Waitall( size(message%mpi_requests_out), message%mpi_requests_out, MPI_STATUSES_IGNORE, ierr ) 
    543     call MPI_Waitall( size(message%mpi_requests_in), message%mpi_requests_in, MPI_STATUSES_IGNORE, ierr ) 
    544     !$omp end master 
    545     call exit_profile(profile_mpi_waitall) 
    546  
    547     call enter_profile(profile_mpi_barrier) 
    548     !$omp barrier 
    549     call exit_profile(profile_mpi_barrier) 
    550  
    551     call enter_profile(profile_mpi_copies) 
    552     !$acc data present(message) async if(message%ondevice) 
    553  
    554     dim4 = size(message%field(1)%rval4d, 3) 
    555     ! Buffer to Halo : copy inbound message to field 
     488     
    556489    !$acc parallel loop default(none) present(message, assigned_domain) async if(message%ondevice) 
    557490    do i = 1, size( message%message_in ) 
     
    563496            !$acc loop 
    564497            do k = 1, message%message_in(i)%npoints  
    565               !print *, "field", message%message_in(i)%ind_loc, message%message_in(i)%displs(k), d3, d4, "-> buff", message%message_in(i)%remote_rank, k+local_displ, message%mpi_buffer_in(message%message_in(i)%remote_rank)%buff(k+local_displ) 
    566498              message%field(message%message_in(i)%ind_loc)%rval4d(message%message_in(i)%displs(k),d3,d4) & 
    567499                  = message%message_in(i)%sign(k)*message%mpi_buffer_in(message%message_in(i)%remote_rank)%buff(k+local_displ) 
     
    571503      endif 
    572504    end do 
    573  
    574     !$acc end data 
     505  end subroutine 
     506 
     507  subroutine send_message(field, message) 
     508    use mpi_mod 
     509    type(t_field),pointer :: field(:) 
     510    type(t_message), target :: message 
     511    integer :: ierr 
     512 
     513    call enter_profile(profile_mpi) 
     514 
     515    ! Needed because rval4d is set in init_message 
     516    if( .not. associated( message%field, field ) ) & 
     517      call dynamico_abort("send_message must be called with the same field used in init_message") 
     518 
     519    !Prepare 'message' for on-device copies if field is on device 
     520    !$omp master 
     521    if( field(1)%ondevice .and. .not. message%ondevice ) call message_create_ondevice(message) 
     522    if( field(1)%ondevice .neqv. message%ondevice ) call dynamico_abort("send_message : internal device/host memory synchronization error") 
     523    ! Check if previous message has been waited 
     524    if(message%send_seq /= message%wait_seq) & 
     525      call dynamico_abort("No matching wait_message before new send_message") 
     526    message%send_seq = message%send_seq + 1 
     527    !$omp end master 
     528 
     529    call enter_profile(profile_mpi_barrier) 
     530    !$omp barrier 
     531    call exit_profile(profile_mpi_barrier) 
     532 
     533    call enter_profile(profile_mpi_copies) 
     534    call copy_HtoB(message) 
     535    call copy_HtoH(message) 
    575536    call exit_profile(profile_mpi_copies) 
    576537 
     538    call enter_profile(profile_mpi_barrier) 
     539    !$acc wait 
     540    !$omp barrier 
     541    call exit_profile(profile_mpi_barrier) 
     542 
     543    !$omp master 
     544    call MPI_Startall( size(message%mpi_requests_out), message%mpi_requests_out, ierr ) 
     545    call MPI_Startall( size(message%mpi_requests_in), message%mpi_requests_in, ierr ) 
     546    !$omp end master 
     547 
     548    call exit_profile(profile_mpi) 
     549  end subroutine 
     550 
     551  subroutine test_message(message) 
     552    use mpi_mod 
     553    type(t_message) :: message 
     554    integer :: ierr 
     555    logical :: completed 
     556 
     557    !$omp master 
     558    call MPI_Testall( size(message%mpi_requests_out), message%mpi_requests_out, completed, MPI_STATUSES_IGNORE, ierr ) 
     559    call MPI_Testall( size(message%mpi_requests_in), message%mpi_requests_in, completed, MPI_STATUSES_IGNORE, ierr ) 
     560    !$omp end master 
     561  end subroutine 
     562 
     563  subroutine wait_message(message) 
     564    use mpi_mod 
     565    type(t_message), target :: message 
     566    integer :: ierr 
     567 
     568    ! Check if message has been sent and not recieved yet 
     569    ! note : barrier needed between this and send_seq increment, and this and wait_seq increment 
     570    ! note : watch out for integer overflow a = b+1 doesn't imply a>b 
     571    if(message%send_seq /= message%wait_seq+1) then 
     572      print*, "WARNING : wait_message called multiple times for one send_message, skipping" 
     573      return ! Don't recieve message if already recieved 
     574    end if 
     575 
     576    call enter_profile(profile_mpi) 
     577 
     578    call enter_profile(profile_mpi_waitall) 
     579    !$omp master 
     580    call MPI_Waitall( size(message%mpi_requests_out), message%mpi_requests_out, MPI_STATUSES_IGNORE, ierr ) 
     581    call MPI_Waitall( size(message%mpi_requests_in), message%mpi_requests_in, MPI_STATUSES_IGNORE, ierr ) 
     582    !$omp end master 
     583    call exit_profile(profile_mpi_waitall) 
     584 
     585    call enter_profile(profile_mpi_barrier) 
     586    !$omp barrier 
     587    call exit_profile(profile_mpi_barrier) 
     588 
     589    call enter_profile(profile_mpi_copies)  
     590    call copy_BtoH(message) 
     591    call exit_profile(profile_mpi_copies) 
     592 
    577593    !$omp master 
    578594    message%wait_seq = message%wait_seq + 1 
     
    582598    !$omp barrier 
    583599    call exit_profile(profile_mpi_barrier) 
     600 
    584601    call exit_profile(profile_mpi) 
    585602  end subroutine 
Note: See TracChangeset for help on using the changeset viewer.