-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathsolve.f
executable file
·7595 lines (7116 loc) · 285 KB
/
solve.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/>.
subroutine adamsb2(ib,ie,jb,je,kb,ke)
c
c Performes time integration with second order
c Adams-Bashforth scheme, i.e
c
c
c n+1 n
c dU U - U n
c ---- = ------------ = 1.5*( -ADV + DIFF + Force) -
c dt dt n-1
c 0.5*( -ADV + DIFF + Force)
c
c This scheme is weakly instabel for pure advection,
c and therefore a very small amount of physical diffusion
c is necessary.
c The timestep is limited with CFL=0.3 (see routine chkdt)
c
c
USE nlist
USE sediment
implicit none
! include 'param.txt'
! include 'common.txt'
integer ib,ie,jb,je,kb,ke,n,t,me
real doldc(nfrac,0:i1,0:j1,0:k1),dnewc(nfrac,0:i1,0:j1,0:k1)
real dold(0:i1,0:j1,0:k1),dnew(0:i1,0:j1,0:k1)
real wsed(nfrac,0:i1,0:j1,0:k1),sumWkm(0:i1,0:j1,0:k1)
real*8 pplus(imax,kmax)
! real dUdt1(0:i1,0:j1,0:k1),dVdt1(0:i1,0:j1,0:k1),dWdt1(0:i1,0:j1,0:k1)
!
! dUdt1=dudt
! dVdt1=dvdt
! dWdt1=dwdt
! real Diffcof(0:i1,0:j1,0:k1)
! Diffcof=ekm/Sc/Rnew
me = momentum_exchange_obstacles !100 or 101 or 111 means no advection momentum exchange with bed-cells with HYB6
c********************************************************************
c CALCULATE slipvelocity
c********************************************************************
if (nfrac>0) then
if (slipvel.eq.1.or.slipvel.eq.2) then
call slipvelocity(cnew,Wnew,wsed,rnew,1,kmax,sumWkm,dt,dz)
!! Set boundary conditions jet in:
do t=1,tmax_inPpunt
i=i_inPpunt(t)
j=j_inPpunt(t)
do n=1,nfrac
do k=kmax,kmax !kmax-kjet,kmax
Wsed(n,i,j,k)=Wnew(i,j,k)
enddo
enddo
if (LOA>0.and.outflow_overflow_down.eq.1) then
do n=1,nfrac
do k=kmax-kjet+1,kmax
Wsed(n,i,j,k)=Wnew(i,j,k)
enddo
enddo
endif
enddo
else
do n=1,nfrac
wsed(n,:,:,:)=wnew
enddo
endif
c********************************************************************
c CALCULATE advection, diffusion Concentration
c********************************************************************
do n=1,nfrac
call advecc_VLE(dnewc(n,:,:,:),cnew(n,:,:,:),Unew,Vnew,Wsed(n,:,:,:),rnew,Ru,Rp,dr,phiv,phipt,dz,
+ i1,j1,k1,ib,ie,jb,je,kb,ke,dt,rank,px,periodicx,periodicy,transporteq_fracs,kbed)
call diffc_CDS2 (dnewc(n,:,:,:),Cnew(n,:,:,:),Diffcof,
+ ib,ie,jb,je,kb,ke)
dcdt(n,:,:,:) =cnew(n,:,:,:) + dt*(1.5*dnewc(n,:,:,:)-0.5*cc(n,:,:,:)) !Adams-Bashford
call bound_c(dcdt(n,:,:,:),frac(n)%c,n,dt) ! bc for intermediate dcdt AB2
cc(n,:,:,:)=dnewc(n,:,:,:)
enddo
call state(dcdt,drdt) ! determine drdt with intermediate dcdt
endif
c********************************************************************
c CALCULATE advection, diffusion and Force U-velocity
c********************************************************************
if (diffusion.eq.'COM4') then
call stress_terms(Unew,Vnew,Wnew,ib,ie,jb,je,kb,ke)
endif
if (convection.eq.'CDS2') then
call advecu_CDS2(dnew,Unew,Vnew,Wnew,Rnew,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px)
elseif(convection.eq.'CDS6') then
call advecu_CDS6(dnew,Unew,Vnew,Wnew,Rnew,rhU,rhV,rhW,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,numdiff,
& periodicx,periodicy,wf,wd,numdiff2)
elseif(convection.eq.'COM4') then
call advecu_COM4(dnew,Unew,Vnew,Wnew,Rnew,Ru,Rp,dr,phivt,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px)
elseif(convection.eq.'CDS4') then
call advecu_CDS4(dnew,Unew,Vnew,Wnew,Rnew,rhU,rhV,rhW,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,numdiff,
& periodicx,periodicy,wf,wd,numdiff2)
elseif(convection.eq.'HYB6') then
call advecu_HYB6(dnew,Unew,Vnew,Wnew,Rnew,rhU,rhV,rhW,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,
& numdiff,periodicx,periodicy,kbed,wf,wd,numdiff2,fc_global,me)
elseif(convection.eq.'HYB4') then
call advecu_HYB4(dnew,Unew,Vnew,Wnew,Rnew,rhU,rhV,rhW,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,
& numdiff,periodicx,periodicy,kbed,wf,wd,numdiff2)
elseif(convection.eq.'C4A6') then
call advecu_C4A6(dnew,Unew,Vnew,Wnew,Rnew,Ru,Rp,dr,phivt,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,numdiff,periodicx,periodicy
&,phiv,wf,wd,numdiff2)
elseif(convection.eq.'uTVD') then
call advecu_TVD(dnew,Unew,Vnew,Wnew,Rnew,Ru,Rp,dr,phip,phiv,phipt,phivt,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,dt
& ,periodicx,periodicy,istep)
endif
if (diffusion.eq.'CDS2') then
call diffu_CDS2 (dnew,Unew,Vnew,Wnew,ib,ie,jb,je,kb,ke)
elseif(diffusion.eq.'COM4') then
call diffu_com4(dnew,Srr,Spr,Szr,Spp,Ru,Rp,dr,phipt,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px)
endif
do k=1,kmax
do j=1,jmax
do i=1,imax
dUdt(i,j,k)=Unew(i,j,k)*0.5*(rnew(i,j,k)+rnew(i+1,j,k))+
1 dt*( 1.5*dnew(i,j,k)-
1 0.5*wx(i,j,k)
1 )
1 -(gx*cos_u(j)+gy*sin_u(j))*dt*(0.75*(rnew(i,j,k)+rnew(i+1,j,k))-0.25*(rold(i,j,k)+rold(i+1,j,k))-rho_b)
1 +Ppropx(i,j,k)*dt
wx(i,j,k)=dnew(i,j,k)
enddo
enddo
enddo
c********************************************************************
c CALCULATE advection, diffusion and Force V-velocity
c********************************************************************
if (convection.eq.'CDS2') then
call advecv_CDS2(dnew,Unew,Vnew,Wnew,Rnew,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px)
elseif(convection.eq.'CDS6') then
call advecv_CDS6(dnew,Unew,Vnew,Wnew,Rnew,rhU,rhV,rhW,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,numdiff,
& periodicx,periodicy,wf,wd)
elseif(convection.eq.'COM4') then
call advecv_COM4(dnew,Unew,Vnew,Wnew,Rnew,Ru,Rp,dr,phivt,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px)
elseif(convection.eq.'CDS4') then
call advecv_CDS4(dnew,Unew,Vnew,Wnew,Rnew,rhU,rhV,rhW,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,numdiff,
& periodicx,periodicy,wf,wd)
elseif(convection.eq.'HYB6') then
call advecv_HYB6(dnew,Unew,Vnew,Wnew,Rnew,rhU,rhV,rhW,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,
& numdiff,periodicx,periodicy,kbed,wf,wd,fc_global,me)
elseif(convection.eq.'HYB4') then
call advecv_HYB4(dnew,Unew,Vnew,Wnew,Rnew,rhU,rhV,rhW,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,
& numdiff,periodicx,periodicy,kbed,wf,wd)
elseif(convection.eq.'C4A6') then
call advecv_C4A6(dnew,Unew,Vnew,Wnew,Rnew,Ru,Rp,dr,phivt,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,numdiff,periodicx,periodicy
&,phiv,wf,wd)
elseif(convection.eq.'uTVD') then
call advecv_TVD(dnew,Unew,Vnew,Wnew,Rnew,Ru,Rp,dr,phip,phiv,phipt,phivt,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,dt
& ,periodicx,periodicy,istep)
endif
if (diffusion.eq.'CDS2') then
call diffv_CDS2 (dnew,Unew,Vnew,Wnew,ib,ie,jb,je,kb,ke)
elseif(diffusion.eq.'COM4') then
call diffv_com4(dnew,Spr,Spp,Spz,Ru,Rp,dr,phipt,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px)
endif
do k=1,kmax
do j=1,jmax
do i=1,imax
dVdt(i,j,k)=Vnew(i,j,k)*0.5*(rnew(i,j,k)+rnew(i,j+1,k))+
1 dt*( 1.5*dnew(i,j,k)-
1 0.5*wy(i,j,k)
1 )
1 -(-gx*sin_v(j)+gy*cos_v(j))*dt*(0.75*(rnew(i,j,k)+rnew(i,j+1,k))-0.25*(rold(i,j,k)+rold(i,j+1,k))-rho_b)
! 1 +Ppropy(i,j,k)/(MAX(1.e-2,ABS(1.5*Vnew(i,j,k)-0.5*Vold(i,j,k))))*dt
1 +Ppropy(i,j,k)*dt
wy(i,j,k)=dnew(i,j,k)
enddo
enddo
enddo
c********************************************************************
c CALCULATE advection, diffusion and Force W-velocity
c********************************************************************
if (convection.eq.'CDS2') then
call advecw_CDS2(dnew,Unew,Vnew,Wnew,Rnew,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px)
elseif(convection.eq.'CDS6') then
call advecw_CDS6(dnew,Unew,Vnew,Wnew,Rnew,rhU,rhV,rhW,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,numdiff,
& periodicx,periodicy,wf,wd)
elseif(convection.eq.'COM4') then
call advecw_COM4(dnew,Unew,Vnew,Wnew,Rnew,Ru,Rp,dr,phivt,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px)
elseif(convection.eq.'CDS4') then
call advecw_CDS4(dnew,Unew,Vnew,Wnew,Rnew,rhU,rhV,rhW,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,numdiff,
& periodicx,periodicy,wf,wd)
elseif(convection.eq.'HYB6') then
call advecw_HYB6(dnew,Unew,Vnew,Wnew,Rnew,rhU,rhV,rhW,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,
& numdiff,periodicx,periodicy,kbed,wf,wd,fc_global,me)
elseif(convection.eq.'HYB4') then
call advecw_HYB4(dnew,Unew,Vnew,Wnew,Rnew,rhU,rhV,rhW,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,
& numdiff,periodicx,periodicy,kbed,wf,wd)
elseif(convection.eq.'C4A6') then
call advecw_C4A6(dnew,Unew,Vnew,Wnew,Rnew,Ru,Rp,dr,phivt,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,numdiff,periodicx,periodicy
&,phiv,wf,wd)
elseif(convection.eq.'uTVD') then
call advecw_TVD(dnew,Unew,Vnew,Wnew,Rnew,Ru,Rp,dr,phip,phiv,phipt,phivt,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,dt
& ,periodicx,periodicy,istep)
endif
if (diffusion.eq.'CDS2') then
call diffw_CDS2 (dnew,Unew,Vnew,Wnew,ib,ie,jb,je,kb,ke)
elseif(diffusion.eq.'COM4') then
call diffw_com4(dnew,Szr,Spz,Szz,Ru,Rp,dr,phipt,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px)
endif
if ((slipvel.eq.1.or.slipvel.eq.2).and.nfrac>0) then
call advecw_driftfluxCDS2(dnew,driftfluxforce_calfac*sumWkm,
& rnew,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke)
endif
do k=1,kmax
do j=1,jmax
do i=1,imax
dWdt(i,j,k)=Wnew(i,j,k)*0.5*(rnew(i,j,k)+rnew(i,j,k+1))+
1 dt*(1.5*dnew(i,j,k)-0.5*wz(i,j,k)) -gz*dt*(0.75*(rnew(i,j,k)+rnew(i,j,k+1))-0.25*(rold(i,j,k)+rold(i,j,k+1))-rho_b)
1 +Ppropz(i,j,k)*dt
wz(i,j,k)=dnew(i,j,k)
enddo
enddo
enddo
pold2=pold !save pressure previous timestep [pold2 = timestep before pold]
pold=p+pold !what is called p here was dp in reality, now p is
call pshiftb(pold,pplus) !,rank,imax,jmax,kmax,px)
c********************************************************************
c CALCULATE pressure-gradient with old pressure:
c********************************************************************
do k=1,kmax
do j=1,jmax
do i=1,imax-1
dUdt(i,j,k) = dUdt(i,j,k) -dt * ( pold(i+1,j,k) - pold(i,j,k) ) /( Rp(i+1) - Rp(i) )
enddo
enddo
enddo
!!! do nothing on imax and i=1, so p(imax+1)=p(imax) --> dp/dn=0 at outflow and inflow
if (periodicx.eq.1) then
!!! !periodic bc in x:
i = imax
do k=1,kmax
do j=1,jmax
dUdt(i,j,k) = dUdt(i,j,k) -dt * ( pold(1,j,k) - pold(i,j,k) ) /( Rp(i+1) - Rp(i) )
enddo
enddo
endif
do k=1,kmax
do j=1,jmax-1
do i=1,imax
dVdt(i,j,k) = dVdt(i,j,k) -dt * ( pold(i,j+1,k) - pold(i,j,k) ) /( Rp(i) * (phip(j+1)-phip(j)) )
enddo
enddo
enddo
if (periodicy.eq.1) then
do k=1,kmax
do i=1,imax
dVdt(i,jmax,k) = dVdt(i,jmax,k) -dt * ( pplus(i,k) - pold(i,jmax,k) ) /( Rp(i) * (phip(jmax+1)-phip(jmax)) )
enddo
enddo
else
if (rank.ne.px-1) then
do k=1,kmax
do i=1,imax
dVdt(i,jmax,k) = dVdt(i,jmax,k) -dt * ( pplus(i,k) - pold(i,jmax,k) ) /( Rp(i) * (phip(jmax+1)-phip(jmax)) )
enddo
enddo
endif
endif
do k=1,kmax-1 ! rij kmax wordt niet meegenomen, dus is p(k+1)=p(k), dus dp/dn=0 bij wateroppervlak
do j=1,jmax
do i=1,imax
dWdt(i,j,k) = dWdt(i,j,k) -dt * ( pold(i,j,k+1) - pold(i,j,k) ) / dz
enddo
enddo
enddo
return
end
subroutine adamsb3(ib,ie,jb,je,kb,ke)
c
c Performes time integration with second order
c Adams-Bashforth-3 scheme, i.e
c
c
c n+1 n
c dU U - U n
c ---- = ------------ = 23/12*( -ADV + DIFF + Force)
c dt dt n-1
c -4/3 *( -ADV + DIFF + Force)
c n-2
c +5/12*( -ADV + DIFF + Force)
c
c The timestep is limited with CFL=0.3-0.6 (see routine chkdt)
c
c
USE nlist
USE sediment
implicit none
! include 'param.txt'
! include 'common.txt'
integer ib,ie,jb,je,kb,ke,n,t,me
real doldc(nfrac,0:i1,0:j1,0:k1),dnewc(nfrac,0:i1,0:j1,0:k1)
real dold(0:i1,0:j1,0:k1),dnew(0:i1,0:j1,0:k1)
real wsed(nfrac,0:i1,0:j1,0:k1),sumWkm(0:i1,0:j1,0:k1)
real*8 pplus(imax,kmax)
real dnewcbot(nfrac,0:i1,0:j1)
! real dUdt1(0:i1,0:j1,0:k1),dVdt1(0:i1,0:j1,0:k1),dWdt1(0:i1,0:j1,0:k1)
!
! dUdt1=dudt
! dVdt1=dvdt
! dWdt1=dwdt
! real Diffcof(0:i1,0:j1,0:k1)
! Diffcof=ekm/Sc/Rnew
me = momentum_exchange_obstacles !100 or 101 or 111 means no advection momentum exchange with bed-cells with HYB6
c********************************************************************
c CALCULATE slipvelocity
c********************************************************************
c********************************************************************
c CALCULATE advection, diffusion Concentration
c********************************************************************
if (nfrac>0) then
if (slipvel.eq.1.or.slipvel.eq.2) then
call slipvelocity(cnew,Wnew,wsed,rnew,1,kmax,sumWkm,dt,dz)
!! Set boundary conditions jet in:
do t=1,tmax_inPpunt
i=i_inPpunt(t)
j=j_inPpunt(t)
do n=1,nfrac
do k=kmax,kmax !-kjet,kmax
Wsed(n,i,j,k)=Wnew(i,j,k)
enddo
enddo
if (LOA>0.and.outflow_overflow_down.eq.1) then
do n=1,nfrac
do k=kmax-kjet+1,kmax
Wsed(n,i,j,k)=Wnew(i,j,k)
enddo
enddo
endif
enddo
else
do n=1,nfrac
wsed(n,:,:,:)=wnew
enddo
endif
do n=1,nfrac
call advecc_VLE(dnewc(n,:,:,:),cnew(n,:,:,:),Unew,Vnew,Wsed(n,:,:,:),rnew,Ru,Rp,dr,phiv,phipt,dz,
+ i1,j1,k1,ib,ie,jb,je,kb,ke,dt,rank,px,periodicx,periodicy,transporteq_fracs,kbed)
call diffc_CDS2 (dnewc(n,:,:,:),Cnew(n,:,:,:),Diffcof,
+ ib,ie,jb,je,kb,ke)
dcdt(n,:,:,:) =cnew(n,:,:,:) + dt*(1.5*dnewc(n,:,:,:)-0.5*cc(n,:,:,:)) !AB2
cc(n,:,:,:)=dnewc(n,:,:,:)
enddo
call slipvelocity_bed(cnew,wnew,wsed,rnew,sumWkm,dt,dz) !driftflux_force must be calculated with settling velocity at bed
if (interaction_bed>0) then
CALL erosion_deposition(dcdt,dcdtbot,unew,vnew,wnew,rnew,cnew,cnewbot,dt,dz) !first two vars are adjusted
do n=1,nfrac
call bound_cbot(dcdtbot(n,:,:))
enddo
!! advec concentration in bed with velocity TSHD:
IF (ABS(U_TSHD)>1e-6) THEN
do n=1,nfrac
call adveccbot_TVD(dnewcbot(n,:,:),dcdtbot(n,:,:),Ubot_TSHD,Vbot_TSHD,Ru,Rp,dr,phiv,phipt,dz,
+ i1,j1,ib,ie,jb,je,dt,rank,px,periodicx,periodicy)
dcdtbot(n,:,:)= dcdtbot(n,:,:) + dt*dnewcbot(n,:,:) ! time update
call bound_cbot(dcdtbot(n,:,:))
enddo
ENDIF
!! for kn1: erosion_deposition must be after advecc_TVD and after dcdt update, because cnew and dcdt are two separate vars
if (interaction_bed.eq.4.or.interaction_bed.eq.6) then
call advec_update_Clivebed(dcdt,dcdtbot,dt)
endif
do n=1,nfrac
call bound_cbot(dcdtbot(n,:,:))
enddo
endif
do n=1,nfrac
call bound_c(dcdt(n,:,:,:),frac(n)%c,n,dt) ! bc after erosion_deposition AB3
enddo
call state(dcdt,drdt) ! determine drdt with intermediate dcdt
endif
if (diffusion.eq.'COM4') then
call stress_terms(Unew,Vnew,Wnew,ib,ie,jb,je,kb,ke)
endif
c********************************************************************
c CALCULATE advection, diffusion and Force U-velocity
c********************************************************************
if (convection.eq.'CDS2') then
call advecu_CDS2(dnew,Unew,Vnew,Wnew,Rnew,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px)
elseif(convection.eq.'CDS6') then
call advecu_CDS6(dnew,Unew,Vnew,Wnew,Rnew,rhU,rhV,rhW,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,numdiff
& ,periodicx,periodicy,wf,wd,numdiff2)
elseif(convection.eq.'COM4') then
call advecu_COM4(dnew,Unew,Vnew,Wnew,Rnew,Ru,Rp,dr,phivt,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px)
elseif(convection.eq.'CDS4') then
call advecu_CDS4(dnew,Unew,Vnew,Wnew,Rnew,rhU,rhV,rhW,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,numdiff
& ,periodicx,periodicy,wf,wd,numdiff2)
elseif(convection.eq.'HYB6') then
call advecu_HYB6(dnew,Unew,Vnew,Wnew,Rnew,rhU,rhV,rhW,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,
& numdiff,periodicx,periodicy,kbed,wf,wd,numdiff2,fc_global,me)
elseif(convection.eq.'HYB4') then
call advecu_HYB4(dnew,Unew,Vnew,Wnew,Rnew,rhU,rhV,rhW,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,
& numdiff,periodicx,periodicy,kbed,wf,wd,numdiff2)
elseif(convection.eq.'C4A6') then
call advecu_C4A6(dnew,Unew,Vnew,Wnew,Rnew,Ru,Rp,dr,phivt,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,numdiff,periodicx,periodicy
&,phiv,wf,wd,numdiff2)
elseif(convection.eq.'uTVD') then
call advecu_TVD(dnew,Unew,Vnew,Wnew,Rnew,Ru,Rp,dr,phip,phiv,phipt,phivt,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,dt
& ,periodicx,periodicy,istep)
endif
if (diffusion.eq.'CDS2') then
call diffu_CDS2 (dnew,Unew,Vnew,Wnew,ib,ie,jb,je,kb,ke)
elseif(diffusion.eq.'COM4') then
call diffu_com4(dnew,Srr,Spr,Szr,Spp,Ru,Rp,dr,phipt,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px)
endif
do k=1,kmax
do j=1,jmax
do i=1,imax
dnew(i,j,k)=dnew(i,j,k)-(gx*cos_u(j)+gy*sin_u(j))*(0.5*(rnew(i,j,k)+rnew(i+1,j,k))-rho_b)
! 1 +Ppropx(i,j,k)/(MAX(1.e-2,(ABS(Unew(i,j,k)))))*dt
1 +Ppropx(i,j,k)
dUdt(i,j,k)=Unew(i,j,k)*0.5*(rnew(i,j,k)+rnew(i+1,j,k))+
1 dt*( 23./12.*dnew(i,j,k)-
1 4./3. *wx(i,j,k)+
1 5./12.*wxold(i,j,k))
! 1 -(gx*cos_u(j)+gy*sin_u(j))*dt*(0.75*(rnew(i,j,k)+rnew(i+1,j,k))-0.25*(rold(i,j,k)+rold(i+1,j,k))-rho_b)
wxold(i,j,k) = wx(i,j,k)
wx(i,j,k) = dnew(i,j,k)
enddo
enddo
enddo
c********************************************************************
c CALCULATE advection, diffusion and Force V-velocity
c********************************************************************
if (convection.eq.'CDS2') then
call advecv_CDS2(dnew,Unew,Vnew,Wnew,Rnew,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px)
elseif(convection.eq.'CDS6') then
call advecv_CDS6(dnew,Unew,Vnew,Wnew,Rnew,rhU,rhV,rhW,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,numdiff,
& periodicx,periodicy,wf,wd)
elseif(convection.eq.'COM4') then
call advecv_COM4(dnew,Unew,Vnew,Wnew,Rnew,Ru,Rp,dr,phivt,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px)
elseif(convection.eq.'CDS4') then
call advecv_CDS4(dnew,Unew,Vnew,Wnew,Rnew,rhU,rhV,rhW,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,numdiff,
& periodicx,periodicy,wf,wd)
elseif(convection.eq.'HYB6') then
call advecv_HYB6(dnew,Unew,Vnew,Wnew,Rnew,rhU,rhV,rhW,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,
& numdiff,periodicx,periodicy,kbed,wf,wd,fc_global,me)
elseif(convection.eq.'HYB4') then
call advecv_HYB4(dnew,Unew,Vnew,Wnew,Rnew,rhU,rhV,rhW,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,
& numdiff,periodicx,periodicy,kbed,wf,wd)
elseif(convection.eq.'C4A6') then
call advecv_C4A6(dnew,Unew,Vnew,Wnew,Rnew,Ru,Rp,dr,phivt,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,numdiff,periodicx,periodicy
&,phiv,wf,wd)
elseif(convection.eq.'uTVD') then
call advecv_TVD(dnew,Unew,Vnew,Wnew,Rnew,Ru,Rp,dr,phip,phiv,phipt,phivt,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,dt
& ,periodicx,periodicy,istep)
endif
if (diffusion.eq.'CDS2') then
call diffv_CDS2 (dnew,Unew,Vnew,Wnew,ib,ie,jb,je,kb,ke)
elseif(diffusion.eq.'COM4') then
call diffv_com4(dnew,Spr,Spp,Spz,Ru,Rp,dr,phipt,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px)
endif
do k=1,kmax
do j=1,jmax
do i=1,imax
dnew(i,j,k)=dnew(i,j,k)-(-gx*sin_v(j)+gy*cos_v(j))*(0.5*(rnew(i,j,k)+rnew(i,j+1,k))-rho_b)
1 +Ppropy(i,j,k)
dVdt(i,j,k)=Vnew(i,j,k)*0.5*(rnew(i,j,k)+rnew(i,j+1,k))+
1 dt*( 23./12.*dnew(i,j,k)-
1 4./3. *wy(i,j,k)+
1 5./12.*wyold(i,j,k))
! 1 -(-gx*sin_v(j)+gy*cos_v(j))*dt*(0.75*(rnew(i,j,k)+rnew(i,j+1,k))-0.25*(rold(i,j,k)+rold(i,j+1,k))-rho_b)
wyold(i,j,k) = wy(i,j,k)
wy(i,j,k) = dnew(i,j,k)
enddo
enddo
enddo
c********************************************************************
c CALCULATE advection, diffusion and Force W-velocity
c********************************************************************
if (convection.eq.'CDS2') then
call advecw_CDS2(dnew,Unew,Vnew,Wnew,Rnew,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px)
elseif(convection.eq.'CDS6') then
call advecw_CDS6(dnew,Unew,Vnew,Wnew,Rnew,rhU,rhV,rhW,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,numdiff,
& periodicx,periodicy,wf,wd)
elseif(convection.eq.'COM4') then
call advecw_COM4(dnew,Unew,Vnew,Wnew,Rnew,Ru,Rp,dr,phivt,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px)
elseif(convection.eq.'CDS4') then
call advecw_CDS4(dnew,Unew,Vnew,Wnew,Rnew,rhU,rhV,rhW,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,numdiff,
& periodicx,periodicy,wf,wd)
elseif(convection.eq.'HYB6') then
call advecw_HYB6(dnew,Unew,Vnew,Wnew,Rnew,rhU,rhV,rhW,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,
& numdiff,periodicx,periodicy,kbed,wf,wd,fc_global,me)
elseif(convection.eq.'HYB4') then
call advecw_HYB4(dnew,Unew,Vnew,Wnew,Rnew,rhU,rhV,rhW,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,
& numdiff,periodicx,periodicy,kbed,wf,wd)
elseif(convection.eq.'C4A6') then
call advecw_C4A6(dnew,Unew,Vnew,Wnew,Rnew,Ru,Rp,dr,phivt,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,numdiff,periodicx,periodicy
&,phiv,wf,wd)
elseif(convection.eq.'uTVD') then
call advecw_TVD(dnew,Unew,Vnew,Wnew,Rnew,Ru,Rp,dr,phip,phiv,phipt,phivt,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px,dt
& ,periodicx,periodicy,istep)
endif
if (diffusion.eq.'CDS2') then
call diffw_CDS2 (dnew,Unew,Vnew,Wnew,ib,ie,jb,je,kb,ke)
elseif(diffusion.eq.'COM4') then
call diffw_com4(dnew,Szr,Spz,Szz,Ru,Rp,dr,phipt,dz,i1,j1,k1,ib,ie,jb,je,kb,ke,rank,px)
endif
if ((slipvel.eq.1.or.slipvel.eq.2).and.nfrac>0) then
call advecw_driftfluxCDS2(dnew,driftfluxforce_calfac*sumWkm,
& rnew,Ru,Rp,dr,phiv,dz,i1,j1,k1,ib,ie,jb,je,kb,ke)
endif
do k=1,kmax
do j=1,jmax
do i=1,imax
dnew(i,j,k)=dnew(i,j,k)-gz*(0.5*(rnew(i,j,k)+rnew(i,j,k+1))-rho_b)
1 +Ppropz(i,j,k)
dWdt(i,j,k)=Wnew(i,j,k)*0.5*(rnew(i,j,k)+rnew(i,j,k+1))+
1 dt*( 23./12.*dnew(i,j,k)-
1 4./3. *wz(i,j,k)+
1 5./12.*wzold(i,j,k))
! 1 dt*(1.5*dnew(i,j,k)-0.5*wz(i,j,k)) -gz*dt*(0.75*(rnew(i,j,k)+rnew(i,j,k+1))-0.25*(rold(i,j,k)+rold(i,j,k+1))-rho_b)
wzold(i,j,k) = wz(i,j,k)
wz(i,j,k) = dnew(i,j,k)
enddo
enddo
enddo
pold2=pold !save pressure previous timestep [pold2 = timestep before pold]
pold=p+pold !what is called p here was dp in reality, now p is
call pshiftb(pold,pplus) !,rank,imax,jmax,kmax,px)
c********************************************************************
c CALCULATE pressure-gradient with old pressure:
c********************************************************************
do k=1,kmax
do j=1,jmax
do i=1,imax-1
dUdt(i,j,k) = dUdt(i,j,k) -dt * ( pold(i+1,j,k) - pold(i,j,k) ) /( Rp(i+1) - Rp(i) )
enddo
enddo
enddo
!!! do nothing on imax and i=1, so p(imax+1)=p(imax) --> dp/dn=0 at outflow and inflow
if (periodicx.eq.1) then
!!! !periodic bc in x:
i = imax
do k=1,kmax
do j=1,jmax
dUdt(i,j,k) = dUdt(i,j,k) -dt * ( pold(1,j,k) - pold(i,j,k) ) /( Rp(i+1) - Rp(i) )
enddo
enddo
endif
do k=1,kmax
do j=1,jmax-1
do i=1,imax
dVdt(i,j,k) = dVdt(i,j,k) -dt * ( pold(i,j+1,k) - pold(i,j,k) ) /( Rp(i) * (phip(j+1)-phip(j)) )
enddo
enddo
enddo
if (periodicy.eq.1) then
do k=1,kmax
do i=1,imax
dVdt(i,jmax,k) = dVdt(i,jmax,k) -dt * ( pplus(i,k) - pold(i,jmax,k) ) /( Rp(i) * (phip(jmax+1)-phip(jmax)) )
enddo
enddo
else
if (rank.ne.px-1) then
do k=1,kmax
do i=1,imax
dVdt(i,jmax,k) = dVdt(i,jmax,k) -dt * ( pplus(i,k) - pold(i,jmax,k) ) /( Rp(i) * (phip(jmax+1)-phip(jmax)) )
enddo
enddo
endif
endif
do k=1,kmax-1 ! rij kmax wordt niet meegenomen, dus is p(k+1)=p(k), dus dp/dn=0 bij wateroppervlak
do j=1,jmax
do i=1,imax
dWdt(i,j,k) = dWdt(i,j,k) -dt * ( pold(i,j,k+1) - pold(i,j,k) ) / dz
enddo
enddo
enddo
return
end
subroutine adamsbv(ib,ie,jb,je,kb,ke)
c
c Performes time integration with second order
c Adams-Bashforth-3 scheme with variable timestep, i.e
c
c
c n+1 n
c dU U - U n
c ---- = ------------ = 23/12*( -ADV + DIFF + Force)
c dt dt n-1
c -4/3 *( -ADV + DIFF + Force)
c n-2
c +5/12*( -ADV + DIFF + Force)
c
c The timestep is limited with CFL=0.7 (see routine chkdt)
c
c
USE nlist
USE sediment
implicit none
include 'mpif.h'
! include 'param.txt'
! include 'common.txt'
integer ib,ie,jb,je,kb,ke,n,t,kpp,kp,km,ip,im,jp,jm,me
real doldc(nfrac,0:i1,0:j1,0:k1),dnewc(nfrac,0:i1,0:j1,0:k1)
real dold(0:i1,0:j1,0:k1),dnew(0:i1,0:j1,0:k1),dnew2(0:i1,0:j1,0:k1)
real wsed(nfrac,0:i1,0:j1,0:k1),sumWkm(0:i1,0:j1,0:k1)
real*8 pplus(imax,kmax)
real fAB3_1,fAB3_2,fAB3_3,fAB3_4,facAB3_1,facAB3_2,facAB3_3,facAB3_4
real dist,timeAB_desired(1:4)
real dnewcbot(nfrac,0:i1,0:j1)
real aaa(0:k1),bbb(0:k1),ccc(0:k1),ekm_min,ekm_plus,rhss(0:k1),rho_min,rho_plus
real aaax(0:i1),bbbx(0:i1),cccx(0:i1),rhssx(0:i1)
real aaay(0:jmax*px+1),bbby(0:jmax*px+1),cccy(0:jmax*px+1),rhssy(0:jmax*px+1)
real ans_T(1:i1,0:jmax*px+1,1:kmax/px+1),ekm_T(1:i1,0:jmax*px+1,1:kmax/px+1)!,ekm2(0:i1,0:j1,0:k1)
real uu_T(1:imax,0:jmax*px+1,1:kmax/px),vv_T(1:imax,0:jmax*px+1,1:kmax/px),ww_T(1:imax,0:jmax*px+1,1:kmax/px)
real rhU_T(1:imax,0:jmax*px+1,1:kmax/px),rhV_T(1:imax,0:jmax*px+1,1:kmax/px),rhW_T(1:imax,0:jmax*px+1,1:kmax/px)
real interface_T(1:i1,0:jmax*px+1)
real utr(0:i1,0:j1,0:k1),vtr(0:i1,0:j1,0:k1),wtr(0:i1,0:j1,0:k1)
real ws(0:i1,0:j1,0:k1),gvector,CNz,dir_order_real
integer ileng,ierr,itag,status(MPI_STATUS_SIZE),perx,pery,dir_order
integer clock
INTEGER, DIMENSION(:), ALLOCATABLE :: seed
!real rhU2(0:i1,0:j1,0:k1),rhV2(0:i1,0:j1,0:k1),rhW2(0:i1,0:j1,0:k1)
!!! real dUdt2(0:i1,0:j1,0:k1),dVdt2(0:i1,0:j1,0:k1),dWdt2(0:i1,0:j1,0:k1)
real ddd1(0:i1,0:j1,0:k1),ddd2(0:i1,0:j1,0:k1),ddd3(0:i1,0:j1,0:k1)
! real eee1(0:i1,0:j1,0:k1),eee2(0:i1,0:j1,0:k1),eee3(0:i1,0:j1,0:k1) !als deze twee regels uitgezet worden dan komen er andere uitkomsten bij CN31?
! real fff1(0:i1,0:j1,0:k1),fff2(0:i1,0:j1,0:k1),fff3(0:i1,0:j1,0:k1)
me = momentum_exchange_obstacles !100 or 101 or 111 means no advection momentum exchange with bed-cells with HYB6
c********************************************************************
c CALCULATE coefficients for AB3 interpolation with variable dt
c********************************************************************
do i=4,1,-1
if (i.eq.1) then
timeAB_real(i)=time_n
else
timeAB_real(i)=timeAB_real(i-1)
endif
timeAB_desired(i)=time_n-(i-1)*dt
enddo
!IF (CNdiffz.eq.11.or.CNdiffz.eq.12) THEN !switch to AB2 with variable time step:
IF (.false.) THEN !always use AB3 with variable time step
! ! !EE1:
! ! facAB3_1=1. !1.5 !23./12.
! ! facAB3_2=0. !-0.5 !-4./3.
! ! facAB3_3=0. !5./12.
! ! facAB3_4=0.
if (istep.lt.5) then !start with standard AB2
facAB3_1=1.5 !23./12.
facAB3_2=-0.5 !-4./3.
facAB3_3=0. !5./12.
facAB3_4=0.
else
dist=timeAB_desired(2)-timeAB_real(2)
if (dist.lt.0) then
fAB3_2=(timeAB_real(2)-timeAB_real(3)+dist)/(timeAB_real(2)-timeAB_real(3))
fAB3_3=1.-(timeAB_real(2)-timeAB_real(3)+dist)/(timeAB_real(2)-timeAB_real(3))
fAB3_1=0.
else
fAB3_2=(timeAB_real(1)-timeAB_real(2)-dist)/(timeAB_real(1)-timeAB_real(2))
fAB3_1=1.-(timeAB_real(1)-timeAB_real(2)-dist)/(timeAB_real(1)-timeAB_real(2))
fAB3_3=0.
endif
facAB3_1=1.5-fAB3_1*0.5 !23./12.-fAB3_1*4./3.
facAB3_2=-fAB3_2*0.5 !-fAB3_2*4./3.
facAB3_3=-fAB3_3*0.5 !-fAB3_3*4./3.
facAB3_4=0.
endif
ELSE !AB3
if (istep.lt.5) then !start with standard AB3
! facAB3_1=1.
! facAB3_2=0.
! facAB3_3=0.
! facAB3_4=0.
facAB3_1=23./12.
facAB3_2=-4./3.
facAB3_3=5./12.
facAB3_4=0.
else
dist=timeAB_desired(2)-timeAB_real(2)
if (dist.lt.0) then
fAB3_2=(timeAB_real(2)-timeAB_real(3)+dist)/(timeAB_real(2)-timeAB_real(3))
fAB3_3=1.-(timeAB_real(2)-timeAB_real(3)+dist)/(timeAB_real(2)-timeAB_real(3))
fAB3_1=0.
else
fAB3_2=(timeAB_real(1)-timeAB_real(2)-dist)/(timeAB_real(1)-timeAB_real(2))
fAB3_1=1.-(timeAB_real(1)-timeAB_real(2)-dist)/(timeAB_real(1)-timeAB_real(2))
fAB3_3=0.
endif
facAB3_1=23./12.-fAB3_1*4./3.
facAB3_2=-fAB3_2*4./3.
facAB3_3=-fAB3_3*4./3.
dist=timeAB_desired(3)-timeAB_real(3)
if (dist.lt.0) then
fAB3_3=(timeAB_real(3)-timeAB_real(4)+dist)/(timeAB_real(3)-timeAB_real(4))
fAB3_4=1.-(timeAB_real(3)-timeAB_real(4)+dist)/(timeAB_real(3)-timeAB_real(4))
fAB3_2=0.
else
fAB3_3=(timeAB_real(2)-timeAB_real(3)-dist)/(timeAB_real(2)-timeAB_real(3))
fAB3_2=1.-(timeAB_real(2)-timeAB_real(3)-dist)/(timeAB_real(2)-timeAB_real(3))
fAB3_4=0.
endif
facAB3_2=facAB3_2+fAB3_2*5./12.
facAB3_3=facAB3_3+fAB3_3*5./12.
facAB3_4=fAB3_4*5./12.
endif
ENDIF
if (applyVOF.eq.1) then ! VOF fraction apply real density for gravity and sediment settling related things
CALL state(cnew,rnew)
endif
c********************************************************************
c CALCULATE tau walls based on divergence free velocity field previous time step
c********************************************************************
dUdt = Unew*rhU ! for some continuity_solver dUdt-end is without rho
dVdt = Vnew*rhV
dWdt = Wnew*rhW
call bound_rhoU(dUdt,dVdt,dWdt,rnew,MIN(3,slip_bot),monopile,time_n,
& Ub1old,Vb1old,Wb1old,Ub2old,Vb2old,Wb2old,Ub3old,Vb3old,Wb3old) ! apply wall shear stress based on UVWnew (valid at time_n and UB1old...)
dcdt = 0.
sumWkm = 0.
dnewc=0.
dnew2=0. !initialize as zero
c********************************************************************
c CALCULATE slipvelocity
c********************************************************************
if (nfrac>0) then
if (slipvel.eq.1.or.slipvel.eq.2) then
call slipvelocity(cnew,Wnew,wsed,rnew,1,kmax,sumWkm,dt,dz)
!! Set boundary conditions jet in:
do t=1,tmax_inPpunt
i=i_inPpunt(t)
j=j_inPpunt(t)
do n=1,nfrac
do k=kmax,kmax !-kjet,kmax
Wsed(n,i,j,k)=Wnew(i,j,k)
enddo
enddo
if (LOA>0.and.outflow_overflow_down.eq.1) then
do n=1,nfrac
do k=kmax-kjet+1,kmax
Wsed(n,i,j,k)=Wnew(i,j,k)
enddo
enddo
endif
enddo
else
do n=1,nfrac
wsed(n,:,:,:)=wnew
enddo
endif
c********************************************************************
c CALCULATE advection, diffusion Concentration
c********************************************************************
do n=1,nfrac
if (transporteq_fracs.eq.'massfrac') then
utr=rhoU
vtr=rhoV
wtr(:,:,0:kmax)=wsed(n,:,:,0:kmax)*0.5*(rnew(:,:,0:kmax)+rnew(:,:,1:k1))
cnew(n,:,:,:)=cnew(n,:,:,:)*frac(n)%rho/rnew(:,:,:) ! go from volume fraction to mass fraction
Diffcof=Diffcof*rnew
elseif (transporteq_fracs.eq.'massfra2') then
cnew(n,:,:,:)=cnew(n,:,:,:)*frac(n)%rho ! go from volume fraction to (mass-fraction*rho_mix)
utr=Unew
vtr=Vnew
wtr=Wsed(n,:,:,:)
else
utr=Unew
vtr=Vnew
wtr=Wsed(n,:,:,:)
endif
if (settling_along_gvector.eq.1) then
gvector=sqrt(gx**2+gy**2+gz**2)
ws = Wsed(n,:,:,:)-Wnew
utr=utr+gx/gvector*ws
vtr=vtr+gy/gvector*ws
wtr=wtr-ws+gz/gvector*ws
endif
call make_UtransC_zeroIBM(utr,vtr,wtr)
if (frac(n)%type.eq.-1) then ! VOF fraction
call advecc_VOF(dnewc(n,:,:,:),cnew(n,:,:,:),utr,vtr,wtr,rnew,Ru,Rp,dr,phiv,phipt,dz,
+ i1,j1,k1,ib,ie,jb,je,kb,ke,dt,rank,px,periodicx,periodicy)
! call advecc_NVD_VOF(dnewc(n,:,:,:),cnew(n,:,:,:),utr,vtr,wtr,rnew,Ru,Rp,dr,phiv,phipt,dz,
! + i1,j1,k1,ib,ie,jb,je,kb,ke,dt,rank,px,periodicx,periodicy)
dcdt(n,:,:,:) =cnew(n,:,:,:) + dt*(dnewc(n,:,:,:)) !update in time with EE1 for TVD, no diffusion
do i=1,imax
do j=1,jmax
do k=1,kmax
dcdt(n,i,j,k)=MAX(0.,dcdt(n,i,j,k))
dcdt(n,i,j,k)=MIN(1.,dcdt(n,i,j,k))
enddo
enddo
enddo
else
if (advec_conc.eq.'NVD') then
call advecc_NVD(dnewc(n,:,:,:),cnew(n,:,:,:),utr,vtr,wtr,rnew,Ru,Rp,dr,phiv,phipt,dz,
+ i1,j1,k1,ib,ie,jb,je,kb,ke,dt,rank,px,periodicx,periodicy)
elseif (advec_conc.eq.'VLE') then
call advecc_VLE(dnewc(n,:,:,:),cnew(n,:,:,:),utr,vtr,wtr,rnew,Ru,Rp,dr,phiv,phipt,dz,
+ i1,j1,k1,ib,ie,jb,je,kb,ke,dt,rank,px,periodicx,periodicy,transporteq_fracs,kbed)
elseif (advec_conc.eq.'SBE') then
call advecc_SBE(dnewc(n,:,:,:),cnew(n,:,:,:),utr,vtr,wtr,rnew,Ru,Rp,dr,phiv,phipt,dz,
+ i1,j1,k1,ib,ie,jb,je,kb,ke,dt,rank,px,periodicx,periodicy,transporteq_fracs,kbed)
elseif (advec_conc.eq.'VL2') then
call advecc_VL2_nocfl(dnewc(n,:,:,:),cnew(n,:,:,:),utr,vtr,wtr,rnew,Ru,Rp,dr,phiv,phipt,dz,
+ i1,j1,k1,ib,ie,jb,je,kb,ke,dt,rank,px,periodicx,periodicy)
elseif (advec_conc.eq.'SB2') then
call advecc_SB2_nocfl(dnewc(n,:,:,:),cnew(n,:,:,:),utr,vtr,wtr,rnew,Ru,Rp,dr,phiv,phipt,dz,
+ i1,j1,k1,ib,ie,jb,je,kb,ke,dt,rank,px,periodicx,periodicy)
elseif (advec_conc.eq.'ARO') then
call advecc_ARO(dnewc(n,:,:,:),cnew(n,:,:,:),utr,vtr,wtr,rnew,Ru,Rp,dr,phiv,phipt,dz,
+ i1,j1,k1,ib,ie,jb,je,kb,ke,dt,rank,px,periodicx,periodicy)
endif
if (Sc<1.e18) then
call diffc_CDS2 (dnewc(n,:,:,:),Cnew(n,:,:,:),Diffcof,
+ ib,ie,jb,je,kb,ke)
endif
dcdt(n,:,:,:) =cnew(n,:,:,:) + dt*(dnewc(n,:,:,:)) !update in time with EE1 for TVD
endif
if (transporteq_fracs.eq.'massfrac') then
dcdt(n,:,:,:)=dcdt(n,:,:,:)+cnew(n,:,:,:)*rnew-cnew(n,:,:,:)
dcdt(n,:,:,:)=dcdt(n,:,:,:)/frac(n)%rho ! from dcdt (mass-fraction*rho_mix) back to volume-fraction
cnew(n,:,:,:)=cnew(n,:,:,:)*rnew/frac(n)%rho ! from mass fraction back to volume fraction
Diffcof=Diffcof/rnew
elseif (transporteq_fracs.eq.'massfra2') then
dcdt(n,:,:,:)=dcdt(n,:,:,:)/frac(n)%rho ! from dcdt (mass-fraction*rho_mix) back to volume-fraction
cnew(n,:,:,:)=cnew(n,:,:,:)/frac(n)%rho ! from (mass-fraction*rho_mix) back to volume fraction
endif
IF (CNdiffz.eq.1.or.CNdiffz.eq.2.or.CNdiffz.eq.11.or.CNdiffz.eq.12) THEN !CN semi implicit treatment diff-z:
IF (CNdiffz.eq.1.or.CNdiffz.eq.11) THEN
CNz=0.5
ELSE
CNz=1.
ENDIF
call bound_c(dcdt(n,:,:,:),frac(n)%c,n,0.) ! bc start CN-diffz ABv
do k=k1,k1 !-kjet,k1
do t=1,tmax_inPpunt
i=i_inPpunt(t)
j=j_inPpunt(t)
dcdt(n,i,j,k) = dcdt(n,i,j,kmax) ! No diffusion over vertical inflow boundary, this line is needed for exact influx
enddo
enddo
!!! dUdt2=dcdt(n,:,:,:) !for non-updated RHS
IF (CNdiffz.eq.11.or.CNdiffz.eq.12) THEN
perx=0
pery=0
IF (periodicx.eq.1) perx=1
IF (periodicy.eq.1) pery=1
call t2np(Diffcof(1:imax,1:jmax,1:kmax),ekm_T(1:imax,1:jmax*px,1:kmax/px))
call t2np(dcdt(n,1:imax,1:jmax,1:kmax),ans_T(1:imax,1:jmax*px,1:kmax/px))
!! also pass over boundaries at j=0 :
IF (rank.eq.0) THEN
do i=1,px-1
call mpi_send(Diffcof(1:imax,0,i*kmax/px+1:(i+1)*kmax/px),imax*kmax/px,MPI_REAL8,i,i,MPI_COMM_WORLD,status,ierr)
call mpi_send(dcdt(n,1:imax,0,i*kmax/px+1:(i+1)*kmax/px),imax*kmax/px,MPI_REAL8,i,i+100000,MPI_COMM_WORLD,status,ierr)
enddo
ekm_T(1:imax,0,1:kmax/px)=Diffcof(1:imax,0,1:kmax/px)
ans_T(1:imax,0,1:kmax/px)=dcdt(n,1:imax,0,1:kmax/px)
ELSE
call mpi_recv(ekm_T(1:imax,0,1:kmax/px),imax*kmax/px,MPI_REAL8,0,rank,MPI_COMM_WORLD,status,ierr)
call mpi_recv(ans_T(1:imax,0,1:kmax/px),imax*kmax/px,MPI_REAL8,0,rank+100000,MPI_COMM_WORLD,status,ierr)
ENDIF
!! also pass over boundaries at j=jmax+1 :
IF (rank.eq.px-1) THEN
do i=0,px-2
call mpi_send(Diffcof(1:imax,jmax+1,i*kmax/px+1:(i+1)*kmax/px),imax*kmax/px,MPI_REAL8,i,i+200000,MPI_COMM_WORLD
& ,status,ierr)
call mpi_send(dcdt(n,1:imax,jmax+1,i*kmax/px+1:(i+1)*kmax/px),imax*kmax/px,MPI_REAL8,i,i+300000,MPI_COMM_WORLD
& ,status,ierr)
enddo
ekm_T(1:imax,jmax*px+1,1:kmax/px)=Diffcof(1:imax,jmax+1,rank*kmax/px+1:(rank+1)*kmax/px)
ans_T(1:imax,jmax*px+1,1:kmax/px)=dcdt(n,1:imax,jmax+1,rank*kmax/px+1:(rank+1)*kmax/px)
ELSE
call mpi_recv(ekm_T(1:imax,jmax*px+1,1:kmax/px),imax*kmax/px,MPI_REAL8,px-1,rank+200000,MPI_COMM_WORLD,status,ierr)
call mpi_recv(ans_T(1:imax,jmax*px+1,1:kmax/px),imax*kmax/px,MPI_REAL8,px-1,rank+300000,MPI_COMM_WORLD,status,ierr)
ENDIF
IF (periodicy.eq.0) THEN !bc defined at both lateral boundaries
do k=1,kmax/px
do i=1,imax
do j=1,px*jmax !0,px*jmax+1
jm=j-1 !MAX(0,j-1)
jp=j+1 !MIN(px*jmax+1,j+1)
ekm_min=0.5* (ekm_T(i,jm,k)+ekm_T(i,j,k))*MIN(fc_global(i,jm,k+rank*kmax/px),fc_global(i,j,k+rank*kmax/px))
ekm_plus=0.5*(ekm_T(i,jp,k)+ekm_T(i,j,k))*MIN(fc_global(i,jp,k+rank*kmax/px),fc_global(i,j,k+rank*kmax/px))
aaay(j)=-CNz*ekm_min*dt/((Rp(i)*(phipt(j)-phipt(jm)))*Rp(i)*dphi2t(j))
cccy(j)=-CNz*ekm_plus*dt/((Rp(i)*(phipt(jp)-phipt(j)))*Rp(i)*dphi2t(j))
bbby(j)=1.-aaay(j)-cccy(j)
enddo
aaay(0)=0.
aaay(px*jmax+1)=0.
cccy(0)=0.
cccy(px*jmax+1)=0.
bbby(0)=1.
bbby(px*jmax+1)=1.
rhssy=ans_T(i,0:px*jmax+1,k)
CALL solve_tridiag_switchperiodic2(ans_T(i,0:px*jmax+1,k),aaay,bbby,cccy,rhssy,px*jmax+2,pery)
enddo
enddo
ELSEIF (periodicy.eq.1) THEN !periodic lateral boundaries
do k=1,kmax/px
do i=1,imax
do j=1,px*jmax !0,px*jmax+1
jm=j-1 !MAX(0,j-1)
jp=j+1 !MIN(px*jmax+1,j+1)
ekm_min=0.5* (ekm_T(i,jm,k)+ekm_T(i,j,k))*MIN(fc_global(i,jm,k+rank*kmax/px),fc_global(i,j,k+rank*kmax/px))
ekm_plus=0.5*(ekm_T(i,jp,k)+ekm_T(i,j,k))*MIN(fc_global(i,jp,k+rank*kmax/px),fc_global(i,j,k+rank*kmax/px))
aaay(j)=-CNz*ekm_min*dt/((Rp(i)*(phipt(j)-phipt(jm)))*Rp(i)*dphi2t(j))
cccy(j)=-CNz*ekm_plus*dt/((Rp(i)*(phipt(jp)-phipt(j)))*Rp(i)*dphi2t(j))
bbby(j)=1.-aaay(j)-cccy(j)
enddo
rhssy(1:px*jmax)=ans_T(i,1:px*jmax,k)
CALL solve_tridiag_switchperiodic2(ans_T(i,1:px*jmax,k),aaay(1:px*jmax),bbby(1:px*jmax),cccy(1:px*jmax),
& rhssy(1:px*jmax),px*jmax,pery)
enddo
enddo
ELSEIF (periodicy.eq.2) THEN !free slip lateral boundaries dUWdn=0 at both ends and V=0
do k=1,kmax/px
do i=1,imax
do j=1,px*jmax !0,px*jmax+1
jm=j-1 !MAX(0,j-1)