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

Minor fixes for scalar transport. #18

Merged
merged 1 commit into from
Jan 12, 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
1 change: 1 addition & 0 deletions src/chkdt.f90
Original file line number Diff line number Diff line change
Expand Up @@ -69,6 +69,7 @@ subroutine chkdt(n,dl,dzci,dzfi,is_solve_ns,mu12,rho12,sigma,gacc,u,v,w,dtmax,ga
!$acc wait(1)
call MPI_ALLREDUCE(MPI_IN_PLACE,dti,1,MPI_REAL_RP,MPI_MAX,MPI_COMM_WORLD,ierr)
dtipsi = gam*seps/dlmin**2
if(dtipsi == 0.) dtipsi = 1.
if(.not.is_solve_ns) then
dtmax = min(dti**(-1),dtipsi**(-1))
else
Expand Down
2 changes: 1 addition & 1 deletion src/initflow.f90
Original file line number Diff line number Diff line change
Expand Up @@ -328,7 +328,7 @@ subroutine initscal(inisca,bcsca,ng,lo,l,dl,dzf,zc,s)
sref = 1.
sref = 0.5*(bcsca(0,3)+bcsca(1,3))
select case(trim(inisca))
case('cou')
case('cou','dhc')
call couette( n(3),zc/l(3),1._rp,s1d)
s1d(:) = s1d(:) + 0.5 ! from 0 to 1
s1d(:) = bcsca(0,3)*(1.-s1d(:)) + bcsca(1,3)*(s1d(:))
Expand Down
14 changes: 7 additions & 7 deletions src/main.f90
Original file line number Diff line number Diff line change
Expand Up @@ -118,7 +118,7 @@ program cans
integer :: savecounter
character(len=7 ) :: fldnum
character(len=4 ) :: chkptnum
character(len=100) :: filename,fexts(5)
character(len=100) :: filename,fexts(6)
integer :: k,kk
logical :: is_done,kill
real(rp), dimension(2) :: tm_coeff
Expand Down Expand Up @@ -308,9 +308,9 @@ program cans
!
write(fldnum,'(i7.7)') istep
!$acc update self(u,v,w,p,psi,kappa)
include 'out1d.h90'
include 'out2d.h90'
include 'out3d.h90'
#include "out1d.h90"
#include "out2d.h90"
#include "out3d.h90"
!
call chkdt(n,dl,dzci,dzfi,is_solve_ns,mu12,rho12,sigma,gacc,u,v,w,dtmax,gam,seps,ka12,cp12)
dt = min(cfl*dtmax,dtmin)
Expand Down Expand Up @@ -435,15 +435,15 @@ program cans
write(fldnum,'(i7.7)') istep
if(mod(istep,iout1d) == 0) then
!$acc update self(u,v,w,p,psi,kappa)
include 'out1d.h90'
#include "out1d.h90"
end if
if(mod(istep,iout2d) == 0) then
!$acc update self(u,v,w,p,psi,kappa)
include 'out2d.h90'
#include "out2d.h90"
end if
if(mod(istep,iout3d) == 0) then
!$acc update self(u,v,w,p,psi,kappa)
include 'out3d.h90'
#include "out3d.h90"
end if
if(mod(istep,isave ) == 0.or.(is_done.and..not.kill)) then
if(is_overwrite_save) then
Expand Down
5 changes: 5 additions & 0 deletions src/out3d.h90
Original file line number Diff line number Diff line change
Expand Up @@ -43,3 +43,8 @@
call write_visu_3d(datadir,'cur_fld_'//fldnum//'.bin','log_visu_3d.out','Kappa', &
[1,1,1],[ng(1),ng(2),ng(3)],[1,1,1],time,istep, &
kappa(1:n(1),1:n(2),1:n(3)))
#if defined(_SCALAR)
call write_visu_3d(datadir,'sca_fld_'//fldnum//'.bin','log_visu_3d.out','Temperature', &
[1,1,1],[ng(1),ng(2),ng(3)],[1,1,1],time,istep, &
s(1:n(1),1:n(2),1:n(3)))
#endif
1 change: 1 addition & 0 deletions src/rk.f90
Original file line number Diff line number Diff line change
Expand Up @@ -144,6 +144,7 @@ subroutine rk_scal(rkpar,n,dli,dzci,dzfi,dt, &
factor1 = rkpar(1)*dt
factor2 = rkpar(2)*dt
factor12 = factor1 + factor2
rhocp12 = rho12(:)*cp12(:)
if(is_first) then ! leverage save attribute to allocate these arrays on the device only once
is_first = .false.
allocate(dsdtrk_t(n(1),n(2),n(3)),dsdtrko_t(n(1),n(2),n(3)))
Expand Down
Loading