program fd include 'mpif.h' integer n,p real tol parameter (n=100,p=4,maxsize=(n-1)/p+2,tol=0.00001) real new(0:n,0:maxsize), old(0:n,0:maxsize),iteration, maxerr, & err, maxerrG integer noprocs, nid, remainder, size, i, j, error character suffix*2 integer status(MPI_STATUS_SIZE), req_send10, req_send20, & req_recv10, req_recv20 call MPI_Init(error) call MPI_Comm_rank(MPI_COMM_WORLD, nid, error) call MPI_Comm_size(MPI_COMM_WORLD, noprocs, error) if (noprocs .ne. p) then if (nid .eq. 0) write(6,*)'Number of processes should be...',p call MPI_Abort(MPI_COMM_WORLD,-1,error) end if remainder = mod ((n-1),noprocs) size = (n - 1 - remainder)/noprocs if (nid .lt. remainder) then size = size + 2 else size = size + 1 end if do 8 i=0,size do 5 j=0,n old(j,i) = 0.0 5 continue 8 continue do 10 i=0,size new(0,i) = 1.0 new(n,i) = 1.0 old(0,i) = 1.0 old(n,i) = 1.0 10 continue if (nid .eq. 0) then do 20 j=1,n-1 new(j,0) = 1.0 old(j,0) = 1.0 20 continue end if if (nid .eq. noprocs-1) then do 30 j=1,n-1 new(j,size) = 1.0 old(j,size) = 1.0 30 continue end if maxerr = iteration(n,size,old,new,1,size) call MPI_Allreduce(maxerr,maxerrG,1,MPI_REAL,MPI_MAX, & MPI_COMM_WORLD,error) 40 continue if (maxerrG .gt. tol) then if (nid .eq. 0) write(6,*)maxerrG do 60 i=0,size do 50 j=0,n old(j,i) = new(j,i) 50 continue 60 continue req_send10 = MPI_REQUEST_NULL req_recv20 = MPI_REQUEST_NULL if (nid .lt. noprocs-1) then call MPI_Isend(old(1,size-1),n-1,MPI_REAL,nid+1,10, & MPI_COMM_WORLD,req_send10,error) call MPI_Irecv(old(1,size),n-1,MPI_REAL,nid+1,20, & MPI_COMM_WORLD,req_recv20,error) end if req_send20 = MPI_REQUEST_NULL req_recv10 = MPI_REQUEST_NULL if (nid .gt. 0) then call MPI_Isend(old(1,1),n-1,MPI_REAL,nid-1,20, & MPI_COMM_WORLD,req_send20,error) call MPI_Irecv(old(1,0),n-1,MPI_REAL,nid-1,10, & MPI_COMM_WORLD,req_recv10,error) end if maxerr = iteration(n,size,old,new,2,size-1) if (nid .lt. noprocs-1) call MPI_Wait(req_recv20,status,error) err = iteration(n,size,old,new,size-1,size) if (err .gt. maxerr) maxerr = err if (nid .gt. 0) call MPI_Wait(req_recv10,status,error) err = iteration(n,size,old,new,1,2) if (err .gt. maxerr) maxerr = err call MPI_Allreduce(maxerr,maxerrG,1,MPI_REAL,MPI_MAX, & MPI_COMM_WORLD,error) goto 40 end if write(suffix,'(i2)')nid if (nid .le. 9) suffix(1:1) = '0' open (10,file='Solution'//suffix//'.Txt',form='formatted') if (nid .eq. 0) then do 70 j=0,n write(10,'(f6.4)')new(j,0) 70 continue end if do 90 i=1,size-1 do 80 j=0,n write(10,'(f6.4)')new(j,i) 80 continue 90 continue if (nid .eq. noprocs-1) then do 100 j=0,n write(10,'(f6.4)')new(j,size) 100 continue end if close(10) call MPI_Finalize(error) stop end ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc real function iteration(n, size, old, new, start, finish) integer n, size, start, finish, i, j real old(0:n,0:size), new(0:n,0:size), diff iteration = 0.0 do 1020 i=start,finish-1 do 1010 j=1,n-1 new(j,i) = 0.25 * ( old(j,i+1) + old(j,i-1) + & old(j+1,i) + old(j-1,i) ) diff = abs(new(j,i) - old(j,i)) if (iteration .lt. diff) iteration = diff 1010 continue 1020 continue return end