WAVEWATCH III  beta 0.0.1
w3src3md Module Reference

The 'WAM4+' source terms based on P.A.E.M. More...

Functions/Subroutines

subroutine w3spr3 (A, CG, WN, EMEAN, FMEAN, FMEANS, WNMEAN, AMAX, U, UDIR, USTAR, USDIR, TAUWX, TAUWY, CD, Z0, CHARN, LLWS, FMEANWS)
 Calculate mean wave parameters for the use in the source term routines. More...
 
subroutine w3sin3 (A, CG, K, U, USTAR, DRAT, AS, USDIR, Z0, CD, TAUWX, TAUWY, TAUWNX, TAUWNY, ICE, S, D, LLWS, IX, IY)
 Calculate diagonal and input source term for WAM4+ approach. More...
 
subroutine insin3
 Initialization for source term routine. More...
 
subroutine tabu_stress
 To generate the friction velocity table, TAUT(TAUW,U10)=SQRT(TAU). More...
 
subroutine tabu_tauhf (FRMAX)
 Tabulation of the high-frequency wave-supported stress. More...
 
subroutine calc_ustar (WINDSPEED, TAUW, USTAR, Z0, CHARN)
 Compute friction velocity based on wind speed U10. More...
 
subroutine w3sds3 (A, K, CG, EMEAN, FMEAN, WNMEAN, USTAR, USDIR, DEPTH, S, D, IX, IY)
 Calculate whitecapping source term and diagonal term of derivative. More...
 

Variables

integer, parameter itaumax =200
 
integer, parameter jumax =200
 
integer, parameter iustar =100
 
integer, parameter ialpha =200
 
integer, parameter ilevtail =50
 
