Ticket #210: sed_fluxes.F

File sed_fluxes.F, 11.4 KB (added by m.hadfield, 16 years ago)
Line 
1#include "cppdefs.h"
2
3 MODULE sed_fluxes_mod
4
5#if defined NONLINEAR && defined SEDIMENT && defined SUSPLOAD
6!
7!svn $Id: sed_fluxes.F 224 2008-08-19 19:39:15Z arango $
8!==================================================== John C. Warner ===
9! Copyright (c) 2002-2008 The ROMS/TOMS Group Hernan G. Arango !
10! Licensed under a MIT/X style license !
11! See License_ROMS.txt !
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_stepping
42 USE mod_bbl
43!
44! Imported variable declarations.
45!
46 integer, intent(in) :: ng, tile
47!
48! Local variable declarations.
49!
50# include "tile.h"
51!
52# ifdef PROFILE
53 CALL wclock_on (ng, iNLM, 16)
54# endif
55 CALL sed_fluxes_tile (ng, tile, &
56 & LBi, UBi, LBj, UBj, &
57 & IminS, ImaxS, JminS, JmaxS, &
58 & nstp(ng), nnew(ng), &
59 & GRID(ng) % Hz, &
60# ifdef WET_DRY
61 & GRID(ng) % rmask_wet, &
62# endif
63# ifdef BBL_MODEL
64 & BBL(ng) % bustrc, &
65 & BBL(ng) % bvstrc, &
66 & BBL(ng) % bustrw, &
67 & BBL(ng) % bvstrw, &
68 & BBL(ng) % bustrcwmax, &
69 & BBL(ng) % bvstrcwmax, &
70# endif
71 & FORCES(ng) % bustr, &
72 & FORCES(ng) % bvstr, &
73 & OCEAN(ng) % t, &
74 & OCEAN(ng) % ero_flux, &
75 & OCEAN(ng) % settling_flux, &
76# if defined SED_MORPH
77 & GRID(ng) % bed_thick, &
78# endif
79 & OCEAN(ng) % bed, &
80 & OCEAN(ng) % bed_frac, &
81 & OCEAN(ng) % bed_mass, &
82 & OCEAN(ng) % bottom)
83# ifdef PROFILE
84 CALL wclock_off (ng, iNLM, 16)
85# endif
86 RETURN
87 END SUBROUTINE sed_fluxes
88!
89!***********************************************************************
90 SUBROUTINE sed_fluxes_tile (ng, tile, &
91 & LBi, UBi, LBj, UBj, &
92 & IminS, ImaxS, JminS, JmaxS, &
93 & nstp, nnew, &
94 & Hz, &
95# ifdef WET_DRY
96 & rmask_wet, &
97# endif
98# ifdef BBL_MODEL
99 & bustrc, bvstrc, &
100 & bustrw, bvstrw, &
101 & bustrcwmax, bvstrcwmax, &
102# endif
103 & bustr, bvstr, &
104 & t, &
105 & ero_flux, settling_flux, &
106# if defined SED_MORPH
107 & bed_thick, &
108# endif
109 & bed, bed_frac, bed_mass, &
110 & bottom)
111!***********************************************************************
112!
113 USE mod_param
114 USE mod_scalars
115 USE mod_sediment
116!
117!
118! Imported variable declarations.
119!
120 integer, intent(in) :: ng, tile
121 integer, intent(in) :: LBi, UBi, LBj, UBj
122 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
123 integer, intent(in) :: nstp, nnew
124!
125# ifdef ASSUMED_SHAPE
126 real(r8), intent(in) :: Hz(LBi:,LBj:,:)
127# ifdef WET_DRY
128 real(r8), intent(in) :: rmask_wet(LBi:,LBj:)
129# endif
130# ifdef BBL_MODEL
131 real(r8), intent(in) :: bustrc(LBi:,LBj:)
132 real(r8), intent(in) :: bvstrc(LBi:,LBj:)
133 real(r8), intent(in) :: bustrw(LBi:,LBj:)
134 real(r8), intent(in) :: bvstrw(LBi:,LBj:)
135 real(r8), intent(in) :: bustrcwmax(LBi:,LBj:)
136 real(r8), intent(in) :: bvstrcwmax(LBi:,LBj:)
137# endif
138 real(r8), intent(in) :: bustr(LBi:,LBj:)
139 real(r8), intent(in) :: bvstr(LBi:,LBj:)
140# if defined SED_MORPH
141 real(r8), intent(inout):: bed_thick(LBi:,LBj:,:)
142# endif
143 real(r8), intent(inout) :: t(LBi:,LBj:,:,:,:)
144 real(r8), intent(inout) :: ero_flux(LBi:,LBj:,:)
145 real(r8), intent(inout) :: settling_flux(LBi:,LBj:,:)
146 real(r8), intent(inout) :: bed(LBi:,LBj:,:,:)
147 real(r8), intent(inout) :: bed_frac(LBi:,LBj:,:,:)
148 real(r8), intent(inout) :: bed_mass(LBi:,LBj:,:,:,:)
149 real(r8), intent(inout) :: bottom(LBi:,LBj:,:)
150# else
151 real(r8), intent(in) :: Hz(LBi:UBi,LBj:UBj,N(ng))
152# ifdef WET_DRY
153 real(r8), intent(in) :: rmask_wet(LBi:UBi,LBj:UBj)
154# endif
155# ifdef BBL_MODEL
156 real(r8), intent(in) :: bustrc(LBi:UBi,LBj:UBj)
157 real(r8), intent(in) :: bvstrc(LBi:UBi,LBj:UBj)
158 real(r8), intent(in) :: bustrw(LBi:UBi,LBj:UBj)
159 real(r8), intent(in) :: bvstrw(LBi:UBi,LBj:UBj)
160 real(r8), intent(in) :: bustrcwmax(LBi:UBi,LBj:UBj)
161 real(r8), intent(in) :: bvstrcwmax(LBi:UBi,LBj:UBj)
162# endif
163 real(r8), intent(in) :: bustr(LBi:UBi,LBj:UBj)
164 real(r8), intent(in) :: bvstr(LBi:UBi,LBj:UBj)
165# if defined SED_MORPH
166 real(r8), intent(inout):: bed_thick(LBi:UBi,LBj:UBj,2)
167# endif
168 real(r8), intent(inout) :: t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
169 real(r8), intent(inout) :: ero_flux(LBi:UBi,LBj:UBj,NST)
170 real(r8), intent(inout) :: settling_flux(LBi:UBi,LBj:UBj,NST)
171 real(r8), intent(inout) :: bed(LBi:UBi,LBj:UBj,Nbed,MBEDP)
172 real(r8), intent(inout) :: bed_frac(LBi:UBi,LBj:UBj,Nbed,NST)
173 real(r8), intent(inout) :: bed_mass(LBi:UBi,LBj:UBj,Nbed,1:2,NST)
174 real(r8), intent(inout) :: bottom(LBi:UBi,LBj:UBj,MBOTP)
175# endif
176!
177! Local variable declarations.
178!
179# ifdef DISTRIBUTE
180# ifdef EW_PERIODIC
181 logical :: EWperiodic=.TRUE.
182# else
183 logical :: EWperiodic=.FALSE.
184# endif
185# ifdef NS_PERIODIC
186 logical :: NSperiodic=.TRUE.
187# else
188 logical :: NSperiodic=.FALSE.
189# endif
190# endif
191 integer :: Ksed, i, indx, ised, j, k, ks
192 integer :: bnew
193
194 real(r8) :: cff, cff1, cff2, cff3, cff4
195
196 real(r8), dimension(IminS:ImaxS,N(ng)) :: Hz_inv
197
198 real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tau_w
199
200# include "set_bounds.h"
201
202# ifdef BEDLOAD
203 bnew=nnew
204# else
205 bnew=nstp
206# endif
207!
208!-----------------------------------------------------------------------
209! Compute sediment deposition, resuspension, and erosion.
210!-----------------------------------------------------------------------
211!
212# if defined BEDLOAD_MPM || defined SUSPLOAD
213# ifdef BBL_MODEL
214 DO j=Jstr-1,Jend+1
215 DO i=Istr-1,Iend+1
216 tau_w(i,j)=SQRT(bustrcwmax(i,j)*bustrcwmax(i,j)+ &
217 & bvstrcwmax(i,j)*bvstrcwmax(i,j))
218# ifdef WET_DRY
219 tau_w(i,j)=tau_w(i,j)*rmask_wet(i,j)
220# endif
221 END DO
222 END DO
223# else
224# ifdef EW_PERIODIC
225# define I_RANGE Istr-1,Iend+1
226# else
227# define I_RANGE MAX(Istr-1,1),MIN(Iend+1,Lm(ng))
228# endif
229# ifdef NS_PERIODIC
230# define J_RANGE Jstr-1,Jend+1
231# else
232# define J_RANGE MAX(Jstr-1,1),MIN(Jend+1,Mm(ng))
233# endif
234 DO i=I_RANGE
235 DO j=J_RANGE
236 tau_w(i,j)=0.5_r8*SQRT((bustr(i,j)+bustr(i+1,j))* &
237 & (bustr(i,j)+bustr(i+1,j))+ &
238 & (bvstr(i,j)+bvstr(i,j+1))* &
239 & (bvstr(i,j)+bvstr(i,j+1)))
240# ifdef WET_DRY
241 tau_w(i,j)=tau_w(i,j)*rmask_wet(i,j)
242# endif
243 END DO
244 END DO
245# undef I_RANGE
246# undef J_RANGE
247# endif
248# endif
249!
250!-----------------------------------------------------------------------
251! Sediment deposition and resuspension near the bottom.
252!-----------------------------------------------------------------------
253!
254! The deposition and resuspension of sediment on the bottom "bed"
255! is due to precepitation settling_flux, already computed, and the
256! resuspension (erosion, hence called ero_flux). The resuspension is
257! applied to the bottom-most grid box value qc(:,1) so the total mass
258! is conserved. Restrict "ero_flux" so that "bed" cannot go negative
259! after both fluxes are applied.
260!
261 J_LOOP : DO j=Jstr,Jend
262 DO k=1,N(ng)
263 DO i=Istr,Iend
264 Hz_inv(i,k)=1.0_r8/Hz(i,j,k)
265 END DO
266 END DO
267!
268 SED_LOOP: DO ised=1,NST
269 indx=idsed(ised)
270 DO i=Istr,Iend
271!
272! Calculate critical shear stress in Pa
273!
274# if defined COHESIVE_BED
275 cff = rho0/bed(i,j,1,ibtcr)
276# elif defined MIXED_BED
277 cff = MAX(bottom(i,j,idprp)*bed(i,j,1,ibtcr)/rho0+ &
278 & (1.0_r8-bottom(i,j,idprp))*tau_ce(ised,ng), &
279 & tau_ce(ised,ng))
280 cff=1.0_r8/cff
281# else
282 cff=1.0_r8/tau_ce(ised,ng)
283# endif
284!
285! Compute erosion, ero_flux (kg/m2).
286!
287 cff1=(1.0_r8-bed(i,j,1,iporo))*bed_frac(i,j,1,ised)
288 cff2=dt(ng)*Erate(ised,ng)*cff1
289 cff3=Srho(ised,ng)*cff1
290 cff4=bed_mass(i,j,1,bnew,ised)
291 ero_flux(i,j,ised)= &
292 & MIN(MAX(0.0_r8,cff2*(cff*tau_w(i,j)-1.0_r8)), &
293 & MIN(cff3*bottom(i,j,iactv),cff4)+ &
294 & settling_flux(i,j,ised))
295!
296! Update global tracer variables (m Tunits for nnew indx, Tuints for 3)
297! for erosive flux.
298!
299 t(i,j,1,nnew,indx)=t(i,j,1,nnew,indx)+ero_flux(i,j,ised)
300# ifdef TS_MPDATA
301 t(i,j,1,3,indx)=t(i,j,1,3,indx)+ &
302 & ero_flux(i,j,ised)*Hz_inv(i,1)
303# endif
304 END DO
305 END DO SED_LOOP
306 END DO J_LOOP
307
308 RETURN
309 END SUBROUTINE sed_fluxes_tile
310#endif
311 END MODULE sed_fluxes_mod