WAVEWATCH III  beta 0.0.1
w3src2md Module Reference

Tolman and Chalikov (1996) input and dissipation source terms. More...

Functions/Subroutines

subroutine w3spr2 (A, CG, WN, DEPTH, FPI, U, USTAR, EMEAN, FMEAN, WNMEAN, AMAX, ALFA, FP)
 Calculate mean wave parameters for the use in the source term routines (Tolman and Chalikov). More...
 
subroutine w3sin2 (A, CG, K, U, UDIR, CD, Z0, FPI, S, D)
 Calculate input source term. More...
 
subroutine w3sds2 (A, CG, K, FPI, USTAR, ALFA, S, D)
 Calculate whitecapping source term and diagonal term of derivative. More...
 
subroutine inptab
 Generate an interpolation table for the air-sea interaction parameter of Chalikov and Belevich (1993). More...
 

Detailed Description

Tolman and Chalikov (1996) input and dissipation source terms.

Bundled with interpolation tables.

Author
H. L. Tolman
Date
29-May-2009

Function/Subroutine Documentation

◆ inptab()

subroutine w3src2md::inptab

Generate an interpolation table for the air-sea interaction parameter of Chalikov and Belevich (1993).

The size of the table is set in parameter statements, the range is set by the input parameters of this routine. The first counter of the table corresponds to the nondimensional frequency

                  SIGMA Ul
        SIGA  =  ----------  COS ( THETA - THETA     )           (1)
                     g                          wind

The second counter of the table represents the drag coefficient. The maximum values of both parameters are passed to the routine through the parameter list.

Author
H. L. Tolman
Date
21-Feb-2004

Definition at line 848 of file w3src2md.F90.

