ROMS
Loading...
Searching...
No Matches
ana_sediment.h
Go to the documentation of this file.
1!!
2 SUBROUTINE ana_sediment (ng, tile, model)
3!
4!! git $Id$
5!!======================================================================
6!! Copyright (c) 2002-2025 The ROMS Group !
7!! Licensed under a MIT/X style license !
8!! See License_ROMS.md !
9!=======================================================================
10! !
11! This routine sets initial conditions for sedimen t tracer fields !
12! concentrations (kg/m3) using analytical expressions for sediment !
13! and/or bottom boundary layer configurations. It also sets initial !
14! bed conditions in each sediment layer. !
15! !
16!=======================================================================
17!
18 USE mod_param
19 USE mod_grid
20 USE mod_ncparam
21 USE mod_ocean
22 USE mod_sedbed
23!
24! Imported variable declarations.
25!
26 integer, intent(in) :: ng, tile, model
27!
28! Local variable declarations.
29!
30 character (len=*), parameter :: MyFile = &
31 & __FILE__
32!
33#include "tile.h"
34!
35 CALL ana_sediment_tile (ng, tile, model, &
36 & lbi, ubi, lbj, ubj, &
37 & imins, imaxs, jmins, jmaxs, &
38 & grid(ng) % pm, &
39 & grid(ng) % pn, &
40 & grid(ng) % xr, &
41 & grid(ng) % yr, &
42#if defined BBL_MODEL && (defined MB_BBL || defined SSW_BBL)
43 & ocean(ng) % rho, &
44#endif
45#ifdef SEDIMENT
46 & ocean(ng) % t, &
47 & sedbed(ng) % bed, &
48 & sedbed(ng) % bed_frac, &
49 & sedbed(ng) % bed_mass, &
50#endif
51 & sedbed(ng) % bottom)
52!
53! Set analytical header file name used.
54!
55#ifdef DISTRIBUTE
56 IF (lanafile) THEN
57#else
58 IF (lanafile.and.(tile.eq.0)) THEN
59#endif
60 ananame(23)=myfile
61 END IF
62!
63 RETURN
64 END SUBROUTINE ana_sediment
65!
66!***********************************************************************
67 SUBROUTINE ana_sediment_tile (ng, tile, model, &
68 & LBi, UBi, LBj, UBj, &
69 & IminS, ImaxS, JminS, JmaxS, &
70 & pm, pn, &
71 & xr, yr, &
72#if defined BBL_MODEL && (defined MB_BBL || defined SSW_BBL)
73 & rho, &
74#endif
75#ifdef SEDIMENT
76 & t, &
77 & bed, bed_frac, bed_mass, &
78#endif
79 & bottom)
80!***********************************************************************
81!
82 USE mod_param
83 USE mod_ncparam
84 USE mod_scalars
85 USE mod_sediment
86!
87! Imported variable declarations.
88!
89 integer, intent(in) :: ng, tile, model
90 integer, intent(in) :: LBi, UBi, LBj, UBj
91 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
92!
93#ifdef ASSUMED_SHAPE
94 real(r8), intent(in) :: pm(LBi:,LBj:)
95 real(r8), intent(in) :: pn(LBi:,LBj:)
96 real(r8), intent(in) :: xr(LBi:,LBj:)
97 real(r8), intent(in) :: yr(LBi:,LBj:)
98# if defined BBL_MODEL && (defined MB_BBL || defined SSW_BBL)
99 real(r8), intent(in) :: rho(LBi:,LBj:,:)
100# endif
101# ifdef SEDIMENT
102 real(r8), intent(inout) :: t(LBi:,LBj:,:,:,:)
103 real(r8), intent(out) :: bed(LBi:,LBj:,:,:)
104 real(r8), intent(out) :: bed_frac(LBi:,LBj:,:,:)
105 real(r8), intent(out) :: bed_mass(LBi:,LBj:,:,:,:)
106# endif
107 real(r8), intent(inout) :: bottom(LBi:,LBj:,:)
108#else
109 real(r8), intent(in) :: pm(LBi:UBi,LBj:UBj)
110 real(r8), intent(in) :: pn(LBi:UBi,LBj:UBj)
111 real(r8), intent(in) :: xr(LBi:UBi,LBj:UBj)
112 real(r8), intent(in) :: yr(LBi:UBi,LBj:UBj)
113# if defined BBL_MODEL && (defined MB_BBL || defined SSW_BBL)
114 real(r8), intent(in) :: rho(LBi:,LBj:,:)
115# endif
116# ifdef SEDIMENT
117 real(r8), intent(inout) :: t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
118 real(r8), intent(out) :: bed(LBi:UBi,LBj:UBj,Nbed,MBEDP)
119 real(r8), intent(out) :: bed_frac(LBi:UBi,LBj:UBj,Nbed,NST)
120 real(r8), intent(out) :: bed_mass(LBi:UBi,LBj:UBj,Nbed,2,NST)
121# endif
122 real(r8), intent(inout) :: bottom(LBi:UBi,LBj:UBj,MBOTP)
123#endif
124!
125! Local variable declarations.
126!
127#ifdef DISTRIBUTE
128 integer :: Tstr, Tend
129#endif
130 integer :: i, ised, j, k
131!
132 real(r8) :: cff1, cff2, cff3, cff4, Kvisc, phinot
133
134#include "set_bounds.h"
135
136#if defined BBL_MODEL && !defined SEDIMENT
137!
138!-----------------------------------------------------------------------
139! If only bottom boundary layer and not sediment model, set bottom
140! sediment grain diameter (m) and density (kg/m3).
141!-----------------------------------------------------------------------
142!
143# if defined BL_TEST || defined NJ_BIGHT
144 DO j=jstrt,jendt
145 DO i=istrt,iendt
146 bottom(i,j,isd50)=0.0005_r8
147 bottom(i,j,idens)=2650.0_r8
148 END DO
149 END DO
150# elif defined LAKE_SIGNELL || defined ADRIA02
151 DO j=jstrt,jendt
152 DO i=istrt,iendt
153 bottom(i,j,isd50)=0.000150_r8 ! 150 microns
154 bottom(i,j,idens)=2650.0_r8
155 END DO
156 END DO
157# elif defined SED_TOY
158 DO j=jstrt,jendt
159 DO i=istrt,iendt
160 bottom(i,j,isd50)=0.0005_r8
161 bottom(i,j,idens)=2650.0_r8
162 END DO
163 END DO
164# else
165 ana_sediment.h: no values provided for bottom(:,:,isd50) and
166 bottom(:,:,idens)
167# endif
168
169# if defined MB_BBL || defined SSW_BBL
170# undef YALIN
171!
172!-----------------------------------------------------------------------
173! If only Blass bottom boundary layer and not sediment model, set
174! set critical (threshold) bedload stress (m2/s2).
175!-----------------------------------------------------------------------
176!
177# ifdef YALIN
178
179! For more accurate estime of critical bedload stress, consider the
180! Yalin method (Miller et. al, 1977).
181!
182 kvisc=0.0013_r8/rho0
183 DO j=jstrt,jendt
184 DO i=istrt,iendt
185 rhowater=rho(i,j,1)+1000.0_r8
186 cff=sqrt((bottom(i,j,idens)-rhowater)* &
187 & g*bottom(i,j,isd50)*bottom(i,j,isd50)* &
188 & bottom(i,j,isd50)/rhowater)/kvisc
189!! D=bottom(i,j,isd50)*g* &
190!! & ((bottom(i,j,idens)/rho0-1.0_r8)/Kvisk)**(1.0_r8/3.0_r8)
191!! theta_cr=0.3_r8./(1.0_r8+1.2_r8*D)+ &
192!! & 0.055_r8*(1.0_r8-EXP(-0.02_r8*D))
193 IF (cff.lt.100.0_r8) THEN
194 theta_cb=0.041_r8*(log(cff)**2)-0.356_r8*log(cff)-0.977_r8
195!! theta_cb=10**theta_cr
196 ELSE IF (cff.gt.3000.0_r8) THEN
197 theta_cb=0.045_r8
198 ELSE
199 theta_cb=0.132_r8*log(cff)-1.804_r8
200!! theta_cb=10.0_r8**theta_cr
201 ENDIF
202 bottom(i,j,itauc)=(bottom(i,j,idens)-rho0)*g* &
203 & bottom(i,j,isd50)*theta_cb/rho0
204 END DO
205 END DO
206# else
207 DO j=jstrt,jendt
208 DO i=istrt,iendt
209 bottom(i,j,itauc)=0.15_r8/rho0
210 END DO
211 END DO
212# endif
213!
214!-----------------------------------------------------------------------
215! If only Blass bottom boundary layer and not sediment model, set
216! sediment settling velocity (m/s).
217!-----------------------------------------------------------------------
218!
219 kvisc=0.0013_r8/rho0
220 DO j=jstrt,jendt
221 DO i=istrt,iendt
222 bottom(i,j,iwsed)=0.02_r8
223!!
224!! Consider Souslby (1997) estimate of settling velocity.
225!!
226!! D=bottom(i,j,isd50)*g* &
227!! & ((bottom(i,j,idens)/rho0-1.0)/Kvisk)**(1.0_r8/3.0_r8)
228!! bottom(i,j,iwsed)=Kvisc*(SQRT(10.36_r8*10.36_r8+
229!! & 1.049_r8*D*D*D)-10.36_r8)/bottom(i,j,isd50)
230 END DO
231 END DO
232# endif
233#endif
234
235#ifdef SEDIMENT
236!
237!-----------------------------------------------------------------------
238! Initial sediment concentrations in the water column.
239!-----------------------------------------------------------------------
240!
241 DO ised=1,nst
242 DO k=1,n(ng)
243 DO j=jstrt,jendt
244 DO i=istrt,iendt
245 t(i,j,k,1,idsed(ised))=csed(ised,ng)
246 END DO
247 END DO
248 END DO
249 END DO
250!
251!-----------------------------------------------------------------------
252! Initial sediment bed layer properties of age, thickness, porosity,
253! and initialize sediment bottom properites of ripple length, ripple
254! height, and default Zob.
255!-----------------------------------------------------------------------
256!
257# if defined LAKE_SIGNELL || defined ADRIA02
258 DO j=jstrt,jendt
259 DO i=istrt,iendt
260!
261! Set bed layer properties.
262!
263 DO k=1,nbed
264 bed(i,j,k,iaged)=time(ng)
265 bed(i,j,k,ithck)=0.10_r8
266 bed(i,j,k,iporo)=0.90_r8
267 DO ised=1,nst
268 bed_frac(i,j,k,ised)=1.0_r8/real(nst,r8)
269 END DO
270 END DO
271!
272! Set exposed sediment layer properties.
273!
274 bottom(i,j,irlen)=0.10_r8
275 bottom(i,j,irhgt)=0.01_r8
276 bottom(i,j,izdef)=zob(ng)
277 END DO
278 END DO
279# elif defined ESTUARY_TEST
280 DO j=jstrt,jendt
281 DO i=istrt,iendt
282!
283! Set bed layer properties.
284!
285 DO k=1,nbed
286 bed(i,j,k,iaged)=time(ng)
287 bed(i,j,k,ithck)=0.001_r8
288 bed(i,j,k,iporo)=0.90_r8
289 DO ised=1,nst
290 bed_frac(i,j,k,ised)=1.0_r8/real(nst,r8)
291 END DO
292 END DO
293!
294! Set exposed sediment layer properties.
295!
296 bottom(i,j,irlen)=0.10_r8
297 bottom(i,j,irhgt)=0.01_r8
298 bottom(i,j,izdef)=zob(ng)
299 END DO
300 END DO
301# elif defined INLET_TEST
302 DO j=jstrt,jendt
303 DO i=istrt,iendt
304!
305! Set bed layer properties.
306!
307 DO k=1,nbed
308 bed(i,j,k,iaged)=time(ng)
309 bed(i,j,k,ithck)=10.0_r8
310 bed(i,j,k,iporo)=0.50_r8
311 DO ised=1,nst
312 bed_frac(i,j,k,ised)=1.0_r8/real(nst,r8)
313 ENDDO
314 END DO
315!
316! Set exposed sediment layer properties.
317!
318 bottom(i,j,irlen)=0.10_r8
319 bottom(i,j,irhgt)=0.01_r8
320 bottom(i,j,izdef)=zob(ng)
321 END DO
322 END DO
323# elif defined SED_TOY
324 DO j=jstrt,jendt
325 DO i=istrt,iendt
326!
327! Set bed layer properties.
328!
329 DO k=1,nbed
330 bed(i,j,k,iaged)=time(ng)
331 bed(i,j,k,ithck)=0.01_r8
332 bed(i,j,k,iporo)=0.30_r8
333!! DO ised=1,NST
334!! bed_frac(i,j,k,ised)=1.0_r8/REAL(NST,r8)
335!! END DO
336 bed_frac(i,j,k,1)=1.0_r8
337 bed_frac(i,j,k,2)=0.0_r8
338 END DO
339!
340! Set exposed sediment layer properties.
341!
342 bottom(i,j,irlen)=0.10_r8
343 bottom(i,j,irhgt)=0.01_r8
344 bottom(i,j,izdef)=zob(ng)
345 END DO
346 END DO
347 DO j=jstrt,jendt
348 DO i=istrt,iendt
349 END DO
350 END DO
351# elif defined SED_TEST1
352 DO j=jstrt,jendt
353 DO i=istrt,iendt
354!
355! Set bed layer properties.
356!
357 DO k=1,nbed
358 bed(i,j,k,iaged)=time(ng)
359 bed(i,j,k,ithck)=15.00_r8
360 bed(i,j,k,iporo)=0.50_r8
361 DO ised=1,nst
362 bed_frac(i,j,k,ised)=1.0_r8/real(nst,r8)
363 END DO
364 END DO
365!
366! Set exposed sediment layer properties.
367!
368 bottom(i,j,irlen)=0.10_r8
369 bottom(i,j,irhgt)=0.01_r8
370 bottom(i,j,izdef)=zob(ng)
371 END DO
372 END DO
373# elif defined SHOREFACE
374 DO j=jstrt,jendt
375 DO i=istrt,iendt
376!
377! Set bed layer properties.
378!
379 DO k=1,nbed
380 bed(i,j,k,iaged)=time(ng)
381 bed(i,j,k,ithck)=5.0_r8
382 bed(i,j,k,iporo)=0.50_r8
383 DO ised=1,nst
384 bed_frac(i,j,k,ised)=1.0_r8/real(nst,r8)
385 END DO
386 END DO
387!
388! Set exposed sediment layer properties.
389!
390 bottom(i,j,irlen)=0.10_r8
391 bottom(i,j,irhgt)=0.01_r8
392 bottom(i,j,izdef)=zob(ng)
393 END DO
394 END DO
395# elif defined TEST_CHAN
396 DO j=jstrt,jendt
397 DO i=istrt,iendt
398!
399! Set bed layer properties.
400!
401 DO k=1,nbed
402 bed(i,j,k,iaged)=time(ng)
403 bed(i,j,k,ithck)=1.0_r8
404 bed(i,j,k,iporo)=0.90_r8
405 DO ised=1,nst
406 bed_frac(i,j,k,ised)=1.0_r8/real(nst,r8)
407 END DO
408 END DO
409!
410! Set exposed sediment layer properties.
411!
412 bottom(i,j,irlen)=0.0_r8
413 bottom(i,j,irhgt)=0.0_r8
414 bottom(i,j,izdef)=zob(ng)
415 END DO
416 END DO
417# else
418 ana_sediment.h: no values provided for bed, bed_mass, bottom.
419# endif
420!
421!-----------------------------------------------------------------------
422! Initial sediment bed_mass and surface layer properties.
423! Same for all applications.
424!-----------------------------------------------------------------------
425!
426 DO k=1,nbed
427 DO j=jstrt,jendt
428 DO i=istrt,iendt
429!
430! Calculate mass so it is consistent with density, thickness, and
431! porosity.
432!
433 DO ised=1,nst
434 bed_mass(i,j,k,1,ised)=bed(i,j,k,ithck)* &
435 & srho(ised,ng)* &
436 & (1.0_r8-bed(i,j,k,iporo))* &
437 & bed_frac(i,j,k,ised)
438 END DO
439 END DO
440 END DO
441 END DO
442!
443! Set exposed sediment layer properties.
444!
445 DO j=jstrt,jendt
446 DO i=istrt,iendt
447 cff1=1.0_r8
448 cff2=1.0_r8
449 cff3=1.0_r8
450 cff4=1.0_r8
451 DO ised=1,nst
452 cff1=cff1*sd50(ised,ng)**bed_frac(i,j,1,ised)
453 cff2=cff2*srho(ised,ng)**bed_frac(i,j,1,ised)
454 cff3=cff3*wsed(ised,ng)**bed_frac(i,j,1,ised)
455 cff4=cff4*tau_ce(ised,ng)**bed_frac(i,j,1,ised)
456 END DO
457 bottom(i,j,isd50)=cff1
458 bottom(i,j,idens)=cff2
459 bottom(i,j,iwsed)=cff3
460 bottom(i,j,itauc)=cff4
461# ifdef SED_BIODIFF
462 bottom(i,j,idoff)=0.0_r8
463 bottom(i,j,idslp)=0.0_r8
464 bottom(i,j,idtim)=0.0_r8
465 bottom(i,j,idbmx)=0.0_r8
466 bottom(i,j,idbmm)=0.0_r8
467 bottom(i,j,idbzs)=0.0_r8
468 bottom(i,j,idbzm)=0.0_r8
469 bottom(i,j,idbzp)=0.0_r8
470# endif
471 END DO
472 END DO
473#endif
474!
475 RETURN
476 END SUBROUTINE ana_sediment_tile
subroutine ana_sediment(ng, tile, model)
Definition ana_sediment.h:3
subroutine ana_sediment_tile(ng, tile, model, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs, pm, pn, xr, yr, rho, t, bed, bed_frac, bed_mass, bottom)
type(t_grid), dimension(:), allocatable grid
Definition mod_grid.F:365
logical lanafile
character(len=256), dimension(39) ananame
type(t_ocean), dimension(:), allocatable ocean
Definition mod_ocean.F:351
real(r8), dimension(:), allocatable zob
real(dp) g
real(dp), dimension(:), allocatable time
real(dp) rho0
type(t_sedbed), dimension(:), allocatable sedbed
Definition sedbed_mod.h:157
integer, parameter iwsed
integer, parameter irlen
integer, parameter idbzs
integer, parameter irhgt
integer, parameter iaged
real(r8), dimension(:,:), allocatable srho
integer, parameter idbmx
integer, dimension(:), allocatable idsed
integer, parameter iporo
integer, parameter idslp
integer, parameter idens
integer, parameter isd50
integer, parameter idtim
integer, parameter idbzp
integer, parameter izdef
real(r8), dimension(:,:), allocatable sd50
real(r8), dimension(:,:), allocatable wsed
integer, parameter idoff
integer, parameter itauc
integer, parameter ithck
integer, parameter idbzm
real(r8), dimension(:,:), allocatable tau_ce
integer, parameter idbmm
real(r8), dimension(:,:), allocatable csed