Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Use newunit in simple I/O tasks in place of hard coded values #48

Merged
merged 1 commit into from
Apr 30, 2024
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
66 changes: 33 additions & 33 deletions src/CalcWriteQ.F90
Original file line number Diff line number Diff line change
Expand Up @@ -20,7 +20,7 @@ subroutine CalcWriteQ
real :: dvzx1,dvzx2,dvzx3
real :: strn, omeg, tprfi

integer :: ic,jc,kc,ip,jp,kp,im,jm,km,itime
integer :: ic,jc,kc,ip,jp,kp,im,jm,km,itime,io

character(30) :: filnam,filnamxdm,path
character(5) :: ipfi
Expand Down Expand Up @@ -98,41 +98,41 @@ subroutine CalcWriteQ
! if(allocated(qtens)) deallocate(qtens)

if (ismaster) then
open(45,file=filnamxdm,status='unknown')
rewind(45)
write(45,'("<?xml version=""1.0"" ?>")')
write(45,'("<!DOCTYPE Xdmf SYSTEM ""Xdmf.dtd"" []>")')
write(45,'("<Xdmf Version=""2.0"">")')
write(45,'("<Domain>")')
write(45,'("<Grid Name=""basegrid"" GridType=""Uniform"">")')

write(45,'("<Topology TopologyType=""3DRectMesh"" NumberOfElements=""",i4," ",i4," ",i4,"""/>")') nzm,nym,nx

write(45,'("<Geometry GeometryType=""VXVYVZ"">")')
write(45,'("<DataItem Dimensions=""",i4,""" NumberType=""Float"" Precision=""4"" Format=""HDF"">")')nx
write(45,'("cordin_info.h5:/xc")')
write(45,'("</DataItem>")')
write(45,'("<DataItem Dimensions=""",i4,""" NumberType=""Float"" Precision=""4"" Format=""HDF"">")')nym
write(45,'("cordin_info.h5:/ym")')
write(45,'("</DataItem>")')
write(45,'("<DataItem Dimensions=""",i4,""" NumberType=""Float"" Precision=""4"" Format=""HDF"">")')nzm
write(45,'("cordin_info.h5:/zm")')
write(45,'("</DataItem>")')
write(45,'("</Geometry>")')

write(45,'("<Attribute Name=""Q"" AttributeType=""Scalar"" Center=""Node"">")')
write(45,'("<DataItem Dimensions=""",i4," ",i4," ",i4,""" NumberType=""Double"" Precision=""2"" Format=""HDF"">")')&
open(newunit=io, file=filnamxdm, status='unknown')
rewind(io)
write(io,'("<?xml version=""1.0"" ?>")')
write(io,'("<!DOCTYPE Xdmf SYSTEM ""Xdmf.dtd"" []>")')
write(io,'("<Xdmf Version=""2.0"">")')
write(io,'("<Domain>")')
write(io,'("<Grid Name=""basegrid"" GridType=""Uniform"">")')

write(io,'("<Topology TopologyType=""3DRectMesh"" NumberOfElements=""",i4," ",i4," ",i4,"""/>")') nzm,nym,nx

write(io,'("<Geometry GeometryType=""VXVYVZ"">")')
write(io,'("<DataItem Dimensions=""",i4,""" NumberType=""Float"" Precision=""4"" Format=""HDF"">")')nx
write(io,'("cordin_info.h5:/xc")')
write(io,'("</DataItem>")')
write(io,'("<DataItem Dimensions=""",i4,""" NumberType=""Float"" Precision=""4"" Format=""HDF"">")')nym
write(io,'("cordin_info.h5:/ym")')
write(io,'("</DataItem>")')
write(io,'("<DataItem Dimensions=""",i4,""" NumberType=""Float"" Precision=""4"" Format=""HDF"">")')nzm
write(io,'("cordin_info.h5:/zm")')
write(io,'("</DataItem>")')
write(io,'("</Geometry>")')

write(io,'("<Attribute Name=""Q"" AttributeType=""Scalar"" Center=""Node"">")')
write(io,'("<DataItem Dimensions=""",i4," ",i4," ",i4,""" NumberType=""Double"" Precision=""2"" Format=""HDF"">")')&
nzm,nym,nx
write(45,*) trim(filnam)//':/var'
write(45,'("</DataItem>")')
write(45,'("</Attribute>")')
write(io,*) trim(filnam)//':/var'
write(io,'("</DataItem>")')
write(io,'("</Attribute>")')

write(45,'("<Time Value=""",e12.5,""" />")') time
write(io,'("<Time Value=""",e12.5,""" />")') time