848  !/
849  !/ +-----------------------------------+
850  !/ | WAVEWATCH III NOAA/NCEP |
851  !/ | H. L. Tolman |
852  !/ | FORTRAN 90 |
853  !/ | Last update : 21-Feb-2004 |
854  !/ +-----------------------------------+
855  !/
856  !/ 03-Jun-1996 : Final version 1.18 / FORTRAN 77 version.
857  !/ 06-Dec-1999 : Upgrade to FORTRAN 90 ( version 2.00 )
858  !/ 21-Feb-2004 : Multiple model version. ( version 3.06 )
859  !/
860  ! 1. Purpose :
861  !
862  ! Generate an interpolation table for the air-sea interaction
863  ! parameter of Chalikov and Belevich (1993).
864  !
865  ! 2. Method :
866  !
867  ! The size of the table is set in parameter statements, the range
868  ! is set by the input parameters of this routine. The first counter
869  ! of the table corresponds to the nondimensional frequency
870  !
871  ! SIGMA Ul
872  ! SIGA = ---------- COS ( THETA - THETA ) (1)
873  ! g wind
874  !
875  ! The second counter of the table represents the drag coefficient.
876  ! The maximum values of both parameters are passed to the routine
877  ! through the parameter list.
878  !
879  ! 3. Parameters :
880  !
881  ! 4. Subroutines used :
882  !
883  ! Name Type Module Description
884  ! ----------------------------------------------------------------
885  ! STRACE Subr. W3SERVMD Subroutine tracing.
886  ! W3BETA Func. Internal Function to calculate the
887  ! interaction parameter.
888  ! ----------------------------------------------------------------
889  !
890  ! 5. Called by :
891  !
892  ! Name Type Module Description
893  ! ----------------------------------------------------------------
894  ! W3IOGR Subr. W3IOGRMD Model definition IO routine.
895  ! ----------------------------------------------------------------
896  !
897  ! 6. Error messages :
898  !
899  ! None.
900  !
901  ! 7. Remarks :
902  !
903  ! 8. Structure :
904  !
905  ! See source code.
906  !
907  ! 9. Switches :
908  !
909  ! !/S Enable subroutine tracing.
910  ! !/T Enable test output.
911  ! !/T0 Print table.
912  ! !/T1 Estimate maximum errors.
913  !
914  ! 10. Source code :
915  !
916  !/ ------------------------------------------------------------------- /
917  USE constants
918  USE w3odatmd, ONLY: ndst
919 #ifdef W3_S
920  USE w3servmd, ONLY: strace
921 #endif
922  !/
923  IMPLICIT NONE
924  !/
925  !/ ------------------------------------------------------------------- /
926  !/ Parameter list
927  !/
928  !/ ------------------------------------------------------------------- /
929  !/ Local parameters
930  !/
931  INTEGER :: ISIGA, IDRAG
932 #ifdef W3_S
933  INTEGER, SAVE :: IENT = 0
934 #endif
935 #ifdef W3_T0
936  INTEGER :: I1
937 #endif
938 #ifdef W3_T1
939  INTEGER :: IE1
940 #endif
941  REAL :: SIGA, DRAG
942 #ifdef W3_T0
943  REAL :: BMIN, BMAX
944 #endif
945 #ifdef W3_T1
946  REAL :: ENORM, ERR(NRDRAG)
947 #endif
948  !/
949  !/ ------------------------------------------------------------------- /
950  !/
951 #ifdef W3_S
952  CALL strace (ient, 'INPTAB')
953 #endif
954  !
955  ! 1. Determine range and increments of table ------------------------ *
956 
957  !
958  dsiga = sigamx / real(nrsiga)
959  ddrag = dragmx / real(nrdrag)
960  !
961 #ifdef W3_T
962  WRITE (ndst,9000) sigamx, dsiga, dragmx, ddrag
963 #endif
964  !
965  ! 2. Fill table ----------------------------------------------------- *
966  !
967  DO isiga=-nrsiga,nrsiga+1
968  siga = real(isiga) * dsiga
969  DO idrag=1, nrdrag+1
970  drag = real(idrag) * ddrag
971  betatb(isiga,idrag) = w3beta( siga, drag , ndst )
972  END DO
973  END DO
974  !
975  ! 3. Test output ---------------------------------------------------- *
976  !
977 #ifdef W3_T0
978  WRITE (ndst,9010)
979  i1 = min(35,nrdrag)
980  DO isiga=-nrsiga,nrsiga
981  siga = real(isiga) * dsiga
982  bmin = 0.
983  bmax = 0.
984  DO idrag=1, nrdrag
985  bmin = min( bmin , betatb(isiga,idrag) )
986  bmax = max( bmax , betatb(isiga,idrag) )
987  END DO
988  bmax = max( bmax , -bmin )
989  WRITE (ndst,9011) isiga, siga, bmax, &
990  (nint(betatb(isiga,idrag)/bmax*100.),idrag=1,i1)
991  IF (i1.LT.nrdrag) WRITE (ndst,9012) &
992  (nint(betatb(isiga,idrag)/bmax*100.),idrag=i1+1,nrdrag)
993  END DO
994 #endif
995  !
996 #ifdef W3_T1
997  WRITE (ndst,9020)
998  ie1 = min(30,nrdrag-1)
999  enorm = 1000. / abs(betatb(0,nrdrag))
1000  DO isiga=-nrsiga,nrsiga
1001  siga = real(isiga) * dsiga
1002  IF ( abs(siga) .LT. 5.01 ) THEN
1003  DO idrag=1, nrdrag-1
1004  drag = ( real(idrag) + 0.5 ) * ddrag
1005  err(idrag) = - w3beta(siga,drag,ndst) + 0.5 * &
1006  ( betatb(isiga,idrag) + betatb(isiga,idrag+1) )
1007  END DO
1008  WRITE (ndst,9021) isiga, siga, &
1009  (nint(enorm*err(idrag)),idrag=1,ie1)
1010  IF (ie1.LT.nrdrag-1) WRITE (ndst,9022) &
1011  (nint(enorm*err(idrag)),idrag=ie1+1,nrdrag-1)
1012  ENDIF
1013  END DO
1014  !
1015  WRITE (ndst,9030)
1016  ie1 = min(30,nrdrag)
1017  enorm = 1000. / abs(betatb(0,nrdrag))
1018  DO isiga=-nrsiga,nrsiga-1
1019  siga = ( real(isiga) + 0.5 ) * dsiga
1020  IF ( abs(siga) .LT. 5.01 ) THEN
1021  DO idrag=1, nrdrag
1022  drag = real(idrag) * ddrag
1023  err(idrag) = - w3beta(siga,drag,ndst) + 0.5 * &
1024  ( betatb(isiga,idrag) + betatb(isiga+1,idrag) )
1025  END DO
1026  WRITE (ndst,9031) isiga, siga, &
1027  (nint(enorm*err(idrag)),idrag=1,ie1)
1028  IF (ie1.LT.nrdrag) WRITE (ndst,9032) &
1029  (nint(enorm*err(idrag)),idrag=ie1+1,nrdrag)
1030  ENDIF
1031  END DO
1032 #endif
1033  !
1034  RETURN
1035  !
1036  ! Formats
1037  !
1038 #ifdef W3_T
1039 9000 FORMAT ( ' TEST INPTAB : SIGAMX, DSIGA : ',f6.2,f8.2/ &
1040  ' DRAGMX, DDRAG : ',f8.4,f9.5)
1041 #endif
1042  !
1043 #ifdef W3_T0
1044 9010 FORMAT (/' TEST INPTAB : TABLE, NORMALIZED WITH ', &
1045  'BETATB(ISIGA,NRDRAG)'/ &
1046  ' ISIGA, SIGA, BETA_MAX, TABLE (x100)')
1047 9011 FORMAT (1x,i4,f7.2,f6.4,1x,35i3)
1048 9012 FORMAT (19x,35i3)
1049 #endif
1050  !
1051 #ifdef W3_T1
1052 9020 FORMAT (/' TEST INPTAB : ERROR DUE TO DRAG, NORMALIZED ', &
1053  'WITH BETATB(ISIGA,NRDRAG)'/ &
1054  ' ISIGA, SIGA, TABLE (x1000)')
1055 9021 FORMAT (1x,i4,f7.2,35i3)
1056 9022 FORMAT (12x,35i3)
1057 9030 FORMAT (/' TEST INPTAB : ERROR DUE TO SIGA, NORMALIZED WITH ', &
1058  'BETATB(ISIGA,NRDRAG)'/ &
1059  ' ISIGA, SIGA, TABLE (x1000)')
1060 9031 FORMAT (1x,i4,f7.2,35i3)
1061 9032 FORMAT (12x,35i3)
1062 #endif
1063  !/
1064  !/ Internal function W3BETA
1065  !/
1066  CONTAINS
1067  !/ ------------------------------------------------------------------- /
1081  REAL FUNCTION W3BETA ( OMA , CL , NDST )
1082  !/
1083  !/ +-----------------------------------+
1084  !/ | WAVEWATCH III NOAA/NCEP |
1085  !/ | H. L. Tolman |
1086  !/ | D.Chalikov |
1087  !/ | FORTRAN 90 |
1088  !/ | Last update : 21-Feb-2004 |
1089  !/ +-----------------------------------+
1090  !/
1091  !/ 06-Dec-1996 : Final version 1.18 / FORTRAN 77 version.
1092  !/ 06-Dec-1999 : Upgrade to FORTRAN 90 ( version 2.00 )
1093  !/ 21-Feb-2004 : Multiple model version. ( version 3.06 )
1094  !/
1095  ! 1. Purpose :
1096  !
1097  ! Calculate wind-wave interaction parameter beta.
1098  !
1099  ! 2. Method :
1100  !
1101  ! Chalikov and Belevich (1992), see also manual.
1102  !
1103  ! 3. Parameters :
1104  !
1105  ! Parameter list
1106  ! ----------------------------------------------------------------
1107  ! W3BETA Real O Wind-wave interaction parameter multiplied
1108  ! by density ratio.
1109  ! OMA Real I Non-dimensional apparent frequency.
1110  !
1111  ! OMA = OMEGA | U | cos(theta-theta ) / g
1112  ! l w
1113  !
1114  ! CL Real I Drag coefficient at height l
1115  ! ----------------------------------------------------------------
1116  !
1117  ! 4. Subroutines used :
1118  !
1119  ! 5. Called by :
1120  !
1121  ! 6. Error messages :
1122  !
1123  ! 7. Remarks :
1124  !
1125  ! 8. Structure :
1126  !
1127  ! See source code.
1128  !
1129  ! 9. Switches :
1130  !
1131  ! !/S Enable subroutine tracing.
1132  ! !/T0 Enable test output.
1133  !
1134  ! 10. Source code :
1135  !
1136  !/ ------------------------------------------------------------------- /
1137  IMPLICIT NONE
1138  !/
1139  !/ ------------------------------------------------------------------- /
1140  !/ Parameter list
1141  !/
1142  INTEGER, INTENT(IN) :: NDST
1143  REAL, INTENT(IN) :: OMA, CL
1144  !/
1145  !/ ------------------------------------------------------------------- /
1146  !/ Local parameters
1147  !/
1148 #ifdef W3_S
1149  INTEGER, SAVE :: IENT = 0
1150 #endif
1151  REAL :: OM1, OM2, A0, A1, A2, A3, A4, A5, &
1152  A6, A7, A8, A9, A10
1153  !/
1154  !/ ------------------------------------------------------------------- /
1155  !/
1156 #ifdef W3_S
1157  CALL strace (ient, 'W3BETA')
1158 #endif
1159  !
1160 #ifdef W3_T0
1161  WRITE (ndst,9000) oma, cl
1162 #endif
1163  !
1164  ! calculate Omegas
1165  !
1166  om1 = 1.075 + 75.*cl
1167  om2 = 1.2 + 300.*cl
1168  !
1169  ! calculate factors a
1170  !
1171  a1 = 0.25 + 395.*cl
1172  a2 = 0.35 + 150.*cl
1173  a4 = 0.3 + 300.*cl
1174  a9 = 0.35 + 240.*cl
1175  a10 = -0.06 + 470.*cl
1176  !
1177  a5 = a4 * om1
1178  a0 = 0.25 * a5**2 / a4
1179  a3 = (a0-a2-a1) / (a0+a4+a5)
1180  a6 = a0 * (1.-a3)
1181  a7 = (a9*(om2-1)**2+a10) / (om2-om1)
1182  a8 = a7 * om1
1183  !
1184 #ifdef W3_T0
1185  WRITE (ndst,9001) om1, om2
1186  WRITE (ndst,9002) a0, a1, a2, a3, a4, a5, a6, a7, a8, a9, a10
1187 #endif
1188  !
1189  ! calculate beta * 1.e4
1190  !
1191  IF ( oma .LT. -1. ) THEN
1192  w3beta = -a1 * oma**2 - a2
1193  ELSE IF (oma .LT. om1/2.) THEN
1194  w3beta = a3 * oma * ( a4 * oma - a5 ) - a6
1195  ELSE IF (oma .LT. om1) THEN
1196  w3beta = oma * ( a4 * oma - a5 )
1197  ELSE IF (oma .LT. om2) THEN
1198  w3beta = a7 * oma - a8
1199  ELSE
1200  w3beta = a9 * (oma-1.)**2 + a10
1201  END IF
1202  !
1203  ! beta * dwat / dair
1204  !
1205  w3beta = w3beta * 1.e-4
1206 #ifdef W3_T0
1207  WRITE (ndst,9003) w3beta
1208 #endif
1209  !
1210  RETURN
1211  !
1212  ! Formats
1213  !
1214 #ifdef W3_T0
1215 9000 FORMAT ( ' TEST W3BETA : INPUT : ',2e10.3)
1216 9001 FORMAT ( ' TEST W3BETA : OM1-2 : ',2e10.3)
1217 9002 FORMAT ( ' TEST W3BETA : A0-10 : ',5e10.3/ &
1218  ' ',6e10.3)
1219 9003 FORMAT ( ' TEST W3BETA : BETA : ',e10.3)
1220 #endif
1221  !/
1222  !/ End of W3BETA ----------------------------------------------------- /
1223  !/
1224  END FUNCTION w3beta
1225  !/
1226  !/ End of INPTAB ----------------------------------------------------- /
1227  !/

