ROMS
Loading...
Searching...
No Matches
sed_fluxes.F
Go to the documentation of this file.
1#include "cppdefs.h"
2
4
5#if defined NONLINEAR && defined SEDIMENT && defined SUSPLOAD
6!
7!git $Id$
8!==================================================== John C. Warner ===
9! Copyright (c) 2002-2025 The ROMS Group Hernan G. Arango !
10! Licensed under a MIT/X style license !
11! See License_ROMS.md !
12!=======================================================================
13! !
14! This computes sediment bed and water column exchanges: deposition, !
15! resuspension, and erosion. !
16! !
17! References: !
18! !
19! Warner, J.C., C.R. Sherwood, R.P. Signell, C.K. Harris, and H.G. !
20! Arango, 2008: Development of a three-dimensional, regional, !
21! coupled wave, current, and sediment-transport model, Computers !
22! & Geosciences, 34, 1284-1306. !
23! !
24!=======================================================================
25!
26 implicit none
27!
28 PRIVATE
29 PUBLIC :: sed_fluxes
30!
31 CONTAINS
32!
33!***********************************************************************
34 SUBROUTINE sed_fluxes (ng, tile)
35!***********************************************************************
36!
37 USE mod_param
38 USE mod_forces
39 USE mod_grid
40 USE mod_ocean
41 USE mod_sedbed
42 USE mod_stepping
43 USE mod_bbl
44!
45! Imported variable declarations.
46!
47 integer, intent(in) :: ng, tile
48!
49! Local variable declarations.
50!
51 character (len=*), parameter :: myfile = &
52 & __FILE__
53!
54# include "tile.h"
55!
56# ifdef PROFILE
57 CALL wclock_on (ng, inlm, 16, __line__, myfile)
58# endif
59 CALL sed_fluxes_tile (ng, tile, &
60 & lbi, ubi, lbj, ubj, &
61 & imins, imaxs, jmins, jmaxs, &
62 & nstp(ng), nnew(ng), &
63 & grid(ng) % Hz, &
64# ifdef WET_DRY
65 & grid(ng) % rmask_wet, &
66# endif
67# ifdef BBL_MODEL
68 & bbl(ng) % bustrc, &
69 & bbl(ng) % bvstrc, &
70 & bbl(ng) % bustrw, &
71 & bbl(ng) % bvstrw, &
72 & bbl(ng) % bustrcwmax, &
73 & bbl(ng) % bvstrcwmax, &
74# endif
75 & forces(ng) % bustr, &
76 & forces(ng) % bvstr, &
77 & ocean(ng) % t, &
78 & sedbed(ng) % ero_flux, &
79 & sedbed(ng) % settling_flux, &
80# if defined SED_MORPH
81 & sedbed(ng) % bed_thick, &
82# endif
83 & sedbed(ng) % bed, &
84 & sedbed(ng) % bed_frac, &
85 & sedbed(ng) % bed_mass, &
86 & sedbed(ng) % bottom)
87# ifdef PROFILE
88 CALL wclock_off (ng, inlm, 16, __line__, myfile)
89# endif
90!
91 RETURN
92 END SUBROUTINE sed_fluxes
93!
94!***********************************************************************
95 SUBROUTINE sed_fluxes_tile (ng, tile, &
96 & LBi, UBi, LBj, UBj, &
97 & IminS, ImaxS, JminS, JmaxS, &
98 & nstp, nnew, &
99 & Hz, &
100# ifdef WET_DRY
101 & rmask_wet, &
102# endif
103# ifdef BBL_MODEL
104 & bustrc, bvstrc, &
105 & bustrw, bvstrw, &
106 & bustrcwmax, bvstrcwmax, &
107# endif
108 & bustr, bvstr, &
109 & t, &
110 & ero_flux, settling_flux, &
111# if defined SED_MORPH
112 & bed_thick, &
113# endif
114 & bed, bed_frac, bed_mass, &
115 & bottom)
116!***********************************************************************
117!
118 USE mod_param
119 USE mod_scalars
120 USE mod_sediment
121!
122!
123! Imported variable declarations.
124!
125 integer, intent(in) :: ng, tile
126 integer, intent(in) :: LBi, UBi, LBj, UBj
127 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
128 integer, intent(in) :: nstp, nnew
129!
130# ifdef ASSUMED_SHAPE
131 real(r8), intent(in) :: Hz(LBi:,LBj:,:)
132# ifdef WET_DRY
133 real(r8), intent(in) :: rmask_wet(LBi:,LBj:)
134# endif
135# ifdef BBL_MODEL
136 real(r8), intent(in) :: bustrc(LBi:,LBj:)
137 real(r8), intent(in) :: bvstrc(LBi:,LBj:)
138 real(r8), intent(in) :: bustrw(LBi:,LBj:)
139 real(r8), intent(in) :: bvstrw(LBi:,LBj:)
140 real(r8), intent(in) :: bustrcwmax(LBi:,LBj:)
141 real(r8), intent(in) :: bvstrcwmax(LBi:,LBj:)
142# endif
143 real(r8), intent(in) :: bustr(LBi:,LBj:)
144 real(r8), intent(in) :: bvstr(LBi:,LBj:)
145# if defined SED_MORPH
146 real(r8), intent(inout):: bed_thick(LBi:,LBj:,:)
147# endif
148 real(r8), intent(inout) :: t(LBi:,LBj:,:,:,:)
149 real(r8), intent(inout) :: ero_flux(LBi:,LBj:,:)
150 real(r8), intent(inout) :: settling_flux(LBi:,LBj:,:)
151 real(r8), intent(inout) :: bed(LBi:,LBj:,:,:)
152 real(r8), intent(inout) :: bed_frac(LBi:,LBj:,:,:)
153 real(r8), intent(inout) :: bed_mass(LBi:,LBj:,:,:,:)
154 real(r8), intent(inout) :: bottom(LBi:,LBj:,:)
155# else
156 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
157# ifdef WET_DRY
158 real(r8), intent(in) :: rmask_wet(LBi:UBi,LBj:UBj)
159# endif
160# ifdef BBL_MODEL
161 real(r8), intent(in) :: bustrc(LBi:UBi,LBj:UBj)
162 real(r8), intent(in) :: bvstrc(LBi:UBi,LBj:UBj)
163 real(r8), intent(in) :: bustrw(LBi:UBi,LBj:UBj)
164 real(r8), intent(in) :: bvstrw(LBi:UBi,LBj:UBj)
165 real(r8), intent(in) :: bustrcwmax(LBi:UBi,LBj:UBj)
166 real(r8), intent(in) :: bvstrcwmax(LBi:UBi,LBj:UBj)
167# endif
168 real(r8), intent(in) :: bustr(LBi:UBi,LBj:UBj)
169 real(r8), intent(in) :: bvstr(LBi:UBi,LBj:UBj)
170# if defined SED_MORPH
171 real(r8), intent(inout):: bed_thick(LBi:UBi,LBj:UBj,3)
172# endif
173 real(r8), intent(inout) :: t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
174 real(r8), intent(inout) :: ero_flux(LBi:UBi,LBj:UBj,NST)
175 real(r8), intent(inout) :: settling_flux(LBi:UBi,LBj:UBj,NST)
176 real(r8), intent(inout) :: bed(LBi:UBi,LBj:UBj,Nbed,MBEDP)
177 real(r8), intent(inout) :: bed_frac(LBi:UBi,LBj:UBj,Nbed,NST)
178 real(r8), intent(inout) :: bed_mass(LBi:UBi,LBj:UBj,Nbed,1:2,NST)
179 real(r8), intent(inout) :: bottom(LBi:UBi,LBj:UBj,MBOTP)
180# endif
181!
182! Local variable declarations.
183!
184 integer :: Ksed, i, indx, ised, j, k, ks
185 integer :: bnew
186
187 real(r8) :: cff, cff1, cff2, cff3, cff4
188
189 real(r8), dimension(IminS:ImaxS,N(ng)) :: Hz_inv
190
191 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tau_w
192
193# include "set_bounds.h"
194
195# ifdef BEDLOAD
196 bnew=nnew
197# else
198 bnew=nstp
199# endif
200!
201!-----------------------------------------------------------------------
202! Compute sediment deposition, resuspension, and erosion.
203!-----------------------------------------------------------------------
204!
205# if defined BEDLOAD_MPM || defined SUSPLOAD
206# ifdef BBL_MODEL
207 DO j=jstr-1,jend+1
208 DO i=istr-1,iend+1
209 tau_w(i,j)=sqrt(bustrcwmax(i,j)*bustrcwmax(i,j)+ &
210 & bvstrcwmax(i,j)*bvstrcwmax(i,j))
211# ifdef WET_DRY
212 tau_w(i,j)=tau_w(i,j)*rmask_wet(i,j)
213# endif
214 END DO
215 END DO
216# else
217 DO j=jstrm1,jendp1
218 DO i=istrm1,iendp1
219 tau_w(i,j)=0.5_r8*sqrt((bustr(i,j)+bustr(i+1,j))* &
220 & (bustr(i,j)+bustr(i+1,j))+ &
221 & (bvstr(i,j)+bvstr(i,j+1))* &
222 & (bvstr(i,j)+bvstr(i,j+1)))
223# ifdef WET_DRY
224 tau_w(i,j)=tau_w(i,j)*rmask_wet(i,j)
225# endif
226 END DO
227 END DO
228# endif
229# endif
230!
231!-----------------------------------------------------------------------
232! Sediment deposition and resuspension near the bottom.
233!-----------------------------------------------------------------------
234!
235! The deposition and resuspension of sediment on the bottom "bed"
236! is due to precepitation settling_flux, already computed, and the
237! resuspension (erosion, hence called ero_flux). The resuspension is
238! applied to the bottom-most grid box value qc(:,1) so the total mass
239! is conserved. Restrict "ero_flux" so that "bed" cannot go negative
240! after both fluxes are applied.
241!
242 j_loop : DO j=jstr,jend
243 DO k=1,n(ng)
244 DO i=istr,iend
245 hz_inv(i,k)=1.0_r8/hz(i,j,k)
246 END DO
247 END DO
248!
249 sed_loop: DO ised=1,nst
250 indx=idsed(ised)
251 DO i=istr,iend
252!
253! Calculate critical shear stress in Pa
254!
255# if defined COHESIVE_BED
256 cff = rho0/bed(i,j,1,ibtcr)
257# elif defined MIXED_BED
258 cff = max(bottom(i,j,idprp)*bed(i,j,1,ibtcr)/rho0+ &
259 & (1.0_r8-bottom(i,j,idprp))*tau_ce(ised,ng), &
260 & tau_ce(ised,ng))
261 cff=1.0_r8/cff
262# else
263 cff=1.0_r8/tau_ce(ised,ng)
264# endif
265!
266! Compute erosion, ero_flux (kg/m2).
267!
268 cff1=(1.0_r8-bed(i,j,1,iporo))*bed_frac(i,j,1,ised)
269 cff2=dt(ng)*erate(ised,ng)*cff1
270 cff3=srho(ised,ng)*cff1
271 cff4=bed_mass(i,j,1,bnew,ised)
272 ero_flux(i,j,ised)= &
273 & min(max(0.0_r8,cff2*(cff*tau_w(i,j)-1.0_r8)), &
274 & min(cff3*bottom(i,j,iactv),cff4)+ &
275 & settling_flux(i,j,ised))
276!
277! Update global tracer variables (m Tunits for nnew indx, Tuints for 3)
278! for erosive flux.
279!
280 t(i,j,1,nnew,indx)=t(i,j,1,nnew,indx)+ero_flux(i,j,ised)
281 END DO
282 END DO sed_loop
283 END DO j_loop
284!
285 RETURN
286 END SUBROUTINE sed_fluxes_tile
287#endif
288 END MODULE sed_fluxes_mod
type(t_bbl), dimension(:), allocatable bbl
Definition mod_bbl.F:62
type(t_forces), dimension(:), allocatable forces
Definition mod_forces.F:554
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
integer, parameter inlm
Definition mod_param.F:662
real(dp), dimension(:), allocatable dt
real(dp) rho0
type(t_sedbed), dimension(:), allocatable sedbed
Definition sedbed_mod.h:157
integer, parameter ibtcr
integer, parameter idprp
real(r8), dimension(:,:), allocatable erate
integer, parameter iactv
real(r8), dimension(:,:), allocatable srho
integer, dimension(:), allocatable idsed
integer, parameter iporo
real(r8), dimension(:,:), allocatable tau_ce
integer, dimension(:), allocatable nnew
integer, dimension(:), allocatable nstp
subroutine sed_fluxes_tile(ng, tile, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, nstp, nnew, hz, rmask_wet, bustrc, bvstrc, bustrw, bvstrw, bustrcwmax, bvstrcwmax, bustr, bvstr, t, ero_flux, settling_flux, bed_thick, bed, bed_frac, bed_mass, bottom)
Definition sed_fluxes.F:116
subroutine, public sed_fluxes(ng, tile)
Definition sed_fluxes.F:35
recursive subroutine wclock_off(ng, model, region, line, routine)
Definition timers.F:148
recursive subroutine wclock_on(ng, model, region, line, routine)
Definition timers.F:3