! part2d_lib_h module for QuickPIC Open Source 1.0 ! update: 04/18/2016 module part2d_lib use mpi implicit none public ! integer :: nproc, lgrp, mreal, mint, mcplx, mdouble, lworld common /PPARMS/ nproc, lgrp, mreal, mint, mcplx, mdouble, lworld ! interface subroutine PPDBLKP2L(part,kpic,npp,noff,nppmx,idimp,npmax,mx,my& &,mx1,mxyp1,irc) implicit none integer, intent(in) :: idimp, npmax, mx, my, mx1, mxyp1, npp integer, intent(in) :: noff integer, intent(inout) :: nppmx, irc real, dimension(idimp,npmax), intent(in) :: part integer, dimension(mxyp1), intent(inout) :: kpic end subroutine end interface ! interface subroutine PPPMOVIN2L(part,ppart,kpic,npp,noff,nppmx,idimp, & &npmax,mx,my,mx1,mxyp1,irc) implicit none integer, intent(in) :: nppmx, idimp, npmax, mx, my, mx1, mxyp1 integer, intent(in) :: npp, noff integer, intent(inout) :: irc real, dimension(idimp,npmax), intent(in) :: part real, dimension(idimp,nppmx,mxyp1), intent(inout) :: ppart integer, dimension(mxyp1), intent(inout) :: kpic end subroutine end interface ! interface subroutine PPPCHECK2L(ppart,kpic,noff,nyp,idimp,nppmx,nx,mx,my,& &mx1,myp1,irc) implicit none integer, intent(in) :: idimp, nppmx, nx, mx, my, mx1, myp1 integer, intent(in) :: noff, nyp integer, intent(inout) :: irc real, dimension(idimp,nppmx,mx1*myp1), intent(in) :: ppart integer, dimension(mx1*myp1), intent(in) :: kpic end subroutine end interface ! interface subroutine PPGPPOST2L(ppart,q,kpic,noff,idimp,nppmx,mx,my, & &nxv,nypmx,mx1,mxyp1) implicit none integer, intent(in) :: idimp, nppmx, mx, my, nxv, nypmx integer, intent(in) :: mx1, mxyp1 integer, intent(in) :: noff real, dimension(idimp,nppmx,mxyp1), intent(in) :: ppart real, dimension(nxv,nypmx), intent(inout) :: q integer, dimension(mxyp1), intent(in) :: kpic end subroutine end interface ! interface subroutine PPGRDCJPPOST2L_QP(ppart,fxy,bxy,psit,cu,dcu,amu,kpic,noff,& &nyp,qbm,dt,ci,idimp,nppmx,nx,mx,my,nxv,nypmx,mx1,mxyp1,dex) implicit none integer, intent(in) :: noff, nyp, idimp, nppmx, nx, mx, my integer, intent(in) :: nxv, nypmx, mx1, mxyp1 real, intent(in) :: qbm, dt, ci real, dimension(idimp,nppmx,mxyp1), intent(inout) :: ppart real, dimension(2,nxv,nypmx), intent(in) :: fxy real, dimension(3,nxv,nypmx), intent(in) :: bxy real, dimension(nxv,nypmx), intent(in) :: psit real, dimension(3,nxv,nypmx), intent(inout) :: cu real, dimension(2,nxv,nypmx), intent(inout) :: dcu real, dimension(3,nxv,nypmx), intent(inout) :: amu integer, dimension(mxyp1), intent(in) :: kpic real, intent(in) :: dex end subroutine end interface ! interface subroutine PPGRBPPUSHF23L_QP(ppart,fxy,bxy,psit,kpic,ncl,ihole,& &noff,nyp,qbm,dt,dtc,ci,ek,idimp,nppmx,nx,ny,mx,my,nxv,nypmx,mx1,mx& &yp1,ntmax,irc,dex) implicit none integer noff, nyp, idimp, nppmx, nx, ny, mx, my, nxv, nypmx integer mx1, mxyp1, ntmax, irc real qbm, dt, dtc, ci, ek, dex real ppart, fxy, bxy, psit integer kpic, ncl, ihole dimension ppart(idimp,nppmx,mxyp1) dimension fxy(2,nxv,nypmx), bxy(3,nxv,nypmx) dimension psit(nxv,nypmx) dimension kpic(mxyp1), ncl(8,mxyp1) dimension ihole(2,ntmax+1,mxyp1) end subroutine end interface ! interface subroutine PPPORDER2LA(ppart,ppbuff,sbufl,sbufr,kpic,ncl,ihole,& &ncll,nclr,noff,nyp,idimp,nppmx,nx,ny,mx,my,mx1,myp1,npbmx,ntmax, & &nbmax,irc) implicit none integer, intent(in) :: idimp, nppmx, nx, ny, mx, my, mx1, myp1 integer, intent(in) :: npbmx, ntmax, nbmax, noff, nyp integer, intent(inout) :: irc real, dimension(idimp,nppmx,mx1*myp1), intent(inout) :: ppart real, dimension(idimp,npbmx,mx1*myp1), intent(inout) :: ppbuff real, dimension(idimp,nbmax), intent(inout) :: sbufl, sbufr integer, dimension(mx1*myp1), intent(in) :: kpic integer, dimension(8,mx1*myp1), intent(inout) :: ncl integer, dimension(2,ntmax+1,mx1*myp1), intent(inout) :: ihole integer, dimension(3,mx1), intent(inout) :: ncll, nclr end subroutine end interface ! interface subroutine PPPORDERF2LA(ppart,ppbuff,sbufl,sbufr,ncl,ihole,ncll& &,nclr,idimp,nppmx,mx1,myp1,npbmx,ntmax,nbmax,irc) implicit none integer, intent(in) :: idimp, nppmx, mx1, myp1, npbmx, ntmax integer, intent(in) :: nbmax integer, intent(inout) :: irc real, dimension(idimp,nppmx,mx1*myp1), intent(inout) :: ppart real, dimension(idimp,npbmx,mx1*myp1), intent(inout) :: ppbuff real, dimension(idimp,nbmax), intent(inout) :: sbufl, sbufr integer, dimension(8,mx1*myp1), intent(inout) :: ncl integer, dimension(2,ntmax+1,mx1*myp1), intent(in) :: ihole integer, dimension(3,mx1), intent(inout) :: ncll, nclr end subroutine end interface ! interface subroutine PPPORDER2LB(ppart,ppbuff,rbufl,rbufr,kpic,ncl,ihole,& &mcll,mclr,idimp,nppmx,mx1,myp1,npbmx,ntmax,nbmax,irc) implicit none integer, intent(in) :: idimp, nppmx, mx1, myp1, npbmx, ntmax integer, intent(in) :: nbmax integer, intent(inout) :: irc real, dimension(idimp,nppmx,mx1*myp1), intent(inout) :: ppart real, dimension(idimp,npbmx,mx1*myp1), intent(in) :: ppbuff real, dimension(idimp,nbmax), intent(in) :: rbufl, rbufr integer, dimension(mx1*myp1), intent(inout) :: kpic integer, dimension(8,mx1*myp1), intent(in) :: ncl integer, dimension(2,ntmax+1,mx1*myp1), intent(in) :: ihole integer, dimension(3,mx1), intent(in) :: mcll, mclr end subroutine end interface ! interface subroutine WPGPSIPOST2L_QP(ppart,psi,kpic,qbm,noff,nyp,idimp,np& &pmx,nx,mx,my,nxv,nypmx,mx1,mxyp1,dex) implicit none integer, intent(in) :: noff, nyp, idimp, nppmx, nx, mx, my integer, intent(in) :: mx1, mxyp1, nxv, nypmx real, intent(in) :: dex,qbm real, dimension(idimp,nppmx,mxyp1), intent(inout) :: ppart real, dimension(nxv,nypmx), intent(in) :: psi integer, dimension(mxyp1), intent(in) :: kpic end subroutine end interface ! interface subroutine PPPCOPYOUT2(part,ppart,kpic,npp,npmax,nppmx,idimp, & &mxyp1,irc) implicit none integer, intent(in) :: npmax, nppmx, idimp, mxyp1 integer, intent(inout) :: npp, irc real, dimension(idimp,npmax), intent(inout) :: part real, dimension(idimp,nppmx,mxyp1), intent(in) :: ppart integer, dimension(mxyp1), intent(in) :: kpic end subroutine end interface ! contains ! subroutine PPPMOVE2(sbufr,sbufl,rbufr,rbufl,ncll,nclr,mcll,mclr, & &kstrt,nvp,idimp,nbmax,mx1) ! this subroutine moves particles into appropriate spatial regions ! for distributed data, with 1d domain decomposition in y. ! tiles are assumed to be arranged in 2D linear memory ! output: rbufr, rbufl, mcll, mclr ! sbufl = buffer for particles being sent to lower processor ! sbufr = buffer for particles being sent to upper processor ! rbufl = buffer for particles being received from lower processor ! rbufr = buffer for particles being received from upper processor ! ncll = particle number being sent to lower processor ! nclr = particle number being sent to upper processor ! mcll = particle number being received from lower processor ! mclr = particle number being received from upper processor ! kstrt = starting data block number ! nvp = number of real or virtual processors ! idimp = size of phase space = 4 or 5 ! nbmax = size of buffers for passing particles between processors ! mx1 = (system length in x direction - 1)/mx + 1 implicit none integer, intent(in) :: kstrt, nvp, idimp, nbmax, mx1 real, dimension(idimp,nbmax), intent(in) :: sbufl, sbufr real, dimension(idimp,nbmax), intent(inout) :: rbufl, rbufr integer, dimension(3,mx1), intent(in) :: ncll, nclr integer, dimension(3,mx1), intent(inout) :: mcll, mclr ! lgrp = current communicator ! mint = default datatype for integers ! mreal = default datatype for reals ! local data integer :: ierr, ks, kl, kr, i, j, jsl, jsr integer :: nbsize, ncsize integer, dimension(8) :: msid integer, dimension(4) :: itg integer, dimension(10) :: istatus data itg /3,4,5,6/ ks = kstrt - 1 nbsize = idimp*nbmax ncsize = 3*mx1 ! copy particle buffers: update rbufl, rbufr, mcll, mclr ! special case for one processor if (nvp==1) then do j = 1, mx1 do i = 1, 3 mcll(i,j) = nclr(i,j) enddo continue enddo do j = 1, mx1 do i = 1, 3 mclr(i,j) = ncll(i,j) enddo enddo do j = 1, nclr(3,mx1) do i = 1, idimp rbufl(i,j) = sbufr(i,j) enddo enddo do j = 1, ncll(3,mx1) do i = 1, idimp rbufr(i,j) = sbufl(i,j) enddo enddo ! this segment is used for mpi computers else ! get particles from below and above kr = ks + 1 if (kr >= nvp) kr = kr - nvp kl = ks - 1 if (kl < 0) kl = kl + nvp ! post receives call MPI_IRECV(mcll,ncsize,mint,kl,itg(1),lgrp,msid(1),ierr) call MPI_IRECV(mclr,ncsize,mint,kr,itg(2),lgrp,msid(2),ierr) call MPI_IRECV(rbufl,nbsize,mreal,kl,itg(3),lgrp,msid(3),ierr) call MPI_IRECV(rbufr,nbsize,mreal,kr,itg(4),lgrp,msid(4),ierr) ! send particle number offsets call MPI_ISEND(nclr,ncsize,mint,kr,itg(1),lgrp,msid(5),ierr) call MPI_ISEND(ncll,ncsize,mint,kl,itg(2),lgrp,msid(6),ierr) call MPI_WAIT(msid(1),istatus,ierr) call MPI_WAIT(msid(2),istatus,ierr) ! send particles jsr = idimp*nclr(3,mx1) call MPI_ISEND(sbufr,jsr,mreal,kr,itg(3),lgrp,msid(7),ierr) jsl = idimp*ncll(3,mx1) call MPI_ISEND(sbufl,jsl,mreal,kl,itg(4),lgrp,msid(8),ierr) call MPI_WAIT(msid(3),istatus,ierr) call MPI_WAIT(msid(4),istatus,ierr) endif ! make sure sbufr, sbufl, ncll, and nclr have been sent if (nvp /= 1) then do i = 1, 4 call MPI_WAIT(msid(i+4),istatus,ierr) enddo endif end subroutine PPPMOVE2 end module part2d_lib