-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathwwm_specparam.F90
1758 lines (1581 loc) · 60.2 KB
/
wwm_specparam.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"
! Last change: 1 9 Jun 2004 1:44 am
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE STOKES_DRIFT_SURFACE_BAROTROPIC(IP,STOKESBOTTX,STOKESBOTTY,STOKESSURFX,STOKESSURFY,STOKESBAROX,STOKESBAROY)
USE DATAPOOL
IMPLICIT NONE
INTEGER, INTENT(IN) :: IP
REAL(rkind), INTENT(OUT) :: STOKESBOTTX, STOKESBOTTY
REAL(rkind), INTENT(OUT) :: STOKESSURFX, STOKESSURFY
REAL(rkind), INTENT(OUT) :: STOKESBAROX, STOKESBAROY
INTEGER :: ID, IS
REAL(rkind) :: eQuot1, eQuot2, eProd1, eProd2
REAL(rkind) :: eMult, eWk, kD, eWkReal
REAL(rkind) :: eSinh2kd, eSinhkd, eSinhkd2
REAL(rkind) :: eSigma, eLoc, eUint, eVint
REAL(rkind) :: eDep
REAL(rkind) :: eProd0, eQuot0
STOKESBOTTX=0
STOKESBOTTY=0
STOKESSURFX=0
STOKESSURFY=0
STOKESBAROX=0
STOKESBAROY=0
eDep=DEP(IP)
DO IS=1,NUMSIG
eMult=SPSIG(IS)*DDIR*DS_INCR(IS)
eWk=WK(IS,IP)
kD=MIN(KDMAX, eWk*eDep)
eWkReal=kD/eDep
eSinh2kd=MySINH(2*kD)
eSinhkd=MySINH(kD)
eSinhkd2=eSinhkd**2
eSigma=SPSIG(IS)
eUint=0
eVint=0
DO ID=1,NUMDIR
eLoc=AC2(IS,ID,IP)*eMult
eUint=eUint + eLoc*COSTH(ID)
eVint=eVint + eLoc*SINTH(ID)
END DO
eQuot0=ONE/eSinhkd2
eProd0=eSigma*eWkReal*eQuot0
eQuot1=MyCOSH(2*kD)/eSinhkd2
eProd1=eSigma*eWkReal*eQuot1
eQuot2=(eSinh2kd/(2*kD))/eSinhkd2
eProd2=eSigma*eWkReal*eQuot2
STOKESBOTTX = STOKESBOTTX + eUint*eProd0
STOKESBOTTY = STOKESBOTTY + eVint*eProd0
STOKESSURFX = STOKESSURFX + eUint*eProd1
STOKESSURFY = STOKESSURFY + eVint*eProd1
STOKESBAROX = STOKESBAROX + eUint*eProd2
STOKESBAROY = STOKESBAROY + eVint*eProd2
END DO
END SUBROUTINE
!**********************************************************************:
!* *
!**********************************************************************
SUBROUTINE STOKES_DRIFT_SURFACE_BAROTROPIC_LOC(WALOC,DEPLOC,WKLOC, &
& STOKESBOTTX,STOKESBOTTY,STOKESSURFX,STOKESSURFY,STOKESBAROX,STOKESBAROY)
USE DATAPOOL
IMPLICIT NONE
REAL(rkind), INTENT(IN) :: WALOC(NUMSIG,NUMDIR)
REAL(rkind), INTENT(IN) :: DEPLOC
REAL(rkind), INTENT(IN) :: WKLOC(NUMSIG)
REAL(rkind), INTENT(OUT) :: STOKESBOTTX, STOKESBOTTY
REAL(rkind), INTENT(OUT) :: STOKESSURFX, STOKESSURFY
REAL(rkind), INTENT(OUT) :: STOKESBAROX, STOKESBAROY
INTEGER :: ID, IS
REAL(rkind) :: eQuot1, eQuot2, eProd1, eProd2
REAL(rkind) :: eMult, eWk, kD, eWkReal
REAL(rkind) :: eSinh2kd, eSinhkd, eSinhkd2
REAL(rkind) :: eSigma, eLoc, eUint, eVint
REAL(rkind) :: eDep
REAL(rkind) :: eProd0, eQuot0
STOKESBOTTX=0
STOKESBOTTY=0
STOKESSURFX=0
STOKESSURFY=0
STOKESBAROX=0
STOKESBAROY=0
eDep=DEPLOC
DO IS=1,NUMSIG
eMult=SPSIG(IS)*DDIR*DS_INCR(IS)
eWk=WKLOC(IS)
kD=MIN(KDMAX, eWk*eDep)
eWkReal=kD/eDep
eSinh2kd=MySINH(2*kD)
eSinhkd=MySINH(kD)
eSinhkd2=eSinhkd**2
eSigma=SPSIG(IS)
eUint=0
eVint=0
DO ID=1,NUMDIR
eLoc=WALOC(IS,ID)*eMult
eUint=eUint + eLoc*COSTH(ID)
eVint=eVint + eLoc*SINTH(ID)
END DO
eQuot0=ONE/eSinhkd2
eProd0=eSigma*eWkReal*eQuot0
eQuot1=MyCOSH(2*kD)/eSinhkd2
eProd1=eSigma*eWkReal*eQuot1
eQuot2=(eSinh2kd/(2*kD))/eSinhkd2
eProd2=eSigma*eWkReal*eQuot2
STOKESBOTTX = STOKESBOTTX + eUint*eProd0
STOKESBOTTY = STOKESBOTTY + eVint*eProd0
STOKESSURFX = STOKESSURFX + eUint*eProd1
STOKESSURFY = STOKESSURFY + eVint*eProd1
STOKESBAROX = STOKESBAROX + eUint*eProd2
STOKESBAROY = STOKESBAROY + eVint*eProd2
END DO
END SUBROUTINE
!**********************************************************************:
!* *
!**********************************************************************
SUBROUTINE MEAN_WAVE_PARAMETER(IP,WALOC,HS,ETOT,SME01,SME10,KME01,KMWAM,KMWAM2)
!AR: This must be replaced by ST4_PRE or the certain WAM routine we are not consistent here this routine needs urgen deletion
USE DATAPOOL
IMPLICIT NONE
INTEGER, INTENT(IN) :: IP
REAL(rkind), INTENT(IN) :: WALOC(NUMSIG,NUMDIR)
REAL(rkind), INTENT(OUT) :: SME01, SME10, KME01
REAL(rkind), INTENT(OUT) :: KMWAM, KMWAM2
REAL(rkind), INTENT(OUT) :: HS
INTEGER :: ID, IS
REAL(rkind) :: ACTOT, ETOT
REAL(rkind) :: ETOT_SPSIG
REAL(rkind) :: ETOT_WK
REAL(rkind) :: ETOT_ISQ_WK
REAL(rkind) :: ETOT_SQ_WK
REAL(rkind) :: Y(NUMSIG), tmp(NUMSIG)
REAL(rkind) :: DS, ATAIL, ETAIL, ESIGTAIL
! REAL(rkind) :: dintspec, dintspec_y
!
! total energy ...
! 2do improve efficiency ...
! 2do check integration style ...
!
!ETOT = DINTSPEC(IP,WALOC)
ETOT = ZERO
do id = 1, NUMDIR
tmp(:) = WALOC(:,id) * spsig
ETOT = ETOT + tmp(1) * ONEHALF * ds_incr(1)*ddir
do is = 2, NUMSIG
ETOT = ETOT + ONEHALF*(tmp(is)+tmp(is-1))*ds_band(is)*ddir
end do
ETOT = ETOT + ONEHALF * tmp(NUMSIG) * ds_incr(NUMSIG)*ddir
end do
!
! if etot too small skip ...
!
if (etot .gt. thr) then
!
! integrals ... inlined ... for speed ...
!
ACTOT = ZERO
y = ONE/SPSIG
do id = 1, NUMDIR
tmp(:) = WALOC(:,id) * spsig * y
ACTOT = ACTOT + tmp(1) * ONEHALF * ds_incr(1)*ddir
do is = 2, NUMSIG
ACTOT = ACTOT + ONEHALF*(tmp(is)+tmp(is-1))*ds_band(is)*ddir
end do
ACTOT = ACTOT + ONEHALF * tmp(NUMSIG) * ds_incr(NUMSIG)*ddir
end do
!tmp = ONE/SPSIG
!ACTOT = DINTSPEC_Y(IP,WALOC,tmp)
ETOT_SPSIG = ZERO
y = SIGPOW(:,1)
do id = 1, NUMDIR
tmp(:) = WALOC(:,id) * spsig * y
ETOT_SPSIG = ETOT_SPSIG + tmp(1) * ONEHALF * ds_incr(1)*ddir
do is = 2, NUMSIG
ETOT_SPSIG = ETOT_SPSIG + ONEHALF*(tmp(is)+tmp(is-1))*ds_band(is)*ddir
end do
ETOT_SPSIG = ETOT_SPSIG + ONEHALF * tmp(NUMSIG) * ds_incr(NUMSIG)*ddir
end do
!tmp = SIGPOW(:,1)
!ETOT_SPSIG = DINTSPEC_Y(IP,WALOC,tmp)
ETOT_WK = ZERO
y = WK(:,IP)
do id = 1, NUMDIR
tmp(:) = WALOC(:,id) * spsig * y
ETOT_WK = ETOT_WK + ONEHALF * tmp(1) * ds_incr(1)*ddir
do is = 2, NUMSIG
ETOT_WK = ETOT_WK + ONEHALF*(tmp(is)+tmp(is-1))*ds_band(is)*ddir
end do
ETOT_WK = ETOT_WK + ONEHALF * tmp(NUMSIG) * ds_incr(NUMSIG)*ddir
end do
!tmp = WK(:,IP)
!ETOT_WK = DINTSPEC_Y(IP,WALOC,tmp)
ETOT_ISQ_WK = ZERO
y = ONE/SQRT(WK(:,IP))
do id = 1, NUMDIR
tmp(:) = WALOC(:,id) * spsig * y
ETOT_ISQ_WK = ETOT_ISQ_WK + ONEHALF * tmp(1) * ds_incr(1)*ddir
do is = 2, NUMSIG
ETOT_ISQ_WK = ETOT_ISQ_WK+ONEHALF*(tmp(is)+tmp(is-1))*ds_band(is)*ddir
end do
ETOT_ISQ_WK = ETOT_ISQ_WK+ONEHALF*tmp(NUMSIG) * ds_incr(NUMSIG)*ddir
end do
!tmp = ONE/SQRT(WK(:,IP))
!ETOT_ISQ_WK = DINTSPEC_Y(IP,WALOC,tmp)
ETOT_SQ_WK = ZERO
y = SQRT(WK(:,IP))
do id = 1, NUMDIR
tmp(:) = WALOC(:,id) * spsig * y
ETOT_SQ_WK = ETOT_SQ_WK + ONEHALF * tmp(1) * ds_incr(1)*ddir
do is = 2, NUMSIG -1
ETOT_SQ_WK = ETOT_SQ_WK + ONEHALF*(tmp(is)+tmp(is-1))*ds_incr(is)*ddir
end do
ETOT_SQ_WK = ETOT_SQ_WK + ONEHALF*tmp(NUMSIG) * ds_incr(NUMSIG)*ddir
end do
!tmp = SQRT(WK(:,IP))
!ETOT_SQ_WK = DINTSPEC_Y(IP,WALOC,tmp)
!
! tail factors ...
!
DS = SPSIG(NUMSIG) - SPSIG(NUMSIG-1)
ATAIL = SUM(WALOC(NUMSIG,:)) * SIGPOW(NUMSIG,1) * DDIR * DS
ETAIL = SUM(WALOC(NUMSIG,:)) * SIGPOW(NUMSIG,2) * DDIR * DS
ESIGTAIL = SUM(WALOC(NUMSIG,:)) * SIGPOW(NUMSIG,3) * DDIR * DS
!
! tail factors ...
!
ACTOT = ACTOT + TAIL_ARR(5) * ATAIL
ETOT = ETOT + TAIL_ARR(6) * ETAIL
ETOT_SPSIG = ETOT_SPSIG + TAIL_ARR(7) * ETAIL
ETOT_ISQ_WK = ETOT_ISQ_WK + TAIL_ARR(5) * ETAIL / (SQRT(WK(NUMSIG,IP)))
ETOT_SQ_WK = ETOT_SQ_WK + TAIL_ARR(5) * ETAIL * (SQRT(WK(NUMSIG,IP)))
ETOT_WK = ETOT_WK + TAIL_ARR(8) * ETAIL * WK(NUMSIG,IP)
!
! integral parameters ...
!
HS = MAX(ZERO,4.*SQRT(ETOT))
IF (LMONO_OUT) HS = HS / SQRT(2.)
SME01 = ETOT_SPSIG / ETOT
KME01 = ETOT_WK / ETOT
SME10 = ETOT / ACTOT
KMWAM = (ETOT/ETOT_ISQ_WK)**2
KMWAM2 = (ETOT_SQ_WK/ETOT)**2
else
!
! no or too less energy ...
!
ETOT = ZERO
HS = ZERO
SME01 = ZERO
KME01 = 10.0_rkind
SME10 = ZERO
KMWAM = 10.0_rkind
KMWAM2 = 10.0_rkind
end if
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
!AR: This must be replaced by ST4_PRE or the certain WAM routine we are not consistent here this routine needs urgen deletion
SUBROUTINE MEAN_PARAMETER_BDCONS(WALOC,HS,TM01,TM02)
USE DATAPOOL
IMPLICIT NONE
REAL(rkind), INTENT(IN) :: WALOC(NUMSIG,NUMDIR)
REAL(rkind), INTENT(OUT) :: HS,TM01,TM02
INTEGER :: ID, IS
REAL(rkind) :: Y(NUMSIG)
REAL(rkind) :: DS, ETAIL
REAL(rkind) :: OMEG2, EAD, ETOT
REAL(rkind) :: EFTAIL,PTAIL_ARR,EFTOT,ETAIL_ARR
REAL(rkind) :: EHFR,AHFR,ATAIL_ARR,EPTOT,APTOT
REAL(rkind) :: tmp(NUMSIG),actmp(NUMSIG)
!
! total energy ...
!
Y = ZERO
tmp = ZERO
actmp = ZERO
ETOT = ZERO
do id = 1, NUMDIR
tmp(:) = WALOC(:,id) * spsig
ETOT = ETOT + tmp(1) * ONEHALF * ds_incr(1)*ddir
do is = 2, NUMSIG
ETOT = ETOT + ONEHALF*(tmp(is)+tmp(is-1))*ds_band(is)*ddir
end do
ETOT = ETOT + ONEHALF * tmp(NUMSIG) * ds_incr(NUMSIG)*ddir
end do
IF (ETOT .GT. THR) THEN
!
! tail ratios
!
DS = SPSIG(NUMSIG) - SPSIG(NUMSIG-1)
ETAIL = SUM(WALOC(NUMSIG,:)) * SIGPOW(NUMSIG,2) * DDIR * DS
ETOT = ETOT + TAIL_ARR(6) * ETAIL
HS = 4*SQRT(ETOT)
APTOT = ZERO
EPTOT = ZERO
PTAIL_ARR = TAIL_ARR(1)
ATAIL_ARR = ONE / (PTAIL_ARR * (ONE + PTAIL_ARR * (FRINTH-ONE)))
PTAIL_ARR = TAIL_ARR(1) - ONE
ETAIL_ARR = ONE / (PTAIL_ARR * (ONE + PTAIL_ARR * (FRINTH-ONE)))
DO ID = 1, NUMDIR
DO IS = 1, NUMSIG
APTOT = APTOT + SPSIG(IS) * WALOC(IS,ID)
EPTOT = EPTOT + SIGPOW(IS,2) * WALOC(IS,ID)
ENDDO
ENDDO
APTOT = APTOT * FRINTF
EPTOT = EPTOT * FRINTF
IF (NUMSIG .GT. 3 .AND. .NOT. LSIGMAX) THEN
DO ID = 1, NUMDIR
AHFR = SPSIG(NUMSIG) * WALOC(NUMSIG,ID)
APTOT = APTOT + ATAIL_ARR * AHFR
EHFR = SPSIG(NUMSIG) * AHFR
EPTOT = EPTOT + ETAIL_ARR * EHFR
ENDDO
ENDIF
IF (EPTOT .GT. ZERO) THEN
TM01 = PI2 * APTOT / EPTOT
ELSE
TM01 = ZERO
END IF
ETOT = ZERO
EFTOT = ZERO
PTAIL_ARR = TAIL_ARR(1) - ONE
ETAIL = ONE / (PTAIL_ARR * (ONE + PTAIL_ARR * (FRINTH-ONE)))
PTAIL_ARR = TAIL_ARR(1) - 3.
EFTAIL = ONE / (PTAIL_ARR * (ONE + PTAIL_ARR * (FRINTH-ONE)))
DO ID=1, NUMDIR
DO IS = 1, NUMSIG
EAD = SIGPOW(IS,2) * WALOC(IS,ID) * FRINTF
OMEG2 = SIGPOW(IS,2)
ETOT = ETOT + EAD
EFTOT = EFTOT + EAD * OMEG2
ENDDO
IF (NUMSIG .GT. 3 .AND. .NOT. LSIGMAX) THEN
EAD = SIGPOW(NUMSIG,2) * WALOC(NUMSIG,ID)
ETOT = ETOT + ETAIL * EAD
EFTOT = EFTOT + EFTAIL * OMEG2 * EAD
ENDIF
ENDDO
IF (EFTOT .GT. sqrt(verysmall)) THEN
TM02 = PI2 * SQRT(ETOT/EFTOT)
ELSE
TM02 = ZERO
END IF
ELSE
HS = ZERO
TM01 = ZERO
TM02 = ZERO
END IF
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
!AR: This must be replaced by ST4_PRE or the certain WAM routine we are not consistent here this routine needs urgen deletion
SUBROUTINE MEAN_PARAMETER(IP,WALOC,ISMAX,HS,TM01,TM02,TM10,KLM,WLM)
USE DATAPOOL
IMPLICIT NONE
!2do ... rewrite this integration ...
INTEGER, INTENT(IN) :: IP,ISMAX
REAL(rkind), INTENT(IN) :: WALOC(NUMSIG,NUMDIR)
REAL(rkind), INTENT(OUT) :: HS,TM01,TM02,KLM,WLM,TM10
INTEGER :: ID, IS
REAL(rkind) :: Y(NUMSIG)
REAL(rkind) :: DS, ETAIL
REAL(rkind) :: OMEG2,OMEG,EAD,UXD, ETOT
REAL(rkind) :: EFTAIL,PTAIL_ARR,EFTOT,ETAIL_ARR
REAL(rkind) :: EHFR,AHFR,ATAIL_ARR,EPTOT,APTOT
REAL(rkind) :: CKTAIL, ETOT1, SIG22, EKTOT, CETAIL
REAL(rkind) :: tmp(NUMSIG),actmp(NUMSIG), SKK
!
! total energy ...
!
Y = ZERO
tmp = ZERO
actmp = ZERO
ETOT = ZERO
do id = 1, NUMDIR
tmp(:) = WALOC(:,id) * spsig
ETOT = ETOT + tmp(1) * ONEHALF * ds_incr(1)*ddir
do is = 2, NUMSIG
ETOT = ETOT + ONEHALF*(tmp(is)+tmp(is-1))*ds_band(is)*ddir
end do
ETOT = ETOT + ONEHALF * tmp(NUMSIG) * ds_incr(NUMSIG)*ddir
end do
IF (ETOT .GT. THR) THEN
!
! tail ratios
!
DS = SPSIG(NUMSIG) - SPSIG(NUMSIG-1)
ETAIL = SUM(WALOC(NUMSIG,:)) * SIGPOW(NUMSIG,2) * DDIR * DS
ETOT = ETOT + TAIL_ARR(6) * ETAIL
HS = 4*SQRT(ETOT)
APTOT = ZERO
EPTOT = ZERO
PTAIL_ARR = TAIL_ARR(1)
ATAIL_ARR = ONE / (PTAIL_ARR * (ONE + PTAIL_ARR * (FRINTH-ONE)))
PTAIL_ARR = TAIL_ARR(1) - ONE
ETAIL_ARR = ONE / (PTAIL_ARR * (ONE + PTAIL_ARR * (FRINTH-ONE)))
DO ID = 1, NUMDIR
DO IS = 1, ISMAX
APTOT = APTOT + SPSIG(IS) * WALOC(IS,ID)
EPTOT = EPTOT + SIGPOW(IS,2) * WALOC(IS,ID)
ENDDO
ENDDO
APTOT = APTOT * FRINTF
EPTOT = EPTOT * FRINTF
IF (NUMSIG .GT. 3 .AND. .NOT. LSIGMAX) THEN
DO ID = 1, NUMDIR
AHFR = SPSIG(NUMSIG) * WALOC(NUMSIG,ID)
APTOT = APTOT + ATAIL_ARR * AHFR
EHFR = SPSIG(NUMSIG) * AHFR
EPTOT = EPTOT + ETAIL_ARR * EHFR
ENDDO
ENDIF
IF (EPTOT .GT. ZERO) THEN
TM01 = PI2 * APTOT / EPTOT
ELSE
TM01 = ZERO
END IF
ETOT = ZERO
EFTOT = ZERO
PTAIL_ARR = TAIL_ARR(1) - ONE
ETAIL = ONE / (PTAIL_ARR * (ONE + PTAIL_ARR * (FRINTH-ONE)))
PTAIL_ARR = TAIL_ARR(1) - 3.
EFTAIL = ONE / (PTAIL_ARR * (ONE + PTAIL_ARR * (FRINTH-ONE)))
DO ID=1, NUMDIR
IF (LSECU .OR. LSTCU) THEN
UXD = CURTXY(IP,1)*COSTH(ID) + CURTXY(IP,2)*SINTH(ID)
ENDIF
DO IS = 1, ISMAX
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
IF (NUMSIG .GT. 3 .AND. .NOT. LSIGMAX) THEN
EAD = SIGPOW(NUMSIG,2) * WALOC(NUMSIG,ID)
ETOT = ETOT + ETAIL * EAD
EFTOT = EFTOT + EFTAIL * OMEG2 * EAD
ENDIF
ENDDO
IF (EFTOT .GT. sqrt(verysmall)) THEN
TM02 = PI2 * SQRT(ETOT/EFTOT)
ELSE
TM02 = ZERO
END IF
ETOT1 = ZERO
EKTOT = ZERO
!
! tail ratios same
!
PTAIL_ARR = TAIL_ARR(1) - ONE
CETAIL = ONE / (PTAIL_ARR * (ONE + PTAIL_ARR * (FRINTH-ONE)))
PTAIL_ARR = TAIL_ARR(1) - ONE - 2*ONE
CKTAIL = ONE / (PTAIL_ARR * (ONE + PTAIL_ARR * (FRINTH-ONE)))
DO IS = 1, ISMAX
SIG22 = SIGPOW(IS,2)
SKK = SIG22 * WK(IS,IP)
DO ID = 1, NUMDIR
ETOT1 = ETOT1 + SIG22 * WALOC(IS,ID)
EKTOT = EKTOT + SKK * WALOC(IS,ID)
ENDDO
ENDDO
ETOT1 = FRINTF * ETOT1
EKTOT = FRINTF * EKTOT
IF (NUMSIG .GT. 3) THEN
DO ID=1,NUMDIR
ETOT1 = ETOT1 + CETAIL * SIG22 * WALOC(NUMSIG,ID)
EKTOT = EKTOT + CKTAIL * SKK * WALOC(NUMSIG,ID)
ENDDO
ENDIF
IF (ETOT1.GT.VERYSMALL.AND.EKTOT.GT.VERYSMALL) THEN
WLM = PI2 * (ETOT1/EKTOT)
KLM = PI2/WLM
ELSE
KLM = 10.0_rkind
WLM = ZERO
ENDIF
APTOT = 0.
EPTOT = 0.
DO ID=1, NUMDIR
DO IS=1,ISMAX
APTOT = APTOT + SPSIG(IS) * WALOC(IS,ID)
EPTOT = EPTOT + SIGPOW(IS,2) * WALOC(IS,ID)
ENDDO
ENDDO
APTOT = APTOT * FRINTF
EPTOT = EPTOT * FRINTF
IF (NUMSIG .GT. 3) THEN
PTAIL_ARR = TAIL_ARR(1)
ATAIL_ARR = 1. / (PTAIL_ARR * (1. + PTAIL_ARR * (FRINTH-1.)))
PTAIL_ARR = TAIL_ARR(1) - 1.
ETAIL_ARR = 1. / (PTAIL_ARR * (1. + PTAIL_ARR * (FRINTH-1.)))
DO ID = 1, NUMDIR
AHFR = SPSIG(NUMSIG) * WALOC(NUMSIG,ID)
APTOT = APTOT + ATAIL_ARR * AHFR
EHFR = SPSIG(NUMSIG) * AHFR
EPTOT = EPTOT + ETAIL_ARR * EHFR
ENDDO
ENDIF
TM10 = 2.*PI * APTOT / EPTOT
ELSE
HS = ZERO
TM01 = ZERO
TM02 = ZERO
TM10 = ZERO
KLM = 10.0_rkind
WLM = ZERO
END IF
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
SUBROUTINE MEAN_PARAMETER_OUTPUT(IP,WALOC,HS,TM01,TM02,TM10,KLM,WLM)
USE DATAPOOL
USE W3SRC4MD
IMPLICIT NONE
INTEGER, INTENT(IN) :: IP
REAL(rkind), INTENT(IN) :: WALOC(NUMSIG,NUMDIR)
REAL(rkind), INTENT(OUT) :: HS,TM01,TM02,KLM,WLM,TM10
INTEGER ISMAX, ID, IS
REAL(rkind) AWW3(NSPEC), FL3(NUMDIR,NUMSIG)
REAL(rkind) WIND10, WINDTH, JAC
REAL(rkind) FMEAN1, WNMEAN, AMAX
REAL(rkind) F1MEAN, AKMEAN, XKMEAN
ISMAX = NUMSIG
CALL MEAN_PARAMETER(IP,WALOC,ISMAX,HS,TM01,TM02,TM10,KLM,WLM)
IF (DEP(IP) .gt. 0) THEN
IF (ISOURCE .eq. 1) THEN
DO IS = 1, NUMSIG
DO ID = 1, NUMDIR
AWW3(ID + (IS-1) * NUMDIR) = WALOC(IS,ID) * CG(IS,IP)
END DO
END DO
LLWS = .TRUE.
CALL SET_WIND( IP, WIND10, WINDTH )
CALL W3SPR4 ( AWW3, CG(:,IP), WK(:,IP), EMEAN(IP), FMEAN(IP), FMEAN1, WNMEAN, AMAX, WIND10, WINDTH, UFRIC(IP), USTDIR(IP), TAUWX(IP), TAUWY(IP), CD(IP), Z0(IP), ALPHA_CH(IP), LLWS, FMEANWS(IP))
HS = 4.0_rkind * SQRT(EMEAN(IP))
ELSE IF (ISOURCE .eq. 2) THEN
DO IS = 1, NUMSIG
JAC = PI2 * SPSIG(IS)
DO ID = 1, NUMDIR
FL3(ID,IS) = WALOC(IS,ID) * JAC
END DO
END DO
CALL FKMEAN_LOCAL(IP, FL3, EMEAN(IP), FMEAN(IP), F1MEAN, AKMEAN, XKMEAN)
HS = 4.0_rkind * SQRT(EMEAN(IP))
END IF
ELSE
HS = ZERO
END IF
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
!AR: This must be replaced by ST4_PRE or the certain WAM routine we are not consistent here this routine needs urgen deletion
SUBROUTINE MEAN_WAVE_PARAMETER_LOC(WALOC,CURTXYLOC,DEPLOC,WKLOC,HS,ETOT,SME01,SME10,KME01,KMWAM,KMWAM2)
USE DATAPOOL
IMPLICIT NONE
REAL(rkind), INTENT(IN) :: WALOC(NUMSIG,NUMDIR)
REAL(rkind), INTENT(IN) :: WKLOC(NUMSIG), DEPLOC
REAL(rkind), INTENT(IN) :: CURTXYLOC(2)
REAL(rkind), INTENT(OUT) :: SME01, SME10, KME01
REAL(rkind), INTENT(OUT) :: KMWAM, KMWAM2
REAL(rkind), INTENT(OUT) :: HS
INTEGER :: ID, IS
REAL(rkind) :: ACTOT, ETOT
REAL(rkind) :: ETOT_SPSIG
REAL(rkind) :: ETOT_WK
REAL(rkind) :: ETOT_ISQ_WK
REAL(rkind) :: ETOT_SQ_WK
REAL(rkind) :: Y(NUMSIG)
REAL(rkind) :: DS, ATAIL, ETAIL, ESIGTAIL
REAL(rkind) :: tmp(NUMSIG)
!
! total energy ...
!
ETOT = ZERO
do id = 1, NUMDIR
tmp(:) = WALOC(:,id) * spsig
ETOT = ETOT + tmp(1) * ONEHALF * ds_incr(1)*ddir
do is = 2, NUMSIG
ETOT = ETOT + ONEHALF*(tmp(is)+tmp(is-1))*ds_band(is)*ddir
end do
ETOT = ETOT + ONEHALF * tmp(NUMSIG) * ds_incr(NUMSIG)*ddir
end do
!
! if etot too small skip ...
!
if (etot .gt. verysmall) then
ACTOT = ZERO
y = ONE/SPSIG
do id = 1, NUMDIR
tmp(:) = WALOC(:,id) * spsig * y
ACTOT = ACTOT + tmp(1) * ONEHALF * ds_incr(1)*ddir
do is = 2, NUMSIG
ACTOT = ACTOT + ONEHALF*(tmp(is)+tmp(is-1))*ds_band(is)*ddir
end do
ACTOT = ACTOT + ONEHALF * tmp(NUMSIG) * ds_incr(NUMSIG)*ddir
end do
ETOT_SPSIG = ZERO
y = SIGPOW(:,1)
do id = 1, NUMDIR
tmp(:) = WALOC(:,id) * spsig * y
ETOT_SPSIG = ETOT_SPSIG + tmp(1) * ONEHALF * ds_incr(1)*ddir
do is = 2, NUMSIG
ETOT_SPSIG = ETOT_SPSIG + ONEHALF*(tmp(is)+tmp(is-1))*ds_band(is)*ddir
end do
ETOT_SPSIG = ETOT_SPSIG + ONEHALF * tmp(NUMSIG) * ds_incr(NUMSIG)*ddir
end do
ETOT_WK = ZERO
y = WKLOC(:)
do id = 1, NUMDIR
tmp(:) = WALOC(:,id) * spsig * y
ETOT_WK = ETOT_WK + ONEHALF * tmp(1) * ds_incr(1)*ddir
do is = 2, NUMSIG
ETOT_WK = ETOT_WK + ONEHALF*(tmp(is)+tmp(is-1))*ds_band(is)*ddir
end do
ETOT_WK = ETOT_WK + ONEHALF * tmp(NUMSIG) * ds_incr(NUMSIG)*ddir
end do
ETOT_ISQ_WK = ZERO
y = ONE/SQRT(WKLOC)
do id = 1, NUMDIR
tmp(:) = WALOC(:,id) * spsig * y
ETOT_ISQ_WK = ETOT_ISQ_WK + ONEHALF * tmp(1) * ds_incr(1)*ddir
do is = 2, NUMSIG -1
ETOT_ISQ_WK = ETOT_ISQ_WK+ONEHALF*(tmp(is)+tmp(is-1))*ds_band(is)*ddir
end do
ETOT_ISQ_WK = ETOT_ISQ_WK+ONEHALF*tmp(NUMSIG) * ds_incr(NUMSIG)*ddir
end do
ETOT_SQ_WK = ZERO
y = SQRT(WKLOC)
do id = 1, NUMDIR
tmp(:) = WALOC(:,id) * spsig * y
ETOT_SQ_WK = ETOT_SQ_WK + ONEHALF * tmp(1) * ds_incr(1)*ddir
do is = 2, NUMSIG -1
ETOT_SQ_WK = ETOT_SQ_WK + ONEHALF*(tmp(is)+tmp(is-1))*ds_incr(is)*ddir
end do
ETOT_SQ_WK = ETOT_SQ_WK + ONEHALF*tmp(NUMSIG) * ds_incr(NUMSIG)*ddir
end do
!
DS = SPSIG(NUMSIG) - SPSIG(NUMSIG-1)
ATAIL = SUM(WALOC(NUMSIG,:)) * SIGPOW(NUMSIG,1) * DDIR * DS
ETAIL = SUM(WALOC(NUMSIG,:)) * SIGPOW(NUMSIG,2) * DDIR * DS
ESIGTAIL = SUM(WALOC(NUMSIG,:)) * SIGPOW(NUMSIG,3) * DDIR * DS
!
! tail factors ...
!
ACTOT = ACTOT + TAIL_ARR(5) * ATAIL
ETOT = ETOT + TAIL_ARR(6) * ETAIL
ETOT_SPSIG = ETOT_SPSIG + TAIL_ARR(7) * ETAIL
ETOT_ISQ_WK = ETOT_ISQ_WK + TAIL_ARR(5) * ETAIL / SQRT(WKLOC(NUMSIG))
ETOT_SQ_WK = ETOT_SQ_WK + TAIL_ARR(5) * ETAIL * SQRT(WKLOC(NUMSIG))
ETOT_WK = ETOT_WK + TAIL_ARR(8) * ETAIL * WKLOC(NUMSIG)
!
! integral parameters ...
!
HS = MAX(ZERO,4.*SQRT(ETOT))
SME01 = ETOT_SPSIG / ETOT
KME01 = ETOT_WK / ETOT
SME10 = ETOT / ACTOT
KMWAM = (ETOT/ETOT_ISQ_WK)**2
KMWAM2 = (ETOT_SQ_WK/ETOT)**2
else
!
! no or too less energy ...
!
HS = ZERO
SME01 = ZERO
KME01 = 10.0_rkind
SME10 = ZERO
KMWAM = 10.0_rkind
KMWAM2 = 10.0_rkind
end if
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
!AR: This must be replaced by ST4_PRE or the certain WAM routine we are not consistent here this routine needs urgen deletion
SUBROUTINE WAVE_CURRENT_PARAMETER(IP,WALOC,UBOT,ORBITAL,BOTEXPER,TMBOT,CALLFROM)
USE DATAPOOL
IMPLICIT NONE
INTEGER, INTENT(IN) :: IP
REAL(rkind), INTENT(IN) :: WALOC(NUMSIG,NUMDIR)
CHARACTER(len=*), INTENT(IN) :: CALLFROM
REAL(rkind), INTENT(OUT) :: UBOT, ORBITAL, BOTEXPER, TMBOT
INTEGER :: ID, IS
REAL(rkind) :: ETOT_SKD
REAL(rkind) :: ETOT_SKDSIG, TMP(NUMSIG), Y(NUMSIG)
!
! integrals ...
!
IF (DEP(IP) .LT. DMIN) RETURN
ETOT_SKD = ZERO
ETOT_SKDSIG = ZERO
y = ONE/SINH(MIN(KDMAX,WK(:,IP)*DEP(IP)))**2
do id = 1, NUMDIR
tmp(:) = WALOC(:,id) * spsig * y
ETOT_SKD = ETOT_SKD + tmp(1) * ONEHALF * ds_incr(1)*ddir
do is = 2, NUMSIG -1
ETOT_SKD = ETOT_SKD + ONEHALF*(tmp(is)+tmp(is-1))*ds_band(is)*ddir
end do
ETOT_SKD = ETOT_SKD + tmp(NUMSIG) * ONEHALF * ds_incr(NUMSIG)*ddir
end do
y = SIGPOW(:,2)*ONE/SINH(MIN(KDMAX,WK(:,IP)*DEP(IP)))**2
do id = 1, NUMDIR
tmp(:) = WALOC(:,id) * spsig * y
ETOT_SKDSIG = ETOT_SKDSIG + tmp(1) * ONEHALF * ds_incr(1)*ddir
do is = 2, NUMSIG -1
ETOT_SKDSIG = ETOT_SKDSIG + ONEHALF*(tmp(is)+tmp(is-1))*ds_band(is)*ddir
end do
ETOT_SKDSIG = ETOT_SKDSIG + tmp(NUMSIG) * ONEHALF * ds_incr(NUMSIG)*ddir
end do
IF (ETOT_SKD .gt. verysmall) THEN
!
! integral parameters ...
!
UBOT = SQRT(ETOT_SKD)
ORBITAL = SQRT(2*ETOT_SKD)
BOTEXPER = SQRT(2*ETOT_SKDSIG)
TMBOT = PI2*SQRT(ETOT_SKDSIG/ETOT_SKD)
ELSE
!
! no or too less energy ...
!
UBOT = ZERO
ORBITAL = ZERO
BOTEXPER = ZERO
TMBOT = ZERO
ENDIF
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
!AR: This must be replaced by ST4_PRE or the certain WAM routine we are not consistent here this routine needs urgen deletion
SUBROUTINE WAVE_CURRENT_PARAMETER_LOC(WALOC,CURTXYLOC,DEPLOC,WKLOC,UBOT,ORBITAL,BOTEXPER,TMBOT)
USE DATAPOOL
IMPLICIT NONE
REAL(rkind), INTENT(IN) :: WALOC(NUMSIG,NUMDIR)
REAL(rkind), INTENT(IN) :: WKLOC(NUMSIG), DEPLOC
REAL(rkind), INTENT(IN) :: CURTXYLOC(2)
REAL(rkind), INTENT(OUT) :: UBOT, ORBITAL, BOTEXPER, TMBOT
INTEGER :: ID, IS
REAL(rkind) :: ETOT_SKD
REAL(rkind) :: ETOT_SKDSIG, TMP(NUMSIG), Y(NUMSIG)
!
! integrals ...
!
IF (DEPLOC .LT. DMIN) RETURN
ETOT_SKD = ZERO
ETOT_SKDSIG = ZERO
y = ONE/SINH(MIN(KDMAX,WKLOC*DEPLOC))**2
do id = 1, NUMDIR
tmp(:) = WALOC(:,id) * spsig * y
ETOT_SKD = ETOT_SKD + tmp(1) * ONEHALF * ds_incr(1)*ddir
do is = 2, NUMSIG -1
ETOT_SKD = ETOT_SKD + ONEHALF*(tmp(is)+tmp(is-1))*ds_band(is)*ddir
end do
ETOT_SKD = ETOT_SKD + tmp(NUMSIG) * ONEHALF * ds_incr(NUMSIG)*ddir
end do
y = SIGPOW(:,2)*ONE/SINH(MIN(KDMAX,WKLOC*DEPLOC))**2
do id = 1, NUMDIR
tmp(:) = WALOC(:,id) * spsig * y
ETOT_SKDSIG = ETOT_SKDSIG + tmp(1) * ONEHALF * ds_incr(1)*ddir
do is = 2, NUMSIG -1
ETOT_SKDSIG = ETOT_SKDSIG + ONEHALF*(tmp(is)+tmp(is-1))*ds_band(is)*ddir
end do
ETOT_SKDSIG = tmp(NUMSIG) * ONEHALF * ds_incr(NUMSIG)*ddir
end do
IF (ETOT_SKD .gt. verysmall) THEN
!
! integral parameters ...
!
UBOT = SQRT(ETOT_SKD)
ORBITAL = SQRT(2*ETOT_SKD)
BOTEXPER = SQRT(2*ETOT_SKDSIG)
TMBOT = PI2*SQRT(ETOT_SKDSIG/ETOT_SKD)
ELSE
!
! no or too less energy ...
!
UBOT = ZERO
ORBITAL = ZERO
BOTEXPER = ZERO
TMBOT = ZERO
ENDIF
!WRITE(*,'(9F15.4)') DEPLOC,SUM(WALOC),CURTXYLOC,SUM(WKLOC), ETOT_SKD, ETOT_SKDSIG
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
!AR: This must be replaced by ST4_PRE or the certain WAM routine we are not consistent here this routine needs urgen deletion
SUBROUTINE URSELL_NUMBER(HS,SME,DEPTH,URSELL)
USE DATAPOOL, ONLY : G9, DMIN, verysmall, rkind, ONE, TWO, ZERO
IMPLICIT NONE
REAL(rkind), INTENT(IN) :: HS, SME, DEPTH
REAL(rkind), INTENT(OUT) :: URSELL
IF (DEPTH .GT. DMIN .AND. SME .GT. verysmall .AND. HS .GT. verysmall) THEN
URSELL = (G9 * HS)/(TWO*SQRT(TWO)*SME**2*DEPTH**2)
ELSE
URSELL = ZERO
END IF
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************
!AR: This must be replaced by ST4_PRE or the certain WAM routine we are not consistent here this routine needs urgen deletion
SUBROUTINE WINDSEASWELLSEP( IP, WALOC, TM_W, CGP_W, CP_W, TP_W, LP_W, HS_W, KP_W )
USE DATAPOOL
IMPLICIT NONE
INTEGER, INTENT(IN) :: IP
REAL(rkind) , INTENT(IN) :: WALOC(NUMSIG,NUMDIR)
REAL(rkind) , INTENT(OUT) :: TM_W, CGP_W, CP_W, TP_W, LP_W, HS_W, KP_W
INTEGER :: ID, IS
REAL(rkind) :: ETOT
REAL(rkind) :: VEC2RAD, EFTOT, OMEG, WINDTH
REAL(rkind) :: EAD, DS, EHFR, EFTAIL, ETAIL, PTAIL_ARR, ACWIND(NUMSIG,NUMDIR)
REAL(rkind) :: UXD, ETOTF3, ETOTF4, FP, WN_W, WVC, WKDEP_W
WINDTH = VEC2RAD(WINDXY(IP,1),WINDXY(IP,2))
ETOT = ZERO
EFTAIL = ONE / (TAIL_ARR(1)-ONE)
DO ID = 1, NUMDIR ! Calculate wind sea energy ... weak criterion
DO IS = 1, NUMSIG
WVC = MyREAL(SPSIG(IS)/WK(IS,IP))
IF ( 1.2*UFRIC(IP)*COS(SPDIR(ID)-WINDTH)*(28./WVC) .LT. ONE) THEN
ACWIND(IS,ID) = ZERO ! Swell
ELSE
ACWIND(IS,ID) = WALOC(IS,ID) ! Wind Sea
END IF
END DO
DO IS = 2, NUMSIG
DS = SPSIG(IS) - SPSIG(IS-1)
EAD = ONEHALF*(SPSIG(IS)*ACWIND(IS,ID)+SPSIG(IS-1)*ACWIND(IS-1,ID))*DS*DDIR
ETOT = ETOT + EAD
END DO
IF (NUMSIG > 3) THEN
EHFR = WALOC(NUMSIG,ID) * SPSIG(NUMSIG)
ETOT = ETOT + DDIR * EHFR * SPSIG(NUMSIG) * EFTAIL
ENDIF
END DO
IF (ETOT .LT. VERYSMALL) RETURN
IF (ETOT > ZERO) THEN
HS_W = 4.0*SQRT(ETOT)
ELSE
HS_W = ZERO
END IF
ETOTF3 = ZERO
ETOTF4 = ZERO
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. VERYSMALL) THEN
FP = ETOTF3/ETOTF4*PI2
TP_W = ONE/FP/PI2
CALL WAVEKCG(DEP(IP), FP, WN_W, CP_W, KP_W, CGP_W)
!CALL ALL_FROM_TABLE(FP,DEP(IP),KP_W,CGP_W,WKDEP_W,WN_W,CP_W)
LP_W = PI2/KP_w
ELSE
FP = ZERO
CP_W = ZERO
KP_W = 10.0_rkind
CGP_W = ZERO
LP_W = ZERO
END IF
ETOT = ZERO
EFTOT = ZERO
PTAIL_ARR = TAIL_ARR(1) - ONE
ETAIL = ONE / (PTAIL_ARR * (ONE + PTAIL_ARR * (FRINTH-ONE)))
PTAIL_ARR = TAIL_ARR(1) - 2.
EFTAIL = ONE / (PTAIL_ARR * (ONE + PTAIL_ARR * (FRINTH-ONE)))
DO ID = 1, NUMDIR
UXD = CURTXY(IP,1)*COSTH(ID) + CURTXY(IP,2)*SINTH(ID)
DO IS = 1, NUMSIG
OMEG = SPSIG(IS) + WK(IS,IP) * UXD
EAD = FRINTF * SIGPOW(IS,2) * ACWIND(IS,ID)
ETOT = ETOT + EAD
EFTOT = EFTOT + EAD * OMEG
ENDDO
IF (NUMSIG .GT. 3) THEN
EAD = SIGPOW(NUMSIG,2) * ACWIND(NUMSIG,ID)
ETOT = ETOT + ETAIL * EAD
EFTOT = EFTOT + EFTAIL * OMEG * EAD
ENDIF
ENDDO
IF (EFTOT.GT.ZERO) THEN
TM_W = PI2 * ETOT / EFTOT
ELSE
TM_W = ZERO
ENDIF
END SUBROUTINE
!**********************************************************************
!* *
!**********************************************************************