plasma.f90 Source File


This file depends on

sourcefile~~plasma.f90~~EfferentGraph sourcefile~plasma.f90 plasma.f90 sourcefile~io.f90 io.f90 sourcefile~plasma.f90->sourcefile~io.f90 sourcefile~reader.f90 reader.f90 sourcefile~plasma.f90->sourcefile~reader.f90 sourcefile~grid.f90 grid.f90 sourcefile~io.f90->sourcefile~grid.f90 sourcefile~mpimod.f90 mpimod.F90 sourcefile~io.f90->sourcefile~mpimod.f90 sourcefile~phys_consts.f90 phys_consts.F90 sourcefile~io.f90->sourcefile~phys_consts.f90 sourcefile~pathlib.f90 pathlib.f90 sourcefile~io.f90->sourcefile~pathlib.f90 sourcefile~reader.f90->sourcefile~phys_consts.f90 sourcefile~grid.f90->sourcefile~reader.f90 sourcefile~grid.f90->sourcefile~mpimod.f90 sourcefile~grid.f90->sourcefile~phys_consts.f90 sourcefile~mesh.f90 mesh.f90 sourcefile~grid.f90->sourcefile~mesh.f90 sourcefile~mpimod.f90->sourcefile~phys_consts.f90 sourcefile~mesh.f90->sourcefile~phys_consts.f90

Files dependent on this one

sourcefile~~plasma.f90~~AfferentGraph sourcefile~plasma.f90 plasma.f90 sourcefile~plasma_output_raw.f90 plasma_output_raw.f90 sourcefile~plasma_output_raw.f90->sourcefile~plasma.f90 sourcefile~plasma_output_ncdf.f90 plasma_output_ncdf.f90 sourcefile~plasma_output_ncdf.f90->sourcefile~plasma.f90 sourcefile~plasma_input_hdf5.f90 plasma_input_hdf5.f90 sourcefile~plasma_input_hdf5.f90->sourcefile~plasma.f90 sourcefile~plasma_input_raw.f90 plasma_input_raw.f90 sourcefile~plasma_input_raw.f90->sourcefile~plasma.f90 sourcefile~plasma_output_hdf5.f90 plasma_output_hdf5.f90 sourcefile~plasma_output_hdf5.f90->sourcefile~plasma.f90

Contents

Source Code


Source Code

submodule (io) plasma
!! plasma.f90 uses submodules in plasma_input_*.f90 and plasma_output_*.f90 for raw, hdf5 or netcdf4 I/O

use reader, only : get_simsize3
implicit none

interface ! plasma_input_*.f90

module subroutine input_root_currents(outdir,flagoutput,ymd,UTsec,J1,J2,J3)
character(*), intent(in) :: outdir
integer, intent(in) :: flagoutput
integer, dimension(3), intent(in) :: ymd
real(wp), intent(in) :: UTsec
real(wp), dimension(:,:,:), intent(out) :: J1,J2,J3
end subroutine input_root_currents


module subroutine input_root_mpi(x1,x2all,x3all,indatsize,ns,vs1,Ts)
real(wp), dimension(-1:), intent(in) :: x1, x2all, x3all
character(*), intent(in) :: indatsize
real(wp), dimension(-1:,-1:,-1:,:), intent(out) :: ns,vs1,Ts
end subroutine input_root_mpi

end interface


interface ! plasma_output_*.f90

module subroutine output_root_stream_mpi(outdir,flagoutput,ymd,UTsec,vs2,vs3,ns,vs1,Ts,Phiall,J1,J2,J3)
character(*), intent(in) :: outdir
integer, intent(in) :: flagoutput

integer, dimension(3), intent(in) :: ymd
real(wp), intent(in) :: UTsec
real(wp), dimension(-1:,-1:,-1:,:), intent(in) :: vs2,vs3,ns,vs1,Ts

real(wp), dimension(:,:,:), intent(in) :: Phiall
real(wp), dimension(:,:,:), intent(in) :: J1,J2,J3
end subroutine output_root_stream_mpi

end interface

contains

module procedure input_plasma
! subroutine input_plasma(x1,x2,x3all,indatsize,ns,vs1,Ts)
!! A BASIC WRAPPER FOR THE ROOT AND WORKER INPUT FUNCTIONS
!! BOTH ROOT AND WORKERS CALL THIS PROCEDURE SO UNALLOCATED
!! VARIABLES MUST BE DECLARED AS ALLOCATABLE, INTENT(INOUT)

if (myid==0) then
  !ROOT FINDS/CALCULATES INITIAL CONDITIONS AND SENDS TO WORKERS
  ! print '(A,/,A,/,A)', 'Assembling initial condition on root with:',indatsize,indatfile
  call input_root_mpi(x1,x2,x3all,indatsize,ns,vs1,Ts)
else
  !WORKERS RECEIVE THE IC DATA FROM ROOT
  call input_workers_mpi(ns,vs1,Ts)
end if

end procedure input_plasma


