-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathsediment.f
executable file
·4460 lines (4287 loc) · 253 KB
/
sediment.f
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
! Part of Dflow3d a 3D Navier Stokes solver with variable density for
! simulations of near field dredge plume mixing
! Copyright (C) 2012 Lynyrd de Wit
! This program is free software: you can redistribute it and/or modify
! it under the terms of the GNU General Public License as published by
! the Free Software Foundation, either version 3 of the License, or
! (at your option) any later version.
! This program is distributed in the hope that it will be useful,
! but WITHOUT ANY WARRANTY; without even the implied warranty of
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
! GNU General Public License for more details.
! You should have received a copy of the GNU General Public License
! along with this program. If not, see <http://www.gnu.org/licenses/>.
MODULE sediment
USE nlist
IMPLICIT NONE
CONTAINS
subroutine init_sediment
USE netcdf
implicit none
INTEGER n,n2
REAL Re_p,atm_pres,rhoair_z,zz,ddzz,z_rks(1:kmax),interpseries,gvector
integer :: ncid, rhVarId, status2, ndims, xtype,natts,status
integer, dimension(nf90_max_var_dims) :: dimids
integer size1,size2,size3,size4
DO n=1,nfrac
Re_p=ABS(frac(n)%ws)*frac(n)%dfloc/(ekm_mol/rho_b)
! hindered_settling = 1 !Hindered settling formula [-] 1=Rowe (1987) (smooth Ri-Za org) (default); 2=Garside (1977); 3=Di Felice (1999)
IF (hindered_settling.eq.1) THEN
! hindered settling function acc. to Rowe (1987),
! is the smoothed representation of the original Richardson and Zaki (1954) relations:
! determined for 0.04<c<0.55 and 0.001<Rep<3e4
frac(n)%n=(4.7+0.41*Re_p**0.75)/(1.+0.175*Re_p**0.75) -1.
ELSEIF (hindered_settling.eq.2) THEN
! hindered settling function acc. to Garside (1977),
! higher value of n than Rowe (1987)
! determined for 0.04<c<0.55 and 0.2<Rep<1e3
frac(n)%n=(5.1+0.27*Re_p**0.9)/(1.+0.1*Re_p**0.9)-1.
ELSEIF (hindered_settling.eq.3) THEN
! hindered settling function acc. to Di Felice (1999),
! much higher value of n than Rowe (1987) and Garside (1977)
! determined for 0<c<0.05 and 0.01<Rep<1e3
frac(n)%n=(6.5+0.3*Re_p**0.74)/(1.+0.1*Re_p**0.74)-1.
ENDIF
IF (rank.eq.0) THEN
write(*,*),'fraction #,ws,Re_p,n-factor hindered settling:',n,frac(n)%ws,Re_p,frac(n)%n+1.
ENDIF
ENDDO
rhocorr_air_z=1.
atm_pres=101325. !standard atmospheric pressure in Pa
gvector=sqrt(gx**2+gy**2+gz**2)
DO n=1,nair !correction for comppressible air
DO k=0,k1
zz=k*dz-0.5*dz
ddzz=zz-(depth-frac(nfrac_air(n))%zair_ref_belowsurf) !positive upward and defined wrt zair_ref_belowsurf
rhoair_z = ((frac(nfrac_air(n))%zair_ref_belowsurf-ddzz)*rho_b*gvector+atm_pres) /
& ((frac(nfrac_air(n))%zair_ref_belowsurf)*rho_b*gvector+atm_pres) * frac(nfrac_air(n))%rho
rhocorr_air_z(nfrac_air(n),k) = frac(nfrac_air(n))%rho / rhoair_z
ENDDO
ENDDO
wscorr_z=1.
DO n=1,nfrac
!simplified correction for startup settling velocity as function of initial vertical level, relevant for coarser rocks, correction Miedema 1981 as mentioned in TU Delft MSc thesis Ravelli 2012
IF (frac(n)%CD>0.) THEN
DO k=0,k1
zz=k*dz-0.5*dz
ddzz=-zz+(depth-frac(n)%zair_ref_belowsurf) !positive downward and defined wrt zair_ref_belowsurf
wscorr_z(n,k) = sqrt(1.-exp(-1.5*frac(n)%CD/frac(n)%dpart*rho_b/frac(n)%rho*MAX(0.,ddzz)))
ENDDO
ENDIF
ENDDO
DO k=1,kmax
z_rks(k)=k*dz-0.5*dz
ENDDO
IF (av_slope_z(1)<0.) THEN
av_slope(1:imax,1:jmax,0:k1)=avalanche_slope(1)
ELSE
IF (av_slope_z(1)>0.) CALL writeerror(141)
n=1
DO WHILE (av_slope_z(n+1)>0.)
n=n+1
IF (av_slope_z(n)-av_slope_z(n-1)<0.) CALL writeerror(141)
END DO
IF (av_slope_z(n)<depth) CALL writeerror(141)
IF (avalanche_slope(n)<0.) CALL writeerror(142)
n2=0
DO k=1,kmax
av_slope(1:imax,1:jmax,k)=interpseries(av_slope_z(1:n),avalanche_slope(1:n),n2,z_rks(k))
ENDDO
av_slope(1:imax,1:jmax,0)=av_slope(1:imax,1:jmax,1)
av_slope(1:imax,1:jmax,k1)=av_slope(1:imax,1:jmax,kmax)
ENDIF
if (avfile.ne.'') then
status2 = nf90_open(avfile, nf90_NoWrite, ncid)
IF (status2/= nf90_noerr) THEN
write(*,*),'file =',avfile
CALL writeerror(703)
ENDIF
status = nf90_inq_varid(ncid, "av_slope",rhVarid)
if (status.eq.nf90_NoErr) then
call check( nf90_inquire_variable(ncid, rhVarid, dimids = dimIDs))
call check( nf90_inquire_dimension(ncid, dimIDs(1), len = size1))
call check( nf90_inquire_dimension(ncid, dimIDs(2), len = size2))
call check( nf90_inquire_dimension(ncid, dimIDs(3), len = size3))
IF(size1.ne.imax.or.size2.ne.jmax*px.or.size3.ne.k1+1) CALL writeerror(704)
call check( nf90_get_var(ncid,rhVarid,av_slope(1:imax,1:jmax,0:k1),
& start=(/1,rank*jmax+1,1/),count=(/imax,jmax,k1+1/)) )
else
write(*,*),'file =',avfile,' variable "av_slope" not found '
CALL writeerror(704)
endif
call check( nf90_close(ncid) )
endif
!IF (rank.eq.0) write(*,*),av_slope,z_rks,av_slope_z(1:n),avalanche_slope(1:n)
END SUBROUTINE init_sediment
SUBROUTINE slipvelocity(csed,Wcfd,wsed,rr,k1b,k1e,sumWkm,dt,dz)
implicit none
INTEGER n,kp ! local variables!
REAL wsed(nfrac,0:i1,0:j1,0:k1)
REAL csed(nfrac,0:i1,0:j1,0:k1),csed2(nfrac,0:i1,0:j1,0:k1)
REAL Wcfd(0:i1,0:j1,0:k1),sumWkm(0:i1,0:j1,0:k1)
REAL rr(0:i1,0:j1,0:k1),rr2(0:i1,0:j1,0:k1)
REAL ctot,sum_c_ws,ws(nfrac),ccc(nfrac)
INTEGER k1b,k1e,kpp,km,iter,kplus,kplus2
REAL dt_dzi,noemer,rrpos,rrneg,limiter,dt,dz,ws_basis(nfrac),W_km_sum
dt_dzi=dt/dz
IF (k1b.eq.1.and.k1e.eq.(k1-1)) THEN ! ws in flow domain is determined
DO n=1,nfrac
ws_basis(n)=frac(n)%ws
ENDDO
ELSE ! ws at bed for ero-depo is determined
DO n=1,nfrac
ws_basis(n)=frac(n)%ws_dep
ENDDO
ENDIF
csed2=cW !csed
IF (applyVOF.eq.1) THEN
call state_edges(cW,rhW)
ENDIF
rr2=rhW
IF (applyVOF.eq.1) THEN
rhW=rho_b
ENDIF
if (nobst>0.or.bedlevelfile.ne.''.or.interaction_bed.ge.4) then ! default C(k=0)=C(k=1); therefore only for cases when bed is not necessarily at k=0 this fix is needed:
DO i=0,i1
DO j=0,j1
kplus = MIN(kbed(i,j)+1,k1)
kplus2 = MIN(kbed(i,j)+2,k1)
!rr2(i,j,kbed(i,j))=rr2(i,j,kplus) ! apply neumann boundary over obstacles to get correct drift flux settling
rr2(i,j,kbed(i,j))=rr(i,j,kplus) ! make rhW(kbed) equal to rho fluid first cell above to get correct drift flux settling
!rr2(i,j,kbed(i,j))=1.5*rr(i,j,kplus)-0.5*rr(i,j,kplus2)
DO n=1,nfrac
!csed2(n,i,j,kbed(i,j))=csed2(n,i,j,kplus) ! apply neumann boundary over obstacles to get correct drift flux settling
csed2(n,i,j,kbed(i,j))=csed(n,i,j,kplus) ! apply neumann boundary over obstacles to get correct drift flux settling
!csed2(n,i,j,kbed(i,j))=1.5*csed(n,i,j,kplus)-0.5*csed(n,i,j,kplus2)
ENDDO
ENDDO
ENDDO
endif
IF (slipvel.eq.2) THEN ! no hindered settling
DO i=0,i1
DO j=0,j1
DO k=k1b,k1e ! at k=0 wsed should be zero at kmax wsed is calculated
ctot=0.
sum_c_ws=0.
kp = MIN(k+1,k1)
DO n=1,nfrac
ccc(n) = csed2(n,i,j,k) !0.5*(csed2(n,i,j,k) + csed2(n,i,j,kp))
ctot=ccc(n)*frac(n)%dfloc/frac(n)%dpart*0.5*(rhocorr_air_z(n,k)+rhocorr_air_z(n,kp))+ctot
ws(n)=ws_basis(n)*wscorr_z(n,k)
! ws is positive downward, wsed is positive upward
sum_c_ws=sum_c_ws+ws(n)*ccc(n)*frac(n)%rho/rr2(i,j,k) !(0.5*(rr2(i,j,k)+rr2(i,j,kp)))
!! ! According to drift velocity literature the drift flux is calculated using the mass-fraction in stead of volume-fraction,
!! ! because sum_c_ws is slipvelocity relative to mixture velocity not relative to volumetric flux!
! therefore an extra frac(n)%rho/rho_mix is included
ENDDO
ctot=MIN(ctot,1.) ! limit on 1, see also Winterwerp 1999 p.46, because Cfloc can be >1
DO n=1,nfrac
wsed(n,i,j,k)=Wcfd(i,j,k)+sum_c_ws-ws(n) ! wsed is positive upward
ENDDO
W_km_sum=0.
do n=1,nfrac
W_km_sum=W_km_sum+ccc(n)*frac(n)%rho*(wsed(n,i,j,k)-Wcfd(i,j,k))
enddo
W_km_sum=W_km_sum+(1.-ctot)*rho_b*sum_c_ws !sum_c_ws=fluid return velocity --> slipvelocity of fluid relative to mixture-velocity
sumWkm(i,j,k)=W_km_sum !NEW 2-10-2018: contains sum of all fractions and fluid phase drift velocity for correct determination driftflux-force
ENDDO
ENDDO
ENDDO
ELSE ! hindered settling
DO i=0,i1
DO j=0,j1
DO k=k1b,k1e ! at k=0 wsed should be zero at kmax wsed is calculated
ctot=0.
sum_c_ws=0.
kp = MIN(k+1,k1)
DO n=1,nfrac
ccc(n) = csed2(n,i,j,k) !0.5*(csed2(n,i,j,k) + csed2(n,i,j,kp))
ctot=ccc(n)*frac(n)%dfloc/frac(n)%dpart*0.5*(rhocorr_air_z(n,k)+rhocorr_air_z(n,kp))+ctot
! for mud Cfloc must be used to calculate Ctot (Cfloc=Ctot*dfloc/dpart)
! ccc(n) is used in drift flux correction sum_c_ws which is calculated with mass concentration based on Ctot (not on Cfloc)
ENDDO
ctot=MIN(ctot,1.) ! limit on 1, see also Winterwerp 1999 p.46, because Cfloc can be >1
DO n=1,nfrac
ws(n)=ws_basis(n)*(1.-ctot)**(frac(n)%n)*wscorr_z(n,k) !frac(n)%n lowered with one already
! ws is positive downward, wsed is positive upward
! Ri_Za is defined with volume concentration ctot
sum_c_ws=sum_c_ws+ws(n)*ccc(n)*frac(n)%rho/rr2(i,j,k) !(0.5*(rr2(i,j,k)+rr2(i,j,kp)))
!! ! According to drift velocity literature the drift flux is calculated using the mass-fraction in stead of volume-fraction,
!! ! because sum_c_ws is slipvelocity relative to mixture velocity not relative to volumetric flux!
! therefore an extra frac(n)%rho/rho_mix is included
ENDDO
DO n=1,nfrac
wsed(n,i,j,k)=Wcfd(i,j,k)+sum_c_ws-ws(n) ! wsed is positive upward
ENDDO
! DO iter=1,1
! ! TVD interpolation csed for correct return current in drift flux:
! ! calculate all again with new TVD interpolation of csed and use first calculation as first guess of direction Wsed
! ctot=0.
! sum_c_ws=0.
! DO n=1,nfrac
! noemer = csed(n,i,j,kp)-csed(n,i,j,k)
! noemer = MAX(ABS(noemer),1.e-6)*sign(1.,noemer)
! rrpos = (csed(n,i,j,k)-csed(n,i,j,km))/noemer
! rrneg = (csed(n,i,j,kpp+1)-csed(n,i,j,kp))/noemer
! if (Wsed(n,i,j,k)>0.) then
! ccc(n)=csed(n,i,j,k )+0.5*limiter(rrpos)*(csed(n,i,j,kp)-csed(n,i,j,k ))*(1.-dt_dzi*ABS(Wsed(n,i,j,k)))
! else
! ccc(n)=csed(n,i,j,kp)+0.5*limiter(rrneg)*(csed(n,i,j,k )-csed(n,i,j,kp))*(1.-dt_dzi*ABS(Wsed(n,i,j,k)))
! endif
! ctot=ccc(n)*frac(n)%dfloc/frac(n)%dpart+ctot
! ENDDO
! ctot=MIN(ctot,1.) ! limit on 1, see also Winterwerp 1999 p.46, because Cfloc can be >1
! DO n=1,nfrac
! ws(n)=frac(n)%ws*(1.-ctot)**(frac(n)%n) !frac(n)%n lowered with one already
! sum_c_ws=sum_c_ws+ws(n)*ccc(n)*frac(n)%rho/(0.5*(rr(i,j,k)+rr(i,j,kp)))
! ENDDO
! DO n=1,nfrac
! wsed(n,i,j,k)=Wcfd(i,j,k)+sum_c_ws-ws(n) ! wsed is positive upward
! ENDDO
! ENDDO
W_km_sum=0.
do n=1,nfrac
W_km_sum=W_km_sum+ccc(n)*frac(n)%rho*(wsed(n,i,j,k)-Wcfd(i,j,k))
enddo
W_km_sum=W_km_sum+(1.-ctot)*rho_b*sum_c_ws !sum_c_ws=fluid return velocity --> slipvelocity of fluid relative to mixture-velocity
sumWkm(i,j,k)=W_km_sum !NEW 2-10-2018: contains sum of all fractions and fluid phase drift velocity for correct determination driftflux-force
ENDDO
ENDDO
ENDDO
ENDIF
IF (k1b.eq.1.and.k1e.eq.(k1-1)) THEN
DO i=0,i1
DO j=0,j1
DO n=1,nfrac
wsed(n,i,j,kbed(i,j))=Wcfd(i,j,kbed(i,j)) ! prevent sediment to flow through the bed or bed-obstacles
wsed(n,i,j,k1)=Wcfd(i,j,k1) ! prevent sediment to flow out of free surface
wsed(n,i,j,k1-1)=-wsed(n,i,j,k1-1)*MIN(0.,frac(n)%ws/(ABS(frac(n)%ws)+1.e-12)) !limit wsed at zero for fractions with downward settling velocity--> no transport of sediment in/out domain at top, but air may escape
ENDDO
ENDDO
ENDDO
Wsed(:,:,:,0)=0.
ENDIF
END SUBROUTINE slipvelocity
SUBROUTINE slipvelocity_bed(csed,Wcfd,wsed,rr,sumWkm,dt,dz)
implicit none
INTEGER n,kp ! local variables!
REAL wsed(nfrac,0:i1,0:j1,0:k1)
REAL csed(nfrac,0:i1,0:j1,0:k1),csed2(nfrac,0:i1,0:j1,0:k1)
REAL Wcfd(0:i1,0:j1,0:k1),sumWkm(0:i1,0:j1,0:k1)
REAL rr(0:i1,0:j1,0:k1),rr2(0:i1,0:j1,0:k1)
REAL ctot,sum_c_ws,ws(nfrac),ccc(nfrac)
INTEGER kpp,km,iter,kplus,kplus2
REAL dt_dzi,noemer,rrpos,rrneg,limiter,dt,dz,ws_basis(nfrac),W_km_sum
dt_dzi=dt/dz
DO n=1,nfrac
ws_basis(n)=frac(n)%ws_dep
ENDDO
csed2=cW !csed
IF (applyVOF.eq.1) THEN
call state_edges(cW,rhW)
ENDIF
rr2=rhW
IF (applyVOF.eq.1) THEN
rhW=rho_b
ENDIF
if (nobst>0.or.bedlevelfile.ne.''.or.interaction_bed.ge.4) then ! default C(k=0)=C(k=1); therefore only for cases when bed is not necessarily at k=0 this fix is needed:
DO i=0,i1
DO j=0,j1
kplus = MIN(kbed(i,j)+1,k1)
kplus2 = MIN(kbed(i,j)+2,k1)
!rr2(i,j,kbed(i,j))=rr2(i,j,kplus) ! apply neumann boundary over obstacles to get correct drift flux settling
rr2(i,j,kbed(i,j))=rr(i,j,kplus) ! make rhW(kbed) equal to rho fluid first cell above to get correct drift flux settling
!rr2(i,j,kbed(i,j))=1.5*rr(i,j,kplus)-0.5*rr(i,j,kplus2)
DO n=1,nfrac
!csed2(n,i,j,kbed(i,j))=csed2(n,i,j,kplus) ! apply neumann boundary over obstacles to get correct drift flux settling
csed2(n,i,j,kbed(i,j))=csed(n,i,j,kplus) ! apply neumann boundary over obstacles to get correct drift flux settling
!csed2(n,i,j,kbed(i,j))=1.5*csed(n,i,j,kplus)-0.5*csed(n,i,j,kplus2)
ENDDO
ENDDO
ENDDO
endif
IF (slipvel.eq.2) THEN ! no hindered settling
DO i=0,i1
DO j=0,j1
k=kbed(i,j)
!DO k=k1b,k1e ! at k=0 wsed should be zero at kmax wsed is calculated
ctot=0.
sum_c_ws=0.
kp = MIN(k+1,k1)
DO n=1,nfrac
ccc(n) = csed2(n,i,j,k) !0.5*(csed2(n,i,j,k) + csed2(n,i,j,kp))
ctot=ccc(n)*frac(n)%dfloc/frac(n)%dpart*0.5*(rhocorr_air_z(n,k)+rhocorr_air_z(n,kp))+ctot
ws(n)=ws_basis(n)*wscorr_z(n,k)
! ws is positive downward, wsed is positive upward
sum_c_ws=sum_c_ws+ws(n)*ccc(n)*frac(n)%rho/rr2(i,j,k) !(0.5*(rr2(i,j,k)+rr2(i,j,kp)))
!! ! According to drift velocity literature the drift flux is calculated using the mass-fraction in stead of volume-fraction,
!! ! because sum_c_ws is slipvelocity relative to mixture velocity not relative to volumetric flux!
! therefore an extra frac(n)%rho/rho_mix is included
ENDDO
ctot=MIN(ctot,1.) ! limit on 1, see also Winterwerp 1999 p.46, because Cfloc can be >1
DO n=1,nfrac
!wsed(n,i,j,k)=Wcfd(i,j,k)+sum_c_ws-ws(n) ! wsed is positive upward
wsed(n,i,j,k)=sum_c_ws-ws(n) ! wsed is positive upward at bed Wcfd is zero
ENDDO
W_km_sum=0.
do n=1,nfrac
!W_km_sum=W_km_sum+ccc(n)*frac(n)%rho*(wsed(n,i,j,k)-Wcfd(i,j,k))
W_km_sum=W_km_sum+ccc(n)*frac(n)%rho*(wsed(n,i,j,k)) !at bed Wcfd is zero
enddo
W_km_sum=W_km_sum+(1.-ctot)*rho_b*sum_c_ws !sum_c_ws=fluid return velocity --> slipvelocity of fluid relative to mixture-velocity
sumWkm(i,j,k)=W_km_sum !NEW 2-10-2018: contains sum of all fractions and fluid phase drift velocity for correct determination driftflux-force
!ENDDO
ENDDO
ENDDO
ELSE ! hindered settling
DO i=0,i1
DO j=0,j1
!DO k=k1b,k1e ! at k=0 wsed should be zero at kmax wsed is calculated
k=kbed(i,j)
ctot=0.
sum_c_ws=0.
kp = MIN(k+1,k1)
DO n=1,nfrac
ccc(n) = csed2(n,i,j,k) !0.5*(csed2(n,i,j,k) + csed2(n,i,j,kp))
ctot=ccc(n)*frac(n)%dfloc/frac(n)%dpart*0.5*(rhocorr_air_z(n,k)+rhocorr_air_z(n,kp))+ctot
! for mud Cfloc must be used to calculate Ctot (Cfloc=Ctot*dfloc/dpart)
! ccc(n) is used in drift flux correction sum_c_ws which is calculated with mass concentration based on Ctot (not on Cfloc)
ENDDO
ctot=MIN(ctot,1.) ! limit on 1, see also Winterwerp 1999 p.46, because Cfloc can be >1
DO n=1,nfrac
ws(n)=ws_basis(n)*(1.-ctot)**(frac(n)%n)*wscorr_z(n,k) !frac(n)%n lowered with one already
! ws is positive downward, wsed is positive upward
! Ri_Za is defined with volume concentration ctot
sum_c_ws=sum_c_ws+ws(n)*ccc(n)*frac(n)%rho/rr2(i,j,k) !(0.5*(rr2(i,j,k)+rr2(i,j,kp)))
!! ! According to drift velocity literature the drift flux is calculated using the mass-fraction in stead of volume-fraction,
!! ! because sum_c_ws is slipvelocity relative to mixture velocity not relative to volumetric flux!
! therefore an extra frac(n)%rho/rho_mix is included
ENDDO
DO n=1,nfrac
!wsed(n,i,j,k)=Wcfd(i,j,k)+sum_c_ws-ws(n) ! wsed is positive upward
wsed(n,i,j,k)=sum_c_ws-ws(n) ! wsed is positive upward at bed Wcfd must be zero
ENDDO
W_km_sum=0.
do n=1,nfrac
!W_km_sum=W_km_sum+ccc(n)*frac(n)%rho*(wsed(n,i,j,k)-Wcfd(i,j,k))
W_km_sum=W_km_sum+ccc(n)*frac(n)%rho*(wsed(n,i,j,k)) !at bed Wcfd is zero
enddo
W_km_sum=W_km_sum+(1.-ctot)*rho_b*sum_c_ws !sum_c_ws=fluid return velocity --> slipvelocity of fluid relative to mixture-velocity
sumWkm(i,j,k)=W_km_sum !NEW 2-10-2018: contains sum of all fractions and fluid phase drift velocity for correct determination driftflux-force
!ENDDO
ENDDO
ENDDO
ENDIF
END SUBROUTINE slipvelocity_bed
subroutine erosion_deposition(ccnew,cbotnew,ucfd,vcfd,wcfd,rcfd,ccfd,cbotcfd,ddt,dz)
implicit none
include 'mpif.h'
integer ierr
real wsed(nfrac,0:i1,0:j1,0:k1),erosion,deposition,uu,vv,absU,z0,ust,yplus,tau !local variables
integer n,tel,kplus,n1,k2,kk,n2
REAL ccnew(nfrac,0:i1,0:j1,0:k1),cbotnew(nfrac,0:i1,0:j1) ! output
REAL ccfd(nfrac,0:i1,0:j1,0:k1),cbotcfd(nfrac,0:i1,0:j1) ! input
REAL ddt,dz ! input
REAL ucfd(0:i1,0:j1,0:k1),vcfd(0:i1,0:j1,0:k1),wcfd(0:i1,0:j1,0:k1),rcfd(0:i1,0:j1,0:k1) ! input
REAL sumWkm(0:i1,0:j1,0:k1) !dummy
REAL cbottot,kn_sed_avg,Mr_avg,tau_e_avg,ctot_firstcel,cbotnewtot,cbotnewtot_pos
REAL PSD_bot_sand_massfrac(nfr_sand),PSD_bed_sand_massfrac(nfr_sand),PSD_sand(0:nfr_sand),factor,d50
REAL cbottot_sand,mbottot_sand,cbedtot,mbedtot_sand,diameter_sand_PSD(0:nfr_sand)
REAL cbedtot_sand,rho_botsand,rho_bedsand,rho_sand,delta,Dstar,Shields_cr,ustc2,TT,phipp,erosion_sand
REAL ws_botsand,ws_bedsand,ws_sand,Re_p,CD,Shields_eff,dubdt,Fd,Fi,Fl,W
REAL erosionf(nfrac),depositionf(nfrac),erosion_avg(nfrac)
REAL zb_all(0:i1,0:j1),maxbedslope(0:i1,0:j1),sl1,sl2,sl3,sl4,sl5,sl6,sl7,sl8
REAL d_cbotnew(nfrac,0:i1,0:j1),dbed,dl,dbed_allowed,dbed_adjust,dz_botlayer,c_adjust,c_adjustA,c_adjustB
REAL*8 cbf(nfrac,0:i1),cbb(nfrac,0:i1),zbf(0:i1),zbb(0:i1),reduced_sed
INTEGER itrgt,jtrgt,nav,n_av,kplus2,kpp
REAL ws_botsand2,rho_botsand2,mbottot_sand2,PSD_bot_sand_massfrac2(nfr_sand),have_avalanched,have_avalanched_tmp,cctot
REAL ccfdtot_firstcel,wsedbed,distance_to_bed,zb_W,gvector,ero_factor
REAL pickup_random(1:imax,1:jmax),vs,ve,Uhor(1:imax,1:jmax,1:kmax),cbed,uu2,vv2,bed_slope,facx,facy,bs_geo
REAL qb,MMM,MME,flux,ucr,qbf(1:nfrac),uuRrel,uuLrel,vvRrel,vvLrel,Shields,absUbl,ve_check,ustbl
integer clock,nnn,k_maxU,itrgt2,jtrgt2
INTEGER, DIMENSION(:), ALLOCATABLE :: seed
INTEGER kbedp(0:i1,0:j1),ibeg,iend,kppE,kppW,kppN,kppS,kbed_new(0:i1,0:j1),k_ust_tau_temp
REAL dzbed_dx,dzbed_dy,dzbed_dn,dzbed_ds,bedslope_angle,bedslope_alpha,Shields_cr_bl,fnorm,dzbed_dl,fcor_slope
REAL ppp(0:i1,0:j1,0:k1),Fix,Fiy,Fix_ref,ddzzE,ddzzW,ddzzS,ddzzN,ustu2,ustv2,ww,c_adjust1,c_adjust2,dz1
REAL dzbed_dl2,dzbed_dl3,vwal_x,vwal_y,dzB,bwal,fluxA,fluxB,fluxC,fluxD,absUU,fcor_pres
REAL zb_avg(0:i1,0:j1),dbedL,dbedR,Shields_cr_den,Shields_cr_num,Shields_cr_den_bl,Shields_cr_num_bl,i_curved
REAL ppp_avg(0:i1,0:j1),ppp_avg2(0:i1,0:j1),absU_rks(1:kmax),pR1,pR2,pL1,pL2,qbx,qby,dz3
erosion=0.
deposition=0.
! IF (nobst>0.or.bedlevelfile.ne.''.or.interaction_bed.eq.4.or.interaction_bed.eq.6) THEN
! call slipvelocity(ccfd,wcfd,wsed,rcfd,0,k1-1,sumWkm,ddt,dz)
! !determine wsed on top of obstacles (actually everywhere in the domain) --> Wsed is not zero on kbed(i,j) in this manner!
! ELSE
! call slipvelocity(ccfd,wcfd,wsed,rcfd,0,0,sumWkm,ddt,dz)
! ENDIF
!call slipvelocity_bed(ccfd,wcfd,wsed,rcfd,sumWkm,ddt,dz)
call slipvelocity_bed(ccfd,0.*wcfd,wsed,rcfd,sumWkm,ddt,dz)
IF (interaction_bed.ge.4) THEN
DO i=1,imax
DO j=1,jmax
zbed(i,j)=REAL(MAX(kbed(i,j)-1,0))*dz+(SUM(cbotcfd(1:nfrac,i,j))+SUM(Clivebed(1:nfrac,i,j,kbed(i,j))))/cfixedbed*dz
ENDDO
ENDDO
call bound_cbot(zbed)
ENDIF
IF (IBMorder.eq.2) THEN
DO i=0,i1
DO j=0,j1
zb_W=zbed(i,j)
! kpp=MIN(CEILING(zb_W/dz+0.5)-1+k_ust_tau,k1) !for k_ust_tau=2 kpp is between 1*dz-2*dz distance from bed
! IF (distance_to_bed<0.5*dz) THEN !first cell too close to bed, therefore use second cell (1-1.5)*dz distance from bed
! kpp=MIN(kpp+1,k1) !kpp is in principle between 1*dz-2*dz distance from bed, but due to this if-statement only 1-1.5 from bed
! distance_to_bed=(REAL(kpp)-0.5)*dz-zb_W
! ENDIF
kpp=MIN(CEILING(zb_W/dz+k_ust_tau),k1) !between 0.5-1.5dz from bed for k_ust_tau=1
!kpp=MIN(CEILING(zb_W/dz+0.5),k1) !kpp is between 0-dz distance from bed
distance_to_bed=(REAL(kpp)-0.5)*dz-zb_W
kbedp(i,j)=kpp
ENDDO
ENDDO
ELSE
DO i=0,i1
DO j=0,j1
kbedp(i,j)=MIN(kbed(i,j)+k_ust_tau,k1)
ENDDO
ENDDO
ENDIF
IF (pres_in_predictor_step.eq.0) THEN
ppp(1:imax,1:jmax,1:kmax)=p
ELSE
ppp(1:imax,1:jmax,1:kmax)=pold+p
ENDIF
call bound_3D(ppp)
pL1=0.
pR1=0.
pL2=0.
pR2=0.
IF (dpdx_ref_j(1)>0) THEN !determine ambient background pressure gradient on left and right edge of CFD domain
IF (rank.eq.0) THEN
do j=dpdx_ref_j(1),dpdx_ref_j(2)
pL1 = pL1+ppp(1,j,MIN(k1,kbed(1,j)+k_ust_tau)) !correction is done with k_ust_tau and not k_ust_tau_temp (which is not known outside do loop, but for Fix_ref this difference is of minor importance)
pR1 = pR1+ppp(imax,j,MIN(k1,kbed(imax,j)+k_ust_tau))
enddo
pL1 = pL1 / DBLE(dpdx_ref_j(2)-dpdx_ref_j(1)+1)
pR1 = pR1 / DBLE(dpdx_ref_j(2)-dpdx_ref_j(1)+1)
ELSEIF (rank.eq.px-1) THEN
do j=dpdx_ref_j(1),dpdx_ref_j(2)
pL2 = pL2+ppp(1,j1-j,MIN(k1,kbed(1,j1-j)+k_ust_tau))
pR2 = pR2+ppp(imax,j1-j,MIN(k1,kbed(imax,j1-j)+k_ust_tau))
enddo
pL2 = pL2 / DBLE(dpdx_ref_j(2)-dpdx_ref_j(1)+1)
pR2 = pR2 / DBLE(dpdx_ref_j(2)-dpdx_ref_j(1)+1)
ENDIF
call mpi_bcast(pL1,1,MPI_REAL8,0,MPI_COMM_WORLD,ierr)
call mpi_bcast(pR1,1,MPI_REAL8,0,MPI_COMM_WORLD,ierr)
call mpi_bcast(pL2,1,MPI_REAL8,px-1,MPI_COMM_WORLD,ierr)
call mpi_bcast(pR2,1,MPI_REAL8,px-1,MPI_COMM_WORLD,ierr)
ENDIF
Fix_ref = 0.5*((pR1+pR2)-(pL1+pL2))/(Rp(imax)-Rp(1))/rho_b
!!! IF(time_n.ge.tstart_morf) THEN
!!! IF (wallmodel_tau_sed.eq.3.or.wallmodel_tau_sed.eq.4.or.wallmodel_tau_sed.eq.8.or.wallmodel_tau_sed.eq.9) THEN
!!!! DO i=1,imax
!!!! DO j=1,jmax
!!!! ppp_avg2(i,j)=0.25*ppp(i-1,j,kbedp(i-1,j))+0.5*ppp(i,j,kbedp(i,j))+0.25*ppp(i+1,j,kbedp(i+1,j))
!!!! ENDDO
!!!! ENDDO
!!!! call bound_cbot(ppp_avg2)
!!!! DO i=1,imax
!!!! DO j=1,jmax
!!!! ppp_avg(i,j)=0.25*ppp_avg2(i,j-1)+0.5*ppp_avg2(i,j)+0.25*ppp_avg2(i,j+1)
!!!! ENDDO
!!!! ENDDO
!!!! call bound_cbot(ppp_avg)
!!! DO i=2,imax-1 !leave inflow and outflow dpdx zero 1,imax
!!! DO j=1,jmax
!!!! Fix = (ppp(i+1,j,kbedp(i+1,j))-ppp(i-1,j,kbedp(i-1,j)))/(Rp(i+1)-Rp(i-1))/rcfd(i,j,kbedp(i,j))
!!!! Fiy = (ppp(i,j+1,kbedp(i,j+1))-ppp(i,j-1,kbedp(i,j-1)))/(Rp(i)*(phip(j+1)-phip(j-1)))/rcfd(i,j,kbedp(i,j))
!!! Fix = (ppp(i+1,j,kbedp(i+1,j))-ppp(i-1,j,kbedp(i-1,j)))/sqrt((Rp(i+1)-Rp(i-1))**2+(zbed(i+1,j)-zbed(i-1,j))**2)
!!! & /rcfd(i,j,kbedp(i,j))
!!! Fiy = (ppp(i,j+1,kbedp(i,j+1))-ppp(i,j-1,kbedp(i,j-1)))/sqrt((Rp(i)*(phip(j+1)-phip(j-1)))**2+(zbed(i,j+1)-zbed(i,j-1))**2)
!!! & /rcfd(i,j,kbedp(i,j))
!!!! Fix = (ppp_avg(i+1,j)-ppp_avg(i-1,j))/sqrt((Rp(i+1)-Rp(i-1))**2+(zbed(i+1,j)-zbed(i-1,j))**2)/rcfd(i,j,kbedp(i,j))
!!!! Fiy = (ppp_avg(i,j+1)-ppp_avg(i,j-1))/sqrt((Rp(i)*(phip(j+1)-phip(j-1)))**2+(zbed(i,j+1)-zbed(i,j-1))**2)/rcfd(i,j,kbedp(i,j))
!!! TBLEsed_dpdx(i,j)=TBLEsed_grad_relax*Fix+(1.-TBLEsed_grad_relax)*TBLEsed_dpdx(i,j)
!!! TBLEsed_dpdy(i,j)=TBLEsed_grad_relax*Fiy+(1.-TBLEsed_grad_relax)*TBLEsed_dpdy(i,j)
!!! ENDDO
!!! ENDDO
!!! ENDIF
!!! ENDIF
IF (wallmodel_tau_sed.eq.11) THEN
IF (mod(istep,ndtbed).eq.0) THEN
telUVWbed=telUVWbed+1
IF (telUVWbed>nrmsbed) telUVWbed=1
ENDIF
ENDIF
d_cbotnew = 0.
IF (interaction_bed.le.3.and.time_n.ge.tstart_morf) THEN ! erosion sedimentation without bed update and for each sediment fraction independently
DO i=1,imax
DO j=1,jmax
IF (k_ust_tau_sed_range(1)>0) THEN
absU_rks = 0.
tel=k_ust_tau_sed_range(1)
DO k=k_ust_tau_sed_range(1),k_ust_tau_sed_range(2)
IF (IBMorder.eq.2) THEN
kppE=MIN(CEILING(0.5*(zbed(i,j)+zbed(i+1,j))/dz+k),k1) !between 0.5-1.5dz from bed for k=1
kppW=MIN(CEILING(0.5*(zbed(i,j)+zbed(i-1,j))/dz+k),k1)
kppN=MIN(CEILING(0.5*(zbed(i,j)+zbed(i,j+1))/dz+k),k1)
kppS=MIN(CEILING(0.5*(zbed(i,j)+zbed(i,j-1))/dz+k),k1)
ddzzE=(REAL(kppE)-0.5)*dz-0.5*(zbed(i,j)+zbed(i+1,j))
ddzzW=(REAL(kppW)-0.5)*dz-0.5*(zbed(i,j)+zbed(i-1,j))
ddzzN=(REAL(kppN)-0.5)*dz-0.5*(zbed(i,j)+zbed(i,j+1))
ddzzS=(REAL(kppS)-0.5)*dz-0.5*(zbed(i,j)+zbed(i,j-1))
!distance_to_bed=MAX(0.5*dz,0.25*(ddzzE+ddzzW+ddzzN+ddzzS))
distance_to_bed=0.25*(ddzzE+ddzzW+ddzzN+ddzzS)
ELSE
kppE = MIN(MAX(kbed(i,j),kbed(i+1,j))+k,k1)
kppW = MIN(MAX(kbed(i,j),kbed(i-1,j))+k,k1)
kppS = MIN(MAX(kbed(i,j),kbed(i,j+1))+k,k1)
kppN = MIN(MAX(kbed(i,j),kbed(i,j-1))+k,k1)
distance_to_bed=(REAL(k)-0.5)*dz
ENDIF
uu=0.5*(ucfd(i,j,kppE)+ucfd(i-1,j,kppW))-Ubot_TSHD(j)
vv=0.5*(vcfd(i,j,kppN)+vcfd(i,j-1,kppS))-Vbot_TSHD(j)
absU_rks(tel)=sqrt(uu**2+vv**2)
tel=tel+1
ENDDO
k_ust_tau_temp = MAXLOC(absU_rks,DIM=1)
ELSE
k_ust_tau_temp = k_ust_tau
ENDIF
IF (wallmodel_tau_sed.eq.3.or.wallmodel_tau_sed.eq.4.or.wallmodel_tau_sed.eq.8.or.wallmodel_tau_sed.eq.9) THEN
Fix = (ppp(i+1,j,MIN(k1,kbed(i+1,j)+k_ust_tau_temp))-ppp(i-1,j,MIN(k1,kbed(i-1,j)+k_ust_tau_temp)))
& /sqrt((Rp(i+1)-Rp(i-1))**2+(zbed(i+1,j)-zbed(i-1,j))**2)/rcfd(i,j,MIN(k1,kbed(i,j)+k_ust_tau_temp))
Fix = Fix - Fix_ref
Fiy = (ppp(i,j+1,MIN(k1,kbed(i,j+1)+k_ust_tau_temp))-ppp(i,j-1,MIN(k1,kbed(i,j-1)+k_ust_tau_temp)))
& /sqrt((Rp(i)*(phip(j+1)-phip(j-1)))**2+(zbed(i,j+1)-zbed(i,j-1))**2)/rcfd(i,j,MIN(k1,kbed(i,j)+k_ust_tau_temp))
TBLEsl_dpdx(i,j)=TBLEsl_grad_relax*Fix+(1.-TBLEsl_grad_relax)*TBLEsl_dpdx(i,j)
TBLEsl_dpdy(i,j)=TBLEsl_grad_relax*Fiy+(1.-TBLEsl_grad_relax)*TBLEsl_dpdy(i,j)
ENDIF
!! First Ubot_TSHD and Vbot_TSHD is subtracted to determine tau
!! only over ambient velocities not over U_TSHD
IF (IBMorder.eq.2) THEN
kppE=MIN(CEILING(0.5*(zbed(i,j)+zbed(i+1,j))/dz+k_ust_tau_temp),k1) !between 0.5-1.5dz from bed
kppW=MIN(CEILING(0.5*(zbed(i,j)+zbed(i-1,j))/dz+k_ust_tau_temp),k1)
kppN=MIN(CEILING(0.5*(zbed(i,j)+zbed(i,j+1))/dz+k_ust_tau_temp),k1)
kppS=MIN(CEILING(0.5*(zbed(i,j)+zbed(i,j-1))/dz+k_ust_tau_temp),k1)
ddzzE=(REAL(kppE)-0.5)*dz-0.5*(zbed(i,j)+zbed(i+1,j))
ddzzW=(REAL(kppW)-0.5)*dz-0.5*(zbed(i,j)+zbed(i-1,j))
ddzzN=(REAL(kppN)-0.5)*dz-0.5*(zbed(i,j)+zbed(i,j+1))
ddzzS=(REAL(kppS)-0.5)*dz-0.5*(zbed(i,j)+zbed(i,j-1))
!distance_to_bed=MAX(0.5*dz,0.25*(ddzzE+ddzzW+ddzzN+ddzzS))
distance_to_bed=0.25*(ddzzE+ddzzW+ddzzN+ddzzS)
! zb_W=zbed(i,j)
! kpp=MIN(CEILING(zb_W/dz+0.5)-1+k_ust_tau_temp,k1) !for k_ust_tau_temp=2 kpp is between 1*dz-2*dz distance from bed
! !kpp=MIN(CEILING(zb_W/dz+0.5),k1) !kpp is between 0-dz distance from bed
! distance_to_bed=(REAL(kpp)-0.5)*dz-zb_W
! !as start use first cell (0.5-1)*dz distance from bed
! IF (distance_to_bed<0.5*dz) THEN !first cell too close to bed, therefore use second cell (1-1.5)*dz distance from bed
! kpp=MIN(kpp+1,k1) !kpp is in principle between 1*dz-2*dz distance from bed, but due to this if-statement only 1-1.5 from bed
! distance_to_bed=(REAL(kpp)-0.5)*dz-zb_W
! ENDIF
ELSE
kppE = MIN(MAX(kbed(i,j),kbed(i+1,j))+k_ust_tau_temp,k1)
kppW = MIN(MAX(kbed(i,j),kbed(i-1,j))+k_ust_tau_temp,k1)
kppS = MIN(MAX(kbed(i,j),kbed(i,j+1))+k_ust_tau_temp,k1)
kppN = MIN(MAX(kbed(i,j),kbed(i,j-1))+k_ust_tau_temp,k1)
distance_to_bed=(REAL(k_ust_tau_temp)-0.5)*dz
!
! kpp = MIN(kbedp(i,j),k1) !for k_ust_tau=2 kpp is 1.5*dz from 0-order ibm bed
! distance_to_bed=(REAL(k_ust_tau)-0.5)*dz
ENDIF
uu=0.5*(ucfd(i,j,kppE)+ucfd(i-1,j,kppW))-Ubot_TSHD(j)
vv=0.5*(vcfd(i,j,kppN)+vcfd(i,j-1,kppS))-Vbot_TSHD(j)
IF (pickup_bedslope_geo.eq.1) THEN
bed_slope = atan((zbed(i+1,j)-zbed(i-1,j))/(Rp(i+1)-Rp(i-1)))
& *MIN(bednotfixed(i+1,j,kbed(i+1,j)),bednotfixed(i-1,j,kbed(i-1,j)))
uu2 = uu*cos(bed_slope)+wcfd(i,j,kbedp(i,j))*sin(bed_slope)
facx = 1./cos(bed_slope)
bed_slope = atan((zbed(i,j+1)-zbed(i,j-1))/(Rp(i)*(phip(j+1)-phip(j-1))))
& *MIN(bednotfixed(i,j+1,kbed(i,j+1)),bednotfixed(i,j-1,kbed(i,j-1)))
vv2 = vv*cos(bed_slope)+wcfd(i,j,kbedp(i,j))*sin(bed_slope)
facy = 1./cos(bed_slope)
bs_geo = facx*facy ! increase in dx and dy (area) over which pickup and deposition take place
absU = sqrt((uu2)**2+(vv2)**2)
absU_sed_relax(i,j) = sl_relax*absU+(1.-sl_relax)*absU_sed_relax(i,j)
absU = absU_sed_relax(i,j)
ELSE
absU=sqrt((uu)**2+(vv)**2)
absU_sed_relax(i,j) = sl_relax*absU+(1.-sl_relax)*absU_sed_relax(i,j)
absU = absU_sed_relax(i,j)
bs_geo = 1.
ENDIF
IF (ABS(z_tau_sed-0.5*dz)>1.e-6) THEN !only if z_tau_sed is user defined (so not 0.5*dz) then do correction below:
ust=0.1*absU
do tel=1,10 ! 10 iter is more than enough
z0=MAX(kn_flow(i,j)/30.+0.11*nu_sediment_pickup/MAX(ust,1.e-9),1e-9)
ust=absU/MAX(1./kappa*log(MAX(distance_to_bed/z0,1.001)),2.) !ust maximal 0.5*absU
enddo
distance_to_bed=z_tau_sed
absU=ust/kappa*log(distance_to_bed/z0) ! replace absU with velocity that is valid at z_tau_sed from bed (user input to make result less dependent of grid resolution)
ENDIF
! when layer above is non-erodible and non-deposition then no sediment tau and no pickup should be calculated, also not for k_ust_tau>1 which might look at flow above mudmat of dz thick
! (which is the case when both bednotfixed and bednotfixed_depo are 0, so when either of those is 1 then it is not an obstacle)
kplus = MIN(k1,kbed(i,j)+1)
absU = absU * MAX(bednotfixed(i,j,kplus),bednotfixed_depo(i,j,kplus))
DO n=1,nfrac
IF (wallmodel_tau_sed.eq.3) THEN
CALL wall_fun_TBL_Ploc(uu/MAX(1.e-6,sqrt(uu**2+vv**2))*absU,vv/MAX(1.e-6,sqrt(uu**2+vv**2))*absU,TBLEsl_dpdx(i,j),
& TBLEsl_dpdy(i,j),rho_b,2.*distance_to_bed,frac(n)%kn_sed,kappa,nu_sediment_pickup,ust_frac_old(n,i,j),ust,nWM)
ELSEIF (wallmodel_tau_sed.eq.4) THEN
CALL wall_fun_TBL2_Ploc(uu/MAX(1.e-6,sqrt(uu**2+vv**2))*absU,vv/MAX(1.e-6,sqrt(uu**2+vv**2))*absU,TBLEsl_dpdx(i,j),
& TBLEsl_dpdy(i,j),rho_b,2.*distance_to_bed,frac(n)%kn_sed,kappa,nu_sediment_pickup,ust_frac_old(n,i,j),ust,nWM)
ELSEIF (wallmodel_tau_sed.eq.5) THEN
CALL wall_fun_VD_Ploc(uu/MAX(1.e-6,sqrt(uu**2+vv**2))*absU,vv/MAX(1.e-6,sqrt(uu**2+vv**2))*absU,
& rho_b,2.*distance_to_bed,frac(n)%kn_sed,kappa,nu_sediment_pickup,ust_frac_old(n,i,j),ust,nWM)
ELSEIF (wallmodel_tau_sed.eq.8) THEN
CALL wall_fun_GWF_Ploc(uu/MAX(1.e-6,sqrt(uu**2+vv**2))*absU,vv/MAX(1.e-6,sqrt(uu**2+vv**2))*absU,TBLEsl_dpdx(i,j),
& TBLEsl_dpdy(i,j),rho_b,2.*distance_to_bed,frac(n)%kn_sed,kappa,nu_sediment_pickup,ust)
ELSEIF (wallmodel_tau_sed.eq.9) THEN
CALL wall_fun_GWF_dpdlfavo_Ploc(uu/MAX(1.e-6,sqrt(uu**2+vv**2))*absU,vv/MAX(1.e-6,sqrt(uu**2+vv**2))*absU,TBLEsl_dpdx(i,j),
& TBLEsl_dpdy(i,j),rho_b,2.*distance_to_bed,frac(n)%kn_sed,kappa,nu_sediment_pickup,ust)
ELSE
ust=0.1*absU
do tel=1,10 ! 10 iter is more than enough
z0=frac(n)%kn_sed/30.+0.11*nu_sediment_pickup/MAX(ust,1.e-9)
! for tau shear on sediment don't use kn (which is result of bed ripples), use frac(n)%kn_sed; it is adviced to use kn_sed=dfloc
ust=absU/MAX(1./kappa*log(distance_to_bed/z0),2.) !ust maximal 0.5*absU
enddo
ENDIF
kplus = MIN(kbed(i,j)+1,k1)
kplus2 = MIN(kbed(i,j)+2,k1)
tau=rho_b*ust*ust
ust_frac_new(n,i,j)=ust
! called before update in correc (in last k3 step of RK3) so Cnewbot and Cnew are used to determine interaction with bed,
! but effect is added to dcdtbot and dcdt to make superposition on other terms already included in dcdt
! dcdtbot contains sediment volume concentration in bed [-] with ghost bed-cells of dz deep
! dcdt contains sediment volume concentration in water-cells [-]
IF (interaction_bed.ge.2) THEN
erosion = frac(n)%M*MAX(0.,(tau/frac(n)%tau_e-1.))*ddt/frac(n)%rho*bs_geo ! m3/m2
ENDIF
! erosion = erosion -MIN(0.,ddt*0.5*(Diffcof(i,j,kplus)+Diffcof(i,j,kplus2))*(ccnew(n,i,j,kplus2)-ccnew(n,i,j,kplus))/dz) !add upward diffusion to erosion flux
IF (interaction_bed.eq.3) THEN
erosion = MIN(erosion,cbotnew(n,i,j)*dz) ! m3/m2, not more material can be eroded than there was in the bed
ENDIF
IF (depo_implicit.eq.1) THEN !determine deposition as sink implicit
kplus = MIN(kbed(i,j)+1,k1)
ccnew(n,i,j,kplus)=(ccnew(n,i,j,kplus)+erosion/dz)/
& (1.-MAX(0.,(1.-tau/frac(n)%tau_d))*MIN(0.,Wsed(n,i,j,kbed(i,j)))*ddt/dz) ! vol conc. [-]
deposition = MAX(0.,(1.-tau/frac(n)%tau_d))*ccnew(n,i,j,kplus)*MIN(0.,Wsed(n,i,j,kbed(i,j)))*ddt ! m !ccfd
cbotnew(n,i,j)=cbotnew(n,i,j)-(erosion+deposition)/(dz) ! vol conc. [-]
ELSE
kplus = MIN(kbed(i,j)+1,k1)
deposition = MAX(0.,(1.-tau/frac(n)%tau_d))*ccnew(n,i,j,kplus)*MIN(0.,Wsed(n,i,j,kbed(i,j)))*ddt ! m !ccfd
ccnew(n,i,j,kplus)=ccnew(n,i,j,kplus)+(erosion+deposition)/(dz) ! vol conc. [-]
cbotnew(n,i,j)=cbotnew(n,i,j)-(erosion+deposition)/(dz) ! vol conc. [-]
ENDIF
ENDDO
ENDDO
ENDDO
ELSEIF((interaction_bed.ge.4).and.time_n.ge.tstart_morf) THEN ! including bedupdate; erosion based on avg sediment properties top layer
IF (pickup_fluctuations.eq.1) THEN
!1 add white noise to pickup
CALL SYSTEM_CLOCK(COUNT=clock)
CALL RANDOM_SEED(size = nnn)
ALLOCATE(seed(nnn))
CALL SYSTEM_CLOCK(COUNT=clock)
seed = clock + 37 * (/ (i - 1, i = 1, nnn) /)
CALL RANDOM_SEED(PUT = seed)
call random_number(pickup_random) ! uniform distribution 0,1
pickup_random=1.+2.*(pickup_random-0.5)*pickup_fluctuations_ampl
ENDIF
IF (cbed_method.eq.2) THEN
Uhor(1:imax,1:jmax,1:kmax)=sqrt((0.5*(ucfd(0:imax-1,1:jmax,1:kmax)+ucfd(1:imax,1:jmax,1:kmax)))**2 +
& (0.5*(vcfd(1:imax,0:jmax-1,1:kmax)+vcfd(1:imax,1:jmax,1:kmax)))**2)
ENDIF
! IF ((periodicx.eq.0.and.monopile<0).and.(interaction_bed.eq.4.or.interaction_bed.eq.6)) THEN
! ibeg=2
! iend=imax-1
! ELSEIF ((periodicx.eq.0.and.monopile>0).and.(interaction_bed.eq.4.or.interaction_bed.eq.6)) THEN
! ibeg=1
! iend=imax-1
! ELSE
ibeg=1
iend=imax
! ENDIF
DO i=ibeg,iend
DO j=1,jmax
erosionf=0.
depositionf=0.
cbotnewtot=0.
cbotnewtot_pos=0.
ctot_firstcel=0.
IF (k_ust_tau_sed_range(1)>0) THEN
absU_rks = 0.
tel=k_ust_tau_sed_range(1)
DO k=k_ust_tau_sed_range(1),k_ust_tau_sed_range(2)
IF (IBMorder.eq.2) THEN
kppE=MIN(CEILING(0.5*(zbed(i,j)+zbed(i+1,j))/dz+k),k1) !between 0.5-1.5dz from bed for k=1
kppW=MIN(CEILING(0.5*(zbed(i,j)+zbed(i-1,j))/dz+k),k1)
kppN=MIN(CEILING(0.5*(zbed(i,j)+zbed(i,j+1))/dz+k),k1)
kppS=MIN(CEILING(0.5*(zbed(i,j)+zbed(i,j-1))/dz+k),k1)
ddzzE=(REAL(kppE)-0.5)*dz-0.5*(zbed(i,j)+zbed(i+1,j))
ddzzW=(REAL(kppW)-0.5)*dz-0.5*(zbed(i,j)+zbed(i-1,j))
ddzzN=(REAL(kppN)-0.5)*dz-0.5*(zbed(i,j)+zbed(i,j+1))
ddzzS=(REAL(kppS)-0.5)*dz-0.5*(zbed(i,j)+zbed(i,j-1))
!distance_to_bed=MAX(0.5*dz,0.25*(ddzzE+ddzzW+ddzzN+ddzzS))
distance_to_bed=0.25*(ddzzE+ddzzW+ddzzN+ddzzS)
ELSE
kppE = MIN(MAX(kbed(i,j),kbed(i+1,j))+k,k1)
kppW = MIN(MAX(kbed(i,j),kbed(i-1,j))+k,k1)
kppS = MIN(MAX(kbed(i,j),kbed(i,j+1))+k,k1)
kppN = MIN(MAX(kbed(i,j),kbed(i,j-1))+k,k1)
distance_to_bed=(REAL(k)-0.5)*dz
ENDIF
uu=0.5*(ucfd(i,j,kppE)+ucfd(i-1,j,kppW))-Ubot_TSHD(j)
vv=0.5*(vcfd(i,j,kppN)+vcfd(i,j-1,kppS))-Vbot_TSHD(j)
absU_rks(tel)=sqrt(uu**2+vv**2)
tel=tel+1
ENDDO
k_ust_tau_temp = MAXLOC(absU_rks,DIM=1)
ELSE
k_ust_tau_temp = k_ust_tau
ENDIF
IF (wallmodel_tau_sed.eq.3.or.wallmodel_tau_sed.eq.4.or.wallmodel_tau_sed.eq.8.or.wallmodel_tau_sed.eq.9) THEN
Fix = (ppp(i+1,j,MIN(k1,kbed(i+1,j)+k_ust_tau_temp))-ppp(i-1,j,MIN(k1,kbed(i-1,j)+k_ust_tau_temp)))
& /sqrt((Rp(i+1)-Rp(i-1))**2+(zbed(i+1,j)-zbed(i-1,j))**2)/rcfd(i,j,MIN(k1,kbed(i,j)+k_ust_tau_temp))
Fix = Fix - Fix_ref
Fiy = (ppp(i,j+1,MIN(k1,kbed(i,j+1)+k_ust_tau_temp))-ppp(i,j-1,MIN(k1,kbed(i,j-1)+k_ust_tau_temp)))
& /sqrt((Rp(i)*(phip(j+1)-phip(j-1)))**2+(zbed(i,j+1)-zbed(i,j-1))**2)/rcfd(i,j,MIN(k1,kbed(i,j)+k_ust_tau_temp))
TBLEsl_dpdx(i,j)=TBLEsl_grad_relax*Fix+(1.-TBLEsl_grad_relax)*TBLEsl_dpdx(i,j)
TBLEsl_dpdy(i,j)=TBLEsl_grad_relax*Fiy+(1.-TBLEsl_grad_relax)*TBLEsl_dpdy(i,j)
TBLEbl_dpdx(i,j)=TBLEbl_grad_relax*Fix+(1.-TBLEbl_grad_relax)*TBLEbl_dpdx(i,j)
TBLEbl_dpdy(i,j)=TBLEbl_grad_relax*Fiy+(1.-TBLEbl_grad_relax)*TBLEbl_dpdy(i,j)
ENDIF
!! First Ubot_TSHD and Vbot_TSHD is subtracted to determine tau
!! only over ambient velocities not over U_TSHD
IF (IBMorder.eq.2) THEN
kppE=MIN(CEILING(0.5*(zbed(i,j)+zbed(i+1,j))/dz+k_ust_tau_temp),k1) !between 0.5-1.5dz from bed
kppW=MIN(CEILING(0.5*(zbed(i,j)+zbed(i-1,j))/dz+k_ust_tau_temp),k1)
kppN=MIN(CEILING(0.5*(zbed(i,j)+zbed(i,j+1))/dz+k_ust_tau_temp),k1)
kppS=MIN(CEILING(0.5*(zbed(i,j)+zbed(i,j-1))/dz+k_ust_tau_temp),k1)
ddzzE=(REAL(kppE)-0.5)*dz-0.5*(zbed(i,j)+zbed(i+1,j))
ddzzW=(REAL(kppW)-0.5)*dz-0.5*(zbed(i,j)+zbed(i-1,j))
ddzzN=(REAL(kppN)-0.5)*dz-0.5*(zbed(i,j)+zbed(i,j+1))
ddzzS=(REAL(kppS)-0.5)*dz-0.5*(zbed(i,j)+zbed(i,j-1))
!distance_to_bed=MAX(0.5*dz,0.25*(ddzzE+ddzzW+ddzzN+ddzzS))
distance_to_bed=0.25*(ddzzE+ddzzW+ddzzN+ddzzS)
! zb_W=zbed(i,j)
! kpp=MIN(CEILING(zb_W/dz+0.5)-1+k_ust_tau_temp,k1) !for k_ust_tau_temp=2 kpp is between 1*dz-2*dz distance from bed
! !kpp=MIN(CEILING(zb_W/dz+0.5),k1) !kpp is between 0-dz distance from bed
! distance_to_bed=(REAL(kpp)-0.5)*dz-zb_W
! !as start use first cell (0.5-1)*dz distance from bed
! IF (distance_to_bed<0.5*dz) THEN !first cell too close to bed, therefore use second cell (1-1.5)*dz distance from bed
! kpp=MIN(kpp+1,k1) !kpp is in principle between 1*dz-2*dz distance from bed, but due to this if-statement only 1-1.5 from bed
! distance_to_bed=(REAL(kpp)-0.5)*dz-zb_W
! ENDIF
ELSE
kppE = MIN(MAX(kbed(i,j),kbed(i+1,j))+k_ust_tau_temp,k1)
kppW = MIN(MAX(kbed(i,j),kbed(i-1,j))+k_ust_tau_temp,k1)
kppS = MIN(MAX(kbed(i,j),kbed(i,j+1))+k_ust_tau_temp,k1)
kppN = MIN(MAX(kbed(i,j),kbed(i,j-1))+k_ust_tau_temp,k1)
distance_to_bed=(REAL(k_ust_tau_temp)-0.5)*dz
!
! kpp = MIN(kbedp(i,j),k1) !for k_ust_tau=2 kpp is 1.5*dz from 0-order ibm bed
! distance_to_bed=(REAL(k_ust_tau)-0.5)*dz
ENDIF
! IF (kbedp(i,j)<kbedp(i+1,j).and.kbedp(i,j)<kbedp(i-1,j)) THEN ! pit
! uu=0.5*(ucfd(i,j,kpp)+ucfd(i-1,j,kpp))-Ubot_TSHD(j)
! ELSE
! uu=0.5*(ucfd(i,j,MAX(kbedp(i,j),kbedp(i+1,j),kpp))+ucfd(i-1,j,MAX(kbedp(i,j),kbedp(i-1,j),kpp)))-Ubot_TSHD(j) !choice to not alter distance_to_bed at slopes
! ENDIF
! IF (kbedp(i,j)<kbedp(i,j+1).and.kbedp(i,j)<kbedp(i,j-1)) THEN ! pit
! vv=0.5*(vcfd(i,j,kpp)+vcfd(i,j-1,kpp))-Vbot_TSHD(j)
! ELSE
! vv=0.5*(vcfd(i,j,MAX(kbedp(i,j),kbedp(i,j+1),kpp))+vcfd(i,j-1,MAX(kbedp(i,j),kbedp(i,j-1),kpp)))-Vbot_TSHD(j) !choice to not alter distance_to_bed at slopes
! ENDIF
uu=0.5*(ucfd(i,j,kppE)+ucfd(i-1,j,kppW))-Ubot_TSHD(j)
vv=0.5*(vcfd(i,j,kppN)+vcfd(i,j-1,kppS))-Vbot_TSHD(j)
IF (pickup_bedslope_geo.eq.1) THEN
! !bedload near bed velocity:
! uuRrel = ucfd(i ,j,MAX(kbedp(i,j),kbedp(i+1,j),kpp))-Ubot_TSHD(j)
! bed_slope = atan((zbed(i+1,j)-zbed(i,j))/(Rp(i+1)-Rp(i)))
! uuRrel = uuRrel*cos(bed_slope)+0.5*(wcfd(i,j,kbedp(i,j))+wcfd(i+1,j,kbedp(i+1,j)))*sin(bed_slope)
! !no correction for pit because uuR from cell i always needs to be same as uuL from i+1 in bedload otherwise interuption and pit may never fill up
! uuLrel = ucfd(i-1,j,MAX(kbedp(i,j),kbedp(i-1,j),kpp))-Ubot_TSHD(j)
! bed_slope = atan((zbed(i,j)-zbed(i-1,j))/(Rp(i)-Rp(i-1)))
! uuLrel = uuLrel*cos(bed_slope)+0.5*(wcfd(i,j,kbedp(i,j))+wcfd(i-1,j,kbedp(i-1,j)))*sin(bed_slope)
! vvRrel = vcfd(i,j ,MAX(kbedp(i,j),kbedp(i,j+1),kpp))-Vbot_TSHD(j)
! bed_slope = atan((zbed(i,j+1)-zbed(i,j))/(Rp(i)*(phip(j+1)-phip(j))))
! vvRrel = vvRrel*cos(bed_slope)+0.5*(wcfd(i,j,kbedp(i,j))+wcfd(i,j+1,kbedp(i,j+1)))*sin(bed_slope)
! vvLrel = vcfd(i,j-1,MAX(kbedp(i,j),kbedp(i,j-1),kpp))-Vbot_TSHD(j)
! bed_slope = atan((zbed(i,j)-zbed(i,j-1))/(Rp(i)*(phip(j)-phip(j-1))))
! vvLrel = vvLrel*cos(bed_slope)+0.5*(wcfd(i,j,kbedp(i,j))+wcfd(i,j-1,kbedp(i,j-1)))*sin(bed_slope)
! uuR_relax(i,j) = bl_relax*uuRrel+(1.-bl_relax)*uuR_relax(i,j) !needed for bedload-fluxes
! uuL_relax(i,j) = bl_relax*uuLrel+(1.-bl_relax)*uuL_relax(i,j)
! vvR_relax(i,j) = bl_relax*vvRrel+(1.-bl_relax)*vvR_relax(i,j)
! vvL_relax(i,j) = bl_relax*vvLrel+(1.-bl_relax)*vvL_relax(i,j)
! absUbl = MAX(sqrt((0.5*(uuR_relax(i,j)+uuL_relax(i,j)))**2+(0.5*(vvR_relax(i,j)+vvL_relax(i,j)))**2),1.e-6)
! uuRrel = uuR_relax(i,j)/absUbl
! uuLrel = uuL_relax(i,j)/absUbl
! vvRrel = vvR_relax(i,j)/absUbl
! vvLrel = vvL_relax(i,j)/absUbl
! !suspension load near bed velocity:
! bed_slope = atan((zbed(i+1,j)-zbed(i-1,j))/(Rp(i+1)-Rp(i-1)))
! uu2 = uu*cos(bed_slope)+wcfd(i,j,kpp)*sin(bed_slope)
! facx = 1./cos(bed_slope)
! bed_slope = atan((zbed(i,j+1)-zbed(i,j-1))/(Rp(i)*(phip(j+1)-phip(j-1))))
! vv2 = vv*cos(bed_slope)+wcfd(i,j,kpp)*sin(bed_slope)
! facy = 1./cos(bed_slope)
! bs_geo = facx*facy ! increase in dx and dy (area) over which pickup and deposition take place
! absU = sqrt((uu2)**2+(vv2)**2)
!bedload near bed velocity:
uuRrel = ucfd(i ,j,kppE)-Ubot_TSHD(j)
bed_slope = atan((zbed(i+1,j)-zbed(i,j))/(Rp(i+1)-Rp(i)))
& *MIN(bednotfixed(i+1,j,kbed(i+1,j)),bednotfixed(i,j,kbed(i,j)))
uuRrel = uuRrel*cos(bed_slope)+0.5*(wcfd(i,j,kbedp(i,j))+wcfd(i+1,j,kbedp(i+1,j)))*sin(bed_slope)
!no correction for pit because uuR from cell i always needs to be same as uuL from i+1 in bedload otherwise interuption and pit may never fill up
uuLrel = ucfd(i-1,j,kppW)-Ubot_TSHD(j)
bed_slope = atan((zbed(i,j)-zbed(i-1,j))/(Rp(i)-Rp(i-1)))
& *MIN(bednotfixed(i,j,kbed(i,j)),bednotfixed(i-1,j,kbed(i-1,j)))
uuLrel = uuLrel*cos(bed_slope)+0.5*(wcfd(i,j,kbedp(i,j))+wcfd(i-1,j,kbedp(i-1,j)))*sin(bed_slope)
vvRrel = vcfd(i,j ,kppN)-Vbot_TSHD(j)
bed_slope = atan((zbed(i,j+1)-zbed(i,j))/(Rp(i)*(phip(j+1)-phip(j))))
& *MIN(bednotfixed(i,j+1,kbed(i,j+1)),bednotfixed(i,j,kbed(i,j)))
vvRrel = vvRrel*cos(bed_slope)+0.5*(wcfd(i,j,kbedp(i,j))+wcfd(i,j+1,kbedp(i,j+1)))*sin(bed_slope)
vvLrel = vcfd(i,j-1,kppS)-Vbot_TSHD(j)
bed_slope = atan((zbed(i,j)-zbed(i,j-1))/(Rp(i)*(phip(j)-phip(j-1))))
& *MIN(bednotfixed(i,j,kbed(i,j)),bednotfixed(i,j-1,kbed(i,j-1)))
vvLrel = vvLrel*cos(bed_slope)+0.5*(wcfd(i,j,kbedp(i,j))+wcfd(i,j-1,kbedp(i,j-1)))*sin(bed_slope)
uuR_relax(i,j) = bl_relax*uuRrel+(1.-bl_relax)*uuR_relax(i,j) !needed for bedload-fluxes
uuL_relax(i,j) = bl_relax*uuLrel+(1.-bl_relax)*uuL_relax(i,j)
vvR_relax(i,j) = bl_relax*vvRrel+(1.-bl_relax)*vvR_relax(i,j)
vvL_relax(i,j) = bl_relax*vvLrel+(1.-bl_relax)*vvL_relax(i,j)
absUbl = MAX(sqrt((0.5*(uuR_relax(i,j)+uuL_relax(i,j)))**2+(0.5*(vvR_relax(i,j)+vvL_relax(i,j)))**2),1.e-6)
uuRrel = uuR_relax(i,j)/absUbl
uuLrel = uuL_relax(i,j)/absUbl
vvRrel = vvR_relax(i,j)/absUbl
vvLrel = vvL_relax(i,j)/absUbl
!suspension load near bed velocity:
bed_slope = atan((zbed(i+1,j)-zbed(i-1,j))/(Rp(i+1)-Rp(i-1)))
& *MIN(bednotfixed(i+1,j,kbed(i+1,j)),bednotfixed(i-1,j,kbed(i-1,j)))
uu2 = uu*cos(bed_slope)+wcfd(i,j,kbedp(i,j))*sin(bed_slope)
facx = 1./cos(bed_slope)
bed_slope = atan((zbed(i,j+1)-zbed(i,j-1))/(Rp(i)*(phip(j+1)-phip(j-1))))
& *MIN(bednotfixed(i,j+1,kbed(i,j+1)),bednotfixed(i,j-1,kbed(i,j-1)))
vv2 = vv*cos(bed_slope)+wcfd(i,j,kbedp(i,j))*sin(bed_slope)
facy = 1./cos(bed_slope)
bs_geo = facx*facy ! increase in dx and dy (area) over which pickup and deposition take place
absU = sqrt((uu2)**2+(vv2)**2)
absU_sed_relax(i,j) = sl_relax*absU+(1.-sl_relax)*absU_sed_relax(i,j)
absU = absU_sed_relax(i,j)
ELSE
!bedload near bed velocity:
! uuRrel = ucfd(i ,j,MAX(kbedp(i,j),kbedp(i+1,j),kpp))-Ubot_TSHD(j) !no correction for pit because uuR from cell i always needs to be same as uuL from i+1 in bedload otherwise interuption and pit may never fill up
! uuLrel = ucfd(i-1,j,MAX(kbedp(i,j),kbedp(i-1,j),kpp))-Ubot_TSHD(j)
! vvRrel = vcfd(i,j ,MAX(kbedp(i,j),kbedp(i,j+1),kpp))-Vbot_TSHD(j)
! vvLrel = vcfd(i,j-1,MAX(kbedp(i,j),kbedp(i,j-1),kpp))-Vbot_TSHD(j)
uuRrel = ucfd(i ,j,kppE)-Ubot_TSHD(j) !no correction for pit because uuR from cell i always needs to be same as uuL from i+1 in bedload otherwise interuption and pit may never fill up
uuLrel = ucfd(i-1,j,kppW)-Ubot_TSHD(j)
vvRrel = vcfd(i,j ,kppN)-Vbot_TSHD(j)
vvLrel = vcfd(i,j-1,kppS)-Vbot_TSHD(j)
uuR_relax(i,j) = bl_relax*uuRrel+(1.-bl_relax)*uuR_relax(i,j) !needed for bedload-fluxes
uuL_relax(i,j) = bl_relax*uuLrel+(1.-bl_relax)*uuL_relax(i,j)
vvR_relax(i,j) = bl_relax*vvRrel+(1.-bl_relax)*vvR_relax(i,j)
vvL_relax(i,j) = bl_relax*vvLrel+(1.-bl_relax)*vvL_relax(i,j)
absUbl = MAX(sqrt((0.5*(uuR_relax(i,j)+uuL_relax(i,j)))**2+(0.5*(vvR_relax(i,j)+vvL_relax(i,j)))**2),1.e-6)
uuRrel = uuR_relax(i,j)/absUbl
uuLrel = uuL_relax(i,j)/absUbl
vvRrel = vvR_relax(i,j)/absUbl
vvLrel = vvL_relax(i,j)/absUbl
!suspension load near bed velocity:
absU=sqrt((uu)**2+(vv)**2)
absU_sed_relax(i,j) = sl_relax*absU+(1.-sl_relax)*absU_sed_relax(i,j)
absU = absU_sed_relax(i,j)
bs_geo = 1.
ENDIF
IF (wallmodel_tau_sed.eq.11) THEN !for 11 use uu and vv not scaled back to z_tau_sed
IF (mod(istep,ndtbed).eq.0) THEN
sigtbed(telUVWbed)=dt
ww=wcfd(i,j,kbedp(i,j)) !need to adjust ww for pickup_bedslope_geo.eq.1 in a later stage when this approach proves to be usefull
sigUWbed(telUVWbed,i,j)=dt*uu*ww
sigVWbed(telUVWbed,i,j)=dt*vv*ww
sigUbed(telUVWbed,i,j)=dt*uu
sigVbed(telUVWbed,i,j)=dt*vv
sigWbed(telUVWbed,i,j)=dt*ww
ENDIF
ENDIF
IF (ABS(z_tau_sed-0.5*dz)>1.e-6) THEN !only if z_tau_sed is user defined (so not 0.5*dz) then do correction below:
ust=0.1*absU
do tel=1,10 ! 10 iter is more than enough
z0=MAX(kn_flow(i,j)/30.+0.11*nu_sediment_pickup/MAX(ust,1.e-9),1e-9)
ust=absU/MAX(1./kappa*log(MAX(distance_to_bed/z0,1.001)),2.) !ust maximal 0.5*absU
enddo
ustbl=0.1*absUbl
do tel=1,10 ! 10 iter is more than enough
z0=MAX(kn_flow(i,j)/30.+0.11*nu_sediment_pickup/MAX(ust,1.e-9),1e-9)
ustbl=absUbl/MAX(1./kappa*log(MAX(distance_to_bed/z0,1.001)),2.) !ust maximal 0.5*absU
enddo
distance_to_bed=z_tau_sed
absU=ust/kappa*log(distance_to_bed/z0) ! replace absU with velocity that is valid at z_tau_sed from bed (user input to make result less dependent of grid resolution)
absUbl=ustbl/kappa*log(distance_to_bed/z0) ! replace absUbl with velocity that is valid at z_tau_sed from bed (user input to make result less dependent of grid resolution)
ENDIF
! when layer above is non-erodible and non-deposition then no sediment tau and no pickup should be calculated, also not for k_ust_tau>1 which might look at flow above mudmat of dz thick
! (which is the case when both bednotfixed and bednotfixed_depo are 0, so when either of those is 1 then it is not an obstacle)
kplus = MIN(k1,kbed(i,j)+1)
absU = absU * MAX(bednotfixed(i,j,kplus),bednotfixed_depo(i,j,kplus))
absUbl = absUbl * MAX(bednotfixed(i,j,kplus),bednotfixed_depo(i,j,kplus))
cbottot=0.
cbedtot=0.
IF (nfr_silt>0) THEN
!! 1 determine erosion/sedimentation of mixture of all silt fractions
DO n1=1,nfr_silt
n=nfrac_silt(n1)
cbottot=cbottot+MAX(cbotnew(n,i,j),0.)