-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathELM.f90
2489 lines (2258 loc) · 165 KB
/
ELM.f90
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
!>\brief Routines related to the Eulerian-Lagragean Method (FU, FV and FW) of TRIM model - Casulli
!>\details
! References:
! [1]
! [2] Cunha, A.H.F.; Fragoso, C.R.; Tavares, M.H.; Cavalcanti, J.R.; Bonnet, M.-P.; Motta-Marques, D. Combined Use of High-Resolution Numerical
! Schemes to Reduce Numerical Diffusion in Coupled Nonhydrostatic Hydrodynamic and Solute Transport Model. Water 2019, 11, 2288.
!>\author Carlos Ruberto Fragoso and Rafael Cavalcanti
Module ELM
private
public:: FuFv
public:: ComputeFW
public:: signa
!-----------------------------------------------------------------
Contains
!> Main subroutine of Eulerian-Lagragean Method
!>\return Fu, Fv and Fw the resultant vectors for each edge
!>\author Carlos Ruberto Fragoso
!>\attention List of Modifications: \n
!> - 15.09.2015: Routine Implementation (Carlos Ruberto Fragoso)
Subroutine FuFv(HydroParam,MeshParam,dt)
Use MeshVars !, Only: Left,Right,nElem,nEdge,EdgeBary,NormalVector,TangentVector,Edge,EdgeLength,Quadri,xNode,yNode,Area,xb,yb,nNode,Kmax,nEgdesatNode,EgdesatNode
Use Hydrodynamic !, Only: Smallm,CapitalM,Z,u,utang,Fu,Fv,H,Pcri,w,Nfut,uxy,uxyback,ElSmallm,ElCapitalM,lat,OMEGA,Pi,Wu,HorViscosity,Fub,FuxyNode
!Use SimulationModel
Implicit None
Integer:: iElem,iEdge,lEdge,iLayer,jlev,l,r,nnel,ndels,Sig,Face,n1,n2,n3,n4,iNode,icount,j,fac,kin,TrajectoryFlag,FuFw_flag,jump, ntals, iEdge0
Real::x0,y0,z0,xt,yt,zt,uuint,vvint,wwint,hhint,wdown,wup,vmag,dtb,aa(4),weit,sinal,dtbCFL
Real:: NearZero = 1e-10
Real::dt, tal
type(MeshGridParam) :: MeshParam
type(HydrodynamicParam) :: HydroParam
! 1.1 Find out uxy velocities components, that are equivalents to normal and tangencial velocities components in particle position after
! the backtrack process:
Do iEdge = 1,MeshParam%nEdge
!Get right/left iElements that share iEdge
l = MeshParam%Left(iEdge)
r = MeshParam%Right(iEdge)
! If is a dry cell, the backtrack velocities is equal to zero.
! Note that in subsurface coupled case %H(iEdge) can takes into account the freatic water level in Cell.
! DZsj represents the Edge freatic component thickness, thus H - DZsj is equivalent only surface water level
! In Cell where no subsufarce layer exist or in only surface flow simulation DZsj is set equal 0.0d0 (ReadHydroIniCond Module)
If (HydroParam%H(iEdge) + HydroParam%sj(iEdge)-HydroParam%hj(iEdge)<=HydroParam%Pcri+NearZero) Then !CAYO
Do iLayer = HydroParam%Smallm(iEdge), HydroParam%CapitalM(iEdge)
HydroParam%uxyback(iLayer,1:2,iEdge) = (/ 0., 0. /)
EndDo
Cycle
EndIf
! In case that Cell is Wet:
Do iLayer = HydroParam%Smallm(iEdge), HydroParam%CapitalM(iEdge)
! The particle is set in iEdge barycenter (x0,y0,z0) on Left side Element (nnel):
nnel = l
jlev = iLayer
x0 = MeshParam%EdgeBary(1,iEdge)
y0 = MeshParam%EdgeBary(2,iEdge)
z0 = 0.5d0*(HydroParam%Z(iLayer+1,iEdge) + HydroParam%Z(iLayer,iEdge) + sum(HydroParam%DZsj(:,iEdge)))
! Velocities in the initial position:
uuint = HydroParam%uxy(iLayer,1,iEdge)
vvint = HydroParam%uxy(iLayer,2,iEdge)
wwint = HydroParam%wfc(iLayer,iEdge)
hhint = HydroParam%H(iEdge)-HydroParam%hj(iEdge)
! Velocity magnitude:
vmag = dsqrt(uuint**2+vvint**2+wwint**2)
! In cases that the velocity magnitude in layer is very low or the iEdge no has neighbour Element (r==0), the backtracking
! velocity in layer is set equal the iEdge:
If(vmag.le.1.e-6) Then ! No activity or Water Level Condition
HydroParam%uxyback(iLayer,1:2,iEdge) = (/ uuint, vvint /)
!Elseif (HydroParam%eta(nnel) < HydroParam%hb(nnel) .and. HydroParam%eta(r) > HydroParam%hb(r)) Then
! HydroParam%uxyback(iLayer,1:2,iEdge) = (/ uuint, vvint /)
Else
! Else, the backtracking process is initialized:
If (r==0) Then !.and.HydroParam%u(iLayer,iEdge)>0
HydroParam%uxyback(iLayer,1:2,iEdge) = (/ uuint, vvint /)
Else
If (iLayer < Hydroparam%ElSmallm(r)) Then
HydroParam%uxyback(iLayer,1,iEdge) = uuint
HydroParam%uxyback(iLayer,2,iEdge) = vvint
Else
! If TrajectoryFlag == 0 System Defined Intervals to Integrate Particle Trajectory,
! else User Defined Intervals to Integrate Particle Trajectory.
TrajectoryFlag = 0
! Find the Local Time Step (dtb) to Integrate Trajectory:
If ( TrajectoryFlag == 0 ) Then ! System Calculation - [7,8]
dtb = MinVal(MeshParam%InCircle)/vmag ! ( InCircle(lElem)/Sqrt( Veloc(1)**2. + Veloc(2)**2. ), InCircle(r)/Sqrt( Veloc(1)**2. + Veloc(2)**2. ) )
!dtbCFL = Max((1/HydroParam%CFL(l)),(1/HydroParam%CFL(r)))
!ndels = Max(Floor(dt/dtb),HydroParam%NFUT,Floor(dt/dtbCFL))
ndels = Max(Floor(dt/dtb),HydroParam%NFUT)
dtb = dt/ndels
ElseIf ( TrajectoryFlag == 1 ) Then
dtbCFL = Max((1/HydroParam%CFL(l)),(1/HydroParam%CFL(r)))
ndels = Max(HydroParam%NFUT,Floor(dt/dtbCFL))
!ndels = HydroParam%NFUT
dtb = dt/ndels
EndIf
FuFw_flag = 0
iEdge0 = iEdge
Call btrack(ndels,dtb,dt,uuint,vvint,wwint,hhint,&
&x0,y0,z0,xt,yt,zt,nnel,jlev,iEdge0,HydroParam,MeshParam, FuFw_flag, iLayer, iEdge)
HydroParam%uxyback(iLayer,1:2,iEdge) = (/ uuint, vvint /)
EndIf
EndIf
Endif
EndDo
EndDo
If (HydroParam%iNonHydro==0.and.MeshParam%Kmax>1) Then
Call ComputeFW(MeshParam,HydroParam,dt)
EndIf
!Explicit terms
Do iEdge=1,MeshParam%nEdge
l = MeshParam%Left(iEdge)
r = MeshParam%Right(iEdge)
Do lEdge=1,4
If (MeshParam%Edge(lEdge,l)==iEdge) Then
Exit
EndIf
EndDo
!If (lEdge == 1) Then
! Do iLayer = HydroParam%Smallm(iEdge), HydroParam%CapitalM(iEdge)
! HydroParam%Fv(iLayer,iEdge) = Sig(l,MeshParam%Right(iEdge),MeshParam%Left(iEdge))*HydroParam%uxyback(iLayer,1,iEdge)*MeshParam%TangentVector(1,lEdge,l) + Sig(l,MeshParam%Right(iEdge),MeshParam%Left(iEdge))*HydroParam%uxyback(iLayer,2,iEdge)*MeshParam%TangentVector(2,lEdge,l)
! EndDo
!ElseIf(lEdge == 2) Then
! Do iLayer = HydroParam%Smallm(iEdge), HydroParam%CapitalM(iEdge)
! HydroParam%Fv(iLayer,iEdge) = Sig(l,MeshParam%Right(iEdge),MeshParam%Left(iEdge))*HydroParam%uxyback(iLayer,1,iEdge)*MeshParam%TangentVector(1,lEdge,l) - Sig(l,MeshParam%Right(iEdge),MeshParam%Left(iEdge))*HydroParam%uxyback(iLayer,2,iEdge)*MeshParam%TangentVector(2,lEdge,l)
! EndDo
!Else ! For iEdge == 3 and 4
! Do iLayer = HydroParam%Smallm(iEdge), HydroParam%CapitalM(iEdge)
! HydroParam%Fv(iLayer,iEdge) = Sig(l,MeshParam%Right(iEdge),MeshParam%Left(iEdge))*(HydroParam%uxyback(iLayer,1,iEdge)*MeshParam%TangentVector(1,lEdge,l) + HydroParam%uxyback(iLayer,2,iEdge)*MeshParam%TangentVector(2,lEdge,l))
! EndDo
!EndIf
If (lEdge == 1) Then
Do iLayer = HydroParam%Smallm(iEdge), HydroParam%CapitalM(iEdge)
HydroParam%Fv(iLayer,iEdge) = HydroParam%uxyback(iLayer,1,iEdge)*MeshParam%TangentVector(1,lEdge,l) + HydroParam%uxyback(iLayer,2,iEdge)*MeshParam%TangentVector(2,lEdge,l)
EndDo
ElseIf(lEdge == 2) Then
Do iLayer = HydroParam%Smallm(iEdge), HydroParam%CapitalM(iEdge)
HydroParam%Fv(iLayer,iEdge) = HydroParam%uxyback(iLayer,1,iEdge)*MeshParam%TangentVector(1,lEdge,l) - HydroParam%uxyback(iLayer,2,iEdge)*MeshParam%TangentVector(2,lEdge,l)
EndDo
Else ! For iEdge == 3 and 4
Do iLayer = HydroParam%Smallm(iEdge), HydroParam%CapitalM(iEdge)
HydroParam%Fv(iLayer,iEdge) = HydroParam%uxyback(iLayer,1,iEdge)*MeshParam%TangentVector(1,lEdge,l) + HydroParam%uxyback(iLayer,2,iEdge)*MeshParam%TangentVector(2,lEdge,l)
EndDo
EndIf
!
! If the face has a boundary condition:
If (HydroParam%IndexInflowEdge(iEdge)>0.or.HydroParam%IndexWaterLevelEdge(iEdge)>0) then
HydroParam%Fu(:,iEdge) = HydroParam%u(:,iEdge)
HydroParam%Fv(:,iEdge) = HydroParam%utang(:,iEdge)
cycle
EndIf
!MeshParam%Neighbor(:,l)
Do iLayer = HydroParam%Smallm(iEdge), HydroParam%CapitalM(iEdge)
If (r==0) Then !Velocity in lateral boundaries is equal to cell center velocity
HydroParam%Fu(iLayer,iEdge) = HydroParam%u(iLayer,iEdge) !Sig(l,MeshParam%Right(iEdge),MeshParam%Left(iEdge))*(HydroParam%Fub(iLayer,1,l))*MeshParam%NormalVector(1,lEdge,l) + Sig(l,MeshParam%Right(iEdge),MeshParam%Left(iEdge))*(HydroParam%Fub(iLayer,2,l))*MeshParam%NormalVector(2,lEdge,l)
Else
!sinal = Sig(l,MeshParam%Right(iEdge),MeshParam%Left(iEdge))
! Do in the Edges - > sig func always return 1, because right and left elements in Edge doesn't change (l == left always) - > check if is necessary
if (iLayer < HydroParam%ElSmallm(r)) Then
HydroParam%Fu(iLayer,iEdge) = HydroParam%u(iLayer,iEdge)
else
HydroParam%Fu(iLayer,iEdge) = Sig(l,MeshParam%Right(iEdge),MeshParam%Left(iEdge))*HydroParam%uxyback(iLayer,1,iEdge)*MeshParam%NormalVector(1,lEdge,l) + Sig(l,MeshParam%Right(iEdge),MeshParam%Left(iEdge))*HydroParam%uxyback(iLayer,2,iEdge)*MeshParam%NormalVector(2,lEdge,l)
!HydroParam%Fu(iLayer,iEdge) = Sig(l,MeshParam%Right(iEdge),MeshParam%Left(iEdge))*HydroParam%uxyback(iLayer,1,iEdge)*MeshParam%NormalVector(1,lEdge,l) + (Sig(l,MeshParam%Right(iEdge),MeshParam%Left(iEdge))*HydroParam%uxyback(iLayer,2,iEdge)*MeshParam%NormalVector(2,lEdge,l)*dt)*((HydroParam%uxyback(iLayer,2,iEdge-3)-2*HydroParam%uxyback(iLayer,2,iEdge)+HydroParam%uxyback(iLayer,2,iEdge+3))/MeshParam%EdgeLength(iEdge))
endif
EndIf
HydroParam%Fv(iLayer,iEdge) = HydroParam%uxyback(iLayer,1,iEdge)*MeshParam%TangentVector(1,lEdge,l) + HydroParam%uxyback(iLayer,2,iEdge)*MeshParam%TangentVector(2,lEdge,l)
EndDo
EndDo
!! When we try to calculating Fu and Fw with the particle starting from the center of the cell
!Do iElem = 1,MeshParam%nElem
! n1=MeshParam%Quadri(1,iElem) + 1
! n2=MeshParam%Quadri(2,iElem) + 1
! n3=MeshParam%Quadri(3,iElem) + 1
! n4=MeshParam%Quadri(4,iElem) + 1
! aa(1)=signa(MeshParam%xNode(n1),MeshParam%xNode(n2),MeshParam%xb(iElem),MeshParam%yNode(n1),MeshParam%yNode(n2),MeshParam%yb(iElem))
! aa(2)=signa(MeshParam%xNode(n2),MeshParam%xNode(n3),MeshParam%xb(iElem),MeshParam%yNode(n2),MeshParam%yNode(n3),MeshParam%yb(iElem))
! aa(3)=signa(MeshParam%xNode(n3),MeshParam%xNode(n4),MeshParam%xb(iElem),MeshParam%yNode(n3),MeshParam%yNode(n4),MeshParam%yb(iElem))
! aa(4)=signa(MeshParam%xNode(n4),MeshParam%xNode(n1),MeshParam%xb(iElem),MeshParam%yNode(n4),MeshParam%yNode(n1),MeshParam%yb(iElem))
! Do iLayer = HydroParam%ElSmallm(iElem), HydroParam%ElCapitalM(iElem)
! HydroParam%Fub(iLayer,1:2,iElem) = 0.
! Do iEdge = 1,4
! HydroParam%Fu(iLayer,MeshParam%Edge(iEdge,iElem)) = Sig(iElem,MeshParam%Right(MeshParam%Edge(iEdge,iElem)),MeshParam%Left(MeshParam%Edge(iEdge,iElem)))*HydroParam%uxyback(iLayer,1,MeshParam%Edge(iEdge,iElem))*MeshParam%NormalVector(1,iEdge,iElem) + Sig(iElem,MeshParam%Right(MeshParam%Edge(iEdge,iElem)),MeshParam%Left(MeshParam%Edge(iEdge,iElem)))*HydroParam%uxyback(iLayer,2,MeshParam%Edge(iEdge,iElem))*MeshParam%NormalVector(2,iEdge,iElem)
! HydroParam%Fv(iLayer,MeshParam%Edge(iEdge,iElem)) = Sig(iElem,MeshParam%Right(MeshParam%Edge(iEdge,iElem)),MeshParam%Left(MeshParam%Edge(iEdge,iElem)))*HydroParam%uxyback(iLayer,1,MeshParam%Edge(iEdge,iElem))*MeshParam%TangentVector(1,iEdge,iElem) + Sig(iElem,MeshParam%Right(MeshParam%Edge(iEdge,iElem)),MeshParam%Left(MeshParam%Edge(iEdge,iElem)))*HydroParam%uxyback(iLayer,2,MeshParam%Edge(iEdge,iElem))*MeshParam%TangentVector(2,iEdge,iElem)
! HydroParam%Fub(iLayer,1,iElem) = HydroParam%ubBack(iLayer,1,iElem) !HydroParam%Fub(iLayer,1,iElem) + HydroParam%uxyback(iLayer,1,MeshParam%Edge(iEdge,iElem))*aa(iEdge)/MeshParam%Area(iElem)
! HydroParam%Fub(iLayer,2,iElem) = HydroParam%ubBack(iLayer,2,iElem) !HydroParam%Fub(iLayer,2,iElem) + HydroParam%uxyback(iLayer,2,MeshParam%Edge(iEdge,iElem))*aa(iEdge)/MeshParam%Area(iElem)
! EndDo
! EndDo
!
!EndDo
HydroParam%FuxyNode = 0.
Do iNode=1,MeshParam%nNode
Do iLayer = 1,MeshParam%Kmax
weit=0
icount=0
Do j=1,MeshParam%nEgdesatNode(iNode)
Face = MeshParam%EgdesatNode(j,iNode)
l = MeshParam%Left(Face)
r = MeshParam%Right(Face)
If (r==0) Then
fac=2
Else
fac=1
EndIf
If(iLayer>=HydroParam%Smallm(Face)) Then
kin=min(iLayer,HydroParam%CapitalM(Face))
HydroParam%FuxyNode(kin,1,iNode)=HydroParam%FuxyNode(kin,1,iNode)+HydroParam%uxyback(kin,1,Face)/MeshParam%EdgeLength(Face)*fac
HydroParam%FuxyNode(kin,2,iNode)=HydroParam%FuxyNode(kin,2,iNode)+HydroParam%uxyback(kin,2,Face)/MeshParam%EdgeLength(Face)*fac
icount=icount+1
EndIf
weit=weit+fac/MeshParam%EdgeLength(Face)
EndDo
If(icount.ne.0) Then
HydroParam%FuxyNode(iLayer,1,iNode) = HydroParam%FuxyNode(iLayer,1,iNode)/weit
HydroParam%FuxyNode(iLayer,2,iNode) = HydroParam%FuxyNode(iLayer,2,iNode)/weit
EndIf
EndDo
EndDo
Return
End Subroutine FuFv
Subroutine btrack(ndelt,dtb,dt,uuint,vvint,wwint,hhint,x0,y0,z0,&
&xt, yt, zt, nnel, jlev, id0, HydroParam, MeshParam,FuFw_flag, iLayer, iEdge)
!> Routine for backtracking.
!> Input: point P0(x0,y0,z0) and nnel that encloses it, and initial vel.,
!> initial level (for quicksearch), and a flag indicating 1st or 2nd tracking.\n
!> Output: destination point Pt(xt,yt,zt), element nnel and level jlev,
!> and vel. there (uuint,vvint,wwint).
!>\param imode mode of backtracking (imode = 1: hydrodynamic)
!>\param ndelt number of sub step time (input)
!>\param dtb sub step time in seg (input)
!>\param uuint x velocity component in the center of the edge (input)
!>\param vvint y velocity component in the center of the edge (input)
!>\param wwint z velocity component in the center of the edge (input)
!>\param nnel number of the element that enclosed the point P0 (input)
!>\param jlev number of the layer that enclosed the point P0 (input)
!>\param x0 x initial position of the point P0 (input)
!>\param y0 y initial position of the point P0(input)
!>\param z0 z initial position of the point P0 (input)
!>\return xt: x final position of the point Pt (output)
!>\return yt: y final position of the point Pt (output)
!>\return zt: z final position of the point Pt (output)
!>\return uuint: x velocity component in the Edge (output)
!>\return vvint: y velocity component in the Edge (output)
!>\return wwint: z velocity component in the Edge (output)
!>\return nnel: number of the element that enclosed the point Pt (output)
!>\return jlev: number of the layer that enclosed the point Pt (output)
!>\author Carlos Ruberto Fragoso
!>\attention List of Modifications: \n
!> - 15.09.2015: Routine Implementation (Carlos Ruberto Fragoso)
Use MeshVars !, Only: Quadri,EdgeDef,xNode,yNode,Area,Edge,EdgeNodes,EdgeBary,Left,Right
Use Hydrodynamic !, Only: Smallm,ElSmallm,ElCapitalM,CapitalM,Z,DZj,Ze,DZi,uNode,uxy,Ze,DZi,w
Implicit none
Real, parameter :: small1=1e-6
integer, intent(in) :: ndelt, FuFw_flag, iLayer
Real, intent(inout) :: dtb, dt
integer, intent(inout) :: nnel,jlev, id0
Real, intent(inout) :: uuint,vvint,wwint,x0,y0,z0, hhint
Real, intent(out) :: xt,yt,zt
Real :: uuNode(3,9), vvNode(3,9), wwNode(3,9), uuNodet(3,9), vvNodet(3,9), wwNodet(3,9), xxNode(3,9), yyNode(3,9), zzNode(3,9), hhNode(3,9)
Real :: uuBtrack, vvBtrack, wwBtrack, hhBTrack, uuBtrack2, vvBtrack2, wwBtrack2, theta, Umax, Umin
Real :: staint(8),t_xi(4),s_xi(4),sig(4),subrat(4)
Integer:: IntersectFlag,nel,nnelIni,idt,iflqs1,i,j,l,nd,nn,lev, jjlev,n1,n2,n3,n4, n5, n6, n7, n8, n9,iNode1,iNode2,iEdge,iNode, Nodes(9), addLayer, N, S, Face, r, ll, nnel0,psi_flag,nel_j
Real:: trat,zup,zrat,aa1,aa2,aa3,aa4,aa,csi,etta,wdown,wup,weit, xtaux, ytaux, ztaux, dtaux
Real:: NearZero = 1e-10 !< Small Number
Real:: talx, taly, talz, tal, timeAcum, dtin,dtbCFL, maxZbed, minZeta
Integer:: i34 = 4
type(MeshGridParam) :: MeshParam
type(HydrodynamicParam) :: HydroParam
Integer:: ELM_flag = 1
Integer :: Interpolate_Flag = 1
Real:: xaux,yaux,zaux
jjlev = jlev !Initial Edge
nnelIni = nnel !Initial Element
nnel0 = nnel !Currently Element
dtaux = dtb
timeAcum = 0.d0 !Integration time for particle tracking
IntersectFlag = 0 !Flag to check if particle intersects bottom
Do While (timeAcum < dt)
!Posição da partícula no primeiro subpasso de tempo
xt=x0-dtb*uuint
yt=y0-dtb*vvint
zt=z0-dtb*wwint
!Call RK4order(uuint, vvint, wwint, dtb, nnel, id0, jlev, xt, yt, zt, x0, y0, z0, dt, timeAcum + dtb, Interpolate_Flag, HydroParam, MeshParam)
dtin = dtb
Call quicksearch(1,nnel,jlev,dtb,dtin,x0,y0,z0,xt,yt,zt,iflqs1,idt,id0,i34,uuint,vvint,wwint,IntersectFlag,nel_j,HydroParam,MeshParam)
dtb = dtin
if(isnan(dtb)) Then
continue
endif
uuBtrack = 0.d0
vvBtrack = 0.d0
wwBtrack = 0.d0
If (IntersectFlag == 0) Then
If (ELM_flag==0) Then !iBilinear Interpolation
timeAcum = timeAcum + dtb
If (FuFw_flag==0) Then !grid for U velocitys
!Call FuVelocities2(uuNode(:,:), vvNode(:,:), wwNode(:,:), xxNode(:,:), yyNode(:,:), zzNode(:,:), nnel,jlev, xt,yt,zt,x0,y0,z0, id0, HydroParam,MeshParam)
Call FuVelocities3(uuNode(:,:), vvNode(:,:), wwNode(:,:), xxNode(:,:), yyNode(:,:), zzNode(:,:), nnel, jlev, xt,yt,zt,x0,y0,z0, id0, HydroParam,MeshParam)
Else !grid for W velocitys
Call FwVelocities2(uuNode(:,:), vvNode(:,:), wwNode(:,:), xxNode(:,:), yyNode(:,:), zzNode(:,:), nnel, jlev, jjlev, xt,yt,zt,x0,y0,z0, HydroParam,MeshParam)
EndIf
Call iBilinear2(uuBtrack, vvBtrack, wwBtrack, uuNode(:,:), vvNode(:,:), wwNode(:,:), xxNode(:,:), yyNode(:,:), zzNode(:,:), xt, yt, zt, x0, y0, z0, id0, nnel, FuFw_flag, MeshParam)
ElseIf (ELM_flag==1) Then !iQuadratic Interpolation
timeAcum = timeAcum + dtb
! Get nodals' velocities and positions:
Call iQuadraticNodes(uuNode(:,:), vvNode(:,:), wwNode(:,:), uuNodet(:,:), vvNodet(:,:), wwNodet(:,:), xxNode(:,:), yyNode(:,:), zzNode(:,:), nnel, jlev, HydroParam, MeshParam)
!Get max value to bed elevation between nine nodes in cell:
maxZbed = maxval(zzNode(1,:))
!Get min value to bed elevation between nine nodes in cell:
minZeta = minval(zzNode(3,:))
!Checking to keep particle z position in the cell limits. If zt is not in the z interval inside cell (out of bounds),
!zt is brought to the z limits:
zt = min(minZeta,max(zt,maxZbed))
If (HydroParam%eta(nnel) - HydroParam%hb(nnel) < HydroParam%PCRI + NearZero) Then
timeAcum = dt
Else
Call iQuadratic(uuBtrack, vvBtrack, wwBtrack, uuNode(:,:), vvNode(:,:), wwNode(:,:), xxNode(:,:), yyNode(:,:), zzNode(:,:), xt, yt, zt)
EndIf
EndIf
Else
timeAcum = dt
EndIf
uuint = uuBtrack
vvint = vvBtrack
wwint = wwBtrack
If(iflqs1.eq.1) exit
x0=xt
y0=yt
z0=zt
!Adaptative sub-time step:
If (nnel /= nnel0 .and. timeAcum < dt ) Then
dtb = dtaux
tal = min(MeshParam%dx/abs(uuint), MeshParam%dy/abs(vvint), (zzNode(3,5) - zzNode(1,5))/abs(wwint))
dtbCFL = 1/HydroParam%CFL(nnel)
if (tal > 0 ) then
dtb = min(tal,dt,dtb)
!dtb = min(tal,dt,dtb,dtbCFL)
endif
nnel0 = nnel
If (dtb + timeAcum > dt) Then
dtb = dt - timeAcum
EndIf
Elseif( uuint**2 + vvint**2 + wwint**2 == 0.0d0 ) Then
timeAcum = dt
EndIf
if (IntersectFlag == 1) Then
timeAcum = dt
uuint = HydroParam%uxy(jlev,1,MeshParam%Edge(nel_j,nel))
vvint = HydroParam%uxy(jlev,2,MeshParam%Edge(nel_j,nel))
EndIf
EndDo
Return
End Subroutine btrack
Subroutine quicksearch(iloc,nnel,jlev,dtb,dtin,x0,y0,z0,&
&xt,yt,zt,nfl,idt,id0,i34,uuint,vvint,wwint,IntersectFlag,nel_j,HydroParam,MeshParam)
!> Straightline search algorithm.
!>\note Initially nnel is an element that encompasses the point P0(x0,y0).\n
!> Input: iloc,nnel,x0,y0,z0,xt,yt,zt,jlev, time, and vt2, ww2 for abnormal cases;\n
!> Output: the updated end pt (xt,yt,zt) (if so), nnel, jlev, a flag nfl.\n
!> Exit btrack if a bnd or dry element is hit and vel. there is small, or death trap is reached.
!>\param iloc nudge initial pt (iloc=0: do not nudge initial pt; iloc=1: nudge)
!>\param dtb sub step time in seg (input)
!>\param nnel number of the element that enclosed the point P0 (input)
!>\param jlev number of the layer that enclosed the point P0 (input)
!>\param x0 x initial position of the point P0 (input)
!>\param y0 y initial position of the point P0(input)
!>\param z0 z initial position of the point P0 (input)
!>\param xt x final position of the point Pt (input)
!>\param yt y final position of the point Pt (input)
!>\param zt z final position of the point Pt (input)
!>\return nfl: The element that enclosed the point Pt was founded (output)
!>\return nnel: number of the element that enclosed the point Pt (output)
!>\return jlev: number of the layer that enclosed the point Pt (output)
!>\author Carlos Ruberto Fragoso
!>\attention List of Modifications: \n
!> - 15.09.2015: Routine Implementation (Carlos Ruberto Fragoso)
Use MeshVars !, Only:Quadri,Edge,EdgeDef,xNode,yNode,Area,xb,yb,Left,Right,Neighbor,EdgeBary
Use Hydrodynamic !, Only: ElSmallm,ElCapitalM,Ze,H,Pcri,uxy,uNode
Implicit none
Real, parameter :: small1=1e-5
Integer, intent(in) :: iloc,idt,i34
Real, intent(in) :: dtb,x0,y0,z0
Integer, intent(out) :: nfl
Integer, intent(inout) :: nnel,jlev,id0
Real:: xpoly(4),ypoly(4)
Real, intent(inout) :: xt,yt,zt
Real:: trm,aa,aa1,ae,xcg,ycg,pathl,xin,yin,zin,tt1,tt2,dist,xvel,yvel,zvel,hvel
Real:: uuint,vvint,wwint,dtin
Integer:: IntersectFlag
Integer:: nel,i,j,k,n1,n2,nel_j,iflag,it,md1,md2,lit,k1,k2,jd1,jd2,r,isd,nel0,INOUT
Real:: NearZero = 1e-10 !< Small Number
Integer :: NWater
type(MeshGridParam) :: MeshParam
type(HydrodynamicParam) :: HydroParam
nel0 = nnel
nfl=0
trm=dtb !time remaining
nel_j = id0
nel = nnel
aa=0
aa1=0
!The ideia is that any polygon area can be represents by sum of multiples triangules.
!Next are calculated the multiples triangules area formed by P0 and iEdges nodes points (acumulated in aa variable)
!and Pt and iEdges points (acumulated in aa1 variable).
Do i=1,i34
!Get the nodes for each edge
n1=MeshParam%Quadri(MeshParam%EdgeDef(1,i),nel) + 1
n2=MeshParam%Quadri(MeshParam%EdgeDef(2,i),nel) + 1
aa=aa+dabs(signa(MeshParam%xNode(n1),MeshParam%xNode(n2),x0,MeshParam%yNode(n1),MeshParam%yNode(n2),y0))
aa1=aa1+dabs(signa(MeshParam%xNode(n1),MeshParam%xNode(n2),xt,MeshParam%yNode(n1),MeshParam%yNode(n2),yt))
EndDo !i
!If P0 is inside on the Element(nnel) iEdge barycenter, then the acumalated area calculate in previous step is equal to %Area(iElement).
!Check if P0(x0,y0,z0) is enclosed in nel:
ae=dabs(aa-MeshParam%Area(nel))/MeshParam%Area(nel)
If(ae.gt.small1) Then
print*,'(x0,y0) not in nnel initially',ae,nnel,id0
pause
stop
EndIf
!If PT is inside the Element(nnel), then the acumalated area calculate in previous step is lower than %Area(iElement).
!In this case, nnel was finded and need to find iLayer inside nnel (go to 400)
ae=dabs(aa1-MeshParam%Area(nel))/MeshParam%Area(nel)
if(ae.lt.small1) Then
nnel=nel
go to 400
endif
!(xt,yt) not in nel, and thus (x0,y0) and (xt,yt) are distinctive
!An interior pt close to (x0,y0) to prevent underflow for iloc >=1.
If(iloc.eq.0) Then
!xcg= MeshParam%EdgeBary(1,id0) !x0
!ycg= MeshParam%EdgeBary(2,id0) !y0
xcg = x0 !MeshParam%EdgeBary(1,id0) !x0
ycg = y0 !MeshParam%EdgeBary(2,id0) !y0
Elseif(iloc.eq.1) Then
!xcg=(1-1.0d-4)*MeshParam%EdgeBary(1,id0)+1.0d-4*MeshParam%xb(nel)
!ycg=(1-1.0d-4)*MeshParam%EdgeBary(2,id0)+1.0d-4*MeshParam%yb(nel)
xcg = (1-1.0d-4)*x0+1.0d-4*MeshParam%xb(nel)
ycg = (1-1.0d-4)*y0+1.0d-4*MeshParam%yb(nel)
! |x0-xcg|/|x0-xctr|=1.0d-3
endif
pathl=dsqrt((xt-xcg)**2+(yt-ycg)**2)
If(xcg.eq.xt.and.ycg.eq.yt.or.pathl.eq.0) Then
print*,'Zero path',x0,y0,xt,yt,xcg,ycg
pause
stop
endif
! Check if the particle position crosses some Edge from Element(nel):
! Starting edge nel_j
Do i=1,i34
!Get the nodes for each edge
n1=MeshParam%Quadri(MeshParam%EdgeDef(1,i),nel) + 1
n2=MeshParam%Quadri(MeshParam%EdgeDef(2,i),nel) + 1
!Check intersection between lines segments, one formed by nodes points (Edge) and another formed by Pt(xt,yt) and Pcg(xcg,ycg):
call intersect2(xcg,xt,MeshParam%xNode(n1),MeshParam%xNode(n2),ycg,yt,MeshParam%yNode(n1),MeshParam%yNode(n2),iflag,xin,yin,tt1,tt2)
!If has a intersection point PI(xin, yin), the Edge i is crosses in particle trajectory:
If(iflag.eq.1) Then
nel_j=i
!Find the moment which the particle crosses the edge:
if(sqrt((x0-xin)**0.5+(y0-yin)**0.5)>small1 .and. abs(sqrt(uuint**0.5+vvint**0.5))>0 )then
dtin = sqrt((x0-xin)**0.5+(y0-yin)**0.5)/sqrt(uuint**0.5+vvint**0.5)
endif
Go to 399
Endif
Enddo !i=1,3
!No one Edge in element (nel) was crossed by particle:
If(iflag.eq.0) Then
xt = x0
yt = y0
zt = z0
return
EndIf
print*,'Found no intersecting edges'
print*,'xcg',xcg
print*,'xt',xt
print*,'ycg',ycg
print*,'yt',yt
print*,'Elem',nel
pause
stop
399 continue
zin=z0 !initialize
it=0
loop4: Do
!----------------------------------------------------------------------------------------
it=it+1
If(it.gt.1000) Then
!write(12,*)'Death trap reached',idt,id0
nfl=1
xt=xin
yt=yin
zt=zin
nnel=nel
Exit loop4
Endif
md1 = MeshParam%Quadri(MeshParam%EdgeDef(1,nel_j),nel) + 1
md2 = MeshParam%Quadri(MeshParam%EdgeDef(2,nel_j),nel) + 1
!Compute z position
dist=dsqrt((xin-xt)**2+(yin-yt)**2)
if(dist/pathl.gt.1+1.0d-4) Then
print*,'Path overshot'
pause
stop
endif
zin=zt-dist/pathl*(zt-zin)
trm=trm*dist/pathl !time remaining
pathl=dsqrt((xin-xt)**2+(yin-yt)**2)
If(pathl.eq.0.or.trm.eq.0) Then
print*,'Target reached'
pause
stop
endif
lit=0 !flag
!For horizontal exit and dry elements, compute tangential vel.,
!update target (xt,yt,zt) and continue.
isd = MeshParam%Edge(nel_j,nel)
id0 = isd
r = MeshParam%Right(isd)
if (nel==179) then
continue
endif
!If particle cross iEdge from iElement(nel), so the new iElement is the neighbour Element which share iEdge(nel_j).
!However, it can a abnormal case which either iEdge(isd) no has neighbour (horizontal exit or wall) or is dry:
If (r == 0 .or. HydroParam%H(isd)-HydroParam%hj(isd)<=HydroParam%Pcri+NearZero) Then
lit=1
!Nudge intersect (xin,yin), and update starting pt
xin=(1-1.0d-4)*xin+1.0d-4*MeshParam%xb(nel)
yin=(1-1.0d-4)*yin+1.0d-4*MeshParam%yb(nel)
xcg=xin
ycg=yin
!Set tang. velocities:
xvel= 0.!uxy(jlev,1,isd)
yvel= 0.!uxy(jlev,2,isd)
zvel=0.5*((HydroParam%uNode(jlev,3,md1)+HydroParam%uNode(jlev,3,md2))/2. + (HydroParam%uNode(jlev+1,3,md1)+HydroParam%uNode(jlev+1,3,md2))/2.)
!Update Pt:
xt=xin-xvel*trm
yt=yin-yvel*trm
zt=zin-zvel*trm
!Horizontal velocity magnitude:
hvel=dsqrt(xvel**2+yvel**2)
If(hvel.lt.1.e-4) Then !Checar essa condi��o, todos os casos entram aqui CAYO
nfl=1
xt=xin
yt=yin
zt=zin
nnel=nel
IntersectFlag = 1
exit loop4
EndIf
pathl=hvel*trm
ElseIf( MeshParam%Neighbor(nel_j,nnel) .ne. 0 ) Then
If (HydroParam%eta(MeshParam%Neighbor(nel_j,nnel)) - HydroParam%hb(MeshParam%Neighbor(nel_j,nnel)) <= HydroParam%PCRI + NearZero) Then !CAYO
lit=1
!Nudge intersect (xin,yin), and update starting pt
xin=(1-1.0d-4)*xin+1.0d-4*MeshParam%xb(nel)
yin=(1-1.0d-4)*yin+1.0d-4*MeshParam%yb(nel)
xcg=xin
ycg=yin
!Set tang. velocities:
xvel = 0.!uxy(jlev,1,isd)
yvel = 0.!uxy(jlev,2,isd)
zvel = 0.5*((HydroParam%uNode(jlev,3,md1)+HydroParam%uNode(jlev,3,md2))/2. + (HydroParam%uNode(jlev+1,3,md1)+HydroParam%uNode(jlev+1,3,md2))/2.)
!Update Pt:
xt=xin-xvel*trm
yt=yin-yvel*trm
zt=zin-zvel*trm
!Horizontal velocity magnitude:
hvel=dsqrt(xvel**2+yvel**2)
If(hvel.lt.1.e-4) Then !Checar essa condi��o, todos os casos entram aqui CAYO
nfl=1
xt=xin
yt=yin
zt=zin
nnel=nel
IntersectFlag = 1
exit loop4
EndIf
pathl=hvel*trm
ElseIf(dmax1(zt,HydroParam%hb(nel0)) < HydroParam%hb(MeshParam%Neighbor(nel_j,nnel)) ) Then
nnel = nel0
xt=(1-1.0d-4)*MeshParam%EdgeBary(1,MeshParam%Edge(nel_j,nel0)) + 1.0d-4*MeshParam%xb(nel0)
yt=(1-1.0d-4)*MeshParam%EdgeBary(2,MeshParam%Edge(nel_j,nel0)) + 1.0d-4*MeshParam%yb(nel0)
nfl=1
IntersectFlag = 1
exit loop4
EndIf
EndIf
!Else in normal cases, we need get the neighbour element which shares the iEdge(nel_j):
If(lit.eq.0) Then
!next front element
nel = MeshParam%Neighbor(nel_j,nel)
EndIf
!With the element updated, we check again if Pt(xt,yt) is inside the element:
aa=0
Do i=1,i34
k1=MeshParam%Quadri(MeshParam%EdgeDef(1,i),nel) + 1
k2=MeshParam%Quadri(MeshParam%EdgeDef(2,i),nel) + 1
aa=aa+dabs(signa(MeshParam%xNode(k1),MeshParam%xNode(k2),xt,MeshParam%yNode(k1),MeshParam%yNode(k2),yt))
EndDo !i
!If is inside, the acumulated area (aa) is lower than %Area(iElement):
ae = dabs(aa-MeshParam%Area(nel))/MeshParam%Area(nel)
If(ae.lt.small1) Then
nnel=nel
! Element find, xt and yt position defined. Go to find zt positon (Go to 400).
Exit loop4
EndIf
!The particle is out the element, find next intersecting edge:
Do j=1,i34
!Get the nodes for each edge:
jd1=MeshParam%Quadri(MeshParam%EdgeDef(1,j),nel) + 1 !nm(nel,nx(i34(nel),j,1))
jd2=MeshParam%Quadri(MeshParam%EdgeDef(2,j),nel) + 1 !nm(nel,nx(i34(nel),j,2))
if(jd1.eq.md1.and.jd2.eq.md2.or.jd2.eq.md1.and.jd1.eq.md2) cycle !iEdge shared with previous element that was checked, skip to next iEdge.
!Check path intersection:
call intersect2(xcg,xt,MeshParam%xNode(jd1),MeshParam%xNode(jd2),ycg,yt,MeshParam%yNode(jd1),MeshParam%yNode(jd2),iflag,xin,yin,tt1,tt2)
!If has a intersection point PI(xin, yin), the Edge i is crosses in particle trajectory:
if(iflag.eq.1) then
nel_j=j !next front edge
if(sqrt((x0-xin)**0.5+(y0-yin)**0.5)>small1 .and. abs(sqrt(uuint**0.5+vvint**0.5))>0 )then
dtin = sqrt((x0-xin)**0.5+(y0-yin)**0.5)/sqrt(uuint**0.5+vvint**0.5)
endif
cycle loop4
endif
EndDo !j
!print*,'Failed to find next edge',lit,xin,yin,xt,yt,nel,&
!&md1,md2,idt
xt = x0
yt = y0
Exit loop4
pause
stop
EndDo loop4
400 Continue
!The element (nnel) was found, now we are looking for the iLayer (jlev) that contains Pt(xt,yt,zt).
!First zt is set as the maximum value between Ze (smallm)lower layer and zt, this ensure that the particle no set in position
!under the bounds of vertical grid discretization. In cases which the zt is set as Ze(smallm) layers, this implies that the
!particle reach to bottom.
!Next zt is set as the minimum between it and Ze(CapitalM+1) upper layer, this ensure that particle not set in a position
!above the bounds of vertical grid discretization.
!zt is set as the minimum value between max[zt calculed and Lower Layer Level in the Element (nnel)] and Upper Layer in Element(nnel)
!We want find the iLayer that includes this value, not specfic point.
!If (nnel /= nel0) Then
! !zt = dmax1(zt,HydroParam%Ze(HydroParam%ElSmallm(nel0),nel0)+sum(HydroParam%DZsi(:,nel0)))
! zt = dmax1(zt,HydroParam%hb(nel0))
! !If (zt < HydroParam%Ze(HydroParam%ElSmallm(nnel),nnel) - HydroParam%hb(nnel) ) then
! If (zt < HydroParam%hb(nnel) ) then
! nnel = nel0
! xt = MeshParam%EdgeBary(1,MeshParam%Edge(nel_j,nel0))
! yt = MeshParam%EdgeBary(2,MeshParam%Edge(nel_j,nel0))
! EndIf
!EndIf
!zt = dmin1(dmax1(zt,HydroParam%Ze(HydroParam%ElSmallm(nnel),nnel)+sum(HydroParam%DZsi(:,nnel))),HydroParam%Ze(HydroParam%ElCapitalM(nnel)+1,nnel))
zt = dmin1(dmax1(zt,HydroParam%Ze(HydroParam%ElSmallm(nnel),nnel)+sum(HydroParam%DZsi(:,nnel))),HydroParam%eta(nnel))
Do k = HydroParam%ElSmallm(nnel), HydroParam%ElCapitalM(nnel)
If (zt.gt.HydroParam%Ze(k,nnel).and.zt.le.HydroParam%Ze(k+1,nnel)) Then
jlev = k
EndIf
EndDo
Return
End Subroutine quicksearch
Subroutine ComputeFW(MeshParam,HydroParam,dt)
! Compute the Explicit Finite Difference Operator Fw
! Based on:
! [1] Wang,B.; Zhao,G.; Fringer,O.B. Reconstruction of vector fields for semi-Lagrangian advection on unstructured, staggered grids.
! Ocean Modelling, 40, p. 52-71, 2011.
! [2] Zhang, Y.; Baptista, A.M. SELFE: A semi-implicit Eulerian�Lagrangian finite-element model for cross-scale ocean circulation.
! Ocean Modelling, 21, p. 71-96, 2008.
! [3] Fringer, O.B.; Gerritsen, M.; Street, R.L. An unstructured-grid, finite-volume, nonhydrostatic, parallel coastal ocean simulator.
! Ocean Modelling, 14, p. 139-173, 2006.
! [4] Cheng, Ralph T.; Casulli, Vincenzo; Gartner, Jeffrey W. Tidal, Residual, Intertidal Mudflat (TRIM) Model and its Applications to San Francisco Bay, California.
! Estuarine, Coastal and Shelf Sciences, 36, p. 235-280. 1993.
! Input:
! Edge-Based Velocity Vector
! Vertical Velocity Field
! Output:
! Advected Vertical Velocity Field
! List of Modifications:
! 28.06.2017: Routine Implementation (J. Rafael Cavalcanti)
! Note:
! -> I'm neglecting the viscosity effects in the Vertical Momentum transport
! Programmer: J. Rafael Cavalcanti
Use Hydrodynamic !, Only: nElem, KMax, ElSmallmOld, ElCapitalMOld, hb, nNode, Ze, nEdge
Use MeshVars
!Use Hydrodynamic, Only: InCircle, Edge, xb, yb, tri, x, y, Edgelength, Area, Dzin
!Use Hydrodynamic, Only: eta, etan, w, u, un, uEdge
!Use Param, Only: Pcri, dt, Theta
Implicit None
Integer:: iElem,iLayer,jlev, jjlev,l,r,nnel,iNode,j,f,FuFw_flag, TrajectoryFlag, ndels, iEdge, Face
Real:: x0,y0,z0,xt,yt,zt,uuint,vvint,wwint,hhint,wdown,wup,vmag,dtb
Real:: NearZero = 1e-10
Real:: dt
type(MeshGridParam) :: MeshParam
type(HydrodynamicParam) :: HydroParam
Integer:: id0=1
HydroParam%Fw = HydroParam%w
Do iElem = 1,MeshParam%nElem
If ( HydroParam%eta(iElem) - HydroParam%hb(iElem) < HydroParam%PCRI/2. ) Then
HydroParam%Fw(:,iElem) = 0.
Cycle
EndIf
do iEdge = 1,4
Face = MeshParam%Edge(iEdge,iElem)
if (HydroParam%IndexInflowEdge(Face)>0.or.HydroParam%IndexInflowEdge(Face)>0) then
HydroParam%Fw(:,iElem) = HydroParam%w(:,iElem)
cycle
endif
enddo
Do iLayer = HydroParam%ElSmallm(iElem)+1, HydroParam%ElCapitalM(iElem)
!Do iEdge = 1,4
! if (MeshParam%Neighbor(iEdge,iElem)/=0) then
! If (iLayer<=HydroParam%ElSmallm(MeshParam%Neighbor(iEdge,iElem))) Then
! HydroParam%Fw(iLayer,iElem) = HydroParam%w(iLayer,iElem)
! continue
! return
! EndIf
! else
! !HydroParam%Fw(iLayer,iElem) = HydroParam%w(iLayer,iElem)
! !continue
! !return
! Endif
!EndDo
nnel = iElem
jlev = iLayer
x0 = MeshParam%xb(iElem)
y0 = MeshParam%yb(iElem)
z0 = HydroParam%Ze(iLayer,iElem)
uuint = HydroParam%uxyL(iLayer,1,iElem)
vvint = HydroParam%uxyL(iLayer,2,iElem)
wwint = HydroParam%w(iLayer,iElem)
vmag = dsqrt(uuint**2+vvint**2+wwint**2)
If(vmag.le.1.e-6) Then ! No activity
HydroParam%Fw(iLayer,iElem) = HydroParam%w(iLayer,iElem)
Else !do backtrack
TrajectoryFlag = 0 ! = 0 -> System Defined Intervals to Integrate Particle Trajectory; = 1 -> User Defined Intervals to Integrate Particle Trajectory
! 3.5 Find the Local Time Step to Integrate Trajectory
If ( TrajectoryFlag == 0 ) Then ! System Calculation - [7,8]
dtb = MinVal(MeshParam%InCircle)/vmag ! ( InCircle(lElem)/Sqrt( Veloc(1)**2. + Veloc(2)**2. ), InCircle(r)/Sqrt( Veloc(1)**2. + Veloc(2)**2. ) )
ndels = Max(Floor(dt/dtb),HydroParam%NFUT)
dtb = dt/ndels !sub-step in backtracking
Else If ( TrajectoryFlag == 1 ) Then ! User Defined sub-steps numbers
ndels = HydroParam%NFUT
dtb = dt/ndels !sub-step in backtracking
EndIf
FuFw_flag = 1
Call btrack(ndels,dtb,dt,uuint,vvint,wwint,hhint,&
&x0,y0,z0,xt,yt,zt,nnel,jlev,id0,HydroParam,MeshParam,FuFw_flag, iLayer, iEdge)
HydroParam%Fw(iLayer,iElem) = wwint
Endif
EndDo
EndDo
Return
End Subroutine ComputeFW
Subroutine intersect2(x1,x2,x3,x4,y1,y2,y3,y4,iflag,xin,yin,tt1,tt2)
!> Subroutine to detect if two segments (1,2) and (3,4) have a common point
!>\note The 4 pts are distinctive.\n
!> The eqs. for the 2 lines are: X=X1+(X2-X1)*tt1 and X=X3+(X4-X3)*tt2.\n
!> Output: iflag: 0: no intersection or colinear; 1: exactly 1 intersection..\n
!> If iflag=1, (xin,yin) is the intersection.
!>\param iflag A common point was founded (iflag=0: no; iflag=1: yes)
!>\param x1 x position of the point P1 (input)
!>\param x2 x position of the point P2 (input)
!>\param x3 x position of the point P3 (input)
!>\param x4 x position of the point P4 (input)
!>\param y1 y position of the point P1 (input)
!>\param y2 y position of the point P2 (input)
!>\param y3 y position of the point P3 (input)
!>\param y4 y position of the point P4 (input)
!>\return xin: x position in the intersection (output)
!>\return yin: y position in the intersection (output)
!>\return tt1: angular coefficient of the segment formmed by P1 and P2
!>\return tt2: angular coefficient of the segment formmed by P3 and P4
!>\author Carlos Ruberto Fragoso
!>\attention List of Modifications: \n
!> - 15.09.2015: Routine Implementation (Carlos Ruberto Fragoso)
Implicit none
Real, parameter :: small2=0.0 !small positive number or 0
Real, intent(in) :: x1,x2,x3,x4,y1,y2,y3,y4
Integer, intent(out) :: iflag
Real, intent(out) :: xin,yin,tt1,tt2
Real:: delta,delta1,delta2
tt1=-1000
tt2=-1000
iflag=0
delta=(x2-x1)*(y3-y4)-(y2-y1)*(x3-x4)
delta1=(x3-x1)*(y3-y4)-(y3-y1)*(x3-x4)
delta2=(x2-x1)*(y3-y1)-(y2-y1)*(x3-x1)
If(delta.ne.0.0d0) Then
tt1=delta1/delta
tt2=delta2/delta
If(tt1.ge.-small2.and.tt1.le.1+small2.and.tt2.ge.-small2.and.tt2.le.1+small2) Then
iflag=1
xin=x1+(x2-x1)*tt1
yin=y1+(y2-y1)*tt1
Endif
Endif
Return
End Subroutine intersect2
Subroutine ibilinear(elem,x1,x2,x3,x4,y1,y2,y3,y4,x,y,xi,eta,shapef)
!> Inverse bilinear mapping for quadrangles
!>\note Convexity of the quad must have been checked, and (x,y) must not be outside the quad.
!>\param elem Number of element (input)
!>\param x1 x position of the point P1 (input)
!>\param x2 x position of the point P2 (input)
!>\param x3 x position of the point P3 (input)
!>\param x4 x position of the point P4 (input)
!>\param y1 y position of the point P1 (input)
!>\param y2 y position of the point P2 (input)
!>\param y3 y position of the point P3 (input)
!>\param y4 y position of the point P4 (input)
!>\param x x position in the interested point (input)
!>\param y y position in the interested point (input)
!>\return shapef: weighted coefficients in each vertices (output)
!>\author Carlos Ruberto Fragoso
!>\attention List of Modifications: \n
!> - 15.09.2015: Routine Implementation (Carlos Ruberto Fragoso)
Implicit None
Real, parameter:: small1=1e-6
Real, parameter:: small3=1.e-5
Real, intent(in) :: elem,x1,x2,x3,x4,y1,y2,y3,y4,x,y !elem for debugging only
Real, intent(out) :: xi,eta,shapef(4)
Real:: x0,y0,axi,aet,bxy,root_xi,root_et,dxi,deta,dd,beta,gamma,delta
Integer:: icaseno,icount,i,j
dimension axi(2),aet(2),bxy(2),root_xi(2),root_et(2)
!Consts.
x0=(x1+x2+x3+x4)/4
y0=(y1+y2+y3+y4)/4
axi(1)=x2-x1+x3-x4
axi(2)=y2-y1+y3-y4
aet(1)=x3+x4-x1-x2
aet(2)=y3+y4-y1-y2
bxy(1)=x1-x2+x3-x4
bxy(2)=y1-y2+y3-y4
dxi=2*((x3-x4)*(y1-y2)-(y3-y4)*(x1-x2))
deta=2*((x4-x1)*(y3-y2)-(y4-y1)*(x3-x2))
!Inverse mapping
If(dabs(bxy(1))<small3.and.dabs(bxy(2))<small3.or.dabs(dxi)<small3.and.dabs(deta)<small3) Then
icaseno=1
!print*, 'Entering case 1'
dd=axi(1)*aet(2)-axi(2)*aet(1)
if(dd==0) then
print*,'Case 1 error:',dd
pause
stop
EndIf
xi=4*(aet(2)*(x-x0)-aet(1)*(y-y0))/dd
eta=4*(axi(1)*(y-y0)-axi(2)*(x-x0))/dd
Elseif(dabs(dxi)<small3.and.dabs(deta)>=small3) Then
icaseno=2
!print*, 'Entering case 2'
eta=4*(bxy(2)*(x-x0)-bxy(1)*(y-y0))/deta
dd=(axi(1)+eta*bxy(1))**2+(axi(2)+eta*bxy(2))**2
If(dd==0) then
print*,'Case 2 error:',dd
pause
stop
EndIf
xi=((4*(x-x0)-eta*aet(1))*(axi(1)+eta*bxy(1))+(4*(y-y0)-eta*aet(2))*(axi(2)+eta*bxy(2)))/dd
Elseif(dabs(dxi)>=small3.and.dabs(deta)<small3) Then
icaseno=3
! print*, 'Entering case 3'
xi=4*(bxy(2)*(x-x0)-bxy(1)*(y-y0))/dxi
dd=(aet(1)+xi*bxy(1))**2+(aet(2)+xi*bxy(2))**2
If(dd==0) Then
print*,'Case 3 error:',dd
pause
stop
EndIf
eta=((4*(x-x0)-xi*axi(1))*(aet(1)+xi*bxy(1))+(4*(y-y0)-xi*axi(2))*(aet(2)+xi*bxy(2)))/dd