Changeset d6a0977 in flexpart.git for src/mpi_mod.f90
- Timestamp:
- Dec 14, 2015, 3:10:04 PM (8 years ago)
- Branches:
- master, 10.4.1_pesei, GFS_025, bugfixes+enhancements, dev, release-10, release-10.4.1, scaling-bug, univie
- Children:
- f75967d
- Parents:
- 88d8c3d
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
src/mpi_mod.f90
ra1f4dd6 rd6a0977 89 89 ! MPI tags/requests for send/receive operation 90 90 integer :: tm1 91 integer, parameter :: nvar_async=2 7 !29 :DBG:91 integer, parameter :: nvar_async=26 !27 !29 :DBG: 92 92 !integer, dimension(:), allocatable :: tags 93 93 integer, dimension(:), allocatable :: reqs … … 119 119 logical, parameter :: mp_dbg_out = .false. 120 120 logical, parameter :: mp_time_barrier=.true. 121 logical, parameter :: mp_measure_time=. true.121 logical, parameter :: mp_measure_time=.false. 122 122 logical, parameter :: mp_exact_numpart=.true. 123 123 … … 207 207 write(*,FMT='(80("#"))') 208 208 ! Force "syncronized" version if all processes will call getfields 209 else if ( .not.lmp_sync.and.mp_np.lt.read_grp_min) then209 else if ((.not.lmp_sync.and.mp_np.lt.read_grp_min).or.(mp_np.eq.1)) then 210 210 if (lroot) then 211 211 write(*,FMT='(80("#"))') … … 963 963 ! cloud water/ice: 964 964 if (readclouds) then 965 call MPI_Bcast(clwc(:,:,:,li:ui),d3s1,mp_pp,id_read,MPI_COMM_WORLD,mp_ierr) 966 if (mp_ierr /= 0) goto 600 967 call MPI_Bcast(ciwc(:,:,:,li:ui),d3s1,mp_pp,id_read,MPI_COMM_WORLD,mp_ierr) 968 if (mp_ierr /= 0) goto 600 965 call MPI_Bcast(icloud_stats(:,:,:,li:ui),d2_size1*5,mp_pp,id_read,MPI_COMM_WORLD,mp_ierr) 966 if (mp_ierr /= 0) goto 600 967 968 ! call MPI_Bcast(clwc(:,:,:,li:ui),d3s1,mp_pp,id_read,MPI_COMM_WORLD,mp_ierr) 969 ! if (mp_ierr /= 0) goto 600 970 ! call MPI_Bcast(ciwc(:,:,:,li:ui),d3s1,mp_pp,id_read,MPI_COMM_WORLD,mp_ierr) 971 ! if (mp_ierr /= 0) goto 600 969 972 end if 970 973 … … 1189 1192 ! 1190 1193 ! TODO 1191 ! Transfer cloud water/ice1192 1194 ! 1193 1195 !******************************************************************************* … … 1334 1336 if (readclouds) then 1335 1337 i=i+1 1336 call MPI_Isend( clwc(:,:,:,mind),d3s1,mp_pp,dest,tm1,&1338 call MPI_Isend(icloud_stats(:,:,:,mind),d2s1*5,mp_pp,dest,tm1,& 1337 1339 &MPI_COMM_WORLD,reqs(i),mp_ierr) 1338 1340 if (mp_ierr /= 0) goto 600 1339 i=i+1 1340 1341 call MPI_Isend(ciwc(:,:,:,mind),d3s1,mp_pp,dest,tm1,& 1342 &MPI_COMM_WORLD,reqs(i),mp_ierr) 1343 if (mp_ierr /= 0) goto 600 1344 ! else 1345 ! i=i+2 1341 1342 ! call MPI_Isend(clwc(:,:,:,mind),d3s1,mp_pp,dest,tm1,& 1343 ! &MPI_COMM_WORLD,reqs(i),mp_ierr) 1344 ! if (mp_ierr /= 0) goto 600 1345 ! i=i+1 1346 1347 ! call MPI_Isend(ciwc(:,:,:,mind),d3s1,mp_pp,dest,tm1,& 1348 ! &MPI_COMM_WORLD,reqs(i),mp_ierr) 1349 ! if (mp_ierr /= 0) goto 600 1350 1346 1351 end if 1347 1348 1352 end do 1349 1353 … … 1371 1375 ! 1372 1376 ! TODO 1373 ! Transfer cloud water/ice1374 1377 ! 1375 1378 !******************************************************************************* … … 1390 1393 ! integer :: d3s1,d3s2,d2s1,d2s2 1391 1394 !******************************************************************************* 1392 1393 ! :TODO: don't need these1394 ! d3s1=d3_size11395 ! d3s2=d3_size21396 ! d2s1=d2_size11397 ! d2s2=d2_size21398 1395 1399 1396 ! At the time this immediate receive is posted, memstat is the state of … … 1545 1542 if (readclouds) then 1546 1543 j=j+1 1547 call MPI_Irecv(clwc(:,:,:,mind),d3s1,mp_pp,id_read,MPI_ANY_TAG,& 1548 &MPI_COMM_WORLD,reqs(j),mp_ierr) 1549 if (mp_ierr /= 0) goto 600 1550 j=j+1 1551 call MPI_Irecv(ciwc(:,:,:,mind),d3s1,mp_pp,id_read,MPI_ANY_TAG,& 1552 &MPI_COMM_WORLD,reqs(j),mp_ierr) 1553 if (mp_ierr /= 0) goto 600 1544 1545 call MPI_Irecv(icloud_stats(:,:,:,mind),d2s1*5,mp_pp,id_read,MPI_ANY_TAG,& 1546 &MPI_COMM_WORLD,reqs(j),mp_ierr) 1547 if (mp_ierr /= 0) goto 600 1548 1549 ! call MPI_Irecv(clwc(:,:,:,mind),d3s1,mp_pp,id_read,MPI_ANY_TAG,& 1550 ! &MPI_COMM_WORLD,reqs(j),mp_ierr) 1551 ! if (mp_ierr /= 0) goto 600 1552 ! j=j+1 1553 ! call MPI_Irecv(ciwc(:,:,:,mind),d3s1,mp_pp,id_read,MPI_ANY_TAG,& 1554 ! &MPI_COMM_WORLD,reqs(j),mp_ierr) 1555 ! if (mp_ierr /= 0) goto 600 1556 1554 1557 end if 1555 1558 … … 2089 2092 implicit none 2090 2093 2091 integer, parameter :: li=1, ui=3 ! wfmem indices (i.e, operate on entire field) 2094 integer :: li=1, ui=2 ! wfmem indices (i.e, operate on entire field) 2095 2096 if (.not.lmp_sync) ui=3 2092 2097 2093 2098
Note: See TracChangeset
for help on using the changeset viewer.