real, dimension(0:itaumax, 0:jumaxtaut
 
real deltauw
 
real delu
 
real, dimension(0:iustar, 0:ialphatauhft
 
real delust
 
real delalp
 
real, parameter umax = 50.
 
real, parameter tauwmax = 2.2361
 

Detailed Description

The 'WAM4+' source terms based on P.A.E.M.

Janssen's work, with extensions by him and by J.-R. Bidlot.

Converted from the original WAM codes by F. Ardhuin, with further extensions to adapt to a saturation-based breaking and observation-based swell dissipation.

Author
F. Ardhuin
H. L. Tolman
Date
02-Sep-2012

Function/Subroutine Documentation

◆ calc_ustar()

subroutine w3src3md::calc_ustar ( real, intent(in)  WINDSPEED,
real, intent(in)  TAUW,
real, intent(out)  USTAR,
real, intent(out)  Z0,
real, intent(out)  CHARN 
)

Compute friction velocity based on wind speed U10.

Computation of u* based on Quasi-linear theory.

Parameters
[in]WINDSPEED10-m wind speed ... should be NEUTRAL
[in]TAUWWave-supported stress.
[out]USTARFriction velocity.
[out]Z0Air-side roughness length.
[out]CHARN
Author
F. Ardhuin
Date
14-Aug-2006

Definition at line 1136 of file w3src3md.F90.

1136  !/
1137  !/ +-----------------------------------+
1138  !/ | WAVEWATCH III NOAA/NCEP |
1139  !/ | F. Ardhuin |
1140  !/ | FORTRAN 90 |
1141  !/ | Last update 2006/08/14 |
1142  !/ +-----------------------------------+
1143  !/
1144  !/ 27-Feb-2004 : Origination in WW3 ( version 2.22-SHOM )
1145  !/ the resulting table was checked to be identical to the original f77 result
1146  !/ 14-Aug-2006 : Modified following Bidlot ( version 2.22-SHOM )
1147  !/ 18-Aug-2006 : Ported to version 3.09
1148  !/ 03-Apr-2010 : Adding output of Charnock parameter ( version 3.14-IFREMER )
1149  !
1150  ! 1. Purpose :
1151  !
1152  ! Compute friction velocity based on wind speed U10
1153  !
1154  ! 2. Method :
1155  !
1156  ! Computation of u* based on Quasi-linear theory
1157  !
1158  ! 3. Parameters :
1159  !
1160  ! Parameter list
1161  ! ----------------------------------------------------------------
1162  ! U10,TAUW,USTAR,Z0
1163  ! ----------------------------------------------------------------
1164  ! WINDSPEED Real I 10-m wind speed ... should be NEUTRAL
1165  ! TAUW Real I Wave-supported stress
1166  ! USTAR Real O Friction velocity.
1167  ! Z0 Real O air-side roughness length
1168  ! ----------------------------------------------------------------
1169  !
1170  ! 4. Subroutines used :
1171  !
1172  ! STRACE Service routine.
1173  !
1174  ! 5. Called by :
1175  !
1176  ! W3SIN3 Wind input Source term routine.
1177  !
1178  ! 6. Error messages :
1179  !
1180  ! 7. Remarks :
1181  !
1182  ! 8. Structure :
1183  !
1184  ! See source code.
1185  !
1186  ! 9. Switches :
1187  !
1188  ! !/S Enable subroutine tracing.
1189  ! !/T Enable test output.
1190  !
1191  ! 10. Source code :
1192  !-----------------------------------------------------------------------------!
1193  USE constants, ONLY: grav
1194  USE w3gdatmd, ONLY: zzwnd, aalpha
1195  IMPLICIT NONE
1196  REAL, intent(in) :: WINDSPEED,TAUW
1197  REAL, intent(out) :: USTAR, Z0, CHARN
1198  ! local variables
1199  REAL SQRTCDM1
1200  REAL X,XI,DELI1,DELI2,XJ,delj1,delj2
1201  REAL UST,DELTOLD,TAUW_LOCAL
1202  INTEGER IND,J
1203  !
1204  tauw_local=max(min(tauw,tauwmax),0.)
1205  xi = sqrt(tauw_local)/deltauw
1206  ind = min( itaumax-1, int(xi)) ! index for stress table
1207  deli1 = min(1.,xi - real(ind)) !interpolation coefficient for stress table
1208  deli2 = 1. - deli1
1209  xj = windspeed/delu
1210  j = min( jumax-1, int(xj) )
1211  delj1 = min(1.,xj - real(j))
1212  delj2 = 1. - delj1
1213  ustar=(taut(ind,j)*deli2+taut(ind+1,j )*deli1)*delj2 &
1214  + (taut(ind,j+1)*deli2+taut(ind+1,j+1)*deli1)*delj1
1215  !
1216  ! Determines roughness length
1217  !
1218  sqrtcdm1 = min(windspeed/ustar,100.0)
1219  z0 = zzwnd*exp(-kappa*sqrtcdm1)
1220  IF (ustar.GT.0.001) THEN
1221  charn = grav*z0/ustar**2
1222  ELSE
1223  charn = aalpha
1224  END IF
1225  !
1226  RETURN

References w3gdatmd::aalpha, deltauw, delu, constants::grav, itaumax, jumax, taut, tauwmax, and w3gdatmd::zzwnd.

Referenced by w3spr3().

◆ insin3()

subroutine w3src3md::insin3

Initialization for source term routine.

Author
F. Ardhuin
Date
23-Jul-2009

Definition at line 727 of file w3src3md.F90.

727  !/
728  !/ +-----------------------------------+
729  !/ | WAVEWATCH III NOAA/NCEP |
730  !/ | SHOM |
731  !/ | F. Ardhuin |
732  !/ | FORTRAN 90 |
733  !/ | Last update : 23-Jul-2009 |
734  !/ +-----------------------------------+
735  !/
736  !/ 23-Jun-2006 : Origination. ( version 3.09 )
737  !/ 23-Jul-2007 : Cleaning up convolutions ( version 3.14-SHOM)
738  !
739  ! 1. Purpose :
740  !
741  ! Initialization for source term routine.
742  !
743  ! 2. Method :
744  !
745  ! 3. Parameters :
746  !
747  ! Parameter list
748  ! ----------------------------------------------------------------
749  ! ----------------------------------------------------------------
750  !
751  ! 4. Subroutines used :
752  !
753  ! Name Type Module Description
754  ! ----------------------------------------------------------------
755  ! STRACE Subr. W3SERVMD Subroutine tracing.
756  ! ----------------------------------------------------------------
757  !
758  ! 5. Called by :
759  !
760  ! Name Type Module Description
761  ! ----------------------------------------------------------------
762  ! W3SIN3 Subr. W3SRC3MD Corresponding source term.
763  ! ----------------------------------------------------------------
764  !
765  ! 6. Error messages :
766  !
767  ! None.
768  !
769  ! 7. Remarks :
770  !
771  ! 8. Structure :
772  !
773  ! See source code.
774  !
775  ! 9. Switches :
776  !
777  ! !/S Enable subroutine tracing.
778  !
779  ! 10. Source code :
780  !
781  !/ ------------------------------------------------------------------- /
782  USE constants, ONLY: tpiinv
783  USE w3gdatmd, ONLY: sig, nk
784 #ifdef W3_S
785  USE w3servmd, ONLY: strace
786 #endif
787  !/
788  IMPLICIT NONE
789  !/
790  !/ ------------------------------------------------------------------- /
791  !/ Parameter list
792  !/
793  !/ ------------------------------------------------------------------- /
794  !/ Local parameters
795  !/
796 #ifdef W3_S
797  INTEGER, SAVE :: IENT = 0
798 #endif
799  !/
800  !/ ------------------------------------------------------------------- /
801  !/
802 #ifdef W3_S
803  CALL strace (ient, 'INSIN3')
804 #endif
805  !
806  ! 1. .... ----------------------------------------------------------- *
807  !
808  CALL tabu_stress
809  CALL tabu_tauhf(sig(nk) * tpiinv) !tabulate high-frequency stress
810  !/
811  !/ End of INSIN3 ----------------------------------------------------- /
812  !/

References w3gdatmd::nk, w3gdatmd::sig, w3servmd::strace(), tabu_stress(), tabu_tauhf(), and constants::tpiinv.

Referenced by w3iogrmd::w3iogr().

◆ tabu_stress()

subroutine w3src3md::tabu_stress

To generate the friction velocity table, TAUT(TAUW,U10)=SQRT(TAU).

Author
F. Ardhuin
Date
17-Oct-2007

Definition at line 823 of file w3src3md.F90.

823  !/
824  !/ +-----------------------------------+
825  !/ | WAVEWATCH III NOAA/NCEP |
826  !/ | F. Ardhuin |
827  !/ | FORTRAN 90 |
828  !/ | Last update : 17-Oct-2007 |
829  !/ +-----------------------------------+
830  !/
831  !/ 23-Jun-2006 : Origination. ( version 3.13 )
832  !/ adapted from WAM, original:P.A.E.M. JANSSEN KNMI AUGUST 1990
833  !/ adapted version (subr. STRESS): J. BIDLOT ECMWF OCTOBER 2004
834  !/ Table values were checkes against the original f90 result and found to
835  !/ be identical (at least at 0.001 m/s accuracy)
836  !/
837  ! 1. Purpose :
838  ! TO GENERATE friction velocity table TAUT(TAUW,U10)=SQRT(TAU).
839  ! METHOD.
840  ! A STEADY STATE WIND PROFILE IS ASSUMED.
841  ! THE WIND STRESS IS COMPUTED USING THE ROUGHNESSLENGTH
842  ! Z1=Z0/SQRT(1-TAUW/TAU)
843  ! WHERE Z0 IS THE CHARNOCK RELATION , TAUW IS THE WAVE-
844  ! INDUCED STRESS AND TAU IS THE TOTAL STRESS.
845  ! WE SEARCH FOR STEADY-STATE SOLUTIONS FOR WHICH TAUW/TAU < 1.
846  ! FOR QUASILINEAR EFFECT SEE PETER A.E.M. JANSSEN,1990.
847  !
848  ! Initialization for source term routine.
849  !
850  ! 2. Method :
851  !
852  ! 3. Parameters :
853  !
854  ! Parameter list
855  ! ----------------------------------------------------------------
856  ! ----------------------------------------------------------------
857  !
858  ! 4. Subroutines used :
859  !
860  ! Name Type Module Description
861  ! ----------------------------------------------------------------
862  ! STRACE Subr. W3SERVMD Subroutine tracing.
863  ! ----------------------------------------------------------------
864  !
865  ! 5. Called by :
866  !
867  ! Name Type Module Description
868  ! ----------------------------------------------------------------
869  ! W3SIN3 Subr. W3SRC3MD Corresponding source term.
870  ! ----------------------------------------------------------------
871  !
872  ! 6. Error messages :
873  !
874  ! None.
875  !
876  ! 7. Remarks :
877  !
878  ! 8. Structure :
879  !
880  ! See source code.
881  !
882  ! 9. Switches :
883  !
884  ! !/S Enable subroutine tracing.
885  !
886  ! 10. Source code :
887  !
888  !/ ------------------------------------------------------------------- /
889  USE constants, ONLY: grav
890  USE w3gdatmd, ONLY: zzwnd, aalpha, zz0max
891  IMPLICIT NONE
892  INTEGER, PARAMETER :: NITER=10
893  REAL , PARAMETER :: XM=0.50, eps1=0.00001
894  ! VARIABLE. TYPE. PURPOSE.
895  ! *XM* REAL POWER OF TAUW/TAU IN ROUGHNESS LENGTH.
896  ! *XNU* REAL KINEMATIC VISCOSITY OF AIR.
897  ! *NITER* INTEGER NUMBER OF ITERATIONS TO OBTAIN TOTAL STRESS
898  ! *EPS1* REAL SMALL NUMBER TO MAKE SURE THAT A SOLUTION
899  ! IS OBTAINED IN ITERATION WITH TAU>TAUW.
900  ! ----------------------------------------------------------------------
901  INTEGER I,J,ITER
902  REAL ZTAUW,UTOP,CDRAG,WCD,USTOLD,TAUOLD
903  REAL X,UST,ZZ0,ZNU,F,DELF,ZZ00
904  !
905  !
906  delu = umax/float(jumax)
907  deltauw = tauwmax/float(itaumax)
908  DO i=0,itaumax
909  ztauw = (real(i)*deltauw)**2
910  DO j=0,jumax
911  utop = float(j)*delu
912  cdrag = 0.0012875
913  wcd = sqrt(cdrag)
914  ustold = utop*wcd
915  tauold = max(ustold**2, ztauw+eps1)
916  DO iter=1,niter
917  x = ztauw/tauold
918  ust = sqrt(tauold)
919  zz00=aalpha*tauold/grav
920  IF (zz0max.NE.0) zz00=min(zz00,zz0max)
921  ! Corrects roughness ZZ00 for quasi-linear effect
922  zz0 = zz00/(1.-x)**xm
923  !ZNU = 0.1*nu_air/UST ! This was removed by Bidlot in 1996
924  !ZZ0 = MAX(ZNU,ZZ0)
925  f = ust-kappa*utop/(alog(zzwnd/zz0))
926  delf= 1.-kappa*utop/(alog(zzwnd/zz0))**2*2./ust &
927  *(1.-(xm+1)*x)/(1.-x)
928  ust = ust-f/delf
929  tauold= max(ust**2., ztauw+eps1)
930  END DO
931  taut(i,j) = sqrt(tauold)
932  END DO
933  END DO
934  i=itaumax
935  j=jumax
936  !
937  ! Force zero wind to have zero stress (Bidlot 1996)
938  !
939  DO i=0,itaumax
940  taut(i,0)=0.0
941  END DO
942  RETURN

References w3gdatmd::aalpha, deltauw, delu, constants::grav, itaumax, jumax, taut, tauwmax, umax, w3gdatmd::zz0max, and w3gdatmd::zzwnd.

Referenced by insin3().

◆ tabu_tauhf()

subroutine w3src3md::tabu_tauhf ( real, intent(in)  FRMAX)

Tabulation of the high-frequency wave-supported stress.

SEE REFERENCE FOR WAVE STRESS CALCULATION. FOR QUASILINEAR EFFECT SEE PETER A.E.M. JANSSEN,1990. See tech. Memo ECMWF 03 december 2003 by Bidlot & Janssen.

Parameters
[in]FRMAXMaximum frequency.
Author
F. Ardhuin
Date
14-Aug-2006

Definition at line 959 of file w3src3md.F90.

959  !/
960  !/ +-----------------------------------+
961  !/ | WAVEWATCH III NOAA/NCEP |
962  !/ | F. Ardhuin |
963  !/ | FORTRAN 90 |
964  !/ | Last update 2006/08/14 |
965  !/ +-----------------------------------+
966  !/
967  !/ 27-Feb-2004 : Origination in WW3 ( version 2.22.SHOM )
968  !/ the resulting table was checked to be identical to the original f77 result
969  !/ 14-Aug-2006 : Modified following Bidlot ( version 2.22.SHOM )
970  !/ 18-Aug-2006 : Ported to version 3.09
971  !
972  ! 1. Purpose :
973  !
974  ! Tabulation of the high-frequency wave-supported stress
975  !
976  ! 2. Method :
977  !
978  ! SEE REFERENCE FOR WAVE STRESS CALCULATION.
979  ! FOR QUASILINEAR EFFECT SEE PETER A.E.M. JANSSEN,1990.
980  ! See tech. Memo ECMWF 03 december 2003 by Bidlot & Janssen
981  !
982  ! 3. Parameters :
983  !
984  ! Parameter list
985  ! ----------------------------------------------------------------
986  ! FRMAX Real I maximum frequency.
987  ! ----------------------------------------------------------------
988  !
989  ! 4. Subroutines used :
990  !
991  ! STRACE Service routine.
992  !
993  ! 5. Called by :
994  !
995  ! W3SIN3 Wind input Source term routine.
996  !
997  ! 6. Error messages :
998  !
999  ! 7. Remarks :
1000  !
1001  ! 8. Structure :
1002  !
1003  ! See source code.
1004  !
1005  ! 9. Switches :
1006  !
1007  ! !/S Enable subroutine tracing.
1008  ! !/T Enable test output.
1009  !
1010  ! 10. Source code :
1011  !
1012  !/ ------------------------------------------------------------------- /
1013  USE constants, ONLY: grav, tpi
1014  USE w3gdatmd, ONLY: aalpha, bbeta, zzalp, xfr, fachfe, zz0max
1015 #ifdef W3_T
1016  USE w3odatmd, ONLY: ndst
1017 #endif
1018 #ifdef W3_S
1019  USE w3servmd, ONLY: strace
1020 #endif
1021  !
1022  IMPLICIT NONE
1023  !/
1024  !/ ------------------------------------------------------------------- /
1025  !/ Parameter list
1026  !/
1027  REAL, intent(in) :: FRMAX ! maximum frequency
1028  !/
1029  !/ ------------------------------------------------------------------- /
1030  !/ Local parameters
1031  !/
1032  ! USTARM R.A. Maximum friction velocity
1033  ! ALPHAM R.A. Maximum Charnock Coefficient
1034  ! WLV R.A. Water levels.
1035  ! UA R.A. Absolute wind speeds.
1036  ! UD R.A. Absolute wind direction.
1037  ! U10 R.A. Wind speed used.
1038  ! U10D R.A. Wind direction used.
1039  ! 10. Source code :
1040  !
1041  !/ ------------------------------------------------------------------- /
1042  REAL :: USTARM, ALPHAM
1043  REAL :: CONST1, OMEGA, OMEGAC
1044  REAL :: UST, ZZ0,OMEGACC, CM
1045  INTEGER, PARAMETER :: JTOT=250
1046  REAL, ALLOCATABLE :: W(:)
1047  REAL :: ZX,ZARG,ZMU,ZLOG,ZZ00,ZBETA
1048  REAL :: Y,YC,DELY
1049  INTEGER :: I,J,K,L
1050  REAL :: X0
1051 #ifdef W3_S
1052  INTEGER, SAVE :: IENT = 0
1053  CALL strace (ient, 'TABU_HF')
1054 #endif
1055  !
1056  ustarm = 5.
1057  alpham = 20.*aalpha
1058  delust = ustarm/real(iustar)
1059  delalp = alpham/real(ialpha)
1060  const1 = bbeta/kappa**2
1061  omegac = tpi*frmax
1062  !
1063  tauhft(0:iustar,0:ialpha)=0. !table initialization
1064  !
1065  ALLOCATE(w(jtot))
1066  w(2:jtot-1)=1.
1067  w(1)=0.5
1068  w(jtot)=0.5
1069  x0 = 0.05
1070  !
1071  DO l=0,ialpha
1072  DO k=0,iustar
1073  ust = max(real(k)*delust,0.000001)
1074  zz00 = ust**2*aalpha/grav
1075  IF (zz0max.NE.0) zz00=min(zz00,zz0max)
1076  zz0 = zz00*(1+float(l)*delalp/aalpha)
1077  omegacc = max(omegac,x0*grav/ust)
1078  yc = omegacc*sqrt(zz0/grav)
1079  dely = max((1.-yc)/real(jtot),0.)
1080  ! For a given value of UST and ALPHA,
1081  ! the wave-supported stress is integrated all the way
1082  ! to 0.05*g/UST
1083  DO j=1,jtot
1084  y = yc+real(j-1)*dely
1085  omega = y*sqrt(grav/zz0)
1086  ! This is the deep water phase speed
1087  cm = grav/omega
1088  !this is the inverse wave age, shifted by ZZALP (tuning)
1089  zx = ust/cm +zzalp
1090  zarg = min(kappa/zx,20.)
1091  zmu = min(grav*zz0/cm**2*exp(zarg),1.)
1092  zlog = min(alog(zmu),0.)
1093  zbeta = const1*zmu*zlog**4
1094  ! Power of Y in denominator should be FACHFE-4
1095  tauhft(k,l) = tauhft(k,l)+w(j)*zbeta/y*dely
1096  END DO
1097  !IF (MOD(K,5).EQ.0.AND.MOD(L,5).EQ.0) &
1098  !WRITE(102,'(2I4,3G16.8)') L,K,UST,AALPHA+FLOAT(L)*DELALP,TAUHFT(K,L)
1099 #ifdef W3_T
1100  WRITE (ndst,9000) l,k,aalpha+float(l)*delalp,ust,tauhft(k,l)
1101 #endif
1102  END DO
1103  END DO
1104  DEALLOCATE(w)
1105  ! DO L=0,IALPHA
1106  ! DO K=0,IUSTAR
1107  !WRITE(101,'(A,2I4,G16.8)') 'L,K,TAUHFT(K,L):',L,K,TAUHFT(K,L)
1108  ! END DO
1109  ! END DO
1110  !WRITE(101,*) 'TAUHFT:',FRMAX,BBETA,AALPHA,CONST1,OMEGAC,TPI
1111  !WRITE(101,'(20G16.8)') TAUHFT
1112  RETURN
1113 #ifdef W3_T
1114 9000 FORMAT (' TABU_HF, L, K :',(2i4,3f8.3)/)
1115 #endif

References w3gdatmd::aalpha, w3gdatmd::bbeta, delalp, delust, w3gdatmd::fachfe, constants::grav, ialpha, iustar, w3odatmd::ndst, w3servmd::strace(), tauhft, constants::tpi, w3gdatmd::xfr, w3gdatmd::zz0max, and w3gdatmd::zzalp.

Referenced by insin3().

◆ w3sds3()

subroutine w3src3md::w3sds3 ( real, dimension(nspec), intent(in)  A,
real, dimension(nk), intent(in)  K,
real, dimension(nk), intent(in)  CG,
real, intent(in)  EMEAN,
real, intent(in)  FMEAN,
real, intent(in)  WNMEAN,
real, intent(in)  USTAR,
real, intent(in)  USDIR,
real, intent(in)  DEPTH,
real, dimension(nspec), intent(out)  S,
real, dimension(nspec), intent(out)  D,
integer, intent(in)  IX,
integer, intent(in)  IY 
)

Calculate whitecapping source term and diagonal term of derivative.

WAM-Cycle 4 and following. The last update (09-May-2005) follows the redefinition of the mean wavenumber as in Bidlot et al. (2005).

Parameters
[in]AAction density spectrum (1-D).
[in]KWavenumber for entire spectrum.
[in]CG
[in]EMEANMean wave energy.
[in]FMEANMean wave frequency.
[in]WNMEANMean wavenumber.
[in]USTARFriction velocity.
[in]USDIRWind stress direction.
[in]DEPTHWater depth.
[out]SSource term (1-D version).
[out]DDiagonal term of derivative.
[in]IX
[in]IY
Author
F. Ardhuin
Date
23-Jul-2009

Definition at line 1255 of file w3src3md.F90.

1255  !/
1256  !/ +-----------------------------------+
1257  !/ | WAVEWATCH III NOAA/NCEP |
1258  !/ ! F. Ardhuin !
1259  !/ | FORTRAN 90 |
1260  !/ | Last update : 23-Jul-2009 |
1261  !/ +-----------------------------------+
1262  !/
1263  !/ 05-Dec-1996 : Final FORTRAN 77 ( version 1.18 )
1264  !/ 08-Dec-1999 : Upgrade to FORTRAN 90 ( version 2.00 )
1265  !/ 14-Aug-2006 : Generic WAM4+ dissipation term ( version 2.22SHOM )
1266  !/ 23-Jul-2009 : Addition of Filipot &al convolution ( version 3.14-SHOM )
1267  !/
1268  ! 1. Purpose :
1269  !
1270  ! Calculate whitecapping source term and diagonal term of derivative.
1271  !
1272  ! 2. Method :
1273  !
1274  ! WAM-Cycle 4 and following.
1275  ! The last update (09-May-2005) follows the redefinition of
1276  ! the mean wavenumber as in Bidlot et al. (2005).
1277  !
1278  ! 3. Parameters :
1279  !
1280  ! Parameter list
1281  ! ----------------------------------------------------------------
1282  ! A R.A. I Action density spectrum (1-D).
1283  ! K R.A. I Wavenumber for entire spectrum. *)
1284  ! EMEAN Real I Mean wave energy.
1285  ! FMEAN Real I Mean wave frequency.
1286  ! WNMEAN Real I Mean wavenumber.
1287  ! USTAR Real I Friction velocity.
1288  ! USDIR Real I wind stress direction.
1289  ! DEPTH Real I Water depth.
1290  ! S R.A. O Source term (1-D version).
1291  ! D R.A. O Diagonal term of derivative. *)
1292  ! ----------------------------------------------------------------
1293  ! *) Stored in 1-D array with dimension NTH*NK
1294  !
1295  ! 4. Subroutines used :
1296  !
1297  ! STRACE Subroutine tracing. ( !/S switch )
1298  ! PRT2DS Print plot of spectrum. ( !/T0 switch )
1299  ! OUTMAT Print out matrix. ( !/T1 switch )
1300  !
1301  ! 5. Called by :
1302  !
1303  ! W3SRCE Source term integration.
1304  ! W3EXPO Point output program.
1305  ! GXEXPO GrADS point output program.
1306  !
1307  ! 6. Error messages :
1308  !
1309  ! 7. Remarks :
1310  !
1311  ! 8. Structure :
1312  !
1313  ! See source code.
1314  !
1315  ! 9. Switches :
1316  !
1317  ! !/S Enable subroutine tracing.
1318  ! !/T Enable general test output.
1319  ! !/T0 2-D print plot of source term.
1320  ! !/T1 Print arrays.
1321  !
1322  ! 10. Source code :
1323  !
1324  !/ ------------------------------------------------------------------- /
1325  USE constants, ONLY: grav, tpi
1326  USE w3gdatmd, ONLY: nspec, nth, nk, ddelta1, ddelta2, &
1327 #ifdef w3_t0
1328  sig, &
1329 #endif
1330  ssdsc1
1331 #ifdef W3_S
1332  USE w3servmd, ONLY: strace
1333 #endif
1334 #ifdef W3_T
1335  USE w3odatmd, ONLY: ndst
1336 #endif
1337 #ifdef W3_T1
1338  USE w3odatmd, ONLY: ndst
1339 #endif
1340 #ifdef W3_T0
1341  USE w3arrymd, ONLY: prt2ds
1342 #endif
1343 #ifdef W3_T1
1344  USE w3arrymd, ONLY: outmat
1345 #endif
1346  !
1347  IMPLICIT NONE
1348  !/
1349  !/ ------------------------------------------------------------------- /
1350  !/ Parameter list
1351  !/
1352  REAL, INTENT(IN) :: A(NSPEC), K(NK), CG(NK), &
1353  DEPTH, USTAR, USDIR, EMEAN, FMEAN, WNMEAN
1354  INTEGER, INTENT(IN) :: IX, IY
1355  REAL, INTENT(OUT) :: S(NSPEC), D(NSPEC)
1356  !/
1357  !/ ------------------------------------------------------------------- /
1358  !/ Local parameters
1359  !/
1360  INTEGER :: IS, IK, ITH
1361 #ifdef W3_S
1362  INTEGER, SAVE :: IENT = 0
1363 #endif
1364  REAL :: FACTOR, FACTOR2
1365  REAL :: ALFAMEAN, WNMEAN2
1366 #ifdef W3_T0
1367  REAL :: DOUT(NK,NTH)
1368 #endif
1369  !/
1370  !/ ------------------------------------------------------------------- /
1371  !/
1372 #ifdef W3_S
1373  CALL strace (ient, 'W3SDS3')
1374 #endif
1375  !
1376  ! 0. Pre-initialization of arrays, should be set before being used
1377  ! but this is helping with bit reproducibility
1378  d=0.
1379 
1380  ! 1. Common factor
1381  !
1382  wnmean2 = max( 1.e-10 , wnmean )
1383  alfamean=wnmean**2*emean
1384  factor = ssdsc1 * tpi*fmean * alfamean**2
1385  !
1386 #ifdef W3_T
1387  WRITE (ndst,9000) ssdsc1, fmean, wnmean, emean, factor
1388 #endif
1389  !
1390  !----------------------------------------------------------------------
1391  !
1392  ! 2. Source term
1393  !
1394  DO ik=1, nk
1395  !
1396  ! Original WAM4/WAM4+ dissipation term
1397  !
1398  factor2=factor*(ddelta1*k(ik)/wnmean2 + ddelta2*(k(ik)/wnmean2)**2)
1399  DO ith=1,nth
1400  is=ith+(ik-1)*nth
1401  d(is)= factor2
1402  END DO
1403  END DO
1404  !
1405  s = d * a
1406  !
1407  ! ... Test output of arrays
1408  !
1409 #ifdef W3_T0
1410  DO ik=1, nk
1411  DO ith=1, nth
1412  dout(ik,ith) = d(ith+(ik-1)*nth)
1413  END DO
1414  END DO
1415  CALL prt2ds (ndst, nk, nk, nth, dout, sig(1:nk), ' ', 1., &
1416  0.0, 0.001, 'Diag Sds', ' ', 'NONAME')
1417 #endif
1418  !
1419 #ifdef W3_T1
1420  CALL outmat (ndst, d, nth, nth, nk, 'diag Sds')
1421 #endif
1422  !
1423  RETURN
1424  !
1425  ! Formats
1426  !
1427 #ifdef W3_T
1428 9000 FORMAT (' TEST W3SDS3 : COMMON FACT.: ',5e10.3)
1429 #endif
1430  !/
1431  !/ End of W3SDS3 ----------------------------------------------------- /
1432  !/

