WAVEWATCH III  beta 0.0.1
ww3_outf.F90
Go to the documentation of this file.
1 
7 
8 #include "w3macros.h"
9 !/ ------------------------------------------------------------------- /
10 
35 PROGRAM w3outf
36  !/
37  !/ +-----------------------------------+
38  !/ | WAVEWATCH III NOAA/NCEP |
39  !/ | H. L. Tolman |
40  !/ | FORTRAN 90 |
41  !/ | Last update : 22-Mar-2021 |
42  !/ +-----------------------------------+
43  !/
44  !/ 19-Oct-1998 : Final FORTRAN 77 ( version 1.18 )
45  !/ 19-Jan-2000 : Upgrade to FORTRAN 90 ( version 2.00 )
46  !/ 24-Jan-2001 : Flat grid version ( version 2.06 )
47  !/ 23-Apr-2002 : Clean-up ( version 2.19 )
48  !/ 29-Apr-2002 : Adding output fields 17-18. ( version 2.20 )
49  !/ 13-Nov-2002 : Add stress vector. ( version 3.00 )
50  !/ 24-Dec-2004 : Multiple grid version. ( version 3.06 )
51  !/ 21-Jul-2005 : Adding output fields 19-21. ( version 3.07 )
52  !/ 28-Jun-2006 : Adding file name preamble. ( version 3.09 )
53  !/ 05-Jul-2006 : Consolidate stress arrays. ( version 3.09 )
54  !/ 28-Mar-2007 : Adding partitioned output. ( version 3.11 )
55  !/ Adding user slots for outputs.
56  !/ 31-Jul-2007 : Fix file extension errors. ( version 3.12 )
57  !/ 29-May-2009 : Preparing distribution version. ( version 3.14 )
58  !/ 30-Oct-2009 : Implement run-time grid selection. ( version 3.14 )
59  !/ (W. E. Rogers & T. J. Campbell, NRL)
60  !/ 30-Oct-2009 : Implement curvilinear grid type. ( version 3.14 )
61  !/ (W. E. Rogers & T. J. Campbell, NRL)
62  !/ 12-Dec-2012 : SMC grid sea-point text output.JG_Li( version 4.08 )
63  !/ 25-Dec-2012 : New structure of output fields. ( version 4.11 )
64  !/ Minor bug fixes and clean up.
65  !/ 11-Nov-2013 : SMC and rotated grid incorporated in the main
66  !/ trunk ( version 4.13 )
67  !/ 27-Aug-2015 : ICEH and ICEF added as output ( version 5.10 )
68  !/ 12-Sep-2018 : Added new partitioned output fields ( version 6.06 )
69  !/ 26-Jan-2021 : Added TP field (derived from FP0) ( version 7.12 )
70  !/ 22-Mar-2021 : New coupling fields output ( version 7.13 )
71  !/
72  !/ Copyright 2009-2012 National Weather Service (NWS),
73  !/ National Oceanic and Atmospheric Administration. All rights
74  !/ reserved. WAVEWATCH III is a trademark of the NWS.
75  !/ No unauthorized use without permission.
76  !/
77  ! 1. Purpose :
78  !
79  ! Post-processing of grid output.
80  !
81  ! 2. Method :
82  !
83  ! Data is read from the grid output file out_grd.ww3 (raw data)
84  ! and from the file ww3_outf.inp ( NDSI, output requests ).
85  ! Model definition and raw data files are read using WAVEWATCH III
86  ! subroutines.
87  !
88  ! Output types :
89  ! 1 : print plots
90  ! 2 : field statistics
91  ! 3 : transfer file
92  ! 4 : text output at sea points (1:NSEA).
93  !
94  ! 3. Parameters :
95  !
96  ! 4. Subroutines used :
97  !
98  ! Name Type Module Description
99  ! ----------------------------------------------------------------
100  ! W3NMOD Subr. W3GDATMD Set number of model.
101  ! W3SETG Subr. Id. Point to selected model.
102  ! W3NDAT Subr. W3WDATMD Set number of model for wave data.
103  ! W3SETW Subr. Id. Point to selected model for wave data.
104  ! W2NAUX Subr. W3ADATMD Set number of model for aux data.
105  ! W3SETA Subr. Id. Point to selected model for aux data.
106  ! ITRACE Subr. W3SERVMD Subroutine tracing initialization.
107  ! STRACE Subr. Id. Subroutine tracing.
108  ! NEXTLN Subr. Id. Get next line from input file.
109  ! EXTCDE Subr. Id. Abort program as graceful as possible.
110  ! STME21 Subr. W3TIMEMD Convert time to string.
111  ! TICK21 Subr. Id. Advance time.
112  ! DSEC21 Func. Id. Difference between times.
113  ! W3IOGR Subr. W3IOGRMD Reading/writing model definition file.
114  ! W3IOGO Subr. W3IOGOMD Reading/writing raw gridded data file.
115  ! W3EXGO Subr. Internal Execute grid output.
116  ! W3TXTS Subr. Internal Text output at sea points only.
117  ! ----------------------------------------------------------------
118  !
119  ! 5. Called by :
120  !
121  ! None, stand-alone program.
122  !
123  ! 6. Error messages :
124  !
125  ! Checks on input, checks in W3IOxx.
126  !
127  ! 7. Remarks :
128  !
129  ! 8. Structure :
130  !
131  ! See source code.
132  !
133  ! 9. Switches :
134  !
135  ! !/S Enable subroutine tracing.
136  !
137  ! 10. Source code :
138  !
139  !/ ------------------------------------------------------------------- /
140  USE constants
141  !/
142  ! USE W3GDATMD, ONLY: W3NMOD, W3SETG
143  USE w3wdatmd, ONLY: w3ndat, w3setw
144  USE w3adatmd, ONLY: w3naux, w3seta
145  USE w3odatmd, ONLY: w3nout, w3seto
146  USE w3servmd, ONLY : itrace, nextln, extcde
147 #ifdef W3_S
148  USE w3servmd, ONLY : strace
149 #endif
150  USE w3timemd
151  USE w3iogrmd, ONLY: w3iogr
152  USE w3iogomd, ONLY: w3iogo, w3readflgrd
153  !/
154  USE w3gdatmd
155  USE w3wdatmd, ONLY: time, wlv, ice, iceh, icef, berg, ust, &
156  ustdir, rhoair
157  USE w3adatmd, ONLY: dw, ua, ud, as, cx, cy, hs, wlm, t0m1, thm, &
158  ths, fp0, thp0, dtdyn, fcut, &
159  aba, abd, uba, ubd, sxx, syy, sxy, usero, &
160  phs, ptp, plp, pdir, psi, pws, pwst, pnr, &
161  ptm1, pt1, pt2, pep, tauocx, tauocy, &
162  pthp0, pqp, psw, ppe, pgw, qp, qkk, &
163  skew, embia1, embia2, &
164  tauox, tauoy, tauwix,bhd, &
166  ussx, ussy, mssx, mssy, mscx, mscy, charn, &
169  cge, t01, hsig, stmaxe, stmaxd, hmaxe, &
170  hcmaxe, hmaxd, hcmaxd, mssd, mscd, wbt, &
172  USE w3odatmd, ONLY: ndso, ndse, ndst, nogrp, ngrpp, idout, &
174  !
175  IMPLICIT NONE
176  !/
177  !/ ------------------------------------------------------------------- /
178  !/ Local parameters
179  !/
180  INTEGER :: ndsi, ndsm, ndsog, ndsdat, ndsdt, &
181  ndstrc, ntrace, ierr, i, j, ifi, ifj,&
182  tout(2), tdum(2), iotest, nout, &
183  itype, ix1, ixn, ixs, iy1, iyn, iys, &
184  idla, idfm, iout, ipart
185 #ifdef W3_S
186  INTEGER, SAVE :: ient = 0
187 #endif
188  REAL :: dtreq, dtest
189  CHARACTER :: comstr*1, idtime*23, iddday*11, &
190  tabnme*9
191  LOGICAL :: flreq(nogrp,ngrpp), flog(nogrp), &
192  scale, vector, ltemp(ngrpp)
193  !/
194  !/ ------------------------------------------------------------------- /
195  !/
196  ! 1. IO set-up.
197  !
198  CALL w3nmod ( 1, 6, 6 )
199  CALL w3setg ( 1, 6, 6 )
200  CALL w3ndat ( 6, 6 )
201  CALL w3setw ( 1, 6, 6 )
202  CALL w3naux ( 6, 6 )
203  CALL w3seta ( 1, 6, 6 )
204  CALL w3nout ( 6, 6 )
205  CALL w3seto ( 1, 6, 6 )
206  !
207  ndsi = 10
208  ndsm = 20
209  ndsog = 20
210  ndsdat = 50
211  !
212  ndstrc = 6
213  ntrace = 10
214  CALL itrace ( ndstrc, ntrace )
215  !
216 #ifdef W3_S
217  CALL strace (ient, 'W3OUTF')
218 #endif
219  !
220  WRITE (ndso,900)
221  !
222  j = len_trim(fnmpre)
223  OPEN (ndsi,file=fnmpre(:j)//'ww3_outf.inp',status='OLD', &
224  err=800,iostat=ierr)
225  READ (ndsi,'(A)',END=801,ERR=802) comstr
226  IF (comstr.EQ.' ') comstr = '$'
227  WRITE (ndso,901) comstr
228  !
229  !--- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
230  ! 2. Read model definition file.
231  !
232  CALL w3iogr ( 'READ', ndsm )
233  WRITE (ndso,920) gname
234  !
235  !--- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
236  ! 3. Read general data and first fields from file
237  !
238  CALL w3iogo ( 'READ', ndsog, iotest )
239  !
240  WRITE (ndso,930)
241  DO ifi=1, nogrp
242  DO ifj=1, ngrpp
243  IF ( flogrd(ifi,ifj) ) WRITE (ndso,931) idout(ifi,ifj)
244  END DO
245  END DO
246  !
247  !--- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
248  ! 4. Read requests from input file.
249  ! Output times
250  !
251  CALL nextln ( comstr , ndsi , ndse )
252  READ (ndsi,*,END=801,ERR=802) TOUT, DTREQ, nout
253  dtreq = max( 0. , dtreq )
254  IF ( dtreq.EQ.0. ) nout = 1
255  nout = max( 1 , nout )
256  !
257  CALL stme21 ( tout , idtime )
258  WRITE (ndso,940) idtime
259  !
260  tdum = 0
261  CALL tick21 ( tdum , dtreq )
262  CALL stme21 ( tdum , idtime )
263  IF ( dtreq .GE. 86400. ) THEN
264  WRITE (iddday,'(I10,1X)') int(dtreq/86400.)
265  ELSE
266  iddday = ' '
267  END IF
268  idtime(1:11) = iddday
269  idtime(21:23) = ' '
270  WRITE (ndso,941) idtime, nout
271  !
272  ! ... Output fields
273  !
274  CALL w3readflgrd ( ndsi, ndso, 9, ndse, comstr, flog, &
275  flreq, 1, 1, ierr )
276  IF (ierr.NE.0) GOTO 800
277 
278  !
279  ! ... Output type
280  !
281  CALL nextln ( comstr , ndsi , ndse )
282  READ (ndsi,*,END=801,ERR=802) ITYPE, ipart
283  !Li IF ( ITYPE.LT.0 .OR. ITYPE.GT.3 ) THEN
284  IF ( itype.LT.0 .OR. itype.GT.4 ) THEN
285  !Li Type 4 for text output at sea points. JGLi12Dec2012
286  WRITE (ndse,1010) itype
287  CALL extcde ( 1 )
288  END IF
289  ipart = max( 0 , min( noswll , ipart ) )
290  !
291  ! ... ITYPE = 0
292  !
293  IF ( itype .EQ. 0 ) THEN
294  WRITE (ndso,942) itype, 'Checking contents of file'
295  DO
296  CALL stme21 ( time , idtime )
297  WRITE (ndso,943) idtime
298  CALL w3iogo ( 'READ', ndsog, iotest )
299  IF ( iotest .EQ. -1 ) THEN
300  WRITE (ndso,944)
301  GOTO 888
302  END IF
303  END DO
304  !
305  ! ... ITYPE = 1
306  !
307  ELSE IF (itype .EQ. 1) THEN
308  WRITE (ndso,942) itype, 'Print plots'
309  CALL nextln ( comstr , ndsi , ndse )
310  READ (ndsi,*,END=801,ERR=802) &
311  ix1, ixn, ixs, iy1, iyn, iys, scale, vector
312  ix1 = max( ix1 , 1 )
313  ixn = min( ixn , nx )
314  ixs = max( ixs , 1 )
315  iy1 = max( iy1 , 1 )
316  iyn = min( iyn , ny )
317  iys = max( iys , 1 )
318  WRITE (ndso,1940) ix1, ixn, ixs, iy1, iyn, iys
319  IF ( scale ) WRITE (ndso,1941)
320  !
321  ! ... ITYPE = 2
322  !
323  ELSE IF (itype .EQ. 2) THEN
324  WRITE (ndso,942) itype, 'Field statistics'
325  ndsdt = ndsdat - 1
326  CALL nextln ( comstr , ndsi , ndse )
327  READ (ndsi,*,END=801,ERR=802) IX1, IXN, IY1, iyn
328  ix1 = max( ix1 , 1 )
329  ixn = min( ixn , nx )
330  iy1 = max( iy1 , 1 )
331  iyn = min( iyn , ny )
332  WRITE (ndso,2940) ix1, ixn, iy1, iyn
333  !
334  ! ... ITYPE = 3
335  !
336  ELSE IF (itype .EQ. 3) THEN
337  WRITE (ndso,942) itype, 'Transfer files'
338  CALL nextln ( comstr , ndsi , ndse )
339  READ (ndsi,*,END=801,ERR=802) &
340  ix1, ixn, iy1, iyn, idla, idfm
341  ix1 = max( ix1 , 1 )
342  ixn = min( ixn , nx )
343  iy1 = max( iy1 , 1 )
344  iyn = min( iyn , ny )
345  IF (idla.LT.1 .OR. idla.GT.5) idla = 1
346  IF (idfm.LT.1 .OR. idfm.GT.3) idfm = 1
347  vector = .true.
348  WRITE (ndso,3940) ix1, ixn, iy1, iyn, idla, idfm
349  !
350  !Li Added sea-point output type 4. JGLi12Dec2012
351  !
352  ! ... ITYPE = 4
353  !
354  ELSE IF (itype .EQ. 4) THEN
355  WRITE (ndso,942) itype, 'Full sea-point output.'
356  CALL nextln ( comstr , ndsi , ndse )
357  READ (ndsi,*,END=801,ERR=802) &
358  ix1, ixn, iy1, iyn, idla, idfm
359  !Li
360  !
361  END IF
362  !
363  ! ... Output of output fields
364  !
365  IF ( itype.NE.2 ) THEN
366  WRITE (ndso,945)
367  ELSE
368  WRITE (ndso,2945)
369  END IF
370  !
371  DO ifi=1, nogrp
372  DO ifj=1, ngrpp
373  IF ( flreq(ifi,ifj) ) THEN
374  IF ( flogrd(ifi,ifj) ) THEN
375  IF ( itype.NE.2 ) THEN
376  WRITE (ndso,946) idout(ifi,ifj), ' '
377  ELSE
378  j = len_trim(fnmpre)
379  ndsdt = ndsdt + 1
380  WRITE (tabnme,'(A3,I2.2,A4)') 'tab', ndsdt, '.ww3'
381  WRITE (ndso,2946) tabnme, idout(ifi,ifj)
382  OPEN (ndsdt,file=fnmpre(:j)//tabnme)
383  WRITE (ndsdt,2947) idout(ifi,ifj)
384  END IF
385  ELSE
386  WRITE (ndso,946) idout(ifi,ifj), '*** NOT AVAILABLE ***'
387  flreq(ifi,ifj) = .false.
388  END IF
389  END IF
390  END DO
391  END DO
392  !
393  IF ( flog(4) ) THEN
394  IF ( ipart .EQ. 0 ) THEN
395  WRITE (ndso,948)
396  ELSE
397  WRITE (ndso,949) ipart
398  END IF
399  END IF
400  !
401  !--- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
402  ! 5. Time management.
403  !
404  iout = 0
405  IF (itype.EQ.3) WRITE (ndso,970)
406  !
407  DO
408  dtest = dsec21( time , tout )
409  IF ( dtest .GT. 0. ) THEN
410  CALL w3iogo ( 'READ', ndsog, iotest )
411  IF ( iotest .EQ. -1 ) THEN
412  WRITE (ndso,944)
413  GOTO 888
414  END IF
415  cycle
416  END IF
417  IF ( dtest .LT. 0. ) THEN
418  CALL tick21 ( tout , dtreq )
419  cycle
420  END IF
421  !
422  iout = iout + 1
423  CALL stme21 ( tout , idtime )
424  IF (itype.EQ.1) THEN
425  WRITE (ndso,950) idtime
426  ELSE IF (itype.EQ.3) THEN
427  WRITE (ndso,971) idtime
428  END IF
429  !
430  CALL w3exgo ( nx, ny, nsea )
431  !
432  CALL tick21 ( tout , dtreq )
433  IF ( iout .GE. nout ) EXIT
434  END DO
435  !
436  IF (itype.EQ.3) WRITE (ndso,972)
437  !
438  GOTO 888
439  !
440  ! Escape locations read errors :
441  !
442 800 CONTINUE
443  WRITE (ndse,1000) ierr
444  CALL extcde ( 10 )
445  !
446 801 CONTINUE
447  WRITE (ndse,1001)
448  CALL extcde ( 11 )
449  !
450 802 CONTINUE
451  WRITE (ndse,1002) ierr
452  CALL extcde ( 12 )
453  !
454 888 CONTINUE
455  WRITE (ndso,999)
456  !
457  ! Formats
458  !
459 900 FORMAT (/15x,' *** WAVEWATCH III Field output postp. *** '/ &
460  15x,'==============================================='/)
461 901 FORMAT ( ' Comment character is ''',a,''''/)
462  !
463 920 FORMAT ( ' Grid name : ',a/)
464  !
465 930 FORMAT ( ' Fields in file : '/ &
466  ' --------------------------')
467 931 FORMAT ( ' ',a)
468  !
469 940 FORMAT (/' Output time data : '/ &
470  ' --------------------------------------------------'/ &
471  ' First time : ',a)
472 941 FORMAT ( ' Interval : ',a/ &
473  ' Number of requests : ',i6)
474 942 FORMAT (/' Output type ',i2,' :'/ &
475  ' --------------------------------------------------'/ &
476  ' ',a/)
477 943 FORMAT ( ' Data for ',a)
478 944 FORMAT (/' End of file reached '/)
479  !
480 945 FORMAT (/' Requested output fields : '/ &
481  ' --------------------------------------------------')
482 2945 FORMAT (/' Output files and fields : '/ &
483  ' --------------------------------------------------')
484 946 FORMAT ( ' ',a,2x,a)
485 2946 FORMAT ( ' ',a,' : ',a)
486 2947 FORMAT ( ' Statitics of ',a/ &
487  ' (time, min, max, avg, std)'/)
488 948 FORMAT (/' Partitioned field data for wind seas')
489 949 FORMAT (/' Partitioned field data for swell field',i2)
490  !
491 1940 FORMAT ( ' X range and interval : ',3i5/ &
492  ' Y range and interval : ',3i5)
493 1941 FORMAT ( ' Data is normalized ')
494  !
495 2940 FORMAT ( ' X range : ',2i5/ &
496  ' Y range : ',2i5)
497  !
498 3940 FORMAT ( ' X range : ',2i5/ &
499  ' Y range : ',2i5/ &
500  ' Layout indicator : ',i5/ &
501  ' Format indicator : ',i5)
502  !
503 950 FORMAT (//' Output for ',a/ &
504  ' --------------------------------------------------')
505  !
506 970 FORMAT (//' Generating files '/ &
507  ' --------------------------------------------------')
508 971 FORMAT ( ' Files for ',a)
509 972 FORMAT ( ' ')
510  !
511 999 FORMAT (/' End of program '/ &
512  ' ========================================='/ &
513  ' WAVEWATCH III Field output '/)
514  !
515 1000 FORMAT (/' *** WAVEWATCH III ERROR IN W3OUTF : '/ &
516  ' ERROR IN OPENING INPUT FILE'/ &
517  ' IOSTAT =',i5/)
518  !
519 1001 FORMAT (/' *** WAVEWATCH III ERROR IN W3OUTF : '/ &
520  ' PREMATURE END OF INPUT FILE'/)
521  !
522 1002 FORMAT (/' *** WAVEWATCH III ERROR IN W3OUTF : '/ &
523  ' ERROR IN READING FROM INPUT FILE'/ &
524  ' IOSTAT =',i5/)
525  !
526 1010 FORMAT (/' *** WAVEWATCH III ERROR IN W3OUTF : '/ &
527  ' ILLEGAL TYPE, ITYPE =',i4/)
528  !/
529  !/ Internal subroutine W3EXGO ---------------------------------------- /
530  !/
531 CONTAINS
532  !/ ------------------------------------------------------------------- /
533 
547  SUBROUTINE w3exgo ( NX, NY, NSEA )
548  !/
549  !/ +-----------------------------------+
550  !/ | WAVEWATCH III NOAA/NCEP |
551  !/ | H. L. Tolman |
552  !/ | FORTRAN 90 |
553  !/ | Last update : 22-Mar-2021 |
554  !/ +-----------------------------------+
555  !/
556  !/ 26-Sep-1997 : Final FORTRAN 77 ( version 1.18 )
557  !/ 19-Jan-2000 : Upgrade to FORTRAN 90 ( version 2.00 )
558  !/ Massive changes to logistics
559  !/ 24-Jan-2001 : Flat grid version ( version 2.06 )
560  !/ 23-Apr-2002 : Clean-up ( version 2.19 )
561  !/ 29-Apr-2002 : Adding output fields 17-18. ( version 2.20 )
562  !/ 16-Oct-2002 : Fix bound. error for stress output. ( version 3.00 )
563  !/ 16-Oct-2002 : Fix statistical output for UNDEF. ( version 3.00 )
564  !/ 13-Nov-2002 : Add stress vector. ( version 3.00 )
565  !/ 24-Dec-2004 : Multiple grid version. ( version 3.06 )
566  !/ 21-Jul-2005 : Adding output fields 19-21. ( version 3.07 )
567  !/ 28-Jun-2006 : Adding file name preamble. ( version 3.09 )
568  !/ 05-Jul-2006 : Consolidate stress arrays. ( version 3.09 )
569  !/ 28-Mar-2007 : Adding partitioned output. ( version 3.11 )
570  !/ Adding user slots for outputs.
571  !/ 31-Jul-2007 : Fix file extension errors. ( version 3.12 )
572  !/ 25-Dec-2012 : New structure of output fields. ( version 4.11 )
573  !/ 25-Jun-2013 : Add type 4 sea point text output. ( version 4.11 )
574  !/ 26-Jan-2021 : Added TP field (derived from FP0) ( version 7.12 )
575  !/ 22-Mar-2021 : New coupling fields output ( version 7.13 )
576  !/
577  ! 1. Purpose :
578  !
579  ! Perform actual grid output.
580  !
581  ! 3. Parameters :
582  !
583  ! Parameter list
584  ! ----------------------------------------------------------------
585  ! NX/Y Int. I Grid dimensions.
586  ! NSEA Int. I Number of sea points.
587  ! ----------------------------------------------------------------
588  !
589  ! Internal parameters
590  ! ----------------------------------------------------------------
591  ! FLONE Log. Flags for one-dimensional field.
592  ! FLTWO Log. Flags for two-dimensional field X Y.
593  ! FLDIR Log. Flags for two-dimensional, directional field.
594  ! FLTRI Log. Flags for three dimensional field.
595  ! X1, X2, XX, XY
596  ! R.A. Output fields
597  ! ----------------------------------------------------------------
598  !
599  ! 4. Subroutines used :
600  !
601  ! Name Type Module Description
602  ! ----------------------------------------------------------------
603  ! STRACE Subr. W3SERVMD Subroutine tracing.
604  ! EXTCDE Subr. Id. Abort program as graceful as possible.
605  ! W3S2XY Subr. Id. Convert from storage to spatial grid.
606  ! PRTBLK Subr. W3ARRYMD Print plot of array.
607  ! OUTA2I Subr. Id. Print array of INTEGERS.
608  ! ----------------------------------------------------------------
609  !
610  ! 5. Called by :
611  !
612  ! Main program in which it is contained.
613  !
614  ! 6. Error messages :
615  !
616  ! None.
617  !
618  ! 7. Remarks :
619  !
620  ! - Note that arrays CX and CY of the main program now contain
621  ! the absolute current speed and direction respectively.
622  !
623  ! 8. Structure :
624  !
625  ! See source code.
626  !
627  ! 9. Switches :
628  !
629  ! !/S Enable subroutine tracing.
630  ! !/T Enable test output.
631  !
632  ! 10. Source code :
633  !
634  !/ ------------------------------------------------------------------- /
635  USE w3servmd, ONLY : w3s2xy
636 #ifdef W3_RTD
637  USE w3servmd, ONLY : w3thrtn, w3xyrtn
638 #endif
639  USE w3arrymd, ONLY : outa2i, prtblk
640  !/
641  !/ ------------------------------------------------------------------- /
642  !/ Parameter list
643  !/
644  INTEGER :: NX, NY, NSEA
645  !/
646  !/ ------------------------------------------------------------------- /
647  !/ Local parameters
648  !/
649  INTEGER :: NXMAX, NXTOT, NBLOK, IH, IM, IS, &
650  MFILL, J, ISEA, IX, IY, IXB, IB, &
651  IXA, NINGRD, JJ, IFI, IFJ
652  INTEGER :: MAP(NX+1,NY), MP2(NX+1,NY), &
653  MX1(NX,NY), MXX(NX,NY), MYY(NX,NY), &
654  MXY(NX,NY)
655  INTEGER, SAVE :: IPASS
656  ! INTEGER, SAVE :: NCOL = 80
657  INTEGER, SAVE :: NCOL = 132
658 #ifdef W3_S
659  INTEGER, SAVE :: IENT = 0
660 #endif
661  REAL :: FSC, CABS, UABS, FSCA, XMIN, XMAX, &
662  XAVG, XSTD, YGBX, XGBX, AABS
663  REAL :: X1(NX+1,NY), X2(NX+1,NY), &
664  XX(NX+1,NY), XY(NX+1,NY), DPTMAX(1)
665  !!Li Type 4 sea point only text output variables. JGLi25Jun2013
666  REAL, Dimension(NSEA) :: XS1, XS2, XS3, XS4, AUX
667  !!Li
668  DOUBLE PRECISION :: XDS, XDSQ
669  LOGICAL :: FLONE, FLTWO, FLDIR, FLTRI
670 #ifdef W3_T
671  LOGICAL :: LTEMP(NGRPP)
672 #endif
673  CHARACTER :: OLDTID*8, FNAME*32, ENAME*7, &
674  FORMG*12, FORMF*11, UNITS*10, FSCS*7
675  CHARACTER, SAVE :: TIMEID*8 = '00000000'
676  CHARACTER, SAVE :: FILEID*13 = 'WAVEWATCH III'
677 #ifdef W3_BT4
678  REAL, PARAMETER :: LOG2=log(2.)
679 #endif
680  !/
681  !/ ------------------------------------------------------------------- /
682  !/
683 #ifdef W3_S
684  CALL strace (ient, 'W3EXGO')
685 #endif
686  !
687 #ifdef W3_T
688  DO ifi=1, nogrp
689  ltemp = flreq(ifi,:)
690  WRITE (ndst,9000) ifi, ltemp
691  END DO
692  WRITE (ndst,9001) itype, ix1, ixn, ixs, iy1, iyn, iys, &
693  scale, vector, ndsdat
694 #endif
695  !
696  !--- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
697  ! 1. Preparations
698  !
699  x1 = undef
700  x2 = undef
701  xx = undef
702  xy = undef
703  !!Li Type 4 sea point only variables
704  xs1 = undef
705  xs2 = undef
706  xs3 = undef
707  xs4 = undef
708  !
709  ! Number of print-plots
710  !
711  IF ( itype .EQ. 1 ) THEN
712  IF ( scale ) THEN
713  nxmax = ( ncol - 10 ) / 2
714  ELSE
715  nxmax = ( ncol - 10 ) / 5
716  END IF
717  nxtot = 1 + (ixn-ix1)/ixs
718  nblok = 1 + (nxtot-1)/nxmax
719 #ifdef W3_T
720  WRITE (ndst,9012) nxmax, nxtot, nblok
721 #endif
722  END IF
723  !
724  ! Output file unit number
725  !
726  IF ( itype .EQ. 2 ) THEN
727  ndsdt = ndsdat - 1
728  ih = time(2) / 10000
729  im = mod( time(2) , 10000 ) / 100
730  is = mod( time(2) , 100 )
731  END IF
732  !
733  ! Set-up transfer files
734  !
735  !!Li Type 4 share filename with type 3 JGLi25Jun2013
736  !! IF ( ITYPE .EQ. 3 ) THEN
737  IF ( itype .EQ. 3 .OR. itype .EQ. 4 ) THEN
738  mfill = -999
739  oldtid = timeid
740  WRITE (timeid,'(I6.6,I2.2)') mod( time(1) , 1000000 ), &
741  time(2)/10000
742  fname(05:12) = timeid
743  fname(13:13) = '.'
744  IF ( timeid .NE. oldtid ) THEN
745  fname(1:4) = 'ww3.'
746  ipass = 1
747  ELSE
748  WRITE (ename,'(A1,I2.2,A1)') 'e', ipass, '.'
749  fname(1:4) = ename
750  ipass = ipass + 1
751  END IF
752 #ifdef W3_T
753  WRITE (ndst,9014) fname(1:13)
754 #endif
755  formg = '((10G12.2))'
756  END IF
757  !
758  !--- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
759  ! 2. Loop over output fields.
760  !
761  DO ifi=1, nogrp
762  DO ifj=1, ngrpp
763  IF ( flreq(ifi,ifj) ) THEN
764  !
765  formf = '(1X,32I4)'
766 #ifdef W3_T
767  WRITE (ndst,9020) idout(ifi,ifj)
768 #endif
769  !
770  ! 2.a Set output arrays and parameters
771 
772  flone = .false.
773  fltwo = .false.
774  fldir = .false.
775  fltri = .false.
776  !
777  IF ( ifi .EQ. 1 .AND. ifj .EQ. 1 ) THEN
778  flone = .true.
779  dptmax = maxval( dw(1:nsea) )
780  fsc = 1.
781  IF ( dptmax(1) .GT. 999. ) THEN
782  fsc = 0.1
783  ELSE IF ( dptmax(1) .GT. 99.9 ) THEN
784  fsc = 0.1
785  ELSE IF ( dptmax(1) .GT. 9.99 ) THEN
786  fsc = 0.01
787  END IF
788  IF ( itype .EQ. 3 ) fsc = 0.01
789  units = 'm'
790  ename = '.dpt'
791  formf = '(1X,17I7)'
792  IF ( itype .EQ. 4 ) THEN
793  xs1 = dw(1:nsea)
794  ELSE
795  CALL w3s2xy ( nsea, nsea, nx+1, ny, dw(1:nsea) &
796  , mapsf, x1 )
797  ENDIF
798  !
799  ELSE IF ( ifi .EQ. 1 .AND. ifj .EQ. 2 ) THEN
800  IF ( vector ) THEN
801  fltwo = .true.
802  ELSE
803  fldir = .true.
804  END IF
805  fsc = 0.01
806  ename = '.cur'
807  units = 'm s-1'
808  formf = '(1X,17I7)'
809 #ifdef W3_RTD
810  ! Rotate x,y vector back to standard pole
811  IF ( flagunr ) CALL w3xyrtn(nsea, cx, cy, angld)
812 #endif
813  IF ( itype .EQ. 4 ) THEN
814  xs1 = cx(1:nsea)
815  xs2 = cy(1:nsea)
816  ELSE
817  CALL w3s2xy ( nsea, nsea, nx+1, ny, cx(1:nsea) &
818  , mapsf, xx )
819  CALL w3s2xy ( nsea, nsea, nx+1, ny, cy(1:nsea) &
820  , mapsf, xy )
821  ENDIF
822  DO isea=1, nsea
823  cabs = sqrt(cx(isea)**2+cy(isea)**2)
824  IF ( cabs .GT. 0.05 ) THEN
825  cy(isea) = mod( 630. - &
826  rade*atan2(cy(isea),cx(isea)) , 360. )
827  ELSE
828  cy(isea) = undef
829  END IF
830  cx(isea) = cabs
831  END DO
832  IF ( itype .EQ. 4 ) THEN
833  xs3 = cx(1:nsea)
834  xs4 = cy(1:nsea)
835  ELSE
836  CALL w3s2xy ( nsea, nsea, nx+1, ny, cx(1:nsea) &
837  , mapsf, x1 )
838  CALL w3s2xy ( nsea, nsea, nx+1, ny, cy(1:nsea) &
839  , mapsf, x2 )
840  ENDIF
841  !
842  ELSE IF ( ifi .EQ. 1 .AND. ifj .EQ. 3 ) THEN
843  IF ( vector ) THEN
844  fltwo = .true.
845  ELSE
846  fldir = .true.
847  END IF
848  fsc = 0.1
849  ename = '.wnd'
850  units = 'm s-1'
851 #ifdef W3_RTD
852  ! Rotate x,y vector back to standard pole
853  IF ( flagunr ) CALL w3xyrtn(nsea, ua, ud, angld)
854 #endif
855  IF ( itype .EQ. 4 ) THEN
856  xs1 = ua(1:nsea)
857  xs2 = ud(1:nsea)
858  ELSE
859  CALL w3s2xy ( nsea, nsea, nx+1, ny, ua(1:nsea) &
860  , mapsf, xx )
861  CALL w3s2xy ( nsea, nsea, nx+1, ny, ud(1:nsea) &
862  , mapsf, xy )
863  ENDIF
864  DO isea=1, nsea
865  uabs = sqrt(ua(isea)**2+ud(isea)**2)
866  IF ( uabs .GT. 1.0 ) THEN
867  ud(isea) = mod( 630. - &
868  rade*atan2(ud(isea),ua(isea)) , 360. )
869  ELSE
870  ud(isea) = undef
871  END IF
872  ua(isea) = uabs
873  END DO
874  IF ( itype .EQ. 4 ) THEN
875  xs3 = ua(1:nsea)
876  xs4 = ud(1:nsea)
877  ELSE
878  CALL w3s2xy ( nsea, nsea, nx+1, ny, ua(1:nsea) &
879  , mapsf, x1 )
880  CALL w3s2xy ( nsea, nsea, nx+1, ny, ud(1:nsea) &
881  , mapsf, x2 )
882  ENDIF
883  !
884  ELSE IF ( ifi .EQ. 1 .AND. ifj .EQ. 4 ) THEN
885  flone = .true.
886  fsc = 0.1
887  ename = '.ast'
888  units = 'K'
889  IF ( itype .EQ. 4 ) THEN
890  xs1 = as(1:nsea)
891  ELSE
892  CALL w3s2xy ( nsea, nsea, nx+1, ny, as(1:nsea) &
893  , mapsf, x1 )
894  ENDIF
895  !
896  ELSE IF ( ifi .EQ. 1 .AND. ifj .EQ. 5 ) THEN
897  flone = .true.
898  fsc = 0.01
899  units = 'm'
900  ename = '.wlv'
901  IF ( itype .EQ. 4 ) THEN
902  xs1 = wlv
903  ELSE
904  CALL w3s2xy ( nsea, nsea, nx+1, ny, wlv , mapsf, x1 )
905  ENDIF
906  !
907  ELSE IF ( ifi .EQ. 1 .AND. ifj .EQ. 6 ) THEN
908  flone = .true.
909  fsc = 0.001
910  units = '1'
911  ename = '.ice'
912  IF ( itype .EQ. 4 ) THEN
913  xs1 = ice
914  ELSE
915  CALL w3s2xy ( nsea, nsea, nx+1, ny, ice , mapsf, x1 )
916  ENDIF
917  !
918  ELSE IF ( ifi .EQ. 1 .AND. ifj .EQ. 7 ) THEN
919  flone = .true.
920  fsc = 0.0002
921  units = 'km-1'
922  ename = '.ibg'
923  WHERE ( berg.NE.undef) berg = berg*0.1
924  IF ( itype .EQ. 4 ) THEN
925  xs1 = berg
926  ELSE
927  CALL w3s2xy ( nsea, nsea, nx+1, ny, berg , mapsf, x1 )
928  ENDIF
929  !
930  ELSE IF ( ifi .EQ. 1 .AND. ifj .EQ. 8 ) THEN
931  !! Note - TAUA and TAUADIR read in from .ww3 file are TAUAX,TAUAY
932  IF ( vector ) THEN
933  fltwo = .true.
934  ELSE
935  fldir = .true.
936  END IF
937  fsc = 0.01
938  units = 'Pa'
939  ename = '.taua'
940 #ifdef W3_RTD
941  ! Rotate x,y vector back to standard pole
942  IF ( flagunr ) THEN
943  CALL w3xyrtn(nsea, taua(1:nsea), tauadir(1:nsea), angld)
944  ENDIF
945 #endif
946  IF ( itype .EQ. 4 ) THEN
947  xs1 = taua(1:nsea)
948  xs2 = tauadir(1:nsea)
949  ELSE
950  CALL w3s2xy ( nsea, nsea, nx+1, ny, taua(1:nsea), mapsf, xx )
951  CALL w3s2xy ( nsea, nsea, nx+1, ny, tauadir(1:nsea), mapsf, xy )
952  ENDIF
953 
954  DO isea=1, nsea
955  uabs = sqrt(taua(isea)**2+tauadir(isea)**2)
956  IF ( uabs .GT. 0.01 ) THEN
957  tauadir(isea) = mod( 630. - &
958  rade*atan2(taua(isea),tauadir(isea)), 360.)
959  ELSE
960  tauadir(isea) = undef
961  END IF
962  ua(isea) = uabs
963  END DO
964  IF ( itype .EQ. 4 ) THEN
965  xs3 = taua(1:nsea)
966  xs4 = tauadir(1:nsea)
967  ELSE
968  CALL w3s2xy ( nsea, nsea, nx+1, ny, taua(1:nsea), mapsf, x1 )
969  CALL w3s2xy ( nsea, nsea, nx+1, ny, tauadir(1:nsea), mapsf, x2)
970  ENDIF
971  !
972  ELSE IF ( ifi .EQ. 1 .AND. ifj .EQ. 9 ) THEN
973  flone = .true.
974  fsc = 0.0001
975  units = 'kg m-3'
976  ename = '.rhoa'
977  IF ( itype .EQ. 4 ) THEN
978  xs1 = rhoair
979  ELSE
980  CALL w3s2xy ( nsea, nsea, nx+1, ny, rhoair, mapsf, x1 )
981  ENDIF
982  !
983 #ifdef W3_BT4
984  ELSE IF ( ifi .EQ. 1 .AND. ifj .EQ. 10 ) THEN
985  flone = .true.
986  fsc = 0.01
987  units = 'Krumbein phi scale'
988  ename = '.d50'
989  WHERE ( sed_d50.NE.undef) sed_d50 = -log(sed_d50/0.001)/log2
990  IF ( itype .EQ. 4 ) THEN
991  xs1 = sed_d50
992  ELSE
993  CALL w3s2xy ( nsea, nsea, nx+1, ny, sed_d50 , mapsf, x1 )
994  ENDIF
995 #endif
996  !
997 #ifdef W3_IS2
998  ELSE IF ( ifi .EQ. 1 .AND. ifj .EQ. 11 ) THEN
999  flone = .true.
1000  fsc = 0.001
1001  units = 'm'
1002  ename = '.ic1'
1003  IF ( itype .EQ. 4) THEN
1004  xs1 = iceh
1005  ELSE
1006  CALL w3s2xy (nsea, nsea, nx+1, ny, iceh, mapsf, x1 )
1007  ENDIF
1008 #endif
1009  !
1010 #ifdef W3_IS2
1011  ELSE IF ( ifi .EQ. 1 .AND. ifj .EQ. 12) THEN
1012  flone = .true.
1013  fsc = 0.001
1014  units = 'm'
1015  ename = '.ic5'
1016  IF ( itype .EQ. 4) THEN
1017  xs1 = icef
1018  ELSE
1019  CALL w3s2xy (nsea, nsea, nx+1, ny, icef, mapsf, x1 )
1020  ENDIF
1021 #endif
1022 
1023  ELSE IF ( ifi .EQ. 2 .AND. ifj .EQ. 1 ) THEN
1024  flone = .true.
1025  fsc = 0.01
1026  units = 'm'
1027  ename = '.hs'
1028  IF ( itype .EQ. 4 ) THEN
1029  xs1 = hs
1030  ELSE
1031  CALL w3s2xy ( nsea, nsea, nx+1, ny, hs , mapsf, x1 )
1032  ENDIF
1033  !
1034  ELSE IF ( ifi .EQ. 2 .AND. ifj .EQ. 2 ) THEN
1035  flone = .true.
1036  fsc = 1.
1037  units = 'm'
1038  ename = '.lm'
1039  IF ( itype .EQ. 4 ) THEN
1040  xs1 = wlm
1041  ELSE
1042  CALL w3s2xy ( nsea, nsea, nx+1, ny, wlm , mapsf, x1 )
1043  ENDIF
1044  !
1045  ELSE IF ( ifi .EQ. 2 .AND. ifj .EQ. 3 ) THEN
1046  flone = .true.
1047  fsc = 0.01
1048  units = 's'
1049  ename = '.t02'
1050  IF ( itype .EQ. 4 ) THEN
1051  xs1 = t02
1052  ELSE
1053  CALL w3s2xy ( nsea, nsea, nx+1, ny, t02 , mapsf, x1 )
1054  ENDIF
1055  !
1056  ELSE IF ( ifi .EQ. 2 .AND. ifj .EQ. 4 ) THEN
1057  flone = .true.
1058  fsc = 0.01
1059  units = 's'
1060  ename = '.t0m1'
1061  IF ( itype .EQ. 4 ) THEN
1062  xs1 = t0m1
1063  ELSE
1064  CALL w3s2xy ( nsea, nsea, nx+1, ny, t0m1 , mapsf, x1 )
1065  ENDIF
1066  !
1067  ELSE IF ( ifi .EQ. 2 .AND. ifj .EQ. 5 ) THEN
1068  flone = .true.
1069  fsc = 0.01
1070  units = 's'
1071  ename = '.t01'
1072  IF ( itype .EQ. 4 ) THEN
1073  xs1 = t01
1074  ELSE
1075  CALL w3s2xy ( nsea, nsea, nx+1, ny, t01 , mapsf, x1 )
1076  ENDIF
1077  !
1078  ELSE IF ( ifi .EQ. 2 .AND. ifj .EQ. 6 ) THEN
1079  flone = .true.
1080  fsc = 0.001
1081  units = 's-1'
1082  ename = '.fp'
1083  IF ( itype .EQ. 4 ) THEN
1084  xs1 = fp0
1085  ELSE
1086  CALL w3s2xy ( nsea, nsea, nx+1, ny, fp0 , mapsf, x1 )
1087  ENDIF
1088  !
1089  ELSE IF ( ifi .EQ. 2 .AND. ifj .EQ. 7 ) THEN
1090  flone = .true.
1091  fsc = 1.
1092  units = 'degree'
1093  ename = '.dir'
1094 #ifdef W3_RTD
1095  ! Rotate direction back to standard pole
1096  IF ( flagunr ) CALL w3thrtn(nsea, thm, angld, .false.)
1097 #endif
1098  DO isea=1, nsea
1099  IF ( thm(isea) .NE. undef ) &
1100  thm(isea) = mod( 630. - rade*thm(isea) , 360. )
1101  END DO
1102  IF ( itype .EQ. 4 ) THEN
1103  xs1 = thm
1104  ELSE
1105  CALL w3s2xy ( nsea, nsea, nx+1, ny, thm , mapsf, x1 )
1106  ENDIF
1107  !
1108  ELSE IF ( ifi .EQ. 2 .AND. ifj .EQ. 8 ) THEN
1109  flone = .true.
1110  fsc = 0.1
1111  units = 'degree'
1112  ename = '.spr'
1113  IF ( itype .EQ. 4 ) THEN
1114  xs1 = ths
1115  ELSE
1116  CALL w3s2xy ( nsea, nsea, nx+1, ny, ths , mapsf, x1 )
1117  ENDIF
1118  !
1119  ELSE IF ( ifi .EQ. 2 .AND. ifj .EQ. 9 ) THEN
1120  flone = .true.
1121  fsc = 1.
1122  units = 'degree'
1123  ename = '.dp'
1124 #ifdef W3_RTD
1125  ! Rotate direction back to standard pole
1126  IF ( flagunr ) CALL w3thrtn(nsea, thp0, angld, .false.)
1127 #endif
1128  DO isea=1, nsea
1129  IF ( thp0(isea) .NE. undef ) THEN
1130  thp0(isea) = mod( 630-rade*thp0(isea) , 360. )
1131  END IF
1132  END DO
1133  IF ( itype .EQ. 4 ) THEN
1134  xs1 = thp0
1135  ELSE
1136  CALL w3s2xy ( nsea, nsea, nx+1, ny, thp0 , mapsf, x1 )
1137  ENDIF
1138  !
1139  ELSE IF ( ifi .EQ. 2 .AND. ifj .EQ. 10 ) THEN
1140  flone = .true.
1141  fsc = 0.001
1142  units = 'm'
1143  ename = '.hig'
1144  IF ( itype .EQ. 4 ) THEN
1145  xs1 = hsig
1146  ELSE
1147  CALL w3s2xy ( nsea, nsea, nx+1, ny, hsig , mapsf, x1 )
1148  END IF
1149  !
1150  ELSE IF ( ifi .EQ. 2 .AND. ifj .EQ. 11 ) THEN
1151  flone = .true.
1152  fsc = 0.002
1153  units = 'm'
1154  ename = '.emc'
1155  IF ( itype .EQ. 4 ) THEN
1156  xs1 = stmaxe
1157  ELSE
1158  CALL w3s2xy ( nsea, nsea, nx+1, ny, stmaxe, mapsf, x1 )
1159  END IF
1160  !
1161  ELSE IF ( ifi .EQ. 2 .AND. ifj .EQ. 12 ) THEN
1162  flone = .true.
1163  fsc = 0.002
1164  units = 'm'
1165  ename = '.smc'
1166  IF ( itype .EQ. 4 ) THEN
1167  xs1 = stmaxd
1168  ELSE
1169  CALL w3s2xy ( nsea, nsea, nx+1, ny, stmaxd, mapsf, x1 )
1170  ENDIF
1171  !
1172  ELSE IF ( ifi .EQ. 2 .AND. ifj .EQ. 13 ) THEN
1173  flone = .true.
1174  fsc = 0.002
1175  units = 'm'
1176  ename = '.emh'
1177  IF ( itype .EQ. 4 ) THEN
1178  xs1 = hmaxe
1179  ELSE
1180  CALL w3s2xy ( nsea, nsea, nx+1, ny, hmaxe, mapsf, x1 )
1181  END IF
1182  !
1183  ELSE IF ( ifi .EQ. 2 .AND. ifj .EQ. 14 ) THEN
1184  flone = .true.
1185  fsc = 0.002
1186  units = 'm'
1187  ename = '.eml'
1188  IF ( itype .EQ. 4 ) THEN
1189  xs1 = hcmaxe
1190  ELSE
1191  CALL w3s2xy ( nsea, nsea, nx+1, ny, hcmaxe, mapsf, x1 )
1192  ENDIF
1193  !
1194  ELSE IF ( ifi .EQ. 2 .AND. ifj .EQ. 15 ) THEN
1195  flone = .true.
1196  fsc = 0.002
1197  units = 'm'
1198  ename = '.smh'
1199  IF ( itype .EQ. 4 ) THEN
1200  xs1 = hmaxd
1201  ELSE
1202  CALL w3s2xy ( nsea, nsea, nx+1, ny, hmaxd, mapsf, x1 )
1203  END IF
1204  !
1205  ELSE IF ( ifi .EQ. 2 .AND. ifj .EQ. 16 ) THEN
1206  flone = .true.
1207  fsc = 0.002
1208  units = 'm'
1209  ename = '.sml'
1210  IF ( itype .EQ. 4 ) THEN
1211  xs1 = hcmaxd
1212  ELSE
1213  CALL w3s2xy ( nsea, nsea, nx+1, ny, hcmaxd, mapsf, x1)
1214  ENDIF
1215  !
1216  ELSE IF ( ifi .EQ. 2 .AND. ifj .EQ. 17 ) THEN
1217  flone = .true.
1218  fsc = 0.001
1219  units = '1'
1220  ename = '.wbt'
1221  IF ( itype .EQ. 4 ) THEN
1222  xs1 = wbt
1223  ELSE
1224  CALL w3s2xy ( nsea, nsea, nx+1, ny, wbt, mapsf, x1)
1225  ENDIF
1226  !
1227  ELSE IF ( ifi .EQ. 2 .AND. ifj .EQ. 18 ) THEN
1228  flone = .true.
1229  fsc = 0.01
1230  units = 's'
1231  ename = '.tp'
1232  DO i=1,nsea
1233  IF(fp0(i) .NE. undef) THEN
1234  aux(i) = 1.0 / fp0(i)
1235  ELSE
1236  aux(i) = undef
1237  ENDIF
1238  ENDDO
1239 
1240  IF ( itype .EQ. 4 ) THEN
1241  xs1 = aux
1242  ELSE
1243  CALL w3s2xy ( nsea, nsea, nx+1, ny, aux, mapsf, x1)
1244  ENDIF
1245  !
1246  ELSE IF ( ifi .EQ. 2 .AND. ifj .EQ. 19 ) THEN
1247  flone = .true.
1248  fsc = 0.001
1249  units = 'm-1'
1250  ename = '.wnm'
1251  IF ( itype .EQ. 4 ) THEN
1252  xs1 = wnmean
1253  ELSE
1254  CALL w3s2xy ( nsea, nsea, nx+1, ny, wnmean, mapsf, x1)
1255  ENDIF
1256  !
1257  ELSE IF ( ifi .EQ. 4 .AND. ifj .EQ. 1 ) THEN
1258  flone = .true.
1259  fsc = 0.01
1260  units = 'm'
1261  ename = '.phs'
1262  IF ( itype .EQ. 4 ) THEN
1263  xs1 = phs(:,ipart)
1264  ELSE
1265  CALL w3s2xy ( nsea, nsea, nx+1, ny, phs(:,ipart) &
1266  , mapsf, x1 )
1267  ENDIF
1268  !
1269  ELSE IF ( ifi .EQ. 4 .AND. ifj .EQ. 2 ) THEN
1270  flone = .true.
1271  fsc = 0.01
1272  units = 's'
1273  ename = '.ptp'
1274  IF ( itype .EQ. 4 ) THEN
1275  xs1 = ptp(:,ipart)
1276  ELSE
1277  CALL w3s2xy ( nsea, nsea, nx+1, ny, ptp(:,ipart) &
1278  , mapsf, x1 )
1279  ENDIF
1280  !
1281  ELSE IF ( ifi .EQ. 4 .AND. ifj .EQ. 3 ) THEN
1282  flone = .true.
1283  fsc = 1.
1284  units = 'm'
1285  ename = '.plp'
1286  IF ( itype .EQ. 4 ) THEN
1287  xs1 = plp(:,ipart)
1288  ELSE
1289  CALL w3s2xy ( nsea, nsea, nx+1, ny, plp(:,ipart) &
1290  , mapsf, x1 )
1291  ENDIF
1292  !
1293  ELSE IF ( ifi .EQ. 4 .AND. ifj .EQ. 4 ) THEN
1294  flone = .true.
1295  fsc = 1.
1296  units = 'degree'
1297  ename = '.pdir'
1298 #ifdef W3_RTD
1299  ! Rotate direction back to standard pole
1300  IF ( flagunr ) CALL w3thrtn(nsea, pdir(:,ipart), angld, .false.)
1301 #endif
1302  DO isea=1, nsea
1303  IF ( pdir(isea,ipart) .NE. undef ) THEN
1304  pdir(isea,ipart) = &
1305  mod( 630-rade*pdir(isea,ipart) , 360. )
1306  END IF
1307  END DO
1308  IF ( itype .EQ. 4 ) THEN
1309  xs1 = pdir(:,ipart)
1310  ELSE
1311  CALL w3s2xy ( nsea, nsea, nx+1, ny, pdir(:,ipart) &
1312  , mapsf, x1 )
1313  ENDIF
1314  !
1315  ELSE IF ( ifi .EQ. 4 .AND. ifj .EQ. 5 ) THEN
1316  flone = .true.
1317  fsc = 0.1
1318  units = 'degree'
1319  ename = '.pspr'
1320  IF ( itype .EQ. 4 ) THEN
1321  xs1 = psi(:,ipart)
1322  ELSE
1323  CALL w3s2xy ( nsea, nsea, nx+1, ny, psi(:,ipart) &
1324  , mapsf, x1 )
1325  ENDIF
1326  !
1327  ELSE IF ( ifi .EQ. 4 .AND. ifj .EQ. 6 ) THEN
1328  flone = .true.
1329  fsc = 0.001
1330  units = '1'
1331  ename = '.pws'
1332  IF ( itype .EQ. 4 ) THEN
1333  xs1 = pws(:,ipart)
1334  ELSE
1335  CALL w3s2xy ( nsea, nsea, nx+1, ny, pws(:,ipart) &
1336  , mapsf, x1 )
1337  ENDIF
1338  !
1339  ELSE IF ( ifi .EQ. 4 .AND. ifj .EQ. 7 ) THEN
1340  flone = .true.
1341  fsc = 1.0
1342  units = 'degree'
1343  ename = '.pdp'
1344 #ifdef W3_RTD
1345  ! Rotate direction back to standard pole
1346  IF ( flagunr ) CALL w3thrtn(nsea, pthp0(:,ipart), angld, .false.)
1347 #endif
1348  DO isea=1, nsea
1349  IF ( pthp0(isea,ipart) .NE. undef ) THEN
1350  pthp0(isea,ipart) = &
1351  mod( 630-rade*pthp0(isea,ipart) , 360. )
1352  END IF
1353  END DO
1354  IF ( itype .EQ. 4 ) THEN
1355  xs1 = pthp0(:,ipart)
1356  ELSE
1357  CALL w3s2xy ( nsea, nsea, nx+1, ny, pthp0(:,ipart), &
1358  mapsf, x1 )
1359  ENDIF
1360  !
1361  ELSE IF ( ifi .EQ. 4 .AND. ifj .EQ. 8 ) THEN
1362  flone = .true.
1363  fsc = 0.01
1364  units = '1'
1365  ename = '.pqp'
1366  IF ( itype .EQ. 4 ) THEN
1367  xs1 = pqp(:,ipart)
1368  ELSE
1369  CALL w3s2xy ( nsea, nsea, nx+1, ny, pqp(:,ipart), mapsf, x1 )
1370  ENDIF
1371  !
1372  ELSE IF ( ifi .EQ. 4 .AND. ifj .EQ. 9 ) THEN
1373  flone = .true.
1374  fsc = 0.01
1375  units = '1'
1376  ename = '.ppe'
1377  IF ( itype .EQ. 4 ) THEN
1378  xs1 = ppe(:,ipart)
1379  ELSE
1380  CALL w3s2xy ( nsea, nsea, nx+1, ny, ppe(:,ipart), mapsf, x1 )
1381  ENDIF
1382  !
1383  ELSE IF ( ifi .EQ. 4 .AND. ifj .EQ. 10 ) THEN
1384  flone = .true.
1385  fsc = 0.0001
1386  units = 's-1'
1387  ename = '.pgw'
1388  IF ( itype .EQ. 4 ) THEN
1389  xs1 = pgw(:,ipart)
1390  ELSE
1391  CALL w3s2xy ( nsea, nsea, nx+1, ny, pgw(:,ipart), mapsf, x1 )
1392  ENDIF
1393  !
1394  ELSE IF ( ifi .EQ. 4 .AND. ifj .EQ. 11 ) THEN
1395  flone = .true.
1396  fsc = 0.0001
1397  units = '1'
1398  ename = '.psw'
1399  IF ( itype .EQ. 4 ) THEN
1400  xs1 = psw(:,ipart)
1401  ELSE
1402  CALL w3s2xy ( nsea, nsea, nx+1, ny, psw(:,ipart), mapsf, x1 )
1403  ENDIF
1404  !
1405  ELSE IF ( ifi .EQ. 4 .AND. ifj .EQ. 12 ) THEN
1406  flone = .true.
1407  fsc = 0.01
1408  units = 's'
1409  ename = '.ptm10'
1410  IF ( itype .EQ. 4 ) THEN
1411  xs1 = ptm1(:,ipart)
1412  ELSE
1413  CALL w3s2xy ( nsea, nsea, nx+1, ny, ptm1(:,ipart), mapsf, x1 )
1414  ENDIF
1415  !
1416  ELSE IF ( ifi .EQ. 4 .AND. ifj .EQ. 13 ) THEN
1417  flone = .true.
1418  fsc = 0.01
1419  units = 's'
1420  ename = '.pt01'
1421  IF ( itype .EQ. 4 ) THEN
1422  xs1 = pt1(:,ipart)
1423  ELSE
1424  CALL w3s2xy ( nsea, nsea, nx+1, ny, pt1(:,ipart), mapsf, x1 )
1425  ENDIF
1426  !
1427  ELSE IF ( ifi .EQ. 4 .AND. ifj .EQ. 14 ) THEN
1428  flone = .true.
1429  fsc = 0.01
1430  units = 's'
1431  ename = '.pt02'
1432  IF ( itype .EQ. 4 ) THEN
1433  xs1 = pt2(:,ipart)
1434  ELSE
1435  CALL w3s2xy ( nsea, nsea, nx+1, ny, pt2(:,ipart), mapsf, x1 )
1436  ENDIF
1437  !
1438  ELSE IF ( ifi .EQ. 4 .AND. ifj .EQ. 15 ) THEN
1439  flone = .true.
1440  fsc = 0.02
1441  units = 'm2 s rad-1'
1442  ename = '.pep'
1443  IF ( itype .EQ. 4 ) THEN
1444  xs1 = pep(:,ipart)
1445  ELSE
1446  CALL w3s2xy ( nsea, nsea, nx+1, ny, pep(:,ipart), mapsf, x1 )
1447  ENDIF
1448  !
1449  ELSE IF ( ifi .EQ. 4 .AND. ifj .EQ. 16 ) THEN
1450  flone = .true.
1451  fsc = 0.001
1452  units = '1'
1453  ename = '.tws'
1454  IF ( itype .EQ. 4 ) THEN
1455  xs1 = pwst(:)
1456  ELSE
1457  CALL w3s2xy ( nsea, nsea, nx+1, ny, pwst(:), mapsf, x1 )
1458  ENDIF
1459  !
1460  ELSE IF ( ifi .EQ. 4 .AND. ifj .EQ. 17 ) THEN
1461  flone = .true.
1462  fsc = 1.
1463  units = '1'
1464  ename = '.pnr'
1465  IF ( itype .EQ. 4 ) THEN
1466  xs1 = pnr(:)
1467  ELSE
1468  CALL w3s2xy ( nsea, nsea, nx+1, ny, pnr(:), mapsf, x1 )
1469  ENDIF
1470  !
1471  ELSE IF ( ifi .EQ. 5 .AND. ifj .EQ. 1 ) THEN
1472  IF ( vector ) THEN
1473  fltwo = .true.
1474  ELSE
1475  fldir = .true.
1476  END IF
1477  fsc = 0.001
1478  ename = '.ust'
1479  formf = '(1X,20I6)'
1480  units = 'm s-1'
1481 #ifdef W3_RTD
1482  ! Rotate x,y vector back to standard pole
1483  IF ( flagunr ) CALL w3xyrtn(nsea, ust, ustdir, angld)
1484 #endif
1485  IF ( itype .EQ. 4 ) THEN
1486  xs1 = ust(1:nsea)
1487  xs2 = ustdir(1:nsea)
1488  ELSE
1489  CALL w3s2xy (nsea,nsea,nx+1,ny, ust(1:nsea) &
1490  , mapsf, xx )
1491  CALL w3s2xy (nsea,nsea,nx+1,ny, ustdir(1:nsea) &
1492  , mapsf, xy )
1493  ENDIF
1494  DO isea=1, nsea
1495  uabs = sqrt(ust(isea)**2+ustdir(isea)**2)
1496  IF ( ust(isea) .EQ. undef ) THEN
1497  ustdir(isea) = undef
1498  uabs = undef
1499  ELSE IF ( uabs .GT. 0.05 ) THEN
1500  ustdir(isea) = mod( 630. - &
1501  rade*atan2(ustdir(isea),ust(isea)) , 360. )
1502  ELSE
1503  ustdir(isea) = undef
1504  END IF
1505  ust(isea) = uabs
1506  END DO
1507  IF ( itype .EQ. 4 ) THEN
1508  xs3 = ust(1:nsea)
1509  xs4 = ustdir(1:nsea)
1510  ELSE
1511  CALL w3s2xy (nsea,nsea,nx+1,ny, ust(1:nsea) &
1512  , mapsf, x1 )
1513  CALL w3s2xy (nsea,nsea,nx+1,ny, ustdir(1:nsea) &
1514  , mapsf, x2 )
1515  ENDIF
1516  !
1517  ELSE IF ( ifi .EQ. 5 .AND. ifj .EQ. 2 ) THEN
1518  flone = .true.
1519  fsc = 1.e-6
1520  units = '1'
1521  ename = '.cha'
1522  IF ( itype .EQ. 4 ) THEN
1523  xs1 = charn(1:nsea)
1524  ELSE
1525  CALL w3s2xy ( nsea, nsea, nx+1, ny, charn(1:nsea) &
1526  , mapsf, x1 )
1527  ENDIF
1528  !
1529  ELSE IF ( ifi .EQ. 5 .AND. ifj .EQ. 3 ) THEN
1530  flone = .true.
1531  fsc = 0.1 !0.01
1532  units = 'kW m-1'
1533  ename = '.cge'
1534  DO isea=1, nsea
1535  IF ( cge(isea) .NE. undef ) &
1536  cge(isea) = 0.001 * cge(isea)
1537  END DO
1538  IF ( itype .EQ. 4 ) THEN
1539  xs1 = cge(1:nsea)
1540  ELSE
1541  CALL w3s2xy ( nsea, nsea, nx+1, ny, cge(1:nsea) &
1542  , mapsf, x1 )
1543  ENDIF
1544  !
1545  ELSE IF ( ifi .EQ. 5 .AND. ifj .EQ. 4 ) THEN
1546  flone = .true.
1547  fsc = 0.01
1548  units = 'W m-2'
1549  ename = '.faw'
1550  DO isea=1, nsea
1551  phiaw(isea)=min(99.98,phiaw(isea))
1552  END DO
1553  IF ( itype .EQ. 4 ) THEN
1554  xs1 = phiaw(1:nsea)
1555  ELSE
1556  CALL w3s2xy ( nsea, nsea, nx+1, ny, phiaw(1:nsea) &
1557  , mapsf, x1 )
1558  ENDIF
1559  !
1560  ELSE IF ( ifi .EQ. 5 .AND. ifj .EQ. 5 ) THEN
1561  IF ( vector ) THEN
1562  fltwo = .true.
1563  ELSE
1564  fldir = .true.
1565  END IF
1566  fsc = 1.e-6
1567  units = 'm2 s-2'
1568  ename = '.taw'
1569 #ifdef W3_RTD
1570  ! Rotate x,y vector back to standard pole
1571  IF ( flagunr ) CALL w3xyrtn(nsea, tauwix, tauwiy, angld)
1572 #endif
1573  IF ( itype .EQ. 4 ) THEN
1574  xs1 = tauwix(1:nsea)
1575  xs2 = tauwiy(1:nsea)
1576  ELSE
1577  CALL w3s2xy ( nsea, nsea, nx+1, ny, tauwix(1:nsea) &
1578  , mapsf, xx )
1579  CALL w3s2xy ( nsea, nsea, nx+1, ny, tauwiy(1:nsea) &
1580  , mapsf, xy )
1581  ENDIF
1582  DO isea=1, nsea
1583  cabs = sqrt(tauwix(isea)**2+tauwiy(isea)**2)
1584  IF ( tauwix(isea) .EQ. undef ) THEN
1585  tauwiy(isea) = undef
1586  cabs = undef
1587  ELSE IF ( tauwix(isea) .EQ. 0. .AND. &
1588  tauwiy(isea) .EQ. 0. ) THEN
1589  tauwiy(isea) = undef
1590  ELSE
1591  tauwiy(isea) = mod( 630. - &
1592  rade*atan2(tauwiy(isea),tauwix(isea)) , 360. )
1593  END IF
1594  tauwix(isea) = cabs
1595  END DO
1596  IF ( itype .EQ. 4 ) THEN
1597  xs3 = tauwix(1:nsea)
1598  xs4 = tauwiy(1:nsea)
1599  ELSE
1600  CALL w3s2xy ( nsea, nsea, nx+1, ny, tauwix(1:nsea) &
1601  , mapsf, x1 )
1602  CALL w3s2xy ( nsea, nsea, nx+1, ny, tauwiy(1:nsea) &
1603  , mapsf, x2 )
1604  ENDIF
1605  !
1606  ELSE IF ( ifi .EQ. 5 .AND. ifj .EQ. 6 ) THEN
1607  IF ( vector ) THEN
1608  fltwo = .true.
1609  ELSE
1610  fldir = .true.
1611  END IF
1612  fsc = 1.e-6
1613  units = 'm2 s-2'
1614  ename = '.twa'
1615 #ifdef W3_RTD
1616  ! Rotate x,y vector back to standard pole
1617  IF ( flagunr ) CALL w3xyrtn(nsea, tauwnx, tauwny, angld)
1618 #endif
1619  IF ( itype .EQ. 4 ) THEN
1620  xs1 = tauwnx(1:nsea)
1621  xs2 = tauwny(1:nsea)
1622  ELSE
1623  CALL w3s2xy ( nsea, nsea, nx+1, ny, tauwnx(1:nsea) &
1624  , mapsf, xx )
1625  CALL w3s2xy ( nsea, nsea, nx+1, ny, tauwny(1:nsea) &
1626  , mapsf, xy )
1627  ENDIF
1628  DO isea=1, nsea
1629  cabs = sqrt(tauwnx(isea)**2+tauwny(isea)**2)
1630  IF ( tauwnx(isea) .EQ. undef ) THEN
1631  tauwny(isea) = undef
1632  cabs = undef
1633  ELSE IF ( tauwnx(isea) .EQ. 0. .AND. &
1634  tauwny(isea) .EQ. 0. ) THEN
1635  tauwny(isea) = undef
1636  ELSE
1637  tauwny(isea) = mod( 630. - &
1638  rade*atan2(tauwny(isea),tauwnx(isea)) , 360. )
1639  END IF
1640  tauwnx(isea) = cabs
1641  END DO
1642  IF ( itype .EQ. 4 ) THEN
1643  xs3 = tauwnx(1:nsea)
1644  xs4 = tauwny(1:nsea)
1645  ELSE
1646  CALL w3s2xy ( nsea, nsea, nx+1, ny, tauwnx(1:nsea) &
1647  , mapsf, x1 )
1648  CALL w3s2xy ( nsea, nsea, nx+1, ny, tauwny(1:nsea) &
1649  , mapsf, x2 )
1650  ENDIF
1651  !
1652  ELSE IF ( ifi .EQ. 5 .AND. ifj .EQ. 7 ) THEN
1653  flone = .true.
1654  fsc = 0.001
1655  units = '1'
1656  ename = '.wcc'
1657  IF ( itype .EQ. 4 ) THEN
1658  xs1 = whitecap(1:nsea,1)
1659  ELSE
1660  CALL w3s2xy ( nsea, nsea, nx+1, ny, whitecap(1:nsea,1) &
1661  , mapsf, x1 )
1662  ENDIF
1663  !
1664  ELSE IF ( ifi .EQ. 5 .AND. ifj .EQ. 8 ) THEN
1665  flone = .true.
1666  fsc = 0.1
1667  units = 'm'
1668  ename = '.wcf'
1669  IF ( itype .EQ. 4 ) THEN
1670  xs1 = whitecap(1:nsea,2)
1671  ELSE
1672  CALL w3s2xy ( nsea, nsea, nx+1, ny, whitecap(1:nsea,2) &
1673  , mapsf, x1 )
1674  ENDIF
1675  !
1676  ELSE IF ( ifi .EQ. 5 .AND. ifj .EQ. 9 ) THEN
1677  flone = .true.
1678  fsc = 0.1
1679  units = 'm'
1680  ename = '.wch'
1681  IF ( itype .EQ. 4 ) THEN
1682  xs1 = whitecap(1:nsea,3)
1683  ELSE
1684  CALL w3s2xy ( nsea, nsea, nx+1, ny, whitecap(1:nsea,3) &
1685  , mapsf, x1 )
1686  ENDIF
1687  !
1688  ELSE IF ( ifi .EQ. 5 .AND. ifj .EQ. 10 ) THEN
1689  flone = .true.
1690  fsc = 0.1
1691  units = '1'
1692  ename = '.wcm'
1693  IF ( itype .EQ. 4 ) THEN
1694  xs1 = whitecap(1:nsea,4)
1695  ELSE
1696  CALL w3s2xy ( nsea, nsea, nx+1, ny, whitecap(1:nsea,4) &
1697  , mapsf, x1 )
1698  ENDIF
1699  !
1700  ELSE IF ( ifi .EQ. 6 .AND. ifj .EQ. 1 ) THEN
1701  fltri = .true.
1702  fsc = 10.
1703  units = 'N m-1'
1704  ename = '.sxy'
1705 #ifdef W3_RTD
1706  ! Radition stress components are always left on rotated pole
1707  ! at present - need to confirm how to de-rotate (A. Saulter)
1708 #endif
1709  IF ( itype .EQ. 4 ) THEN
1710  xs1 = sxx(1:nsea)
1711  xs2 = syy(1:nsea)
1712  xs3 = sxy(1:nsea)
1713  ELSE
1714  CALL w3s2xy ( nsea, nsea, nx+1, ny, sxx(1:nsea) &
1715  , mapsf, x1 )
1716  CALL w3s2xy ( nsea, nsea, nx+1, ny, syy(1:nsea) &
1717  , mapsf, x2 )
1718  CALL w3s2xy ( nsea, nsea, nx+1, ny, sxy(1:nsea) &
1719  , mapsf, xy )
1720  ENDIF
1721  !
1722  ELSE IF ( ifi .EQ. 6 .AND. ifj .EQ. 2 ) THEN
1723  IF ( vector ) THEN
1724  fltwo = .true.
1725  ELSE
1726  fldir = .true.
1727  END IF
1728  fsc = 1.e-6
1729  units = 'm2 s-2'
1730  ename = '.two'
1731 #ifdef W3_RTD
1732  ! Rotate x,y vector back to standard pole
1733  IF ( flagunr ) CALL w3xyrtn(nsea, tauox, tauoy, angld)
1734 #endif
1735  IF ( itype .EQ. 4 ) THEN
1736  xs1 = tauox(1:nsea)
1737  xs2 = tauoy(1:nsea)
1738  ELSE
1739  CALL w3s2xy ( nsea, nsea, nx+1, ny, tauox(1:nsea) &
1740  , mapsf, xx )
1741  CALL w3s2xy ( nsea, nsea, nx+1, ny, tauoy(1:nsea) &
1742  , mapsf, xy )
1743  ENDIF
1744  DO isea=1, nsea
1745  uabs = sqrt(tauox(isea)**2+tauoy(isea)**2)
1746  IF ( tauox(isea) .EQ. undef ) THEN
1747  tauoy(isea) = undef
1748  uabs = undef
1749  ELSE IF ( uabs .GT. 1.e-8 ) THEN
1750  tauoy(isea) = mod( 630. - &
1751  rade*atan2(tauoy(isea),tauox(isea)) , 360. )
1752  ELSE
1753  tauoy(isea) = undef
1754  END IF
1755  tauox(isea) = uabs
1756  END DO
1757  IF ( itype .EQ. 4 ) THEN
1758  xs3 = tauox(1:nsea)
1759  xs4 = tauoy(1:nsea)
1760  ELSE
1761  CALL w3s2xy (nsea,nsea,nx+1,ny, tauox(1:nsea) &
1762  , mapsf, x1 )
1763  CALL w3s2xy (nsea,nsea,nx+1,ny, tauoy(1:nsea) &
1764  , mapsf, x2 )
1765  ENDIF
1766  !
1767  ELSE IF ( ifi .EQ. 6 .AND. ifj .EQ.3 ) THEN
1768  flone = .true.
1769  fsc = 0.001
1770  units = 'N m-1'
1771  ename = '.bhd'
1772  IF ( itype .EQ. 4 ) THEN
1773  xs1 = bhd(1:nsea)
1774  ELSE
1775  CALL w3s2xy ( nsea, nsea, nx+1, ny, bhd(1:nsea) &
1776  , mapsf, x1 )
1777  ENDIF
1778  !
1779  ELSE IF ( ifi .EQ. 6 .AND. ifj .EQ. 4 ) THEN
1780  flone = .true.
1781  fsc = 0.1
1782  units = 'W m-2'
1783  ename = '.foc'
1784  DO isea=1, nsea
1785  phioc(isea)=min(99.98,phioc(isea))
1786  END DO
1787  IF ( itype .EQ. 4 ) THEN
1788  xs1 = phioc(1:nsea)
1789  ELSE
1790  CALL w3s2xy ( nsea, nsea, nx+1, ny, phioc(1:nsea) &
1791  , mapsf, x1 )
1792  ENDIF
1793  !
1794  ELSE IF ( ifi .EQ. 6 .AND. ifj .EQ. 5 ) THEN
1795  IF ( vector ) THEN
1796  fltwo = .true.
1797  ELSE
1798  fldir = .true.
1799  END IF
1800  fsc = 0.001
1801  units = 'm2 s-1'
1802  ename = '.tus'
1803 #ifdef W3_RTD
1804  ! Rotate x,y vector back to standard pole
1805  IF ( flagunr ) CALL w3xyrtn(nsea, tusx, tusy, angld)
1806 #endif
1807  IF ( itype .EQ. 4 ) THEN
1808  xs1 = tusx(1:nsea)
1809  xs2 = tusy(1:nsea)
1810  ELSE
1811  CALL w3s2xy ( nsea, nsea, nx+1, ny, tusx(1:nsea) &
1812  , mapsf, xx )
1813  CALL w3s2xy ( nsea, nsea, nx+1, ny, tusy(1:nsea) &
1814  , mapsf, xy )
1815  ENDIF
1816  DO isea=1, nsea
1817  cabs = sqrt(tusx(isea)**2+tusy(isea)**2)
1818  IF ( tusx(isea) .NE. undef ) THEN
1819  tusy(isea) = mod( 630. - &
1820  rade*atan2(tusy(isea),tusx(isea)) , 360. )
1821  ELSE
1822  tusy(isea) = undef
1823  cabs = undef
1824  END IF
1825  tusx(isea) = cabs
1826  END DO
1827  IF ( itype .EQ. 4 ) THEN
1828  xs3 = tusx(1:nsea)
1829  xs4 = tusy(1:nsea)
1830  ELSE
1831  CALL w3s2xy ( nsea, nsea, nx+1, ny,tusx,mapsf, x1 )
1832  CALL w3s2xy ( nsea, nsea, nx+1, ny,tusy,mapsf, x2 )
1833  ENDIF
1834  !
1835  ELSE IF ( ifi .EQ. 6 .AND. ifj .EQ. 6 ) THEN
1836  IF ( vector ) THEN
1837  fltwo = .true.
1838  ELSE
1839  fldir = .true.
1840  END IF
1841  fsc = 0.001
1842  units = 'm s-1'
1843  ename = '.uss'
1844  DO isea=1, nsea
1845  IF (ussx(isea) .NE. undef ) THEN
1846  ussx(isea)=max(-0.9998,min(0.9998,ussx(isea)))
1847  ussy(isea)=max(-0.9998,min(0.9998,ussy(isea)))
1848  END IF
1849  END DO
1850 #ifdef W3_RTD
1851  ! Rotate x,y vector back to standard pole
1852  IF ( flagunr ) CALL w3xyrtn(nsea, ussx, ussy, angld)
1853 #endif
1854  IF ( itype .EQ. 4 ) THEN
1855  xs1 = ussx(1:nsea)
1856  xs2 = ussy(1:nsea)
1857  ELSE
1858  CALL w3s2xy ( nsea, nsea, nx+1, ny, ussx(1:nsea) &
1859  , mapsf, xx )
1860  CALL w3s2xy ( nsea, nsea, nx+1, ny, ussy(1:nsea) &
1861  , mapsf, xy )
1862  ENDIF
1863  DO isea=1, nsea
1864  cabs = sqrt(ussx(isea)**2+ussy(isea)**2)
1865  IF ( ussx(isea) .NE. undef ) THEN
1866  ussy(isea) = mod( 630. - &
1867  rade*atan2(ussy(isea),ussx(isea)) , 360. )
1868  ELSE
1869  ussy(isea) = undef
1870  cabs = undef
1871  END IF
1872  ussx(isea) = cabs
1873  END DO
1874  IF ( itype .EQ. 4 ) THEN
1875  xs3 = ussx(1:nsea)
1876  xs4 = ussy(1:nsea)
1877  ELSE
1878  CALL w3s2xy ( nsea, nsea, nx+1, ny, ussx(1:nsea), &
1879  mapsf, x1 )
1880  CALL w3s2xy ( nsea, nsea, nx+1, ny, ussy(1:nsea), &
1881  mapsf, x2 )
1882  ENDIF
1883  !
1884  ELSE IF ( ifi .EQ. 6 .AND. ifj .EQ. 7 ) THEN
1885  fltwo = .true.
1886  fsc = 0.01
1887  ename = '.p2s'
1888  units = 'm4'
1889  DO isea=1, nsea
1890  prms(isea)=prms(isea)
1891  END DO
1892  IF ( itype .EQ. 4 ) THEN
1893  xs1 = prms(1:nsea)
1894  xs2 = tpms(1:nsea)
1895  ELSE
1896  CALL w3s2xy ( nsea, nsea, nx+1,ny,prms,mapsf, x1 )
1897  CALL w3s2xy ( nsea, nsea, nx+1,ny,tpms,mapsf, x2 )
1898  ENDIF
1899  !
1900  ELSE IF ( ifi .EQ. 6 .AND. ifj .EQ. 13 ) THEN
1901  IF ( vector ) THEN
1902  fltwo = .true.
1903  ELSE
1904  fldir = .true.
1905  END IF
1906  fsc = 1.e-6
1907  units = 'm2 s-2'
1908  ename = '.toc'
1909 #ifdef W3_RTD
1910  ! Rotate x,y vector back to standard pole
1911  IF ( flagunr ) CALL w3xyrtn(nsea, tauocx, tauocy, angld)
1912 #endif
1913  IF ( itype .EQ. 4 ) THEN
1914  xs1 = tauocx(1:nsea)
1915  xs2 = tauocy(1:nsea)
1916  ELSE
1917  CALL w3s2xy ( nsea, nsea, nx+1, ny, tauocx(1:nsea) &
1918  , mapsf, xx )
1919  CALL w3s2xy ( nsea, nsea, nx+1, ny, tauocy(1:nsea) &
1920  , mapsf, xy )
1921  ENDIF
1922  DO isea=1, nsea
1923  uabs = sqrt(tauocx(isea)**2+tauocy(isea)**2)
1924  IF ( tauocx(isea) .EQ. undef ) THEN
1925  tauocy(isea) = undef
1926  uabs = undef
1927  ELSE IF ( uabs .GT. 1.e-8 ) THEN
1928  tauocy(isea) = mod( 630. - &
1929  rade*atan2(tauocy(isea),tauocx(isea)) , 360. )
1930  ELSE
1931  tauocy(isea) = undef
1932  END IF
1933  tauocx(isea) = uabs
1934  END DO
1935  IF ( itype .EQ. 4 ) THEN
1936  xs3 = tauocx(1:nsea)
1937  xs4 = tauocy(1:nsea)
1938  ELSE
1939  CALL w3s2xy (nsea,nsea,nx+1,ny, tauocx(1:nsea) &
1940  , mapsf, x1 )
1941  CALL w3s2xy (nsea,nsea,nx+1,ny, tauocy(1:nsea) &
1942  , mapsf, x2 )
1943  ENDIF
1944  !
1945  ELSE IF ( ifi .EQ. 7 .AND. ifj .EQ. 1 ) THEN
1946  IF ( vector ) THEN
1947  fltwo = .true.
1948  ELSE
1949  fldir = .true.
1950  END IF
1951  fsc = 0.01
1952  ename = '.abr'
1953  units = 'm'
1954 #ifdef W3_RTD
1955  ! Rotate x,y vector back to standard pole
1956  IF ( flagunr ) CALL w3xyrtn(nsea, aba, abd, angld)
1957 #endif
1958  IF ( itype .EQ. 4 ) THEN
1959  xs1 = aba(1:nsea)
1960  xs2 = abd(1:nsea)
1961  ELSE
1962  CALL w3s2xy ( nsea, nsea, nx+1, ny, aba(1:nsea) &
1963  , mapsf, xx )
1964  CALL w3s2xy ( nsea, nsea, nx+1, ny, abd(1:nsea) &
1965  , mapsf, xy )
1966  ENDIF
1967  DO isea=1, nsea
1968  IF ( aba(isea) .NE. undef ) THEN
1969  aabs = sqrt(aba(isea)**2+abd(isea)**2)
1970  IF ( aabs .GT. 0.005 ) THEN
1971  abd(isea) = mod( 630. - &
1972  rade*atan2(abd(isea),aba(isea)) , 360. )
1973  ELSE
1974  abd(isea) = undef
1975  END IF
1976  aba(isea) = aabs
1977  END IF
1978  END DO
1979  IF ( itype .EQ. 4 ) THEN
1980  xs3 = aba(1:nsea)
1981  xs4 = abd(1:nsea)
1982  ELSE
1983  CALL w3s2xy ( nsea, nsea, nx+1, ny, aba(1:nsea) &
1984  , mapsf, x1 )
1985  CALL w3s2xy ( nsea, nsea, nx+1, ny, abd(1:nsea) &
1986  , mapsf, x2 )
1987  ENDIF
1988  !
1989  ELSE IF ( ifi .EQ. 7 .AND. ifj .EQ. 2 ) THEN
1990  IF ( vector ) THEN
1991  fltwo = .true.
1992  ELSE
1993  fldir = .true.
1994  END IF
1995  fsc = 0.01
1996  ename = '.ubr'
1997  units = 'm s-1'
1998 #ifdef W3_RTD
1999  ! Rotate x,y vector back to standard pole
2000  IF ( flagunr ) CALL w3xyrtn(nsea, uba, ubd, angld)
2001 #endif
2002  IF ( itype .EQ. 4 ) THEN
2003  xs1 = uba(1:nsea)
2004  xs2 = ubd(1:nsea)
2005  ELSE
2006  CALL w3s2xy ( nsea, nsea, nx+1, ny, uba(1:nsea) &
2007  , mapsf, xx )
2008  CALL w3s2xy ( nsea, nsea, nx+1, ny, ubd(1:nsea) &
2009  , mapsf, xy )
2010  ENDIF
2011  DO isea=1, nsea
2012  IF ( uba(isea) .NE. undef ) THEN
2013  uabs = sqrt(uba(isea)**2+ubd(isea)**2)
2014  IF ( uabs .GT. 0.005 ) THEN
2015  ubd(isea) = mod( 630. - &
2016  rade*atan2(ubd(isea),uba(isea)) , 360. )
2017  ELSE
2018  ubd(isea) = undef
2019  END IF
2020  uba(isea) = uabs
2021  END IF
2022  END DO
2023  IF ( itype .EQ. 4 ) THEN
2024  xs3 = uba(1:nsea)
2025  xs4 = ubd(1:nsea)
2026  ELSE
2027  CALL w3s2xy ( nsea, nsea, nx+1, ny, uba(1:nsea) &
2028  , mapsf, x1 )
2029  CALL w3s2xy ( nsea, nsea, nx+1, ny, ubd(1:nsea) &
2030  , mapsf, x2 )
2031  ENDIF
2032  !
2033  ELSE IF ( ifi .EQ. 7 .AND. ifj .EQ. 3 ) THEN
2034  fltri = .true.
2035  fsc = 1.e-2
2036  units = 'm'
2037  ename = '.bed'
2038 #ifdef W3_RTD
2039  ! Rotate x,y vector back to standard pole
2040  IF ( flagunr ) CALL w3xyrtn(nsea, bedforms(1:nsea,2), &
2041  bedforms(1:nsea,3), angld)
2042 #endif
2043  IF ( itype .EQ. 4 ) THEN
2044  xs1 = bedforms(1:nsea,1)
2045  xs2 = bedforms(1:nsea,2)
2046  xs3 = bedforms(1:nsea,3)
2047  ELSE
2048  CALL w3s2xy ( nsea, nsea, nx+1, ny, bedforms(1:nsea,1) &
2049  , mapsf, x1 )
2050  CALL w3s2xy ( nsea, nsea, nx+1, ny, bedforms(1:nsea,2) &
2051  , mapsf, x2 )
2052  CALL w3s2xy ( nsea, nsea, nx+1, ny, bedforms(1:nsea,3) &
2053  , mapsf, xy )
2054  ENDIF
2055  !
2056  ELSE IF ( ifi .EQ. 7 .AND. ifj .EQ. 4 ) THEN
2057  flone = .true.
2058  fsc = 0.1
2059  units = 'W m-2'
2060  ename = '.fbb'
2061  IF ( itype .EQ. 4 ) THEN
2062  xs1 = phibbl(1:nsea)
2063  ELSE
2064  CALL w3s2xy ( nsea, nsea, nx+1, ny, phibbl(1:nsea) &
2065  , mapsf, x1 )
2066  ENDIF
2067  !
2068  ELSE IF ( ifi .EQ. 7 .AND. ifj .EQ. 5 ) THEN
2069  fltwo = .true.
2070  fsc = 1.e-6
2071  units = 'm2 s-2'
2072  ename = '.tbb'
2073 #ifdef W3_RTD
2074  ! Rotate x,y vector back to standard pole
2075  IF ( flagunr ) CALL w3xyrtn(nsea, taubbl(1:nsea,1), &
2076  taubbl(1:nsea,2), angld)
2077 #endif
2078  IF ( itype .EQ. 4 ) THEN
2079  xs1 = taubbl(1:nsea,1)
2080  xs2 = taubbl(1:nsea,2)
2081  ELSE
2082  CALL w3s2xy ( nsea, nsea, nx+1, ny, taubbl(1:nsea,1) &
2083  , mapsf, xx )
2084  CALL w3s2xy ( nsea, nsea, nx+1, ny, taubbl(1:nsea,2) &
2085  , mapsf, xy )
2086  ENDIF
2087  !
2088  ELSE IF ( ifi .EQ. 8 .AND. ifj .EQ. 1 ) THEN
2089  IF ( vector ) THEN
2090  fltwo = .true.
2091  ELSE
2092  fldir = .true.
2093  END IF
2094  fsc = 1.e-6
2095  ename = '.mss'
2096  formf = '(1X,20I6)'
2097  units = '1'
2098 #ifdef W3_RTD
2099  ! Rotate x,y vector back to standard pole
2100  IF ( flagunr ) CALL w3xyrtn(nsea, mssx, mssy, angld)
2101 #endif
2102  IF ( itype .EQ. 4 ) THEN
2103  xs1 = mssx(1:nsea)
2104  xs2 = mssy(1:nsea)
2105  ELSE
2106  CALL w3s2xy ( nsea, nsea, nx+1, ny, mssx(1:nsea), &
2107  mapsf, xx )
2108  CALL w3s2xy ( nsea, nsea, nx+1, ny ,mssy(1:nsea), &
2109  mapsf, xy )
2110  ENDIF
2111  !
2112  ELSE IF ( ifi .EQ. 8 .AND. ifj .EQ. 2 ) THEN
2113  IF ( vector ) THEN
2114  fltwo = .true.
2115  ELSE
2116  fldir = .true.
2117  END IF
2118  fsc = 0.00001
2119  ename = '.msc'
2120  units = '1'
2121 #ifdef W3_RTD
2122  ! Rotate x,y vector back to standard pole
2123  IF ( flagunr ) CALL w3xyrtn(nsea, mscx, mscy, angld)
2124 #endif
2125  IF ( itype .EQ. 4 ) THEN
2126  xs1 = mscx(1:nsea)
2127  xs2 = mscy(1:nsea)
2128  ELSE
2129  CALL w3s2xy ( nsea, nsea, nx+1, ny, mscx(1:nsea), &
2130  mapsf, xx )
2131  CALL w3s2xy ( nsea, nsea, nx+1, ny, mscy(1:nsea), &
2132  mapsf, xy )
2133  ENDIF
2134  DO isea=1, nsea
2135  cabs = sqrt(mscx(isea)**2+mscy(isea)**2)
2136  IF ( mscx(isea) .EQ. undef ) THEN
2137  mscy(isea) = undef
2138  cabs = undef
2139  ELSE IF ( mscx(isea) .EQ. 0. .AND. &
2140  mscy(isea) .EQ. 0. ) THEN
2141  mscy(isea) = undef
2142  ELSE
2143  mscy(isea) = mod( 630. - &
2144  rade*atan2(mscy(isea),mscx(isea)) , 360. )
2145  END IF
2146  mscx(isea) = cabs
2147  END DO
2148  IF ( itype .EQ. 4 ) THEN
2149  xs3 = mscx(1:nsea)
2150  xs4 = mscy(1:nsea)
2151  ELSE
2152  CALL w3s2xy ( nsea, nsea, nx+1, ny, mscx(1:nsea), &
2153  mapsf, x1 )
2154  CALL w3s2xy ( nsea, nsea, nx+1, ny, mscy(1:nsea), &
2155  mapsf, x2 )
2156  ENDIF
2157  !
2158  ELSE IF ( ifi .EQ. 8 .AND. ifj .EQ. 3 ) THEN
2159  flone = .true.
2160  fsc = 0.1
2161  units = 'degree'
2162  ename = '.msd'
2163 #ifdef W3_RTD
2164  ! Rotate direction back to standard pole
2165  IF ( flagunr ) CALL w3thrtn(nsea, mssd, angld, .false.)
2166 #endif
2167  DO isea=1, nsea
2168  IF ( mssd(isea) .NE. undef ) THEN
2169  mssd(isea) = mod( 630. - rade*mssd(isea) , 180. )
2170  END IF
2171  END DO
2172  IF ( itype .EQ. 4 ) THEN
2173  xs1 = mssd(1:nsea)
2174  ELSE
2175  CALL w3s2xy ( nsea, nsea, nx+1, ny, mssd(1:nsea), mapsf, x1 )
2176  ENDIF
2177  !
2178  ELSE IF ( ifi .EQ. 8 .AND. ifj .EQ. 4 ) THEN
2179  flone = .true.
2180  fsc = 0.1
2181  units = 'degree'
2182  ename = '.mcd'
2183 #ifdef W3_RTD
2184  ! Rotate direction back to standard pole
2185  IF ( flagunr ) CALL w3thrtn(nsea, mscd, angld, .false.)
2186 #endif
2187  DO isea=1, nsea
2188  IF ( mscd(isea) .NE. undef ) THEN
2189  mscd(isea) = mod( 630. - rade*mscd(isea) , 180. )
2190  END IF
2191  END DO
2192  IF ( itype .EQ. 4 ) THEN
2193  xs1 = mscd(1:nsea)
2194  ELSE
2195  CALL w3s2xy ( nsea, nsea, nx+1, ny, mscd(1:nsea), mapsf, x1 )
2196  ENDIF
2197  !
2198  ELSE IF ( ifi .EQ. 8 .AND. ifj .EQ. 5 ) THEN
2199  flone = .true.
2200  fsc = 0.001
2201  units = '1'
2202  ename = '.qp'
2203  IF ( itype .EQ. 4 ) THEN
2204  xs1 = qp
2205  ELSE
2206  CALL w3s2xy ( nsea, nsea, nx+1, ny, qp, mapsf, x1 )
2207  ENDIF
2208  !
2209  ELSE IF ( ifi .EQ. 8 .AND. ifj .EQ. 6 ) THEN
2210  flone = .true.
2211  fsc = 0.05
2212  units = '1'
2213  ename = '.qkk'
2214  IF ( itype .EQ. 4 ) THEN
2215  xs1 = qkk
2216  ELSE
2217  CALL w3s2xy ( nsea, nsea, nx+1, ny, qkk, mapsf, x1 )
2218  ENDIF
2219  !
2220  ELSE IF ( ifi .EQ. 8 .AND. ifj .EQ. 7 ) THEN
2221  flone = .true.
2222  fsc = 0.01
2223  units = '1'
2224  ename = '.skw'
2225  IF ( itype .EQ. 4 ) THEN
2226  xs1 = skew
2227  ELSE
2228  CALL w3s2xy ( nsea, nsea, nx+1, ny, skew, mapsf, x1 )
2229  ENDIF
2230  !
2231  ELSE IF ( ifi .EQ. 8 .AND. ifj .EQ. 8 ) THEN
2232  flone = .true.
2233  fsc = 0.0001
2234  units = '1'
2235  ename = '.emb'
2236  IF ( itype .EQ. 4 ) THEN
2237  xs1 = embia1
2238  ELSE
2239  CALL w3s2xy ( nsea, nsea, nx+1, ny, embia1, mapsf, x1 )
2240  ENDIF
2241  !
2242  ELSE IF ( ifi .EQ. 8 .AND. ifj .EQ. 9 ) THEN
2243  flone = .true.
2244  fsc = 0.0001
2245  units = '1'
2246  ename = '.emc'
2247  IF ( itype .EQ. 4 ) THEN
2248  xs1 = embia2
2249  ELSE
2250  CALL w3s2xy ( nsea, nsea, nx+1, ny, embia2, mapsf, x1 )
2251  ENDIF
2252  !
2253  ELSE IF ( ifi .EQ. 9 .AND. ifj .EQ. 1 ) THEN
2254  flone = .true.
2255  fsc = 0.1
2256  units = 'min.'
2257  ename = '.dtd'
2258  DO isea=1, nsea
2259  IF ( dtdyn(isea) .NE. undef ) &
2260  dtdyn(isea) = dtdyn(isea) / 60.
2261  END DO
2262  IF ( itype .EQ. 4 ) THEN
2263  xs1 = dtdyn(1:nsea)
2264  ELSE
2265  CALL w3s2xy ( nsea, nsea, nx+1, ny, dtdyn , mapsf, x1 )
2266  ENDIF
2267  !
2268  ELSE IF ( ifi .EQ. 9 .AND. ifj .EQ. 2 ) THEN
2269  flone = .true.
2270  fsc = 0.001
2271  units = 's-1'
2272  ename = '.fc'
2273  IF ( itype .EQ. 4 ) THEN
2274  xs1 = fcut
2275  ELSE
2276  CALL w3s2xy ( nsea, nsea, nx+1, ny, fcut , mapsf, x1 )
2277  ENDIF
2278  !
2279  ELSE IF ( ifi .EQ. 9 .AND. ifj .EQ. 3 ) THEN
2280  flone = .true.
2281  fsc = 0.001
2282  fsc = 1.
2283  units = '1'
2284  ename = '.cfx'
2285  IF ( itype .EQ. 4 ) THEN
2286  xs1 = cflxymax
2287  ELSE
2288  CALL w3s2xy ( nsea, nsea, nx+1, ny, cflxymax, mapsf, x1 )
2289  ENDIF
2290  !
2291  ELSE IF ( ifi .EQ. 9 .AND. ifj .EQ. 4 ) THEN
2292  flone = .true.
2293  fsc = 0.001
2294  units = '1'
2295  ename = '.cfd'
2296  IF ( itype .EQ. 4 ) THEN
2297  xs1 = cflthmax
2298  ELSE
2299  CALL w3s2xy ( nsea, nsea, nx+1, ny, cflthmax, mapsf, x1 )
2300  ENDIF
2301  !
2302  ELSE IF ( ifi .EQ. 9 .AND. ifj .EQ. 5 ) THEN
2303  flone = .true.
2304  fsc = 0.001
2305  units = '1'
2306  ename = '.cfk'
2307  IF ( itype .EQ. 4 ) THEN
2308  xs1 = cflkmax
2309  ELSE
2310  CALL w3s2xy ( nsea, nsea, nx+1, ny, cflkmax, mapsf, x1 )
2311  ENDIF
2312  !
2313  ELSE IF ( ifi .EQ. 10 ) THEN
2314  flone = .true.
2315  fsc = 1.
2316  units = 'TBD'
2317  WRITE (ename,'(A2,I2.2)') '.u', ifj
2318  IF ( itype .EQ. 4 ) THEN
2319  xs1 = usero(:,ifj)
2320  ELSE
2321  CALL w3s2xy ( nsea, nsea, nx+1, ny, usero(:,ifj) &
2322  , mapsf, x1 )
2323  ENDIF
2324  !
2325  ELSE
2326  WRITE (ndse,990) ifi,ifj
2327  WRITE (ndse,999)
2328  CALL extcde ( 1 )
2329  !
2330  END IF
2331  !
2332  ! 2.b Make map
2333  !
2334  DO ix=1, nx
2335  DO iy=1, ny
2336  IF ( mapsta(iy,ix) .EQ. 0 ) THEN
2337  x1(ix,iy) = undef
2338  x2(ix,iy) = undef
2339  xx(ix,iy) = undef
2340  xy(ix,iy) = undef
2341  END IF
2342  IF ( x1(ix,iy) .EQ. undef ) THEN
2343  map(ix,iy) = 0
2344  ELSE
2345  map(ix,iy) = 1
2346  END IF
2347  IF ( x2(ix,iy) .EQ. undef ) THEN
2348  mp2(ix,iy) = 0
2349  ELSE
2350  mp2(ix,iy) = 1
2351  END IF
2352  END DO
2353  END DO
2354  !
2355  ! 2.c Perform output type 1 ( print plots )
2356  !
2357  IF ( itype .EQ. 1 ) THEN
2358  !
2359  IF ( scale ) THEN
2360  fsc = 0.
2361  fsca = 0.
2362  ELSE
2363  fsca = 1.
2364  END IF
2365  ixb = ix1 - ixs
2366  !
2367  DO ib=1, nblok
2368  ixa = ixb + ixs
2369  ixb = ixa + (nxmax-1)*ixs
2370  ixb = min( ixb , ixn )
2371  IF ( fltri ) THEN
2372  CALL prtblk (ndso, nx, ny, nx+1, x1, map, 0, &
2373  fsc, ixa, ixb, ixs, iy1, iyn, iys, &
2374  idout(ifi,ifj), units)
2375  CALL prtblk (ndso, nx, ny, nx+1, x2, map, 0, &
2376  fsc, ixa, ixb, ixs, iy1, iyn, iys, &
2377  idout(ifi,ifj), units)
2378  CALL prtblk (ndso, nx, ny, nx+1, xy, map, 0, &
2379  fsc, ixa, ixb, ixs, iy1, iyn, iys, &
2380  idout(ifi,ifj), units)
2381  ELSE IF ( flone ) THEN
2382  CALL prtblk (ndso, nx, ny, nx+1, x1, map, 0, &
2383  fsc, ixa, ixb, ixs, iy1, iyn, iys, &
2384  idout(ifi,ifj), units)
2385  ELSE IF ( fltwo ) THEN
2386  CALL prtblk (ndso, nx, ny, nx+1, xx, map, 0, &
2387  fsc, ixa, ixb, ixs, iy1, iyn, iys, &
2388  idout(ifi,ifj), units)
2389  CALL prtblk (ndso, nx, ny, nx+1, xy, map, 0, &
2390  fsc, ixa, ixb, ixs, iy1, iyn, iys, &
2391  idout(ifi,ifj), units)
2392  ELSE IF ( fldir ) THEN
2393  CALL prtblk (ndso, nx, ny, nx+1, x1, map, 0, &
2394  fsc, ixa, ixb, ixs, iy1, iyn, iys, &
2395  idout(ifi,ifj), units)
2396  CALL prtblk (ndso, nx, ny, nx+1, x2, mp2, 0, &
2397  fsca, ixa, ixb, ixs, iy1, iyn, iys, &
2398  idout(ifi,ifj), 'Deg.')
2399  END IF
2400  END DO
2401  !
2402  ! 2.d Perform output type 2 ( statistics )
2403  !
2404  ELSE IF ( itype .EQ. 2 ) THEN
2405  xmin = 1.e20
2406  xmax = -1.e20
2407  xds = 0.d0
2408  xdsq = 0.d0
2409  ningrd = 0
2410  !
2411  DO ix=ix1, ixn
2412  DO iy=iy1, iyn
2413  IF ( x1(ix,iy) .NE. undef ) THEN
2414  ningrd = ningrd + 1
2415  xmin = min( xmin , x1(ix,iy) )
2416  xmax = max( xmax , x1(ix,iy) )
2417  xds = xds + dble(x1(ix,iy))
2418  xdsq = xdsq + dble(x1(ix,iy))**2
2419  END IF
2420  END DO
2421  END DO
2422  !
2423  ndsdt = ndsdt + 1
2424  !
2425  IF ( ningrd .EQ. 0 ) THEN
2426  WRITE (ndsdt,940) time(1), ih, im, is
2427  ELSE IF ( ningrd .LE. 2 ) THEN
2428  xavg = real( xds / dble(ningrd) )
2429  WRITE (ndsdt,940) time(1), ih, im, is, &
2430  xmin, xmax
2431  ELSE
2432  xavg = real( xds / dble(ningrd) )
2433  xstd = real( ( xdsq - xds**2/dble(ningrd) ) &
2434  / dble(ningrd-1) )
2435  xstd = sqrt( max( xstd , 0. ) )
2436  WRITE (ndsdt,940) time(1), ih, im, is, &
2437  xmin, xmax, xavg, xstd
2438  END IF
2439  !
2440  ! 2.e Perform output type 3 ( file )
2441  !
2442  ELSE IF ( itype .EQ. 3 ) THEN
2443  !
2444  fname(13:) = ename
2445  IF ( idfm .EQ. 3 ) THEN
2446  IF(gtype .NE. ungtype) THEN
2447  jj = len_trim(fnmpre)
2448  OPEN (ndsdat,file=fnmpre(:jj)//fname, &
2449  form='UNFORMATTED', convert=file_endian,err=800,iostat=ierr)
2450  WRITE (ndsdat) fileid, time, &
2451  minval(xgrd(iy1:iyn,ix1:ixn)), &
2452  maxval(xgrd(iy1:iyn,ix1:ixn)), ixn-ix1+1, &
2453  minval(ygrd(iy1:iyn,ix1:ixn)), &
2454  maxval(ygrd(iy1:iyn,ix1:ixn)), iyn-iy1+1, &
2455  ename, fsc, units, idla, idfm, formf, mfill
2456  ELSE
2457  OPEN (ndsdat,file=fname, &
2458  form='UNFORMATTED', convert=file_endian,err=800,iostat=ierr)
2459  WRITE (ndsdat) fileid, time, &
2460  x0,maxx,nx, &
2461  y0,maxy,ny, &
2462  ename, fsc, units, idla, idfm, formf, mfill
2463  ENDIF
2464  ELSE
2465  IF(gtype .NE. ungtype) THEN
2466  jj = len_trim(fnmpre)
2467  OPEN (ndsdat,file=fnmpre(:jj)//fname,err=800, &
2468  iostat=ierr)
2469  IF (fsc.LT.1e-4) THEN
2470  WRITE(fscs,'(G8.1)') fsc
2471  ELSE
2472  WRITE(fscs,'(F7.4)') fsc
2473  END IF
2474  IF ( flagll ) THEN
2475  WRITE (ndsdat,950) fileid, time, &
2476  minval(xgrd(iy1:iyn,ix1:ixn)), &
2477  maxval(xgrd(iy1:iyn,ix1:ixn)), ixn-ix1+1, &
2478  minval(ygrd(iy1:iyn,ix1:ixn)), &
2479  maxval(ygrd(iy1:iyn,ix1:ixn)), iyn-iy1+1, &
2480  ename, fscs, units, idla, idfm, formf, mfill
2481  ELSE
2482  WRITE (ndsdat,960) fileid, time, &
2483  minval(xgrd(iy1:iyn,ix1:ixn)), &
2484  maxval(xgrd(iy1:iyn,ix1:ixn)), ixn-ix1+1, &
2485  minval(ygrd(iy1:iyn,ix1:ixn)), &
2486  maxval(ygrd(iy1:iyn,ix1:ixn)), iyn-iy1+1, &
2487  ename, fscs, units, idla, idfm, formf, mfill
2488  END IF
2489  ELSE
2490  OPEN (ndsdat,file=fname, &
2491  err=800,iostat=ierr)
2492  WRITE (ndsdat, 949) fileid, time, &
2493  x0,maxx,nx, &
2494  y0,maxy,ny, &
2495  ename, fsc, units, idla, idfm, formf, mfill
2496  ENDIF
2497  END IF
2498  !
2499  IF ( fltri ) THEN
2500  DO ix=ix1, ixn
2501  DO iy=iy1, iyn
2502  IF ( xx(ix,iy) .EQ. undef ) THEN
2503  mxx(ix,iy) = mfill
2504  myy(ix,iy) = mfill
2505  mxy(ix,iy) = mfill
2506  ELSE
2507  mxx(ix,iy) = nint(x1(ix,iy)/fsc)
2508  myy(ix,iy) = nint(x2(ix,iy)/fsc)
2509  mxy(ix,iy) = nint(xy(ix,iy)/fsc)
2510  END IF
2511  END DO
2512  END DO
2513  IF ( idla .NE. 5 ) THEN
2514  CALL outa2i ( mxx, nx, ny, ix1, ixn, iy1, iyn, &
2515  ndsdat, ndst, ndse, idfm, formf, idla, 1, 0 )
2516  CALL outa2i ( myy, nx, ny, ix1, ixn, iy1, iyn, &
2517  ndsdat, ndst, ndse, idfm, formf, idla, 1, 0 )
2518  CALL outa2i ( mxy, nx, ny, ix1, ixn, iy1, iyn, &
2519  ndsdat, ndst, ndse, idfm, formf, idla, 1, 0 )
2520  ELSE
2521  DO iy=iy1,iyn
2522  ygbx = y0 + real(iy-1)*sy
2523  DO ix=ix1, ixn
2524  xgbx = x0 + real(ix-1)*sx
2525  IF ( mxx(ix,iy) .NE. mfill ) THEN
2526  IF ( idfm .EQ. 3 ) THEN
2527  WRITE (ndsdat) &
2528  xgbx, ygbx, mxx(ix,iy), myy(ix,iy)
2529  ELSE
2530  WRITE (ndsdat,951) &
2531  xgbx, ygbx, mxx(ix,iy), myy(ix,iy)
2532  END IF
2533  END IF
2534  END DO
2535  END DO
2536  END IF
2537  ELSE
2538  IF ( fltwo .OR. fldir ) THEN
2539  DO ix=ix1, ixn
2540  DO iy=iy1, iyn
2541  IF ( xx(ix,iy) .EQ. undef ) THEN
2542  mxx(ix,iy) = mfill
2543  myy(ix,iy) = mfill
2544  ELSE
2545  mxx(ix,iy) = nint(xx(ix,iy)/fsc)
2546  myy(ix,iy) = nint(xy(ix,iy)/fsc)
2547  END IF
2548  END DO
2549  END DO
2550  IF ( idla .NE. 5 ) THEN
2551  CALL outa2i ( mxx, nx, ny, ix1, ixn, iy1, iyn, &
2552  ndsdat, ndst, ndse, idfm, formf, idla, 1,0)
2553  CALL outa2i ( myy, nx, ny, ix1, ixn, iy1, iyn, &
2554  ndsdat, ndst, ndse, idfm, formf, idla, 1,0)
2555  ELSE
2556  DO iy=iy1,iyn
2557  DO ix=ix1, ixn
2558  ygbx = ygrd(iy,ix)
2559  xgbx = xgrd(iy,ix)
2560  IF ( mxx(ix,iy) .NE. mfill ) THEN
2561  IF ( idfm .EQ. 3 ) THEN
2562  WRITE (ndsdat) &
2563  xgbx, ygbx, mxx(ix,iy), myy(ix,iy)
2564  ELSE
2565  IF ( flagll ) THEN
2566  WRITE (ndsdat,951) xgbx, ygbx, &
2567  mxx(ix,iy), myy(ix,iy)
2568  ELSE
2569  WRITE (ndsdat,961) xgbx, ygbx, &
2570  mxx(ix,iy), myy(ix,iy)
2571  END IF
2572  END IF
2573  END IF
2574  END DO
2575  END DO
2576  END IF
2577  ELSE
2578  DO ix=ix1, ixn
2579  DO iy=iy1, iyn
2580  IF ( x1(ix,iy) .EQ. undef ) THEN
2581  mx1(ix,iy) = mfill
2582  ELSE
2583  mx1(ix,iy) = nint(x1(ix,iy)/fsc)
2584  END IF
2585  END DO
2586  END DO
2587  IF ( idla .NE. 5 ) THEN
2588  CALL outa2i ( mx1, nx, ny, ix1, ixn, iy1, iyn, &
2589  ndsdat, ndst, ndse, idfm, formf, idla, 1,0)
2590  ELSE
2591  DO iy=iy1,iyn
2592  DO ix=ix1, ixn
2593  ygbx = ygrd(iy,ix)
2594  xgbx = xgrd(iy,ix)
2595  IF ( mx1(ix,iy) .NE. mfill ) THEN
2596  IF ( idfm .EQ. 3 ) THEN
2597  WRITE (ndsdat) &
2598  xgbx, ygbx, mx1(ix,iy)
2599  ELSE
2600  IF ( flagll ) THEN
2601  WRITE (ndsdat,951) xgbx, ygbx, &
2602  mx1(ix,iy)
2603  ELSE
2604  WRITE (ndsdat,961) xgbx, ygbx, &
2605  mx1(ix,iy)
2606  END IF
2607  END IF
2608  END IF
2609  END DO
2610  END DO
2611  END IF
2612  END IF
2613  END IF
2614  !
2615  CLOSE (ndsdat)
2616  !
2617  ELSE IF ( itype .EQ. 4 ) THEN
2618  !
2619  fname(13:) = ename
2620  jj = len_trim(fnmpre)
2621  OPEN (ndsdat,file=fnmpre(:jj)//fname,err=800, &
2622  iostat=ierr)
2623  WRITE (6,*) fname(1:16)
2624  !
2625  IF ( fltri ) THEN
2626  WRITE (ndsdat,980) fileid, time, nsea, 3, &
2627  fsc, ename, units, gname
2628  WRITE(ndsdat, 113) xs1
2629  WRITE(ndsdat, 113) xs2
2630  WRITE(ndsdat, 113) xs3
2631  ENDIF
2632  IF ( fltwo .OR. fldir ) THEN
2633  WRITE (ndsdat,980) fileid, time, nsea, 2, &
2634  fsc, ename, units, gname
2635  WRITE(ndsdat, 113) xs1
2636  WRITE(ndsdat, 113) xs2
2637  ENDIF
2638  IF ( flone ) THEN
2639  WRITE (ndsdat,980) fileid, time, nsea, 1, &
2640  fsc, ename, units, gname
2641  WRITE(ndsdat, 113) xs1
2642  ENDIF
2643  !
2644  CLOSE (ndsdat)
2645  !
2646  END IF
2647  !
2648  ! ... End of fields loop
2649  !
2650  END IF
2651  END DO
2652  END DO
2653  !
2654  RETURN
2655  !
2656  ! Error escape locations
2657  !
2658 800 CONTINUE
2659  WRITE (ndse,1000) ierr
2660  CALL extcde (2)
2661  !
2662  ! Formats
2663  !
2664 113 FORMAT ((10es11.3))
2665 980 FORMAT (1x,a13,i9.8,i7.6,i9,i3,es10.2,1x,a4,1x,a10,1x,a30)
2666 
2667 940 FORMAT (1x,i8,3i3.2,2x,4e12.4)
2668 949 FORMAT (1x,a13,i9.8,i7.6,2(2f8.2,i8), &
2669  1x,a4,f8.4,1x,a10,2i2,1x,a11,i4)
2670 950 FORMAT (1x,a13,1x,i9.8,1x,i7.6,2(1x,2f8.2,1x,i4), &
2671  1x,a4,1x,a7,1x,a10,1x,2i2,1x,a11,1x,i4)
2672 951 FORMAT (1x,2f10.5,2i8)
2673 960 FORMAT (1x,a13,i9.8,i7.6,2(2e11.3,i4), &
2674  1x,a4,1x,a7,1x,a10,2i2,1x,a11,i4)
2675 961 FORMAT (1x,2e12.4,2i8)
2676  !
2677 990 FORMAT (/' *** WAVEWATCH III ERROR IN W3EXGO :'/ &
2678  ' GROUP',i2,' PARAMETER',i3,' NOT LISTED ' )
2679 999 FORMAT (/' *** WAVEWATCH III ERROR IN W3EXGO :'/ &
2680  ' PLEASE UPDATE FIELDS !!! '/)
2681  !
2682 1000 FORMAT (/' *** WAVEWATCH III ERROR IN W3EXGO : '/ &
2683  ' ERROR IN OPENING OUTPUT FILE'/ &
2684  ' IOSTAT =',i5/)
2685  !
2686 #ifdef W3_T
2687 9000 FORMAT (' TEST W3EXGO : FLAGS :',i3,2x,20l2)
2688 9001 FORMAT (' TEST W3EXGO : ITPYE :',i4/ &
2689  ' IX1/N/S :',3i4/ &
2690  ' IY1/N/S :',3i4/ &
2691  ' SCALE, VECTOR :',2l2/ &
2692  ' NDSDAT :',i4)
2693 #endif
2694  !
2695 #ifdef W3_T
2696 9012 FORMAT (' TEST W3EXGO : BLOK PARS : ',3i4)
2697 9014 FORMAT (' BASE NAME : ',a)
2698 #endif
2699  !
2700 #ifdef W3_T
2701 9020 FORMAT (' TEST W3EXGO : OUTPUT FIELD : ',a)
2702 #endif
2703  !/
2704  !/ End of W3EXGO ----------------------------------------------------- /
2705  !/
2706  END SUBROUTINE w3exgo
2707  !/
2708  !
2709  !/ End of W3OUTF ----------------------------------------------------- /
2710  !/
2711 END PROGRAM w3outf
w3adatmd::pt2
real, dimension(:,:), pointer pt2
Definition: w3adatmd.F90:597
w3adatmd::psw
real, dimension(:,:), pointer psw
Definition: w3adatmd.F90:597
w3adatmd::hcmaxe
real, dimension(:), pointer hcmaxe
Definition: w3adatmd.F90:587
w3servmd::nextln
subroutine nextln(CHCKC, NDSI, NDSE)
Definition: w3servmd.F90:222
w3wdatmd::iceh
real, dimension(:), pointer iceh
Definition: w3wdatmd.F90:183
w3adatmd::charn
real, dimension(:), pointer charn
Definition: w3adatmd.F90:603
w3adatmd::dtdyn
real, dimension(:), pointer dtdyn
Definition: w3adatmd.F90:620
w3adatmd::as
real, dimension(:), pointer as
Definition: w3adatmd.F90:584
w3adatmd
Define data structures to set up wave model auxiliary data for several models simultaneously.
Definition: w3adatmd.F90:26
w3adatmd::hcmaxd
real, dimension(:), pointer hcmaxd
Definition: w3adatmd.F90:587
w3arrymd::outa2i
subroutine outa2i(ARRAY, MX, MY, LX, HX, LY, HY, NDS, NDST, NDSE, IDFM, RFORM, IDLA, VSC, VOF)
Definition: w3arrymd.F90:627
w3adatmd::ussy
real, dimension(:), pointer ussy
Definition: w3adatmd.F90:607
w3adatmd::pep
real, dimension(:,:), pointer pep
Definition: w3adatmd.F90:597
w3adatmd::mscd
real, dimension(:), pointer mscd
Definition: w3adatmd.F90:617
w3adatmd::abd
real, dimension(:), pointer abd
Definition: w3adatmd.F90:614
w3wdatmd
Define data structures to set up wave model dynamic data for several models simultaneously.
Definition: w3wdatmd.F90:18
w3adatmd::stmaxe
real, dimension(:), pointer stmaxe
Definition: w3adatmd.F90:587
w3adatmd::t02
real, dimension(:), pointer t02
Definition: w3adatmd.F90:587
w3adatmd::cflxymax
real, dimension(:), pointer cflxymax
Definition: w3adatmd.F90:620
w3adatmd::tusx
real, dimension(:), pointer tusx
Definition: w3adatmd.F90:607
w3adatmd::fcut
real, dimension(:), pointer fcut
Definition: w3adatmd.F90:620
w3adatmd::tusy
real, dimension(:), pointer tusy
Definition: w3adatmd.F90:607
w3adatmd::ptp
real, dimension(:,:), pointer ptp
Definition: w3adatmd.F90:597
w3adatmd::dw
real, dimension(:), pointer dw
Definition: w3adatmd.F90:584
w3wdatmd::icef
real, dimension(:), pointer icef
Definition: w3wdatmd.F90:183
w3adatmd::t01
real, dimension(:), pointer t01
Definition: w3adatmd.F90:587
w3adatmd::cge
real, dimension(:), pointer cge
Definition: w3adatmd.F90:603
w3wdatmd::wlv
real, dimension(:), pointer wlv
Definition: w3wdatmd.F90:183
w3odatmd::ngrpp
integer, parameter ngrpp
Definition: w3odatmd.F90:324
w3adatmd::pdir
real, dimension(:,:), pointer pdir
Definition: w3adatmd.F90:597
w3wdatmd::time
integer, dimension(:), pointer time
Definition: w3wdatmd.F90:172
w3adatmd::thp0
real, dimension(:), pointer thp0
Definition: w3adatmd.F90:587
w3adatmd::tauocy
real, dimension(:), pointer tauocy
Definition: w3adatmd.F90:607
constants::rade
real, parameter rade
RADE Conversion factor from radians to degrees.
Definition: constants.F90:76
w3adatmd::phs
real, dimension(:,:), pointer phs
Definition: w3adatmd.F90:597
w3adatmd::hs
real, dimension(:), pointer hs
Definition: w3adatmd.F90:587
w3adatmd::uba
real, dimension(:), pointer uba
Definition: w3adatmd.F90:614
w3odatmd::fnmpre
character(len=80) fnmpre
Definition: w3odatmd.F90:330
w3adatmd::pqp
real, dimension(:,:), pointer pqp
Definition: w3adatmd.F90:597
w3adatmd::tauwix
real, dimension(:), pointer tauwix
Definition: w3adatmd.F90:603
w3adatmd::pthp0
real, dimension(:,:), pointer pthp0
Definition: w3adatmd.F90:597
w3gdatmd::w3setg
subroutine w3setg(IMOD, NDSE, NDST)
Definition: w3gdatmd.F90:2152
w3adatmd::tauwiy
real, dimension(:), pointer tauwiy
Definition: w3adatmd.F90:603
w3adatmd::taubbl
real, dimension(:,:), pointer taubbl
Definition: w3adatmd.F90:614
w3odatmd::ndse
integer, pointer ndse
Definition: w3odatmd.F90:456
w3wdatmd::berg
real, dimension(:), pointer berg
Definition: w3wdatmd.F90:183
w3adatmd::w3seta
subroutine w3seta(IMOD, NDSE, NDST)
Select one of the WAVEWATCH III grids / models.
Definition: w3adatmd.F90:2645
w3adatmd::tauwny
real, dimension(:), pointer tauwny
Definition: w3adatmd.F90:603
w3servmd::w3thrtn
subroutine w3thrtn(NSEA, THETA, AnglD, Degrees)
Definition: w3servmd.F90:1333
w3wdatmd::w3ndat
subroutine w3ndat(NDSE, NDST)
Set up the number of grids to be used.
Definition: w3wdatmd.F90:210
w3adatmd::phiaw
real, dimension(:), pointer phiaw
Definition: w3adatmd.F90:603
w3adatmd::pws
real, dimension(:,:), pointer pws
Definition: w3adatmd.F90:597
w3adatmd::plp
real, dimension(:,:), pointer plp
Definition: w3adatmd.F90:597
w3adatmd::phibbl
real, dimension(:), pointer phibbl
Definition: w3adatmd.F90:614
w3adatmd::cflthmax
real, dimension(:), pointer cflthmax
Definition: w3adatmd.F90:620
w3adatmd::skew
real, dimension(:), pointer skew
Definition: w3adatmd.F90:617
w3adatmd::psi
real, dimension(:,:), pointer psi
Definition: w3adatmd.F90:597
w3adatmd::tpms
real, dimension(:), pointer tpms
Definition: w3adatmd.F90:607
w3servmd
Definition: w3servmd.F90:3
w3adatmd::embia1
real, dimension(:), pointer embia1
Definition: w3adatmd.F90:617
w3odatmd::flogrd
logical, dimension(:,:), pointer flogrd
Definition: w3odatmd.F90:478
w3adatmd::ths
real, dimension(:), pointer ths
Definition: w3adatmd.F90:587
w3adatmd::bedforms
real, dimension(:,:), pointer bedforms
Definition: w3adatmd.F90:614
w3adatmd::ud
real, dimension(:), pointer ud
Definition: w3adatmd.F90:584
w3adatmd::hmaxd
real, dimension(:), pointer hmaxd
Definition: w3adatmd.F90:587
w3adatmd::pwst
real, dimension(:), pointer pwst
Definition: w3adatmd.F90:597
w3wdatmd::w3setw
subroutine w3setw(IMOD, NDSE, NDST)
Select one of the WAVEWATCH III grids / models.
Definition: w3wdatmd.F90:660
w3adatmd::qkk
real, dimension(:), pointer qkk
Definition: w3adatmd.F90:617
w3odatmd::noge
integer, dimension(nogrp) noge
Definition: w3odatmd.F90:326
w3odatmd::w3seto
subroutine w3seto(IMOD, NDSERR, NDSTST)
Definition: w3odatmd.F90:1523
w3adatmd::sxy
real, dimension(:), pointer sxy
Definition: w3adatmd.F90:607
w3adatmd::tauwnx
real, dimension(:), pointer tauwnx
Definition: w3adatmd.F90:603
w3odatmd
Definition: w3odatmd.F90:3
w3adatmd::wbt
real, dimension(:), pointer wbt
Definition: w3adatmd.F90:587
w3adatmd::bhd
real, dimension(:), pointer bhd
Definition: w3adatmd.F90:607
w3adatmd::w3naux
subroutine w3naux(NDSE, NDST)
Set up the number of grids to be used.
Definition: w3adatmd.F90:704
w3adatmd::cy
real, dimension(:), pointer cy
Definition: w3adatmd.F90:584
w3adatmd::pnr
real, dimension(:), pointer pnr
Definition: w3adatmd.F90:597
w3adatmd::taua
real, dimension(:), pointer taua
Definition: w3adatmd.F90:584
w3adatmd::hmaxe
real, dimension(:), pointer hmaxe
Definition: w3adatmd.F90:587
w3adatmd::pt1
real, dimension(:,:), pointer pt1
Definition: w3adatmd.F90:597
w3adatmd::wlm
real, dimension(:), pointer wlm
Definition: w3adatmd.F90:587
w3iogrmd::w3iogr
subroutine w3iogr(INXOUT, NDSM, IMOD, FEXT ifdef W3_ASCII
Reading and writing of the model definition file.
Definition: w3iogrmd.F90:117
w3adatmd::wnmean
real, dimension(:), pointer wnmean
Definition: w3adatmd.F90:587
w3adatmd::cflkmax
real, dimension(:), pointer cflkmax
Definition: w3adatmd.F90:620
w3iogomd::w3iogo
subroutine w3iogo(INXOUT, NDSOG, IOTST, IMOD ifdef W3_ASCII
Read/write gridded output.
Definition: w3iogomd.F90:2396
w3adatmd::phioc
real, dimension(:), pointer phioc
Definition: w3adatmd.F90:607
file
file(STRINGS ${CMAKE_BINARY_DIR}/switch switch_strings) separate_arguments(switches UNIX_COMMAND $
Definition: CMakeLists.txt:3
w3iogrmd
Reading/writing of model definition file.
Definition: w3iogrmd.F90:20
w3adatmd::qp
real, dimension(:), pointer qp
Definition: w3adatmd.F90:587
w3iogomd
Gridded output of mean wave parameters.
Definition: w3iogomd.F90:15
w3adatmd::stmaxd
real, dimension(:), pointer stmaxd
Definition: w3adatmd.F90:587
w3servmd::strace
subroutine strace(IENT, SNAME)
Definition: w3servmd.F90:148
w3odatmd::ndso
integer, pointer ndso
Definition: w3odatmd.F90:456
w3wdatmd::ice
real, dimension(:), pointer ice
Definition: w3wdatmd.F90:183
w3gdatmd::w3nmod
subroutine w3nmod(NUMBER, NDSE, NDST, NAUX)
Definition: w3gdatmd.F90:1433
w3arrymd
Definition: w3arrymd.F90:3
w3adatmd::whitecap
real, dimension(:,:), pointer whitecap
Definition: w3adatmd.F90:603
w3adatmd::tauox
real, dimension(:), pointer tauox
Definition: w3adatmd.F90:607
w3odatmd::idout
character(len=20), dimension(nogrp, ngrpp) idout
Definition: w3odatmd.F90:329
w3adatmd::prms
real, dimension(:), pointer prms
Definition: w3adatmd.F90:607
w3adatmd::usero
real, dimension(:,:), pointer usero
Definition: w3adatmd.F90:623
w3iogomd::w3readflgrd
subroutine w3readflgrd(NDSI, NDSO, NDSS, NDSEN, COMSTR, FLG1D, FLG2D, IAPROC, NAPOUT, IERR)
Fills in FLG1D and FLG2D arrays from ASCII input file.
Definition: w3iogomd.F90:336
w3adatmd::mssy
real, dimension(:), pointer mssy
Definition: w3adatmd.F90:617
w3odatmd::ndst
integer, pointer ndst
Definition: w3odatmd.F90:456
w3wdatmd::ust
real, dimension(:), pointer ust
Definition: w3wdatmd.F90:183
constants
Define some much-used constants for global use (all defined as PARAMETER).
Definition: constants.F90:20
w3adatmd::ua
real, dimension(:), pointer ua
Definition: w3adatmd.F90:584
w3adatmd::tauoy
real, dimension(:), pointer tauoy
Definition: w3adatmd.F90:607
w3odatmd::nogrp
integer, parameter nogrp
Definition: w3odatmd.F90:323
w3adatmd::fp0
real, dimension(:), pointer fp0
Definition: w3adatmd.F90:587
w3gdatmd
Definition: w3gdatmd.F90:16
w3adatmd::hsig
real, dimension(:), pointer hsig
Definition: w3adatmd.F90:587
w3servmd::w3xyrtn
subroutine w3xyrtn(NSEA, XVEC, YVEC, AnglD)
Definition: w3servmd.F90:1387
w3adatmd::ussx
real, dimension(:), pointer ussx
Definition: w3adatmd.F90:607
w3adatmd::embia2
real, dimension(:), pointer embia2
Definition: w3adatmd.F90:617
w3adatmd::mssx
real, dimension(:), pointer mssx
Definition: w3adatmd.F90:617
constants::file_endian
character(*), parameter file_endian
FILE_ENDIAN Filled by preprocessor with 'big_endian', 'little_endian', or 'native'.
Definition: constants.F90:86
w3adatmd::tauadir
real, dimension(:), pointer tauadir
Definition: w3adatmd.F90:584
w3servmd::extcde
subroutine extcde(IEXIT, UNIT, MSG, FILE, LINE, COMM)
Definition: w3servmd.F90:736
w3odatmd::w3nout
subroutine w3nout(NDSERR, NDSTST)
Definition: w3odatmd.F90:561
w3outf
program w3outf
Post-processing of grid output.
Definition: ww3_outf.F90:35
w3adatmd::mscy
real, dimension(:), pointer mscy
Definition: w3adatmd.F90:617
w3wdatmd::rhoair
real, dimension(:), pointer rhoair
Definition: w3wdatmd.F90:183
w3servmd::itrace
subroutine itrace(NDS, NMAX)
Definition: w3servmd.F90:91
w3wdatmd::ustdir
real, dimension(:), pointer ustdir
Definition: w3wdatmd.F90:183
w3adatmd::ptm1
real, dimension(:,:), pointer ptm1
Definition: w3adatmd.F90:597
w3odatmd::noswll
integer, pointer noswll
Definition: w3odatmd.F90:460
w3adatmd::thm
real, dimension(:), pointer thm
Definition: w3adatmd.F90:587
w3adatmd::mscx
real, dimension(:), pointer mscx
Definition: w3adatmd.F90:617
w3adatmd::ppe
real, dimension(:,:), pointer ppe
Definition: w3adatmd.F90:597
w3adatmd::cx
real, dimension(:), pointer cx
Definition: w3adatmd.F90:584
w3adatmd::pgw
real, dimension(:,:), pointer pgw
Definition: w3adatmd.F90:597
w3timemd
Definition: w3timemd.F90:3
constants::undef
real undef
UNDEF Value for undefined variable in output.
Definition: constants.F90:84
w3exgo
subroutine w3exgo(NX, NY, NSEA)
Perform actual grid output.
Definition: ww3_outf.F90:548
w3adatmd::mssd
real, dimension(:), pointer mssd
Definition: w3adatmd.F90:617
w3servmd::w3s2xy
subroutine w3s2xy(NSEA, MSEA, MX, MY, S, MAPSF, XY)
Definition: w3servmd.F90:337
w3adatmd::t0m1
real, dimension(:), pointer t0m1
Definition: w3adatmd.F90:587
w3arrymd::prtblk
subroutine prtblk(NDS, NX, NY, MX, F, MAP, MAP0, FSC, IX1, IX2, IX3, IY1, IY2, IY3, PRVAR, PRUNIT)
Definition: w3arrymd.F90:1112
w3adatmd::tauocx
real, dimension(:), pointer tauocx
Definition: w3adatmd.F90:607
w3adatmd::aba
real, dimension(:), pointer aba
Definition: w3adatmd.F90:614
w3adatmd::syy
real, dimension(:), pointer syy
Definition: w3adatmd.F90:607
w3adatmd::ubd
real, dimension(:), pointer ubd
Definition: w3adatmd.F90:614
w3adatmd::sxx
real, dimension(:), pointer sxx
Definition: w3adatmd.F90:607