WAVEWATCH III  beta 0.0.1
w3fld2md Module Reference

Functions/Subroutines

subroutine w3fld2 (ASPC, FPI, WNDX, WNDY, ZWND, DEPTH, RIB, DAIR, UST, USTD, Z0, TAUNUX, TAUNUY, CHARN)
 
subroutine wnd2sat (WND10, SAT)
 

Function/Subroutine Documentation

◆ w3fld2()

subroutine w3fld2md::w3fld2 ( real, dimension(nspec), intent(in)  ASPC,
real, intent(in)  FPI,
real, intent(in)  WNDX,
real, intent(in)  WNDY,
real, intent(in)  ZWND,
real, intent(in)  DEPTH,
real, intent(in)  RIB,
real, intent(in)  DAIR,
real, intent(out)  UST,
real, intent(out)  USTD,
real, intent(out)  Z0,
real, intent(inout)  TAUNUX,
real, intent(inout)  TAUNUY,
real, intent(out), optional  CHARN 
)

Definition at line 67 of file w3fld2md.F90.

67  !/
68  !/ +-----------------------------------+
69  !/ | WAVEWATCH III NOAA/NCEP/NOPP |
70  !/ | B. G. Reichl |
71  !/ | FORTRAN 90 |
72  !/ | Last update : 22-Mar-2021 |
73  !/ +-----------------------------------+
74  !/
75  !/ 01-Jul-2013 : Origination (version 3.14)
76  !/ 19-Mar-2015 : Clean-up for submission (version 5.12)
77  !/ 22-Mar-2021 : Consider DAIR a variable ( version 7.13 )
78  !/
79  ! 1. Purpose :
80  !
81  ! Wind stress vector calculation from wave spectrum and
82  ! n-meter wind speed vector.
83  !
84  ! 2. Method :
85  ! See Donelan et al. (2012).
86  !
87  ! 3. Parameters :
88  !
89  ! Parameter list
90  ! ----------------------------------------------------------------
91  ! ASPC Real I 1-D Wave action spectrum.
92  ! FPI Real I Peak input frequency.
93  ! WNDX Real I X-dir wind (assumed referenced to current)
94  ! WNDY Real I Y-dir wind (assumed referenced to current)
95  ! ZWND Real I Wind height.
96  ! DEPTH Real I Water depth.
97  ! RIB Real I Bulk Richardson number in lower atm
98  ! DAIR Real I Air density
99  ! TAUNUX Real 0 X-dir viscous stress (guessed from prev.)
100  ! TAUNUY Real 0 Y-dir viscous stress (guessed from prev.)
101  ! UST Real O Friction velocity.
102  ! USTD Real O Direction of friction velocity.
103  ! Z0 Real O Surface roughness length
104  ! CHARN Real O,optional Charnock parameter
105  ! ----------------------------------------------------------------
106  !
107  ! 4. Subroutines used :
108  !
109  ! Name Type Module Description
110  ! ----------------------------------------------------------------
111  ! STRACE Subr. W3SERVMD Subroutine tracing.
112  ! APPENDTAIL Subr. W3FLD1MD Modification of tail for calculation
113  ! SIG2WN Subr. W3FLD1MD Depth-dependent dispersion relation
114  ! MFLUX Subr. W3FLD1MD MO stability correction
115  ! WND2Z0M Subr. W3FLD1MD Bulk Z0 from wind
116  ! CALC_FPI Subr. W3FLD1MD Calculate peak frequency
117  ! ----------------------------------------------------------------
118  !
119  ! 5. Called by :
120  !
121  ! Name Type Module Description
122  ! ----------------------------------------------------------------
123  ! W3ASIM Subr. W3ASIMMD Air-sea interface module.
124  ! W3EXPO Subr. N/A Point output post-processor.
125  ! GXEXPO Subr. N/A GrADS point output post-processor.
126  ! ----------------------------------------------------------------
127  !
128  ! 6. Error messages :
129  !
130  ! None.
131  !
132  ! 7. Remarks :
133  !
134  ! 8. Structure :
135  !
136  ! See source code.
137  !
138  ! 9. Switches :
139  !
140  ! !/S Enable subroutine tracing.
141  !
142  ! 10. Source code :
143  !
144  !/ ------------------------------------------------------------------- /
145  USE constants, ONLY: dwat, grav, tpi, pi, kappa
146  USE w3gdatmd, ONLY: nk, nth, nspec, sig, dth, xfr, th
147  USE w3odatmd, ONLY: ndse
148  USE w3servmd, ONLY: extcde
152 #ifdef W3_S
153  USE w3servmd, ONLY: strace
154 #endif
155  !/
156  IMPLICIT NONE
157  !/
158  !/ ------------------------------------------------------------------- /
159  !/ Parameter list
160  !/
161  REAL, INTENT(IN) :: ASPC(NSPEC), WNDX, WNDY, &
162  ZWND, DEPTH, RIB, DAIR, FPI
163  REAL, INTENT(OUT) :: UST, USTD, Z0
164  REAL, INTENT(OUT),OPTIONAL :: CHARN
165  REAL, INTENT(INOUT) :: TAUNUX, TAUNUY
166  !/
167  !/ ------------------------------------------------------------------- /
168  !/ Local parameters
169  !/
170  !-Parameters
171  REAL, PARAMETER :: NU=0.105/10000.0
172  !-Commonly used values
173  REAL :: UREF, UREFD
174  !-Tail
175  REAL :: SAT
176  REAL :: KMAX, KTAILA, KTAILB, KTAILC
177  INTEGER :: KA1, KA2, KA3, NKT
178  !-Extended spectrum
179  REAL, ALLOCATABLE, DIMENSION(:) :: WN, DWN, CP, sig2,TAUINTX, TAUINTY
180  REAL, ALLOCATABLE, DIMENSION(:,:) :: SPC2
181  !-Stress Calculation
182  INTEGER :: K, T, ITS
183  REAL :: TAUXW, TAUYW, TAUX, TAUY
184  REAL :: USTRA, USTRB, USTSM
185  REAL :: A1, SCIN
186  REAL :: CD, CDF, CDS
187  real :: wnd_z, wnd_z_mag, wnd_z_proj, wnd_effect
188  ! Stress iteration
189  REAL :: B1, B2
190  REAL :: USTRI1, USTRF1, USTRI2, USTRF2
191  REAL :: USTGRA, SLO
192  LOGICAL :: UST_IT_FLG(2)
193  !-Z0 iteration
194  REAL :: z01,z02
195  !-Wind iteration
196  real :: wnd_10_x, wnd_10_y, wnd_10_mag, wnd_10_dir
197  real :: u35_1, v35_1, u35_2, v35_2, u35_3, v35_3
198  REAL :: DIFU10xx, DIFU10yx, DIFU10xy, DIFU10yy
199  REAL :: fd_a, fd_b, fd_c, fd_d
200  REAL :: DU, DV, UITV, VITV, CH
201  REAL :: APAR, DTX(3), DTY(3), DT
202  LOGICAL :: WIFLG, WND_IT_FLG
203  !-MO stability correction
204  LOGICAL :: HEIGHTFLG
205  integer :: wi_count, wi
206  real :: wnd_ref_al,wnd_ref_ax
207  real :: wndpa, wndpax, wndpay, wndpe,wndpex, wndpey
208  LOGICAL :: NO_ERR
209  LOGICAL :: ITERFLAG
210  INTEGER :: ITTOT
211  INTEGER :: COUNT
212 #ifdef W3_S
213  INTEGER, SAVE :: IENT = 0
214 #endif
215  LOGICAL, SAVE :: FIRST = .true.
216 #ifdef W3_OMPG
217  !$omp threadprivate( FIRST )
218 #endif
219  !/
220  !/ ------------------------------------------------------------------- /
221  !/
222 #ifdef W3_S
223  CALL strace (ient, 'C3FLD2')
224 #endif
225  !
226  ! 0. Initializations ------------------------------------------------ *
227  !
228  ! **********************************************************
229  ! *** The initialization routine should include all ***
230  ! *** initialization, including reading data from files. ***
231  ! **********************************************************
232  !
233  IF ( first ) THEN
234  CALL infld
235  first = .false.
236  END IF
237  !----------------------------|
238  ! Calculate Reference height |
239  ! wind magnitude |
240  !----------------------------|
241  uref=sqrt(wndx**2+wndy**2)
242  urefd=atan2(wndy,wndx)
243  !----------------------------------------------|
244  ! Check if wind height not equal to 10 m |
245  !----------------------------------------------|
246  !HeightFLG = (abs(zwnd-10.).GT.0.1) ! True if not 10m
247  !----------------------------------------------|
248  ! Assume bulk and calculate 10 m wind guess for|
249  ! defining tail level |
250  !----------------------------------------------|
251  CALL wnd2z0m(uref,z01) ! first guess at z0
252  wnd_10_mag=uref
253  ittot=1
254  ! If input wind is not 10m, solve for approx 10 m
255  !-------------------------------------------------|
256  ! If height != 10 m, then iterate to get 10 m wind|
257  ! (assuming neutral -- this is just a guess) |
258  !-------------------------------------------------|
259  !IF (HeightFLG) THEN
260  ! IterFLAG=.true.
261  ! COUNT = 1 !COUNT is now counting iteration over z0
262  ! do while(IterFLAG)
263  ! wnd_10_mag=UREF*log(10./z01)/log(zwnd/z01)
264  ! CALL wnd2z0m(wnd_10_mag,z02)
265  ! if ( (abs(z01/z02-1.).GT.0.001) .AND. &
266  ! (COUNT.LT.10))THEN
267  ! z01 = z02
268  ! else
269  ! IterFLAG = .false.
270  ! endif
271  ! COUNT = COUNT + 1
272  ! enddo
273  ! ITTOT = 3 !extra iterations for 10m wind
274  !ELSE
275  ! wnd_10_mag = uref
276  ! ITTOT = 1 !no iteration needed
277  !ENDIF
278  if (tail_choice.eq.0) then
279  sat=tail_level
280  elseif (tail_choice.eq.1) then
281  CALL wnd2sat(wnd_10_mag,sat)
282  endif
283  ! now you have the guess at 10 m wind mag. and z01
284  !/
285  !--------------------------|
286  ! Get first guess at ustar |
287  !--------------------------|
288  ustra = uref*kappa/log(zwnd/z01)
289  ustd = urefd
290  wnd_10_dir = urefd
291  wnd_10_x=wnd_10_mag*cos(wnd_10_dir)
292  wnd_10_y=wnd_10_mag*sin(wnd_10_dir)
293  !
294  ! 1. Attach Tail ---------------------------------------------------- *
295  !
296  call sig2wn ( sig(nk),depth,kmax)
297  nkt = nk
298  DO WHILE ( kmax .LT. 366. )
299  nkt = nkt + 1
300  kmax = ( kmax * xfr**2 )
301  ENDDO
302  ALLOCATE( wn(nkt), dwn(nkt), cp(nkt),sig2(nkt), spc2(nkt,nth), &
303  tauintx(nkt),tauinty(nkt))
304  !|--------------------------------------------------------------------|
305  !|----Build Discrete Wavenumbers for defining spectrum on-------------|
306  !|--------------------------------------------------------------------|
307  DO k = 1, nk
308  call sig2wn(sig(k),depth,wn(k))
309  cp(k) = sig(k) / wn(k)
310  sig2(k) = sig(k)
311  ENDDO
312  DO k = ( nk + 1 ), ( nkt)
313  sig2(k)=sig2(k-1)*xfr
314  call sig2wn(sig2(k),depth,wn(k))
315  cp(k)=sig2(k)/wn(k)
316  ENDDO
317  DO k = 2, nkt-1
318  dwn(k) = (wn(k+1) - wn(k-1)) / 2.0
319  ENDDO
320  dwn(1) = ( wn(2)- ( wn(1) / (xfr **2.0) ) ) / 2.0
321  dwn(nkt) = ( wn(nkt)*(xfr**2.0) - wn(nkt-1)) / 2.0
322  !|---------------------------------------------------------------------|
323  !|---Attach initial tail-----------------------------------------------|
324  !|---------------------------------------------------------------------|
325  count=0 !Count is now counting step through 1-d spectrum
326  DO k=1, nk
327  DO t=1, nth
328  count = count + 1
329  spc2(k,t) = aspc(count) * sig(k)
330  ENDDO
331  ENDDO
332  DO k=nk+1, nkt
333  DO t=1, nth
334  spc2(k,t)=spc2(nk,t)*wn(nk)**3.0/wn(k)**(3.0)
335  ENDDO
336  ENDDO
337  !
338  ! 1c. Calculate transitions for new (constant saturation ) tail ------ *
339  !
340  !-----Wavenumber for beginning of (spectrum level) transition to tail- *
341  call sig2wn (fpi*tpi*tail_transition_ratio1,depth,ktaila )
342  !-----Wavenumber for end of (spectrum level) transition to tail------- *
343  call sig2wn (fpi*tpi*tail_transition_ratio2,depth,ktailb )
344  !-----Wavenumber for end of (spectrum direction) transition to tail--- *
345  ktailc= ktailb * 2.0
346  ka1 = 2 ! Do not modify 1st wavenumber bin
347  DO WHILE ( ( ktaila .GE. wn(ka1) ) .AND. (ka1 .LT. nkt-6) )
348  ka1 = ka1 + 1
349  ENDDO
350  ka2 = ka1+2
351  DO WHILE ( ( ktailb .GE. wn(ka2) ) .AND. (ka2 .LT. nkt-4) )
352  ka2 = ka2 + 1
353  ENDDO
354  ka3 = ka2+2
355  DO WHILE ( ( ktailc .GE. wn(ka3)) .AND. (ka3 .LT. nkt-2) )
356  ka3 = ka3 + 1
357  ENDDO
358  CALL appendtail(spc2,wn,nkt,ka1,ka2,ka3,atan2(wndy,wndx),sat)
359  ! Now the spectrum is set w/ tail level SAT
360  !
361  ! 2. Enter iteration ------------------------------------------------- *
362  !
363  ! Add new iteration for wind
364  !
365  ! Wind perturbations for iteration
366  !DT = 1.E-04
367  !DTX = (/ 1. , -1. , 0. /)
368  !DTY = (/ 0. , 1. , -1. /)
369  !/
370  heightflg=.false.!Not set-up for non-10 m winds
371  wiflg = .true. !This kicks out when wind iteration complete
372  no_err = .true. !This kicks out when there is an error
373  wi_count = 1 !Count is now counting wind iterations
374  ! - start of wind iteration (if applicable)
375  DO WHILE ( wiflg .AND. no_err ) !Wind iteration
376  !/
377  DO wi = 1, ittot !Newton-Raphson solve for derivatives if zwnd not 10 m.
378  ! If iterating over 10 m wind need to adjust guesses to get slopes
379  IF (heightflg) THEN
380  wnd_10_x = wnd_10_x + dtx(wi)*dt
381  wnd_10_y = wnd_10_y + dty(wi)*dt
382  wnd_10_mag = sqrt(wnd_10_x**2+wnd_10_y**2)
383  wnd_10_dir = atan2(wnd_10_y,wnd_10_x)
384  ENDIF
385  !
386  ! Stress iteration (inside wind iteration solve for stress)
387  its = 1 !ITS is counting stress iteration
388  ust_it_flg(1)=.true.
389  ust_it_flg(2)=.true.
390  DO WHILE ((ust_it_flg(1) .AND. ust_it_flg(2)) .AND. no_err)
391  !Get z0 from (guessed) stress and wind magnitude
392  z0 = 10. / ( exp( kappa * wnd_10_mag / ustra ) )
393  tauintx(1:nkt) = 0.0
394  tauinty(1:nkt) = 0.0
395  DO k = 1, nkt
396  !Waves 'feel' wind at height related to wavelength
397  wnd_z = min( pi / wn(k), 20.0 )
398  wnd_z_mag = ( ustra / kappa ) * (log(wnd_z/z0))
399  DO t = 1, nth
400  !projected component of wind in wave direction
401  wnd_z_proj = wnd_z_mag * cos( wnd_10_dir-th(t) )
402  IF (wnd_z_proj .GT. cp(k)) THEN
403  !Waves slower than wind
404  a1 = 0.11
405  ELSEIF (( wnd_z_proj .GE. 0 ) .AND. ( wnd_z_proj .LE. cp(k) )) THEN
406  !Wave faster than wind
407  a1 = 0.01
408  ELSEIF (wnd_z_proj .LT. 0) THEN
409  !Waves opposed to wind
410  a1 = 0.1
411  ENDIF
412  wnd_effect = wnd_z_proj - cp(k)
413  scin = a1 * wnd_effect * abs( wnd_effect ) * dair / dwat * &
414  wn(k) / cp(k)
415  ! -- Original version assumed g/Cp = sig,(a.k.a in deep water.)
416  ! TAUINTX(K) = TAUINTX(K) + SPC2(K,T) * SCIN &
417  ! * COS( TH(T) ) / CP(K) * DTH
418  ! TAUINTY(K) = TAUINTY(K) + SPC2(K,T) * SCIN &
419  ! * SIN( TH(T) ) / CP(K) * DTH
420 
421  tauintx(k) = tauintx(k) + spc2(k,t) * scin &
422  * cos( th(t) ) *sig2(k) * dth
423  tauinty(k) = tauinty(k) + spc2(k,t) * scin &
424  * sin( th(t) ) *sig2(k) * dth
425  ENDDO
426  ENDDO
427  tauyw = 0.0
428  tauxw = 0.0
429  DO k = 1, nkt
430  ! TAUXW = TAUXW + DWAT * GRAV * DWN(K) * TAUINTX(K)
431  ! TAUYW = TAUYW + DWAT * GRAV * DWN(K) * TAUINTY(K)
432  tauxw = tauxw + dwat * dwn(k) * tauintx(k)
433  tauyw = tauyw + dwat * dwn(k) * tauinty(k)
434  ENDDO
435  cdf = ( sqrt(tauxw**2.0+tauyw**2.0) / dair ) / wnd_10_mag**2.0
436  !|---------------------------------------------------------------------|
437  !|----Solve for the smooth drag coefficient to use as initial guess----|
438  !|----for the viscous stress-------------------------------------------|
439  !|---------------------------------------------------------------------|
440  IF (uref .LT. 0.01) THEN
441  ustsm = 0.0
442  iterflag = .false.
443  ELSE
444  z02 = 0.001
445  iterflag = .true.
446  ENDIF
447  count = 1
448  ! Finding smooth z0 to get smooth drag
449  DO WHILE( (iterflag ) .AND. (count .LT. 10) )
450  z01 = z02
451  ustsm = kappa * wnd_10_mag / ( log( 10. / z01 ) )
452  z02 = 0.132 * nu / ustsm
453  IF (abs( z02/z01-1.0) .LT. 10.0**(-4)) THEN
454  iterflag = .false.
455  ELSE
456  iterflag = .true.
457  ENDIF
458  count = count + 1
459  ENDDO
460  cds = ustsm**2.0 / wnd_10_mag**2.0
461  ! smooth drag adjustment based on full drag
462  cds = cds / 3.0 * ( 1.0 + 2.0 * cds / ( cds + cdf ) )
463  !-----Solve for viscous stress from smooth Cd
464  taunux = dair * cds * wnd_10_mag**2.0 * cos( wnd_10_dir )
465  taunuy = dair * cds * wnd_10_mag**2.0 * sin( wnd_10_dir )
466  !-----Sum drag components
467  taux = taunux + tauxw
468  tauy = taunuy + tauyw
469  !-----Calculate USTAR
470  ustrb = sqrt( sqrt( tauy**2.0 + taux**2.0) / dair )
471  !-----Calculate stress direction
472  ustd = atan2(tauy,taux)
473  !Checking ustar. ustra=guess. ustrb=found.
474  b1 = ( ustra - ustrb )
475  b2 = ( ustra + ustrb ) / 2.0
476  its = its + 1
477  !Check for convergence
478  ust_it_flg(1)=( abs(b1*100.0/b2) .GE. 0.01)
479  !If not converged after 20 iterations, quit.
480  ust_it_flg(2)=( its .LT. 20 )
481  IF ( ust_it_flg(1) .AND. ust_it_flg(2)) THEN
482  ! Toyed with methods for improving iteration.
483  ! ultimately this was sufficient.
484  ! May be imporved upon in future...
485  ustra = ustrb*.5 + ustrb*.5
486  ELSEIF (abs(b1*100.0/b2) .GE. 5.) THEN
487  !After 20 iterations, >5% from converged
488  ust_it_flg(1) = .false.
489  ust_it_flg(2) = .false.
490  print*,'Attn: Stress not converged for windspeed: ',uref
491  ust = -999.
492  no_err = .false.
493  ENDIF
494  ENDDO
495  !IF (HeightFLG) THEN
496  ! ! Get along stress wind at top wave boundary layer (10m)
497  ! WNDPA=WND_10_mag*COS(WND_10_dir-USTD)
498  ! WNDPE=WND_10_mag*SIN(WND_10_dir-USTD)
499  ! ! Calculate Cartesian of across wind
500  ! WNDPEx=WNDPE*cos(ustd+pi/2.)! add pi/2 since referenced
501  ! WNDPEy=WNDPE*sin(ustd+pi/2.)! to right of stress angle
502  ! !Approx as neutral inside 10 m (WBL) calculate z0
503  ! wnd_ref_al=uref*cos(urefd-ustd)
504  ! z0=10. / exp( wnd_10_mag * KAPPA / ustra )
505  ! ! Use that z0 to calculate stability
506  ! ! Cd to ref height (based on input wind)
507  ! Below is subroutine for computing stability effects
508  ! call mflux(wnd_ref_al,zwnd,z0,rib,cd)
509  ! 2. Get CD with stability
510  ! 3. Get New 35-m wind based on calculated stress Cd from MO
511  ! WNDPA=ustra/sqrt(cd)
512  ! WNDPAX=WNDPA*cos(ustd)
513  ! WNDPAY=WNDPA*sin(ustd)
514  ! if (wi.eq.3) then
515  ! u35_1=WNDPAX+WNDPEX
516  ! v35_1=WNDPAY+WNDPEY
517  ! elseif (wi.eq.2) then
518  ! u35_2=WNDPAX+WNDPEX
519  ! v35_2=WNDPAY+WNDPEY
520  ! elseif (wi.eq.1) then
521  ! u35_3=WNDPAX+WNDPEX
522  ! v35_3=WNDPAY+WNDPEY
523  ! endif
524  !ENDIF
525  ENDDO
526  !IF (HeightFLG) THEN
527  ! DIFU10xx= u35_3-u35_1
528  ! DIFU10yx= v35_3-v35_1
529  ! DIFU10xy= u35_2-u35_1
530  ! DIFU10yy= v35_2-v35_1
531  ! fd_a = difu10xx / DT
532  ! fd_b = difu10xy / DT
533  ! fd_c = difu10yx / DT
534  ! fd_d = difu10yy / DT
535  ! du = -wndx+u35_1
536  ! dv = -wndy+v35_1
537  ! uitv= abs(du)
538  ! vitv=abs(dv)
539  ! ch=sqrt(uitv*uitv+vitv*vitv)
540  ! if (ch.gt.10) then
541  ! apar=0.5/(fd_a*fd_d-fd_b*fd_c)
542  ! else
543  ! apar=1.0/(fd_a*fd_d-fd_b*fd_c)
544  ! endif
545 
546  ! Check for wind convergence
547  ! WND_IT_FLG = (((du**2+dv**2)/(wndx**2+wndy**2)).GT.0.001)
548  !
549  ! IF ( WND_IT_FLG .AND. WI_COUNT.LT.10 ) THEN
550  ! ! New guesses
551  ! wnd_10_x=wnd_10_x-apar*( FD_D * DU - FD_B * DV )
552  ! wnd_10_y=wnd_10_y-apar*( -FD_C * DU +FD_A * DV )
553  ! ELSE
554  ! wiflg = .FALSE.
555  ! IF (WI_COUNT .gt. 10 .AND. &
556  ! ((du**2+dv**2)/(wndx**2+wndy**2)).GT.0.05 ) then
557  ! print*,'Attn: W3FLD2 Wind error gt 5%'
558  ! !print*,' Wind y/error: ',wndy,dv
559  ! !print*,' Wind x/error: ',wndx,du
560  ! NO_ERR = .false.
561  ! endif
562  ! ENDIF
563  ! WI_COUNT = WI_COUNT + 1
564  !ELSE
565  wiflg=.false. ! If already 10 m wind then complete.
566  !ENDIF
567  enddo
568 
569  ust = ustrb
570  ustd = atan2(tauy, taux)
571  cd = ust**2 / uref**2
572  IF (PRESENT(charn)) THEN
573  charn = 0.01/sqrt(sqrt( taunux**2 + taunuy**2 )/(ust**2))
574  ENDIF
575  IF (.not.((cd .LT. 0.01).AND.(cd .GT. 0.0005)).or. .not.(no_err)) THEN
576  !Fail safe to bulk
577  print*,'Attn: W3FLD2 failed, using bulk stress'
578  print*,'Calculated Wind/Cd: ',uref,cd,ust
579  call wnd2z0m(uref,z0)
580  ust = uref * kappa / log(zwnd/z0)
581  cd = ust**2/uref**2
582  ustd = urefd
583  ENDIF
584  DEALLOCATE( tauinty , tauintx , &
585  spc2, sig2, cp , dwn , wn )
586  !STOP
587 
588  !/ End of W3FLD2 ----------------------------------------------------- /
589  !/
590  RETURN
591  !