References w3gdatmd::ddelta1, w3gdatmd::ddelta2, constants::grav, w3odatmd::ndst, w3gdatmd::nk, w3gdatmd::nspec, w3gdatmd::nth, w3arrymd::outmat(), w3arrymd::prt2ds(), w3servmd::strace(), and constants::tpi.

Referenced by gxexpo(), w3exnc(), w3expo(), and w3srcemd::w3srce().

◆ w3sin3()

subroutine w3src3md::w3sin3 ( real, dimension(nspec), intent(in)  A,
real, dimension(nk), intent(in)  CG,
real, dimension(nspec), intent(in)  K,
real, intent(in)  U,
real, intent(in)  USTAR,
real, intent(in)  DRAT,
real, intent(in)  AS,
real, intent(in)  USDIR,
real, intent(in)  Z0,
real, intent(in)  CD,
real, intent(out)  TAUWX,
real, intent(out)  TAUWY,
real, intent(out)  TAUWNX,
real, intent(out)  TAUWNY,
real, intent(in)  ICE,
real, dimension(nspec), intent(out)  S,
real, dimension(nspec), intent(out)  D,
logical, dimension(nspec), intent(out)  LLWS,
integer, intent(in)  IX,
integer, intent(in)  IY 
)

Calculate diagonal and input source term for WAM4+ approach.

       WAM-4 : Janssen et al.
       WAM-"4.5" : gustiness effect (Cavaleri et al. )
       SWELLF: damping coefficient (=1) for Janssen (2004) theory