References w3odatmd::ndst, w3servmd::strace(), and w3beta().

Referenced by w3iogrmd::w3iogr().

◆ w3sds2()

subroutine w3src2md::w3sds2 ( real, dimension(nth,nk), intent(in)  A,
real, dimension(nk), intent(in)  CG,
real, dimension(nk), intent(in)  K,
real, intent(in)  FPI,
real, intent(in)  USTAR,
real, dimension(nk), intent(in)  ALFA,
real, dimension(nth,nk), intent(out)  S,
real, dimension(nth,nk), intent(out)  D 
)

Calculate whitecapping source term and diagonal term of derivative.

Parameters
[in]AInput action density spectrum.
[in]CGGroup velocity array.
[in]KWavenumber array.
[in]FPI'Peak frequency' of input (rad/s).
[in]USTARFriction velocity (m/s).
[in]ALFAPhillips' constant.
[out]SSource term (1-D version).
[out]DDiagonal term of derivative (1-D version).
Author
H. L. Tolman
Date
21-Feb-2004

Definition at line 583 of file w3src2md.F90.

583  !/
584  !/ +-----------------------------------+
585  !/ | WAVEWATCH III NOAA/NCEP |
586  !/ | H. L. Tolman |
587  !/ | FORTRAN 90 |
588  !/ | Last update : 21-Feb-2004 |
589  !/ +-----------------------------------+
590  !/
591  !/ 12-Jun-1996 : Final FORTRAN 77 ( version 1.18 )
592  !/ 04-Feb-2000 : Upgrade to FORTRAN 90 ( version 2.00 )
593  !/ 23-Apr-2002 : Erick Rogers' fix ( version 2.19 )
594  !/ 21-Feb-2004 : Multiple model version. ( version 3.06 )
595  !/
596  ! 1. Purpose :
597  !
598  ! Calculate whitecapping source term and diagonal term of der.
599  !
600  ! 2. Method :
601  !
602  ! Tolman and Chalikov (1995).
603  !
604  ! 3. Parameters :
605  !
606  ! Parameter list
607  ! ----------------------------------------------------------------
608  ! A R.A. I Input action density spectrum.
609  ! CG R.A. I Group velocity array.
610  ! K R.A. I Wavenumber array.
611  ! FPI Real I 'Peak frequency' of input (rad/s).
612  ! USTAR Real I Friction velocity (m/s).
613  ! ALFA R.A. I Phillips' constant.
614  ! S R.A. O Source term (1-D version).
615  ! D R.A. O Diagonal term of derivative (1-D version).
616  ! ----------------------------------------------------------------
617  !
618  ! 4. Subroutines used :
619  !
620  ! Name Type Module Description
621  ! ----------------------------------------------------------------
622  ! STRACE Subr. W3SERVMD Subroutine tracing.
623  ! PRT2DS Subr. W3ARRYMD Print plot of spectra.
624  ! OUTMAT Subr. W3WRRYMD Print out 2D matrix.
625  ! ----------------------------------------------------------------
626  !
627  ! 5. Called by :
628  !
629  ! Name Type Module Description
630  ! ----------------------------------------------------------------
631  ! W3SRCE Subr. W3SRCEMD Source term integration.
632  ! W3EXPO Subr. N/A Point output post-processor.
633  ! GXEXPO Subr. N/A GrADS point output post-processor.
634  ! ----------------------------------------------------------------
635  !
636  ! 6. Error messages :
637  !
638  ! 7. Remarks :
639  !
640  ! 8. Structure :
641  !
642  ! See source code.
643  !
644  ! 9. Switches :
645  !
646  ! !/S Enable subroutine tracing.
647  ! !/T Enable general test output.
648  ! !/T0 Print arrays.
649  ! !/T1 Print filter and constituents.
650  !
651  ! 10. Source code :
652  !
653  !/ ------------------------------------------------------------------- /
654  USE constants
655  USE w3gdatmd, ONLY: nk, nth, sig, dden, dth, fte, fpimin, &
656  facti1, facti2, xf1, xf2, xfh, sdsaln, &
657  cdsa0, cdsa1, cdsa2, cdsb0, cdsb1, cdsb2, &
658  cdsb3
659 #ifdef W3_T
660  USE w3odatmd, ONLY: ndst
661 #endif
662 #ifdef W3_S
663  USE w3servmd, ONLY: strace
664 #endif
665 #ifdef W3_T0
666  USE w3arrymd, ONLY: prt2ds
667 #endif
668 #ifdef W3_T1
669  USE w3arrymd, ONLY: outmat
670 #endif
671  !
672  IMPLICIT NONE
673  !/
674  !/ ------------------------------------------------------------------- /
675  !/ Parameter list
676  !/
677  REAL, INTENT(IN) :: A(NTH,NK), CG(NK), K(NK), FPI, &
678  USTAR, ALFA(NK)
679  REAL, INTENT(OUT) :: S(NTH,NK), D(NTH,NK)
680  !/
681  !/ ------------------------------------------------------------------- /
682  !/ Local parameters
683  !/
684  INTEGER :: IK, ITH, IKHW
685 #ifdef W3_S
686  INTEGER, SAVE :: IENT = 0
687 #endif
688  REAL :: FHW, XHW, FPIT, PHI, AF1, AF2, &
689  AFILT, BFILT, CDIST, FILT, POW, &
690  CDISH, CDISP, HW, EHIGH, EBD(NK)
691 #ifdef W3_T
692  REAL POWMAX
693 #endif
694 #ifdef W3_T0
695  REAL DOUT(NK,NTH)
696 #endif
697  !/
698  !/ ------------------------------------------------------------------- /
699  !/
700 #ifdef W3_S
701  CALL strace (ient, 'W3SDS2')
702 #endif
703  !
704 #ifdef W3_T
705  WRITE (ndst,9000) fpi, ustar
706 #endif
707  !
708  ! 1. Preparations
709  ! 1.a HW
710  !
711  fhw = xfh*fpi
712  xhw = facti2 + facti1*log(fhw)
713  ikhw = min( nk , int( xhw + 0.5 ) )
714  DO ik=ikhw, nk
715  ebd(ik) = 0.
716  DO ith=1, nth
717  ebd(ik) = ebd(ik) + a(ith,ik)
718  END DO
719  END DO
720  !
721  IF ( fhw .LT. sig(nk+1) ) THEN
722  xhw = 1. - mod( xhw + 0.5 , 1. )
723  IF ( ikhw .EQ. nk ) xhw = max( 0. , xhw - 0.5 )
724  hw = xhw * ebd(ikhw)*dden(ikhw)/cg(ikhw)
725  DO ik=ikhw+1, nk
726  hw = hw + ebd(ik)*dden(ik)/cg(ik)
727  END DO
728  hw = 4. * sqrt( hw + ebd(nk)/cg(nk)*fte )
729  ELSE
730  ehigh = ebd(nk)/cg(nk) * sig(nk)*dth * (sig(nk)/fhw)**5
731  hw = 4. * sqrt( 0.25 * fhw * ehigh )
732  END IF
733  !
734  ! 1.b PHI
735  !
736  fpit = max( fpimin , fpi*tpiinv*ustar/grav )
737  phi = cdsb0 + cdsb1*fpit + cdsb2/fpit**cdsb3
738  !
739  ! 1.c Set-up filter
740  !
741  af2 = xf2*fpi
742  af1 = xf1*fpi
743  bfilt = 1. / ( af2 - af1 )
744  afilt = - bfilt * af1
745  !
746  ! 1.d Constants
747  !
748  cdist = - 2. * ustar * hw * phi
749  cdish = g2pi3i * ustar**2
750  cdisp = g1pi1i * ustar
751  !
752  ! 2. Combined diagonal factor
753  !
754 #ifdef W3_T2
755  WRITE (ndst,9020)
756 #endif
757 #ifdef W3_T
758  powmax = 0.
759 #endif
760  DO ik=1, nk
761  filt = min( 1., max( 0. , afilt + bfilt*sig(ik) ))
762  pow = min( 25. , cdsa1 / ( cdisp*sig(ik) )**cdsa2 )
763  IF ( filt .GT. 0. ) THEN
764  d(1,ik) = (1.-filt) * cdist * k(ik)**2 &
765  - filt * cdsa0 * cdish * sig(ik)**3 &
766  * (alfa(ik)/sdsaln)**pow
767  ELSE
768  d(1,ik) = (1.-filt) * cdist * k(ik)**2
769  END IF
770 #ifdef W3_T
771  powmax = max(pow*filt,powmax)
772 #endif
773 #ifdef W3_T2
774  WRITE (ndst,9021) ik, filt, alfa(ik)/sdsaln, &
775  cdist*phi*k(ik)**2, cdsa0*cdish*sig(ik)**3 &
776  * (alfa(ik)/sdsaln)**pow, d(1,ik)
777 #endif
778  END DO
779  !
780 #ifdef W3_T
781  WRITE (ndst,9010) af1, af2, afilt, bfilt, powmax
782 #endif
783  !
784  ! 3. 2-D diagonal array
785  !
786  DO ik=1, nk
787  DO ith=2, nth
788  d(ith,ik) = d(1,ik)
789  END DO
790  END DO
791  !
792  s = d * a
793  !
794  ! ... Test output of arrays
795  !
796 #ifdef W3_T0
797  DO ik=1, nk
798  DO ith=1, nth
799  dout(ik,ith) = d(ith,ik)
800  END DO
801  END DO
802  CALL prt2ds (ndst, nk, nk, nth, dout, sig(1:), ' ', 1., &
803  0.0, 0.001, 'Diag Sds', ' ', 'NONAME')
804 #endif
805  !
806 #ifdef W3_T1
807  CALL outmat (ndst, d, nth, nth, nk, 'diag Sds')
808 #endif
809  !
810  RETURN
811  !
812  ! Formats
813  !
814 #ifdef W3_T
815 9000 FORMAT (' TEST W3SDS2 : FPI, USTAR : ',2f8.3)
816 9010 FORMAT (' TEST W3SDS2 : AF1-2, A-BFILT, PMAX : ',4f7.3,e10.3)
817 #endif
818 #ifdef W3_T2
819 9020 FORMAT (' TEST W3SDS2 : IK, FILT, ALFA, DDST, DDSH, DDS')
820 9021 FORMAT (' ',i6,2f7.3,3e11.3)
821 #endif
822  !/
823  !/ End of W3SDS2 ----------------------------------------------------- /
824  !/