write(45,'("</Grid>")')
write(45,'("</Domain>")')
write(45,'("</Xdmf>")')
close(45)
write(io,'("</Grid>")')
write(io,'("</Domain>")')
write(io,'("</Xdmf>")')
close(io)
end if

end subroutine CalcWriteQ
Expand Down
14 changes: 7 additions & 7 deletions src/GlobalQuantities.F90
Original file line number Diff line number Diff line change
Expand Up @@ -14,7 +14,7 @@ subroutine GlobalQuantities
use decomp_2d, only: xstart,xend
use mpih
implicit none
integer :: jc,kc,kp,ic
integer :: jc,kc,kp,ic,io
real :: anusin,vol,vxcen,fac2,tempcen
real :: vy_rms_vol,vz_rms_vol
real :: vx_rms_vol,vzvyvx_rms_vol,rradpr
Expand Down Expand Up @@ -78,10 +78,10 @@ subroutine GlobalQuantities

anusin=1.d0 + dsqrt(prat*abs(rayt))*anusin*vol

open(95,file='outputdir/nu_vol.out',status='unknown',access='sequential', &
open(newunit=io,file='outputdir/nu_vol.out',status='unknown',access='sequential', &
position='append')
write(95,*) time, anusin
close(95)
write(io,*) time, anusin
close(io)

rradpr=dsqrt(abs(rayt)/prat)
tempm=tempm*vol
Expand All @@ -90,11 +90,11 @@ subroutine GlobalQuantities
vz_rms_vol=dsqrt(vz_rms_vol*vol)*rradpr
vzvyvx_rms_vol=dsqrt(vzvyvx_rms_vol*vol)*rradpr

open(94,file='outputdir/rms_vel.out',status='unknown',position='append', &
open(newunit=io,file='outputdir/rms_vel.out',status='unknown',position='append', &
access='sequential')
write(94,*) time,vz_rms_vol,vy_rms_vol,vx_rms_vol, &
write(io,*) time,vz_rms_vol,vy_rms_vol,vx_rms_vol, &
& vzvyvx_rms_vol
close(94)
close(io)

endif

Expand Down
96 changes: 48 additions & 48 deletions src/flow_solver/2decomp/decomp_2d.F90
Original file line number Diff line number Diff line change
Expand Up @@ -279,7 +279,7 @@ subroutine decomp_2d_init(nx,ny,nz,nxr,nyr,nzr,p_row,p_col,periodic_bc)
integer, intent(INOUT) :: p_row,p_col
logical, dimension(3), intent(IN), optional :: periodic_bc

integer :: errorcode, ierror, row, col
integer :: errorcode, ierror, row, col, io

#ifdef SHM_DEBUG
character(len=80) fname
Expand Down Expand Up @@ -389,55 +389,55 @@ subroutine decomp_2d_init(nx,ny,nz,nxr,nyr,nzr,p_row,p_col,periodic_bc)
! print out shared-memory information
write(fname,99) nrank
99 format('log',I2.2)
open(10,file=fname)
write(10,*)'I am mpi rank ', nrank, 'Total ranks ', nproc
write(10,*)' '
write(10,*)'Global data size:'
write(10,*)'nx*ny*nz', nx,ny,nz
write(10,*)' '
write(10,*)'2D processor grid:'
write(10,*)'p_row*p_col:', dims(1), dims(2)
write(10,*)' '
write(10,*)'Portion of global data held locally:'
write(10,*)'xsize:',xsize
write(10,*)'ysize:',ysize
write(10,*)'zsize:',zsize
write(10,*)' '
write(10,*)'How pensils are to be divided and sent in alltoallv:'
write(10,*)'x1dist:',decomp_main%x1dist
write(10,*)'y1dist:',decomp_main%y1dist
write(10,*)'y2dist:',decomp_main%y2dist
write(10,*)'z2dist:',decomp_main%z2dist
write(10,*)' '
write(10,*)'######Shared buffer set up after this point######'
write(10,*)' '
write(10,*) 'col communicator detais:'
open(newunit=io,file=fname)
write(io,*)'I am mpi rank ', nrank, 'Total ranks ', nproc
write(io,*)' '
write(io,*)'Global data size:'
write(io,*)'nx*ny*nz', nx,ny,nz
write(io,*)' '
write(io,*)'2D processor grid:'
write(io,*)'p_row*p_col:', dims(1), dims(2)
write(io,*)' '
write(io,*)'Portion of global data held locally:'
write(io,*)'xsize:',xsize
write(io,*)'ysize:',ysize
write(io,*)'zsize:',zsize
write(io,*)' '
write(io,*)'How pensils are to be divided and sent in alltoallv:'
write(io,*)'x1dist:',decomp_main%x1dist
write(io,*)'y1dist:',decomp_main%y1dist
write(io,*)'y2dist:',decomp_main%y2dist
write(io,*)'z2dist:',decomp_main%z2dist
write(io,*)' '
write(io,*)'######Shared buffer set up after this point######'
write(io,*)' '
write(io,*) 'col communicator detais:'
call print_smp_info(decomp_main%COL_INFO)
write(10,*)' '
write(10,*) 'row communicator detais:'
write(io,*)' '
write(io,*) 'row communicator detais:'
call print_smp_info(decomp_main%ROW_INFO)
write(10,*)' '
write(10,*)'Buffer count and dispalcement of per-core buffers'
write(10,*)'x1cnts:',decomp_main%x1cnts
write(10,*)'y1cnts:',decomp_main%y1cnts
write(10,*)'y2cnts:',decomp_main%y2cnts
write(10,*)'z2cnts:',decomp_main%z2cnts
write(10,*)'x1disp:',decomp_main%x1disp
write(10,*)'y1disp:',decomp_main%y1disp
write(10,*)'y2disp:',decomp_main%y2disp
write(10,*)'z2disp:',decomp_main%z2disp
write(10,*)' '
write(10,*)'Buffer count and dispalcement of shared buffers'
write(10,*)'x1cnts:',decomp_main%x1cnts_s
write(10,*)'y1cnts:',decomp_main%y1cnts_s
write(10,*)'y2cnts:',decomp_main%y2cnts_s
write(10,*)'z2cnts:',decomp_main%z2cnts_s
write(10,*)'x1disp:',decomp_main%x1disp_s
write(10,*)'y1disp:',decomp_main%y1disp_s
write(10,*)'y2disp:',decomp_main%y2disp_s
write(10,*)'z2disp:',decomp_main%z2disp_s
write(10,*)' '
close(10)
write(io,*)' '
write(io,*)'Buffer count and dispalcement of per-core buffers'
write(io,*)'x1cnts:',decomp_main%x1cnts
write(io,*)'y1cnts:',decomp_main%y1cnts
write(io,*)'y2cnts:',decomp_main%y2cnts
write(io,*)'z2cnts:',decomp_main%z2cnts
write(io,*)'x1disp:',decomp_main%x1disp
write(io,*)'y1disp:',decomp_main%y1disp
write(io,*)'y2disp:',decomp_main%y2disp
write(io,*)'z2disp:',decomp_main%z2disp
write(io,*)' '
write(io,*)'Buffer count and dispalcement of shared buffers'
write(io,*)'x1cnts:',decomp_main%x1cnts_s
write(io,*)'y1cnts:',decomp_main%y1cnts_s
write(io,*)'y2cnts:',decomp_main%y2cnts_s
write(io,*)'z2cnts:',decomp_main%z2cnts_s
write(io,*)'x1disp:',decomp_main%x1disp_s
write(io,*)'y1disp:',decomp_main%y1disp_s
write(io,*)'y2disp:',decomp_main%y2disp_s
write(io,*)'z2disp:',decomp_main%z2disp_s
write(io,*)' '
close(io)
#endif

! determine the number of bytes per float number
Expand Down
16 changes: 8 additions & 8 deletions src/flow_solver/CreateGrid.F90
Original file line number Diff line number Diff line change
Expand Up @@ -14,7 +14,7 @@ subroutine CreateGrid
use GridModule
implicit none

integer :: kc
integer :: kc, io

do kc=1,nxm
kmv(kc)=kc-1
Expand Down Expand Up @@ -138,22 +138,22 @@ subroutine CreateGrid
udx3c(nx) = dx/g3rc(nx)
!m====================================================
if(ismaster) then
open(unit=78,file='outputdir/axicor.out',status='unknown')
open(newunit=io,file='outputdir/axicor.out',status='unknown')
do kc=1,nx
write(78,345) kc,xc(kc),xm(kc),g3rc(kc),g3rm(kc)
write(io,345) kc,xc(kc),xm(kc),g3rc(kc),g3rm(kc)
end do
close(78)
close(io)
345 format(i4,4(2x,e23.15))
!m===================================================
!
! QUANTITIES FOR DERIVATIVES
!
open(unit=78,file='outputdir/fact3.out',status='unknown')
open(newunit=io,file='outputdir/fact3.out',status='unknown')
do kc=1,nxm
write(78,*) kc,udx3m(kc),udx3c(kc)
write(io,*) kc,udx3m(kc),udx3c(kc)
end do
write(78,*) nx,udx3m(nxm),udx3c(nx)
close(78)
write(io,*) nx,udx3m(nxm),udx3c(nx)
close(io)
end if

!
Expand Down
116 changes: 58 additions & 58 deletions src/flow_solver/ReadInputFile.F90
Original file line number Diff line number Diff line change
Expand Up @@ -13,68 +13,68 @@ subroutine ReadInputFile
use afid_phasefield, only: pf_A, pf_D, pf_eps, pf_Lambda, pf_S, pf_Tm
implicit none
character(len=4) :: dummy
integer :: flagmelt
integer :: flagmelt, io
integer :: flagMR, flagsal, flagPF
integer :: FFscaleS

open(unit=15,file='bou.in',status='old')
read(15,301) dummy
read(15,301) dummy
read(15,*) nxm, nym, nzm, nsst
read(15,301) dummy
read(15,301) dummy
read(15,301) dummy
read(15,*) flagMR, nxmr, nymr, nzmr
read(15,301) dummy
read(15,301) dummy
read(15,301) dummy
read(15,*) flagsal, flagPF
read(15,301) dummy
read(15,301) dummy
read(15,301) dummy
read(15,*) nread, ireset
read(15,301) dummy
read(15,301) dummy
read(15,301) dummy
read(15,*) ntst, walltimemax, tmax
read(15,301) dummy
read(15,301) dummy
read(15,301) dummy
read(15,*) tout, tframe, save_3D
read(15,301) dummy
read(15,301) dummy
read(15,301) dummy
read(15,*) alx3, ylen, zlen
read(15,301) dummy
read(15,301) dummy
read(15,301) dummy
read(15,*) istr3, str3, istr3r
read(15,301) dummy
read(15,301) dummy
read(15,301) dummy
read(15,*) rayt, prat, rays, pras, FFscaleS
read(15,301) dummy
read(15,301) dummy
read(15,301) dummy
read(15,*) idtv, dt, resid, limitCFL, dtmin, dtmax
read(15,301) dummy
read(15,301) dummy
read(15,301) dummy
read(15,*) inslws, inslwn, TfixS, TfixN, SfixS, SfixN
read(15,301) dummy
read(15,301) dummy
read(15,301) dummy
read(15,*) active_T, active_S, gAxis
read(15,301) dummy
read(15,301) dummy
read(15,301) dummy
read(15,*) xplusU, xminusU, dPdy, dPdz, flagmelt
read(15,301) dummy
read(15,301) dummy
read(15,301) dummy
read(15,*) pf_D, pf_A, pf_S, pf_Tm, solidtype, pf_IC
open(newunit=io,file='bou.in',status='old')
read(io,301) dummy
read(io,301) dummy
read(io,*) nxm, nym, nzm, nsst
read(io,301) dummy
read(io,301) dummy
read(io,301) dummy
read(io,*) flagMR, nxmr, nymr, nzmr
read(io,301) dummy
read(io,301) dummy
read(io,301) dummy
read(io,*) flagsal, flagPF
read(io,301) dummy
read(io,301) dummy
read(io,301) dummy
read(io,*) nread, ireset
read(io,301) dummy
read(io,301) dummy
read(io,301) dummy
read(io,*) ntst, walltimemax, tmax
read(io,301) dummy
read(io,301) dummy
read(io,301) dummy
read(io,*) tout, tframe, save_3D
read(io,301) dummy
read(io,301) dummy
read(io,301) dummy
read(io,*) alx3, ylen, zlen
read(io,301) dummy
read(io,301) dummy
read(io,301) dummy
read(io,*) istr3, str3, istr3r
read(io,301) dummy
read(io,301) dummy
read(io,301) dummy
read(io,*) rayt, prat, rays, pras, FFscaleS
read(io,301) dummy
read(io,301) dummy
read(io,301) dummy
read(io,*) idtv, dt, resid, limitCFL, dtmin, dtmax
read(io,301) dummy
read(io,301) dummy
read(io,301) dummy
read(io,*) inslws, inslwn, TfixS, TfixN, SfixS, SfixN
read(io,301) dummy
read(io,301) dummy
read(io,301) dummy
read(io,*) active_T, active_S, gAxis
read(io,301) dummy
read(io,301) dummy
read(io,301) dummy
read(io,*) xplusU, xminusU, dPdy, dPdz, flagmelt
read(io,301) dummy
read(io,301) dummy
read(io,301) dummy
read(io,*) pf_D, pf_A, pf_S, pf_Tm, solidtype, pf_IC
301 format(a4)
close(15)
close(io)

nx=nxm+1
ny=nym+1
Expand Down
Loading