WAVEWATCH III  beta 0.0.1
w3sis2md Module Reference

Floe-size dependant scattering of waves in the marginal ice zone. More...

Functions/Subroutines

subroutine, public insis2
 Fill tables used for scattering. More...
 
subroutine, public w3sis2 (A, DEPTH, CICE, ICEH, ICEF, ICEDMAX, IX, IY, S, D, DISSIP, WN, CG, WN_R, CG_ICE, R)
 Wave scattering in the MIZ, adapted from Dumont et al. More...
 
subroutine, public w3rpwnice (ICEH, WN_I, DAMPING, CG_I)
 NA. More...
 

Variables

integer, parameter nthick = 20
 
integer, parameter niced = 500
 
real, parameter fragility = 0.9
 
real thick1 = 0.1
 
real dthick = 0.25
 
real icedmin
 
real, dimension(nicedicedavetab
 
double precision, dimension(:,:), allocatable, public is2eigvec
 
double precision, dimension(:,:), allocatable is2scatmat
 
double precision, dimension(:), allocatable, public is2eigval
 

Detailed Description

Floe-size dependant scattering of waves in the marginal ice zone.

based on tabulated scattering coefficients for a semi-infinite ice sheet. See papers by Dumont et al. (JGR 2011) and Williams et al. (OM 2013) combined with flexural dissipation and ice break-up.

Author
F. Ardhuin
P. Nicot
C. Sevigny
G. Boutin
Date
21-Jan-2018

Function/Subroutine Documentation

◆ insis2()

subroutine, public w3sis2md::insis2

Fill tables used for scattering.

Linear interpolation.

Author
P. Nicot
F. Ardhuin
Date
21-Jan-2018

Definition at line 104 of file w3sis2md.F90.

104  !/
105  !/ +-----------------------------------+
106  !/ | WAVEWATCH III NOAA/NCEP |
107  !/ | P. Nicot & F. Ardhuin |
108  !/ | FORTRAN 90 |
109  !/ | Last update : 21-Jan-2018 |
110  !/ +-----------------------------------+
111  !/
112  !/ 01-Apr-2014 : Creation ( version 4.18 )
113  !/ 13-Dec-2015 : Adds diagonalization of scat. matrix( version 5.10 )
114  !/ 21-Jan-2018 : Implements non-isotropic example ( version 6.04 )
115  !/
116  ! 1. Purpose :
117  !
118  ! Fill tables used for scattering
119  !
120  ! 2. Method :
121  !
122  ! Linear interpolation
123  !
124  ! 3. Parameters :
125  !
126  ! See module documentation.
127  !
128  ! 4. Error messages :
129  !
130  ! - None.
131  !
132  ! 5. Called by :
133  !
134  ! - W3IOGR (initialization after reading mod_def.ww3)
135  !
136  ! 6. Subroutines used :
137  !
138  ! - None
139  !
140  ! 7. Remarks :
141  !
142  !
143  ! 8. Structure
144  !
145  ! 9. Switches :
146  !
147  ! !/S Enable subroutine tracing.
148  !
149  ! 10. Source code :
150  !
151  !/ ------------------------------------------------------------------- /
152  !/
153  USE w3gdatmd, ONLY: sig, dsip, nk, nth, is2pars, &
154  ec2, es2, esc, esin, ecos
155  USE constants, ONLY: tpi, tpiinv
156  USE w3servmd, ONLY: diagonalize
157 #ifdef W3_S
158  USE w3servmd, ONLY: strace
159 #endif
160  !
161  IMPLICIT NONE
162  !/
163  !/ ------------------------------------------------------------------- /
164  !/ Local parameters
165  !/
166  INTEGER :: I, J, K, IND, NFTAB, NROT
167  REAL :: SIS1HTABLE(20), SIS1FTABLE(25)
168  REAL :: SIS1ALPHATABLE(NTHICK,25), X
169  REAL :: SIS1ALPHATABLE2(NTHICK,25)
170 #ifdef W3_S
171  INTEGER, SAVE :: IENT = 0
172 #endif
173  !/
174  !/ ------------------------------------------------------------------- /
175  !/
176 #ifdef W3_S
177  CALL strace (ient, 'SIS2ALPHATAB')
178 #endif
179  !
180  ! -------------------------------------------------------------------- /
181  ! 1. Fills array of reflection as a function of frequency and ice thickness
182  !
183  ALLOCATE(sis2alphas(nthick,nk))
184  ALLOCATE(sis2alpha2(nthick,nk))
185  !
186  ! Table of ice thickness for which the reflection was computed
187  !
188  sis1htable = (/ 0.1, 0.35, 0.6, 0.85, 1.1, 1.35, 1.6, 1.85, 2.1, 2.35, &
189  2.6, 2.85, 3.1, 3.35, 3.6, 3.85, 4.1, 4.35, 4.6, 4.85 /)
190  nftab = 25
191 
192  !
193  ! Table of frequencies for which the reflection was computed
194  !
195  sis1ftable = (/ 0.0420, 0.04620, 0.050820, 0.0559020, 0.06149220, 0.06764142,0.0744055620000000, &
196  0.0818461182, 0.09003073002, 0.099033803022, 0.1089371833242, 0.11983090, 0.13181399, &
197  0.144995391, 0.15949493 , 0.175444423115458, 0.192988865427003, 0.212287751969704, &
198  0.233516527, 0.256868, 0.28255499787, 0.310810497658843, 0.341891547424728, 0.376080702167200, 0.413688772383920 /)
199  IF (is2pars(18).LT.0.5) THEN
200  sis1alphatable = reshape((/ &
201  1.78e-007, 2.21e-006, 6.57e-006, 1.31e-005, 2.28e-005, 3.60e-005, 5.23e-005, 7.19e-005, 9.60e-005, 0.0001260665, &
202  0.0001621645, 0.0002032202, 0.0002483457, 0.0002987263, 0.0003571903, &
203  0.0004247358, 0.00049714, 0.0005689017, 0.0006469729, 0.0007470985, &
204  2.74e-007, 3.40e-006, 1.04e-005, 2.10e-005, 3.67e-005, 5.82e-005, &
205  8.52e-005, 0.0001179843, 0.0001581719, 0.0002075378, &
206  0.0002663323, 0.0003330095, 0.0004063158, 0.000487845, 0.0005815828, &
207  0.0006893214, 0.0008061048, 0.0009252626, 0.0010565619, 0.0012218782, &
208  4.16e-007, 5.21e-006, 1.62e-005, 3.32e-005, 5.91e-005, 9.56e-005, 0.0001423817, &
209  0.0002001783, 0.0002716622, 0.0003599526, 0.0004656353, &
210  0.0005861828, 0.0007194785, 0.0008683691, 0.0010401368, 0.0012386623, 0.0014562639, &
211  0.0016819935, 0.0019335276, 0.0022507523, &
212  6.25e-007, 7.95e-006, 2.50e-005, 5.28e-005, 9.71e-005, 0.0001615969, 0.0002467311, &
213  0.0003540015, 0.0004891301, 0.0006588358, 0.0008648602, &
214  0.0011024698, 0.0013675698, 0.0016664349, 0.0020150216, 0.0024225002, 0.0028732348, &
215  0.0033437743, 0.0038725363, 0.0045480433, &
216  9.32e-007, 1.21e-005, 3.91e-005, 8.60e-005, 0.000165062, 0.0002848805, 0.0004473356, &
217  0.0006565411, 0.0009257279, 0.0012707025, 0.0016966237, &
218  0.0021939892, 0.002754291, 0.0033926056, 0.0041466541, 0.0050385513, 0.0060311467, &
219  0.0070669538, 0.0082363302, 0.0097527074, &
220  1.38e-006, 1.86e-005, 6.30e-005, 0.0001464051, 0.0002950074, 0.0005287398, 0.0008535728, &
221  0.0012801516, 0.0018391215, 0.0025674934, 0.0034790152,&
222  0.0045539812, 0.0057742581, 0.0071757314, 0.0088467499, 0.0108382999, &
223  0.0130588683, 0.0153657347, 0.0179697664, 0.0213785773, &
224  2.05e-006, 2.91e-005, 0.0001064725, 0.0002646727, 0.0005601397, 0.0010380532, &
225  0.0017148531, 0.0026159282, 0.0038087163, 0.0053744112, &
226  0.0073439386, 0.0096752268, 0.0123295445, 0.0153861547, 0.0190373907, 0.0233881203, &
227  0.0282220195, 0.0332084504, 0.0388053869, 0.0461221517, &
228  3.03e-006, 4.71e-005, 0.0001918323, 0.000513363, 0.0011329164, &
229  0.00214921, 0.0036015483, 0.0055437802, 0.0081104965, 0.0114597841, 0.0156453149, &
230  0.0205751599, 0.0261681604, 0.0325784506, 0.040175062, 0.0491364185, &
231  0.0589976693, 0.0690831851, 0.0802738154, 0.0946821772, &
232  4.53e-006, 8.02e-005, 0.0003720589, 0.001067721, 0.0024215478, &
233  0.0046335423, 0.0077820302, 0.0119609444, 0.0173982353, 0.0243482166, &
234  0.0328666235, 0.0427569128, 0.0538584883, 0.0664227317, 0.0810452442, &
235  0.0979544205, 0.1162816559, 0.1348568677, 0.1551243874, 0.1804646861, &
236  6.84e-006, 0.0001465372, 0.0007768361, 0.0023492877, 0.0053624751, &
237  0.0101694228, 0.0168827981, 0.0256087065, 0.0366354159, 0.0502611381, &
238  0.0664627541, 0.0848577616, 0.1051596929, 0.1276991046, 0.1532607056, &
239  0.1820303259, 0.2126339282, 0.2433779188, 0.276250064, 0.3157667291, &
240  1.06e-005, 0.0002908265, 0.0017246725, 0.0053239059, 0.0119340853, &
241  0.0220417832, 0.035687347, 0.0528302869, 0.073649217, 0.098327901, &
242  0.1266340957, 0.1579155214, 0.1917144685, 0.2283923343, 0.2688390962, &
243  0.3131039168, 0.359282913, 0.4052034867, 0.4533028057, 0.5090275342, &
244  1.68e-005, 0.0006285069, 0.0039687874, 0.0119764661, 0.0257390989, &
245  0.0456577018, 0.0713257381, 0.1021551099, 0.1379813612, 0.1787934975, &
246  0.2240723105, 0.2727534775, 0.3241066558, 0.3786049119, 0.4374575502, &
247  0.5006384418, 0.5653760852, 0.6286093824, 0.6937583927, 0.7684424911, &
248  2.82e-005, 0.0014590465, 0.0091291531, 0.0257006983, 0.051998859, &
249  0.0879406127, 0.1317806874, 0.181756664, 0.2375759087, 0.2995488732, &
250  0.3668698493, 0.4374715634, 0.5099661328, 0.5855233305, 0.6668652222, &
251  0.754201342, 0.8422038741, 0.925009729, 1.0096608603, 1.1112763605, &
252  5.06e-005, 0.0035218608, 0.0201532593, 0.0509283911, 0.0962607638, &
253  0.1553638291, 0.2233707746, 0.2967358894, 0.3765436929, 0.4650074273, &
254  0.5607827923, 0.6590228472, 0.756685343, 0.8572097524, 0.9677298862,&
255  1.0892101997, 1.2094355521, 1.3151711317, 1.4236588573, 1.5698199322, &
256  9.99e-005, 0.0084186275, 0.0412091805, 0.0918067318, 0.1629311719, &
257  0.2530998556, 0.3511295005, 0.4513397894, 0.5593550205, 0.6817564883, &
258  0.815653431, 0.9499632903, 1.0782538468, 1.2090808019, 1.3588213156, &
259  1.5298071338, 1.6953922745, 1.8273584757, 1.9646959304, 2.1824935098, &
260  0.0002197702, 0.0188788997, 0.076522568, 0.1517352519, 0.2567566667, &
261  0.3885311336, 0.5247441812, 0.6573850155, 0.8007012438, 0.9681737964, &
262  1.1534672992, 1.3345908371, 1.4996817772, 1.666222509, 1.8655897771, &
263  2.1018122262, 2.3245234997, 2.4819251275, 2.6492217053, 2.9647668115, &
264  0.0005376597, 0.0381101752, 0.1293683129, 0.2359311434, 0.3874916537, &
265  0.576671054, 0.7645781239, 0.9404856491, 1.1305780655, 1.356723506, &
266  1.6068566622, 1.844140891, 2.0502598744, 2.2550088512, 2.5085040228, &
267  2.815733869, 3.096422752, 3.2715199709, 3.462031842, 3.8799520309, &
268  0.0014123358, 0.0685490628, 0.203577783, 0.354894301, 0.5714749709, &
269  0.8380065625, 1.0957796855, 1.3307089678, 1.5810360789, 1.8762692254,&
270  2.197628874, 2.4940515596, 2.7425260367, 2.9852393578, 3.2886471617, &
271  3.6568112069, 3.9833963075, 4.1693103781, 4.3750188868, 4.868451483, &
272  0.0037199175, 0.1126937266, 0.3066036571, 0.5225755275, 0.8226663225, &
273  1.1819668289, 1.5247740065, 1.8332140124, 2.152798111, 2.5177593097, &
274  2.9057350396, 3.2594397982, 3.5553142762, 3.8420910788, 4.1929464318,&
275  4.6097641365, 4.975238806, 5.1850270882, 5.4153321104, 5.9516719058, &
276  0.009037082, 0.1761561171, 0.4486599083, 0.747307698, 1.1406376531, &
277  1.5966776005, 2.0299434584, 2.4186982894, 2.8096987013, 3.2399866752, &
278  3.6899433755, 4.105747098, 4.4658248514, 4.8172532951, 5.2300901821, &
279  5.7046450426, 6.1270633187, 6.3975380116, 6.6851833028, 7.2647998139, &
280  0.0191118494, 0.2660897811, 0.6359174284, 1.0316197503, 1.5313672569, &
281  2.0961171527, 2.6294074335, 3.1054178673, 3.5736395046, 4.0762156975, &
282  4.5979660806, 5.0883063988, 5.5262800229, 5.9567260476, 6.4472029911, &
283  6.9971641707, 7.4947435133, 7.841785986, 8.2023388524, 8.8466621683, &
284  0.0352335589, 0.3840118688, 0.8727819679, 1.3979077636, 2.0496249898, &
285  2.7683726012, 3.4325758193, 4.0133546648, 4.572769913, 5.1624475081, &
286  5.7663362962, 6.3283500032, 6.8264493734, 7.3114532863, 7.8576036504, &
287  8.4626783134, 9.0041860776, 9.378581111, 9.7670030485, 10.46036883, &
288  0.0588123086, 0.5407217302, 1.1866725247, 1.8744101474, 2.7029832853, &
289  3.586993579, 4.3802119773, 5.0573372156, 5.6991092734, 6.3694458919, &
290  7.0521552001, 7.6857035561, 8.2469182713, 8.7926955021, 9.4042992778, &
291  10.0780089999, 10.6798375161, 11.0993158393, 11.5359791079, 12.3054058274, &
292  0.0912925253, 0.7587053794, 1.5933132286, 2.4467380656, 3.4467575836, &
293  4.4845019618, 5.3910117536, 6.1520041613, 6.8749329856, 7.6395860278, &
294  8.4275368684, 9.1660844666, 9.8268448897, 10.4746437466, 11.2031746957, &
295  12.0078539018, 12.7326579385, 13.2484048884, 13.7868789237, 14.7164075718, &
296  0.1320035456, 1.258347597, 2.2697363962, 2.8626691529, 3.4800285532, &
297  4.1591450976, 4.8036627318, 5.4044338329, 6.0517560386, 6.8415021436, &
298  7.8080659794, 8.9257124242, 10.1640694041, 11.5476843252, 13.1494718506, &
299  14.9965430946, 16.972215146, 18.903468102, 20.9125176371, 23.3776351255 &
300  /) ,(/nthick,nftab/))
301 
302  ELSE
303  ! May be changed, but according to T. Williams, wim1 is okay from 0.25 only
304  sis1htable = (/ 0.25, 0.35, 0.6, 0.85, 1.1, 1.35, 1.6, 1.85, 2.1, 2.35, &
305  2.6, 2.85, 3.1, 3.35, 3.6, 3.85, 4.1, 4.35, 4.6, 4.85 /)
306  sis1alphatable = reshape((/ &
307  3.80373e-06 , 6.02822e-06 , 1.12121e-05 , 2.24588e-05 , &
308  3.05165e-05 , 3.85142e-05 , 5.29712e-05 , 7.55453e-05 , &
309  1.05531e-04 , 1.38952e-04 , 1.73378e-04 , 2.09421e-04 , &
310  2.47299e-04 , 2.87669e-04 , 3.34304e-04 , 3.91094e-04 , &
311  4.57710e-04 , 5.31901e-04 , 6.14490e-04 , 7.08299e-04 , &
312  4.76980e-06 , 6.95433e-06 , 1.54179e-05 , 2.46138e-05 , &
313  3.79955e-05 , 5.96725e-05 , 8.81668e-05 , 1.20881e-04 , &
314  1.57364e-04 , 1.99666e-04 , 2.50723e-04 , 3.13346e-04 , &
315  3.89907e-04 , 4.81606e-04 , 5.88164e-04 , 7.08707e-04 , &
316  8.42536e-04 , 9.88622e-04 , 1.14498e-03 , 1.30931e-03 , &
317  5.84556e-06 , 8.58738e-06 , 2.00947e-05 , 3.18410e-05 , &
318  5.67385e-05 , 9.37456e-05 , 1.37304e-04 , 1.87832e-04 , &
319  2.49218e-04 , 3.28900e-04 , 4.31988e-04 , 5.59364e-04 , &
320  7.11302e-04 , 8.86736e-04 , 1.08074e-03 , 1.28808e-03 , &
321  1.50809e-03 , 1.74230e-03 , 1.98930e-03 , 2.24607e-03 , &
322  6.90305e-06 , 1.12155e-05 , 2.59514e-05 , 5.15063e-05 , &
323  9.18902e-05 , 1.45703e-04 , 2.18304e-04 , 3.18476e-04 , &
324  4.51992e-04 , 6.21037e-04 , 8.24428e-04 , 1.05881e-03 , &
325  1.32112e-03 , 1.60970e-03 , 1.92396e-03 , 2.26484e-03 , &
326  2.63577e-03 , 3.04238e-03 , 3.49116e-03 , 3.98911e-03 , &
327  8.46004e-06 , 1.46791e-05 , 3.69837e-05 , 8.32870e-05 , &
328  1.49151e-04 , 2.46175e-04 , 3.92165e-04 , 5.95871e-04 , &
329  8.57938e-04 , 1.17280e-03 , 1.53668e-03 , 1.95109e-03 , &
330  2.42048e-03 , 2.95341e-03 , 3.56468e-03 , 4.27150e-03 , &
331  5.08843e-03 , 6.02927e-03 , 7.11153e-03 , 8.35459e-03 , &
332  1.11513e-05 , 1.95090e-05 , 5.78215e-05 , 1.32267e-04 , &
333  2.57932e-04 , 4.61051e-04 , 7.55044e-04 , 1.14193e-03 , &
334  1.62478e-03 , 2.21208e-03 , 2.92101e-03 , 3.77684e-03 , &
335  4.80845e-03 , 6.04635e-03 , 7.52256e-03 , 9.26818e-03 , &
336  1.13101e-02 , 1.36712e-02 , 1.63732e-02 , 1.94360e-02 , &
337  1.52571e-05 , 2.78128e-05 , 9.43665e-05 , 2.30636e-04 , &
338  4.90057e-04 , 9.02105e-04 , 1.48012e-03 , 2.24822e-03 , &
339  3.24754e-03 , 4.53444e-03 , 6.16979e-03 , 8.21171e-03 , &
340  1.07145e-02 , 1.37253e-02 , 1.72805e-02 , 2.14091e-02 , &
341  2.61367e-02 , 3.14840e-02 , 3.74632e-02 , 4.40798e-02 , &
342  2.13731e-05 , 4.29771e-05 , 1.64050e-04 , 4.49510e-04 , &
343  9.84576e-04 , 1.82921e-03 , 3.06425e-03 , 4.79802e-03 , &
344  7.14520e-03 , 1.02139e-02 , 1.40940e-02 , 1.88529e-02 , &
345  2.45396e-02 , 3.11875e-02 , 3.88141e-02 , 4.74255e-02 , &
346  5.70233e-02 , 6.76049e-02 , 7.91611e-02 , 9.16776e-02 , &
347  3.17688e-05 , 7.05749e-05 , 3.13449e-04 , 9.26652e-04 , &
348  2.08059e-03 , 3.98338e-03 , 6.87070e-03 , 1.09528e-02 , &
349  1.63868e-02 , 2.32726e-02 , 3.16649e-02 , 4.15851e-02 , &
350  5.30293e-02 , 6.59789e-02 , 8.04085e-02 , 9.62863e-02 , &
351  1.13575e-01 , 1.32236e-01 , 1.52233e-01 , 1.73533e-01 , &
352  5.14019e-05 , 1.24975e-04 , 6.49063e-04 , 2.01553e-03 , &
353  4.70028e-03 , 9.17913e-03 , 1.58010e-02 , 2.47576e-02 , &
354  3.61186e-02 , 4.98708e-02 , 6.59568e-02 , 8.42994e-02 , &
355  1.04811e-01 , 1.27404e-01 , 1.51996e-01 , 1.78511e-01 , &
356  2.06873e-01 , 2.37011e-01 , 2.68862e-01 , 3.02365e-01 , &
357  9.05394e-05 , 2.44433e-04 , 1.43884e-03 , 4.66981e-03 , &
358  1.09558e-02 , 2.08947e-02 , 3.46681e-02 , 5.22109e-02 , &
359  7.33489e-02 , 9.78738e-02 , 1.25579e-01 , 1.56274e-01 , &
360  1.89790e-01 , 2.25973e-01 , 2.64687e-01 , 3.05808e-01 , &
361  3.49223e-01 , 3.94830e-01 , 4.42530e-01 , 4.92227e-01 , &
362  1.75074e-04 , 5.25080e-04 , 3.39573e-03 , 1.10036e-02 , &
363  2.46312e-02 , 4.44096e-02 , 6.99511e-02 , 1.00742e-01 , &
364  1.36302e-01 , 1.76219e-01 , 2.20147e-01 , 2.67802e-01 , &
365  3.18937e-01 , 3.73339e-01 , 4.30813e-01 , 4.91183e-01 , &
366  5.54282e-01 , 6.19954e-01 , 6.88048e-01 , 7.58418e-01 , &
367  3.74908e-04 , 1.22443e-03 , 8.17252e-03 , 2.46994e-02 , &
368  5.11207e-02 , 8.63422e-02 , 1.29173e-01 , 1.78641e-01 , &
369  2.34006e-01 , 2.94691e-01 , 3.60235e-01 , 4.30258e-01 , &
370  5.04426e-01 , 5.82438e-01 , 6.64019e-01 , 7.48914e-01 , &
371  8.36872e-01 , 9.27649e-01 , 1.02102e+00 , 1.11676e+00 , &
372  8.82800e-04 , 3.03043e-03 , 1.88248e-02 , 5.09125e-02 , &
373  9.68861e-02 , 1.54064e-01 , 2.20524e-01 , 2.94949e-01 , &
374  3.76394e-01 , 4.64138e-01 , 5.57588e-01 , 6.56231e-01 , &
375  7.59605e-01 , 8.67281e-01 , 9.78849e-01 , 1.09393e+00 , &
376  1.21214e+00 , 1.33315e+00 , 1.45662e+00 , 1.58224e+00 , &
377  2.22521e-03 , 7.51676e-03 , 3.98357e-02 , 9.57376e-02 , &
378  1.68983e-01 , 2.55527e-01 , 3.52939e-01 , 4.59644e-01 , &
379  5.74490e-01 , 6.96566e-01 , 8.25086e-01 , 9.59331e-01 , &
380  1.09864e+00 , 1.24238e+00 , 1.38997e+00 , 1.54085e+00 , &
381  1.69451e+00 , 1.85049e+00 , 2.00833e+00 , 2.16765e+00 , &
382  5.66833e-03 , 1.75464e-02 , 7.66453e-02 , 1.65727e-01 , &
383  2.74951e-01 , 3.99380e-01 , 5.36249e-01 , 6.83721e-01 , &
384  8.40370e-01 , 1.00497e+00 , 1.17642e+00 , 1.35368e+00 , &
385  1.53580e+00 , 1.72190e+00 , 1.91118e+00 , 2.10291e+00 , &
386  2.29645e+00 , 2.49123e+00 , 2.68675e+00 , 2.88259e+00 , &
387  1.36539e-02 , 3.72456e-02 , 1.35140e-01 , 2.67780e-01 , &
388  4.22901e-01 , 5.95111e-01 , 7.81333e-01 , 9.79307e-01 , &
389  1.18712e+00 , 1.40304e+00 , 1.62548e+00 , 1.85300e+00 , &
390  2.08433e+00 , 2.31832e+00 , 2.55400e+00 , 2.79055e+00 , &
391  3.02727e+00 , 3.26361e+00 , 3.49910e+00 , 3.73337e+00 , &
392  2.98722e-02 , 7.15314e-02 , 2.21423e-01 , 4.09271e-01 , &
393  6.21810e-01 , 8.53369e-01 , 1.10028e+00 , 1.35949e+00 , &
394  1.62829e+00 , 1.90425e+00 , 2.18522e+00 , 2.46937e+00 , &
395  2.75516e+00 , 3.04136e+00 , 3.32695e+00 , 3.61118e+00 , &
396  3.89345e+00 , 4.17331e+00 , 4.45047e+00 , 4.72471e+00 , &
397  5.88578e-02 , 1.25612e-01 , 3.41821e-01 , 5.98292e-01 , &
398  8.81793e-01 , 1.18598e+00 , 1.50596e+00 , 1.83748e+00 , &
399  2.17675e+00 , 2.52062e+00 , 2.86649e+00 , 3.21233e+00 , &
400  3.55662e+00 , 3.89823e+00 , 4.23636e+00 , 4.57049e+00 , &
401  4.90029e+00 , 5.22559e+00 , 5.54632e+00 , 5.86250e+00 , &
402  1.05468e-01 , 2.04724e-01 , 5.03158e-01 , 8.44026e-01 , &
403  1.21416e+00 , 1.60557e+00 , 2.01139e+00 , 2.42573e+00 , &
404  2.84383e+00 , 3.26205e+00 , 3.67776e+00 , 4.08914e+00 , &
405  4.49503e+00 , 4.89472e+00 , 5.28788e+00 , 5.67440e+00 , &
406  6.05432e+00 , 6.42782e+00 , 6.79514e+00 , 7.15656e+00 , &
407  1.74545e-01 , 3.14137e-01 , 7.13063e-01 , 1.15691e+00 , &
408  1.63126e+00 , 2.12499e+00 , 2.62866e+00 , 3.13493e+00 , &
409  3.63860e+00 , 4.13623e+00 , 4.62576e+00 , 5.10610e+00 , &
410  5.57680e+00 , 6.03787e+00 , 6.48955e+00 , 6.93227e+00 , &
411  7.36652e+00 , 7.79283e+00 , 8.21174e+00 , 8.62378e+00 , &
412  2.70834e-01 , 4.59350e-01 , 9.80330e-01 , 1.54856e+00 , &
413  2.14570e+00 , 2.75617e+00 , 3.36813e+00 , 3.97376e+00 , &
414  4.56847e+00 , 5.14997e+00 , 5.71743e+00 , 6.27093e+00 , &
415  6.81101e+00 , 7.33850e+00 , 7.85432e+00 , 8.35942e+00 , &
416  8.85471e+00 , 9.34106e+00 , 9.81925e+00 , 1.02900e+01 , &
417  3.99120e-01 , 6.46435e-01 , 1.31507e+00 , 2.03136e+00 , &
418  2.76957e+00 , 3.50946e+00 , 4.23832e+00 , 4.94958e+00 , &
419  5.64069e+00 , 6.31137e+00 , 6.96253e+00 , 7.59564e+00 , &
420  8.21234e+00 , 8.81426e+00 , 9.40295e+00 , 9.97983e+00 , &
421  1.05461e+01 , 1.11030e+01 , 1.16515e+01 , 1.21924e+01 , &
422  5.64538e-01 , 8.82405e-01 , 1.72870e+00 , 2.61770e+00 , &
423  3.51364e+00 , 4.39362e+00 , 5.24695e+00 , 6.07052e+00 , &
424  6.86499e+00 , 7.63269e+00 , 8.37647e+00 , 9.09919e+00 , &
425  9.80351e+00 , 1.04918e+01 , 1.11660e+01 , 1.18280e+01 , &
426  1.24792e+01 , 1.31209e+01 , 1.37543e+01 , 1.43803e+01 , &
427  7.72889e-01 , 1.17552e+00 , 2.23335e+00 , 3.31894e+00 , &
428  4.38704e+00 , 5.41672e+00 , 6.40295e+00 , 7.34788e+00 , &
429  8.25618e+00 , 9.13289e+00 , 9.98272e+00 , 1.08097e+01 , &
430  1.16174e+01 , 1.24086e+01 , 1.31857e+01 , 1.39506e+01 , &
431  1.47049e+01 , 1.54501e+01 , 1.61872e+01 , 1.69173e+01 &
432  /) ,(/nthick,nftab/))
433  END IF
434  DO i=1,nk
435  DO j=1,nthick
436  IF (sig(i)*tpiinv.LT.sis1ftable(1)) THEN
437  sis2alphas(j,i) = sis1alphatable(j,1)
438  ELSE IF (sig(i)*tpiinv.GT. sis1ftable(nftab)) THEN
439  sis2alphas(j,i) = sis1alphatable(j,nftab)
440  ELSE
441  ind = 1
442  DO k = 1, nftab-1
443  IF (sis1ftable(k).LT.sig(i)*tpiinv) ind = k
444  END DO
445  x=(sig(i)*tpiinv-sis1ftable(ind))/(sis1ftable(ind+1)-sis1ftable(ind))
446  sis2alphas(j,i)=sis1alphatable(j,ind)*(1-x)+sis1alphatable(j,ind+1)*x
447  END IF
448  ! WRITE(998,*) I, J, SIG(I)*TPIINV,SIS1FTABLE(NFTAB), X, IND, SIS2ALPHAS(J,I),SIS1ALPHATABLE(J,NFTAB)
449  END DO
450  END DO
451  !
452 
453  !
454  sis1alphatable2 = reshape((/ &
455  0.000001693306, 0.000001694535, 0.000001709985, 0.000001718371, 0.000001718689, &
456  0.000001715291, 0.000001711765, 0.000001709562, 0.000001708573, 0.000001708047, &
457  0.000001707216, 0.000001705508, 0.000001702489, 0.000001697757, 0.000001690926, &
458  0.000001681715, 0.000001670039, 0.000001655970, 0.000001639540, 0.000001620666, &
459  0.000002345944, 0.000002373942, 0.000002352788, 0.000002336058, 0.000002340157, &
460  0.000002361578, 0.000002389818, 0.000002414910, 0.000002430787, 0.000002435510, &
461  0.000002429671, 0.000002414451, 0.000002390516, 0.000002358094, 0.000002317788, &
462  0.000002271221, 0.000002220674, 0.000002167641, 0.000002111561, 0.000002051905, &
463  0.000003223057, 0.000003267049, 0.000003224977, 0.000003204685, 0.000003222416, &
464  0.000003261346, 0.000003299552, 0.000003321243, 0.000003319410, 0.000003293875, &
465  0.000003247777, 0.000003184750, 0.000003107742, 0.000003019246, 0.000002922032, &
466  0.000002819435, 0.000002714739, 0.000002609966, 0.000002505271, 0.000002400993, &
467  0.000004390462, 0.000004433012, 0.000004403102, 0.000004405652, 0.000004433604, &
468  0.000004455373, 0.000004447247, 0.000004400212, 0.000004316137, 0.000004201667, &
469  0.000004064075, 0.000003909816, 0.000003744650, 0.000003573915, 0.000003402281, &
470  0.000003233187, 0.000003068613, 0.000002909686, 0.000002757781, 0.000002614371, &
471  0.000005932907, 0.000005957402, 0.000005982223, 0.000006017129, 0.000006011597, &
472  0.000005931834, 0.000005775734, 0.000005561133, 0.000005309860, 0.000005038831, &
473  0.000004758982, 0.000004477929, 0.000004202314, 0.000003938073, 0.000003689194, &
474  0.000003456677, 0.000003239230, 0.000003035853, 0.000002848219, 0.000002677864, &
475  0.000007962387, 0.000007965182, 0.000008068497, 0.000008081519, 0.000007918056, &
476  0.000007586175, 0.000007143591, 0.000006652694, 0.000006157274, 0.000005680186, &
477  0.000005230854, 0.000004813199, 0.000004429817, 0.000004082375, 0.000003770294, &
478  0.000003490093, 0.000003236746, 0.000003006720, 0.000002800108, 0.000002616762, &
479  0.000010630745, 0.000010631503, 0.000010754695, 0.000010555056, 0.000009993969, &
480  0.000009206201, 0.000008341755, 0.000007500204, 0.000006728370, 0.000006039825, &
481  0.000005432987, 0.000004901489, 0.000004438118, 0.000004035401, 0.000003685333, &
482  0.000003379621, 0.000003110711, 0.000002873062, 0.000002663457, 0.000002478682, &
483  0.000014147503, 0.000014179502, 0.000014068022, 0.000013255509, 0.000011955952, &
484  0.000010517681, 0.000009161184, 0.000007970672, 0.000006956434, 0.000006102184, &
485  0.000005386202, 0.000004787198, 0.000004285040, 0.000003861250, 0.000003500220, &
486  0.000003190247, 0.000002923212, 0.000002692594, 0.000002491216, 0.000002312444, &
487  0.000018803447, 0.000018843681, 0.000017888469, 0.000015846997, 0.000013469781, &
488  0.000011291688, 0.000009484334, 0.000008030183, 0.000006862685, 0.000005921675, &
489  0.000005161810, 0.000004546889, 0.000004045234, 0.000003629478, 0.000003278778, &
490  0.000002980118, 0.000002726294, 0.000002510734, 0.000002323048, 0.000002154027, &
491  0.000024999293, 0.000024767047, 0.000021866213, 0.000017907275, 0.000014286707, &
492  0.000011453674, 0.000009339677, 0.000007754652, 0.000006540404, 0.000005593636, &
493  0.000004848620, 0.000004258144, 0.000003783667, 0.000003393700, 0.000003065741, &
494  0.000002787091, 0.000002551635, 0.000002353114, 0.000002179924, 0.000002022176, &
495  0.000033276027, 0.000031812088, 0.000025413911, 0.000019081783, 0.000014360014, &
496  0.000011113013, 0.000008878307, 0.000007287536, 0.000006111165, 0.000005216269, &
497  0.000004523259, 0.000003978625, 0.000003542278, 0.000003183888, 0.000002882748, &
498  0.000002627060, 0.000002410531, 0.000002226889, 0.000002066156, 0.000001920400, &
499  0.000044336466, 0.000039340417, 0.000027862847, 0.000019248910, 0.000013854771, &
500  0.000010491793, 0.000008293280, 0.000006775936, 0.000005680752, 0.000004862822, &
501  0.000004234499, 0.000003739067, 0.000003338786, 0.000003008719, 0.000002732485, &
502  0.000002498817, 0.000002299071, 0.000002126144, 0.000001974716, 0.000001841606, &
503  0.000059030780, 0.000046157443, 0.000028771759, 0.000018580063, 0.000013047269, &
504  0.000009812168, 0.000007739586, 0.000006322672, 0.000005313968, 0.000004570787, &
505  0.000004001034, 0.000003546046, 0.000003172357, 0.000002862735, 0.000002606289, &
506  0.000002391441, 0.000002205202, 0.000002038837, 0.000001893870, 0.000001773775, &
507  0.000078242681, 0.000050872198, 0.000028192943, 0.000017445281, 0.000012194039, &
508  0.000009215299, 0.000007293336, 0.000005968851, 0.000005030058, 0.000004344187, &
509  0.000003817214, 0.000003389423, 0.000003031921, 0.000002735020, 0.000002492733, &
510  0.000002292508, 0.000002116550, 0.000001954358, 0.000001814339, 0.000001706341, &
511  0.000102553944, 0.000052628108, 0.000026654291, 0.000016235235, 0.000011448810, &
512  0.000008738110, 0.000006951009, 0.000005700488, 0.000004811678, 0.000004163305, &
513  0.000003662892, 0.000003251855, 0.000002905015, 0.000002617308, 0.000002385271, &
514  0.000002195475, 0.000002027393, 0.000001869747, 0.000001734845, 0.000001635748, &
515  0.000131528803, 0.000051648919, 0.000024847696, 0.000015214276, 0.000010847813, &
516  0.000008336757, 0.000006656766, 0.000005471182, 0.000004622914, 0.000003999733, &
517  0.000003516266, 0.000003119051, 0.000002785020, 0.000002508708, 0.000002285704, &
518  0.000002102851, 0.000001941273, 0.000001790864, 0.000001662390, 0.000001566768, &
519  0.000162640998, 0.000049040647, 0.000023269162, 0.000014460197, 0.000010345625, &
520  0.000007943173, 0.000006352063, 0.000005238097, 0.000004434583, 0.000003834908, &
521  0.000003367217, 0.000002987629, 0.000002673539, 0.000002414062, 0.000002200682, &
522  0.000002022244, 0.000001866432, 0.000001726382, 0.000001605605, 0.000001507441, &
523  0.000190593962, 0.000046018523, 0.000022036249, 0.000013884361, 0.000009879926, &
524  0.000007534440, 0.000006029990, 0.000004998037, 0.000004246561, 0.000003673278, &
525  0.000003223583, 0.000002865466, 0.000002575963, 0.000002336404, 0.000002132787, &
526  0.000001957123, 0.000001806342, 0.000001677844, 0.000001564855, 0.000001460377, &
527  0.000208777708, 0.000043254960, 0.000020953050, 0.000013327645, 0.000009433705, &
528  0.000007161759, 0.000005746018, 0.000004790209, 0.000004085598, 0.000003536085, &
529  0.000003102390, 0.000002762488, 0.000002493229, 0.000002269788, 0.000002073862, &
530  0.000001900259, 0.000001753553, 0.000001634519, 0.000001527844, 0.000001418102, &
531  0.000213643676, 0.000040801760, 0.000019820045, 0.000012719494, 0.000009042103, &
532  0.000006888517, 0.000005551025, 0.000004645394, 0.000003968666, 0.000003433279, &
533  0.000003008733, 0.000002677964, 0.000002418184, 0.000002202481, 0.000002011120, &
534  0.000001840100, 0.000001696853, 0.000001583290, 0.000001480832, 0.000001370962, &
535  0.000208202335, 0.000038538552, 0.000018774435, 0.000012168047, 0.000008723021, &
536  0.000006681551, 0.000005399851, 0.000004523630, 0.000003863808, 0.000003339344, &
537  0.000002922742, 0.000002598268, 0.000002343776, 0.000002132781, 0.000001945837, &
538  0.000001779030, 0.000001639702, 0.000001529688, 0.000001430567, 0.000001324016, &
539  0.000198860865, 0.000036596393, 0.000018094381, 0.000011782298, 0.000008423914, &
540  0.000006432368, 0.000005194822, 0.000004356699, 0.000003726655, 0.000003224518, &
541  0.000002825485, 0.000002515982, 0.000002274472, 0.000002073961, 0.000001894595, &
542  0.000001733161, 0.000001598541, 0.000001493403, 0.000001398055, 0.000001292778, &
543  0.000188086831, 0.000035005432, 0.000017449034, 0.000011378985, 0.000008128438, &
544  0.000006206255, 0.000005018705, 0.000004216819, 0.000003612753, 0.000003129463, &
545  0.000002744757, 0.000002446723, 0.000002214548, 0.000002021421, 0.000001847651, &
546  0.000001690537, 0.000001559668, 0.000001458051, 0.000001365539, 0.000001262003, &
547  0.000176803511, 0.000033479878, 0.000016713766, 0.000010962586, 0.000007891418, &
548  0.000006057305, 0.000004907363, 0.000004122061, 0.000003528516, 0.000003054204, &
549  0.000002676573, 0.000002382977, 0.000002153342, 0.000001962633, 0.000001792438, &
550  0.000001639642, 0.000001512208, 0.000001412404, 0.000001322051, 0.000001223019, &
551  0.000157501278, 0.000031061181, 0.000013781133, 0.000009601049, 0.000007947466, &
552  0.000006760411, 0.000005655833, 0.000004680347, 0.000003892156, 0.000003277589, &
553  0.000002790517, 0.000002392986, 0.000002067558, 0.000001808761, 0.000001609460, &
554  0.000001454171, 0.000001323348, 0.000001206327, 0.000001111024, 0.000001045460, &
555  0.000007782240, 0.000002755736, 0.000000005024, 0.000000101927, 0.000088216814, &
556  0.009792270789, 0.030135654570, 0.009293514348, 0.001549599316, 0.000298014001, &
557  0.000055965309, 0.000007063938, 0.000000637549, 0.000000072256, 0.000000017854, &
558  0.000000008901, 0.000000004238, 0.000000001169, 0.000000000492, 0.000000001113 /), &
559  (/nthick,nftab/))
560 
561  DO i=1,nk
562  DO j=1,nthick
563  IF (sig(i)*tpiinv.LT. sis1ftable(1)) THEN
564  sis2alpha2(j,i) = sis1alphatable2(j,1)
565  ELSE IF (sig(i)*tpiinv.GT. sis1ftable(nftab)) THEN
566  sis2alpha2(j,i) = sis1alphatable2(j,nftab)
567  ELSE
568  ind = 1
569  DO k = 1, nftab-1
570  IF (sis1ftable(k).LT.sig(i)*tpiinv) ind = k
571  END DO
572  x=(sig(i)*tpiinv-sis1ftable(ind))/(sis1ftable(ind+1)-sis1ftable(ind))
573  sis2alpha2(j,i)=sis1alphatable2(j,ind)*(1-x)+sis1alphatable2(j,ind+1)*x
574  END IF
575  END DO
576  END DO
577  !
578  ! -------------------------------------------------------------------- /
579  ! 2. Fills array of ICEDMAX to ICEDAVE
580  !
581  DO i=1,niced
582  icedavetab(i) = w3fsd_dave(is2pars(9),real(i),is2pars(8))
583  ENDDO
584  !
585  ! -------------------------------------------------------------------- /
586  ! 2. Defines and diagonalizes the scattering matrix
587  !
588  ALLOCATE(is2scatmat(nth,nth))
589  ALLOCATE(is2eigvec(nth,nth))
590  ALLOCATE(is2eigval(nth))
591  !
592  DO i=1,nth
593  DO j=1,nth
594  ! This is for isotropic back-scatter
595  ! IS2SCATMAT(I,J)=-1./DBLE(NTH)
596  ! Other example that looks like figure 12 in Masson & LeBlond
597  is2scatmat(i,j)=-1./dble(nth)*(2.*ec2(abs(i-j)+1)**2+0.8*ecos(abs(i-j)+1)**3)
598  IF (ecos(abs(i-j)+1).LT.0.001) is2scatmat(i,j)=is2scatmat(i,j)-1./dble(nth)*0.8*es2(abs(i-j)+1)
599  END DO
600  !WRITE(997,'(36G16.8)') IS2SCATMAT(I,:)
601  ! Now removes sum from diagonal to enforce energy conservation ...
602  is2scatmat(i,i)=is2scatmat(i,i)-sum(is2scatmat(i,1:nth))
603  END DO
604  CALL diagonalize(is2scatmat,is2eigval,is2eigvec,nrot)
605  DO i=1,nth
606  is2eigval(i)=max(0.d0,is2eigval(i))
607  !WRITE(994,'(36G16.8)') I,IS2EIGVAL(I)
608  !WRITE(995,'(36G16.8)') IS2EIGVEC(I,:)
609  !WRITE(996,'(36G16.8)') IS2SCATMAT(I,:)
610  END DO
611  !CLOSE(994)
612  !CLOSE(995)
613  !CLOSE(996)
614  !CLOSE(997)
615 
616  RETURN
617  !/
618  !/ End of INSIS2 ----------------------------------------------------- /
619  !/

References w3servmd::diagonalize(), w3gdatmd::dsip, w3gdatmd::ec2, w3gdatmd::ecos, w3gdatmd::es2, w3gdatmd::esc, w3gdatmd::esin, icedavetab, is2eigval, is2eigvec, w3gdatmd::is2pars, is2scatmat, niced, w3gdatmd::nk, w3gdatmd::nth, nthick, w3gdatmd::sig, w3servmd::strace(), constants::tpi, and constants::tpiinv.

Referenced by w3iogrmd::w3iogr().

◆ w3rpwnice()

subroutine, public w3sis2md::w3rpwnice ( real, intent(in)  ICEH,
real, dimension(:), intent(inout)  WN_I,
real, dimension(:), intent(inout)  DAMPING,
real, dimension(:), intent(inout)  CG_I 
)

NA.

Parameters
[in]ICEH
[in,out]WN_I
[in,out]DAMPING
[in,out]CG_I
Author
C. Sevigny
Date
27-Aug-2015

Definition at line 1100 of file w3sis2md.F90.

1100  !/
1101  !/ +-----------------------------------+
1102  !/ | WAVEWATCH III NOAA/NCEP |
1103  !/ | C. Sevigny |
1104  !/ | FORTRAN 90 |
1105  !/ | Last update : 27-Aug-2015 |
1106  !/ +-----------------------------------+
1107  !/
1108  !/ 27-Aug-2015 : Origination. ( version 5.10 )
1109  !
1110  ! 1. Purpose :
1111  !
1112  !/ ------------------------------------------------------------------- /
1113  !
1114  ! 2. Method :
1115  !
1116  ! 3. Parameters :
1117  !
1118  ! Parameter list
1119  ! ----------------------------------------------------------------
1120  ! ----------------------------------------------------------------
1121  !
1122  ! 4. Subroutines used :
1123  !
1124  ! Name Type Module Description
1125  ! ----------------------------------------------------------------
1126  ! ----------------------------------------------------------------
1127  !
1128  ! 5. Called by :
1129  !
1130  ! Name Type Module Description
1131  ! ----------------------------------------------------------------
1132  ! ----------------------------------------------------------------
1133  !
1134  ! 6. Error messages :
1135  !
1136  ! None.
1137  !
1138  ! 7. Remarks :
1139  !
1140  ! 8. Structure :
1141  !
1142  ! See source code.
1143  !
1144  ! 9. Switches :
1145  !
1146  ! 10. Source code :
1147  !
1148  !/ ------------------------------------------------------------------- /
1149  USE constants, ONLY: grav
1150  USE w3gdatmd, ONLY: nk, sig
1151 
1152  IMPLICIT NONE
1153  !/
1154  !/ ------------------------------------------------------------------- /
1155  !/ Parameter list
1156  REAL, INTENT(IN) :: ICEH
1157  REAL, INTENT(INOUT) :: WN_I(:), DAMPING(:), CG_I(:)
1158 
1159  !/
1160  !/ ------------------------------------------------------------------- /
1161  !/ Local parameters
1162  !/
1163  DOUBLE COMPLEX :: WN_ROOT, GS1
1164  INTEGER :: IK
1165  REAL :: FLEX_RIGID
1166  REAL, PARAMETER :: VISC_RP = 10
1167  REAL, PARAMETER :: DENS = 1025.0
1168  REAL, PARAMETER :: DENS_ICE = 922.5
1169  REAL, PARAMETER :: POISSON = 0.3
1170  REAL, PARAMETER :: YOUNG = 5.49e+9
1171 
1172 
1173  flex_rigid = young * iceh**3 /(12 *(1-poisson**2) )
1174  ! Guess value for roots
1175  gs1 = cmplx(sig(1)**2/grav,0.)
1176 
1177  DO ik=1,nk
1178  CALL findroots_nr(gs1,0.,wn_root,iceh,sig(ik))
1179  wn_i(ik) = real(wn_root)
1180  CALL findroots_nr(gs1,visc_rp,wn_root,iceh,sig(ik))
1181  damping(ik) = aimag(wn_root)
1182  gs1 = wn_i(ik)
1183  cg_i(ik) = (5* flex_rigid*wn_i(ik)**4 + dens*grav - dens_ice*iceh*sig(ik)**2) &
1184  /(2*sig(ik)*(dens+dens_ice*wn_i(ik)*iceh))
1185  END DO
1186 
1187  !/

References w3servmd::extcde(), fragility, constants::grav, icedmin, w3odatmd::ndse, w3gdatmd::nk, and w3gdatmd::sig.

Referenced by w3sis2().

◆ w3sis2()

subroutine, public w3sis2md::w3sis2 ( real, dimension(nspec), intent(in)  A,
real, intent(in)  DEPTH,
real, intent(in)  CICE,
real, intent(in)  ICEH,
real, intent(inout)  ICEF,
real, intent(in)  ICEDMAX,
integer, intent(in)  IX,
integer, intent(in)  IY,
real, dimension(nspec), intent(out)  S,
real, dimension(nspec), intent(out)  D,
real, dimension(nspec), intent(out)  DISSIP,
real, dimension(nk), intent(in)  WN,
real, dimension(nk), intent(in)  CG,
real, dimension(nk), intent(in)  WN_R,
real, dimension(nk), intent(in)  CG_ICE,
real, dimension(nk), intent(out)  R 
)

Wave scattering in the MIZ, adapted from Dumont et al.

This scattering routine allows the estimation of the maximum floe size and an estimate of the creep-induced dissipation. For the scattering, it is based on the normal incidence results of Kohout and Meylan which are provided in a table.

Parameters
[in]AAction density spectrum (1-D).
[in]DEPTHWater depth.
[in]CICESea ice concentration.
[in]ICEHIce thickness.
[in,out]ICEFMaximum floe size (updated).
[in]ICEDMAXMaximum floe size.
[in]IXNot used.
[in]IYNot used.
[out]SSource term (1-D version).
[out]DDiagonal part of scattering (1-D version).
[out]DISSIPDiagonal dissipation term (1-D version)
[in]WNWave number.
[in]CGGroup speed.
[in]WN_RWave number in ice.
[in]CG_ICEGroup speed in ice.
[out]RRatio of energy to wave energy without ice.
Author
P. Nicot
F. Ardhuin
G. Boutin
Date
04-May-2016

Definition at line 654 of file w3sis2md.F90.

654  !/
655  !/ +-----------------------------------+
656  !/ | WAVEWATCH III NOAA/NCEP |
657  !/ | P. Nicot & F. Ardhuin & G. Boutin |
658  !/ | FORTRAN 90 |
659  !/ | Last update : 04-May-2016 |
660  !/ +-----------------------------------+
661  !/
662  !/ 16-Mar-2014 : Origination. ( version 4.18 )
663  !/ 19-Sep-2014 : Correcting group speed factor ( version 5.03 )
664  !/ 20-Sep-2014 : Adding back-scattered energy ( version 5.03 )
665  !/ 27-Aug-2015 : Add breaking criterion, WIM1d ( version 5.05 )
666  !/ (ref. Williams, 2012)
667  !/ 02-Nov-2015 : Integration of strain over bandwidth( version 5.05 )
668  !/ 13-Jan-2016 : Changed initialization of ICEDMAX ( version 5.10 )
669  !/ 06-Feb-2016 : Added IICEHMIN and creep dissipation( version 5.10 )
670  !/ 10-Mar-2016 : Added depth and call to Liu disp. ( version 5.10 )
671  !/ 02-May-2016 : Call to Liu disp moved to w3srce ( version 5.10 )
672  !
673  ! 1. Purpose :
674  ! Wave scattering in the MIZ, adapted from Dumont et al.
675  !
676  !/ ------------------------------------------------------------------- /
677  !
678  ! 2. Method :
679  ! This scattering routine allows the estimation of the maximum floe
680  ! size and an estimate of the creep-induced dissipation.
681  ! For the scattering, it is based on the normal incidence results of
682  ! Kohout and Meylan which are provided in a table.
683  !
684  ! 3. Parameters :
685  !
686  ! Parameter list
687  ! ----------------------------------------------------------------
688  ! A R.A. I Action density spectrum (1-D)
689  ! DEPTH Real I Water depth
690  ! CICE Real I Sea ice concentration
691  ! ICEH Real I ice thickness
692  ! ICEF Real I/O Maximum floe size (updated)
693  ! ICEDMAX Real I Maximum floe size
694  ! IX,IY Int I Not used
695  ! S R.A. O Source term (1-D version)
696  ! D R.A. O Diagonal part of scattering (1-D version)
697  ! DISSIP R.A. O Diagonal dissipation term (1-D version)
698  ! WN R.A. I Wave number
699  ! CG R.A. I Group speed
700  ! WN_R R.A. I Wave number in ice
701  ! CG_ICE R.A. I Group speed in ice
702  ! R R.A. O Ratio of energy to wave energy without ice
703  ! ----------------------------------------------------------------
704  !
705  ! 4. Subroutines used :
706  !
707  ! Name Type Module Description
708  ! ----------------------------------------------------------------
709  ! ----------------------------------------------------------------
710  !
711  ! 5. Called by :
712  !
713  ! Name Type Module Description
714  ! ----------------------------------------------------------------
715  ! W3SRCE Subr. W3SRCEMD Source term integration.
716  ! W3EXPO Subr. N/A ASCII Point output post-processor.
717  ! W3EXNC Subr. N/A NetCDF Point output post-processor.
718  ! GXEXPO Subr. N/A GrADS point output post-processor.
719  ! ----------------------------------------------------------------
720  !
721  ! 6. Error messages :
722  !
723  ! None.
724  !
725  ! 7. Remarks :
726  !
727  ! If ice concentration is zero, no calculations are made.
728  !
729  ! 8. Structure :
730  !
731  ! See source code.
732  !
733  ! 9. Switches :
734  !
735  ! !/S Enable subroutine tracing.
736  ! !/T Enable general test output.
737  ! 2-D print plot of source term.
738  !
739  ! 10. Source code :
740  !
741  !/ ------------------------------------------------------------------- /
742  USE w3odatmd, ONLY: ndse
743  USE constants, ONLY: tpiinv, pi, tpi, grav, dwat
744  USE w3servmd, ONLY: extcde
745  USE w3gdatmd, ONLY: nk, nth, nspec, sig, sig2, dden, is2pars, xfr, &
747 #ifdef W3_T
748  USE w3odatmd, ONLY: ndst
749 #endif
750 #ifdef W3_S
751  USE w3servmd, ONLY: strace
752 #endif
753 #ifdef W3_T
754  USE w3arrymd, ONLY: prt2ds
755 #endif
756  USE w3dispmd
757  !
758  IMPLICIT NONE
759  !/
760  !/ ------------------------------------------------------------------- /
761  !/ Parameter list
762  REAL, INTENT(IN) :: A(NSPEC), DEPTH, CICE, ICEH, ICEDMAX
763  INTEGER, INTENT(IN) :: IX, IY
764  REAL, INTENT(OUT) :: S(NSPEC), D(NSPEC), DISSIP(NSPEC), R(NK)
765  REAL, INTENT(INOUT) :: ICEF
766  REAL, INTENT(IN) :: WN(NK), CG(NK), WN_R(NK), CG_ICE(NK)
767 
768  !/
769  !/ ------------------------------------------------------------------- /
770  !/ Local parameters
771  !/
772 #ifdef W3_S
773  INTEGER, SAVE :: IENT = 0
774 #endif
775  INTEGER :: IK, IKP1, IKM1, ITH, ITH2, IS, IS2, IND1, IND2
776  REAL :: W
777  INTEGER :: IKBREAK, NSUM
778  LOGICAL :: BRK_CRIT_W(NK), BRK_CRIT
779  REAL :: ALPHA, STRAIN_C, WAMP(NK), D_FLEX_FAIL, &
780  SMOOTHD, TAU_D, S_D(NK), ALPHA_D, DELTA_D,B_COLE, &
781  DMAX, S_ATT, FACTOR, BETA
782  REAL :: ICEDAVE(NK), CURVTOSTRAIN, CREEPFAC, MP2, B, ICEF_CREEP
783  REAL :: SUMALLDIR, SUMA, SUME, CURVSPEC(NK), ESPEC(NK),STRAIN
784  REAL, PARAMETER :: YOUNG = 5.49e+9 ! Young modulus
785  REAL, PARAMETER :: POISSON = 0.3 ! Poisson Ratio
786  REAL :: SIGMA_C
787  REAL, PARAMETER :: DENS = 1025.0 ! ice density
788  REAL :: GAMMA_TOY
789  REAL, DIMENSION(NK) :: WN_I, WN_RP, WSQ, WLG, WLG_I, CG_I, &
790  CURV, CGRATIO, CG_EFF, DUMMY, ALPHA_DISP
791 #ifdef W3_T
792  REAL :: SOUT(NK,NTH)
793 #endif
794  !/
795  !/ ------------------------------------------------------------------- /
796  !/
797 #ifdef W3_S
798  CALL strace (ient, 'W3SIS1')
799 #endif
800  !
801  ! 0. Initializations ------------------------------------------------ *
802  !
803  s = 0.
804  d = 0.
805  dissip = 0.
806  dummy = wn
807  wn_i = 0.
808  wn_rp = 0.
809  cg_i = 0.
810  cg_eff = 0.
811  sigma_c = is2pars(19)
812  mp2=(1-poisson**2)
813  gamma_toy = 2 + log(0.9)/log(2.)
814  ! Variables from Cole et al. 1995
815  alpha_d=0.54
816  b_cole=1.205e-9 * exp(is2pars(24)*1.60218e-19/(1.38064852e-23*268.15))
817  tau_d=b_cole/0.07
818  s_d=log(sig(1:nk)*tau_d)
819  delta_d=is2pars(21)
820  IF (is2pars(9).GT.0) icedmin = is2pars(9)
821  IF (is2pars(12).GT.0) THEN
822  b=is2pars(12)
823  ! 2 is the ratio of <H/2^4> to (Hrms/2)^4 for a Rayleigh distribution
824  ! 0.375 is the average of cos^4
825  ! 0.4 is 2/5
826  creepfac = -2*(0.25/(is2pars(15)+2))*0.375*b*iceh**(is2pars(15)+2) &
827  *(young/(2*b*mp2))**(is2pars(15)+1)/(dwat*grav)
828  ELSE
829  creepfac=1.
830  ENDIF
831 
832  wlg = tpi/wn
833  icef_creep=icef
834  dmax = icef
835  brk_crit = .false.
836  nsum = nint(0.3/(xfr-1.))
837  !
838  strain_c = sigma_c*mp2/young
839 
840  ! Minimum floe size that can break
841  d_flex_fail = 0.5* ( (pi**4*young*iceh**3)/( 48*dens*grav*mp2 ) )**.25
842  ! Estimates mean floe diameter from max floe diameter
843  IF (dmax.GT.niced) THEN
844  icedave=dmax
845  ELSE
846  IF (is2pars(9).GT.0) THEN
847  icedave= icedavetab(max(1,nint(dmax)))
848  ELSE
849  DO ik=1,nk
850  IF (is2pars(14)*(tpi/wn_r(ik)).LT.dmax) THEN
851  icedave(ik)=w3fsd_dave( (tpi/wn_r(ik))*is2pars(14),dmax,is2pars(8))
852  ELSE
853  icedave(ik)=dmax
854  ENDIF
855  END DO
856  ENDIF
857  ENDIF
858  !
859 #ifdef W3_T
860  sout = 0.
861 #endif
862  !
863  IF (cice .GT. 0) THEN
864  !
865  ! 1. Calculate wavelength, Robinson-Palmer dispersion relation
866  ! (should be tabulated )
867  ! Warning: it only applies to unbroken ice.
868  !
869  IF (is2pars(6).GT.0.5) THEN
870  CALL w3rpwnice(iceh,wn_i,wn_rp,cg_i)
871  ELSE
872  wn_i = wn_r
873  cg_i = cg_ice
874  wn_rp=wn*0.
875  END IF
876  !
877  wlg_i = tpi/wn_i ! Ice wavelength
878  wsq = (wlg/wlg_i)
879  !
880  IF (is2pars(16).GT.0.5) THEN
881  cgratio(:)=cg_i(:)/cg(:)
882  ELSE
883  cgratio(:)=1
884  END IF
885  !
886  ! 2. gets reflection coefficient from table
887  !
888  ind1 = 1+floor((iceh-thick1)/dthick)
889  ind2 = ind1+1
890  ! defines weight for interpolation of ice thickness
891  w = (iceh-thick1)/dthick - (ind1-1)
892  IF (ind1.LT.1) THEN
893  ind1 = 1
894  w = 0
895  ELSE IF (ind2.GT.nthick) THEN
896  ind2 = nthick
897  ind1 = nthick
898  w = 0
899  END IF
900  !
901  DO ik = 1,nk
902  !
903  ! Spatial decay scale taken from table. This corresponds the the values shown
904  ! in Dumont et al. (JGR 2011, fig 3).
905  !
906  cg_eff(ik) = cice*cg_i(ik) + (1-cice)*cg(ik)
907 
908  ! Note by FA: dissipation should be done by ICx not by ISx
909  ! the RP damping is thus defined by an optional IS2PARS(7), which is 0 by default
910  alpha = -1.*( (sis2alphas(ind1,ik)*(1-w)+sis2alphas(ind2,ik)*w)/icedave(ik) &
911  +2.*is2pars(7)*wn_rp(ik) )*cice
912  !
913  ! Additional scattering for pack ice defined by IS2PARS(4:5) (see Squire et al. GRL 2009)
914  !
915  alpha = is2pars(1) * alpha -2.*is2pars(4)*exp(-1.*is2pars(5)/sig(ik))
916  IF (is2pars(11).GT.0) THEN
917  IF (cice.LT.0.2) alpha = 0.
918  IF (cice.GT.0.8) alpha = 0.
919  IF (cice.GE.0.2.AND.cice.LE.0.8) alpha = alpha*(cice-0.2)*(0.8-cice)
920  END IF
921  !
922  ! time decay
923  !
924  beta = alpha * cice * cg_eff(ik)
925  !
926  ! 3. attenuation due to scattering for all spectral components
927  ! with added backscattering for energy conservation ( if IS2PARS(2).EQ.1)
928  !
929  sumalldir= 0.
930  suma = 0.
931  curvtostrain = (0.25*max(iceh,iicehmin)**2)
932  DO ith = 1,nth
933  is = ith+(ik-1)*nth
934  d(is) = beta
935  s(is) = beta * a(is)
936  sumalldir = sumalldir + s(is)
937  suma = suma + a(is)
938  END DO ! loop over directions
939  !
940  ! R is the ratio of energy (including bending of ice) to wave energy without ice
941  ! Wadhams 1973 eq. 34, warning, his ice thickness is 2*h
942  ! Warning : R uses DMAX=ICEF, even if IS2DUPDATE=F
943  !
944  IF (iicesmooth) THEN
945  IF (is2pars(14)*wlg_i(ik).LT.dmax) THEN
946  smoothd=tanh((dmax-is2pars(14)*wlg_i(ik))/(dmax*is2pars(13)))
947  ELSE
948  smoothd=0.
949  END IF
950  ELSE
951  smoothd=1.
952  END IF
953  !
954  r(ik) =1+is2pars(16)*smoothd*4*young*iceh**3*(pi/wlg_i(ik))**4/(3*dwat*grav*mp2)
955  !
956  ! Converting action to surface elevation variance SUME with units m^2
957  !
958  sume = suma*dden(ik) / cg(ik) / (r(ik)*cgratio(ik))
959  !
960  ! CURVSPEC is the curvature variance = elevation variance * k^4
961  !
962  curvspec(ik) = sume * (2*pi/ wlg_i(ik))**4
963  espec(ik) = sume
964  sumalldir = sumalldir / real(nth)
965  !
966  ! Adds the scattered energy isotropically to conserve the energy
967  ! This may not be a very good scheme numerically. Another possible
968  ! approach is the matrix inversion used for bottom scattering (w3sbs1md.ftn)
969  !
970  s(1+(ik-1)*nth:ik*nth)=s(1+(ik-1)*nth:ik*nth)-sumalldir*is2pars(2)
971  END DO ! loop over wavenumbers IK
972  !
973  ! 4. update of floe size
974  !
975  IF (is2pars(10).LT.0.5) THEN ! resets max floe size to the last forcing or initial value
976  dmax = icedmax
977  icef = icedmax
978  END IF
979  !
980  DO ik = 1, nk
981  ! CURV is the variance of the curvature integrated over a finite bandwidth
982  curv(ik) = sum(curvspec(max(1,ik-nsum):min(nk,ik+nsum)))
983  ! Now converts curvature variance to strain variance
984  END DO ! end of loop on IK
985  !
986  ! If IS2PARS(3)=IS2BREAK is set to true in ww3_grid, then activates ice break-up
987  !
988  IF (is2pars(3).GT.0.5) THEN
989  ikbreak=0
990  DO ik = 1, nk
991  strain = curv(ik)*curvtostrain
992  IF (d_flex_fail .LT. dmax) THEN
993  ! Note that Williams et al. used IS2PARS(17)=SQRT(2), Here our default is 3.6
994  wamp(ik)= is2pars(17)*sqrt(strain)
995  !
996  IF (is2pars(9).EQ.0) THEN
997  icedmin=(tpi/wn_r(ik))*is2pars(14)
998  END IF
999  brk_crit_w(ik) = wamp(ik) .GT. strain_c .AND. wlg_i(ik)/2 .GT. icedmin .AND. wlg_i(ik)/2 .LT. dmax &
1000  .AND. wlg_i(ik)/2 .GT. d_flex_fail
1001 
1002  !
1003  IF (brk_crit_w(ik)) THEN
1004  ikbreak=ik
1005  brk_crit = .true.
1006  END IF
1007  END IF
1008  END DO ! end of loop on IK
1009  !
1010  ! 4.b Correction for bias introduced by the finite bandwidth sum .. .
1011  !
1012  IF (brk_crit) THEN
1013  DO ik=max(ikbreak-nsum,1),ikbreak,1
1014  ! Modified by F.A. on Jan. 31, 2017: uses the maximum of CURVSPEC instead of CURV.
1015  ! this is better for very narrow spectra.
1016  IF (curvspec(ik).GE.curvspec(ikbreak).AND.dmax.GE.(wlg_i(ik)/2)) THEN
1017  ikbreak = ik
1018  END IF
1019  END DO
1020  !
1021  dmax = wlg_i(ikbreak)/2
1022  !
1023  ! Uses a weighting by CURVSPEC to have a continuous shift of DMAX ..
1024  !
1025  ikp1=min(ikbreak+1,nk)
1026  ikm1=max(ikbreak-1,1)
1027  IF (brk_crit_w(ikp1).AND.brk_crit_w(ikm1)) THEN
1028  dmax = (wlg_i(ikbreak)*curvspec(ikbreak) &
1029  +wlg_i(ikp1)*curvspec(ikp1)+wlg_i(ikm1)*curvspec(ikm1)) &
1030  /(2.*(curvspec(ikbreak)+curvspec(ikm1)+curvspec(ikp1)))
1031  END IF
1032  !
1033  icef = dmax
1034  END IF
1035  END IF !end of test (IS2PARS(3).GT.0.5)
1036  !
1037  ! 5. inelastic or anelastic dissipation
1038  !
1039  IF (is2pars(12).GT.0) THEN
1040  DO ik = 1, nk
1041  !
1042  ! The TANH((DMAX-D*WLG_I(IK))/DMAX*C)
1043  ! is an ad hoc factor that goes to zero for WLG << DMAX and 1 for WLG >> DMAX
1044  ! this should probably be adjusted.
1045  !
1046  IF (is2pars(14)*wlg_i(ik).LT.dmax) THEN
1047  smoothd=tanh((dmax-is2pars(14)*wlg_i(ik))/(dmax*is2pars(13)))
1048  IF (is2pars(23).LE.0.5) THEN ! this is the inelastic option
1049  dissip(1+(ik-1)*nth:ik*nth)=creepfac*4*curv(ik) &
1050  *((2*pi)/wlg_i(ik))**(is2pars(15)+1) &
1051  /(cgratio(ik)**1*r(ik)**2) &
1052  *smoothd
1053  ELSE ! this is the inelastic option
1054  dissip(1+(ik-1)*nth:ik*nth) =-4*4/3*sig(ik)* delta_d*alpha_d *wn_i(ik)**4 * (young/mp2)**2 &
1055  * (iceh/2)**3/3 * 1/( exp(alpha_d*s_d(ik)) + exp(-alpha_d*s_d(ik))) &
1056  * smoothd /(r(ik)**2*cgratio(ik)) / (dwat*grav) *tpiinv
1057  END IF
1058  END IF
1059  s=s+dissip*cice*a
1060  END DO ! end of loop on IK
1061  ENDIF ! end of test (IS2PARS(12).GT.0)
1062  !
1063  ! 6. Case of no scattering nor dissipation
1064  !
1065  ELSE
1066  dmax = 0.
1067  icef = 0.
1068  END IF ! end of test (CICE .GT. 0 .AND. ICEDAVE .GT. 0)
1069  !
1070 #ifdef W3_T
1071  DO ik = 1, nk
1072  DO ith = 1, nth
1073  is = ith+(ik-1)*nth
1074  sout(ik,ith) = s(is)
1075  END DO
1076  END DO
1077  CALL prt2ds (ndst, nk, nk, nth, sout, sig(1:nk), ' ', 1., &
1078  0.0, 0.001, 'Diag Sir1', ' ', 'NONAME')
1079 #endif
1080  !
1081  ! Formats
1082 8000 FORMAT (' TEST W3SIS2 : ALPHA :',e10.3)
1083  !
1084  !/
1085  !/ End of W3SIS2 ----------------------------------------------------- /
1086  !/

References w3gdatmd::dden, dthick, constants::dwat, w3servmd::extcde(), constants::grav, icedavetab, icedmin, w3gdatmd::iicehmin, w3gdatmd::iicesmooth, w3gdatmd::is2pars, w3odatmd::ndse, w3odatmd::ndst, niced, w3gdatmd::nk, w3gdatmd::nspec, w3gdatmd::nth, nthick, constants::pi, w3arrymd::prt2ds(), w3gdatmd::sig, w3gdatmd::sig2, w3servmd::strace(), thick1, constants::tpi, constants::tpiinv, w3rpwnice(), and w3gdatmd::xfr.

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

Variable Documentation

◆ dthick

real w3sis2md::dthick = 0.25

Definition at line 80 of file w3sis2md.F90.

Referenced by w3sis2().

◆ fragility

real, parameter w3sis2md::fragility = 0.9

Definition at line 79 of file w3sis2md.F90.

79  REAL, PARAMETER :: FRAGILITY = 0.9

Referenced by w3rpwnice().

◆ icedavetab

real, dimension(niced) w3sis2md::icedavetab

Definition at line 82 of file w3sis2md.F90.

82  REAL :: ICEDAVETAB(NICED)

Referenced by insis2(), and w3sis2().

◆ icedmin

real w3sis2md::icedmin

Definition at line 81 of file w3sis2md.F90.

81  REAL :: ICEDMIN ! minimum floe diameter

Referenced by w3rpwnice(), and w3sis2().

◆ is2eigval

double precision, dimension(:), allocatable, public w3sis2md::is2eigval

Definition at line 85 of file w3sis2md.F90.

85  DOUBLE PRECISION , ALLOCATABLE,DIMENSION(:) :: IS2EIGVAL

Referenced by insis2().

◆ is2eigvec

double precision, dimension(:,:), allocatable, public w3sis2md::is2eigvec

Definition at line 84 of file w3sis2md.F90.

84  DOUBLE PRECISION, ALLOCATABLE,DIMENSION(:,:) :: IS2EIGVEC, IS2SCATMAT

Referenced by insis2().

◆ is2scatmat

double precision, dimension(:,:), allocatable w3sis2md::is2scatmat

Definition at line 84 of file w3sis2md.F90.

Referenced by insis2().

◆ niced

integer, parameter w3sis2md::niced = 500

Definition at line 78 of file w3sis2md.F90.

Referenced by insis2(), and w3sis2().

◆ nthick

integer, parameter w3sis2md::nthick = 20

Definition at line 78 of file w3sis2md.F90.

78  INTEGER , PARAMETER :: NTHICK = 20, niced = 500

Referenced by insis2(), and w3sis2().

◆ thick1

real w3sis2md::thick1 = 0.1

Definition at line 80 of file w3sis2md.F90.

80  REAL :: THICK1 = 0.1, dthick = 0.25

Referenced by w3sis2().

w3gdatmd::nk
integer, pointer nk
Definition: w3gdatmd.F90:1230
w3gdatmd::esc
real, dimension(:), pointer esc
Definition: w3gdatmd.F90:1234
constants::pi
real, parameter pi
PI Value of Pi.
Definition: constants.F90:71
w3wdatmd::iceh
real, dimension(:), pointer iceh
Definition: w3wdatmd.F90:183
w3gdatmd::nspec
integer, pointer nspec
Definition: w3gdatmd.F90:1230
w3gdatmd::sig
real, dimension(:), pointer sig
Definition: w3gdatmd.F90:1234
w3wdatmd::icef
real, dimension(:), pointer icef
Definition: w3wdatmd.F90:183
w3gdatmd::ecos
real, dimension(:), pointer ecos
Definition: w3gdatmd.F90:1234
w3gdatmd::dsip
real, dimension(:), pointer dsip
Definition: w3gdatmd.F90:1234
w3odatmd::ndse
integer, pointer ndse
Definition: w3odatmd.F90:456
w3gdatmd::es2
real, dimension(:), pointer es2
Definition: w3gdatmd.F90:1234
w3gdatmd::esin
real, dimension(:), pointer esin
Definition: w3gdatmd.F90:1234
w3gdatmd::iicesmooth
logical, pointer iicesmooth
Definition: w3gdatmd.F90:1217
w3servmd
Definition: w3servmd.F90:3
w3gdatmd::iicehmin
real, pointer iicehmin
Definition: w3gdatmd.F90:1183
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
constants::dwat
real, parameter dwat
DWAT Density of water (kg/m3).
Definition: constants.F90:62
w3gdatmd::sig2
real, dimension(:), pointer sig2
Definition: w3gdatmd.F90:1234
w3servmd::diagonalize
subroutine diagonalize(a1, d, v, nrot)
Definition: w3servmd.F90:1840
constants::tpi
real, parameter tpi
TPI 2*Pi.
Definition: constants.F90:72
w3servmd::strace
subroutine strace(IENT, SNAME)
Definition: w3servmd.F90:148
w3wdatmd::icedmax
real, dimension(:), pointer icedmax
Definition: w3wdatmd.F90:183
w3arrymd
Definition: w3arrymd.F90:3
w3gdatmd::xfr
real, pointer xfr
Definition: w3gdatmd.F90:1232
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
w3servmd::extcde
subroutine extcde(IEXIT, UNIT, MSG, FILE, LINE, COMM)
Definition: w3servmd.F90:736
w3arrymd::prt2ds
subroutine prt2ds(NDS, NFR0, NFR, NTH, E, FR, UFR, FACSP, FSC, RRCUT, PRVAR, PRUNIT, PNTNME)
Definition: w3arrymd.F90:1943
w3gdatmd::is2pars
real, dimension(:), pointer is2pars
Definition: w3gdatmd.F90:1161
w3gdatmd::ec2
real, dimension(:), pointer ec2
Definition: w3gdatmd.F90:1234
w3dispmd
Definition: w3dispmd.F90:3
constants::grav
real, parameter grav
GRAV Acc.
Definition: constants.F90:61