Skip to content

Commit db4cf92

Browse files
committed
Add comment and change flag name
1 parent 0a74b34 commit db4cf92

File tree

1 file changed

+15
-12
lines changed

1 file changed

+15
-12
lines changed

src/Marine.f90

+15-12
Original file line numberDiff line numberDiff line change
@@ -10,17 +10,17 @@ subroutine Marine()
1010

1111
double precision, dimension(:), allocatable :: flux,shelfdepth,ht,Fs,dh,dh1,dh2,Fmixt,mwater
1212
double precision, dimension(:), allocatable :: dhs, dhs1, F1, F2, zi, zo
13-
integer, dimension(:), allocatable :: flag,mmnrec,mmstack
13+
integer, dimension(:), allocatable :: COTflag,mmnrec,mmstack
1414
integer, dimension(:,:), allocatable :: mmrec
1515
double precision, dimension(:,:), allocatable :: mmwrec,mmlrec
1616
double precision shelfslope,ratio1,ratio2,dx,dy
1717
integer ij,ijr,ijk,k
1818

19-
allocate (flux(nn),shelfdepth(nn),ht(nn),Fs(nn),dh(nn),dh1(nn),dh2(nn),Fmixt(nn),flag(nn))
19+
allocate (flux(nn),shelfdepth(nn),ht(nn),Fs(nn),dh(nn),dh1(nn),dh2(nn),Fmixt(nn),COTflag(nn))
2020
allocate (dhs(nn),dhs1(nn),F1(nn),F2(nn),zi(nn),zo(nn))
2121

2222
! set nodes at transition between ocean and continent
23-
flag=0
23+
COTflag=0
2424

2525
dx=xl/(nx-1)
2626
dy=yl/(ny-1)
@@ -40,10 +40,13 @@ subroutine Marine()
4040
where (h.gt.sealevel) flux=0.d0
4141

4242
! set nodes at transition between ocean and continent
43-
!where (flux.gt.tiny(flux)) flag=1
43+
!where (flux.gt.tiny(flux)) COTflag=1
44+
! use single flow direction to set the flag that marks
45+
! the recieving node below/at sealevel of a node above/at
46+
! sealevel as a continent-ocean transition node.
4447
do ij=1,nn
4548
ijr=rec(ij)
46-
if (h(ij).ge.sealevel.and.h(ijr).le.sealevel) flag(ijr)=1
49+
if (h(ij).ge.sealevel.and.h(ijr).le.sealevel) COTflag(ijr)=1
4750
enddo
4851

4952
! decompact volume of pure solid phase (silt and sand) from onshore
@@ -114,7 +117,7 @@ subroutine Marine()
114117

115118
! silt and sand coupling diffusion in ocean
116119
call SiltSandCouplingDiffusion (h,Fmix,flux*Fs,flux*(1.d0-Fs), &
117-
nx,ny,dx,dy,dt,sealevel,layer,kdsea1,kdsea2,nGSMarine,flag,bounds_ibc)
120+
nx,ny,dx,dy,dt,sealevel,layer,kdsea1,kdsea2,nGSMarine,COTflag,bounds_ibc)
118121

119122
! pure silt and sand during deposition/erosion
120123
dh1=((h-ht)*Fmix+layer*(Fmix-Fmixt))*(1.d0-poro1)
@@ -173,7 +176,7 @@ end subroutine Marine
173176
!----------------------------------------------------------------------------------
174177

175178
subroutine SiltSandCouplingDiffusion (h,f,Q1,Q2,nx,ny,dx,dy,dt, &
176-
sealevel,L,kdsea1,kdsea2,niter,flag,ibc)
179+
sealevel,L,kdsea1,kdsea2,niter,COTflag,ibc)
177180

178181
implicit none
179182

@@ -182,7 +185,7 @@ subroutine SiltSandCouplingDiffusion (h,f,Q1,Q2,nx,ny,dx,dy,dt, &
182185
double precision h(nx*ny),f(nx*ny),Q1(nx*ny),Q2(nx*ny)
183186
double precision, dimension(:), allocatable :: hp,fp,ht,ft,hhalf,fhalf,fhalfp
184187
double precision, dimension(:), allocatable :: diag,sup,inf,rhs,res,tint
185-
integer flag(nx*ny)
188+
integer COTflag(nx*ny)
186189

187190
double precision dx,dy,dt,sealevel,L,kdsea1,kdsea2
188191
double precision K1,K2,tol,err1,err2
@@ -228,7 +231,7 @@ subroutine SiltSandCouplingDiffusion (h,f,Q1,Q2,nx,ny,dx,dy,dt, &
228231
ijp=(j)*nx+i
229232
ijm=(j-2)*nx+i
230233
! in ocean and not at ocean-continent transition
231-
if (ht(ij).le.sealevel.and.flag(ij).eq.0) then
234+
if (ht(ij).le.sealevel.and.COTflag(ij).eq.0) then
232235
if (i.eq.1) then
233236
if (cbc(4:4).eq.'1') then
234237
diag(i)=1.d0
@@ -302,7 +305,7 @@ subroutine SiltSandCouplingDiffusion (h,f,Q1,Q2,nx,ny,dx,dy,dt, &
302305
ijp=(j)*nx+i
303306
ijm=(j-2)*nx+i
304307
! in ocean and not at ocean-continent transition
305-
if (ht(ij).le.sealevel.and.flag(ij).eq.0) then
308+
if (ht(ij).le.sealevel.and.COTflag(ij).eq.0) then
306309
! deposition
307310
if (hhalf(ij).ge.(1.d0+1.d-6)*ht(ij)) then
308311
Dp=(hhalf(ij)-ht(ij))/dt
@@ -359,7 +362,7 @@ subroutine SiltSandCouplingDiffusion (h,f,Q1,Q2,nx,ny,dx,dy,dt, &
359362
ijp=(j)*nx+i
360363
ijm=(j-2)*nx+i
361364
! in ocean and not at ocean-continent transition
362-
if (ht(ij).le.sealevel.and.flag(ij).eq.0) then
365+
if (ht(ij).le.sealevel.and.COTflag(ij).eq.0) then
363366
if (j.eq.1) then
364367
if (cbc(1:1).eq.'1') then
365368
diag(j)=1.d0
@@ -433,7 +436,7 @@ subroutine SiltSandCouplingDiffusion (h,f,Q1,Q2,nx,ny,dx,dy,dt, &
433436
ijp=(j)*nx+i
434437
ijm=(j-2)*nx+i
435438
! in ocean and not at ocean-continent transition
436-
if (ht(ij).le.sealevel.and.flag(ij).eq.0) then
439+
if (ht(ij).le.sealevel.and.COTflag(ij).eq.0) then
437440
! deposition
438441
if (h(ij).ge.(1.d0+1.d-6)*hhalf(ij)) then
439442
Dp=(h(ij)-hhalf(ij))/dt

0 commit comments

Comments
 (0)