References w3fld1md::appendtail(), w3gdatmd::dth, constants::dwat, w3servmd::extcde(), constants::grav, w3fld1md::infld(), constants::kappa, w3odatmd::ndse, w3gdatmd::nk, w3gdatmd::nspec, w3gdatmd::nth, constants::pi, w3gdatmd::sig, w3fld1md::sig2wn(), w3servmd::strace(), w3fld1md::tail_choice, w3fld1md::tail_level, w3fld1md::tail_transition_ratio1, w3fld1md::tail_transition_ratio2, w3gdatmd::th, constants::tpi, wnd2sat(), w3fld1md::wnd2z0m(), and w3gdatmd::xfr.

Referenced by w3srcemd::w3srce().

◆ wnd2sat()

subroutine w3fld2md::wnd2sat ( real, intent(in)  WND10,
real, intent(out)  SAT 
)

Definition at line 595 of file w3fld2md.F90.

595  !/ +-----------------------------------+
596  !/ | WAVEWATCH III NOAA/NCEP |
597  !/ | B. G. Reichl |
598  !/ | FORTRAN 90 |
599  !/ | Last update : 04-Aug-2016 |
600  !/ +-----------------------------------+
601  !/
602  !/ 15-Jan-2016 : Origination. ( version 5.12 )
603  !/ 05-Aug-2016 : Updated for 2016 HWRF CD/U10 curve ( J. Meixner )
604  !/
605  ! 1. Purpose :
606  !
607  ! Gives level of saturation spectrum to produce
608  ! equivalent Cd as in wnd2z0m (for neutral 10m wind)
609  ! tuned for method of Donelan et al. 2012
610  !
611  ! 2. Method :
612  !
613  ! 3. Parameters :
614  !
615  ! Parameter list
616  ! ----------------------------------------------------------------
617  !Input: WND10 - 10 m neutral wind [m/s]
618  !Output: SAT - Level 1-d saturation spectrum in tail [non-dim]
619  ! ----------------------------------------------------------------
620  ! 4. Subroutines used :
621  !
622  ! Name Type Module Description
623  ! ----------------------------------------------------------------
624  ! STRACE Subr. W3SERVMD Subroutine tracing.
625  ! ----------------------------------------------------------------
626  !
627  ! 5. Called by :
628  !
629  ! Name Type Module Description
630  ! ----------------------------------------------------------------
631  ! W3FLD2 Subr. W3FLD2MD Corresponding source term.
632  ! ----------------------------------------------------------------
633  !
634  ! 6. Error messages :
635  !
636  ! None.
637  !
638  ! 7. Remarks :
639  !
640  ! 8. Structure :
641  !
642  ! See source code.
643  !
644  ! 9. Switches :
645  !
646  ! !/S Enable subroutine tracing.
647  !
648  ! 10. Source code :
649  !
650  !/ ------------------------------------------------------------------- /
651 #ifdef W3_S
652  USE w3servmd, ONLY: strace
653 #endif
654  !
655  IMPLICIT NONE
656  REAL, INTENT(IN) :: WND10
657  REAL, INTENT(OUT) :: SAT
658 #ifdef W3_S
659  INTEGER, SAVE :: IENT = 0
660  CALL strace (ient, 'WND2SAT')
661 #endif
662  !
663  ! ST2, previous HWRF relationship:
664  ! SAT =0.000000000000349* WND10**6 +&
665  ! -0.000000000250547* WND10**5 +&
666  ! 0.000000039543565* WND10**4 +&
667  ! -0.000002229206185* WND10**3 +&
668  ! 0.000034922624204* WND10**2 +&
669  ! 0.000339117617027* WND10**1 +&
670  ! 0.003521314236550
671  ! SAT values based on
672  ! HWRF 2016 CD curve, created with fetch limited cases ST4 physics
673  IF (wnd10<20) THEN
674  sat = -0.022919753482426e-3* wnd10**2 &
675  +0.960758623686446e-3* wnd10 &
676  -0.084461041915030e-3
677  ELSEIF (wnd10<45) THEN
678  sat = -0.000000006585745* wnd10**4 &
679  +0.000001058147954* wnd10**3 &
680  -0.000065829151883* wnd10**2 &
681  +0.001587028483595* wnd10 &
682  -0.002857112191889
683  ELSE
684  sat = -0.000178498197241* wnd10 &
685  +0.012706067280674
686  ENDIF
687  sat = min(max(sat,0.002),0.014)
688 

