| 1 | #include "cppdefs.h"
|
|---|
| 2 |
|
|---|
| 3 | MODULE sed_bed_mod
|
|---|
| 4 |
|
|---|
| 5 | #if defined NONLINEAR && defined SEDIMENT && !defined COHESIVE_BED
|
|---|
| 6 | !
|
|---|
| 7 | !svn $Id: sed_bed.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 routine computes sediment bed layer stratigraphy. !
|
|---|
| 15 | ! !
|
|---|
| 16 | ! Warner, J.C., C.R. Sherwood, R.P. Signell, C.K. Harris, and H.G. !
|
|---|
| 17 | ! Arango, 2008: Development of a three-dimensional, regional, !
|
|---|
| 18 | ! coupled wave, current, and sediment-transport model, Computers !
|
|---|
| 19 | ! & Geosciences, 34, 1284-1306. !
|
|---|
| 20 | ! !
|
|---|
| 21 | !=======================================================================
|
|---|
| 22 | !
|
|---|
| 23 | implicit none
|
|---|
| 24 |
|
|---|
| 25 | PRIVATE
|
|---|
| 26 | PUBLIC :: sed_bed
|
|---|
| 27 |
|
|---|
| 28 | CONTAINS
|
|---|
| 29 | !
|
|---|
| 30 | !***********************************************************************
|
|---|
| 31 | SUBROUTINE sed_bed (ng, tile)
|
|---|
| 32 | !***********************************************************************
|
|---|
| 33 | !
|
|---|
| 34 | USE mod_param
|
|---|
| 35 | USE mod_forces
|
|---|
| 36 | USE mod_grid
|
|---|
| 37 | USE mod_ocean
|
|---|
| 38 | USE mod_stepping
|
|---|
| 39 | # ifdef BBL_MODEL
|
|---|
| 40 | USE mod_bbl
|
|---|
| 41 | # endif
|
|---|
| 42 | !
|
|---|
| 43 | ! Imported variable declarations.
|
|---|
| 44 | !
|
|---|
| 45 | integer, intent(in) :: ng, tile
|
|---|
| 46 | !
|
|---|
| 47 | ! Local variable declarations.
|
|---|
| 48 | !
|
|---|
| 49 | # include "tile.h"
|
|---|
| 50 | !
|
|---|
| 51 | # ifdef PROFILE
|
|---|
| 52 | CALL wclock_on (ng, iNLM, 16)
|
|---|
| 53 | # endif
|
|---|
| 54 | CALL sed_bed_tile (ng, tile, &
|
|---|
| 55 | & LBi, UBi, LBj, UBj, &
|
|---|
| 56 | & IminS, ImaxS, JminS, JmaxS, &
|
|---|
| 57 | & nstp(ng), nnew(ng), &
|
|---|
| 58 | # ifdef WET_DRY
|
|---|
| 59 | & GRID(ng) % rmask_wet, &
|
|---|
| 60 | # endif
|
|---|
| 61 | # ifdef BBL_MODEL
|
|---|
| 62 | & BBL(ng) % bustrc, &
|
|---|
| 63 | & BBL(ng) % bvstrc, &
|
|---|
| 64 | & BBL(ng) % bustrw, &
|
|---|
| 65 | & BBL(ng) % bvstrw, &
|
|---|
| 66 | & BBL(ng) % bustrcwmax, &
|
|---|
| 67 | & BBL(ng) % bvstrcwmax, &
|
|---|
| 68 | # else
|
|---|
| 69 | & FORCES(ng) % bustr, &
|
|---|
| 70 | & FORCES(ng) % bvstr, &
|
|---|
| 71 | # endif
|
|---|
| 72 | & OCEAN(ng) % t, &
|
|---|
| 73 | & OCEAN(ng) % ero_flux, &
|
|---|
| 74 | & OCEAN(ng) % settling_flux, &
|
|---|
| 75 | # if defined SED_MORPH
|
|---|
| 76 | & GRID(ng) % bed_thick, &
|
|---|
| 77 | # endif
|
|---|
| 78 | & OCEAN(ng) % bed, &
|
|---|
| 79 | & OCEAN(ng) % bed_frac, &
|
|---|
| 80 | & OCEAN(ng) % bed_mass, &
|
|---|
| 81 | & OCEAN(ng) % bottom)
|
|---|
| 82 | # ifdef PROFILE
|
|---|
| 83 | CALL wclock_off (ng, iNLM, 16)
|
|---|
| 84 | # endif
|
|---|
| 85 | RETURN
|
|---|
| 86 | END SUBROUTINE sed_bed
|
|---|
| 87 | !
|
|---|
| 88 | !***********************************************************************
|
|---|
| 89 | SUBROUTINE sed_bed_tile (ng, tile, &
|
|---|
| 90 | & LBi, UBi, LBj, UBj, &
|
|---|
| 91 | & IminS, ImaxS, JminS, JmaxS, &
|
|---|
| 92 | & nstp, nnew, &
|
|---|
| 93 | # ifdef WET_DRY
|
|---|
| 94 | & rmask_wet, &
|
|---|
| 95 | # endif
|
|---|
| 96 | # ifdef BBL_MODEL
|
|---|
| 97 | & bustrc, bvstrc, &
|
|---|
| 98 | & bustrw, bvstrw, &
|
|---|
| 99 | & bustrcwmax, bvstrcwmax, &
|
|---|
| 100 | # else
|
|---|
| 101 | & bustr, bvstr, &
|
|---|
| 102 | # endif
|
|---|
| 103 | & t, &
|
|---|
| 104 | & ero_flux, settling_flux, &
|
|---|
| 105 | # if defined SED_MORPH
|
|---|
| 106 | & bed_thick, &
|
|---|
| 107 | # endif
|
|---|
| 108 | & bed, bed_frac, bed_mass, &
|
|---|
| 109 | & bottom)
|
|---|
| 110 | !***********************************************************************
|
|---|
| 111 | !
|
|---|
| 112 | USE mod_param
|
|---|
| 113 | USE mod_scalars
|
|---|
| 114 | USE mod_sediment
|
|---|
| 115 | !
|
|---|
| 116 | USE bc_3d_mod, ONLY : bc_r3d_tile
|
|---|
| 117 | # if defined EW_PERIODIC || defined NS_PERIODIC
|
|---|
| 118 | USE exchange_2d_mod, ONLY : exchange_r2d_tile
|
|---|
| 119 | # endif
|
|---|
| 120 | # ifdef DISTRIBUTE
|
|---|
| 121 | USE mp_exchange_mod, ONLY : mp_exchange3d, mp_exchange4d
|
|---|
| 122 | # endif
|
|---|
| 123 | !
|
|---|
| 124 | ! Imported variable declarations.
|
|---|
| 125 | !
|
|---|
| 126 | integer, intent(in) :: ng, tile
|
|---|
| 127 | integer, intent(in) :: LBi, UBi, LBj, UBj
|
|---|
| 128 | integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
|
|---|
| 129 | integer, intent(in) :: nstp, nnew
|
|---|
| 130 | !
|
|---|
| 131 | # ifdef ASSUMED_SHAPE
|
|---|
| 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 | # else
|
|---|
| 143 | real(r8), intent(in) :: bustr(LBi:,LBj:)
|
|---|
| 144 | real(r8), intent(in) :: bvstr(LBi:,LBj:)
|
|---|
| 145 | # endif
|
|---|
| 146 | # if defined SED_MORPH
|
|---|
| 147 | real(r8), intent(inout):: bed_thick(LBi:,LBj:,:)
|
|---|
| 148 | # endif
|
|---|
| 149 | real(r8), intent(inout) :: t(LBi:,LBj:,:,:,:)
|
|---|
| 150 | real(r8), intent(inout) :: ero_flux(LBi:,LBj:,:)
|
|---|
| 151 | real(r8), intent(inout) :: settling_flux(LBi:,LBj:,:)
|
|---|
| 152 | real(r8), intent(inout) :: bed(LBi:,LBj:,:,:)
|
|---|
| 153 | real(r8), intent(inout) :: bed_frac(LBi:,LBj:,:,:)
|
|---|
| 154 | real(r8), intent(inout) :: bed_mass(LBi:,LBj:,:,:,:)
|
|---|
| 155 | real(r8), intent(inout) :: bottom(LBi:,LBj:,:)
|
|---|
| 156 | # else
|
|---|
| 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 | # else
|
|---|
| 168 | real(r8), intent(in) :: bustr(LBi:UBi,LBj:UBj)
|
|---|
| 169 | real(r8), intent(in) :: bvstr(LBi:UBi,LBj:UBj)
|
|---|
| 170 | # endif
|
|---|
| 171 | # if defined SED_MORPH
|
|---|
| 172 | real(r8), intent(inout):: bed_thick(LBi:UBi,LBj:UBj,2)
|
|---|
| 173 | # endif
|
|---|
| 174 | real(r8), intent(inout) :: t(LBi:UBi,LBj:UBj,N(ng),3,NT(ng))
|
|---|
| 175 | real(r8), intent(inout) :: ero_flux(LBi:UBi,LBj:UBj,NST)
|
|---|
| 176 | real(r8), intent(inout) :: settling_flux(LBi:UBi,LBj:UBj,NST)
|
|---|
| 177 | real(r8), intent(inout) :: bed(LBi:UBi,LBj:UBj,Nbed,MBEDP)
|
|---|
| 178 | real(r8), intent(inout) :: bed_frac(LBi:UBi,LBj:UBj,Nbed,NST)
|
|---|
| 179 | real(r8), intent(inout) :: bed_mass(LBi:UBi,LBj:UBj,Nbed,1:2,NST)
|
|---|
| 180 | real(r8), intent(inout) :: bottom(LBi:UBi,LBj:UBj,MBOTP)
|
|---|
| 181 | # endif
|
|---|
| 182 | !
|
|---|
| 183 | ! Local variable declarations.
|
|---|
| 184 | !
|
|---|
| 185 | # ifdef DISTRIBUTE
|
|---|
| 186 | # ifdef EW_PERIODIC
|
|---|
| 187 | logical :: EWperiodic=.TRUE.
|
|---|
| 188 | # else
|
|---|
| 189 | logical :: EWperiodic=.FALSE.
|
|---|
| 190 | # endif
|
|---|
| 191 | # ifdef NS_PERIODIC
|
|---|
| 192 | logical :: NSperiodic=.TRUE.
|
|---|
| 193 | # else
|
|---|
| 194 | logical :: NSperiodic=.FALSE.
|
|---|
| 195 | # endif
|
|---|
| 196 | # endif
|
|---|
| 197 | integer :: Ksed, i, ised, j, k, ks
|
|---|
| 198 | integer :: bnew
|
|---|
| 199 |
|
|---|
| 200 | real(r8), parameter :: eps = 1.0E-14_r8
|
|---|
| 201 |
|
|---|
| 202 | real(r8) :: cff, cff1, cff2, cff3
|
|---|
| 203 | real(r8) :: thck_avail, thck_to_add
|
|---|
| 204 |
|
|---|
| 205 | real(r8), dimension(IminS:ImaxS,NST) :: dep_mass
|
|---|
| 206 | real(r8), dimension(IminS:ImaxS,JminS:JmaxS) :: tau_w
|
|---|
| 207 |
|
|---|
| 208 | # include "set_bounds.h"
|
|---|
| 209 |
|
|---|
| 210 | # ifdef BEDLOAD
|
|---|
| 211 | bnew=nnew
|
|---|
| 212 | # else
|
|---|
| 213 | bnew=nstp
|
|---|
| 214 | # endif
|
|---|
| 215 | !
|
|---|
| 216 | !-----------------------------------------------------------------------
|
|---|
| 217 | ! Compute sediment bed layer stratigraphy.
|
|---|
| 218 | !-----------------------------------------------------------------------
|
|---|
| 219 | !
|
|---|
| 220 | # if defined BEDLOAD_MPM || defined SUSPLOAD
|
|---|
| 221 | # ifdef BBL_MODEL
|
|---|
| 222 | DO j=Jstr-1,Jend+1
|
|---|
| 223 | DO i=Istr-1,Iend+1
|
|---|
| 224 | tau_w(i,j)=SQRT(bustrcwmax(i,j)*bustrcwmax(i,j)+ &
|
|---|
| 225 | & bvstrcwmax(i,j)*bvstrcwmax(i,j))
|
|---|
| 226 | # ifdef WET_DRY
|
|---|
| 227 | tau_w(i,j)=tau_w(i,j)*rmask_wet(i,j)
|
|---|
| 228 | # endif
|
|---|
| 229 | END DO
|
|---|
| 230 | END DO
|
|---|
| 231 | # else
|
|---|
| 232 | # ifdef EW_PERIODIC
|
|---|
| 233 | # define I_RANGE Istr-1,Iend+1
|
|---|
| 234 | # else
|
|---|
| 235 | # define I_RANGE MAX(Istr-1,1),MIN(Iend+1,Lm(ng))
|
|---|
| 236 | # endif
|
|---|
| 237 | # ifdef NS_PERIODIC
|
|---|
| 238 | # define J_RANGE Jstr-1,Jend+1
|
|---|
| 239 | # else
|
|---|
| 240 | # define J_RANGE MAX(Jstr-1,1),MIN(Jend+1,Mm(ng))
|
|---|
| 241 | # endif
|
|---|
| 242 | DO i=I_RANGE
|
|---|
| 243 | DO j=J_RANGE
|
|---|
| 244 | tau_w(i,j)=0.5_r8*SQRT((bustr(i,j)+bustr(i+1,j))* &
|
|---|
| 245 | & (bustr(i,j)+bustr(i+1,j))+ &
|
|---|
| 246 | & (bvstr(i,j)+bvstr(i,j+1))* &
|
|---|
| 247 | & (bvstr(i,j)+bvstr(i,j+1)))
|
|---|
| 248 | # ifdef WET_DRY
|
|---|
| 249 | tau_w(i,j)=tau_w(i,j)*rmask_wet(i,j)
|
|---|
| 250 | # endif
|
|---|
| 251 | END DO
|
|---|
| 252 | END DO
|
|---|
| 253 | # undef I_RANGE
|
|---|
| 254 | # undef J_RANGE
|
|---|
| 255 | # endif
|
|---|
| 256 | # endif
|
|---|
| 257 | !
|
|---|
| 258 | !-----------------------------------------------------------------------
|
|---|
| 259 | ! Update bed properties according to ero_flux and dep_flux.
|
|---|
| 260 | !-----------------------------------------------------------------------
|
|---|
| 261 | !
|
|---|
| 262 | # ifdef SUSPLOAD
|
|---|
| 263 | J_LOOP : DO j=Jstr,Jend
|
|---|
| 264 | SED_LOOP: DO ised=1,NST
|
|---|
| 265 | !
|
|---|
| 266 | ! The deposition and resuspension of sediment on the bottom "bed"
|
|---|
| 267 | ! is due to precepitation flux FC(:,0), already computed, and the
|
|---|
| 268 | ! resuspension (erosion, hence called ero_flux). The resuspension is
|
|---|
| 269 | ! applied to the bottom-most grid box value qc(:,1) so the total mass
|
|---|
| 270 | ! is conserved. Restrict "ero_flux" so that "bed" cannot go negative
|
|---|
| 271 | ! after both fluxes are applied.
|
|---|
| 272 | !
|
|---|
| 273 | DO i=Istr,Iend
|
|---|
| 274 | dep_mass(i,ised)=0.0_r8
|
|---|
| 275 |
|
|---|
| 276 | # ifdef SED_MORPH
|
|---|
| 277 | !
|
|---|
| 278 | ! Apply morphology factor.
|
|---|
| 279 | !
|
|---|
| 280 | ero_flux(i,j,ised)=ero_flux(i,j,ised)*morph_fac(ised,ng)
|
|---|
| 281 | settling_flux(i,j,ised)=settling_flux(i,j,ised)* &
|
|---|
| 282 | & morph_fac(ised,ng)
|
|---|
| 283 | !
|
|---|
| 284 | # endif
|
|---|
| 285 | IF ((ero_flux(i,j,ised)-settling_flux(i,j,ised)).lt. &
|
|---|
| 286 | & 0.0_r8) THEN
|
|---|
| 287 | !
|
|---|
| 288 | ! If first time step of deposit, then store deposit material in
|
|---|
| 289 | ! temporary array, dep_mass.
|
|---|
| 290 | !
|
|---|
| 291 | IF ((time(ng).gt.(bed(i,j,1,iaged)+1.1_r8*dt(ng))).and. &
|
|---|
| 292 | & (bed(i,j,1,ithck).gt.newlayer_thick(ng))) THEN
|
|---|
| 293 | dep_mass(i,ised)=settling_flux(i,j,ised)- &
|
|---|
| 294 | & ero_flux(i,j,ised)
|
|---|
| 295 | END IF
|
|---|
| 296 | bed(i,j,1,iaged)=time(ng)
|
|---|
| 297 | END IF
|
|---|
| 298 | !
|
|---|
| 299 | ! Update bed mass arrays.
|
|---|
| 300 | !
|
|---|
| 301 | bed_mass(i,j,1,nnew,ised)=MAX(bed_mass(i,j,1,bnew,ised)- &
|
|---|
| 302 | & (ero_flux(i,j,ised)- &
|
|---|
| 303 | & settling_flux(i,j,ised)), &
|
|---|
| 304 | & 0.0_r8)
|
|---|
| 305 | DO k=2,Nbed
|
|---|
| 306 | bed_mass(i,j,k,nnew,ised)=bed_mass(i,j,k,nstp,ised)
|
|---|
| 307 | END DO
|
|---|
| 308 | END DO
|
|---|
| 309 | END DO SED_LOOP
|
|---|
| 310 | !
|
|---|
| 311 | ! If first time step of deposit, create new layer and combine bottom
|
|---|
| 312 | ! two bed layers.
|
|---|
| 313 | !
|
|---|
| 314 | DO i=Istr,Iend
|
|---|
| 315 | cff=0.0_r8
|
|---|
| 316 | !
|
|---|
| 317 | ! Determine if deposition ocurred here.
|
|---|
| 318 | !
|
|---|
| 319 | IF (Nbed.gt.1) THEN
|
|---|
| 320 | DO ised=1,NST
|
|---|
| 321 | cff=cff+dep_mass(i,ised)
|
|---|
| 322 | END DO
|
|---|
| 323 | IF (cff.gt.0.0_r8) THEN
|
|---|
| 324 | !
|
|---|
| 325 | ! Combine bottom layers.
|
|---|
| 326 | !
|
|---|
| 327 | bed(i,j,Nbed,iporo)=0.5_r8*(bed(i,j,Nbed-1,iporo)+ &
|
|---|
| 328 | & bed(i,j,Nbed,iporo))
|
|---|
| 329 | bed(i,j,Nbed,iaged)=0.5_r8*(bed(i,j,Nbed-1,iaged)+ &
|
|---|
| 330 | & bed(i,j,Nbed,iaged))
|
|---|
| 331 | DO ised=1,NST
|
|---|
| 332 | bed_mass(i,j,Nbed,nnew,ised)= &
|
|---|
| 333 | & bed_mass(i,j,Nbed-1,nnew,ised)+ &
|
|---|
| 334 | & bed_mass(i,j,Nbed ,nnew,ised)
|
|---|
| 335 | END DO
|
|---|
| 336 | !
|
|---|
| 337 | ! Push layers down.
|
|---|
| 338 | !
|
|---|
| 339 | DO k=Nbed-1,2,-1
|
|---|
| 340 | bed(i,j,k,iporo)=bed(i,j,k-1,iporo)
|
|---|
| 341 | bed(i,j,k,iaged)=bed(i,j,k-1,iaged)
|
|---|
| 342 | DO ised =1,NST
|
|---|
| 343 | bed_mass(i,j,k,nnew,ised)=bed_mass(i,j,k-1,nnew,ised)
|
|---|
| 344 | END DO
|
|---|
| 345 | END DO
|
|---|
| 346 | !
|
|---|
| 347 | ! Set new top layer parameters.
|
|---|
| 348 | !
|
|---|
| 349 | DO ised=1,NST
|
|---|
| 350 | bed_mass(i,j,2,nnew,ised)=MAX(bed_mass(i,j,2,nnew,ised)- &
|
|---|
| 351 | & dep_mass(i,ised),0.0_r8)
|
|---|
| 352 | bed_mass(i,j,1,nnew,ised)=dep_mass(i,ised)
|
|---|
| 353 | END DO
|
|---|
| 354 | END IF
|
|---|
| 355 | END IF !NBED=1
|
|---|
| 356 | !
|
|---|
| 357 | ! Recalculate thickness and fractions for all layers.
|
|---|
| 358 | !
|
|---|
| 359 | DO k=1,Nbed
|
|---|
| 360 | cff3=0.0_r8
|
|---|
| 361 | DO ised=1,NST
|
|---|
| 362 | cff3=cff3+bed_mass(i,j,k,nnew,ised)
|
|---|
| 363 | END DO
|
|---|
| 364 | IF (cff3.eq.0.0_r8) THEN
|
|---|
| 365 | cff3=eps
|
|---|
| 366 | END IF
|
|---|
| 367 | bed(i,j,k,ithck)=0.0_r8
|
|---|
| 368 | DO ised=1,NST
|
|---|
| 369 | bed_frac(i,j,k,ised)=bed_mass(i,j,k,nnew,ised)/cff3
|
|---|
| 370 | bed(i,j,k,ithck)=MAX(bed(i,j,k,ithck)+ &
|
|---|
| 371 | & bed_mass(i,j,k,nnew,ised)/ &
|
|---|
| 372 | & (Srho(ised,ng)* &
|
|---|
| 373 | & (1.0_r8-bed(i,j,k,iporo))),0.0_r8)
|
|---|
| 374 | END DO
|
|---|
| 375 | END DO
|
|---|
| 376 | END DO
|
|---|
| 377 | END DO J_LOOP
|
|---|
| 378 | !
|
|---|
| 379 | ! End of Suspended Sediment only section.
|
|---|
| 380 | !
|
|---|
| 381 | # endif
|
|---|
| 382 | !
|
|---|
| 383 | ! Ensure top bed layer thickness is greater or equal than active layer
|
|---|
| 384 | ! thickness. If need to add sed to top layer, then entrain from lower
|
|---|
| 385 | ! levels. Create new layers at bottom to maintain Nbed.
|
|---|
| 386 | !
|
|---|
| 387 | J_LOOP2 : DO j=Jstr,Jend
|
|---|
| 388 | DO i=Istr,Iend
|
|---|
| 389 | !
|
|---|
| 390 | ! Calculate active layer thickness, bottom(i,j,iactv).
|
|---|
| 391 | !
|
|---|
| 392 | bottom(i,j,iactv)=MAX(0.0_r8, &
|
|---|
| 393 | & 0.007_r8* &
|
|---|
| 394 | & (tau_w(i,j)-bottom(i,j,itauc))*rho0)+ &
|
|---|
| 395 | & 6.0_r8*bottom(i,j,isd50)
|
|---|
| 396 | !
|
|---|
| 397 | # ifdef SED_MORPH
|
|---|
| 398 | !
|
|---|
| 399 | ! Apply morphology factor.
|
|---|
| 400 | !
|
|---|
| 401 | bottom(i,j,iactv)=MAX(bottom(i,j,iactv)*morph_fac(1,ng), &
|
|---|
| 402 | & bottom(i,j,iactv))
|
|---|
| 403 | # endif
|
|---|
| 404 | !
|
|---|
| 405 | IF (bottom(i,j,iactv).gt.bed(i,j,1,ithck)) THEN
|
|---|
| 406 | IF (Nbed.eq.1) THEN
|
|---|
| 407 | bottom(i,j,iactv)=bed(i,j,1,ithck)
|
|---|
| 408 | ELSE
|
|---|
| 409 | thck_to_add=bottom(i,j,iactv)-bed(i,j,1,ithck)
|
|---|
| 410 | thck_avail=0.0_r8
|
|---|
| 411 | Ksed=1 ! initialize
|
|---|
| 412 | DO k=2,Nbed
|
|---|
| 413 | IF (thck_avail.lt.thck_to_add) THEN
|
|---|
| 414 | thck_avail=thck_avail+bed(i,j,k,ithck)
|
|---|
| 415 | Ksed=k
|
|---|
| 416 | END IF
|
|---|
| 417 | END DO
|
|---|
| 418 | !
|
|---|
| 419 | ! Catch here if there was not enough bed material.
|
|---|
| 420 | !
|
|---|
| 421 | IF (thck_avail.lt.thck_to_add) THEN
|
|---|
| 422 | bottom(i,j,iactv)=bed(i,j,1,ithck)+thck_avail
|
|---|
| 423 | thck_to_add=thck_avail
|
|---|
| 424 | END IF
|
|---|
| 425 | !
|
|---|
| 426 | ! Update bed mass of top layer and fractional layer.
|
|---|
| 427 | !
|
|---|
| 428 | cff2=MAX(thck_avail-thck_to_add,0.0_r8)/ &
|
|---|
| 429 | & MAX(bed(i,j,Ksed,ithck),eps)
|
|---|
| 430 | DO ised=1,NST
|
|---|
| 431 | cff1=0.0_r8
|
|---|
| 432 | DO k=1,Ksed
|
|---|
| 433 | cff1=cff1+bed_mass(i,j,k,nnew,ised)
|
|---|
| 434 | END DO
|
|---|
| 435 | cff3=cff2*bed_mass(i,j,Ksed,nnew,ised)
|
|---|
| 436 | bed_mass(i,j,1 ,nnew,ised)=cff1-cff3
|
|---|
| 437 | bed_mass(i,j,Ksed,nnew,ised)=cff3
|
|---|
| 438 | END DO
|
|---|
| 439 | !
|
|---|
| 440 | ! Update thickness of fractional layer ksource_sed.
|
|---|
| 441 | !
|
|---|
| 442 | bed(i,j,Ksed,ithck)=MAX(thck_avail-thck_to_add,0.0_r8)
|
|---|
| 443 | !
|
|---|
| 444 | ! Upate bed fraction of top layer.
|
|---|
| 445 | !
|
|---|
| 446 | cff3=0.0_r8
|
|---|
| 447 | DO ised=1,NST
|
|---|
| 448 | cff3=cff3+bed_mass(i,j,1,nnew,ised)
|
|---|
| 449 | END DO
|
|---|
| 450 | IF (cff3.eq.0.0_r8) THEN
|
|---|
| 451 | cff3=eps
|
|---|
| 452 | END IF
|
|---|
| 453 | DO ised=1,NST
|
|---|
| 454 | bed_frac(i,j,1,ised)=bed_mass(i,j,1,nnew,ised)/cff3
|
|---|
| 455 | END DO
|
|---|
| 456 | !
|
|---|
| 457 | ! Upate bed thickness of top layer.
|
|---|
| 458 | !
|
|---|
| 459 | bed(i,j,1,ithck)=bottom(i,j,iactv)
|
|---|
| 460 | !
|
|---|
| 461 | ! Pull all layers closer to the surface.
|
|---|
| 462 | !
|
|---|
| 463 | DO k=Ksed,Nbed
|
|---|
| 464 | ks=Ksed-2
|
|---|
| 465 | bed(i,j,k-ks,ithck)=bed(i,j,k,ithck)
|
|---|
| 466 | bed(i,j,k-ks,iporo)=bed(i,j,k,iporo)
|
|---|
| 467 | bed(i,j,k-ks,iaged)=bed(i,j,k,iaged)
|
|---|
| 468 | DO ised=1,NST
|
|---|
| 469 | bed_frac(i,j,k-ks,ised)=bed_frac(i,j,k,ised)
|
|---|
| 470 | bed_mass(i,j,k-ks,nnew,ised)=bed_mass(i,j,k,nnew,ised)
|
|---|
| 471 | END DO
|
|---|
| 472 | END DO
|
|---|
| 473 | !
|
|---|
| 474 | ! Add new layers onto the bottom. Split what was in the bottom layer to
|
|---|
| 475 | ! fill these new empty cells. ("ks" is the number of new layers).
|
|---|
| 476 | !
|
|---|
| 477 | ks=Ksed-2
|
|---|
| 478 | cff=1.0_r8/REAL(ks+1,r8)
|
|---|
| 479 | DO k=Nbed,Nbed-ks,-1
|
|---|
| 480 | bed(i,j,k,ithck)=bed(i,j,Nbed-ks,ithck)*cff
|
|---|
| 481 | bed(i,j,k,iaged)=bed(i,j,Nbed-ks,iaged)
|
|---|
| 482 | DO ised=1,NST
|
|---|
| 483 | bed_frac(i,j,k,ised)=bed_frac(i,j,Nbed-ks,ised)
|
|---|
| 484 | bed_mass(i,j,k,nnew,ised)= &
|
|---|
| 485 | & bed_mass(i,j,Nbed-ks,nnew,ised)*cff
|
|---|
| 486 | END DO
|
|---|
| 487 | END DO
|
|---|
| 488 | END IF ! Nbed > 1
|
|---|
| 489 | END IF ! increase top bed layer
|
|---|
| 490 | END DO
|
|---|
| 491 | END DO J_LOOP2
|
|---|
| 492 | !
|
|---|
| 493 | !-----------------------------------------------------------------------
|
|---|
| 494 | ! Store old bed thickness.
|
|---|
| 495 | !-----------------------------------------------------------------------
|
|---|
| 496 | !
|
|---|
| 497 | # if defined SED_MORPH
|
|---|
| 498 | DO j=JstrR,JendR
|
|---|
| 499 | DO i=IstrR,IendR
|
|---|
| 500 | bed_thick(i,j,nnew)=0.0_r8
|
|---|
| 501 | DO k=1,Nbed
|
|---|
| 502 | bed_thick(i,j,nnew)=bed_thick(i,j,nnew)+ &
|
|---|
| 503 | & bed(i,j,k,ithck)
|
|---|
| 504 | END DO
|
|---|
| 505 | END DO
|
|---|
| 506 | END DO
|
|---|
| 507 | # if defined EW_PERIODIC || defined NS_PERIODIC
|
|---|
| 508 | CALL exchange_r2d_tile (ng, tile, &
|
|---|
| 509 | & LBi, UBi, LBj, UBj, &
|
|---|
| 510 | & bed_thick(:,:,nnew))
|
|---|
| 511 | # endif
|
|---|
| 512 | # endif
|
|---|
| 513 | !
|
|---|
| 514 | !-----------------------------------------------------------------------
|
|---|
| 515 | ! Apply periodic or gradient boundary conditions to property arrays.
|
|---|
| 516 | !-----------------------------------------------------------------------
|
|---|
| 517 | !
|
|---|
| 518 | DO ised=1,NST
|
|---|
| 519 | CALL bc_r3d_tile (ng, tile, &
|
|---|
| 520 | & LBi, UBi, LBj, UBj, 1, Nbed, &
|
|---|
| 521 | & bed_frac(:,:,:,ised))
|
|---|
| 522 | CALL bc_r3d_tile (ng, tile, &
|
|---|
| 523 | & LBi, UBi, LBj, UBj, 1, Nbed, &
|
|---|
| 524 | & bed_mass(:,:,:,nnew,ised))
|
|---|
| 525 | END DO
|
|---|
| 526 | # ifdef DISTRIBUTE
|
|---|
| 527 | CALL mp_exchange4d (ng, tile, iNLM, 2, &
|
|---|
| 528 | & LBi, UBi, LBj, UBj, 1, Nbed, 1, NST, &
|
|---|
| 529 | & NghostPoints, EWperiodic, NSperiodic, &
|
|---|
| 530 | & bed_frac, &
|
|---|
| 531 | & bed_mass(:,:,:,nnew,:))
|
|---|
| 532 | # endif
|
|---|
| 533 |
|
|---|
| 534 | DO i=1,MBEDP
|
|---|
| 535 | CALL bc_r3d_tile (ng, tile, &
|
|---|
| 536 | & LBi, UBi, LBj, UBj, 1, Nbed, &
|
|---|
| 537 | & bed(:,:,:,i))
|
|---|
| 538 | END DO
|
|---|
| 539 | # ifdef DISTRIBUTE
|
|---|
| 540 | CALL mp_exchange4d (ng, tile, iNLM, 1, &
|
|---|
| 541 | & LBi, UBi, LBj, UBj, 1, Nbed, 1, MBEDP, &
|
|---|
| 542 | & NghostPoints, EWperiodic, NSperiodic, &
|
|---|
| 543 | & bed)
|
|---|
| 544 | # endif
|
|---|
| 545 |
|
|---|
| 546 | RETURN
|
|---|
| 547 | END SUBROUTINE sed_bed_tile
|
|---|
| 548 | #endif
|
|---|
| 549 | END MODULE sed_bed_mod
|
|---|