Parameters
[in]AAction density spectrum (1-D).
[in]CGGroup speed.
[in]KWavenumber for entire spectrum.
[in]UWind speed.
[in]USTARFriction velocity.
[in]DRATAir/water density ratio.
[in]ASAir-sea temperature difference.
[in]USDIRWind stress direction.
[in]Z0Air-sea roughness length.
[in]CDWind drag coefficient.
[out]TAUWXComponent of the wave-supported stress.
[out]TAUWYComponent of the wave-supported stress.
[out]TAUWNXComponent of the negative wave-supported stress.
[out]TAUWNYComponent of the negative wave-supported stress.
[in]ICESea ice fraction.
[out]SSource term (1-D version).
[out]DDiagonal term of derivative.
[out]LLWSWind sea true/false array for each component.
[in]IX
[in]IY
Author
F. Ardhuin
H. L. Tolman
Date
02-Sep-2012

Definition at line 386 of file w3src3md.F90.

386  !/
387  !/ +-----------------------------------+
388  !/ | WAVEWATCH III SHOM |
389  !/ ! F. Ardhuin !
390  !/ | H. L. Tolman |
391  !/ | FORTRAN 90 |
392  !/ | Last update : 02-Sep-2012 |
393  !/ +-----------------------------------+
394  !/
395  !/ 09-Oct-2007 : Origination. ( version 3.13 )
396  !/ 16-May-2010 : Adding sea ice ( version 3.14_Ifremer )
397  !/ 02-Sep-2012 : Clean up test output T, T0, T1 ( version 4.07 )
398  !/
399  ! 1. Purpose :
400  !
401  ! Calculate diagonal and input source term for WAM4+ approach.
402  !
403  ! 2. Method :
404  !
405  ! WAM-4 : Janssen et al.
406  ! WAM-"4.5" : gustiness effect (Cavaleri et al. )
407  ! SWELLF: damping coefficient (=1) for Janssen (2004) theory
408  !
409  ! 3. Parameters :
410  !
411  ! Parameter list
412  ! ----------------------------------------------------------------
413  ! A R.A. I Action density spectrum (1-D).
414  ! CG R.A. I Group speed *)
415  ! K R.A. I Wavenumber for entire spectrum. *)
416  ! U Real I WIND SPEED
417  ! USTAR Real I Friction velocity.
418  ! DRAT Real I Air/water density ratio.
419  ! AS Real I Air-sea temperature difference
420  ! USDIR Real I wind stress direction
421  ! Z0 Real I Air-side roughness lengh.
422  ! CD Real I Wind drag coefficient.
423  ! USDIR Real I Direction of friction velocity
424  ! TAUWX-Y Real I Components of the wave-supported stress.
425  ! TAUWNX Real I Component of the negative wave-supported stress.
426  ! TAUWNY Real I Component of the negative wave-supported stress.
427  ! ICE Real I Sea ice fraction. !/Stefan: ICE is DUMMY argument; remove later.
428  ! S R.A. O Source term (1-D version).
429  ! D R.A. O Diagonal term of derivative. *)
430  ! LLWS L.A. O Wind sea true/false array for each component
431  ! ----------------------------------------------------------------
432  ! *) Stored as 1-D array with dimension NTH*NK
433  !
434  ! 4. Subroutines used :
435  !
436  ! STRACE Subroutine tracing. ( !/S switch )
437  ! PRT2DS Print plot of spectrum. ( !/T0 switch )
438  ! OUTMAT Print out matrix. ( !/T1 switch )
439  !
440  ! 5. Called by :
441  !
442  ! W3SRCE Source term integration.
443  ! W3EXPO Point output program.
444  ! GXEXPO GrADS point output program.
445  !
446  ! 6. Error messages :
447  !
448  ! 7. Remarks :
449  !
450  ! 8. Structure :
451  !
452  ! See source code.
453  !
454  ! 9. Switches :
455  !
456  ! !/S Enable subroutine tracing.
457  ! !/T Enable general test output.
458  ! !/T0 2-D print plot of source term.
459  ! !/T1 Print arrays.
460  !
461  ! 10. Source code :
462  !
463  !/ ------------------------------------------------------------------- /
464  USE constants, ONLY: grav, tpi
465 #ifdef W3_T
466  USE constants, ONLY: rade
467 #endif
468  USE w3gdatmd, ONLY: nk, nth, nspec, xfr, dden, sig, sig2, th, &
470  sswellf, &
472 #ifdef W3_S
473  USE w3servmd, ONLY: strace
474 #endif
475 #ifdef W3_T
476  USE w3odatmd, ONLY: ndst
477 #endif
478 #ifdef W3_T0
479  USE w3odatmd, ONLY: ndst
480 #endif
481 #ifdef W3_T1
482  USE w3odatmd, ONLY: ndst
483 #endif
484 #ifdef W3_T0
485  USE w3arrymd, ONLY: prt2ds
486 #endif
487 #ifdef W3_T1
488  USE w3arrymd, ONLY: outmat
489 #endif
490  !
491  IMPLICIT NONE
492  !/
493  !/ ------------------------------------------------------------------- /
494  !/ Parameter list
495  !/
496  REAL, INTENT(IN) :: A(NSPEC)
497  REAL, INTENT(IN) :: CG(NK), K(NSPEC),Z0,U, CD
498  REAL, INTENT(IN) :: USTAR, USDIR, AS, DRAT, ICE !/Stefan: ICE is DUMMY argument; remove later.
499  REAL, INTENT(OUT) :: S(NSPEC), D(NSPEC), TAUWX, TAUWY, TAUWNX, TAUWNY
500  LOGICAL, INTENT(OUT) :: LLWS(NSPEC)
501  INTEGER, INTENT(IN) :: IX, IY
502  !/
503  !/ ------------------------------------------------------------------- /
504  !/ Local parameters
505  !/
506  INTEGER :: IS,IK,ITH
507 #ifdef W3_S
508  INTEGER, SAVE :: IENT = 0
509 #endif
510  REAL :: COSU, SINU, TAUX, TAUY
511  REAL :: UST2, TAUW, TAUWB
512  REAL , PARAMETER :: EPS1 = 0.00001, eps2 = 0.000001
513 #ifdef W3_STAB3
514  REAL :: Usigma !standard deviation of U due to gustiness
515  REAL :: USTARsigma !standard deviation of USTAR due to gustiness
516 #endif
517  REAL :: CM,UCO,UCN,ZCN, &
518  Z0VISC
519  REAL XI,DELI1,DELI2
520  REAL XJ,DELJ1,DELJ2
521  REAL :: CONST, CONST0, CONST2, CONST3, TAU1
522  REAL :: X,ZARG,ZLOG,ZBETA,UST
523  REAL COSWIND,XSTRESS,YSTRESS,TAUHF
524  REAL TEMP, TEMP2
525  INTEGER IND,J,ISTAB
526  REAL DSTAB(3,NSPEC)
527  REAL STRESSSTAB(3,2),STRESSSTABN(3,2)
528 #ifdef W3_T0
529  REAL :: DOUT(NK,NTH)
530 #endif
531  !/
532  !/ ------------------------------------------------------------------- /
533  !/
534 #ifdef W3_S
535  CALL strace (ient, 'W3SIN3')
536 #endif
537  !
538 #ifdef W3_T
539  WRITE (ndst,9000) bbeta, ustar, usdir*rade
540 #endif
541  !
542  ! 1. Preparations
543  !
544  ! JDM: initializing arrays (shouldn't change answers)
545  dstab=0.; stressstab=0. ;stressstabn=0.
546  !
547  ! 1.a estimation of surface roughness parameters
548  !
549  z0visc = 0.1*nu_air/max(ustar,0.0001)
550  !
551  ! 2. Diagonal
552  !
553  ! Here AS is the air-sea temperature difference in degrees. Expression given by
554  ! Abdalla & Cavaleri, JGR 2002 for Usigma. For USTARsigma ... I do not see where
555  ! I got it from, maybe just made up from drag law ...
556  !
557 #ifdef W3_STAB3
558  usigma=max(0.,-0.025*as)
559  ustarsigma=(1.0+u/(10.+u))*usigma
560 #endif
561  ust=ustar
562  istab=3
563 #ifdef W3_STAB3
564  DO istab=1,2
565  IF (istab.EQ.1) ust=ustar*(1.-ustarsigma)
566  IF (istab.EQ.2) ust=ustar*(1.+ustarsigma)
567 #endif
568  taux = ust**2* cos(usdir)
569  tauy = ust**2* sin(usdir)
570  !
571  ! Loop over the resolved part of the spectrum
572  !
573  stressstab(istab,:)=0.
574  stressstabn(istab,:)=0.
575  const0=bbeta*drat/(kappa**2)
576  ! SSWELLF(1) is IDAMPING in ECMWAM
577  const3 = sswellf(1)*2.*kappa*drat
578  !
579  cosu = cos(usdir)
580  sinu = sin(usdir)
581  DO ik=1, nk
582  is=1+(ik-1)*nth
583  cm=k(is)/sig2(is) !inverse of phase speed
584  ucn=ust*cm+zzalp !this is the inverse wave age
585  !
586  ! the stress is the real stress (N/m^2) divided by
587  ! rho_a, and thus comparable to USTAR**2
588  ! it is the integral of rho_w g Sin/C /rho_a
589  ! (air-> waves momentum flux)
590  !
591  const2=dden2(is)/cg(ik) & !Jacobian to get energy in band
592  *grav/(sig(ik)/k(is)*drat) ! coefficient to get momentum
593  const=sig2(is)*const0
594  ! this CM parameter is 1 / C_phi
595  ! this is the "correct" shallow-water expression
596  ! here Z0 corresponds to Z0+Z1 of the Janssen eq. 14
597  zcn=alog(k(is)*z0)
598  ! commented below is the original WAM version (OK for deep water)
599  ! ZCN=ALOG(G*Z0b(I)*CM(I)**2)
600  DO ith=1,nth
601  is=ith+(ik-1)*nth
602  coswind=(ecos(is)*cosu+esin(is)*sinu)
603  IF (coswind.GT.0.01) THEN
604  x=coswind*ucn
605  ! this ZARG term is the argument of the exponential
606  ! in Janssen 1991 eq. 16.
607  zarg=kappa/x
608  ! ZLOG is ALOG(MU) where MU is defined by Janssen 1991 eq. 15
609  zlog=zcn+zarg
610  zbeta = const3*(coswind+kappa/(zcn*ust*cm))*ucn**2
611  IF (zlog.LT.0.) THEN
612  ! The source term Sp is beta * omega * X**2
613  ! as given by Janssen 1991 eq. 19
614  dstab(istab,is) = const*exp(zlog)*zlog**4*ucn**2*coswind**ssinthp
615  llws(is)=.true.
616  ELSE
617  dstab(istab,is) = zbeta
618  llws(is)=.true.
619  END IF
620  ELSE
621  zbeta = const3*(coswind+kappa/(zcn*ust*cm))*ucn**2
622  dstab(istab,is) = zbeta
623  llws(is)=.false.
624  END IF
625  !
626  temp2=const2*dstab(istab,is)*a(is)
627  IF (dstab(istab,is).LT.0) THEN
628  stressstabn(istab,1)=stressstabn(istab,1)+temp2*ecos(is)
629  stressstabn(istab,2)=stressstabn(istab,2)+temp2*esin(is)
630  END IF
631  stressstab(istab,1)=stressstab(istab,1)+temp2*ecos(is)
632  stressstab(istab,2)=stressstab(istab,2)+temp2*esin(is)
633  END DO
634  END DO
635  !
636  d(:)=dstab(3,:)
637  xstress=stressstab(3,1)
638  ystress=stressstab(3,2)
639  tauwnx =stressstabn(3,1)
640  tauwny =stressstabn(3,2)
641  ! WRITE(995,'(A,11G14.5)') 'NEGSTRESS: ',TAUWNX,TAUWNY,FW*UORB**3
642 #ifdef W3_STAB3
643  END DO
644  d(:)=0.5*(dstab(1,:)+dstab(2,:))
645  xstress=0.5*(stressstab(1,1)+stressstab(2,1))
646  ystress=0.5*(stressstab(1,2)+stressstab(2,2))
647  tauwnx=0.5*(stressstabn(1,1)+stressstabn(2,1))
648  tauwny=0.5*(stressstabn(1,2)+stressstabn(2,2))
649 #endif
650  s = d * a
651  !
652  ! ... Test output of arrays
653  !
654 #ifdef W3_T0
655  DO ik=1, nk
656  DO ith=1, nth
657  dout(ik,ith) = d(ith+(ik-1)*nth)
658  END DO
659  END DO
660  CALL prt2ds (ndst, nk, nk, nth, dout, sig(1:nk), ' ', 1., &
661  0.0, 0.001, 'Diag Sin', ' ', 'NONAME')
662 #endif
663  !
664 #ifdef W3_T1
665  CALL outmat (ndst, d, nth, nth, nk, 'diag Sin')
666 #endif
667  !
668  ! Computes the high-frequency contribution
669  ! the difference in spectal density (kx,ky) to (f,theta)
670  ! is integrated in this modified CONST0
671  const0=dth*sig(nk)**5/((grav**2)*tpi) &
672  *tpi*sig(nk) / cg(nk) !conversion WAM (E(f,theta) to WW3 A(k,theta)
673  temp=0.
674  DO ith=1,nth
675  is=ith+(nk-1)*nth
676  coswind=(ecos(is)*cosu+esin(is)*sinu)
677  temp=temp+a(is)*(max(coswind,0.))**3
678  END DO
679  !
680  ! finds the values in the tabulated stress TAUHFT
681  !
682  xi=ust/delust
683  ind = max(1,min(iustar-1, int(xi)))
684  deli1= max(min(1. ,xi-float(ind)),0.)
685  deli2= 1. - deli1
686  xj=max(0.,(grav*z0/max(ust,0.00001)**2-aalpha) / delalp)
687  j = max(1 ,min(ialpha-1, int(xj)))
688  delj1= max(0.,min(1. , xj-float(j)))
689  delj2=1. - delj1
690  tau1 =(tauhft(ind,j)*deli2+tauhft(ind+1,j)*deli1 )*delj2 &
691  +(tauhft(ind,j+1)*deli2+tauhft(ind+1,j+1)*deli1)*delj1
692  tauhf = const0*temp*ust**2*tau1
693  tauwx = xstress+tauhf*cos(usdir)
694  tauwy = ystress+tauhf*sin(usdir)
695  !
696  ! Reduces tail effect to make sure that wave-supported stress
697  ! is less than total stress, this is borrowed from ECWAM Stresso.F
698  !
699  tauw = sqrt(tauwx**2+tauwy**2)
700  ust2 = max(ustar,eps2)**2
701  tauwb = min(tauw,max(ust2-eps1,eps2**2))
702  IF (tauwb.LT.tauw) THEN
703  tauwx=tauwx*tauwb/tauw
704  tauwy=tauwy*tauwb/tauw
705  END IF
706 
707  RETURN
708  !
709  ! Formats
710  !
711 #ifdef W3_T
712 9000 FORMAT (' TEST W3SIN3 : COMMON FACT.: ',3e10.3)
713 #endif
714  !/
715  !/ End of W3SIN3 ----------------------------------------------------- /
716  !/

References w3gdatmd::aalpha, w3gdatmd::bbeta, w3gdatmd::dden, w3gdatmd::dden2, delalp, delust, w3gdatmd::dth, w3gdatmd::ec2, w3gdatmd::ecos, w3gdatmd::esin, constants::grav, ialpha, iustar, w3odatmd::ndst, w3gdatmd::nk, w3gdatmd::nspec, w3gdatmd::nth, w3arrymd::outmat(), w3arrymd::prt2ds(), constants::rade, w3gdatmd::sig, w3gdatmd::sig2, w3gdatmd::ssinthp, w3gdatmd::sswellf, w3servmd::strace(), tauhft, w3gdatmd::th, constants::tpi, w3gdatmd::xfr, w3gdatmd::zz0rat, w3gdatmd::zzalp, and w3gdatmd::zzwnd.

Referenced by gxexpo(), w3exnc(), w3expo(), and w3srcemd::w3srce().

◆ w3spr3()

subroutine w3src3md::w3spr3 ( real, dimension(nth,nk), intent(in)  A,
real, dimension(nk), intent(in)  CG,
real, dimension(nk), intent(in)  WN,
real, intent(out)  EMEAN,
real, intent(out)  FMEAN,
real, intent(out)  FMEANS,
real, intent(out)  WNMEAN,
real, intent(out)  AMAX,
real, intent(in)  U,
real, intent(in)  UDIR,
real, intent(inout)  USTAR,
real, intent(inout)  USDIR,
real, intent(in)  TAUWX,
real, intent(in)  TAUWY,
real, intent(out)  CD,
real, intent(out)  Z0,
real, intent(out)  CHARN,
logical, dimension(nspec), intent(in)  LLWS,
real, intent(out)  FMEANWS 
)

Calculate mean wave parameters for the use in the source term routines.

Parameters
[in]AAction density spectrum.
[in]CGGroup velocities.
[in]WNWavenumbers.
[out]EMEANEnergy.
[out]FMEANMean frequency for determination of tail.
[out]FMEANSMean frequency for dissipation source term.
[out]WNMEANMean wavenumber.
[out]AMAXMaximum of action spectrum.
[in]UWind speed.
[in]UDIRWind direction.
[in,out]USTARFriction velocity
[in,out]USDIRWind stress direction.
[in]TAUWXComponent of wave-supported stress.
[in]TAUWYComponent of wave-supported stress.
[out]CDDrag coefficient at wind level ZWND.
[out]Z0Corresponding z0.
[out]CHARN
[in]LLWSWind sea true/false array for each component.
[out]FMEANWSMean frequency of wind sea, used for tail.
Author
F. Ardhuin
H. L. Tolman
Date
17-Oct-2007

Definition at line 137 of file w3src3md.F90.

137  !/
138  !/ +-----------------------------------+
139  !/ | WAVEWATCH III SHOM |
140  !/ ! F. Ardhuin !
141  !/ | H. L. Tolman |
142  !/ | FORTRAN 90 |
143  !/ | Last update : 17-Oct-2007 |
144  !/ +-----------------------------------+
145  !/
146  !/ 03-Oct-2007 : Origination. ( version 3.13 )
147  !/
148  ! 1. Purpose :
149  !
150  ! Calculate mean wave parameters for the use in the source term
151  ! routines.
152  !
153  ! 2. Method :
154  !
155  ! See source term routines.
156  !
157  ! 3. Parameters :
158  !
159  ! Parameter list
160  ! ----------------------------------------------------------------
161  ! A R.A. I Action density spectrum.
162  ! CG R.A. I Group velocities.
163  ! WN R.A. I Wavenumbers.
164  ! EMEAN Real O Energy
165  ! FMEAN Real O Mean frequency for determination of tail
166  ! FMEANS Real O Mean frequency for dissipation source term
167  ! WNMEAN Real O Mean wavenumber.
168  ! AMAX Real O Maximum of action spectrum.
169  ! U Real I Wind speed.
170  ! UDIR Real I Wind direction.
171  ! USTAR Real I/O Friction velocity.
172  ! USDIR Real I/O wind stress direction.
173  ! TAUWX-Y Real I Components of wave-supported stress.
174  ! CD Real O Drag coefficient at wind level ZWND.
175  ! Z0 Real O Corresponding z0.
176  ! LLWS L.A. I Wind sea true/false array for each component
177  ! FMEANWS Real O Mean frequency of wind sea, used for tail
178  ! ----------------------------------------------------------------
179  !
180  ! 4. Subroutines used :
181  !
182  ! STRACE Service routine.
183  !
184  ! 5. Called by :
185  !
186  ! W3SRCE Source term integration routine.
187  ! W3OUTP Point output program.
188  ! GXEXPO GrADS point output program.
189  !
190  ! 6. Error messages :
191  !
192  ! 7. Remarks :
193  !
194  ! 8. Structure :
195  !
196  ! See source code.
197  !
198  ! 9. Switches :
199  !
200  ! !/S Enable subroutine tracing.
201  ! !/T Enable test output.
202  !
203  ! 10. Source code :
204  !
205  !/ ------------------------------------------------------------------- /
206  USE constants, ONLY: tpiinv
207  USE w3gdatmd, ONLY: nk, nth, nspec, sig, dth, dden, wwnmeanp, &
210 #ifdef W3_T
211  USE w3odatmd, ONLY: ndst
212 #endif
213 #ifdef W3_S
214  USE w3servmd, ONLY: strace
215 #endif
216  !
217  IMPLICIT NONE
218  !/
219  !/ ------------------------------------------------------------------- /
220  !/ Parameter list
221  !/
222  REAL, INTENT(IN) :: A(NTH,NK), CG(NK), WN(NK), U, UDIR
223  REAL, INTENT(IN) :: TAUWX, TAUWY
224  LOGICAL, INTENT(IN) :: LLWS(NSPEC)
225  REAL, INTENT(INOUT) :: USTAR ,USDIR
226  REAL, INTENT(OUT) :: EMEAN, FMEAN, FMEANS, WNMEAN, AMAX, &
227  CD, Z0, CHARN, FMEANWS
228  !/
229  !/ ------------------------------------------------------------------- /
230  !/ Local parameters
231  !/
232  INTEGER :: IS, IK, ITH
233 #ifdef W3_S
234  INTEGER, SAVE :: IENT = 0
235 #endif
236 
237  REAL :: TAUW, EBAND, EMEANWS, UNZ, &
238  EB(NK),EB2(NK),ALFA(NK)
239  !/
240  !/ ------------------------------------------------------------------- /
241  !/
242 #ifdef W3_S
243  CALL strace (ient, 'W3SPR3')
244 #endif
245  !
246  unz = max( 0.01 , u )
247  ustar = max( 0.0001 , ustar )
248  !
249  emean = 0.
250  emeanws= 0.
251  fmeanws= 0.
252  fmean = 0.
253  fmeans = 0.
254  wnmean = 0.
255  amax = 0.
256  !
257  ! 1. Integral over directions and maximum --------------------------- *
258  !
259  DO ik=1, nk
260  eb(ik) = 0.
261  eb2(ik) = 0.
262  DO ith=1, nth
263  is=ith+(ik-1)*nth
264  eb(ik) = eb(ik) + a(ith,ik)
265  IF (llws(is)) eb2(ik) = eb2(ik) + a(ith,ik)
266  amax = max( amax , a(ith,ik) )
267  END DO
268  END DO
269  !
270  ! 2. Integrate over directions -------------------------------------- *
271  !
272  DO ik=1, nk
273  alfa(ik) = 2. * dth * sig(ik) * eb(ik) * wn(ik)**3
274  eb(ik) = eb(ik) * dden(ik) / cg(ik)
275  eb2(ik) = eb2(ik) * dden(ik) / cg(ik)
276  emean = emean + eb(ik)
277  fmean = fmean + eb(ik) *(sig(ik)**(2.*wwnmeanptail))
278  fmeans = fmeans + eb(ik) *(sig(ik)**(2.*wwnmeanp))
279  wnmean = wnmean + eb(ik) *(wn(ik)**wwnmeanp)
280  emeanws = emeanws+ eb2(ik)
281  fmeanws = fmeanws+ eb2(ik)*(sig(ik)**(2.*wwnmeanptail))
282  END DO
283  !
284  ! 3. Add tail beyond discrete spectrum and get mean pars ------------ *
285  ! ( DTH * SIG absorbed in FTxx )
286  !
287  eband = eb(nk) / dden(nk)
288  emean = emean + eband * fte
289  fmean = fmean + eband * sstxftftail
290  fmeans = fmeans + eband * sstxftf
291  wnmean = wnmean + eband * sstxftwn
292  eband = eb2(nk) / dden(nk)
293  emeanws = emeanws + eband * fte
294  fmeanws = fmeanws + eband * sstxftftail
295  !
296  ! 4. Final processing
297  !
298 
299  IF (fmean.LT.1.e-7) THEN
300  fmean=tpiinv * sig(nk)
301  ELSE
302  fmean = tpiinv *( max( 1.e-7 , fmean ) &
303  / max( 1.e-7 , emean ))**(1/(2.*wwnmeanptail))
304  END IF
305  IF (fmeans.LT.1.e-7) THEN
306  fmeans=tpiinv * sig(nk)
307  ELSE
308  fmeans = tpiinv *( max( 1.e-7 , fmeans ) &
309  / max( 1.e-7 , emean ))**(1/(2.*wwnmeanp))
310  END IF
311  wnmean = ( max( 1.e-7 , wnmean ) &
312  / max( 1.e-7 , emean ) )**(1/wwnmeanp)
313  IF (fmeanws.LT.1.e-7.OR.emeanws.LT.1.e-7) THEN
314  fmeanws=tpiinv * sig(nk)
315  ELSE
316  fmeanws = tpiinv *( max( 1.e-7 , fmeanws ) &
317  / max( 1.e-7 , emeanws ))**(1/(2.*wwnmeanptail))
318  END IF
319  !
320  ! 5. Cd and z0 ----------------------------------------------- *
321  !
322  tauw = sqrt(tauwx**2+tauwy**2)
323 
324  z0=0.
325  CALL calc_ustar(u,tauw,ustar,z0,charn)
326  unz = max( 0.01 , u )
327  cd = (ustar/unz)**2
328  usdir = udir
329  !
330  ! 6. Final test output ---------------------------------------------- *
331  !
332 #ifdef W3_T
333  WRITE (ndst,9060) emean, wnmean, tpiinv, ustar, cd, z0
334 #endif
335  !
336  RETURN
337  !
338  ! Formats
339  !
340 #ifdef W3_T
341 9060 FORMAT (' TEST W3SPR3 : E,WN MN :',f8.3,f8.4/ &
342  ' TPIINV, USTAR, CD, Z0 :',f8.3,f7.2,1x,2f9.5)
343 #endif
344  !/
345  !/ End of W3SPR3 ----------------------------------------------------- /
346  !/

References calc_ustar(), w3gdatmd::dden, w3gdatmd::dth, w3gdatmd::fte, w3gdatmd::ftf, w3odatmd::ndst, w3gdatmd::nk, w3gdatmd::nspec, w3gdatmd::nth, w3gdatmd::sig, w3gdatmd::sstxftf, w3gdatmd::sstxftftail, w3gdatmd::sstxftwn, w3gdatmd::sswellf, w3servmd::strace(), constants::tpiinv, w3gdatmd::wwnmeanp, and w3gdatmd::wwnmeanptail.

Referenced by wmesmfmd::fieldindex(), gxexpo(), w3exnc(), w3expo(), and w3srcemd::w3srce().

Variable Documentation

◆ delalp

real w3src3md::delalp

Definition at line 99 of file w3src3md.F90.

Referenced by tabu_tauhf(), and w3sin3().

◆ deltauw

real w3src3md::deltauw

Definition at line 97 of file w3src3md.F90.

Referenced by calc_ustar(), and tabu_stress().

◆ delu

real w3src3md::delu

Definition at line 97 of file w3src3md.F90.

Referenced by calc_ustar(), and tabu_stress().

◆ delust

real w3src3md::delust

Definition at line 99 of file w3src3md.F90.

Referenced by tabu_tauhf(), and w3sin3().

◆ ialpha

integer, parameter w3src3md::ialpha =200

Definition at line 96 of file w3src3md.F90.

Referenced by tabu_tauhf(), and w3sin3().

◆ ilevtail

integer, parameter w3src3md::ilevtail =50

Definition at line 96 of file w3src3md.F90.

◆ itaumax

integer, parameter w3src3md::itaumax =200

Definition at line 95 of file w3src3md.F90.

95  INTEGER, PARAMETER :: ITAUMAX=200,jumax=200

Referenced by calc_ustar(), and tabu_stress().

◆ iustar

integer, parameter w3src3md::iustar =100

Definition at line 96 of file w3src3md.F90.

96  INTEGER, PARAMETER :: IUSTAR=100,ialpha=200, ilevtail=50

Referenced by tabu_tauhf(), and w3sin3().

◆ jumax

integer, parameter w3src3md::jumax =200

Definition at line 95 of file w3src3md.F90.

Referenced by calc_ustar(), and tabu_stress().

◆ tauhft

real, dimension(0:iustar,0:ialpha) w3src3md::tauhft

Definition at line 99 of file w3src3md.F90.

99  REAL :: TAUHFT(0:IUSTAR,0:IALPHA), DELUST, DELALP

Referenced by tabu_tauhf(), and w3sin3().

◆ taut

real, dimension(0:itaumax,0:jumax) w3src3md::taut

Definition at line 97 of file w3src3md.F90.

97  REAL :: TAUT(0:ITAUMAX,0:JUMAX), DELTAUW, DELU

Referenced by calc_ustar(), and tabu_stress().

◆ tauwmax

real, parameter w3src3md::tauwmax = 2.2361

Definition at line 101 of file w3src3md.F90.

101  REAL, PARAMETER :: TAUWMAX = 2.2361 !SQRT(5.)

Referenced by calc_ustar(), and tabu_stress().

◆ umax

real, parameter w3src3md::umax = 50.

Definition at line 100 of file w3src3md.F90.

100  REAL, PARAMETER :: UMAX = 50.

Referenced by tabu_stress().

w3gdatmd::nk
integer, pointer nk
Definition: w3gdatmd.F90:1230
w3gdatmd::ddelta2
real, pointer ddelta2
Definition: w3gdatmd.F90:1311
w3gdatmd::dth
real, pointer dth
Definition: w3gdatmd.F90:1232
w3gdatmd::sstxftftail
real, pointer sstxftftail
Definition: w3gdatmd.F90:1340
w3gdatmd::nspec
integer, pointer nspec
Definition: w3gdatmd.F90:1230
w3gdatmd::wwnmeanp
real, pointer wwnmeanp
Definition: w3gdatmd.F90:1311
w3gdatmd::ddelta1
real, pointer ddelta1
Definition: w3gdatmd.F90:1311
w3gdatmd::sstxftwn
real, pointer sstxftwn
Definition: w3gdatmd.F90:1311
w3src4md::delust
real delust
Definition: w3src4md.F90:96
w3gdatmd::sig
real, dimension(:), pointer sig
Definition: w3gdatmd.F90:1234
w3gdatmd::wwnmeanptail
real, pointer wwnmeanptail
Definition: w3gdatmd.F90:1340
w3gdatmd::fachfe
real, pointer fachfe
Definition: w3gdatmd.F90:1232
w3gdatmd::ecos
real, dimension(:), pointer ecos
Definition: w3gdatmd.F90:1234
constants::rade
real, parameter rade
RADE Conversion factor from radians to degrees.
Definition: constants.F90:76
w3gdatmd::dden2
real, dimension(:), pointer dden2
Definition: w3gdatmd.F90:1234
w3gdatmd::zzwnd
real, pointer zzwnd
Definition: w3gdatmd.F90:1311
w3arrymd::outmat
subroutine outmat(NDS, A, MX, NX, NY, MNAME)
Definition: w3arrymd.F90:988
w3gdatmd::th
real, dimension(:), pointer th
Definition: w3gdatmd.F90:1234
w3src4md::taut
real, dimension(:,:), allocatable taut
Definition: w3src4md.F90:95
w3gdatmd::zzalp
real, pointer zzalp
Definition: w3gdatmd.F90:1311
w3gdatmd::zz0rat
real, pointer zz0rat
Definition: w3gdatmd.F90:1311
w3gdatmd::sstxftf
real, pointer sstxftf
Definition: w3gdatmd.F90:1311
w3gdatmd::esin
real, dimension(:), pointer esin
Definition: w3gdatmd.F90:1234
w3gdatmd::bbeta
real, pointer bbeta
Definition: w3gdatmd.F90:1311
w3src4md::delalp
real delalp
Definition: w3src4md.F90:96
w3servmd
Definition: w3servmd.F90:3
constants::tpiinv
real, parameter tpiinv
TPIINV Inverse of 2*Pi.
Definition: constants.F90:74
w3gdatmd::nth
integer, pointer nth
Definition: w3gdatmd.F90:1230
w3odatmd
Definition: w3odatmd.F90:3
w3gdatmd::sig2
real, dimension(:), pointer sig2
Definition: w3gdatmd.F90:1234
constants::tpi
real, parameter tpi
TPI 2*Pi.
Definition: constants.F90:72
w3servmd::strace
subroutine strace(IENT, SNAME)
Definition: w3servmd.F90:148
w3arrymd
Definition: w3arrymd.F90:3
w3gdatmd::xfr
real, pointer xfr
Definition: w3gdatmd.F90:1232
w3gdatmd::ssinthp
real, pointer ssinthp
Definition: w3gdatmd.F90:1311
w3gdatmd::fte
real, pointer fte
Definition: w3gdatmd.F90:1232
w3src4md::tauhft
real, dimension(:,:), allocatable tauhft
Definition: w3src4md.F90:95
w3src4md::deltauw
real deltauw
Definition: w3src4md.F90:96
w3odatmd::ndst
integer, pointer ndst
Definition: w3odatmd.F90:456
constants
Define some much-used constants for global use (all defined as PARAMETER).
Definition: constants.F90:20
w3gdatmd::dden
real, dimension(:), pointer dden
Definition: w3gdatmd.F90:1234
w3gdatmd
Definition: w3gdatmd.F90:16
w3gdatmd::sswellf
real, dimension(:), pointer sswellf
Definition: w3gdatmd.F90:1311
w3arrymd::prt2ds
subroutine prt2ds(NDS, NFR0, NFR, NTH, E, FR, UFR, FACSP, FSC, RRCUT, PRVAR, PRUNIT, PNTNME)
Definition: w3arrymd.F90:1943
w3gdatmd::aalpha
real, pointer aalpha
Definition: w3gdatmd.F90:1311
w3gdatmd::ftf
real, pointer ftf
Definition: w3gdatmd.F90:1232
w3src4md::delu
real delu
Definition: w3src4md.F90:96
w3gdatmd::ec2
real, dimension(:), pointer ec2
Definition: w3gdatmd.F90:1234
w3gdatmd::zz0max
real, pointer zz0max
Definition: w3gdatmd.F90:1311
constants::grav
real, parameter grav
GRAV Acc.
Definition: constants.F90:61