References w3gdatmd::cdsa0, w3gdatmd::cdsa1, w3gdatmd::cdsa2, w3gdatmd::cdsb0, w3gdatmd::cdsb1, w3gdatmd::cdsb2, w3gdatmd::cdsb3, w3gdatmd::dden, w3gdatmd::dth, w3gdatmd::facti1, w3gdatmd::facti2, w3gdatmd::fpimin, w3gdatmd::fte, constants::g1pi1i, constants::g2pi3i, constants::grav, w3odatmd::ndst, w3gdatmd::nk, w3gdatmd::nth, w3arrymd::outmat(), w3arrymd::prt2ds(), w3gdatmd::sdsaln, w3gdatmd::sig, w3servmd::strace(), constants::tpiinv, w3gdatmd::xf1, w3gdatmd::xf2, and w3gdatmd::xfh.

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

◆ w3sin2()

subroutine w3src2md::w3sin2 ( 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)  UDIR,
real, intent(in)  CD,
real, intent(in)  Z0,
real, intent(out)  FPI,
real, dimension(nspec), intent(out)  S,
real, dimension(nspec), intent(out)  D 
)

Calculate input source term.

Parameters
[in]AAction density spectrum (1-D).
[in]CGGroup velocities for k-axis of spectrum.
[in]KWavenumber for entire spectrum (1-D).
[in]UWind speed at reference height.
[in]UDIRDirection of U.
[in]CDDrag coefficient at wind level ZWIND.
[in]Z0Corresponding z0.
[out]FPIInput 'peak' frequency.
[out]SSource term (1-D version).
[out]DDiagonal term of derivative (1-D version).
Author
H. L. Tolman
D. Chalikov
Date
21-Feb-2004

