From eede0633d33bf6acb1c9cc7e1a4a2d8e945fea05 Mon Sep 17 00:00:00 2001 From: Atanas Trayanov Date: Mon, 13 Jan 2025 13:56:41 -0500 Subject: [PATCH] Fixes #3307. Fixed more MPI communicator assignments --- base/FileIOShared.F90 | 20 ++++++++++---------- 1 file changed, 10 insertions(+), 10 deletions(-) diff --git a/base/FileIOShared.F90 b/base/FileIOShared.F90 index fba72dc619f8..0d638b1c78cd 100644 --- a/base/FileIOShared.F90 +++ b/base/FileIOShared.F90 @@ -544,12 +544,12 @@ subroutine ArrDescrInit(ArrDes,comm,im_world,jm_world,lm_world,nx,ny,num_readers ArrDes%jm_world=jm_world ArrDes%lm_world=lm_world - ArrDes%readers_comm = readers_comm - ArrDes%ioscattercomm = ioscattercomm - ArrDes%writers_comm = writers_comm - ArrDes%iogathercomm = iogathercomm - ArrDes%xcomm = xcomm - ArrDes%ycomm = ycomm + call MPI_Comm_Dup(readers_comm, ArrDes%readers_comm, status) + call MPI_Comm_Dup(ioscattercomm, ArrDes%IOscattercomm, status) + call MPI_Comm_Dup(writers_comm, ArrDes%writers_comm, status) + call MPI_Comm_Dup(iogathercomm, ArrDes%IOgathercomm, status) + call MPI_Comm_Dup(xcomm, ArrDes%Xcomm, status) + call MPI_Comm_Dup(ycomm, ArrDes%Ycomm, status) call mpi_comm_rank(arrdes%ycomm,arrdes%myrow,status) _VERIFY(status) @@ -598,13 +598,13 @@ subroutine ArrDescrSet(ArrDes, offset, & call MPI_Comm_Dup(readers_comm, ArrDes%readers_comm, status) end if if(present(ioscattercomm)) then - call MPI_Comm_Dup(ioscattercomm, ArrDes%ioscattercomm, status) + call MPI_Comm_Dup(IOscattercomm, ArrDes%IOscattercomm, status) end if if(present(writers_comm )) then call MPI_Comm_Dup(writers_comm, ArrDes%writers_comm, status) end if if(present(iogathercomm)) then - call MPI_Comm_Dup(iogathercomm, ArrDes%iogathercomm, status) + call MPI_Comm_Dup(IOgathercomm, ArrDes%IOgathercomm, status) end if if(present(i1 )) ArrDes%i1 => i1 if(present(in )) ArrDes%in => in @@ -645,7 +645,7 @@ subroutine ArrDescrCreateWriterComm(arrdes, full_comm, num_writers, rc) call MPI_COMM_SPLIT(full_comm, color, myid, arrdes%writers_comm, status) _VERIFY(status) if (num_writers==ny) then - arrdes%IOgathercomm = arrdes%Xcomm + call MPI_Comm_Dup(arrdes%Xcomm, ArrDes%IOgathercomm, status) else j = arrdes%NY0 - mod(arrdes%NY0-1,ny_by_writers) call MPI_COMM_SPLIT(full_comm, j, myid, arrdes%IOgathercomm, status) @@ -693,7 +693,7 @@ subroutine ArrDescrCreateReaderComm(arrdes, full_comm, num_readers, rc) call MPI_COMM_SPLIT(full_comm, color, MYID, arrdes%readers_comm, status) _VERIFY(status) if (num_readers==ny) then - arrdes%IOscattercomm = arrdes%Xcomm + call MPI_Comm_Dup(arrdes%Xcomm, ArrDes%IOscattercomm, status) else j = arrdes%NY0 - mod(arrdes%NY0-1,ny_by_readers) call MPI_COMM_SPLIT(full_comm, j, MYID, arrdes%IOscattercomm, status)