References w3servmd::strace().

Referenced by w3fld2().

w3gdatmd::nk
integer, pointer nk
Definition: w3gdatmd.F90:1230
constants::pi
real, parameter pi
PI Value of Pi.
Definition: constants.F90:71
w3adatmd::charn
real, dimension(:), pointer charn
Definition: w3adatmd.F90:603
w3gdatmd::dth
real, pointer dth
Definition: w3gdatmd.F90:1232
w3gdatmd::nspec
integer, pointer nspec
Definition: w3gdatmd.F90:1230
w3fld1md::tail_choice
integer, save tail_choice
Definition: w3fld1md.F90:77
w3gdatmd::sig
real, dimension(:), pointer sig
Definition: w3gdatmd.F90:1234
w3fld1md::tail_transition_ratio2
real, save tail_transition_ratio2
Definition: w3fld1md.F90:84
w3fld1md
Definition: w3fld1md.F90:2
m_constants::nu
real nu
kinematic viscosity of water
Definition: mod_constants.f90:22
w3gdatmd::th
real, dimension(:), pointer th
Definition: w3gdatmd.F90:1234
w3fld1md::sig2wn
subroutine sig2wn(SIG, DEPTH, WN)
Definition: w3fld1md.F90:1184
w3odatmd::ndse
integer, pointer ndse
Definition: w3odatmd.F90:456
w3fld1md::tail_transition_ratio1
real, save tail_transition_ratio1
Definition: w3fld1md.F90:82
constants::kappa
real, parameter kappa
KAPPA von Karman's constant (N.D.).
Definition: constants.F90:69
w3servmd
Definition: w3servmd.F90:3
w3gdatmd::nth
integer, pointer nth
Definition: w3gdatmd.F90:1230
w3odatmd
Definition: w3odatmd.F90:3
constants::dwat
real, parameter dwat
DWAT Density of water (kg/m3).
Definition: constants.F90:62
constants::tpi
real, parameter tpi
TPI 2*Pi.
Definition: constants.F90:72
w3servmd::strace
subroutine strace(IENT, SNAME)
Definition: w3servmd.F90:148
w3gdatmd::xfr
real, pointer xfr
Definition: w3gdatmd.F90:1232
constants
Define some much-used constants for global use (all defined as PARAMETER).
Definition: constants.F90:20
w3gdatmd
Definition: w3gdatmd.F90:16
w3fld1md::appendtail
subroutine appendtail(INSPC, WN2, NKT, KA1, KA2, KA3, WNDDIR, SAT)
Definition: w3fld1md.F90:886
w3servmd::extcde
subroutine extcde(IEXIT, UNIT, MSG, FILE, LINE, COMM)
Definition: w3servmd.F90:736
w3fld1md::tail_level
real, save tail_level
Definition: w3fld1md.F90:81
constants::grav
real, parameter grav
GRAV Acc.
Definition: constants.F90:61
w3fld1md::wnd2z0m
subroutine wnd2z0m(W10M, ZNOTM)
Definition: w3fld1md.F90:1261
w3fld1md::infld
subroutine infld
Definition: w3fld1md.F90:787