WAVEWATCH III  beta 0.0.1
w3sic1md Module Reference

Calculate ice source term S_{ice} according to simple methods. More...

Functions/Subroutines

subroutine, public w3sic1 (A, DEPTH, CG, IX, IY, S, D)
 S_{ice} source term using 5 parameters read from input files. More...
 

Detailed Description

Calculate ice source term S_{ice} according to simple methods.

Author
E. Rogers
S. Zieger
Date
11-Oct-2013

Function/Subroutine Documentation

◆ w3sic1()

subroutine, public w3sic1md::w3sic1 ( real, dimension(nspec), intent(in)  A,
real, intent(in)  DEPTH,
real, dimension(nk), intent(in)  CG,
integer, intent(in)  IX,
integer, intent(in)  IY,
real, dimension(nspec), intent(out)  S,
real, dimension(nspec), intent(out)  D 
)

S_{ice} source term using 5 parameters read from input files.

Parameters
[in]AAction density spectrum (1-D).
[in]DEPTHLocal water depth.
[in]CGGroup velocities.
[in]IXGrid index.
[in]IYGrid index.
[out]SSource term (1-D version).
[out]DDiagonal term of derivative (1-D version).
Author
E. Rogers
S. Zieger
Date
11-Oct-2013

Definition at line 94 of file w3sic1md.F90.

