Ignore:
Timestamp:
07/15/19 12:29:31 (5 years ago)
Author:
adurocher
Message:

trunk : GPU implementation with OpenACC ( merge from glcp.idris.fr )

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

Legend:

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

    r892 r953  
    1414  INTEGER,SAVE :: id_mpi ! id for profiling 
    1515 
     16  INTEGER, SAVE :: device_id 
     17 
    1618  INTERFACE allocate_mpi_buffer 
    1719    MODULE PROCEDURE allocate_mpi_buffer_r2, allocate_mpi_buffer_r3,allocate_mpi_buffer_r4 
     
    4648  USE xios 
    4749#endif 
     50  USE abort_mod 
    4851  IMPLICIT NONE 
    4952    CHARACTER(LEN=256) :: required_mode_str 
     
    7578      END SELECT 
    7679       
     80      IF (required_mode==MPI_THREAD_SERIALIZED .OR. required_mode==MPI_THREAD_MULTIPLE) THEN 
     81        CALL abort_acc("mpi_threading_mode /= 'single' .AND. mpi_threading_mode /= 'funneled'") 
     82      ENDIF 
    7783 
    7884      IF (required_mode==MPI_THREAD_SINGLE)     PRINT*,'MPI_INIT_THREAD : MPI_SINGLE_THREAD required' 
     
    120126    ENDIF 
    121127     
     128 
     129#ifdef _OPENACC 
     130    device_id = setDevice(mpi_size, mpi_rank) 
     131    PRINT *, 'GPU device ', device_id 
     132#else 
     133     device_id = -1 
     134#endif 
     135 
    122136  END SUBROUTINE  init_mpipara 
    123137 
     
    233247  END SUBROUTINE free_mpi_buffer_r4 
    234248    
     249#ifdef _OPENACC 
     250  FUNCTION setDevice(nprocs, myrank) 
     251    use iso_c_binding 
     252    use openacc 
     253    USE mpi_mod 
     254    implicit none 
     255     
     256    interface  
     257       function gethostid() bind(C) 
     258       use iso_c_binding 
     259       integer(C_INT) :: gethostid 
     260     end function gethostid 
     261  end interface 
     262   
     263  integer, intent(in) :: nprocs, myrank 
     264  integer :: hostids(nprocs), localprocs(nprocs) 
     265  integer :: hostid, ierr, numdev, mydev, i, numlocal 
     266  integer :: setDevice 
     267 
     268  ! get the hostids so we can determine what other processes are on this node 
     269  hostid = gethostid() 
     270  call mpi_allgather(hostid,1,MPI_INTEGER,hostids,1,MPI_INTEGER, MPI_COMM_WORLD, ierr) 
     271 
     272  ! determine which processors are on this node 
     273  numlocal = 0 
     274  localprocs(:) = 0 
     275  do i = 1, nprocs 
     276     if (hostid == hostids(i)) then 
     277        localprocs(i) = numlocal 
     278        numlocal = numlocal + 1 
     279     end if 
     280  end do 
     281 
     282  ! get the number of device on this node 
     283  numdev = acc_get_num_devices(ACC_DEVICE_NVIDIA) 
     284 
     285  if (numdev < 1) then 
     286     print *, "Error: there are no devices available on this host. ABORTING", myrank 
     287     stop 
     288  end if 
     289 
     290  ! print a warning if the number of devices is less than the number of processes on this node. Having multiple processes share a devices is not recommended 
     291  if (numdev < numlocal) then 
     292     if (localprocs(myrank+1) == 1) then 
     293        ! print warning message only once per node 
     294        print *, "WARNING: the number of process is greater than the number of GPUs.", myrank 
     295     end if 
     296     mydev = mod(localprocs(myrank+1), numdev) 
     297  else 
     298     mydev = localprocs(myrank+1) 
     299  end if 
     300 
     301  call acc_set_device_num(mydev,ACC_DEVICE_NVIDIA) 
     302  call acc_init(ACC_DEVICE_NVIDIA) 
     303  setDevice = acc_get_device_num(ACC_DEVICE_NVIDIA) 
     304END FUNCTION setDevice 
     305 
     306#endif 
     307 
     308      
     309      
    235310END MODULE mpipara 
  • codes/icosagcm/trunk/src/parallel/transfert_mpi.f90

    r901 r953  
    55   
    66  TYPE array 
    7     INTEGER,POINTER :: value(:) 
    8     INTEGER,POINTER :: sign(:) 
     7    INTEGER,POINTER :: value(:)=>null() 
     8    INTEGER,POINTER :: sign(:)=>null() 
    99    INTEGER         :: domain 
    1010    INTEGER         :: rank 
     
    1313    INTEGER         :: offset 
    1414    INTEGER         :: ireq 
    15     INTEGER,POINTER :: buffer(:) 
    16     REAL,POINTER    :: buffer_r(:) 
    17     INTEGER,POINTER :: src_value(:) 
     15    INTEGER,POINTER :: buffer(:)=>null() 
     16    REAL,POINTER    :: buffer_r(:)=>null() 
     17    INTEGER,POINTER :: src_value(:)=>null() 
    1818  END TYPE array 
    1919   
     
    7777    LOGICAL :: open      ! for debug 
    7878    INTEGER :: number 
     79    LOGICAL :: ondevice=.false. !<Ready to transfer ondevice field 
    7980  END TYPE t_message 
    8081 
     
    942943!$OMP BARRIER 
    943944!$OMP MASTER 
     945    !init off-device 
     946    message%ondevice=.false. 
     947 
    944948    IF(PRESENT(name)) THEN 
    945949       message%name = TRIM(name) 
     
    10581062    TYPE(t_message) :: message 
    10591063 
    1060     INTEGER :: ireq 
     1064    INTEGER :: ireq, ibuff 
    10611065 
    10621066!$OMP BARRIER 
     
    10841088        ENDDO 
    10851089 
    1086       ENDIF       
     1090      ENDIF 
    10871091    ENDIF 
    1088      
     1092 
     1093    !deallocate device data if ondevice 
     1094    if(message%ondevice) then 
     1095      do ireq=1, ndomain 
     1096        do ibuff=1,message%request(ireq)%nsend 
     1097          !$acc exit data delete(message%request(ireq)%send(ibuff)%buffer(:)) 
     1098          !$acc exit data delete(message%request(ireq)%send(ibuff)%buffer_r(:)) 
     1099          !$acc exit data delete(message%request(ireq)%send(ibuff)%sign(:)) 
     1100          !$acc exit data delete(message%request(ireq)%send(ibuff)%src_value(:)) 
     1101          !$acc exit data delete(message%request(ireq)%send(ibuff)%value(:)) 
     1102        end do 
     1103        do ibuff=1,message%request(ireq)%nrecv 
     1104          !$acc exit data delete(message%request(ireq)%recv(ibuff)%buffer(:)) 
     1105          !$acc exit data delete(message%request(ireq)%recv(ibuff)%buffer_r(:)) 
     1106          !$acc exit data delete(message%request(ireq)%recv(ibuff)%sign(:)) 
     1107          !$acc exit data delete(message%request(ireq)%recv(ibuff)%src_value(:)) 
     1108          !$acc exit data delete(message%request(ireq)%recv(ibuff)%value(:)) 
     1109        end do 
     1110      end do 
     1111      DO ireq=1,message%nreq 
     1112        !$acc exit data delete(message%buffers(ireq)%r) 
     1113      ENDDO 
     1114      message%ondevice=.false. 
     1115    end if 
     1116 
    10891117    DEALLOCATE(message%mpi_req) 
    10901118    DEALLOCATE(message%buffers) 
     
    11201148 
    11211149 
     1150  !!!Update buffers on device for 'message' 
     1151  !!! does create_device_message when not already ondevice 
     1152  SUBROUTINE update_device_message_mpi(message) 
     1153    USE domain_mod 
     1154    IMPLICIT NONE 
     1155    TYPE(t_message), intent(inout) :: message 
     1156    INTEGER :: ireq, ibuff 
     1157 
     1158    !if(.not. message%ondevice) call create_device_message_mpi(message) 
     1159 
     1160    do ireq=1, ndomain 
     1161      do ibuff=1,message%request(ireq)%nsend 
     1162        !device buffers updated even if pointers not attached : 
     1163        !non allocated buffers in 'message' must be set to NULL() 
     1164        !$acc enter data copyin(message%request(ireq)%send(ibuff)%buffer(:)) async 
     1165        !$acc enter data copyin(message%request(ireq)%send(ibuff)%buffer_r(:)) async 
     1166        !$acc enter data copyin(message%request(ireq)%send(ibuff)%sign(:)) async 
     1167        !$acc enter data copyin(message%request(ireq)%send(ibuff)%src_value(:)) async 
     1168        !$acc enter data copyin(message%request(ireq)%send(ibuff)%value(:)) async 
     1169      end do 
     1170      do ibuff=1,message%request(ireq)%nrecv 
     1171        !$acc enter data copyin(message%request(ireq)%recv(ibuff)%buffer(:)) async 
     1172        !$acc enter data copyin(message%request(ireq)%recv(ibuff)%buffer_r(:)) async 
     1173        !$acc enter data copyin(message%request(ireq)%recv(ibuff)%sign(:)) async 
     1174        !$acc enter data copyin(message%request(ireq)%recv(ibuff)%src_value(:)) async 
     1175        !$acc enter data copyin(message%request(ireq)%recv(ibuff)%value(:)) async 
     1176      end do 
     1177    end do 
     1178    DO ireq=1,message%nreq 
     1179      !$acc enter data copyin(message%buffers(ireq)%r) async 
     1180    ENDDO 
     1181    message%ondevice=.true. 
     1182  END SUBROUTINE 
     1183 
     1184  !TODO : add openacc with multiple process 
    11221185  SUBROUTINE send_message_mpi(field,message) 
    11231186  USE abort_mod 
     
    11291192  USE omp_para 
    11301193  USE trace 
     1194  USE abort_mod 
    11311195  IMPLICIT NONE 
    11321196    TYPE(t_field),POINTER :: field(:) 
     
    11451209    INTEGER :: dim3,dim4,d3,d4 
    11461210    INTEGER,POINTER :: src_value(:) 
    1147     INTEGER :: offset,msize,rank 
     1211    INTEGER :: offset,msize,moffset,rank 
    11481212    INTEGER :: lbegin, lend 
    11491213    INTEGER :: max_req 
     
    11521216 
    11531217    CALL enter_profile(id_mpi) 
     1218 
     1219    !Prepare 'message' for on-device copies if field is on device 
     1220    if( field(1)%ondevice .AND. .NOT. message%ondevice ) call update_device_message_mpi(message) 
    11541221 
    11551222!$OMP BARRIER 
     
    11951262              buffer_r=>message%buffers(ireq)%r 
    11961263              offset=send%offset 
    1197                                              
    1198               DO n=1,send%size 
     1264              msize=send%size 
     1265              !$acc parallel loop default(present) async if (field(ind)%ondevice) 
     1266              DO n=1,msize 
    11991267                buffer_r(offset+n)=rval2d(value(n)) 
    12001268              ENDDO 
    12011269 
    12021270              IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN 
     1271                CALL abort_acc("mpi_threading_mode==MPI_THREAD_SERIALIZED") 
    12031272                !$OMP CRITICAL             
    12041273                CALL MPI_ISEND(buffer_r,send%size,MPI_REAL8,send%rank,               & 
     
    12061275                !$OMP END CRITICAL 
    12071276              ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN 
     1277                CALL abort_acc("mpi_threading_mode==MPI_THREAD_MULTIPLE") 
    12081278                CALL MPI_ISEND(buffer_r,send%size,MPI_REAL8,send%rank,               & 
    12091279                  send%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr) 
     
    12281298              src_rval2d=>field(recv%domain)%rval2d 
    12291299              sgn=>recv%sign 
    1230  
    1231               DO n=1,recv%size 
     1300              msize=recv%size 
     1301              !$acc parallel loop default(present) async if (field(ind)%ondevice) 
     1302              DO n=1,msize 
    12321303                rval2d(value(n))=src_rval2d(src_value(n))*sgn(n) 
    12331304              ENDDO 
     
    12391310              buffer_r=>message%buffers(ireq)%r 
    12401311              IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN 
     1312                CALL abort_acc("mpi_threading_mode==MPI_THREAD_SERIALIZED") 
    12411313               !$OMP CRITICAL             
    12421314                CALL MPI_IRECV(buffer_r,recv%size,MPI_REAL8,recv%rank,               & 
     
    12441316               !$OMP END CRITICAL 
    12451317              ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN 
     1318                 CALL abort_acc("mpi_threading_mode==MPI_THREAD_MULTIPLE") 
    12461319                 CALL MPI_IRECV(buffer_r,recv%size,MPI_REAL8,recv%rank,              & 
    12471320                   recv%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr) 
     
    12781351              buffer_r=>message%buffers(ireq)%r 
    12791352 
    1280               offset=send%offset*dim3 + (lbegin-1)*send%size 
    1281            
     1353              msize=send%size 
     1354              moffset=send%offset 
    12821355              CALL trace_start("copy_to_buffer") 
    12831356 
     1357              !$acc parallel loop default(present) async if (field(ind)%ondevice) 
    12841358              DO d3=lbegin,lend 
    1285                 DO n=1,send%size 
     1359                offset=moffset*dim3 + (d3-1)*msize 
     1360                !$acc loop 
     1361                DO n=1,msize 
    12861362                  buffer_r(n+offset)=rval3d(value(n),d3) 
    12871363                ENDDO 
    1288                 offset=offset+send%size 
    12891364              ENDDO 
    12901365              CALL trace_end("copy_to_buffer") 
     
    12961371              IF (is_omp_level_master) THEN 
    12971372                IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN 
     1373                  CALL abort_acc("mpi_threading_mode==MPI_THREAD_SERIALIZED") 
    12981374                  !$OMP CRITICAL    
    12991375                  CALL MPI_ISEND(buffer_r,send%size*dim3,MPI_REAL8,send%rank,        & 
     
    13011377                  !$OMP END CRITICAL 
    13021378                ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN 
     1379                  CALL abort_acc("mpi_threading_mode==MPI_THREAD_MULTIPLE") 
    13031380                  CALL MPI_ISEND(buffer_r,send%size*dim3,MPI_REAL8,send%rank,        & 
    13041381                    send%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr) 
     
    13351412              src_rval3d=>field(recv%domain)%rval3d 
    13361413              sgn=>recv%sign 
     1414              msize=recv%size 
    13371415 
    13381416              CALL trace_start("copy_data") 
    1339  
     1417              !$acc parallel loop collapse(2) default(present) async if (field(ind)%ondevice) 
    13401418              DO d3=lbegin,lend 
    1341                 DO n=1,recv%size 
     1419                DO n=1,msize 
    13421420                  rval3d(value(n),d3)=src_rval3d(src_value(n),d3)*sgn(n) 
    13431421                ENDDO 
     
    13511429              IF (is_omp_level_master) THEN 
    13521430                IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN 
     1431                  CALL abort_acc("mpi_threading_mode==MPI_THREAD_SERIALIZED") 
    13531432                  !$OMP CRITICAL 
    13541433                  CALL MPI_IRECV(buffer_r,recv%size*dim3,MPI_REAL8,recv%rank,        & 
     
    13561435                  !$OMP END CRITICAL 
    13571436                ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN 
     1437                  CALL abort_acc("mpi_threading_mode==MPI_THREAD_MULTIPLE") 
    13581438                  CALL MPI_IRECV(buffer_r,recv%size*dim3,MPI_REAL8,recv%rank,        & 
    13591439                    recv%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr) 
     
    13901470              ireq=send%ireq 
    13911471              buffer_r=>message%buffers(ireq)%r 
     1472              msize=send%size 
     1473              moffset=send%offset 
    13921474 
    13931475              CALL trace_start("copy_to_buffer") 
     1476              !$acc parallel loop default(present) collapse(2) async if (field(ind)%ondevice) 
    13941477              DO d4=1,dim4 
    1395                 offset=send%offset*dim3*dim4 + dim3*send%size*(d4-1) +               & 
    1396                   (lbegin-1)*send%size 
    1397  
    13981478                DO d3=lbegin,lend 
    1399                   DO n=1,send%size 
     1479                  offset=moffset*dim3*dim4 + dim3*msize*(d4-1) +               & 
     1480                    (d3-1)*msize 
     1481                  !$acc loop 
     1482                  DO n=1,msize 
    14001483                    buffer_r(n+offset)=rval4d(value(n),d3,d4) 
    14011484                  ENDDO 
    1402                   offset=offset+send%size 
    14031485                ENDDO 
    14041486              ENDDO 
     
    14111493              IF (is_omp_level_master) THEN 
    14121494                IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN 
     1495                  CALL abort_acc("mpi_threading_mode==MPI_THREAD_SERIALIZED") 
    14131496                  !$OMP CRITICAL 
    14141497                  CALL MPI_ISEND(buffer_r,send%size*dim3*dim4,MPI_REAL8,send%rank,   & 
     
    14161499                  !$OMP END CRITICAL 
    14171500                ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN 
     1501                  CALL abort_acc("mpi_threading_mode==MPI_THREAD_MULTIPLE") 
    14181502                  CALL MPI_ISEND(buffer_r,send%size*dim3*dim4,MPI_REAL8,send%rank,   & 
    14191503                    send%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr) 
     
    14511535              src_rval4d=>field(recv%domain)%rval4d 
    14521536              sgn=>recv%sign 
    1453  
     1537              msize=recv%size 
    14541538              CALL trace_start("copy_data") 
     1539 
     1540              !$acc parallel loop collapse(3) default(present) async if (field(ind)%ondevice) 
    14551541              DO d4=1,dim4 
    14561542                DO d3=lbegin,lend 
    1457                   DO n=1,recv%size 
     1543                  DO n=1,msize 
    14581544                    rval4d(value(n),d3,d4)=src_rval4d(src_value(n),d3,d4)*sgn(n) 
    14591545                  ENDDO 
     
    14691555              IF (is_omp_level_master) THEN 
    14701556                IF (mpi_threading_mode==MPI_THREAD_SERIALIZED) THEN 
     1557                  CALL abort_acc("mpi_threading_mode==MPI_THREAD_SERIALIZED") 
    14711558                 !$OMP CRITICAL            
    14721559                  CALL MPI_IRECV(buffer_r,recv%size*dim3*dim4,MPI_REAL8,recv%rank,   & 
     
    14741561                  !$OMP END CRITICAL 
    14751562                ELSE IF (mpi_threading_mode==MPI_THREAD_MULTIPLE) THEN 
     1563                  CALL abort_acc("mpi_threading_mode==MPI_THREAD_MULTIPLE") 
    14761564                  CALL MPI_IRECV(buffer_r,recv%size*dim3*dim4,MPI_REAL8,recv%rank,   & 
    14771565                    recv%tag+1000*message%number,comm_icosa, message%mpi_req(ireq),ierr)  
     
    14851573 
    14861574      IF (mpi_threading_mode==MPI_THREAD_FUNNELED .OR. mpi_threading_mode==MPI_THREAD_SINGLE) THEN 
     1575!$acc wait 
    14871576!$OMP BARRIER 
    14881577!$OMP MASTER         
     
    14921581          msize=message%buffers(ireq)%size 
    14931582          rank=message%buffers(ireq)%rank 
    1494           CALL MPI_ISEND(buffer_r,msize,MPI_REAL8,rank,1000*message%number,          & 
    1495             comm_icosa, message%mpi_req(ireq),ierr) 
     1583          ! Using the "if" clause on the "host_data" directive seems to cause a crash at compilation 
     1584          IF (message%ondevice) THEN 
     1585            !$acc host_data use_device(buffer_r) ! if (message%ondevice) 
     1586            CALL MPI_ISEND(buffer_r,msize,MPI_REAL8,rank,1000*message%number,          & 
     1587              comm_icosa, message%mpi_req(ireq),ierr) 
     1588            !$acc end host_data 
     1589          ELSE 
     1590            CALL MPI_ISEND(buffer_r,msize,MPI_REAL8,rank,1000*message%number,          & 
     1591              comm_icosa, message%mpi_req(ireq),ierr) 
     1592          ENDIF 
    14961593        ENDDO 
    14971594 
     
    15001597          msize=message%buffers(ireq)%size 
    15011598          rank=message%buffers(ireq)%rank 
    1502           CALL MPI_IRECV(buffer_r,msize,MPI_REAL8,rank,1000*message%number,          & 
    1503             comm_icosa, message%mpi_req(ireq),ierr) 
     1599          ! Using the "if" clause on the "host_data" directive seems to cause a crash at compilation 
     1600          IF (message%ondevice) THEN 
     1601            !$acc host_data use_device(buffer_r) ! if (message%ondevice) 
     1602            CALL MPI_IRECV(buffer_r,msize,MPI_REAL8,rank,1000*message%number,          & 
     1603              comm_icosa, message%mpi_req(ireq),ierr) 
     1604            !$acc end host_data 
     1605          ELSE 
     1606            CALL MPI_IRECV(buffer_r,msize,MPI_REAL8,rank,1000*message%number,          & 
     1607              comm_icosa, message%mpi_req(ireq),ierr) 
     1608          ENDIF 
    15041609        ENDDO 
    15051610 
     
    15521657    INTEGER :: ind,n 
    15531658    INTEGER :: dim3,dim4,d3,d4,lbegin,lend 
    1554     INTEGER :: offset 
     1659    INTEGER :: offset, msize, moffset 
    15551660 
    15561661    message%open=.FALSE. 
     
    15861691              sgn=>recv%sign 
    15871692              offset=recv%offset 
    1588               DO n=1,recv%size 
     1693              msize=recv%size 
     1694              !$acc parallel loop default(present) async if (field(ind)%ondevice) 
     1695              DO n=1,msize 
    15891696                rval2d(value(n))=buffer_r(n+offset)*sgn(n)   
    15901697              ENDDO 
     
    16211728     
    16221729              CALL distrib_level(1,dim3, lbegin,lend) 
    1623               offset=recv%offset*dim3 + (lbegin-1)*recv%size 
     1730              msize=recv%size 
     1731              moffset=recv%offset 
    16241732              CALL trace_start("copy_from_buffer") 
    16251733               
    16261734              IF (req%vector) THEN 
     1735                !$acc parallel loop default(present) async if (field(ind)%ondevice) 
    16271736                DO d3=lbegin,lend 
    1628                   DO n=1,recv%size 
     1737                  offset=moffset*dim3 + (d3-1)*msize 
     1738                  !$acc loop 
     1739                  DO n=1,msize 
    16291740                    rval3d(value(n),d3)=buffer_r(n+offset)*sgn(n)   
    16301741                  ENDDO 
    1631                   offset=offset+recv%size 
    16321742                ENDDO 
    16331743              ELSE 
     1744                !$acc parallel loop default(present) async if (field(ind)%ondevice) 
    16341745                DO d3=lbegin,lend 
    1635                   DO n=1,recv%size 
     1746                  offset=moffset*dim3 + (d3-1)*msize 
     1747                  !$acc loop 
     1748                  DO n=1,msize 
    16361749                    rval3d(value(n),d3)=buffer_r(n+offset)   
    16371750                  ENDDO 
    1638                   offset=offset+recv%size 
    16391751                ENDDO 
    16401752              ENDIF 
     
    16701782              CALL distrib_level(1,dim3, lbegin,lend) 
    16711783              dim4=size(rval4d,3) 
     1784              msize=recv%size 
     1785              moffset=recv%offset 
    16721786              CALL trace_start("copy_from_buffer") 
     1787              !$acc parallel loop default(present) collapse(2) async if (field(ind)%ondevice) 
    16731788              DO d4=1,dim4 
    1674                 offset=recv%offset*dim3*dim4 + dim3*recv%size*(d4-1) +               & 
    1675                   (lbegin-1)*recv%size 
    16761789                DO d3=lbegin,lend 
    1677                   DO n=1,recv%size 
     1790                  offset=moffset*dim3*dim4 + dim3*msize*(d4-1) +               & 
     1791                    (d3-1)*msize 
     1792                  !$acc loop 
     1793                  DO n=1,msize 
    16781794                    rval4d(value(n),d3,d4)=buffer_r(n+offset)*sgn(n)  
    16791795                  ENDDO 
    1680                   offset=offset+recv%size 
    16811796                ENDDO 
    16821797              ENDDO 
     
    16971812!    CALL trace_end("wait_message_mpi") 
    16981813!$OMP BARRIER 
    1699      
     1814   
    17001815    CALL exit_profile(id_mpi) 
    17011816 
Note: See TracChangeset for help on using the changeset viewer.