-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathwwm_aux.F90
2991 lines (2763 loc) · 98.9 KB
/
wwm_aux.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
#include "wwm_functions.h"
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE GRADDEP()
USE DATAPOOL
IMPLICIT NONE
INTEGER :: IP
REAL(rkind) :: GDL, GDD
DDEP(:,:) = 0.0
SELECT CASE (DIMMODE)
CASE (1)
#ifdef MPI_PARALL_GRID
call parallel_abort('WWM - 1d mode cannot work with SCHISM ')
#endif
CALL DIFFERENTIATE_XDIR(DEP,DDEP(:,1))
IF (LSPHE) THEN
DO IP = 1, MNP
DDEP(IP,1) = DDEP(IP,1)/( DEGRAD*REARTH*COS(YP(IP)*DEGRAD) )
END DO
END IF
IF (LTEST .AND. ITEST > 100) THEN
WRITE(STAT%FHNDL,*) 'Gradients of depth '
WRITE(STAT%FHNDL,*) '@D/@X '
DO IP = 1, MNP
WRITE(STAT%FHNDL,'(1X,I5,3F10.5)') IP, DDEP(IP,1)
END DO
FLUSH(STAT%FHNDL)
END IF
CASE (2)
CALL DIFFERENTIATE_XYDIR(DEP,DDEP(:,1),DDEP(:,2))
IF (LSPHE) THEN
DO IP = 1, MNP
DDEP(IP,1) = DDEP(IP,1)/( DEGRAD*REARTH*COS(YP(IP)*DEGRAD) )
DDEP(IP,2) = DDEP(IP,2)/( DEGRAD*REARTH )
END DO
END IF
IF (LSLOP) THEN
DO IP = 1, MNP
GDL = SQRT(DDEP(IP,1)**2 + DDEP(IP,2)**2) ! Achtung Gradient mit SQRT berechnet ...
GDD = MyATAN2(DDEP(IP,2), DDEP(IP,1))
IF (GDL < 1.0E-8) CYCLE
GDL = SQRT(GDL)
IF (GDL > SLMAX) THEN
DDEP(IP,1) = SLMAX*COS(GDD)
DDEP(IP,2) = SLMAX*SIN(GDD)
WRITE(STAT%FHNDL,*) IP, SLMAX, GDL, GDD , 'MAXSLOPE'
END IF
END DO
FLUSH(STAT%FHNDL)
END IF
CASE DEFAULT
END SELECT
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
#ifdef TIMINGS
SUBROUTINE WAV_MY_WTIME(wtime)
USE DATAPOOL, only : rkind
implicit none
real(rkind), intent(inout) :: wtime
# ifdef MPI_PARALL_GRID
real(8) mpi_wtime
wtime=mpi_wtime()
# else
CALL cpu_time(wtime)
# endif
END SUBROUTINE
#endif
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE GRADCURT()
USE DATAPOOL
IMPLICIT NONE
INTEGER :: IP
DCUX(:,:) = 0.0
DCUY(:,:) = 0.0
SELECT CASE (DIMMODE)
CASE (1)
CALL DIFFERENTIATE_XDIR(CURTXY(:,1),DCUX(:,1))
CALL DIFFERENTIATE_XDIR(CURTXY(:,2),DCUY(:,1))
IF (LSPHE) THEN
DO IP = 1, MNP
DCUX(IP,1) = DCUX(IP,1)/( DEGRAD*REARTH*COS(YP(IP)*DEGRAD) )
DCUY(IP,1) = DCUY(IP,1)/( DEGRAD*REARTH*COS(YP(IP)*DEGRAD) )
END DO
END IF
IF (LTEST .AND. ITEST > 100) THEN
WRITE(STAT%FHNDL,*) 'Gradients of depth and current'
WRITE(STAT%FHNDL,*) '@U/@X @V/@X'
DO IP = 1, MNP
WRITE(STAT%FHNDL,'(1X,I5,3F10.5)') IP, DCUX, DCUY
END DO
FLUSH(STAT%FHNDL)
END IF
CASE (2)
CALL DIFFERENTIATE_XYDIR(CURTXY(:,1),DCUX(:,1),DCUX(:,2))
CALL DIFFERENTIATE_XYDIR(CURTXY(:,2),DCUY(:,1),DCUY(:,2))
IF (LSPHE) THEN
DO IP = 1, MNP
DCUX(IP,1) = DCUX(IP,1)/( DEGRAD*REARTH*COS(YP(IP)*DEGRAD) )
DCUY(IP,1) = DCUY(IP,1)/( DEGRAD*REARTH*COS(YP(IP)*DEGRAD) )
DCUX(IP,2) = DCUX(IP,2)/( DEGRAD*REARTH )
DCUY(IP,2) = DCUY(IP,2)/( DEGRAD*REARTH )
END DO
END IF
IF (LTEST .AND. ITEST > 100) THEN
WRITE(STAT%FHNDL,*) 'The Gradient of Depth and Current'
WRITE(STAT%FHNDL,*) ' @U/@X @U/@Y @V/@X @V/@Y'
DO IP = 1, MNP
WRITE(STAT%FHNDL,'(1X,I5,4F15.7)') IP, DCUX(IP,1), DCUX(IP,2), DCUY(IP,1), DCUY(IP,2)
END DO
FLUSH(STAT%FHNDL)
END IF
CASE DEFAULT
END SELECT
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE GRAD_CG_K()
USE DATAPOOL
IMPLICIT NONE
INTEGER :: IS
DCGDX(:,:) = 0.0
DCGDY(:,:) = 0.0
DWKDX(:,:) = 0.0
DWKDY(:,:) = 0.0
SELECT CASE (DIMMODE)
CASE (1)
DO IS = 1, NUMSIG
CALL DIFFERENTIATE_XDIR(CG(IS,:),DCGDX(:,IS))
CALL DIFFERENTIATE_XDIR(WK(IS,:),DWKDX(:,IS))
ENDDO
CASE (2)
DO IS = 1, NUMSIG
CALL DIFFERENTIATE_XYDIR(CG(IS,:),DCGDX(:,IS),DCGDY(:,IS))
CALL DIFFERENTIATE_XYDIR(WK(IS,:),DWKDX(:,IS),DWKDY(:,IS))
ENDDO
CASE DEFAULT
END SELECT
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE SMOOTH( BETA, MNP, XP, YP, VAR )
USE DATAPOOL, ONLY : RKIND, ZERO, ONE, TWO
IMPLICIT NONE
INTEGER :: MNP
REAL(rkind), INTENT(IN) :: XP(MNP), YP(MNP)
REAL(rkind), INTENT(INOUT) :: VAR(MNP)
REAL(rkind), INTENT(IN) :: BETA
REAL(rkind) :: VART(MNP)
REAL(rkind) :: SW, SWQ, DISX, DISY, DIST, DIS
INTEGER :: I, J
DO I = 1, MNP
SW = ZERO
SWQ = ZERO
DO J = 1, MNP
DISX = (XP(I) - XP(J))**2
DISY = (YP(I) - YP(J))**2
DIST = DISX + DISY
IF (DIST > TINY(1.)) THEN
DIS = SQRT(DIST)**BETA
ELSE
DIS = ONE
END IF
SW = SW + DIS
SWQ = SWQ + DIS*VAR(J)
END DO
VART(I) = SWQ / SW
END DO
VAR(:) = VART(:)
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE SMOOTH_V2(VAR_IN, VAR_OUT)
USE DATAPOOL
IMPLICIT NONE
REAL(rkind), INTENT(IN) :: VAR_IN(MNP)
REAL(rkind), INTENT(OUT) :: VAR_OUT(MNP)
INTEGER IP, IADJ, IP_ADJ
REAL(rkind) :: SumVAR, SumSI, eVal
DO IP = 1, NP_RES
SumVAR=SI(IP) * VAR_IN(IP)
SumSI =SI(IP)
DO IADJ=1,VERT_DEG(IP)
IP_ADJ=LIST_ADJ_VERT(IADJ,IP)
SumVAR=SumVAR + SI(IP_ADJ)*VAR_IN(IP_ADJ)
SumSI =SumSI + SI(IP_ADJ)
END DO
eVal=SumVAR/SumSI
VAR_OUT(IP)=eVal
END DO
#ifdef MPI_PARALL_GRID
CALL exchange_p2d(VAR_OUT)
#endif
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE SMOOTH_ON_TRIANGLE(VAR_IN, VAR_OUT, nbIter)
USE DATAPOOL
IMPLICIT NONE
integer, intent(in) :: nbIter
REAL(rkind), INTENT(IN) :: VAR_IN(MNE)
REAL(rkind), INTENT(OUT) :: VAR_OUT(MNE)
INTEGER IE, IEadj, I, nb, iIter
REAL(rkind) :: eVal, eValB
REAL(rkind) :: VARextent(MNEextent)
VAR_OUT=VAR_IN
DO iIter=1,nbIter
VARextent(1:MNE)=VAR_OUT
#ifdef MPI_PARALL_GRID
CALL TRIG_SYNCHRONIZATION(VARextent)
#endif
DO IE=1,MNE
eVal=ZERO
nb=0
DO I=1,3
IEadj=IEneighbor(I,IE)
IF (IEadj .gt. 0) THEN
eVal=eVal + VARextent(IEadj)
nb=nb+1
END IF
END DO
eValB=(MyREAL(6-nb)*VARextent(IE) + eVal)/6.0_rkind
VAR_OUT(IE)=eValB
END DO
END DO
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE DIFFERENTIATE_XDIR(VAR, DVDX)
USE DATAPOOL
IMPLICIT NONE
REAL(rkind), INTENT(IN) :: VAR(MNP)
REAL(rkind), INTENT(OUT) :: DVDX(MNP)
INTEGER :: IP
REAL(rkind) :: TMP1, TMP2
DVDX(1) = (VAR(2)-VAR(1))/(XP(2)-XP(1))
DVDX(MNP) = (VAR(MNP)-VAR(MNP-1))/(XP(MNP)-XP(MNP-1))
DO IP = 2, MNP-1
! DVDX(IP) = ( VAR(IP)-VAR(IP-1) ) / ( XP(IP)-XP(IP-1) )
! DVDX(IP) = (VAR(IP+1)-VAR(IP))/(XP(IP+1)-XP(IP))
! TMP1 = (VAR(IP)-VAR(IP-1))/(XP(IP)-XP(IP-1))
! TMP2 = (VAR(IP+1)-VAR(IP))/(XP(IP+1)-XP(IP))
! DVDX(IP) = 0.5 * (TMP1 + TMP2)
TMP1 = VAR(IP+1)-VAR(IP-1)
TMP2 = XP(IP+1)-XP(IP-1)
DVDX(IP) = TMP1/TMP2
END DO
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE DIFFERENTIATE_XYDIR(VAR, DVDX, DVDY)
USE DATAPOOL
IMPLICIT NONE
REAL(rkind), INTENT(IN) :: VAR(MNP)
REAL(rkind), INTENT(OUT) :: DVDX(MNP), DVDY(MNP)
INTEGER :: NI(3)
INTEGER :: IE, IP
REAL(rkind) :: DEDY(3),DEDX(3)
REAL(rkind) :: DVDXIE, DVDYIE
REAL(rkind) :: WEI(MNP)
WEI(:) = 0.0_rkind
DVDX(:) = 0.0_rkind
DVDY(:) = 0.0_rkind
#ifdef DEBUG
WRITE(STAT%FHNDL,*) 'DIFFERENTIATE_XYDIR'
WRITE(STAT%FHNDL,*) 'sum(VAR ) = ', sum(VAR)
WRITE(STAT%FHNDL,*) 'sum(IEN ) = ', sum(IEN)
WRITE(STAT%FHNDL,*) 'sum(TRIA) = ', sum(TRIA)
#endif
DO IE = 1, MNE
NI = INE(:,IE)
WEI(NI) = WEI(NI) + 2.*TRIA(IE)
DEDX(1) = IEN(1,IE)
DEDX(2) = IEN(3,IE)
DEDX(3) = IEN(5,IE)
DEDY(1) = IEN(2,IE)
DEDY(2) = IEN(4,IE)
DEDY(3) = IEN(6,IE)
DVDXIE = DOT_PRODUCT( VAR(NI),DEDX)
DVDYIE = DOT_PRODUCT( VAR(NI),DEDY)
DVDX(NI) = DVDX(NI) + DVDXIE
DVDY(NI) = DVDY(NI) + DVDYIE
END DO
DVDX = DVDX/WEI
DVDY = DVDY/WEI
DO IP = 1, MNP
IF (DEP(IP) .LT. DMIN) THEN
DVDX(IP) = 0.
DVDY(IP) = 0.
END IF
END DO
#ifdef MPI_PARALL_GRID
CALL exchange_p2d(DVDX)
CALL exchange_p2d(DVDY)
#endif
#ifdef DEBUG
WRITE(STAT%FHNDL,*) 'sum(DVDX) = ', sum(DVDX)
WRITE(STAT%FHNDL,*) 'sum(DVDY) = ', sum(DVDY)
#endif
IF (.FALSE.) THEN
OPEN(2305, FILE = 'erggrad.bin' , FORM = 'UNFORMATTED')
WRITE(2305) 1.
WRITE(2305) (DVDX(IP), DVDY(IP), SQRT(DVDY(IP)**2+DVDY(IP)**2), IP = 1, MNP)
ENDIF
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE DIFFERENTIATE_XYDIR_NEIGHADJ(VAR, DVDX, DVDY)
USE DATAPOOL
IMPLICIT NONE
REAL(rkind), INTENT(IN) :: VAR(MNP)
REAL(rkind), INTENT(OUT) :: DVDX(MNP), DVDY(MNP)
INTEGER :: IP, IADJ, IP_ADJ
REAL(rkind) :: xpdiff, ypdiff, dist, eSum
REAL(rkind) :: XPcomp, YPcomp, eXP, eYP, eCoeff
DO IP=1,NP_RES
eSum=ZERO
XPcomp=ZERO
YPcomp=ZERO
DO IADJ=1,VERT_DEG(IP)
IP_ADJ=LIST_ADJ_VERT(IADJ,IP)
xpdiff=XP(IP_ADJ) - XP(IP)
ypdiff=YP(IP_ADJ) - YP(IP)
dist=sqrt(xpdiff*xpdiff+ypdiff*ypdiff)
eSum=eSum + ONE/dist
eCoeff=(VAR(IP_ADJ) - VAR(IP))/(dist**3)
XPcomp = XPcomp + xpdiff*eCoeff
YPcomp = YPcomp + ypdiff*eCoeff
END DO
eXP=XPcomp / eSum
eYP=YPcomp / eSum
DVDX(IP)=eXP
DVDY(IP)=eYP
END DO
#ifdef MPI_PARALL_GRID
CALL exchange_p2d(DVDX)
CALL exchange_p2d(DVDY)
#endif
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE DIFFERENTIATE_XYDIR_AND_SMOOTH(VAR, DVDX, DVDY, nbIter)
USE DATAPOOL
IMPLICIT NONE
integer, intent(in) :: nbIter
REAL(rkind), INTENT(IN) :: VAR(MNP)
REAL(rkind), INTENT(OUT) :: DVDX(MNP), DVDY(MNP)
INTEGER :: NI(3)
INTEGER :: IE, IP
REAL(rkind) :: DEDY(3),DEDX(3)
REAL(rkind) :: DVDXIE, DVDYIE
REAL(rkind) :: WEI(MNP)
REAL(rkind) :: IE_DY(MNE),IE_DX(MNE)
REAL(rkind) :: IE_DY_SMO(MNE),IE_DX_SMO(MNE)
WEI(:) = 0.0_rkind
DO IE = 1, MNE
NI = INE(:,IE)
WEI(NI) = WEI(NI) + 2.*TRIA(IE)
DEDX(1) = IEN(1,IE)
DEDX(2) = IEN(3,IE)
DEDX(3) = IEN(5,IE)
DEDY(1) = IEN(2,IE)
DEDY(2) = IEN(4,IE)
DEDY(3) = IEN(6,IE)
DVDXIE = DOT_PRODUCT( VAR(NI),DEDX)
DVDYIE = DOT_PRODUCT( VAR(NI),DEDY)
IE_DX(IE)=DVDXIE
IE_DY(IE)=DVDYIE
END DO
!
! Now smoothing
!
CALL SMOOTH_ON_TRIANGLE(IE_DX, IE_DX_SMO, nbIter)
CALL SMOOTH_ON_TRIANGLE(IE_DY, IE_DY_SMO, nbIter)
!
! Mapping to nodes
!
DVDX(:) = 0.0_rkind
DVDY(:) = 0.0_rkind
DO IE = 1, MNE
NI = INE(:,IE)
DVDX(NI) = DVDX(NI) + IE_DX_SMO(IE)
DVDY(NI) = DVDY(NI) + IE_DY_SMO(IE)
END DO
DVDX(:) = DVDX(:)/WEI(:)
DVDY(:) = DVDY(:)/WEI(:)
DO IP = 1, MNP
IF (DEP(IP) .LT. DMIN) THEN
DVDX(IP) = 0.
DVDY(IP) = 0.
END IF
END DO
#ifdef MPI_PARALL_GRID
CALL exchange_p2d(DVDX)
CALL exchange_p2d(DVDY)
#endif
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE WAVEKCG(DEPLOC, SIGIN, WN, WVC, WVK, WVCG2)
USE DATAPOOL
IMPLICIT NONE
REAL(rkind), INTENT(IN) :: DEPLOC, SIGIN
REAL(rkind), INTENT(OUT) :: WVC, WVK, WVCG2, WN
REAL(rkind) :: SGDLS , AUX1, AUX2
REAL(rkind) :: WKDEP
!
WVC = 0.
WVK = 0.
WVCG2 = 0.
WN = 0.
IF (SIGIN .LT. VERYSMALL) THEN
WN = 0.
WVK=10.
WVCG2=0.
RETURN
END IF
IF (DEPLOC .GT. DMIN) THEN
SGDLS = SIGIN*SIGIN*DEPLOC/G9
AUX1 = 1.0+0.6522*SGDLS+0.4622*(SGDLS**2.0)+0.0864*(SGDLS**4.0)+0.0675*(SGDLS**5.0)
AUX2 = 1.0/(SGDLS+1.0/AUX1)
WVC = SQRT(AUX2*G9*DEPLOC)
WVK = SIGIN/WVC
WKDEP = WVK*DEPLOC
IF (WKDEP > 13.0_rkind) THEN
WN = 0.5_rkind
ELSE
WN = 0.5_rkind*(ONE+TWO*WKDEP/SINH(MIN(TWO*KDMAX,TWO*WKDEP)))
END IF
WVCG2 = WN*WVC
ELSE
WVC = 0.
WVK = 10.
WVCG2 = 0.
END IF
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE MAKE_WAVE_TABLE()
USE DATAPOOL
IMPLICIT NONE
INTEGER :: IS
REAL(rkind) :: SIGIN, WVN, WVC, WVK, WVCG
REAL(rkind) :: DEPTH, SIGMAMAX
NMAX = IDISPTAB - 1
DEPTH = 1.
SIGMAMAX = SQRT (G9 * DEPFAC)
DSIGTAB = SIGMAMAX / MyREAL(NMAX)
TABK(0) = 0.
TABCG(0) = SQRT(G9)
DO IS = 1, NMAX
SIGIN = MyREAL(IS)*DSIGTAB
CALL WAVEKCG(DEPTH, SIGIN, WVN, WVC, WVK, WVCG)
TABK(IS) = WVK
TABCG(IS) = WVCG
END DO
IS = NMAX + 1
SIGIN = MyREAL(IS)*DSIGTAB
CALL WAVEKCG(DEPTH, SIGIN, WVN, WVC, WVK, WVCG)
TABK(IS) = WVK
TABCG(IS) = WVCG
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE ALL_FROM_TABLE(SIGIN,DEPTH,WVK,WVCG,WVKDEP,WVN,WVC)
USE DATAPOOL
IMPLICIT NONE
REAL(rkind), INTENT(IN) :: SIGIN, DEPTH
REAL(rkind), INTENT(OUT) :: WVK, WVCG, WVKDEP, WVN, WVC
REAL(rkind) :: SQRTDEP, SIGSQDEP, R1, R2, DEPLOC
INTEGER :: ITAB, ITAB2
DEPLOC = MAX(DMIN,DEPTH)
SQRTDEP = MySQRT(DEPLOC)
SIGSQDEP = SIGIN * SQRTDEP
ITAB = INT(SIGSQDEP/DSIGTAB)
!
IF (ITAB.LE.NMAX.AND.ITAB.GE.0) THEN
ITAB2 = ITAB + 1
R1 = SIGSQDEP/DSIGTAB - MyREAL(ITAB)
R2 = ONE - R1
WVK = ( R2*TABK(ITAB) + R1*TABK(ITAB2) ) / DEPLOC
WVCG = ( R2*TABCG(ITAB) + R1*TABCG(ITAB2) ) * SQRTDEP
WVKDEP = WVK * DEPLOC
WVN = ZEROFIVE*(ONE+TWO*WVKDEP/MySINH(MIN(KDMAX,TWO*WVKDEP)))
WVC = WVCG/WVN
ELSE
WVK = SIGIN*SIGIN/G9
WVCG = ZEROFIVE*G9/SIGIN
WVKDEP = WVK * DEPLOC
WVN = ZEROFIVE
WVC = TWO*WVCG
END IF
!
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE CG_FROM_TABLE(SIGIN,DEPTH,WVCG)
USE DATAPOOL
IMPLICIT NONE
REAL(rkind), INTENT(IN) :: SIGIN, DEPTH
REAL(rkind), INTENT(OUT) :: WVCG
REAL(rkind) :: SQRTDEP, SIGSQDEP, R1, DEPLOC
INTEGER :: ITAB
DEPLOC = MAX(DMIN,DEPTH)
SQRTDEP = SQRT(DEPLOC)
SIGSQDEP = SIGIN*SQRTDEP
ITAB = INT(SIGSQDEP/DSIGTAB)
IF (ITAB.LE.NMAX.AND.ITAB.GE.0) THEN
R1 = SIGSQDEP/DSIGTAB - MyREAL(ITAB)
WVCG = ((1.-R1)*TABCG(ITAB)+R1*TABCG(ITAB+1))*SQRTDEP
ELSE
WVCG = .5*G9/SIGIN
END IF
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE K_FROM_TABLE(SIGIN,DEPTH,WVK)
USE DATAPOOL
IMPLICIT NONE
REAL(rkind), INTENT(IN) :: SIGIN, DEPTH
REAL(rkind), INTENT(OUT) :: WVK
REAL(rkind) :: SQRTDEP, SIGSQDEP, R1, DEPLOC
INTEGER :: ITAB
DEPLOC = MAX(DMIN,DEPTH)
SQRTDEP = SQRT(DEPLOC)
SIGSQDEP = SIGIN * SQRTDEP
ITAB = INT(SIGSQDEP/DSIGTAB)
IF (ITAB.LE.NMAX.AND.ITAB.GE.0) THEN
R1 = SIGSQDEP/DSIGTAB - MyREAL(ITAB)
WVK = ((1.-R1)*TABK(ITAB)+R1*TABK(ITAB + 1))/DEPLOC
ELSE
WVK = SIGIN*SIGIN/G9
END IF
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE WAVE_K_C_CG
USE DATAPOOL
IMPLICIT NONE
INTEGER :: IP, IS
REAL(rkind) :: DEPLOC
REAL(rkind) :: WVK,WVCG,WVN,WVC,SPSIGLOC
DO IP = 1, MNP
DEPLOC = MAX(DMIN,DEP(IP))
! DEPLOC = DEP(IP)
! WRITE(740+myrank,*) 'IP=', IP, ' DEPLOC=', DEPLOC
DO IS = 1, NUMSIG
SPSIGLOC = SPSIG(IS)
! CALL ALL_FROM_TABLE(SPSIGLOC,DEPLOC,WVK,WVCG,WVKDEP,WVN,WVC)
CALL WAVEKCG(DEPLOC,SPSIGLOC,WVN,WVC,WVK,WVCG)
! WRITE(740+myrank,*) ' IS=', IS, 'WVK/WVCG/WVC=', WVK, WVCG, WVC
WK(IS,IP) = WVK
CG(IS,IP) = WVCG
WC(IP,IS) = WVC
END DO
END DO
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE CHECK_FOR_BLOW
USE DATAPOOL
IMPLICIT NONE
INTEGER IP, ID, IS
REAL(rkind) :: WALOC(NUMSIG,NUMDIR)
REAL(rkind) :: ETOT, DS, EAD, HS2
DO IP=1,MNP
WALOC(:,:) = AC2(:,:,IP)
ETOT = 0.0
DO ID = 1, NUMDIR
DO IS = 2, NUMSIG
DS = SPSIG(IS) - SPSIG(IS-1)
EAD = 0.5*(SPSIG(IS)*WALOC(IS,ID)+SPSIG(IS-1)*WALOC(IS-1,ID))*DS*DDIR
ETOT = ETOT + EAD
END DO
END DO
HS2 = 4.0_rkind*SQRT(ETOT)
IF (HS2 > LEVEL_HS_BLOW) THEN
CALL WWM_ABORT('Terminate due to wrong HS')
END IF
END DO
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE CHECK_STEADY(TIME, CONV1, CONV2, CONV3, CONV4, CONV5)
USE DATAPOOL
IMPLICIT NONE
!
REAL(rkind), INTENT(IN) :: TIME
REAL(rkind), INTENT(OUT) :: CONV1, CONV2, CONV3, CONV4, CONV5
INTEGER :: I, IP, IE, IS, ID, NI(3)
INTEGER :: IPCONV1, IPCONV2, IPCONV3, IPCONV4, IPCONV5, ISCONV(MNP)
REAL(rkind) :: SUMAC, WALOC(NUMSIG,NUMDIR)
REAL(rkind) :: ETOT, EAD, DS, HS2
REAL(rkind) :: ETOTF3, ETOTF4, TP, KHS2, EFTOT, TM02
REAL(rkind) :: FP, CP, KPP, CGP, WNP, UXD, OMEG, OMEG2
REAL(rkind) :: CONVK1, CONVK2, CONVK3, CONVK4, CONVK5
INTEGER :: ITMP
IPCONV1 = 0
IPCONV2 = 0
IPCONV3 = 0
IPCONV4 = 0
IPCONV5 = 0
#ifndef MPI_PARALL_GRID
DO IP = 1, MNP
#else
DO IP = 1, NP_RES
IF(ASSOCIATED(IPGL(IPLG(IP))%NEXT)) THEN !interface node
IF(IPGL(IPLG(ip))%NEXT%RANK < MYRANK) CYCLE !already in the sum so skip
ENDIF
#endif
WALOC(:,:) = AC2(:,:,IP)
SUMAC = SUM(WALOC)
ETOT = 0.0
DO ID = 1, NUMDIR
DO IS = 2, NUMSIG
DS = SPSIG(IS) - SPSIG(IS-1)
EAD = 0.5*(SPSIG(IS)*WALOC(IS,ID)+SPSIG(IS-1)*WALOC(IS-1,ID))*DS*DDIR
ETOT = ETOT + EAD
END DO
END DO
HS2 = 4.0_rkind*SQRT(ETOT)
ETOTF3 = 0.
ETOTF4 = 0.
DO IS = 1, NUMSIG
DO ID = 1, NUMDIR
ETOTF3 = ETOTF3 + SPSIG(IS) * WALOC(IS,ID)**4 * DDIR * DS_BAND(IS)
ETOTF4 = ETOTF4 + WALOC(IS,ID)**4 * DDIR * DS_BAND(IS)
END DO
END DO
IF(ETOTF4 .GT. THR8 .AND. ETOTF3 .GT. THR8) THEN
FP = ETOTF3/ETOTF4*PI2
CALL WAVEKCG(DEP(IP), FP, WNP, CP, KPP, CGP)
! CALL ALL_FROM_TABLE(FP,DEP(IP),KPP,CGP,KD,WNP,CP)
TP = 1.0_rkind/FP/PI2
KHS2 = HS2 * KPP
ELSE
KHS2 = 0.0_rkind
END IF
ETOT = 0.
EFTOT = 0.
DO ID=1, NUMDIR
IF (LSECU .OR. LSTCU) THEN
UXD = CURTXY(IP,1)*COSTH(ID) + CURTXY(IP,2)*SINTH(ID)
ENDIF
DO IS=1,NUMSIG
EAD = SIGPOW(IS,2) * WALOC(IS,ID) * FRINTF
IF (LSECU .OR. LSTCU) THEN
OMEG = SPSIG(IS) + WK(IS,IP) * UXD
OMEG2 = OMEG**2
ELSE
OMEG2 = SIGPOW(IS,2)
ENDIF
ETOT = ETOT + EAD
EFTOT = EFTOT + EAD * OMEG2
ENDDO
ENDDO
IF (ETOT/EFTOT .GT. THR) THEN
TM02 = PI2 * SQRT(ETOT/EFTOT)
ELSE
TM02 = 0.
END IF
IF (DEP(IP) .LT. DMIN .OR. SUMAC .LT. THR .OR. IOBP(IP) .EQ. 2 .OR. HS2 .LT. THR) THEN
IPCONV1 = IPCONV1 + 1 ! Summation of the converged grid points ...
IPCONV2 = IPCONV2 + 1
IPCONV3 = IPCONV3 + 1
IPCONV4 = IPCONV4 + 1
IPCONV5 = IPCONV5 + 1
ISCONV(IP) = 1
!write(dbg%fhndl,*) 'boundary -------1------', TIME, IP, IP_IS_STEADY(IP)
CYCLE
ELSE
CONVK1 = ABS(HSOLD(IP)-HS2)/HS2
CONVK2 = ABS(HS2-HSOLD(IP))
CONVK3 = ABS(SUMACOLD(IP)-SUMAC)/SUMAC
CONVK4 = ABS(KHS2-KHSOLD(IP))/KHSOLD(IP)
CONVK5 = ABS(TM02-TM02OLD(IP))/TM02OLD(IP)
IF (CONVK1 .LT. EPSH1) IPCONV1 = IPCONV1 + 1
IF (CONVK2 .LT. EPSH2) IPCONV2 = IPCONV2 + 1
IF (CONVK3 .LT. EPSH3) IPCONV3 = IPCONV3 + 1
IF (CONVK4 .LT. EPSH4) IPCONV4 = IPCONV4 + 1
IF (CONVK5 .LT. EPSH5) IPCONV5 = IPCONV5 + 1
IF (CONVK1 .LT. EPSH1 .AND. CONVK2 .LT. EPSH2 .AND. CONVK3 .LT. EPSH3 .AND. CONVK4 .LT. EPSH4 .AND. CONVK5 .LT. EPSH5) THEN
ISCONV(IP) = 1
!write(dbg%fhndl,*) 'converged -------2------', TIME, IP, IP_IS_STEADY(IP)
ELSE
ISCONV(IP) = 0
!write(dbg%fhndl,*) 'not converged -------3------', TIME, IP, IP_IS_STEADY(IP)
ENDIF
END IF
HSOLD(IP) = HS2
SUMACOLD(IP) = SUMAC
KHSOLD(IP) = KHS2
TM02OLD(IP) = TM02
END DO ! IP
!write(*,*) time, maxval(IP_IS_STEADY), minval(IP_IS_STEADY)
DO IE = 1, MNE
NI = INE(:,IE)
IF (SUM(ISCONV(NI)) .EQ. 3) THEN
IE_IS_STEADY(IE) = IE_IS_STEADY(IE) + 1
ELSE
IE_IS_STEADY(IE) = 0
ENDIF
!WRITE(*,*) IE, IE_IS_STEADY(IE)
ENDDO
DO IP = 1, MNP
ITMP = 0
DO I = 1, CCON(IP)
IF (IE_IS_STEADY(IE_CELL2(IP,I)) .GE. 1) ITMP = ITMP + 1
ENDDO
IF (ITMP .EQ. CCON(IP)) THEN
IP_IS_STEADY(IP) = IP_IS_STEADY(IP) + 1
!IF (IP_IS_STEADY(IP) .GT. 2) WRITE(*,*) TIME, IP, IP_IS_STEADY(IP)
ELSE
IP_IS_STEADY(IP) = 0
ENDIF
ENDDO
#ifdef MPI_PARALL_GRID
CALL MPI_ALLREDUCE(IPCONV1, itmp, 1, itype, MPI_SUM, COMM, ierr)
IPCONV1 = itmp
CALL MPI_ALLREDUCE(IPCONV2, itmp, 1, itype, MPI_SUM, COMM, ierr)
IPCONV2 = itmp
CALL MPI_ALLREDUCE(IPCONV3, itmp, 1, itype, MPI_SUM, COMM, ierr)
IPCONV3 = itmp
CALL MPI_ALLREDUCE(IPCONV4, itmp, 1, itype, MPI_SUM, COMM, ierr)
IPCONV4 = itmp
CALL MPI_ALLREDUCE(IPCONV5, itmp, 1, itype, MPI_SUM, COMM, ierr)
IPCONV5 = itmp
CONV1 = MyREAL(IPCONV1)/MyREAL(NP_GLOBAL)*100.0
CONV2 = MyREAL(IPCONV2)/MyREAL(NP_GLOBAL)*100.0
CONV3 = MyREAL(IPCONV3)/MyREAL(NP_GLOBAL)*100.0
CONV4 = MyREAL(IPCONV4)/MyREAL(NP_GLOBAL)*100.0
CONV5 = MyREAL(IPCONV5)/MyREAL(NP_GLOBAL)*100.0
#else
CONV1 = MyREAL(IPCONV1)/MyREAL(MNP)*100.0
CONV2 = MyREAL(IPCONV2)/MyREAL(MNP)*100.0
CONV3 = MyREAL(IPCONV3)/MyREAL(MNP)*100.0
CONV4 = MyREAL(IPCONV4)/MyREAL(MNP)*100.0
CONV5 = MyREAL(IPCONV5)/MyREAL(MNP)*100.0
#endif
#ifdef MPI_PARALL_GRID
IF (myrank == 0) THEN
WRITE(STAT%FHNDL,*) 'CONVERGENCE CRIT. 1 REACHED IN', CONV1, '% GRIDPOINTS'
WRITE(STAT%FHNDL,*) 'CONVERGENCE CRIT. 2 REACHED IN', CONV2, '% GRIDPOINTS'
WRITE(STAT%FHNDL,*) 'CONVERGENCE CRIT. 3 REACHED IN', CONV3, '% GRIDPOINTS'
WRITE(STAT%FHNDL,*) 'CONVERGENCE CRIT. 4 REACHED IN', CONV4, '% GRIDPOINTS'
WRITE(STAT%FHNDL,*) 'CONVERGENCE CRIT. 5 REACHED IN', CONV5, '% GRIDPOINTS'
END IF
#else
WRITE(STAT%FHNDL,*) 'CONVERGENCE CRIT. 1 REACHED IN', CONV1, '% GRIDPOINTS'
WRITE(STAT%FHNDL,*) 'CONVERGENCE CRIT. 2 REACHED IN', CONV2, '% GRIDPOINTS'
WRITE(STAT%FHNDL,*) 'CONVERGENCE CRIT. 3 REACHED IN', CONV3, '% GRIDPOINTS'
WRITE(STAT%FHNDL,*) 'CONVERGENCE CRIT. 4 REACHED IN', CONV4, '% GRIDPOINTS'
WRITE(STAT%FHNDL,*) 'CONVERGENCE CRIT. 5 REACHED IN', CONV5, '% GRIDPOINTS'
#endif
FLUSH(STAT%FHNDL)
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE TWOD2ONED(WALOC,AC1D)
USE DATAPOOL
IMPLICIT NONE
REAL(rkind), INTENT(IN) :: WALOC(NUMSIG,NUMDIR)
REAL(rkind), INTENT(OUT) :: AC1D(NUMSIG*NUMDIR)
INTEGER :: IS, ID
DO IS = 1, NUMSIG
DO ID = 1, NUMDIR
AC1D(ID + (IS-1) * NUMDIR) = WALOC(IS,ID)
END DO
END DO
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE CONVERT_VS_VD_WWM(IP, VS, VD, PHI_R, DPHIDN_R)
USE DATAPOOL
IMPLICIT NONE
INTEGER, intent(in) :: IP
REAL(rkind), intent(in) :: VS(NSPEC), VD(NSPEC)
REAL(rkind), intent(out) :: PHI_R(NUMSIG,NUMDIR), DPHIDN_R(NUMSIG,NUMDIR)
INTEGER :: IS, ID
DO IS=1,NUMSIG
DO ID=1,NUMDIR
PHI_R(IS,ID) = VS(ID + (IS-1) * NUMDIR) / CG(IS,IP)
DPHIDN_R(IS,ID) = VD(ID + (IS-1) * NUMDIR)
END DO
END DO
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE CONVERT_WWM_VS_VD(IP, VS, VD, PHI_R, DPHIDN_R)
USE DATAPOOL
IMPLICIT NONE
INTEGER, intent(in) :: IP
REAL(rkind), intent(out) :: VS(NSPEC), VD(NSPEC)
REAL(rkind), intent(in) :: PHI_R(NUMSIG,NUMDIR), DPHIDN_R(NUMSIG,NUMDIR)
INTEGER :: IS, ID
DO IS=1,NUMSIG
DO ID=1,NUMDIR
VS(ID + (IS-1) * NUMDIR) = PHI_R(IS,ID) * CG(IS,IP)
VD(ID + (IS-1) * NUMDIR) = DPHIDN_R(IS,ID)
END DO
END DO
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE ONED2TWOD(AC1D,WALOC)
USE DATAPOOL
IMPLICIT NONE
REAL(rkind), INTENT(OUT) :: WALOC(NUMSIG,NUMDIR)
REAL(rkind), INTENT(IN) :: AC1D(NUMSIG*NUMDIR)
INTEGER :: IS, ID
DO IS = 1, NUMSIG
DO ID = 1, NUMDIR
WALOC(IS,ID) = AC1D(ID + (IS-1) * NUMDIR)
END DO
END DO
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
! REAL(rkind) FUNCTION GAMMA_FUNC(XX)
FUNCTION GAMMA_FUNC(XX)
USE DATAPOOL, ONLY: rkind
IMPLICIT NONE
REAL(rkind) :: GAMMA_FUNC
!
! Purpose:
! Compute the transcendental function Gamma
!
! Subroutines used
! GAMMLN (Numerical Recipes)
!
real(rkind) GAMMLN
REAL(rkind) XX, YY, ABIG
SAVE ABIG
DATA ABIG /30./
YY = GAMMLN(XX)
IF (YY > ABIG) YY = ABIG
IF (YY < -ABIG) YY = -ABIG
GAMMA_FUNC = EXP(YY)
END FUNCTION
!**********************************************************************
!* *
!**********************************************************************
FUNCTION GAMMLN(XX)
USE DATAPOOL, ONLY: rkind, ONEHALF, ONE
IMPLICIT NONE
real(rkind) :: GAMMLN
!
! Method:
! function is copied from: Press et al., "Numerical Recipes"
!
real(rkind) XX
INTEGER J
real(rkind) COF(6),STP,FPF,X,TMP,SER
DATA COF,STP/76.18009173_rkind,-86.50532033_rkind, &
& 24.01409822_rkind,-1.231739516_rkind, &
& .120858003E-2_rkind,-.536382E-5_rkind, &
& 2.50662827465_rkind/
DATA FPF/5.5_rkind/
X = XX-ONE
TMP = X+FPF
TMP = (X+ONEHALF)*LOG(TMP)-TMP
SER = ONE
DO J = 1, 6
X = X+ONE
SER = SER+COF(J)/X
END DO
GAMMLN = TMP+LOG(STP*SER)
END FUNCTION
!**********************************************************************
!* *
!**********************************************************************
FUNCTION VEC2RAD(U,V)
USE DATAPOOL, ONLY: rkind, pi
IMPLICIT NONE
REAL(RKIND) :: VEC2RAD
REAL(rkind) :: U,V
VEC2RAD = MyATAN2(V,U) * 180/PI
IF (VEC2RAD < 0.0) VEC2RAD = VEC2RAD + 360.0_rkind
VEC2RAD = VEC2RAD * PI/180.
END FUNCTION
!**********************************************************************
!* *
!**********************************************************************
FUNCTION VEC2DEG(U,V)
USE DATAPOOL, ONLY: rkind, pi
IMPLICIT NONE
REAL(RKIND) :: VEC2DEG
REAL(rkind), intent(in) :: U,V
VEC2DEG = MyATAN2(V,U) * 180./PI
IF (VEC2DEG < 0.0) VEC2DEG = VEC2DEG + 360.0_rkind
END FUNCTION
!**********************************************************************
!* *
!**********************************************************************
FUNCTION DVEC2RAD(U,V)
USE DATAPOOL, ONLY: rkind, pi