WAVEWATCH III  beta 0.0.1
w3snl2md.F90
Go to the documentation of this file.
1 
8 
9 #include "w3macros.h"
10 !/ ------------------------------------------------------------------- /
23 MODULE w3snl2md
24  !/
25  !/ +-----------------------------------+
26  !/ | WAVEWATCH III NOAA/NCEP |
27  !/ | H. L. Tolman |
28  !/ | G. Ph. van Vledder |
29  !/ | FORTRAN 90 |
30  !/ | Last update : 29-May-2009 |
31  !/ +-----------------------------------+
32  !/
33  !/ 14-Feb-2000 : Origination. ( version 2.01 )
34  !/ 02-Feb-2001 : Exact-NL version 3.0 ( version 2.07 )
35  !/ 26-Aug-2002 : Exact-NL version 4.0 ( version 2.22 )
36  !/ 11-Nov-2002 : Interface fix. ( version 3.00 )
37  !/ 25-Sep-2003 : Exact-NL version 5.0 ( version 3.05 )
38  !/ 24-Dec-2004 : Multiple model version. ( version 3.06 )
39  !/ 29-May-2009 : Preparing distribution version. ( version 3.14 )
40  !/
41  !/ Copyright 2009 National Weather Service (NWS),
42  !/ National Oceanic and Atmospheric Administration. All rights
43  !/ reserved. WAVEWATCH III is a trademark of the NWS.
44  !/ No unauthorized use without permission.
45  !/
46  ! 1. Purpose :
47  !
48  ! Interface module to exact nonlinear interactions.
49  !
50  ! 2. Variables and types :
51  !
52  ! 3. Subroutines and functions :
53  !
54  ! Name Type Scope Description
55  ! ----------------------------------------------------------------
56  ! W3SNL2 Subr. Public Interface to Xnl calculation routines.
57  ! INSNL2 Subr. Public Initialization of Xnl routines.
58  ! ----------------------------------------------------------------
59  !
60  ! 4. Subroutines and functions used :
61  !
62  ! See subroutine.
63  !
64  ! 5. Remarks :
65  !
66  ! 6. Switches :
67  !
68  ! !/S Enable subroutine tracing.
69  ! !/T Enable general test output.
70  ! !/T0 2-D print plot of source term.
71  ! !/T1 Print arrays.
72  !
73  ! 7. Source code :
74  !/
75  !/ ------------------------------------------------------------------- /
76  !/
77  PUBLIC
78  !/
79 CONTAINS
80 !/ ------------------------------------------------------------------- /
95  SUBROUTINE w3snl2 ( A, CG, DEPTH, S, D )
96  !/
97  !/ +-----------------------------------+
98  !/ | WAVEWATCH III NOAA/NCEP |
99  !/ | H. L. Tolman |
100  !/ | G. Ph. van Vledder |
101  !/ | FORTRAN 90 |
102  !/ | Last update : 24-Dec-2004 |
103  !/ +-----------------------------------+
104  !/
105  !/ 14-Feb-2000 : Origination ( version 2.01 )
106  !/ 02-Feb-2001 : Exact-NL version 3.0 ( version 2.07 )
107  !/ 26-Aug-2002 : Exact-NL version 4.0 ( version 2.22 )
108  !/ 11-Nov-2002 : Interface fix ( version 3.00 )
109  !/ 25-Sep-2003 : Exact-NL version 5.0 ( version 3.05 )
110  !/ 24-Dec-2004 : Multiple model version. ( version 3.06 )
111  !/
112  ! 1. Purpose :
113  !
114  ! Interface to exact interactions
115  !
116  ! 2. Method :
117  !
118  ! 3. Parameters :
119  !
120  ! Parameter list
121  ! ----------------------------------------------------------------
122  ! A R.A. I Action spectrum A(ITH,IK) as a function of
123  ! direction (rad) and wavenumber.
124  ! CG R.A. I Group velocities (dimension NK).
125  ! DEPTH Real I Water depth in meters.
126  ! S R.A. O Source term.
127  ! D R.A. O Diagonal term of derivative.
128  ! ----------------------------------------------------------------
129  !
130  ! 4. Subroutines used :
131  !
132  ! Name Type Module Description
133  ! ----------------------------------------------------------------
134  ! STRACE Subr. W3SERVMD Subroutine tracing.
135  ! xnl_main Subr. m_xnldata Main Xnl routine.
136  ! ----------------------------------------------------------------
137  !
138  ! 5. Called by :
139  !
140  ! Name Type Module Description
141  ! ----------------------------------------------------------------
142  ! W3SRCE Subr. W3SRCEMD Source term integration.
143  ! W3EXPO Subr. N/A Point output post-processor.
144  ! GXEXPO Subr. N/A GrADS point output post-processor.
145  ! ----------------------------------------------------------------
146  !
147  ! 6. Error messages :
148  !
149  ! None.
150  !
151  ! 7. Remarks :
152  !
153  ! - The following settings are hardwired into the xnl_init routine
154  ! of Gerbrant van Vledder.
155  !
156  ! iufind = 0
157  ! iq_prt = 0
158  ! iq_test = 0
159  ! iq_trace = 0
160  ! iq_log = 0
161  !
162  ! 8. Structure :
163  !
164  ! See source code.
165  !
166  ! 9. Switches :
167  !
168  ! !/S Enable subroutine tracing.
169  ! !/T Enable general test output.
170  ! !/T0 2-D print plot of source term.
171  ! !/T1 Print arrays.
172  !
173  ! 10. Source code :
174  !
175  !/ ------------------------------------------------------------------- /
176  !/
177  USE constants
178  USE w3gdatmd, ONLY: nk, nth, sig, th, iqtpe
179  USE w3odatmd, ONLY: ndse, ndst, iaproc, naperr
180  USE w3servmd, ONLY: extcde
181 #ifdef W3_S
182  USE w3servmd, ONLY: strace
183 #endif
184 #ifdef W3_T0
185  USE w3arrymd, ONLY: prt2ds
186 #endif
187 #ifdef W3_T1
188  USE w3arrymd, ONLY: outmat
189 #endif
190  USE m_xnldata, ONLY: xnl_main
191  !/
192  IMPLICIT NONE
193  !/
194  !/ ------------------------------------------------------------------- /
195  !/ Parameter list
196  !/
197  REAL, INTENT(IN) :: A(NTH,NK), CG(NK), DEPTH
198  REAL, INTENT(OUT) :: S(NTH,NK), D(NTH,NK)
199  !/
200  !/ ------------------------------------------------------------------- /
201  !/ Local parameters
202  !/
203  INTEGER :: IK, ITH, IERR = 0
204 #ifdef W3_S
205  INTEGER, SAVE :: IENT = 0
206 #endif
207  REAL :: A2(NK,NTH), S2(NK,NTH), D2(NK,NTH)
208 #ifdef W3_T0
209  REAL :: SOUT(NK,NK), DOUT(NK,NK)
210 #endif
211  !/
212  !/ ------------------------------------------------------------------- /
213  !/
214 #ifdef W3_S
215  CALL strace (ient, 'W3SNL2')
216 #endif
217 #ifdef W3_T
218  WRITE (ndst,9000) iqtpe
219 #endif
220  !
221  ! 1. Convert input spectrum ----------------------------------------- *
222  ! (Action sigma spectrum, reversed indices)
223  !
224  DO ik=1, nk
225  DO ith=1, nth
226  a2(ik,ith) = a(ith,ik) / cg(ik)
227  END DO
228  END DO
229  !
230  ! 2. Call exact interaction routines -------------------------------- *
231  !
232  CALL xnl_main ( a2, sig(1:nk), th, nk, nth, depth, iqtpe, &
233  s2, d2, iaproc, ierr )
234  !
235  IF ( ierr .NE. 0 ) GOTO 800
236  !
237  ! 3. Pack results in proper format ---------------------------------- *
238  !
239  DO ik=1, nk
240  DO ith=1, nth
241  s(ith,ik) = s2(ik,ith) * cg(ik)
242  d(ith,ik) = d2(ik,ith)
243  END DO
244  END DO
245  !
246  ! ... Test output :
247  !
248 #ifdef W3_T0
249  DO ik=1, nk
250  DO ith=1, nth
251  sout(ik,ith) = s(ik,ith) * tpi * sig(ik) / cg(ik)
252  dout(ik,ith) = d(ik,ith)
253  END DO
254  END DO
255  CALL prt2ds (ndst, nk, nk, nth, sout, sig(1:nk), ' ', 1., &
256  0.0, 0.001, 'Snl(f,t)', ' ', 'NONAME')
257  CALL prt2ds (ndst, nk, nk, nth, dout, sig(1:nk), ' ', 1., &
258  0.0, 0.001, 'Diag Snl', ' ', 'NONAME')
259 #endif
260  !
261 #ifdef W3_T1
262  CALL outmat (ndst, s, nth, nth, nk, 'Snl')
263  CALL outmat (ndst, d, nth, nth, nk, 'Diag Snl')
264 #endif
265  !
266  RETURN
267  !
268  ! Error escape locations
269  !
270 800 CONTINUE
271  IF ( iaproc .EQ. naperr ) WRITE (ndse,1000) ierr
272  CALL extcde ( 1 )
273  !
274  ! Format statements
275  !
276 1000 FORMAT (/' *** WAVEWATCH III ERROR IN W3SNL2 :'/ &
277  ' xnl_main RETURN CODE NON ZERO : ',i4,' ***'/)
278  !
279 #ifdef W3_T
280 9000 FORMAT (' TEST W3SNL2 : IQTPE :',i4)
281 #endif
282  !/
283  !/ End of W3SNL2 ----------------------------------------------------- /
284  !/
285  END SUBROUTINE w3snl2
286  !/ ------------------------------------------------------------------- /
294  SUBROUTINE insnl2
295  !/
296  !/ +-----------------------------------+
297  !/ | WAVEWATCH III NOAA/NCEP |
298  !/ | H. L. Tolman |
299  !/ | G. Ph. van Vledder |
300  !/ | FORTRAN 90 |
301  !/ | Last update : 24-Dec-2004 |
302  !/ +-----------------------------------+
303  !/
304  !/ 02-Feb-2001 : Origination. ( version 2.07 )
305  !/ 25-Sep-2003 : Exact-NL version 5.0 ( version 3.05 )
306  !/ 24-Dec-2004 : Multiple model version. ( version 3.06 )
307  !/
308  ! 1. Purpose :
309  !
310  ! Preprocessing for nonlinear interactions (Xnl).
311  !
312  ! 2. Method :
313  !
314  ! See Xnl documentation.
315  !
316  ! 3. Parameters :
317  !
318  ! 4. Subroutines used :
319  !
320  ! Name Type Module Description
321  ! ----------------------------------------------------------------
322  ! STRACE Subr. W3SERVMD Subroutine tracing.
323  ! init_constants
324  ! Subr. m_xnldata Xnl initialization routine.
325  ! xnl_init Subr. m_constants Xnl initialization routine.
326  ! ----------------------------------------------------------------
327  !
328  ! 5. Called by :
329  !
330  ! Name Type Module Description
331  ! ----------------------------------------------------------------
332  ! W3IOGR Subr. W3IOGRMD Model definition file management.
333  ! ----------------------------------------------------------------
334  !
335  ! 6. Error messages :
336  !
337  ! 7. Remarks :
338  !
339  ! 8. Structure :
340  !
341  ! - See source code.
342  !
343  ! 9. Switches :
344  !
345  ! !/S Enable subroutine tracing.
346  !
347  ! 10. Source code :
348  !
349  !/ ------------------------------------------------------------------- /
350  USE constants
351  USE w3gdatmd, ONLY: nk, nth, sig, th, &
353  USE w3odatmd, ONLY: ndse, ndst, iaproc, naperr
354  USE w3servmd, ONLY: extcde
355 #ifdef W3_S
356  USE w3servmd, ONLY: strace
357 #endif
358  USE m_xnldata
359  USE m_constants, ONLY: init_constants
360  !/
361  IMPLICIT NONE
362  !/
363  !/ Local parameters
364  !/
365  INTEGER :: IGRD, IERR
366 #ifdef W3_S
367  INTEGER, SAVE :: IENT = 0
368 #endif
369  REAL :: XGRAV
370  !/
371  !/ ------------------------------------------------------------------- /
372  !/
373 #ifdef W3_S
374  CALL strace (ient, 'INSNL2')
375 #endif
376  !
377  ! 1. Set necessary values : ----------------------------------------- *
378  !
379  xgrav = grav
380  igrd = 3
381  !
382 #ifdef W3_T
383  WRITE (ndst,9000) nltail, xgrav, iqtpe, igrd, ndpths
384  WRITE (ndst,9001) dpthnl
385  WRITE (ndst,9002) sig(1)*tpiinv, sig(nk)*tpiinv, &
386  th(1)*rade, th(nth)*rade
387 #endif
388  !
389  ! 2. Call initialization routines : --------------------------------- *
390  !
391  CALL init_constants
392  !
393  CALL xnl_init ( sig(1:nk), th, nk, nth, nltail, xgrav, &
394  dpthnl, ndpths, iqtpe, igrd, iaproc, ierr )
395  !
396  IF ( ierr .NE. 0 ) GOTO 800
397  !
398  RETURN
399  !
400  ! Error escape locations
401  !
402 800 CONTINUE
403  IF ( iaproc .EQ. naperr ) WRITE (ndse,1000) ierr
404  CALL extcde ( 1 )
405  !
406  ! Format statements
407  !
408 1000 FORMAT (/' *** WAVEWATCH III ERROR IN INSNL2 :'/ &
409  ' xnl_init RETURN CODE NON ZERO : ',i8/)
410  !
411 #ifdef W3_T
412 9000 FORMAT (' TEST INSNL2 : NLTAIL :',f6.1/ &
413  ' XGRAV :',f8.3/ &
414  ' IQTPE :',i4/ &
415  ' IGRD :',i4/ &
416  ' NDPTHS :',i4,' (depths follow)')
417 9001 FORMAT (' ',5e10.3)
418 9002 FORMAT (' FREQS :',2f8.3/ &
419  ' DIRS :',2f6.1)
420 #endif
421  !/
422  !/ End of INSNL2 ----------------------------------------------------- /
423  !/
424  END SUBROUTINE insnl2
425  !/
426  !/ End of module W3SNL2MD -------------------------------------------- /
427  !/
428 END MODULE w3snl2md
w3gdatmd::nk
integer, pointer nk
Definition: w3gdatmd.F90:1230
w3snl2md::insnl2
subroutine insnl2
Preprocessing for nonlinear interactions (Xnl).
Definition: w3snl2md.F90:295
w3gdatmd::sig
real, dimension(:), pointer sig
Definition: w3gdatmd.F90:1234
w3odatmd::iaproc
integer, pointer iaproc
Definition: w3odatmd.F90:457
constants::rade
real, parameter rade
RADE Conversion factor from radians to degrees.
Definition: constants.F90:76
w3gdatmd::nltail
real, pointer nltail
Definition: w3gdatmd.F90:1346
w3arrymd::outmat
subroutine outmat(NDS, A, MX, NX, NY, MNAME)
Definition: w3arrymd.F90:988
w3gdatmd::th
real, dimension(:), pointer th
Definition: w3gdatmd.F90:1234
w3odatmd::ndse
integer, pointer ndse
Definition: w3odatmd.F90:456
w3odatmd::naperr
integer, pointer naperr
Definition: w3odatmd.F90:457
w3servmd
Definition: w3servmd.F90:3
m_xnldata::xnl_init
subroutine xnl_init(sigma, dird, nsigma, ndir, pftail, x_grav, depth, ndepth, iquad, iqgrid, iproc, ierr)
Initialize coefficients, integration space, file i/o for computation nonlinear quadruplet wave-wave i...
Definition: mod_xnl4v5.f90:465
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
m_constants::init_constants
subroutine init_constants
Subroutine init_constants sets constant values.
Definition: mod_constants.f90:49
w3gdatmd::dpthnl
real, dimension(:), pointer dpthnl
Definition: w3gdatmd.F90:1353
constants::tpi
real, parameter tpi
TPI 2*Pi.
Definition: constants.F90:72
w3servmd::strace
subroutine strace(IENT, SNAME)
Definition: w3servmd.F90:148
w3snl2md::w3snl2
subroutine w3snl2(A, CG, DEPTH, S, D)
Interface to exact interactions.
Definition: w3snl2md.F90:96
w3arrymd
Definition: w3arrymd.F90:3
m_constants
Module for m_constants.
Definition: mod_constants.f90:14
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
w3gdatmd::ndpths
integer, pointer ndpths
Definition: w3gdatmd.F90:1351
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
w3snl2md
Interface module to exact nonlinear interactions.
Definition: w3snl2md.F90:23
m_xnldata
Module for computing the quadruplet interaction.
Definition: mod_xnl4v5.f90:14
m_xnldata::xnl_main
subroutine xnl_main(aspec, sigma, angle, nsig, ndir, depth, iquad, xnl, diag, iproc, ierr)
Compute nonlinear transfer for a given action density spectrum on a given sigma and direction grid (W...
Definition: mod_xnl4v5.f90:810
w3gdatmd::iqtpe
integer, pointer iqtpe
Definition: w3gdatmd.F90:1345
constants::grav
real, parameter grav
GRAV Acc.
Definition: constants.F90:61