Changeset 965


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

trunk : Implemented req_z1_scal in new transfert_mpi

Location:
codes/icosagcm/trunk/src/parallel
Files:
3 edited

Legend:

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

    r964 r965  
    1515                                wait_message, & 
    1616                                test_message 
    17    
     17 
    1818#else 
    1919  use transfert_mpi_legacy_mod, only :  t_message, t_request, & 
     
    4949          req_i1, req_e1_scal, req_e1_vect, & 
    5050          req_i0, req_e0_scal, req_e0_vect, & 
    51 #if !defined(CPP_USING_MPI_NEW) 
    5251          req_z1_scal, & 
    53 #endif 
    5452          init_transfert, & 
    5553          init_message, & 
     
    9997  IMPLICIT NONE 
    10098    CHARACTER(LEN=*),INTENT(INOUT) :: Var 
    101     
    102 !$OMP MASTER 
    103     CALL bcast_mpi(Var) 
    104 !$OMP END MASTER 
    105     CALL bcast_omp(Var) 
    106      
     99 
     100!$OMP MASTER 
     101    CALL bcast_mpi(Var) 
     102!$OMP END MASTER 
     103    CALL bcast_omp(Var) 
     104 
    107105  END SUBROUTINE bcast_c 
    108106 
    109107!! -- Les entiers -- !! 
    110    
     108 
    111109  SUBROUTINE bcast_i(var) 
    112110  IMPLICIT NONE 
     
    116114!$OMP END MASTER 
    117115    CALL bcast_omp(Var) 
    118      
     116 
    119117  END SUBROUTINE bcast_i 
    120118 
     
    122120  IMPLICIT NONE 
    123121    INTEGER,INTENT(INOUT) :: Var(:) 
    124     
    125 !$OMP MASTER 
    126     CALL bcast_mpi(Var) 
    127 !$OMP END MASTER 
    128     CALL bcast_omp(Var) 
    129      
     122 
     123!$OMP MASTER 
     124    CALL bcast_mpi(Var) 
     125!$OMP END MASTER 
     126    CALL bcast_omp(Var) 
     127 
    130128  END SUBROUTINE bcast_i1 
    131129 
     
    134132  IMPLICIT NONE 
    135133    INTEGER,INTENT(INOUT) :: Var(:,:) 
    136     
    137 !$OMP MASTER 
    138     CALL bcast_mpi(Var) 
    139 !$OMP END MASTER 
    140     CALL bcast_omp(Var) 
    141      
     134 
     135!$OMP MASTER 
     136    CALL bcast_mpi(Var) 
     137!$OMP END MASTER 
     138    CALL bcast_omp(Var) 
     139 
    142140  END SUBROUTINE bcast_i2 
    143141 
     
    146144  IMPLICIT NONE 
    147145    INTEGER,INTENT(INOUT) :: Var(:,:,:) 
    148     
    149 !$OMP MASTER 
    150     CALL bcast_mpi(Var) 
    151 !$OMP END MASTER 
    152     CALL bcast_omp(Var) 
    153      
     146 
     147!$OMP MASTER 
     148    CALL bcast_mpi(Var) 
     149!$OMP END MASTER 
     150    CALL bcast_omp(Var) 
     151 
    154152  END SUBROUTINE bcast_i3 
    155153 
     
    158156  IMPLICIT NONE 
    159157    INTEGER,INTENT(INOUT) :: Var(:,:,:,:) 
    160     
    161 !$OMP MASTER 
    162     CALL bcast_mpi(Var) 
    163 !$OMP END MASTER 
    164     CALL bcast_omp(Var) 
    165      
     158 
     159!$OMP MASTER 
     160    CALL bcast_mpi(Var) 
     161!$OMP END MASTER 
     162    CALL bcast_omp(Var) 
     163 
    166164  END SUBROUTINE bcast_i4 
    167165 
    168   
     166 
    169167!! -- Les reels -- !! 
    170    
     168 
    171169  SUBROUTINE bcast_r(var) 
    172170  IMPLICIT NONE 
     
    177175!$OMP END MASTER 
    178176    CALL bcast_omp(Var) 
    179      
     177 
    180178  END SUBROUTINE bcast_r 
    181179 
     
    183181  IMPLICIT NONE 
    184182    REAL,INTENT(INOUT) :: Var(:) 
    185     
    186 !$OMP MASTER 
    187     CALL bcast_mpi(Var) 
    188 !$OMP END MASTER 
    189     CALL bcast_omp(Var) 
    190      
     183 
     184!$OMP MASTER 
     185    CALL bcast_mpi(Var) 
     186!$OMP END MASTER 
     187    CALL bcast_omp(Var) 
     188 
    191189  END SUBROUTINE bcast_r1 
    192190 
     
    195193  IMPLICIT NONE 
    196194    REAL,INTENT(INOUT) :: Var(:,:) 
    197     
    198 !$OMP MASTER 
    199     CALL bcast_mpi(Var) 
    200 !$OMP END MASTER 
    201     CALL bcast_omp(Var) 
    202      
     195 
     196!$OMP MASTER 
     197    CALL bcast_mpi(Var) 
     198!$OMP END MASTER 
     199    CALL bcast_omp(Var) 
     200 
    203201  END SUBROUTINE bcast_r2 
    204202 
     
    207205  IMPLICIT NONE 
    208206    REAL,INTENT(INOUT) :: Var(:,:,:) 
    209     
    210 !$OMP MASTER 
    211     CALL bcast_mpi(Var) 
    212 !$OMP END MASTER 
    213     CALL bcast_omp(Var) 
    214      
     207 
     208!$OMP MASTER 
     209    CALL bcast_mpi(Var) 
     210!$OMP END MASTER 
     211    CALL bcast_omp(Var) 
     212 
    215213  END SUBROUTINE bcast_r3 
    216214 
     
    219217  IMPLICIT NONE 
    220218    REAL,INTENT(INOUT) :: Var(:,:,:,:) 
    221     
    222 !$OMP MASTER 
    223     CALL bcast_mpi(Var) 
    224 !$OMP END MASTER 
    225     CALL bcast_omp(Var) 
    226      
    227   END SUBROUTINE bcast_r4  
     219 
     220!$OMP MASTER 
     221    CALL bcast_mpi(Var) 
     222!$OMP END MASTER 
     223    CALL bcast_omp(Var) 
     224 
     225  END SUBROUTINE bcast_r4 
    228226 
    229227 
    230228!! -- Les booleens -- !! 
    231    
     229 
    232230  SUBROUTINE bcast_l(var) 
    233231  IMPLICIT NONE 
     
    237235!$OMP END MASTER 
    238236    CALL bcast_omp(Var) 
    239      
     237 
    240238  END SUBROUTINE bcast_l 
    241239 
     
    243241  IMPLICIT NONE 
    244242    LOGICAL,INTENT(INOUT) :: Var(:) 
    245     
    246 !$OMP MASTER 
    247     CALL bcast_mpi(Var) 
    248 !$OMP END MASTER 
    249     CALL bcast_omp(Var) 
    250      
     243 
     244!$OMP MASTER 
     245    CALL bcast_mpi(Var) 
     246!$OMP END MASTER 
     247    CALL bcast_omp(Var) 
     248 
    251249  END SUBROUTINE bcast_l1 
    252250 
     
    255253  IMPLICIT NONE 
    256254    LOGICAL,INTENT(INOUT) :: Var(:,:) 
    257     
    258 !$OMP MASTER 
    259     CALL bcast_mpi(Var) 
    260 !$OMP END MASTER 
    261     CALL bcast_omp(Var) 
    262      
     255 
     256!$OMP MASTER 
     257    CALL bcast_mpi(Var) 
     258!$OMP END MASTER 
     259    CALL bcast_omp(Var) 
     260 
    263261  END SUBROUTINE bcast_l2 
    264262 
     
    267265  IMPLICIT NONE 
    268266    LOGICAL,INTENT(INOUT) :: Var(:,:,:) 
    269     
    270 !$OMP MASTER 
    271     CALL bcast_mpi(Var) 
    272 !$OMP END MASTER 
    273     CALL bcast_omp(Var) 
    274      
     267 
     268!$OMP MASTER 
     269    CALL bcast_mpi(Var) 
     270!$OMP END MASTER 
     271    CALL bcast_omp(Var) 
     272 
    275273  END SUBROUTINE bcast_l3 
    276274 
     
    279277  IMPLICIT NONE 
    280278    LOGICAL,INTENT(INOUT) :: Var(:,:,:,:) 
    281     
    282 !$OMP MASTER 
    283     CALL bcast_mpi(Var) 
    284 !$OMP END MASTER 
    285     CALL bcast_omp(Var) 
    286      
     279 
     280!$OMP MASTER 
     281    CALL bcast_mpi(Var) 
     282!$OMP END MASTER 
     283    CALL bcast_omp(Var) 
     284 
    287285  END SUBROUTINE bcast_l4 
    288286 
    289    
     287 
    290288END MODULE transfert_mod 
  • codes/icosagcm/trunk/src/parallel/transfert_mpi.f90

    r963 r965  
    55  use profiling_mod, only : enter_profile, exit_profile, register_id 
    66  use domain_mod, only : ndomain, ndomain_glo, domain, domain_glo, domloc_glo_ind, domglo_rank, domglo_loc_ind 
    7   use field_mod, only : t_field, field_T, field_U 
     7  use field_mod, only : t_field, field_T, field_U, field_Z 
    88  use transfert_request_mod 
    99  implicit none 
     
    103103    type(t_submessage), allocatable :: message_in_tmp(:), message_out_tmp(:) 
    104104    type(t_submessage), pointer :: submessage 
     105    integer :: field_type 
    105106 
    106107    !$omp barrier 
     
    110111    message%send_seq = 0 
    111112    message%wait_seq = 0 
     113 
     114    if( request(1)%field_type /= field(1)%field_type ) call dynamico_abort( "init_message : field_type/request mismatch" ) 
     115    field_type = request(1)%field_type 
     116 
    112117    ! Set field%rval4d pointer to always use 4d array 
    113118    do ind = 1, ndomain 
     
    134139          if(request(ind_loc)%points_HtoB(remote_ind_glo)%npoints > 0 ) then ! Add only non-empty messages 
    135140            ! Add local message ind_loc -> remote_ind_glo, aggregarting submessage_in and submessage_out into submessage_local 
    136             submessage_out = make_submessage( request(ind_loc)%points_HtoB(remote_ind_glo), & 
     141            submessage_out = make_submessage( field_type, request(ind_loc)%points_HtoB(remote_ind_glo), & 
    137142                                              ind_loc, remote_ind_glo, dim3, dim4, request(1)%vector ) 
    138             submessage_in = make_submessage(request(domglo_loc_ind(remote_ind_glo))%points_BtoH(domloc_glo_ind(ind_loc)), & 
     143            submessage_in = make_submessage( field_type, request(domglo_loc_ind(remote_ind_glo))%points_BtoH(domloc_glo_ind(ind_loc)), & 
    139144                                            domglo_loc_ind(remote_ind_glo), domloc_glo_ind(ind_loc), dim3, dim4, request(1)%vector) 
    140145            submessage_local%src_ind_loc = ind_loc 
     
    149154          ! When data to send to remote_domain, add submessage in message%message_out 
    150155          if( request(ind_loc)%points_HtoB(remote_ind_glo)%npoints > 0 ) then 
    151             submessage_out = make_submessage( request(ind_loc)%points_HtoB(remote_ind_glo), & 
     156            submessage_out = make_submessage( field_type, request(ind_loc)%points_HtoB(remote_ind_glo), & 
    152157                                              ind_loc, remote_ind_glo, dim3, dim4, request(1)%vector ) 
    153158            call array_append_submessage( message_out_tmp, message_out_size, submessage_out ) 
    154159          end if 
    155160          if( request(ind_loc)%points_BtoH(remote_ind_glo)%npoints > 0 ) then 
    156             submessage_in = make_submessage( request(ind_loc)%points_BtoH(remote_ind_glo), & 
     161            submessage_in = make_submessage( field_type, request(ind_loc)%points_BtoH(remote_ind_glo), & 
    157162                                             ind_loc, remote_ind_glo, dim3, dim4, request(1)%vector ) 
    158163            call array_append_submessage( message_in_tmp, message_in_size, submessage_in ) 
     
    185190  contains 
    186191    ! Generate submessage from points 
    187     function make_submessage(points, ind_loc, remote_ind_glo, dim3, dim4, vector) result(submessage) 
    188       use dimensions, only : swap_dimensions, iim, u_pos 
     192    function make_submessage(field_type, points, ind_loc, remote_ind_glo, dim3, dim4, vector) result(submessage) 
     193      use dimensions, only : swap_dimensions, iim, u_pos, z_pos 
     194      integer, intent(in) :: field_type 
    189195      type(t_points), intent(in) :: points 
    190196      integer, intent(in) :: ind_loc, remote_ind_glo, dim3, dim4 
     
    199205      allocate( submessage%displs( points%npoints ) ) 
    200206      submessage%displs(:) = points%i + (points%j-1)*iim 
    201       if(allocated(points%edge)) submessage%displs = submessage%displs + u_pos( points%edge ) 
     207      if(field_type == field_U) submessage%displs = submessage%displs + u_pos( points%elt ) 
     208      if(field_type == field_Z) submessage%displs = submessage%displs + z_pos( points%elt ) 
    202209      allocate(submessage%sign( points%npoints )) 
    203       if( vector ) then 
    204         submessage%sign(:) = (/( domain(ind_loc)%edge_assign_sign(points%edge(k)-1, points%i(k), points%j(k)) ,k=1,points%npoints)/) 
     210      if( vector ) then ! For U fields only 
     211        submessage%sign(:) = (/( domain(ind_loc)%edge_assign_sign(points%elt(k)-1, points%i(k), points%j(k)) ,k=1,points%npoints)/) 
    205212      else 
    206213        submessage%sign(:) = 1 
     
    380387    message%send_seq = message%send_seq + 1 
    381388    !$omp end master 
    382      
     389 
    383390    call enter_profile(profile_mpi_barrier) 
    384391    !$omp barrier 
     
    411418    do i = 1, size( message%message_local ) 
    412419      if( assigned_domain( message%message_local(i)%dest_ind_loc ) ) then 
    413         !$acc loop collapse(2)  
     420        !$acc loop collapse(2) 
    414421        do d4 = 1, dim4 
    415422          do d3 = d3_begin, d3_end 
  • 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.