94  !/
95  !/ +-----------------------------------+
96  !/ | WAVEWATCH III NOAA/NCEP |
97  !/ | E. Rogers |
98  !/ | S. Zieger |
99  !/ | FORTRAN 90 |
100  !/ | Last update : 11-Oct-2013 |
101  !/ +-----------------------------------+
102  !/
103  !/ 16-Oct-2012 : Origination. ( version 4.04 )
104  !/ (E. Rogers)
105  !/ 09-Oct-2013 : W3SIC1 SUBTYPE=2 outsourced to W3SIC2 (S. Zieger)
106  !/
107  !/ FIXME : Move field input to W3SRCE and provide
108  !/ (S.Zieger) input parameter to W3SIC1 to make the subroutine
109  !/ : versatile for point output processors ww3_outp
110  !/ and ww3_ounp.
111  !/
112  !/ Copyright 2009 National Weather Service (NWS),
113  !/ National Oceanic and Atmospheric Administration. All rights
114  !/ reserved. WAVEWATCH III is a trademark of the NWS.
115  !/ No unauthorized use without permission.
116  !/
117  ! 1. Purpose :
118  !
119  ! S_{ice} source term using 5 parameters read from input files.
120  ! These parameters are allowed to vary in space and time.
121  ! The parameters control the exponential decay rate k_i
122  ! Since there are 5 parameters, this permits description of
123  ! dependence of k_i on frequency or wavenumber.
124  !
125  !/ ------------------------------------------------------------------- /
126  !
127  ! 2. Method :
128  !
129  ! Regarding i/o (general to all Sice modules): S_{ice} source term
130  ! is calculated using up to 5 parameters read from input files.
131  ! These parameters are allowed to vary in space and time.
132  ! The parameters control the exponential decay rate k_i
133  ! Since there are 5 parameters, this permits description of
134  ! dependence of k_i on frequency or wavenumber.
135  !
136  ! Sea ice affects the wavenumber k of wind-generated ocean waves.
137  ! The ice-modified wavenumber can be expressed as a complex number
138  ! k = k_r + i*k_i, with the real part k_r representing impact of
139  ! the sea ice on the physical wavelength and propagation speeds,
140  ! producing something analogous to shoaling and refraction by
141  ! bathymetry, whereas the imaginary part of the complex
142  ! wavenumber, k_i, is an exponential decay coefficient
143  ! k_i(x,y,t,sigma) (depending on location, time and frequency,
144  ! respectively), representing wave attenuation, and can be
145  ! introduced in a wave model such as WW3 as S_ice/E=-2*Cg*k_i,
146  ! where S_ice is one of several dissipation mechanisms, along
147  ! with whitecapping, for example, S_ds=S_wc+S_ice+⋯. The k_r -
148  ! modified by ice would enter the model via the C calculations
149  ! on the left-hand side of the governing equation.The fundamentals
150  ! are straightforward, e.g. Rogers and Holland (2009 and
151  ! subsequent unpublished work) modified a similar model, SWAN
152  ! (Booij et al. 1999) to include the effects of a viscous mud
153  ! layer using the same approach (k = k_r + i*k_i) previously.
154  !
155  ! General approach is analogous to Rogers and Holland (2009)
156  ! approach for mud.
157  ! See text near their eq. 1 :
158  ! k = k_r + i * k_i
159  ! eta(x,t) = Real( a * exp( i * ( k * x - sigma * t ) ) )
160  ! a = a0 * exp( -k_i * x )
161  ! S / E = -2 * Cg * k_i (see also Komen et al. (1994, pg. 170)
162  !
163  ! Following W3SBT1 as a guide, equation 1 of W3SBT1 says:
164  ! S = D * E
165  ! However, the code of W3SBT1 has
166  ! S = D * A
167  ! This leads me to believe that the calling routine is
168  ! expecting "S/sigma" not "S"
169  ! Thus we will use D = S/E = -2 * Cg * k_i
170  !
171  ! Notes regarding numerics:
172  ! Experiments with constant k_i values suggest that :
173  ! for dx=20.0 km, k_i should not exceed 3.5e-6
174  ! (assumes 2.7% Hs error in my particular test case is intolerable)
175  ! for dx=5.0 km, k_i should not exceed 2.0e-5
176  ! for dx=2.5 km, k_i should not exceed 5.0e-5
177  ! for dx=1.0 km, k_i should not exceed 2.0e-4
178  ! for dx=0.35 km, error is less than 2.1% for all k_i tested
179  ! for dx=0.10 km, error is less than 1.3% for all k_i tested
180  ! "Ground truth" used for this is an exponential decay profile.
181  !
182  ! For reference, ACNFS is 1/12th deg, so delta_latitude=9.25 km.
183  !
184  ! {put more equations here}
185  !
186  ! 3. Parameters :
187  !
188  ! Parameter list
189  ! ----------------------------------------------------------------
190  ! A R.A. I Action density spectrum (1-D)
191  ! DEPTH Real I Local water depth
192  ! CG R.A. I Group velocities.
193  ! IX,IY I.S. I Grid indices.
194  ! S R.A. O Source term (1-D version).
195  ! D R.A. O Diagonal term of derivative (1-D version).
196  ! ----------------------------------------------------------------
197  !
198  ! 4. Subroutines used :
199  !
200  ! Name Type Module Description
201  ! ----------------------------------------------------------------
202  ! STRACE Subr. W3SERVMD Subroutine tracing (!/S switch).
203  ! PRT2DS Subr. W3ARRYMD Print plot output (!/T1 switch).
204  ! OUTMAT Subr. W3ARRYMD Matrix output (!/T2 switch).
205  ! ----------------------------------------------------------------
206  !
207  ! 5. Called by :
208  !
209  ! Name Type Module Description
210  ! ----------------------------------------------------------------
211  ! W3SRCE Subr. W3SRCEMD Source term integration.
212  ! W3EXPO Subr. N/A ASCII Point output post-processor.
213  ! W3EXNC Subr. N/A NetCDF Point output post-processor.
214  ! GXEXPO Subr. N/A GrADS point output post-processor.
215  ! ----------------------------------------------------------------
216  !
217  ! 6. Error messages :
218  !
219  ! None.
220  !
221  ! 7. Remarks :
222  !
223  ! If ice parameter 1 is zero, no calculations are made.
224  !
225  ! 8. Structure :
226  !
227  ! See source code.
228  !
229  ! 9. Switches :
230  !
231  ! !/S Enable subroutine tracing.
232  ! !/T Enable general test output.
233  ! !/T0 2-D print plot of source term.
234  ! !/T1 Print arrays.
235  !
236  ! 10. Source code :
237  !
238  !/ ------------------------------------------------------------------- /
239  USE constants, ONLY: tpi
240  USE w3odatmd, ONLY: ndse
241  USE w3servmd, ONLY: extcde
242  USE w3gdatmd, ONLY: nk, nth, nspec, sig, mapwn
243  USE w3idatmd, ONLY: icep1, icep2, icep3, icep4, icep5, inflags2
244 #ifdef W3_T
245  USE w3odatmd, ONLY: ndst
246 #endif
247 #ifdef W3_S
248  USE w3servmd, ONLY: strace
249 #endif
250 #ifdef W3_T0
251  USE w3arrymd, ONLY: prt2ds
252 #endif
253 #ifdef W3_T1
254  USE w3arrymd, ONLY: outmat
255 #endif
256  !
257  IMPLICIT NONE
258  !/
259  !/ ------------------------------------------------------------------- /
260  !/ Parameter list
261  REAL, INTENT(IN) :: CG(NK), A(NSPEC), DEPTH
262  REAL, INTENT(OUT) :: S(NSPEC), D(NSPEC)
263  INTEGER, INTENT(IN) :: IX, IY
264  !/
265  !/ ------------------------------------------------------------------- /
266  !/ Local parameters
267  !/
268 #ifdef W3_S
269  INTEGER, SAVE :: IENT = 0
270 #endif
271 #ifdef W3_T0
272  INTEGER :: ITH
273  REAL :: DOUT(NK,NTH)
274 #endif
275  INTEGER :: IKTH, IK
276  REAL :: D1D(NK) !In SBT1: D1D was named "CBETA"
277  REAL :: ICECOEF1, ICECOEF2, ICECOEF3, &
278  ICECOEF4, ICECOEF5
279  REAL, ALLOCATABLE :: WN_I(:) ! exponential decay rate for amplitude
280  !/
281  !/ ------------------------------------------------------------------- /
282  !/
283 #ifdef W3_S
284  CALL strace (ient, 'W3SIC1')
285 #endif
286  !
287  ! 0. Initializations ------------------------------------------------ *
288  !
289  d = 0.0
290  !
291  ALLOCATE(wn_i(nk))
292  wn_i = 0.0
293  icecoef1 = 0.0
294  icecoef2 = 0.0
295  icecoef3 = 0.0
296  icecoef4 = 0.0
297  icecoef5 = 0.0
298  !
299  IF (.NOT.inflags2(-7))THEN
300  WRITE (ndse,1001) 'ICE PARAMETER 1'
301  CALL extcde(2)
302  ENDIF
303  !
304  icecoef1 = icep1(ix,iy)
305 
306  !
307  ! 1. No ice --------------------------------------------------------- /
308  !
309  IF ( icecoef1==0. ) THEN
310  d = 0.
311  !
312  ! 2. Ice ------------------------------------------------------------ /
313  ELSE
314  !
315  ! 2.a Set constant(s) and write test output -------------------------- /
316  !
317  ! (none)
318  !
319 #ifdef W3_T38
320  WRITE (ndst,9000) depth,icecoef1,icecoef2,icecoef3,icecoef4
321 #endif
322  !
323  ! 2.b Make calculations ---------------------------------------------- /
324  wn_i = icecoef1 ! uniform in k
325 
326  DO ik=1, nk
327  ! SBT1 has: D1D(IK) = FACTOR * MAX(0., (CG(IK)*WN(IK)/SIG(IK)-0.5) )
328  ! recall that D=S/E=-2*Cg*k_i
329  d1d(ik) = -2. * cg(ik) * wn_i(ik)
330  END DO
331  !
332  ! 2.c Fill diagional matrix
333  !
334  DO ikth=1, nspec
335  d(ikth) = d1d(mapwn(ikth))
336  END DO
337  !
338  END IF
339  !
340  s = d * a
341  !
342  ! ... Test output of arrays
343  !
344 #ifdef W3_T0
345  DO ik=1, nk
346  DO ith=1, nth
347  dout(ik,ith) = d(ith+(ik-1)*nth)
348  END DO
349  END DO
350  CALL prt2ds (ndst, nk, nk, nth, dout, sig(1:), ' ', 1., &
351  0.0, 0.001, 'Diag Sice', ' ', 'NONAME')
352 #endif
353  !
354 #ifdef W3_T1
355  CALL outmat (ndst, d, nth, nth, nk, 'diag Sice')
356 #endif
357  !
358  ! Formats
359  !
360 1001 FORMAT (/' *** WAVEWATCH III ERROR IN W3SIC1 : '/ &
361  ' ',a,' REQUIRED BUT NOT SELECTED'/)
362  !
363 #ifdef W3_T
364 9000 FORMAT (' TEST W3SIC1 : DEPTH,ICECOEF1 : ',2e10.3)
365 #endif
366  !/
367  !/ End of W3SIC1 ----------------------------------------------------- /
368  !/

References w3servmd::extcde(), w3idatmd::inflags2, w3gdatmd::mapwn, w3odatmd::ndse, w3odatmd::ndst, w3gdatmd::nk, w3gdatmd::nspec, w3gdatmd::nth, w3arrymd::outmat(), w3arrymd::prt2ds(), w3gdatmd::sig, w3servmd::strace(), and constants::tpi.

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

w3gdatmd::nk
integer, pointer nk
Definition: w3gdatmd.F90:1230
w3gdatmd::nspec
integer, pointer nspec
Definition: w3gdatmd.F90:1230
w3idatmd::inflags2
logical, dimension(:), pointer inflags2
Definition: w3idatmd.F90:260
w3gdatmd::sig
real, dimension(:), pointer sig
Definition: w3gdatmd.F90:1234
w3arrymd::outmat
subroutine outmat(NDS, A, MX, NX, NY, MNAME)
Definition: w3arrymd.F90:988
w3odatmd::ndse
integer, pointer ndse
Definition: w3odatmd.F90:456
w3servmd
Definition: w3servmd.F90:3
w3gdatmd::nth
integer, pointer nth
Definition: w3gdatmd.F90:1230
w3odatmd
Definition: w3odatmd.F90:3
constants::tpi
real, parameter tpi
TPI 2*Pi.
Definition: constants.F90:72
w3servmd::strace
subroutine strace(IENT, SNAME)
Definition: w3servmd.F90:148
w3idatmd
Define data structures to set up wave model input data for several models simultaneously.
Definition: w3idatmd.F90:16
w3gdatmd::mapwn
integer, dimension(:), pointer mapwn
Definition: w3gdatmd.F90:1231
w3arrymd
Definition: w3arrymd.F90:3
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
Definition: w3gdatmd.F90:16
w3arrymd::prt2ds
subroutine prt2ds(NDS, NFR0, NFR, NTH, E, FR, UFR, FACSP, FSC, RRCUT, PRVAR, PRUNIT, PNTNME)
Definition: w3arrymd.F90:1943