module procedure input_plasma_currents
! module subroutine input_plasma_currents(outdir,flagoutput,ymd,UTsec,J1,J2,J3)
!! READS, AS INPUT, A FILE GENERATED BY THE GEMINI.F90 PROGRAM.
!! THIS SUBROUTINE IS A WRAPPER FOR SEPARATE ROOT/WORKER CALLS

if (myid==0) then
  !> ROOT FINDS/CALCULATES INITIAL CONDITIONS AND SENDS TO WORKERS
  print *, 'Assembling current density data on root...  '
  call input_root_currents(outdir,flagoutput,ymd,UTsec,J1,J2,J3)
else
  !> WORKERS RECEIVE THE IC DATA FROM ROOT
  call input_workers_currents(J1,J2,J3)
end if

end procedure input_plasma_currents


subroutine input_workers_currents(J1,J2,J3)
!! WORKER INPUT FUNCTIONS FOR GETTING CURRENT DENSITIES

real(wp), dimension(:,:,:), intent(out) :: J1,J2,J3


!> ALL WE HAVE TO DO IS WAIT TO RECEIVE OUR PIECE OF DATA FROM ROOT
call bcast_recv(J1,tagJ1)
call bcast_recv(J2,tagJ2)
call bcast_recv(J3,tagJ3)

end subroutine input_workers_currents

subroutine input_workers_mpi(ns,vs1,Ts)

!------------------------------------------------------------
!-------RECEIVE INITIAL CONDITIONS FROM ROOT PROCESS
!------------------------------------------------------------

real(wp), dimension(-1:,-1:,-1:,:), intent(out) :: ns,vs1,Ts

call bcast_recv(ns,tagns)
call bcast_recv(vs1,tagvs1)
call bcast_recv(Ts,tagTs)


if (.false.) then
  print*, myid
  print *, 'Min/max input density:  ',     minval(ns(:,:,:,7)),  maxval(ns(:,:,:,7))
  print *, 'Min/max input velocity:  ',    minval(vs1(:,:,:,:)), maxval(vs1(:,:,:,:))
  print *, 'Min/max input temperature:  ', minval(Ts(:,:,:,:)),  maxval(Ts(:,:,:,:))
endif

end subroutine input_workers_mpi


subroutine output_workers_mpi(vs2,vs3,ns,vs1,Ts,J1,J2,J3)

!------------------------------------------------------------
!-------SEND COMPLETE DATA FROM WORKERS TO ROOT PROCESS FOR OUTPUT.
!-------STATE VARS ARE EXPECTED TO INCLUDE GHOST CELLS
!------------------------------------------------------------

real(wp), dimension(-1:,-1:,-1:,:), intent(in) :: vs2,vs3,ns,vs1,Ts
real(wp), dimension(:,:,:), intent(in) :: J1,J2,J3

integer :: lx1,lx2,lx3,lx3all,isp
real(wp), dimension(1:size(ns,1)-4,1:size(ns,2)-4,1:size(ns,3)-4) :: v2avg,v3avg


!SYSTEM SIZES (W/O GHOST CELLS)
lx1=size(ns,1)-4
lx2=size(ns,2)-4
lx3=size(ns,3)-4


!ONLY AVERAGE DRIFTS PERP TO B NEEDED FOR OUTPUT
v2avg=sum(ns(1:lx1,1:lx2,1:lx3,1:lsp-1)*vs2(1:lx1,1:lx2,1:lx3,1:lsp-1),4)
v2avg=v2avg/ns(1:lx1,1:lx2,1:lx3,lsp)    !compute averages for output.
v3avg=sum(ns(1:lx1,1:lx2,1:lx3,1:lsp-1)*vs3(1:lx1,1:lx2,1:lx3,1:lsp-1),4)
v3avg=v3avg/ns(1:lx1,1:lx2,1:lx3,lsp)


!SEND MY GRID DATA TO THE ROOT PROCESS
call gather_send(v2avg,tagv2)
call gather_send(v3avg,tagv3)
call gather_send(ns,tagns)
call gather_send(vs1,tagvs1)
call gather_send(Ts,tagTs)


!------- SEND ELECTRODYNAMIC PARAMETERS TO ROOT
call gather_send(J1,tagJ1)
call gather_send(J2,tagJ2)
call gather_send(J3,tagJ3)

end subroutine output_workers_mpi


module procedure output_plasma
! subroutine output_plasma(outdir,flagoutput,ymd,UTsec,vs2,vs3,ns,vs1,Ts,Phiall,J1,J2,J3)
!! A BASIC WRAPPER FOR THE ROOT AND WORKER OUTPUT FUNCTIONS
!! BOTH ROOT AND WORKERS CALL THIS PROCEDURE SO UNALLOCATED
!! VARIABLES MUST BE DECLARED AS ALLOCATABLE, INTENT(INOUT)

if (myid/=0) then
  call output_workers_mpi(vs2,vs3,ns,vs1,Ts,J1,J2,J3)
else
  call output_root_stream_mpi(outdir,flagoutput,ymd,UTsec,vs2,vs3,ns,vs1,Ts,Phiall,J1,J2,J3)    !only option that works with >2GB files
end if

end procedure output_plasma

end submodule plasma