-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathw3fld1md.ftn
1431 lines (1413 loc) · 53 KB
/
w3fld1md.ftn
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
!/ ------------------------------------------------------------------- /
Module W3FLD1MD
!/
!/ +-----------------------------------+
!/ | WAVEWATCH III NOAA/NCEP/NOPP |
!/ | B. G. Reichl |
!/ | FORTRAN 90 |
!/ | Last update : 18-Mar-2015 |
!/ +-----------------------------------+
!/
!/ 01-Jul-2013 : Origination. ( version 4.xx )
!/ 18-Mar-2015 : Clean-up/prepare for distribution ( version 5.12 )
!/ 15-Jan-2016 : Updates ( version 5.12 )
!/ ( B. G. Reichl )
!/ 27-Jul-2016 : Added Charnock output (J.Meixner) ( version 5.12 )
!/ 22-Jun-2018 : updated SIG2WN subroutine (X.Chen) ( version 6.06 )
!/ modified the range of wind profile computation;
!/ corrected direction of the shortest waves
!/
!/
!/ Copyright 2009 National Weather Service (NWS),
!/ National Oceanic and Atmospheric Administration. All rights
!/ reserved. WAVEWATCH III is a trademark of the NWS.
!/ No unauthorized use without permission.
!/
! 1. Purpose :
!
! This Module contains routines to compute the wind stress vector
! from the wave spectrum, the wind vector, and the lower atmospheric
! stability (the form included here is for neutral conditions, but
! the structure needed to include stability is included in comments).
! The stress calculated via this subroutine is
! intended for coupling to serve as the boundary condition
! between the ocean and atmosphere, and (for now)
! and has no impact on the wave spectrum calculated.
! The calculation in w3fld1 is based on the method
! presented in Reichl, Hara, and Ginis (2014), "Sea State Dependence
! of the Wind Stress under Hurricane Winds."
!
! 2. Variables and types :
!
! Not applicable.
!
! 3. Subroutines and functions :
!
! Name Type Scope Description
! ----------------------------------------------------------------
! W3FLD1 Subr. Public Reichl et al. 2014 stress calculation
! INFLD1 Subr. Public Corresponding initialization routine.
! APPENDTAIL Subr. Public Modification of tail for calculation
! SIG2WN Subr. Public Depth-dependent dispersion relation
! WND2Z0M Subr. Public Wind to roughness length
! ----------------------------------------------------------------
!
! 4. Subroutines and functions used :
!
! Name Type Module Description
! ----------------------------------------------------------------
! STRACE Subr. W3SERVMD Subroutine tracing.
! ----------------------------------------------------------------
!
! 5. Remarks :
!
! 6. Switches :
!
! !/S Enable subroutine tracing.
! !/
!
! 7. Source code :
!/
!/ ------------------------------------------------------------------- /
!/
!
PUBLIC
! Tail_Choice: Chose the method to determine the level of the tail
INTEGER, SAVE :: Tail_Choice
REAL, SAVE :: Tail_Level !if Tail_Choice=0, tail is constant
REAL, SAVE :: Tail_transition_ratio1! freq/fpi where tail
! adjustment begins
REAL, SAVE :: Tail_transition_ratio2! freq/fpi where tail
! adjustment ends
!/
CONTAINS
!/ ------------------------------------------------------------------- /
SUBROUTINE W3FLD1( ASPC, FPI, WNDX,WNDY, ZWND, &
DEPTH, RIB, UST, USTD, Z0, TAUNUX,TAUNUY, CHARN)
!/
!/ +-----------------------------------+
!/ | WAVEWATCH III NOAA/NCEP/NOPP |
!/ | B. G. Reichl |
!/ | FORTRAN 90 |
!/ | Last update : 18-Mar-2015 |
!/ +-----------------------------------+
!/
!/ 01-Jul-2013 : Origination. ( version 4.xx )
!/ 18-Mar-2015 : Prepare for submission ( version 5.12 )
!/
! 1. Purpose :
!
! Diagnostic stress vector calculation from wave spectrum, lower
! atmosphere stability, and wind vector (at some given height).
! The height of wind vector is assumed to be within the constant
! stress layer. These parameterizations are meant to be performed
! at wind speeds > 10 m/s, and may not converge for extremely young
! seas (i.e. starting from flat sea conditions).
!
! 2. Method :
! See Reichl et al. (2014).
!
! 3. Parameters :
!
! Parameter list
! ----------------------------------------------------------------
! ASPC Real I 1-D Wave action spectrum.
! FPI Real I Peak input frequency.
! WNDX Real I X-dir wind (assumed referenced to current)
! WNDY Real I Y-dir wind (assumed referenced to current)
! ZWND Real I Wind height.
! DEPTH Real I Water depth.
! RIB REAL I Bulk Richardson in lower atmosphere
! (for determining stability in ABL to get
! 10 m neutral wind)
! TAUNUX Real 0 X-dir viscous stress (guessed from prev.)
! TAUNUY Real 0 Y-dir viscous stress (guessed from prev.)
! UST Real O Friction velocity.
! USTD Real O Direction of friction velocity.
! Z0 Real O Surface roughness length
! CHARN Real O,optional Charnock parameter
! ----------------------------------------------------------------
!
! 4. Subroutines used :
!
! Name Type Module Description
! ----------------------------------------------------------------
! STRACE Subr. W3SERVMD Subroutine tracing.
! ----------------------------------------------------------------
!
! 5. Called by :
!
! Name Type Module Description
! ----------------------------------------------------------------
! W3ASIM Subr. W3ASIMMD Air-sea interface module.
! W3EXPO Subr. N/A Point output post-processor.
! GXEXPO Subr. N/A GrADS point output post-processor.
! ----------------------------------------------------------------
!
! 6. Error messages :
!
! None.
!
! 7. Remarks :
!
! 8. Structure :
!
! See source code.
!
! 9. Switches :
!
! !/S Enable subroutine tracing.
!
! 10. Source code :
!
!/ ------------------------------------------------------------------- /
USE CONSTANTS, ONLY: GRAV, DWAT, DAIR, TPI, PI, KAPPA
USE W3GDATMD, ONLY: NK, NTH, NSPEC, SIG, DTH, XFR, TH
USE W3ODATMD, ONLY: NDSE
USE W3SERVMD, ONLY: EXTCDE
!/S USE W3SERVMD, ONLY: STRACE
!/
IMPLICIT NONE
!/
!/ ------------------------------------------------------------------- /
!/ Parameter list
!/
REAL, INTENT(IN) :: ASPC(NSPEC), WNDX, WNDY, &
ZWND, DEPTH, RIB, FPI
REAL, INTENT(OUT) :: UST, USTD, Z0
REAL, INTENT(OUT), OPTIONAL :: CHARN
REAL, INTENT(INOUT) :: TAUNUX, TAUNUY
!/
!/ ------------------------------------------------------------------- /
!/ Local parameters
!/
REAL, PARAMETER :: NU=0.105/10000.0
REAL, PARAMETER :: DELTA=0.03
! Commonly used parameters
REAL :: wnd_in_mag, wnd_in_dir
!For Calculating Tail
REAL :: KMAX, KTAILA, KTAILB, KTAILC
REAL :: SAT, z01, z02, u10
LOGICAL :: ITERFLAG
INTEGER :: COUNT
!For Iterations
REAL :: DTX, DTY, iter_thresh, &
USTSM, Z0SM, Z1
!For stress calculation
REAL :: WAGE, CBETA, BP, CD, &
USTRB, ANGDIF, USTAR, ZNU, &
TAUT, TAUX, TAUY, BETAG, TAUDIR, &
TAUDIRB
!For wind profile calculation
REAL :: UPROFV, VPROFV
!For wind profile iteration
REAL :: WND_1X, WND_1Y, &
WND_2X, WND_2Y, &
WND_3X, WND_3Y, &
DIFU10XX, DIFU10YX, DIFU10XY, DIFU10YY, &
FD_A, FD_B, FD_C, FD_D, &
DWNDX, DWNDY, &
APAR, CH,UITV, VITV,USTL,&
CK
!For adding stability to wind profile
REAL :: WND_TOP, ANG_TOP, WND_PA, WND_PE, &
WND_PEx, WND_PEy, WND_PAx, WND_PAy, &
CDM
INTEGER :: NKT, K, T, Z2, ITER, ZI, ZII, &
I, CTR, ITERATION, KA1, KA2, &
KA3, KB
! For defining extended spectrum with appended tail.
REAL, ALLOCATABLE, DIMENSION(:) :: WN, DWN, CP,SIG2
REAL, ALLOCATABLE, DIMENSION(:,:) :: SPC2
REAL, ALLOCATABLE, DIMENSION(:) :: TLTN, TLTE, TAUD, &
TLTND, &
TLTED, ZOFK, UPROF, VPROF, &
FTILDE, UP1, VP1, UP, VP, &
TLTNA, TLTEA
!/S INTEGER, SAVE :: IENT = 0
LOGICAL :: FSFL1,FSFL2, CRIT1, CRIT2
LOGICAL :: IT_FLAG1, IT_FLAG2
LOGICAL, SAVE :: FIRST = .TRUE.
!/
!/ ------------------------------------------------------------------- /
!/
!/S CALL STRACE (IENT, 'W3FLD1')
!
! 0. Initializations ------------------------------------------------ *
!
! **********************************************************
! *** The initialization routine should include all ***
! *** initialization, including reading data from files. ***
! **********************************************************
!
IF ( FIRST ) THEN
CALL INFLD
FIRST = .FALSE.
END IF
wnd_in_mag = sqrt( wndx**2 + wndy**2 )
wnd_in_dir = atan2(wndy, wndx)
!----------------------------------------------------------+
! Assume wind input is neutral 10 m wind. If wind input +
! is not 10 m, tail level will need to be calculated based +
! on esimation of 10 m wind. +
!----------------------------------------------------------+
u10 = wnd_in_mag
! - Get tail level
if (Tail_Choice.eq.0) then
SAT=Tail_Level
elseif (Tail_Choice.eq.1) then
CALL WND2SAT(U10,SAT)
endif
!
! 1. Attach Tail ---------------------------------------------------- *
!
! If the depth remains constant, the allocation could be limited to the
! first time step. Since this code is designed for coupled
! implementation where the water level can change, I keep it the
! allocation on every time step. When computational efficiency is
! important, this process may be rethought.
!
! i. Find maximum wavenumber of input spectrum
call sig2wn(sig(nk),depth,kmax)
NKT = NK
! ii. Find additional wavenumber bins to extended to cm scale waves
DO WHILE ( KMAX .LT. 366.0 )
NKT = NKT + 1
KMAX = ( KMAX * XFR**2 )
ENDDO!K<366
! iii. Allocate new "extended" spectrum
ALLOCATE( WN(NKT), DWN(NKT), CP(NKT), SIG2(NKT),SPC2(NKT,NTH), &
TLTN(NKT), TLTE(NKT), TAUD(NKT), &
TLTND(NKT), TLTED(NKT), ZOFK(NKT), UPROF(NKT+1),&
VPROF(NKT+1), FTILDE(NKT), UP1(NKT+1),VP1(NKT+1), &
UP(NKT+1), VP(NKT+1), TLTNA(NKT),TLTEA(NKT))
!
! 1a. Build Discrete Wavenumbers for defining extended spectrum on---- *
!
!i. Copy existing sig to extended sig2, calculate phase speed.
DO K = 1, NK !existing spectrum
call sig2wn(sig(k),depth,wn(k))
CP(K) = ( SIG(K) / WN(K) )
sig2(k) = sig(k)
ENDDO!K
!ii. Calculate extended sig2 and phase speed.
DO K = ( NK + 1 ), ( NKT) !extension
sig2(k) = sig2(k-1) *XFR
call sig2wn(sig2(k),depth,wn(k))
CP(K) = SIG2(K) / WN(K)
ENDDO!K
!iii. Calculate dk's for integrations.
DO K = 1, NKT-1
DWN(K) = WN(K+1) - WN(K)
ENDDO
DWN(NKT) = WN(NKT)*XFR**2 - WN(NKT)
!
! 1b. Attach initial tail--------------------------------------------- *
!
!i. Convert action spectrum to variance spectrum
! SPC(k,theta) = A(k,theta) * sig(k)
! This could be redone for computational efficiency
I=0
DO K=1, NK
DO T=1, NTH
I = I + 1
SPC2(K,T) = ASPC(I) * SIG(K)
ENDDO!T
ENDDO!K
!ii. Extend k^-3 tail to extended spectrum
DO K=NK+1, NKT
DO T=1, NTH
SPC2(K,T)=SPC2(NK,T)*WN(NK)**3.0/WN(K)**(3.0)
ENDDO!T
ENDDO!K
!
! 1c. Calculate transitions for new (constant saturation ) tail ------ *
!
!
!i. Find wavenumber for beginning spc level transition to tail
call sig2wn (FPI*TPI*tail_transition_ratio1,depth,ktaila )
!ii. Find wavenumber for ending spc level transition to tail
call sig2wn (FPI*TPI*tail_transition_ratio2,depth,ktailb )
!iii. Find wavenumber for ending spc direction transition to tail
KTAILC= KTAILB * 2.0
!iv. Find corresponding indices of wavenumber transitions
KA1 = 2 ! Do not modify 1st wavenumber bin
DO WHILE ( ( KTAILA .GE. WN(KA1) ) .AND. (KA1 .LT. NKT-6) )
KA1 = KA1 + 1
ENDDO
KA2 = KA1+2
DO WHILE ( ( KTAILB .GE. WN(KA2) ) .AND. (KA2 .LT. NKT-4) )
KA2 = KA2 + 1
ENDDO
KA3 = KA2+2
DO WHILE ( ( KTAILC .GE. WN(KA3)) .AND. (KA3 .LT. NKT-2) )
KA3 = KA3 + 1
ENDDO
!v. Call subroutine to perform actually tail truncation
! only if there is some energy in spectrum
CALL APPENDTAIL(SPC2, WN, NKT, KA1, KA2, KA3,&
wnd_in_dir, SAT)
! Spectrum is now set for stress integration
!
! 2. Prepare for iterative calculation of wave-form stress----------- *
!
DTX = 0.00005
DTY = 0.00005
iter_thresh = 0.001
!
! 2a. Calculate initial guess for viscous stress from smooth-law------ *
! (Would be preferable to use prev. step)
!
Z0SM = 0.001 !Guess
IT_FLAG1 = .true.
ITERATION = 0
DO WHILE( IT_FLAG1 )
ITERATION = ITERATION + 1
Z1 = Z0SM
USTSM = KAPPA * wnd_in_mag / ( LOG( ZWND / Z1 ) )
Z0SM = 0.132 * NU / USTSM
IF ( (ABS( Z0SM - Z1 ) .LT. 10.0**(-6)) .OR.&
( ITERATION .GT. 5 )) THEN
IT_FLAG1 = .false.
ENDIF
ENDDO
ITERATION = 1
! Guessed values of viscous stress
TAUNUX = USTSM**2 * DAIR * wndx / wnd_in_mag
TAUNUY = USTSM**2 * DAIR * wndy / wnd_in_mag
!
! 3. Enter iterative calculation of wave form/skin stress---------- *
!
IT_FLAG1 = .true.
DO WHILE (IT_FLAG1)
DO ITER=1, 3 !3 loops for TAUNU iteration
Z2 = NKT
! First : TAUNUX + DX
IF (ITER .EQ. 1) THEN
TAUNUX = TAUNUX + DTX
! Second : TAUNUY + DY
ELSEIF (ITER .EQ. 2) THEN
TAUNUX = TAUNUX - DTX
TAUNUY = TAUNUY + DTY
! Third : unmodified
ELSEIF (ITER .EQ. 3) THEN
TAUNUY = TAUNUY - DTY
ENDIF
! Near surface turbulent stress = taunu
TLTN(1) = TAUNUY
TLTE(1) = TAUNUX
CALL APPENDTAIL(SPC2, WN, NKT, KA1, KA2, KA3,&
atan2(TAUNUY,TAUNUX), SAT)
!|---------------------------------------------------------------------|
!|-----Calculate first guess at growth rate and local turbulent stress-|
!|-----for integration as a function of wavedirection------------------|
!|---------------------------------------------------------------------|
DO ZI = 2, NKT
USTL=0.0
TLTND(zi)=0.0
TLTED(zi)=0.0
Z2 = Z2 - 1
! Use value of prev. wavenumber/height
TAUD(ZI) = ATAN2( TLTN(ZI-1), TLTE(ZI-1))
USTL = SQRT( SQRT( TLTN(ZI-1)**2 + TLTE(ZI-1)**2 )/ DAIR )
DO T = 1, NTH
ANGDIF=TAUD(ZI)-TH(T) !stress/wave angle
IF ( COS( ANGDIF ) .GE. 0.0 ) THEN !Waves aligned
WAGE = CP(Z2) / (USTL)
! First, waves much slower than wind.
IF ( WAGE .LT. 10. ) THEN
CBETA = 25.0
! Transition from waves slower than wind to faster
ELSEIF ( ( WAGE .GE. 10.0 ) .AND. &
( WAGE .LE. 25.0 ) ) THEN
CBETA = 10.0 + 15.0 * COS( PI * ( WAGE - 10.0 ) &
/ 15.0 )
! Waves faster than wind
ELSEIF ( WAGE .GT. 25.0 ) THEN
CBETA = -5.0
ENDIF
! Waves opposing wind
ELSE
CBETA = -25.0
ENDIF
!Integrate turbulent stress
TLTND(ZI) =TLTND(ZI)+( SIN( TH(T) ) * COS( ANGDIF )**2)&
* CBETA * SPC2(Z2,T) * &
SQRT( TLTE(ZI-1)**2 + TLTN(ZI-1)**2.0 ) &
* ( WN(Z2)**2.0 )*DTH
TLTED(ZI) = TLTED(ZI)+(COS( TH(T) ) * COS( ANGDIF )**2)&
* CBETA * SPC2(Z2,T) * &
SQRT( TLTE(ZI-1)**2 + TLTN(ZI-1)**2.0 ) &
* ( WN(Z2)**2.0 )*DTH
ENDDO
!|---------------------------------------------------------------------|
!|-----Complete the integrations---------------------------------------|
!|---------------------------------------------------------------------|
IF (ZI .EQ. 2) THEN
!First turbulent stress bin above taunu
TLTNA(ZI) = TLTND(ZI) * DWN(Z2) * 0.5
TLTEA(ZI) = TLTED(ZI) * DWN(Z2) * 0.5
ELSE
TLTNA(ZI)=(TLTND(ZI)+TLTND(ZI-1))*0.5*DWN(Z2)
TLTEA(ZI)=(TLTED(ZI)+TLTED(ZI-1))*0.5*DWN(Z2)
ENDIF
TLTN(ZI)=TLTN(ZI-1)+TLTNA(ZI)
TLTE(ZI)=TLTE(ZI-1)+TLTEA(ZI)
ENDDO
TAUY=TLTN(NKT)
TAUX=TLTE(NKT)
! This is the first guess at the stress.
!|---------------------------------------------------------------------|
!|----Iterate til convergence------------------------------------------|
!|---------------------------------------------------------------------|
USTRB=SQRT(SQRT(TAUY**2.0+TAUX**2.0)/DAIR)
TAUDIRB=atan2(TAUY,TAUX)
IT_FLAG2 = .TRUE.
CTR=1
DO WHILE ( (IT_FLAG2) .AND. ( CTR .LT. 10 ) )
Z2=NKT+1
DO ZI=1, NKT
Z2=Z2-1
USTL=0.0
TLTED(zi)=0.0
TLTND(zi)=0.0
FTILDE(z2)=0.0
TAUD(ZI) = ATAN2(TLTN(ZI),TLTE(ZI))
USTL = SQRT(SQRT(TLTN(ZI)**2+TLTE(ZI)**2)/DAIR)
DO T=1, NTH
BETAG=0.0
ANGDIF = TAUD(ZI)-TH(T)
IF ( COS( ANGDIF ) .GE. 0.0 ) THEN
WAGE = CP(Z2) / (USTL)
IF ( WAGE .LT. 10 ) THEN
CBETA = 25.0
ELSEIF ( ( WAGE .GE. 10.0 ) .AND. &
( WAGE .LE. 25.0 ) ) THEN
CBETA = 10.0 + 15.0 * COS( PI * ( WAGE - 10.0 ) &
/ 15.0 )
ELSEIF ( WAGE .GT. 25.0 ) THEN
CBETA = -5.0
ENDIF
ELSE
CBETA = -25.0
ENDIF
BP = SQRT( (COS( TH(T) ) * COS( ANGDIF )**2.0)**2.0 &
+ (SIN( TH(T) ) * COS( ANGDIF )**2.0)**2.0 )
BETAG=BP*CBETA*SQRT(TLTE(ZI)**2.0+TLTN(ZI)**2.0) &
/(DWAT)*SIG2(Z2)/CP(Z2)**2
FTILDE(Z2) = FTILDE(Z2) + BETAG * DWAT * GRAV &
* SPC2(Z2,T) * DTH
TLTND(zi) =tltnd(zi)+ (SIN( TH(T) ) * COS( ANGDIF )**2.0)&
* CBETA * SPC2(Z2,T) * SQRT( &
TLTE(ZI)**2.0 + TLTN(ZI)**2.0 ) * &
( WN(Z2)**2.0 )*dth
TLTED(zi) = tlted(zi)+(COS( TH(T) ) * COS( ANGDIF )**2.0)&
* CBETA * SPC2(Z2,T) * SQRT( &
TLTE(ZI)**2.0 + TLTN(ZI)**2.0 ) * &
( WN(Z2)**2.0 )*dth
ENDDO
IF (ZI .EQ. 1) THEN
TLTNA(ZI)=TLTND(ZI)*DWN(Z2)*0.5
TLTEA(ZI)=TLTED(ZI)*DWN(Z2)*0.5
ELSE
TLTNA(ZI)=(TLTND(ZI)+TLTND(ZI-1))*0.5*DWN(Z2)
TLTEA(ZI)=(TLTED(ZI)+TLTED(ZI-1))*0.5*DWN(Z2)
ENDIF
IF (ZI.GT.1) then
TLTN(ZI)=TLTN(ZI-1)+TLTNA(ZI)
TLTE(ZI)=TLTE(ZI-1)+TLTEA(ZI)
else
TLTN(ZI)=TAUNUY+TLTNA(ZI)
TLTE(ZI)=TAUNUX+TLTEA(ZI)
endif
ENDDO
TAUY=TLTN(NKT) !by NKT full stress is entirely
TAUX=TLTE(NKT) !from turbulent stress
TAUT=SQRT(TAUY**2.0+TAUX**2.0)
USTAR=SQRT(SQRT(TAUY**2.0+TAUX**2.0)/DAIR)
TAUDIR=atan2(TAUY, TAUX)
! Note: add another criterion (stress direction) for iteration.
CRIT1=(ABS(USTAR-USTRB)*100.0)/((USTAR+USTRB)*0.5) .GT. 0.1
CRIT2=(ABS(TAUDIR-TAUDIRB)*100.0/(TAUDIR+TAUDIRB)*0.5) .GT. 0.1
IF (CRIT1 .OR. CRIT2) THEN
! IF ((ABS(USTAR-USTRB)*100.0)/((USTAR+USTRB)*0.5) .GT. 0.1) THEN
USTRB=USTAR
TAUDIRB=TAUDIR
CTR=CTR+1
ELSE
IT_FLAG2 = .FALSE.
ENDIF
ENDDO
! Note: search for the top of WBL from top to bottom (avoid problems
! caused by for very long swell)
KB=NKT
DO WHILE(((TLTN(KB)**2+TLTE(KB)**2)/(TAUX**2+TAUY**2)).GT. &
.99)
KB=KB-1
ENDDO
KB=KB+1
!|---------------------------------------------------------------------|
!|----Now begin work on wind profile-----------------------------------|
!|---------------------------------------------------------------------|
DO I=1,NKT
ZOFK(I)=DELTA/WN(I)
ENDDO
ZNU=0.1 * 1.45E-5 / SQRT(SQRT(TAUNUX**2.0+TAUNUY**2.0)/DAIR)
UPROF(1:NKT+1)=0.0
VPROF(1:NKT+1)=0.0
UPROFV=0.0
VPROFV=0.0
ZI=1
Z2=NKT
UP1(ZI) = ( ( ( WN(Z2)**2 / DELTA ) * FTILDE(z2) ) + &
( DAIR / ( ZOFK(Z2) * KAPPA ) ) * ( SQRT( &
TLTN(ZI)**2 + TLTE(ZI)**2 ) / DAIR )**(3/2) ) &
* ( TLTE(ZI) ) / ( TLTE(ZI) * TAUX &
+ TLTN(ZI) * TAUY )
VP1(ZI) = ( ( ( WN(Z2)**2 / DELTA ) * FTILDE(z2) ) + &
( DAIR / ( ZOFK(Z2) * KAPPA ) ) * ( SQRT ( &
TLTN(ZI)**2 + TLTE(ZI)**2 ) / DAIR )**(3/2) ) &
* ( TLTN(ZI) ) / ( TLTE(ZI) * TAUX &
+ TLTN(ZI) * TAUY )
UP(ZI) = UP1(ZI)
VP(ZI) = VP1(ZI)
UPROF(ZI) = DAIR / KAPPA * ( SQRT( TAUNUX**2.0 + TAUNUY**2.0 ) &
/ DAIR )**(1.5) * ( TAUNUX / ( TAUX * &
TAUNUX + TAUY * TAUNUY ) ) * LOG( &
ZOFK(Z2) / ZNU )
VPROF(ZI) = DAIR / KAPPA * ( SQRT( TAUNUX**2.0 + TAUNUY**2.0 ) &
/ DAIR )**(1.5) * ( TAUNUY / ( TAUX * &
TAUNUX + TAUY * TAUNUY ) ) * LOG( &
ZOFK(Z2) / ZNU )
!Noted: wind profile computed till the inner layer height of the longest
!wave, not just to the top of wave boundary layer (previous)
DO ZI=2, NKT
Z2 = Z2 - 1
UP1(ZI) = ( ( ( WN(Z2)**2.0 / DELTA ) * FTILDE(Z2) ) + &
( DAIR / ( ZOFK(Z2) * KAPPA ) ) * ( SQRT( &
TLTN(ZI)**2.0 + TLTE(ZI)**2.0 ) / DAIR )**(1.5) ) &
* ( TLTE(ZI) ) / ( TLTE(ZI) * TAUX + &
TLTN(ZI) * TAUY )
VP1(ZI) = ( ( ( WN(Z2)**2.0 / DELTA ) * FTILDE(Z2) ) + &
( DAIR / ( ZOFK(Z2) * KAPPA ) ) * ( SQRT( &
TLTN(ZI)**2.0 + TLTE(ZI)**2.0 ) / DAIR )**(1.5) ) &
* ( TLTN(ZI) ) / ( TLTE(ZI) * TAUX + &
TLTN(ZI) * TAUY )
UP(ZI) = UP1(ZI) * 0.5 + UP1(ZI-1) * 0.5
VP(ZI) = VP1(ZI) * 0.5 + VP1(ZI-1) * 0.5
UPROF(ZI) = UPROF(ZI-1) + UP(ZI) * ( ZOFK(Z2) - ZOFK(Z2+1) )
VPROF(ZI) = VPROF(ZI-1) + VP(ZI) * ( ZOFK(Z2) - ZOFK(Z2+1) )
ENDDO
!|---------------------------------------------------------------------|
!|----Iteration completion/checks--------------------------------------|
!|---------------------------------------------------------------------|
!ZI = ( KB + 1 )
! Now solving for 'ZWND' height wind
UPROF(NKT+1) = UPROF(NKT) + ( SQRT( SQRT( TAUY**2.0 + &
TAUX**2.0 ) / DAIR ) ) / KAPPA * TAUX &
/ SQRT( TAUY**2.0 +TAUX**2.0 ) * LOG( ZWND &
/ ZOFK(Z2) )
VPROF(NKT+1) = VPROF(NKT) + ( SQRT( SQRT( TAUY**2.0 + &
TAUX**2.0 ) / DAIR ) ) / KAPPA * TAUY &
/ SQRT( TAUY**2.0 +TAUX**2.0 ) * LOG( ZWND &
/ ZOFK(Z2) )
IF (ITER .EQ. 3) THEN
WND_1X = UPROF(NKT+1)
WND_1Y = VPROF(NKT+1)
ELSEIF (ITER .EQ. 2) THEN
WND_2X = UPROF(NKT+1)
WND_2Y = VPROF(NKT+1)
ELSEIF (ITER .EQ. 1) THEN
WND_3X = UPROF(NKT+1)
WND_3Y = VPROF(NKT+1)
ENDIF
!-------------------------------------+
! Guide for adding stability effects +
!-------------------------------------+
!Get Wind at top of wave boundary layer
! WND_TOP=SQRT(UPROF(KB)**2+VPROF(KB)**2)
! Get Wind Angle at top of wave boundary layer
! ANG_TOP=ATAN2(VPROF(KB),UPROF(KB))
! Stress and direction
! USTD = ATAN2(TAUY,TAUX)
! UST = SQRT( SQRT( TAUX**2 + TAUY**2 ) / DAIR)
! Calclate along (PA) and across (PE) wind components
! WND_PA=WND_TOP*COS(ANG_TOP-USTD)
! WND_PE=WND_TOP*SIN(ANG_TOP-USTD)
! Calculate cartesian aligned wind
! WND_PAx=WND_PA*cos(ustd)
! WND_PAy=WND_PA*sin(USTd)
!Calculate cartesion across wind
! WND_PEx=WND_PE*cos(ustd+pi/2.)
! WND_PEy=WND_PE*sin(ustd+pi/2.)
!----------------------------------------------------+
! If a non-neutral profile is used the effective z0 +
! should be computed. This z0 can then be used +
! with stability information to derive a Cd, which +
! can be used to project the along-stress wind to +
! the given height. +
! i.e.: Assume neutral inside WBL calculate Z0 +
! Z0=ZOFK(Z2)*EXP(-WND_PA*kappa/UST) +
! WND_PA=UST/SQRT(CDM) +
!----------------------------------------------------+
! WND_PAx=WND_PA*cos(ustd)
! WND_PAy=WND_PA*sin(USTd)
! IF (ITER .EQ. 3) THEN
! WND_1X = WND_PAx+WND_PEx
! WND_1Y = WND_PAy+WND_PEy
! ELSEIF (ITER .EQ. 2) THEN
! WND_2X = WND_PAx+WND_PEx
! WND_2Y = WND_PAy+WND_PEy
! ELSEIF (ITER .EQ. 1) THEN
! WND_3X = WND_PAx+WND_PEx
! WND_3Y = WND_PAy+WND_PEy
! ENDIF
ENDDO
ITERATION = ITERATION + 1
DIFU10XX = WND_3X - WND_1X
DIFU10YX = WND_3Y - WND_1Y
DIFU10XY = WND_2X - WND_1X
DIFU10YY = WND_2Y - WND_1Y
FD_A = DIFU10XX / DTX
FD_B = DIFU10XY / DTY
FD_C = DIFU10YX / DTX
FD_D = DIFU10YY / DTY
DWNDX = - WNDX + WND_1X
DWNDY = - WNDY + WND_1Y
UITV = ABS( DWNDX )
VITV = ABS( DWNDY )
CH = SQRT( UITV**2.0 + VITV**2.0 )
IF (CH .GT. 15.) THEN
APAR = 0.5 / ( FD_A * FD_D - FD_B * FD_C )
ELSE
APAR = 1.0 / ( FD_A * FD_D - FD_B * FD_C )
ENDIF
CK=4.
IF (((VITV/MAX(ABS(WNDY),CK) .GT. iter_thresh) .OR. &
(UITV/MAX(ABS(WNDX),CK) .GT. iter_thresh)) .AND. &
(ITERATION .LT. 2)) THEN
TAUNUX = TAUNUX - APAR * ( FD_D * DWNDX - FD_B * DWNDY )
TAUNUY = TAUNUY - APAR * ( -FD_C * DWNDX +FD_A * DWNDY )
ELSEIF (((VITV/MAX(ABS(WNDY),CK) .GT. iter_thresh) .OR. &
(UITV/MAX(ABS(WNDX),CK) .GT. iter_thresh)) .AND. &
(ITERATION .LT. 24)) THEN
iter_thresh = 0.001
TAUNUX = TAUNUX - APAR * ( FD_D * DWNDX - FD_B * DWNDY )
TAUNUY = TAUNUY - APAR * ( -FD_C * DWNDX +FD_A * DWNDY )
ELSEIF (((VITV/MAX(ABS(WNDY),CK) .GT. iter_thresh) .OR. &
(UITV/MAX(ABS(WNDX),CK) .GT. iter_thresh)) .AND. &
(ITERATION .LT. 26)) THEN
iter_thresh = 0.01
TAUNUX = TAUNUX - APAR * ( FD_D * DWNDX - FD_B * DWNDY )
TAUNUY = TAUNUY - APAR * ( -FD_C * DWNDX +FD_A * DWNDY )
ELSEIF (((VITV/MAX(ABS(WNDY),CK) .GT. iter_thresh) .OR. &
(UITV/MAX(ABS(WNDX),CK) .GT. iter_thresh)) .AND. &
(ITERATION .LT. 30)) THEN
iter_thresh = 0.05
TAUNUX = TAUNUX - APAR * ( FD_D * DWNDX - FD_B * DWNDY )
TAUNUY = TAUNUY - APAR * ( -FD_C * DWNDX +FD_A * DWNDY )
ELSEIF (ITERATION .GE. 30) THEN
write(*,*)'Attn: W3FLD1 not converged.'
write(*,*)' Wind (X/Y): ',WNDX,WNDY
IT_FLAG1 = .FALSE.
UST=-999
TAUNUX=0.
TAUNUY=0.
ELSEIF (((VITV/MAX(ABS(WNDY),CK) .LT. iter_thresh) .AND.&
(UITV/MAX(ABS(WNDX),CK) .LT. iter_thresh)) .AND. &
(ITERATION .GE. 2)) THEN
IT_FLAG1 = .FALSE.
ENDIF
! if taunu iteration is unstable try to reset with new guess...
if (.not.(cos(wnd_in_dir-atan2(taunuy,taunux)).ge.0.0)) then
TAUNUX = USTSM**2 * DAIR * wndx / wnd_in_mag*.95
TAUNUY = USTSM**2 * DAIR * wndy / wnd_in_mag*.95
endif
ENDDO
!|---------------------------------------------------------------------|
!|----Finish-----------------------------------------------------------|
!|---------------------------------------------------------------------|
USTD = ATAN2(TAUY,TAUX)
UST = SQRT( SQRT( TAUX**2 + TAUY**2 ) / DAIR)
! Get Z0 from aligned wind
WND_PA=wnd_in_mag*COS(wnd_in_dir-USTD)
Z0 = ZWND/exp(wnd_pa*kappa/ust)
CD = UST**2 / wnd_in_mag**2
IF (PRESENT(CHARN)) THEN
CHARN = 0.01/SQRT(SQRT( TAUNUX**2 + TAUNUY**2 )/(UST**2))
ENDIF
FSFL1=.not.((CD .LT. 0.01).AND.(CD .GT. 0.0001))
FSFL2=.not.(cos(wnd_in_dir-ustd).GT.0.9)
IF (FSFL1 .or. FSFL2) THEN
!Fail safe to bulk
write(*,*)'Attn: W3FLD1 failed, will output bulk...'
CALL wnd2z0m(wnd_in_mag,z0)
UST = wnd_in_mag*kappa/log(zwnd/z0)
USTD = wnd_in_dir
CD = UST**2 / wnd_in_mag**2
ENDIF
DEALLOCATE(WN, DWN, CP,SIG2, SPC2, TLTN, TLTE, TAUD, &
TLTND, TLTED, ZOFK, UPROF, &
VPROF, FTILDE, UP1, VP1, UP, VP, TLTNA, TLTEA)
!/ End of W3FLD1 ----------------------------------------------------- /
!/
RETURN
!
END SUBROUTINE W3FLD1
!/ ------------------------------------------------------------------- /
SUBROUTINE INFLD
!/
!/ +-----------------------------------+
!/ | WAVEWATCH III NOAA/NCEP |
!/ | B. G. Reichl |
!/ | FORTRAN 90 |
!/ | Last update : 15-Jan-2016 |
!/ +-----------------------------------+
!/
!/ 15-Jan-2016 : Origination. ( version 5.12 )
!/
! 1. Purpose :
!
! Initialization for w3fld1 (also used by w3fld2)
!
! 2. Method :
!
! 3. Parameters :
!
! Parameter list
! ----------------------------------------------------------------
! ----------------------------------------------------------------
!
! 4. Subroutines used :
!
! Name Type Module Description
! ----------------------------------------------------------------
! STRACE Subr. W3SERVMD Subroutine tracing.
! ----------------------------------------------------------------
!
! 5. Called by :
!
! Name Type Module Description
! ----------------------------------------------------------------
! W3FLDX Subr. W3FLDXMD Corresponding source term.
! ----------------------------------------------------------------
!
! 6. Error messages :
!
! None.
!
! 7. Remarks :
!
! 8. Structure :
!
! See source code.
!
! 9. Switches :
!
! !/S Enable subroutine tracing.
!
! 10. Source code :
!
!/ ------------------------------------------------------------------- /
USE W3ODATMD, ONLY: NDSE
USE W3GDATMD, ONLY: TAIL_ID, TAIL_LEV, TAIL_TRAN1, TAIL_TRAN2
USE W3SERVMD, ONLY: EXTCDE
!/S USE W3SERVMD, ONLY: STRACE
!/
IMPLICIT NONE
!/
!/ ------------------------------------------------------------------- /
!/ Parameter list
!/
!/
!/ ------------------------------------------------------------------- /
!/ Local parameters
!/
!/S INTEGER, SAVE :: IENT = 0
!/
!/ ------------------------------------------------------------------- /
!/
!/S CALL STRACE (IENT, 'INFLD')
!
! 1. .... ----------------------------------------------------------- *
!
Tail_Choice=Tail_ID
Tail_Level=TAIL_Lev
Tail_transition_ratio1 = TAIL_TRAN1
Tail_transition_ratio2 = TAIL_TRAN2
!
RETURN
!
! Formats
!
!/
!/ End of INFLD1 ----------------------------------------------------- /
!/
END SUBROUTINE INFLD
!/
!/ ------------------------------------------------------------------- /
SUBROUTINE APPENDTAIL(INSPC, WN2, NKT, KA1, KA2, KA3, WNDDIR,SAT)
!/
!/ +-----------------------------------+
!/ | WAVEWATCH III NOAA/NCEP |
!/ | B. G. Reichl |
!/ | FORTRAN 90 |
!/ | Last update : 15-Jan-2016 |
!/ +-----------------------------------+
!/
!/ 15-Jan-2016 : Origination. ( version 5.12 )
!/
! 1. Purpose :
!
! Set tail for stress calculation.
!
! 2. Method :
!
! 3. Parameters :
!
! Parameter list
! ----------------------------------------------------------------
! ----------------------------------------------------------------
!
! 4. Subroutines used :
!
! Name Type Module Description
! ----------------------------------------------------------------
! STRACE Subr. W3SERVMD Subroutine tracing.
! ----------------------------------------------------------------
!
! 5. Called by :
!
! Name Type Module Description
! ----------------------------------------------------------------
! W3FLD1 Subr. W3FLD1MD Corresponding source term.
! ----------------------------------------------------------------
!
! 6. Error messages :
!
! None.
!
! 7. Remarks :
!
! 8. Structure :
!
! See source code.
!
! 9. Switches :
!
! !/S Enable subroutine tracing.
!
! 10. Source code :
!
!/ ------------------------------------------------------------------- /
USE CONSTANTS, ONLY: TPI, PI
USE W3GDATMD, ONLY: NTH, TH, DTH
USE W3ODATMD, ONLY: NDSE
USE W3SERVMD, ONLY: EXTCDE
!/S USE W3SERVMD, ONLY: STRACE
!/
IMPLICIT NONE
!/
!/ ------------------------------------------------------------------- /
!/ Parameter list
!/
INTEGER, INTENT(IN) :: NKT, KA1, KA2, KA3
REAL, INTENT(IN) :: WN2(NKT), WNDDIR,SAT
REAL, INTENT(INOUT) :: INSPC(NKT,NTH)
!/
!/ ------------------------------------------------------------------- /
!/ Local parameters
!/
!/S INTEGER, SAVE :: IENT = 0
REAL :: BT(NKT), IC, ANGLE2, ANG(NKT),&
NORMSPC(NTH), AVG, ANGDIF, M, MAXANG, &
MAXAN, MINAN
INTEGER :: MAI, I, K, T
REAL, ALLOCATABLE, DIMENSION(:) :: ANGLE1
!/
!/ ------------------------------------------------------------------- /
!/
!/S CALL STRACE (IENT, 'APPENDTAIL')
!
! 1. .... ----------------------------------------------------------- *
!
!|###############################################################|
!|##1. Get the level of the saturation spectrum in transition
!|## region A
!|###############################################################|
!-------------------------------------------
! 1a, get saturation level at KA1 (1.25xFPI)
!-------------------------------------------
BT(KA1) = 0
ANG = 0.0
DO T=1, NTH
BT(KA1)=BT(KA1)+INSPC(KA1,T)*WN2(KA1)**3.0*DTH
ENDDO
!-----------------------------------------------
! 1b, Set saturation level at KA2 (3xFPI) to SAT
!-----------------------------------------------
BT(KA2) = SAT
!-------------------------------------------------------------
! 1c, Find slope of saturation spectrum in transition region A
!-------------------------------------------------------------
M = ( BT(KA2) - BT(KA1) ) / ( WN2(KA2) - WN2(KA1) )
!----------------------------------------------------------------
! 1d, Find intercept of saturation spectrum in transition region
! A
!----------------------------------------------------------------
IC = BT(KA1) - M * WN2(KA1)
!------------------------------------------------------
! 1e, Calculate saturation level for all wavenumbers in
! transition region A
!------------------------------------------------------
DO K=KA1,KA2
BT(K)=M*WN2(K)+IC
ENDDO
!|###############################################################|
!|##2. Determine the directionality at each wavenumber in
!|## transition region B
!|###############################################################|
!-----------------------------------------------
! 2a, Find angle of spectral peak at KA2 (3xFPI)
!-----------------------------------------------
MAXANG = 0.0
DO T=1, NTH
IF (INSPC(KA2,T) .GT. MAXANG) THEN
MAXANG=INSPC(KA2,T)
ENDIF
ENDDO
!-------------------------------
! 2b, Check if peak spans 2 bins
!-------------------------------
!MAI = total number of angles of peak (if it spans more than 1)
MAI = 0
DO T=1, NTH
IF (MAXANG .EQ. INSPC(KA2,T)) THEN
MAI = MAI+1
ENDIF
ENDDO
!ANGLE1 = angles that correspond to peak (array)
MAI = MAX(1,MAI)
ALLOCATE(ANGLE1(MAI))
!-----------------------------------------------------
! 2c, If peak spans 2 or more bins it must be averaged
!-----------------------------------------------------