Definition at line 309 of file w3src2md.F90.

309  !/
310  !/ +-----------------------------------+
311  !/ | WAVEWATCH III NOAA/NCEP |
312  !/ | H. L. Tolman |
313  !/ | D.Chalikov |
314  !/ | FORTRAN 90 |
315  !/ | Last update : 21-Feb-2004 |
316  !/ +-----------------------------------+
317  !/
318  !/ 14-Jan-1997 : Final FORTRAN 77 ( version 1.18 )
319  !/ 04-Feb-2000 : Upgrade to FORTRAN 90 ( version 2.00 )
320  !/ 21-Feb-2004 : Multiple model version. ( version 3.06 )
321  !/
322  ! 1. Purpose :
323  !
324  ! Calculate input source term.
325  !
326  ! 2. Method :
327  !
328  ! Tolman and Chalikov (1996), see manual.
329  !
330  ! 3. Parameters :
331  !
332  ! Parameter list
333  ! ----------------------------------------------------------------
334  ! A R.A. I Action density spectrum (1-D).
335  ! CG R.A. I Group velocities for k-axis of spectrum.
336  ! K R.A. I Wavenumber for entire spectrum (1-D).
337  ! U Real I Wind speed at reference height.
338  ! UDIR Real I Direction of U.
339  ! CD Real I Drag coefficient at wind level ZWIND.
340  ! Z0 Real I Corresponding z0.
341  ! FPI R.A. O Input 'peak' frequency.
342  ! S R.A. O Source term (1-D version).
343  ! D R.A. O Diagonal term of derivative (1-D version).
344  ! ----------------------------------------------------------------
345  !
346  ! 4. Subroutines used :
347  !
348  ! Name Type Module Description
349  ! ----------------------------------------------------------------
350  ! STRACE Subr. W3SERVMD Subroutine tracing.
351  ! PRT2DS Subr. W3ARRYMD Print plot of spectra.
352  ! OUTMAT Subr. W3WRRYMD Print out 2D matrix.
353  ! ----------------------------------------------------------------
354  !
355  ! 5. Called by :
356  !
357  ! Name Type Module Description
358  ! ----------------------------------------------------------------
359  ! W3SRCE Subr. W3SRCEMD Source term integration.
360  ! W3EXPO Subr. N/A Point output post-processor.
361  ! GXEXPO Subr. N/A GrADS point output post-processor.
362  ! ----------------------------------------------------------------
363  !
364  ! 6. Error messages :
365  !
366  ! 7. Remarks :
367  !
368  ! - Actual height of wind speed does not need to be 10 m, but is
369  ! given by ZWIND.
370  ! - Abs(cos) > 0.0087 to asure continuity in beta. Corresponds
371  ! to shift of up to half a degree.
372  !
373  ! 8. Structure :
374  !
375  ! See source code.
376  !
377  ! 9. Switches :
378  !
379  ! !/S Enable subroutine tracing.
380  ! !/T Enable general test output.
381  ! !/T0 Print arrays.
382  ! !/T1 Calculation of diagonal within spectrum
383  !
384  ! 10. Source code :
385  !
386  !/ ------------------------------------------------------------------- /
387  USE constants
388  USE w3gdatmd, ONLY: nk, nth, nspec, xfr, dden, sig, sig2, &
389  esin, ecos, fte, fttr, fpimin, zwind, &
391 #ifdef W3_T
392  USE w3odatmd, ONLY: ndst
393 #endif
394 #ifdef W3_S
395  USE w3servmd, ONLY: strace
396 #endif
397 #ifdef W3_T0
398  USE w3arrymd, ONLY: prt2ds
399 #endif
400 #ifdef W3_T1
401  USE w3arrymd, ONLY: outmat
402 #endif
403  !/
404  IMPLICIT NONE
405  !/
406  !/ ------------------------------------------------------------------- /
407  !/ Parameter list
408  !/
409  REAL, INTENT(IN) :: A(NSPEC), CG(NK), K(NSPEC), U, UDIR, &
410  CD, Z0
411  REAL, INTENT(OUT) :: S(NSPEC), D(NSPEC), FPI
412  !/
413  !/ ------------------------------------------------------------------- /
414  !/ Local parameters
415  !/
416  INTEGER :: IS, IK, IOMA, ICL, NKFILT, NKFIL2
417 #ifdef W3_S
418  INTEGER, SAVE :: IENT = 0
419 #endif
420 #ifdef W3_T0
421  INTEGER ITH
422 #endif
423  REAL :: COSU, SINU, COSFAC, LAMBDA, ULAM, &
424  CLAM, OMA, M0, M1, RD1, RD2, BETA, &
425  FACLN1, FACLN2, USTAR, TRANS, FPISTR,&
426  FP1STR, FP1, SIN1A(NK)
427  REAL, PARAMETER :: TRANSF = 0.75
428  REAL, PARAMETER :: PEAKFC = 0.8
429 #ifdef W3_T0
430  REAL :: DOUT(NK,NTH)
431 #endif
432  !/
433  !/ ------------------------------------------------------------------- /
434  !/
435 #ifdef W3_S
436  CALL strace (ient, 'W3SIN2')
437 #endif
438  !
439 #ifdef W3_T
440  WRITE (ndst,9000) dsiga, ddrag, u, udir*rade, cd, z0
441 #endif
442  !
443  ! 1. Preparations
444  !
445  cosu = cos(udir)
446  sinu = sin(udir)
447  !
448  ! 2. Loop over spectrum
449  !
450 #ifdef W3_T2
451  WRITE (ndst,9020)
452 #endif
453  !
454  facln1 = u / log(zwind/z0)
455  facln2 = log(z0)
456  !
457  DO is=1, nspec
458  cosfac = ecos(is)*cosu + esin(is)*sinu
459  cosfac = sign( max( 0.0087 , abs(cosfac) ) , cosfac )
460  lambda = tpi / ( k(is) * abs(cosfac) )
461  ulam = facln1 * ( log(lambda) - facln2 )
462  clam = cd * ( u / ulam )**2
463  oma = k(is) * ulam * cosfac / sig2(is)
464  ioma = int( oma/dsiga ) + &
465  min( 0 , int( sign( -1.1 , oma ) ) )
466  icl = int( clam/ddrag )
467  rd1 = oma/dsiga - real(ioma)
468  rd2 = clam/ddrag - real(icl)
469  ioma = max( -nrsiga , min( nrsiga , ioma ) )
470  icl = max( 1 , min( nrdrag , icl ) )
471  beta = (1.-rd1) * (1.-rd2) * betatb( ioma , icl ) &
472  + rd1 * (1.-rd2) * betatb(ioma+1, icl ) &
473  + (1.-rd1) * rd2 * betatb( ioma ,icl+1) &
474  + rd1 * rd2 * betatb(ioma+1,icl+1)
475  d(is) = beta * sig2(is)
476  s(is) = a(is) * d(is)
477 #ifdef W3_T2
478  WRITE (ndst,9021) is, cosfac, lambda, ulam, clam*1.e3, &
479  oma, beta*1.e4
480 #endif
481  END DO
482  !
483  ! 3. Calculate FPI
484  !
485  DO ik=1, nk
486  sin1a(ik) = 0.
487  DO is=(ik-1)*nth+1, ik*nth
488  sin1a(ik) = sin1a(ik) + max( 0. , s(is) )
489  END DO
490  END DO
491  !
492  m0 = 0.
493  m1 = 0.
494  DO ik=1, nk
495  sin1a(ik) = sin1a(ik) * dden(ik) / ( cg(ik) * sig(ik)**3 )
496  m0 = m0 + sin1a(ik)
497  m1 = m1 + sin1a(ik)/sig(ik)
498  END DO
499  !
500  sin1a(nk) = sin1a(nk) / dden(nk)
501  m0 = m0 + sin1a(nk) * fte
502  m1 = m1 + sin1a(nk) * fttr
503  IF ( m1 .LT. 1e-20 ) THEN
504  fpi = xfr * sig(nk)
505  ELSE
506  fpi = m0 / m1
507  END IF
508  !
509  ! 4. Filter for swell
510  !
511  ustar = u * sqrt(cd)
512  fpistr = max( fpimin , fpi * ustar / grav )
513  fp1str = 3.6e-4 + 0.92*fpistr - 6.3e-10/fpistr**3
514  fp1 = peakfc * fp1str * grav / ustar
515  !
516  nkfilt = min( nk , int(facti2+facti1*log(fp1)) )
517  nkfil2 = min( nk , int(facti2+facti1*log(transf*fp1)) )
518  nkfil2 = max( 0 , nkfil2 )
519  !
520  DO is=1, nkfil2*nth
521  d(is) = max( d(is) , fswell*d(is) )
522  s(is) = a(is) * d(is)
523  END DO
524  !
525  DO ik=nkfil2+1, nkfilt
526  trans = ( sig(ik)/fp1 - transf ) / (1.-transf)
527  DO is=(ik-1)*nth+1, ik*nth
528  d(is) = (1.-trans)*max(d(is),fswell*d(is)) + trans*d(is)
529  s(is) = a(is) * d(is)
530  END DO
531  END DO
532  !
533  ! ... Test output of arrays
534  !
535 #ifdef W3_T0
536  DO ik=1, nk
537  DO ith=1, nth
538  dout(ik,ith) = d(ith+(ik-1)*nth)
539  END DO
540  END DO
541  CALL prt2ds (ndst, nk, nk, nth, dout, sig(1:), ' ', 1., &
542  0.0, 0.001, 'Diag Sin', ' ', 'NONAME')
543 #endif
544  !
545 #ifdef W3_T1
546  CALL outmat (ndst, d, nth, nth, nk, 'diag Sin')
547 #endif
548  !
549  RETURN
550  !
551  ! Formats
552  !
553 #ifdef W3_T
554 9000 FORMAT (' TEST W3SIN2 : DSIGA,DDRAG,U,UDIR,CD,Z0(IN) : '/ &
555  ' ',f8.4,f9.6,f7.2,f6.1,f8.5,f8.5)
556 #endif
557  !
558 #ifdef W3_T2
559 9020 FORMAT (' TEST W3SIN2 : IS, COS, LAMBDA, ULAM, CLAM*1E3, ', &
560  'OMA, BETA*1E4')
561 9021 FORMAT (6x,i6,f7.2,1x,f6.1,2(1x,f5.2),2(1x,f6.2))
562 #endif
563  !/
564  !/ End of W3SIN2 ----------------------------------------------------- /
565  !/

