Ignore:
Timestamp:
07/25/19 19:53:31 (5 years ago)
Author:
adurocher
Message:

trunk : Implemented req_z1_scal in new transfert_mpi

File:
1 edited

Legend:

Unmodified
Added
Removed
  • codes/icosagcm/trunk/src/parallel/transfert_requests.F90

    r963 r965  
    22module transfert_request_mod 
    33  use abort_mod 
    4   use field_mod, only : t_field, field_T, field_U 
     4  use field_mod, only : t_field, field_T, field_U, field_Z 
    55  use domain_mod, only : ndomain, ndomain_glo, domain, domain_glo, domloc_glo_ind, domglo_rank, domglo_loc_ind 
    66  implicit none 
     
    99  ! Points for t_request : list of points coordinates in a domain 
    1010  type t_points 
    11     integer, allocatable :: i(:), j(:), edge(:) 
     11    integer, allocatable :: i(:), j(:) 
     12    integer, allocatable :: elt(:) ! Element : cell edge or cell point 
    1213    integer :: npoints 
    1314  end type 
     
    2627  type(t_request),save,pointer :: req_e1_scal(:) ! Halos for U fields 
    2728  type(t_request),save,pointer :: req_e1_vect(:) ! Halos for U fields (with sign change) 
    28   type(t_request),save,pointer :: req_z1_scal(:) => null() ! DO NOT USE : not initialized 
     29  type(t_request),save,pointer :: req_z1_scal(:) ! Halos for Z fields 
    2930 
    3031  type(t_request),save,pointer :: req_i0(:) ! Duplicated cells for T fields 
     
    4243  subroutine init_all_requests 
    4344    use dimensions, only : swap_dimensions, ii_begin, ii_end, jj_begin, jj_end 
    44     use metric, only : right, rdown, ldown, left, lup, rup 
     45    use metric, only : right, rdown, ldown, left, lup, rup, vdown, vlup, vrdown, vrup, vup 
    4546    integer :: ind, i, j 
    4647 
     
    131132    req_e1_vect = req_e1_scal 
    132133    req_e1_vect%vector = .true. 
     134 
     135    call request_init(req_z1_scal, field_z) 
     136    do ind=1,ndomain 
     137      call swap_dimensions(ind) 
     138      do i=ii_begin,ii_end 
     139        call request_add_point(req_z1_scal,ind,i,jj_begin-1,vrup) 
     140        call request_add_point(req_z1_scal,ind,i+1,jj_begin-1,vup) 
     141      enddo 
     142      do j=jj_begin,jj_end 
     143        call request_add_point(req_z1_scal,ind,ii_end+1,j,vlup) 
     144      enddo 
     145      do j=jj_begin,jj_end 
     146        call request_add_point(req_z1_scal,ind,ii_end+1,j-1,vup) 
     147      enddo 
     148      do i=ii_begin,ii_end 
     149        call request_add_point(req_z1_scal,ind,i,jj_end+1,vdown) 
     150        call request_add_point(req_z1_scal,ind,i-1,jj_end+1,vrdown) 
     151      enddo 
     152      do j=jj_begin,jj_end 
     153        call request_add_point(req_z1_scal,ind,ii_begin-1,j,vrup) 
     154      enddo 
     155      do j=jj_begin,jj_end 
     156        call request_add_point(req_z1_scal,ind,ii_begin-1,j,vrdown) 
     157      enddo 
     158    enddo 
     159    call request_exchange(req_z1_scal) 
    133160  contains 
    134161 
    135162    ! ---- Points ---- 
    136163    ! Initialize (allocate) points 
    137     subroutine points_init(points, has_edges) 
     164    subroutine points_init(points, has_elt) 
    138165      type(t_points):: points 
    139       logical :: has_edges 
     166      logical :: has_elt 
    140167      integer, parameter :: INITIAL_ALLOC_SIZE = 10 
    141168 
    142169      allocate(points%i(INITIAL_ALLOC_SIZE)) 
    143170      allocate(points%j(INITIAL_ALLOC_SIZE)) 
    144       if(has_edges) allocate(points%edge(INITIAL_ALLOC_SIZE)) 
     171      if(has_elt) allocate(points%elt(INITIAL_ALLOC_SIZE)) 
    145172      points%npoints = 0 
    146173    end subroutine 
    147174 
    148     ! Append point (i, j [,edge]) to point list 
    149     subroutine points_append(points, i, j, edge) 
     175    ! Append point (i, j [,elt]) to point list 
     176    subroutine points_append(points, i, j, elt) 
    150177      type(t_points), intent(inout):: points 
    151178      integer, intent(in) :: i, j 
    152       integer, intent(in), optional :: edge 
     179      integer, intent(in), optional :: elt 
    153180      integer :: begin_size 
    154181 
     
    157184      begin_size=points%npoints 
    158185      call array_append( points%j, begin_size, j ) 
    159       if(present(edge)) then 
     186      if(present(elt)) then 
    160187        begin_size=points%npoints 
    161         call array_append( points%edge, begin_size, edge ) 
     188        call array_append( points%elt, begin_size, elt ) 
    162189      end if 
    163190 
    164191      points%npoints = points%npoints + 1 
    165192    end subroutine 
    166      
     193 
    167194    ! Add element to array, and reallocate if necessary 
    168195    subroutine array_append( a, a_size, elt ) 
     
    172199      integer, allocatable :: a_tmp(:) 
    173200      integer, parameter :: GROW_FACTOR = 2 
    174        
     201 
    175202      if( size( a ) <= a_size ) then 
    176203        allocate( a_tmp ( a_size * GROW_FACTOR ) ) 
     
    192219 
    193220      displs = (points%i-1) + (points%j-1)*iim 
    194       if(allocated(points%edge)) displs = displs + (points%edge-1)*iim*jjm 
     221      if(allocated(points%elt)) displs = displs + (points%elt-1)*iim*jjm 
    195222 
    196223      last=0 
     
    210237        points%j(k) = mod((displs(k)-(points%i(k)-1))/iim, jjm)+1 
    211238      end do 
    212       if(allocated(points%edge)) then 
    213         deallocate(points%edge); allocate(points%edge(last)) 
    214         points%edge(:) = displs(1:last)/(iim*jjm) + 1 
     239      if(allocated(points%elt)) then 
     240        deallocate(points%elt); allocate(points%elt(last)) 
     241        points%elt(:) = displs(1:last)/(iim*jjm) + 1 
    215242      end if 
    216243 
     
    231258        allocate(request(ind)%points_HtoB(ndomain_glo)) 
    232259        do remote_ind_glo = 1, ndomain_glo 
    233           call points_init(request(ind)%points_BtoH(remote_ind_glo), field_type == field_U) 
     260          call points_init(request(ind)%points_BtoH(remote_ind_glo), field_type /= field_T) 
    234261          !call points_init(request(ind)%points_HtoB(remote_ind_glo)) ! Initialized before MPI communication in request_exchange 
    235262        end do 
     
    238265 
    239266    ! Append local point (ind_glo,i,j) to the point list of domain ind, where ind_glo is the domain owning the point 
    240     subroutine request_add_point(req, ind, i, j, edge) 
     267    subroutine request_add_point(req, ind, i, j, elt) 
    241268      type(t_request), intent(inout) :: req(:) 
    242269      integer :: ind, i, j 
    243       integer, optional :: edge 
    244  
    245       if( present(edge) ) then ! U field 
    246         if(domain(ind)%edge_assign_domain(edge-1,i,j) /= domloc_glo_ind(ind) ) & ! Add only if remote domain != local domain 
    247         call points_append(req(ind)%points_BtoH(domain(ind)%edge_assign_domain(edge-1,i,j)), i, j, edge ) 
     270      integer, optional :: elt 
     271 
     272      if( req(ind)%field_type == field_U ) then 
     273        if(domain(ind)%edge_assign_domain(elt-1,i,j) /= domloc_glo_ind(ind) ) & ! Add only if remote domain != local domain 
     274        call points_append(req(ind)%points_BtoH(domain(ind)%edge_assign_domain(elt-1,i,j)), i, j, elt ) 
     275      else if( req(ind)%field_type == field_Z ) then 
     276        if(domain(ind)%vertex_assign_domain(elt-1,i,j) /= domloc_glo_ind(ind) ) & ! Add only if remote domain != local domain 
     277        call points_append(req(ind)%points_BtoH(domain(ind)%vertex_assign_domain(elt-1,i,j)), i, j, elt ) 
    248278      else ! T field 
    249279        if(domain(ind)%assign_domain(i,j) /= domloc_glo_ind(ind) ) & 
     
    262292 
    263293      integer :: ind_loc, remote_ind_glo, k 
    264       integer :: i, j, edge 
     294      integer :: i, j, elt 
    265295      integer :: reqs(ndomain, ndomain_glo, 2), reqs2(ndomain, ndomain_glo, 2), reqs3(ndomain, ndomain_glo, 3), ierr 
    266296 
     
    271301          points => req(ind_loc)%points_BtoH(remote_ind_glo) 
    272302          call points_sort(points) 
    273           call points_init( points_send(ind_loc,remote_ind_glo), req(ind_loc)%field_type == field_U ) 
     303          call points_init( points_send(ind_loc,remote_ind_glo), req(ind_loc)%field_type /= field_T ) 
    274304          do k = 1, points%npoints 
    275305            i = points%i(k) 
    276306            j = points%j(k) 
    277             if(req(ind_loc)%field_type == field_T) then 
     307            if( req(ind_loc)%field_type == field_U ) then 
     308              elt = points%elt(k) - 1 ! WARNING : domain%edge_assign_* use index in 0:5 instead of 1:6 everywhere else 
     309              call points_append(points_send(ind_loc,remote_ind_glo), domain(ind_loc)%edge_assign_i(elt,i,j), & 
     310                domain(ind_loc)%edge_assign_j(elt,i,j), & 
     311                domain(ind_loc)%edge_assign_pos(elt,i,j)+1) 
     312            else if( req(ind_loc)%field_type == field_Z ) then 
     313              elt = points%elt(k) - 1 ! WARNING : domain%vertex_assign_* use index in 0:5 instead of 1:6 everywhere else 
     314              call points_append(points_send(ind_loc,remote_ind_glo), domain(ind_loc)%vertex_assign_i(elt,i,j), & 
     315                domain(ind_loc)%vertex_assign_j(elt,i,j), & 
     316                domain(ind_loc)%vertex_assign_pos(elt,i,j)+1) 
     317            else ! T field 
    278318              call points_append(points_send(ind_loc,remote_ind_glo), domain(ind_loc)%assign_i(i,j), & 
    279319                domain(ind_loc)%assign_j(i,j)) 
    280             else 
    281               edge = points%edge(k) - 1 ! WARNING : domain%edge_assign_* use index in 0:5 instead of 1:6 everywhere else 
    282               call points_append(points_send(ind_loc,remote_ind_glo), domain(ind_loc)%edge_assign_i(edge,i,j), & 
    283                 domain(ind_loc)%edge_assign_j(edge,i,j), & 
    284                 domain(ind_loc)%edge_assign_pos(edge,i,j)+1) 
    285             end if 
     320            endif 
    286321          end do 
    287322        end do 
     
    317352            domloc_glo_ind(ind_loc)+ndomain_glo*remote_ind_glo, comm_icosa, reqs2(ind_loc,remote_ind_glo,2), ierr) 
    318353 
    319           if(req(1)%field_type == field_U) then 
    320             call MPI_Isend( points_send(ind_loc,remote_ind_glo)%edge, points_send(ind_loc,remote_ind_glo)%npoints, MPI_INT, & 
     354          if(req(1)%field_type /= field_T) then 
     355            call MPI_Isend( points_send(ind_loc,remote_ind_glo)%elt, points_send(ind_loc,remote_ind_glo)%npoints, MPI_INT, & 
    321356              domglo_rank(remote_ind_glo), remote_ind_glo+ndomain_glo*domloc_glo_ind(ind_loc), & 
    322357              comm_icosa, reqs3(ind_loc,remote_ind_glo,1), ierr  ) 
    323             allocate( points%edge( points%npoints ) ) 
    324             call MPI_Irecv( points%edge, points%npoints, MPI_INT, & 
     358            allocate( points%elt( points%npoints ) ) 
     359            call MPI_Irecv( points%elt, points%npoints, MPI_INT, & 
    325360              domglo_rank(remote_ind_glo), domloc_glo_ind(ind_loc)+ndomain_glo*remote_ind_glo, & 
    326361              comm_icosa, reqs3(ind_loc,remote_ind_glo,2), ierr ) 
     
    330365      call waitall(reqs,ndomain*ndomain_glo*2) 
    331366      call waitall(reqs2,ndomain*ndomain_glo*2) 
    332       if(req(1)%field_type == field_U) call waitall(reqs3,ndomain*ndomain_glo*2) 
     367      if(req(1)%field_type /= field_T) call waitall(reqs3,ndomain*ndomain_glo*2) 
    333368    end subroutine 
    334369 
Note: See TracChangeset for help on using the changeset viewer.