References w3gdatmd::dden, w3gdatmd::ecos, w3gdatmd::esin, w3gdatmd::facti1, w3gdatmd::facti2, w3gdatmd::fpimin, w3gdatmd::fswell, w3gdatmd::fte, w3gdatmd::fttr, constants::grav, w3odatmd::ndst, w3gdatmd::nk, w3gdatmd::nspec, w3gdatmd::nth, w3arrymd::outmat(), w3arrymd::prt2ds(), constants::rade, w3gdatmd::sig, w3gdatmd::sig2, w3servmd::strace(), constants::tpi, w3gdatmd::xfr, and w3gdatmd::zwind.

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

◆ w3spr2()

subroutine w3src2md::w3spr2 ( real, dimension(nth,nk), intent(in)  A,
real, dimension(nk), intent(in)  CG,
real, dimension(nk), intent(in)  WN,
real, intent(in)  DEPTH,
real, intent(in)  FPI,
real, intent(in)  U,
real, intent(in)  USTAR,
real, intent(out)  EMEAN,
real, intent(out)  FMEAN,
real, intent(out)  WNMEAN,
real, intent(out)  AMAX,
real, dimension(nk), intent(out)  ALFA,
real, intent(out)  FP 
)

Calculate mean wave parameters for the use in the source term routines (Tolman and Chalikov).

Parameters
[in]AAction density spectrum.
[in]CGGroup velocities.
[in]WNWavenumbers.
[in]DEPTHWater depth.
[in]FPIPeak input frequency.
[in]UWind speed.
[in]USTARFriction velocity.
[out]EMEANTotal energy (variance).
[out]FMEANMean frequency.
[out]WNMEANMean wavenumber.
[out]AMAXMaximum of action spectrum.
[out]ALFAPhillips' constant.
[out]FPPeak frequency.
Author
H. L. Tolman
D. Chalikov
Date
13-Apr-2007

Definition at line 127 of file w3src2md.F90.

127  !/
128  !/ +-----------------------------------+
129  !/ | WAVEWATCH III NOAA/NCEP |
130  !/ | H. L. Tolman |
131  !/ | D.Chalikov |
132  !/ | FORTRAN 90 |
133  !/ | Last update : 13-Apr-2007 |
134  !/ +-----------------------------------+
135  !/
136  !/ 06-Dec-1996 : Final version 1.18 / FORTRAN 77 version.
137  !/ 16-Nov-1999 : Add itteration to section 5. for removal of W3APR2.
138  !/ 04-Feb-2000 : Upgrade to FORTRAN 90 ( version 2.00 )
139  !/ 21-Dec-2004 : Multiple model version. ( version 3.06 )
140  !/ 03-Jul-2006 : Extract stress computation. ( version 3.09 )
141  !/ 13-Apr-2007 : EMEAN in parameter list. ( version 3.11 )
142  !
143  ! 1. Purpose :
144  !
145  ! Calculate mean wave parameters for the use in the source term
146  ! routines. (Tolman and Chalikov)
147  !
148  ! 2. Method :
149  !
150  ! See source term routines.
151  !
152  ! 3. Parameters :
153  !
154  ! Parameter list
155  ! ----------------------------------------------------------------
156  ! A R.A. I Action density spectrum.
157  ! CG R.A. I Group velocities.
158  ! WN R.A. I Wavenumbers.
159  ! DEPTH Real I Water depth.
160  ! FPI Real I Peak input frequency.
161  ! U Real I Wind speed.
162  ! USTAR Real I Friction velocity.
163  ! EMEAN Real O Total energy (variance).
164  ! FMEAN Real O Mean frequency.
165  ! WNMEAN Real O Mean wavenumber.
166  ! AMAX Real O Maximum of action spectrum.
167  ! ALFA R.A. O Phillips' constant.
168  ! FP Real O Peak frequency.
169  ! ----------------------------------------------------------------
170  !
171  ! 4. Subroutines used :
172  !
173  ! Name Type Module Description
174  ! ----------------------------------------------------------------
175  ! STRACE Subr. W3SERVMD Subroutine tracing.
176  ! ----------------------------------------------------------------
177  !
178  ! 5. Called by :
179  !
180  ! Name Type Module Description
181  ! ----------------------------------------------------------------
182  ! W3SRCE Subr. W3SRCEMD Source term integration.
183  ! W3EXPO Subr. N/A Point output post-processor.
184  ! GXEXPO Subr. N/A GrADS point output post-processor.
185  ! ----------------------------------------------------------------
186  !
187  ! 6. Error messages :
188  !
189  ! 7. Remarks :
190  !
191  ! 8. Structure :
192  !
193  ! See source code.
194  !
195  ! 9. Switches :
196  !
197  ! !/S Enable subroutine tracing.
198  ! !/T Enable test output.
199  !
200  ! 10. Source code :
201  !
202  !/ ------------------------------------------------------------------- /
203  USE constants
204  USE w3gdatmd, ONLY: nk, nth, dth, sig, dden, fte, ftf, ftwn, &
206 #ifdef W3_T
207  USE w3odatmd, ONLY: ndst
208 #endif
209 
210 #ifdef W3_S
211  USE w3servmd, ONLY: strace
212 #endif
213  USE w3dispmd, ONLY: nar1d, dfac, n1max, ecg1, ewn1, dsie
214  !/
215  IMPLICIT NONE
216  !/
217  !/ ------------------------------------------------------------------- /
218  !/ Parameter list
219  !/
220  REAL, INTENT(IN) :: A(NTH,NK), CG(NK), WN(NK), DEPTH, &
221  FPI, U, USTAR
222  REAL, INTENT(OUT) :: EMEAN, FMEAN, WNMEAN, AMAX, &
223  ALFA(NK), FP
224  !/
225  !/ ------------------------------------------------------------------- /
226  !/ Local parameters
227  !/
228  INTEGER :: IK, ITH, I1, ITT
229 #ifdef W3_S
230  INTEGER, SAVE :: IENT = 0
231 #endif
232  REAL :: EBAND, FPISTR, EB(NK), UST
233  !/
234  !/ ------------------------------------------------------------------- /
235  !/
236 #ifdef W3_S
237  CALL strace (ient, 'W3SPR2')
238 #endif
239  !
240  ust = max( 0.0001 , ustar )
241  !
242  emean = 0.
243  fmean = 0.
244  wnmean = 0.
245  amax = 0.
246  !
247  ! 1. Integral over directions and maximum --------------------------- *
248  !
249  DO ik=1, nk
250  eb(ik) = 0.
251  DO ith=1, nth
252  eb(ik) = eb(ik) + a(ith,ik)
253  amax = max( amax , a(ith,ik) )
254  END DO
255  END DO
256  !
257  ! 2. Integrate over directions -------------------------------------- *
258  !
259  DO ik=1, nk
260  alfa(ik) = 2. * dth * sig(ik) * eb(ik) * wn(ik)**3
261  eb(ik) = eb(ik) * dden(ik) / cg(ik)
262  emean = emean + eb(ik)
263  fmean = fmean + eb(ik) / sig(ik)
264  wnmean = wnmean + eb(ik) / sqrt(wn(ik))
265  END DO
266  !
267  ! 3. Add tail beyond discrete spectrum and get mean pars ------------ *
268  ! ( DTH * SIG absorbed in FTxx )
269  !
270  eband = eb(nk) / dden(nk)
271  emean = emean + eband * fte
272  fmean = fmean + eband * ftf
273  wnmean = wnmean + eband * ftwn
274  !
275  fmean = tpiinv * emean / max( 1.e-7 , fmean )
276  wnmean = ( emean / max( 1.e-7 , wnmean ) )**2
277  !
278  ! 4. Estimate peak frequency from FPI ------------------------------- *
279  !
280  fpistr = max( 0.008 , fpi * ust / grav )
281  fp = ( 3.6e-4 + 0.92*fpistr - 6.3e-10/fpistr**3 )/ust*grav
282  fp = fp * tpiinv
283  !
284  RETURN
285  !/
286  !/ End of W3SPR2 ----------------------------------------------------- /
287  !/

References w3gdatmd::cinxsi, w3gdatmd::dden, w3dispmd::dfac, w3dispmd::dsie, w3gdatmd::dth, w3dispmd::ecg1, w3dispmd::ewn1, w3gdatmd::fte, w3gdatmd::ftf, w3gdatmd::ftwn, constants::grav, w3dispmd::n1max, w3dispmd::nar1d, w3odatmd::ndst, w3gdatmd::nittin, w3gdatmd::nk, w3gdatmd::nth, w3gdatmd::sig, w3servmd::strace(), constants::tpiinv, and w3gdatmd::zwind.

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

w3dispmd::dfac
real, parameter dfac
Definition: w3dispmd.F90:75
w3gdatmd::nk
integer, pointer nk
Definition: w3gdatmd.F90:1230
w3gdatmd::dth
real, pointer dth
Definition: w3gdatmd.F90:1232
w3gdatmd::nspec
integer, pointer nspec
Definition: w3gdatmd.F90:1230
constants::g1pi1i
real, parameter g1pi1i
G1PI1I Inverse of gravity * 2 * Pi.
Definition: constants.F90:82
w3snl4md::oma
real, dimension(:), allocatable oma
Definition: w3snl4md.F90:143
w3gdatmd::cdsa2
real, pointer cdsa2
Definition: w3gdatmd.F90:1304
w3gdatmd::sig
real, dimension(:), pointer sig
Definition: w3gdatmd.F90:1234
w3gdatmd::cdsb3
real, pointer cdsb3
Definition: w3gdatmd.F90:1304
constants::g2pi3i
real, parameter g2pi3i
G2PI3I Inverse of gravity^2 * (2*Pi)^3.
Definition: constants.F90:81
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::xf1
real, pointer xf1
Definition: w3gdatmd.F90:1304
w3gdatmd::fswell
real, pointer fswell
Definition: w3gdatmd.F90:1304
w3gdatmd::fpimin
real, pointer fpimin
Definition: w3gdatmd.F90:1304
w3arrymd::outmat
subroutine outmat(NDS, A, MX, NX, NY, MNAME)
Definition: w3arrymd.F90:988
w3dispmd::ewn1
real, dimension(0:nar1d) ewn1
Definition: w3dispmd.F90:78
w3gdatmd::cdsb0
real, pointer cdsb0
Definition: w3gdatmd.F90:1304
w3gdatmd::esin
real, dimension(:), pointer esin
Definition: w3gdatmd.F90:1234
w3gdatmd::cdsa1
real, pointer cdsa1
Definition: w3gdatmd.F90:1304
w3gdatmd::sdsaln
real, pointer sdsaln
Definition: w3gdatmd.F90:1304
w3servmd
Definition: w3servmd.F90:3
constants::tpiinv
real, parameter tpiinv
TPIINV Inverse of 2*Pi.
Definition: constants.F90:74
w3gdatmd::cdsb1
real, pointer cdsb1
Definition: w3gdatmd.F90:1304
w3gdatmd::nth
integer, pointer nth
Definition: w3gdatmd.F90:1230
w3gdatmd::xf2
real, pointer xf2
Definition: w3gdatmd.F90:1304
w3gdatmd::fttr
real, pointer fttr
Definition: w3gdatmd.F90:1232
w3odatmd
Definition: w3odatmd.F90:3
w3gdatmd::cdsa0
real, pointer cdsa0
Definition: w3gdatmd.F90:1304
w3dispmd::ecg1
real, dimension(0:nar1d) ecg1
Definition: w3dispmd.F90:78
w3gdatmd::cinxsi
real, pointer cinxsi
Definition: w3gdatmd.F90:1282
w3gdatmd::cdsb2
real, pointer cdsb2
Definition: w3gdatmd.F90:1304
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
w3gdatmd::facti2
real, pointer facti2
Definition: w3gdatmd.F90:1232
w3arrymd
Definition: w3arrymd.F90:3
w3gdatmd::xfr
real, pointer xfr
Definition: w3gdatmd.F90:1232
w3gdatmd::fte
real, pointer fte
Definition: w3gdatmd.F90:1232
w3odatmd::ndst
integer, pointer ndst
Definition: w3odatmd.F90:456
w3beta
real function w3beta(OMA, CL, NDST)
Calculate wind-wave interaction parameter beta.
Definition: w3src2md.F90:1082
constants
Define some much-used constants for global use (all defined as PARAMETER).
Definition: constants.F90:20
w3gdatmd::nittin
integer, pointer nittin
Definition: w3gdatmd.F90:1281
w3gdatmd::dden
real, dimension(:), pointer dden
Definition: w3gdatmd.F90:1234
w3gdatmd
Definition: w3gdatmd.F90:16
w3dispmd::nar1d
integer, parameter nar1d
Definition: w3dispmd.F90:74
w3arrymd::prt2ds
subroutine prt2ds(NDS, NFR0, NFR, NTH, E, FR, UFR, FACSP, FSC, RRCUT, PRVAR, PRUNIT, PNTNME)
Definition: w3arrymd.F90:1943
w3gdatmd::zwind
real, pointer zwind
Definition: w3gdatmd.F90:1304
w3gdatmd::ftf
real, pointer ftf
Definition: w3gdatmd.F90:1232
w3gdatmd::facti1
real, pointer facti1
Definition: w3gdatmd.F90:1232
w3dispmd
Definition: w3dispmd.F90:3
w3dispmd::n1max
integer n1max
Definition: w3dispmd.F90:77
w3gdatmd::ftwn
real, pointer ftwn
Definition: w3gdatmd.F90:1232
w3gdatmd::xfh
real, pointer xfh
Definition: w3gdatmd.F90:1304
constants::grav
real, parameter grav
GRAV Acc.
Definition: constants.F90:61
w3dispmd::dsie
real dsie
Definition: w3dispmd.F90:78