WAVEWATCH III  beta 0.0.1
ww3_trnc.F90
Go to the documentation of this file.
1 
5 
6 #include "w3macros.h"
7 !/ ------------------------------------------------------------------- /
15 PROGRAM w3trnc
16  !/
17  !/ +-----------------------------------+
18  !/ | WAVEWATCH III NOAA/NCEP |
19  !/ | M. Accensi |
20  !/ | FORTRAN 90 |
21  !/ | Last update : 15-May-2018 |
22  !/ +-----------------------------------+
23  !/
24  !/ 17-Feb-2016 : Creation ( version 5.11 )
25  !/ 11-Apr-2016 : Adapted to use more options ( version 5.11 )
26  !/ 15-May-2018 : Add namelist feature ( version 6.05 )
27  !/ 18-Jun-2020 : Support for 360-day calendar. ( version 7.08 )
28  !/
29  !/ Copyright 2014 National Weather Service (NWS),
30  !/ National Oceanic and Atmospheric Administration. All rights
31  !/ reserved. WAVEWATCH III is a trademark of the NWS.
32  !/ No unauthorized use without permission.
33  !/
34  ! 1. Purpose :
35  !
36  ! Convert direct access track output file to netCDF file.
37  !
38  ! 2. Method :
39  !
40  ! Info read from track_o.ww3, written to track.nc
41  !
42  ! 3. Parameters :
43  !
44  ! 4. Subroutines used :
45  !
46  ! Name Type Module Description
47  ! ----------------------------------------------------------------
48  ! W3NMOD Subr. W3GDATMD Set number of model.
49  ! W3NOUT Subr. W3ODATMD Set number of model for output.
50  ! W3IOGR Subr. W3IOGRMD Reading/writing model definition file.
51  ! ----------------------------------------------------------------
52  !
53  ! 5. Called by :
54  !
55  ! None, stand-alone program.
56  !
57  ! 6. Error messages :
58  !
59  ! 7. Remarks :
60  !
61  ! 8. Structure :
62  !
63  ! See source code.
64  !
65  ! 9. Switches :
66  !
67  ! !/S Enable subroutine tracing.
68  !
69  ! 10. Source code :
70  !
71  !/ ------------------------------------------------------------------- /
72  USE constants
73 
74 #ifdef W3_NL1
75  USE w3adatmd, ONLY : w3naux, w3seta
76 #endif
77  USE w3gdatmd, ONLY : w3nmod, w3setg, flagll, xfr, gname
78  USE w3odatmd, ONLY : w3nout, w3seto, fnmpre
79  USE w3servmd, ONLY : itrace, nextln, extcde
80 #ifdef W3_S
81  USE w3servmd, ONLY : strace
82 #endif
83  USE w3timemd
84  USE w3iogrmd, ONLY: w3iogr
85  !
86  USE w3odatmd, ONLY: ndso, ndse
87  !
88  USE w3nmltrncmd
89  USE netcdf
90  !
91  IMPLICIT NONE
92  !/
93  !/ ------------------------------------------------------------------- /
94  !/ Local parameters
95  !/
96  TYPE(nml_track_t) :: nml_track
97  TYPE(nml_file_t) :: nml_file
98  !
99  INTEGER :: ndsi, ndsinp, ndsm, &
100  ndsout, ndstrc, ntrace, &
101  nspec, ierr, mk, mth, it, &
102  iloc, ispec, s3, iout, &
103  iret, nctype,ncid, ith
104 
105  INTEGER :: time(2), tout(2), nout, tdum(2), &
106  dimid(4), varid(18), dimln(4), &
107  stopdate(8)
108 #ifdef W3_S
109  INTEGER, SAVE :: ient = 0
110 #endif
111  !
112  REAL :: th1, dth, x, y, dw, cx, cy, cao, cdo,&
113  wx, wy, wao, wdo, ust, as, dtest, &
114  dtreq, dthd, rth0, m2km
115  !
116  REAL, ALLOCATABLE :: freq(:), freq1(:), freq2(:), dsip(:),&
117  spec(:,:), e(:,:), thd(:), dir(:)
118  !
119  CHARACTER*34, PARAMETER :: &
120  idtst = 'WAVEWATCH III TRACK OUTPUT SPECTRA'
121  CHARACTER*30 :: fileprefix, strstopdate
122  CHARACTER*20 :: format1
123  CHARACTER :: idtime*23, iddday*11, trckid*32, &
124  comstr*1, idstr*34, tststr*3, stime*23
125  !
126  LOGICAL :: flgnml
127 
128 
129  !/
130  !/ ------------------------------------------------------------------- /
131  !/
132  !
133  ! 0. Initialize data structure
134  !
135  CALL w3nmod ( 1, 6, 6 )
136  CALL w3setg ( 1, 6, 6 )
137 #ifdef W3_NL1
138  CALL w3naux ( 6, 6 )
139  CALL w3seta ( 1, 6, 6 )
140 #endif
141  CALL w3nout ( 6, 6 )
142  CALL w3seto ( 1, 6, 6 )
143  !
144  ! 1. IO set-up.
145  !
146  ndsi = 10
147  ndsm = 20
148  ndsinp = 11
149  ndsout = 51
150  !
151  ndstrc = 6
152  ntrace = 10
153  CALL itrace ( ndstrc, ntrace )
154  !
155 #ifdef W3_S
156  CALL strace ( ient, 'W3TRNC' )
157 #endif
158  !
159  WRITE (ndso,900)
160  !
161  !--- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
162  ! 2. Read model definition file.
163  !
164  CALL w3iogr ( 'READ', ndsm )
165  WRITE (ndso,920) gname
166 
167 
168  !
169  !--- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
170  ! 3. Read requests from input file.
171  !
172 
173  !
174  ! process ww3_trnc namelist
175  !
176  INQUIRE(file=trim(fnmpre)//"ww3_trnc.nml", exist=flgnml)
177  IF (flgnml) THEN
178  ! Read namelist
179  CALL w3nmltrnc (ndsi, trim(fnmpre)//'ww3_trnc.nml', nml_track, nml_file, ierr)
180 
181  ! 3.1 Time setup IDTIME, DTREQ, NOUT
182  READ(nml_track%TIMESTRIDE, *) dtreq
183  READ(nml_track%TIMECOUNT, *) nout
184  READ(nml_track%TIMESTART, *) tout(1), tout(2)
185 
186 
187  ! 3.2 Output type
188  nctype = nml_file%NETCDF
189  fileprefix = nml_file%PREFIX
190  s3 = nml_track%TIMESPLIT
191 
192 
193  END IF ! FLGNML
194 
195  !
196  ! process old ww3_trnc.inp format
197  !
198  IF (.NOT. flgnml) THEN
199  OPEN (ndsi,file=trim(fnmpre)//'ww3_trnc.inp',status='OLD',err=805,iostat=ierr)
200  rewind(ndsi)
201 
202  READ (ndsi,'(A)',END=806,ERR=807,IOSTAT=IERR) comstr
203  IF (comstr.EQ.' ') comstr = '$'
204  WRITE (ndso,901) comstr
205 
206 
207  ! 3.1 Time setup IDTIME, DTREQ, NOUT
208  CALL nextln ( comstr , ndsi , ndse )
209  READ (ndsi,*,END=806,ERR=807) TOUT, DTREQ, nout
210 
211 
212  ! 3.2 Output type
213  CALL nextln ( comstr , ndsi , ndse )
214  READ (ndsi,*,END=806,ERR=807) nctype
215  CALL nextln ( comstr , ndsi , ndse )
216  fileprefix= 'ww3.'
217  READ (ndsi,*,END=806,ERR=807) fileprefix
218  CALL nextln ( comstr , ndsi , ndse )
219  READ (ndsi,*,END=806,ERR=807) s3
220 
221  END IF ! .NOT. FLGNML
222 
223 
224 
225 
226  ! 3.3 Time setup IDTIME, DTREQ, NOUT
227  dtreq = max( 0. , dtreq )
228  IF ( dtreq.EQ.0. ) nout = 1
229  nout = max( 1 , nout )
230  CALL stme21 ( tout , idtime )
231  WRITE (ndso,940) idtime
232  tdum = 0
233  CALL tick21 ( tdum , dtreq )
234  CALL stme21 ( tdum , idtime )
235  IF ( dtreq .GE. 86400. ) THEN
236  WRITE (iddday,'(I10,1X)') int(dtreq/86400.)
237  ELSE
238  iddday = ' '
239  END IF
240  idtime(1:11) = iddday
241  idtime(21:23) = ' '
242  WRITE (ndso,941) idtime, nout
243 
244 
245  ! 3.4 Output type
246  IF ( nctype.LT.3 .OR. nctype.GT.4 ) THEN
247  WRITE (ndse,1010) nctype
248  CALL extcde ( 1 )
249  END IF
250  ! S3 defines the number of characters in the date for the filename
251  ! S3=4-> YYYY, S3=6 -> YYYYMM, S3=10 -> YYYYMMDDTHHZ ...
252 
253 
254 
255  !
256  !--- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
257  ! 4. Check consistency with input file and track_o.ww3
258  !
259  OPEN (ndsinp,file=trim(fnmpre)//'track_o.ww3',form='UNFORMATTED', convert=file_endian, &
260  status='OLD',err=800,iostat=ierr)
261  READ (ndsinp,err=801,iostat=ierr) idstr, flagll, mk, mth, xfr
262  !
263  IF ( flagll ) THEN
264  m2km = 1.
265  ELSE
266  m2km = 1.e-3
267  END IF
268  !
269  IF ( idstr .NE. idtst ) GOTO 810
270 
271  WRITE (ndso,902) mk, mth
272  nspec = mk * mth
273  ALLOCATE ( freq(mk), freq1(mk), freq2(mk), dsip(mk), &
274  spec(mk,mth), e(mk,mth), thd(mth), dir(mth) )
275  !
276  READ (ndsinp,err=801,iostat=ierr) th1, dth, freq, dsip
277 
278  !
279  !--- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
280  ! 5. Time management.
281  !
282  iout = 0
283  ncid = 0
284  WRITE (ndso,970)
285  READ (ndsinp,END=444, ERR=801,IOSTAT=IERR) time
286  backspace(ndsinp)
287 
288 
289  ! 5.1 Loops on track_o.ww3 to read the time and data
290  DO
291  dtest = dsec21( time , tout )
292 
293  ! cycle to reach the start time of input file
294  IF ( dtest .LT. 0. ) THEN
295  CALL tick21 ( tout , dtreq )
296  cycle
297  END IF
298 
299  IF ( dtest .GE. 0. ) THEN
300  trckid=''
301  READ (ndsinp,END=444, ERR=801,IOSTAT=IERR) TIME, X, Y, TSTSTR, trckid
302  IF ( tststr .EQ. 'SEA' ) THEN
303  READ (ndsinp,err=801,iostat=ierr) dw, cx, cy, wx, wy, ust, &
304  as, spec
305  END IF
306  IF ( ierr .EQ. -1 ) THEN
307  WRITE (ndso,944)
308  EXIT
309  END IF
310 
311 
312  IF ( time(1).EQ.tout(1) .AND. time(2).EQ.tout(2) ) THEN
313  iloc = iloc + 1
314  IF ( tststr .EQ. 'SEA' ) ispec = ispec + 1
315  ENDIF
316  IF ( time(1).GT.tout(1) .OR. time(2).GT.tout(2) ) THEN
317  CALL stme21 ( time , stime )
318  WRITE (ndso,945) stime, iloc, ispec
319  iloc = 1
320  ispec = 0
321  IF ( tststr .EQ. 'SEA' ) ispec = ispec + 1
322  tout(1) = time(1)
323  tout(2) = time(2)
324  ENDIF
325  END IF
326 
327 
328  ! 5.1.1 Increments the global time counter IOUT
329  iout = iout + 1
330  CALL stme21 ( tout , idtime )
331  WRITE (ndso,971) idtime
332 
333 
334  ! 5.1.2 Processes the variable value for the time step IOUT
335  CALL w3exnc ( fileprefix, nctype, ncid, s3, strstopdate, mk, mth )
336 
337 
338  ! 5.1.3 Defines the stop date
339  CALL t2d(tout,stopdate,ierr)
340  WRITE(strstopdate,'(I4.4,A,4(I2.2,A),I2.2)') stopdate(1),'-',stopdate(2), &
341  '-',stopdate(3),' ',stopdate(5),':',stopdate(6),':',stopdate(7)
342 
343  IF ( iout .GE. nout ) EXIT
344  END DO
345 
346 
347 444 CONTINUE
348 
349  ! 5.2 Closes the netCDF file
350  IF (ncid.NE.0) THEN
351  iret = nf90_redef(ncid)
352  CALL check_err(iret)
353  iret=nf90_put_att(ncid,nf90_global,'stop_date',strstopdate)
354  CALL check_err(iret)
355  iret=nf90_close(ncid)
356  CALL check_err(iret)
357  END IF
358 
359  !
360  GOTO 888
361  !
362  ! Escape locations read errors :
363  !
364 800 CONTINUE
365  WRITE (ndse,1000) ierr
366  CALL extcde ( 10 )
367  !
368 801 CONTINUE
369  WRITE (ndse,1001)
370  CALL extcde ( 11 )
371  !
372 805 CONTINUE
373  WRITE (ndse,1004) ierr
374  CALL extcde ( 14 )
375  !
376 806 CONTINUE
377  WRITE (ndse,1005) ierr
378  CALL extcde ( 15 )
379  !
380 807 CONTINUE
381  WRITE (ndse,1006) ierr
382  CALL extcde ( 16 )
383  !
384 
385 810 CONTINUE
386  WRITE (ndse,1010) idstr, idtst
387  CALL extcde ( 20 )
388  !
389 888 CONTINUE
390  WRITE (ndso,999)
391  !
392  ! Formats
393  !
394 900 FORMAT (/15x,' *** WAVEWATCH III Track output postp. *** '/ &
395  15x,'==============================================='/)
396 901 FORMAT ( ' Comment character is ''',a,''''/)
397  !
398 902 FORMAT ( ' Spectral grid size : ',i3,' by ',i3// &
399  ' Opening file : '/ &
400  ' -----------------------------------------------')
401 920 FORMAT ( ' Grid name : ',a/)
402  !
403 940 FORMAT (/' Output time data : '/ &
404  ' --------------------------------------------------'/ &
405  ' First time : ',a)
406 941 FORMAT ( ' Interval : ',a/ &
407  ' Number of requests : ',i10)
408  !
409 944 FORMAT (/' End of file reached '/)
410  !
411 945 FORMAT ( ' ',a,' :',i6,' points and',i6,' spectra.')
412  !
413 970 FORMAT (//' Generating files '/ &
414  ' --------------------------------------------------')
415 971 FORMAT ( ' Files for ',a)
416  !
417 999 FORMAT (/' End of program '/ &
418  ' ========================================='/ &
419  ' WAVEWATCH III Track output '/)
420  !
421 1000 FORMAT (/' *** WAVEWATCH III ERROR IN W3TRNC : '/ &
422  ' ERROR IN OPENING INPUT FILE'/ &
423  ' IOSTAT =',i5/)
424  !
425 1001 FORMAT (/' *** WAVEWATCH III ERROR IN W3TRNC : '/ &
426  ' PREMATURE END OF INPUT FILE'/)
427  !
428 1004 FORMAT (/' *** WAVEWATCH III ERROR IN W3TRNC : '/ &
429  ' ERROR IN OPENING INPUT FILE'/ &
430  ' IOSTAT =',i5/)
431  !
432 1005 FORMAT (/' *** WAVEWATCH III ERROR IN W3TRNC : '/ &
433  ' ERROR IN READING FROM INPUT FILE'/ &
434  ' IOSTAT =',i5/)
435  !
436 1006 FORMAT (/' *** WAVEWATCH III ERROR IN W3TRNC : '/ &
437  ' ERROR IN OPENING OUTPUT FILE'/ &
438  ' IOSTAT =',i5/)
439  !
440 1010 FORMAT (/' *** WAVEWATCH III ERROR IN W3TRNC : '/ &
441  ' ILLEGAL TYPE, NCTYPE =',i4/)
442  !
443  !/
444  !/ Internal subroutine W3EXNC ---------------------------------------- /
445  !/
446 CONTAINS
447  !/ ------------------------------------------------------------------- /
448 
461  SUBROUTINE w3exnc ( FILEPREFIX, NCTYPE, NCID, S3, STRSTOPDATE, MK, MTH )
462  !/
463  !/ +-----------------------------------+
464  !/ | M. Accensi |
465  !/ | FORTRAN 90 |
466  !/ | Last update : 8-Apr-2016 |
467  !/ +-----------------------------------+
468  !/
469  !/ 8-apr-2016 : Creation ( version 5.11 )
470  !/
471  ! 1. Purpose :
472  !
473  ! Perform actual track output in NetCDF file.
474  !
475  ! 3. Parameters :
476  !
477  ! Parameter list
478  ! ----------------------------------------------------------------
479  ! ----------------------------------------------------------------
480  !
481  ! Internal parameters
482  ! ----------------------------------------------------------------
483  ! ----------------------------------------------------------------
484  !
485  ! 4. Subroutines used :
486  !
487  ! Name Type Module Description
488  ! ----------------------------------------------------------------
489  ! STRACE Subr. W3SERVMD Subroutine tracing.
490  ! EXTCDE Subr. Id. Abort program as graceful as possible.
491  ! ----------------------------------------------------------------
492  !
493  ! 5. Called by :
494  !
495  ! Main program in which it is contained.
496  !
497  ! 6. Error messages :
498  !
499  ! None.
500  !
501  ! 7. Remarks :
502  !
503  ! None.
504  !
505  ! 8. Structure :
506  !
507  ! See source code.
508  !
509  ! 9. Switches :
510  !
511  ! !/S Enable subroutine tracing.
512  ! !/T Enable test output.
513  !
514  ! 10. Source code :
515  !
516  !/ ------------------------------------------------------------------- /
517  USE netcdf
518  IMPLICIT NONE
519 
520  !/
521  !/ ------------------------------------------------------------------- /
522  !/ Parameter list
523  !/
524  INTEGER, INTENT(IN) :: NCTYPE, MK, MTH
525  CHARACTER(30), INTENT(IN) :: FILEPREFIX, STRSTOPDATE
526  INTEGER, INTENT(INOUT) :: NCID, S3
527  !/
528  !/ ------------------------------------------------------------------- /
529  !/ Local parameters
530  !/
531  INTEGER :: S1, S2, S4, S5, NDSDAT, IRET
532  INTEGER :: STARTDATE(8), CURDATE(8), REFDATE(8)
533  INTEGER :: DEFLATE=1
534 #ifdef W3_S
535  INTEGER, SAVE :: IENT = 0
536 #endif
537  !
538  DOUBLE PRECISION :: OUTJULDAY
539  !
540  CHARACTER*30 :: STRSTARTDATE
541  CHARACTER :: FNAMENC*50, ENAME*6
542  CHARACTER, SAVE :: OLDTIMEID*16 = '0000000000000000'
543  CHARACTER, SAVE :: TIMEID*16 = '0000000000000000'
544 
545  !/
546  !/ ------------------------------------------------------------------- /
547  !/
548  !
549 #ifdef W3_S
550  CALL strace (ient, 'W3EXNC')
551 #endif
552  !
553  CALL u2d('days since 1990-01-01 00:00:00',refdate,ierr)
554 
555  ! 1.1 Sets the date as ISO8601 convention
556  ! S3 defines the number of characters in the date for the filename
557  ! S3=4-> YYYY, S3=6 -> YYYYMM, S3=10 -> YYYYMMDDHH
558  ! Setups min and max date format
559  IF (s3.LT.4) s3=4
560  IF (s3.GT.10) s3=10
561  !
562  ! Defines the format of FILETIME
563  s5=s3-8
564  s4=s3
565  oldtimeid=timeid
566  ! if S3=>YYYYMMDDHH then filetime='YYYYMMDDTHHMMSSZ'
567  IF (s3.EQ.10) THEN
568  s4=s4+2 ! add chars for ISO8601 : day T hours Z
569  WRITE(format1,'(A,I1,A,I1,A)') '(I8.8,A1,I',s5,'.',s5,',A1)'
570  WRITE (timeid,format1) time(1), 'T', &
571  floor(real(time(2))/nint(10.**(6-s5))), 'Z'
572  ! if S3=>YYYYMMDD then filetime='YYYYMMDD'
573  ELSE IF (s3.EQ.8) THEN
574  WRITE(format1,'(A,I1,A,I1,A)') '(I',s3,'.',s3,')'
575  WRITE (timeid,format1) time(1)
576  ! if S3=>YYYYMM then filetime='YYYYMM'
577  ! or S3=>YYYY then filetime='YYYY'
578  ELSE
579  WRITE(format1,'(A,I1,A,I1,A)') '(I',s3,'.',s3,')'
580  WRITE (timeid,format1) floor(real(time(1))/nint(10.**(8-s3)))
581  END IF
582  ! redefines filename with updated date format
583  s1=len_trim(fileprefix)
584  fnamenc=''
585  fnamenc(1:s1)=fileprefix(1:s1)
586  fnamenc(s1+1:s1+s4) = timeid(1:s4)
587 
588 
589  ! 1.2 Setups the output type 4 ( NetCDF file )
590 
591  ename='.trck'
592  s2=len_trim(ename)
593  s1=len_trim(fileprefix)+s4
594  fnamenc(s1+1:50)=' '
595  fnamenc(s1+1:s1+1) = '_'
596 
597  ! add variable name in file name
598  fnamenc(s1+2:s1+s2) = ename(2:s2)
599 
600  ! Defines the netcdf extension
601  fnamenc(s1+s2+1:s1+s2+3) = '.nc'
602  fnamenc(s1+s2+4:s1+s2+6) = ' '
603 
604  ! Defines the dimensions
605  dimln(1)=nf90_unlimited ! time
606  dimln(2)=mk ! frequency
607  dimln(3)=mth ! direction
608  dimln(4)=32 ! string track name length
609 
610 
611  ! 1.3 Gets the netcdf id
612 
613  ndsdat=30
614  OPEN (ndsdat,file=fnamenc,status='new',iostat=iret)
615  IF (iret.EQ.0) THEN
616  ! CLOSE old file
617  IF (index('0000000000000000',oldtimeid).EQ.0 .AND. index(timeid,oldtimeid).EQ.0) THEN
618  iret = nf90_redef(ncid)
619  CALL check_err(iret)
620  iret=nf90_put_att(ncid,nf90_global,'stop_date',strstopdate)
621  CALL check_err(iret)
622  iret=nf90_close(ncid)
623  CALL check_err(iret)
624  END IF
625  ncid=0
626  ELSE
627  ncid=ncid
628  END IF
629 
630 
631  ! 1.4 Creates the netcdf file
632 
633  IF (ncid.EQ.0) THEN
634 
635  ! Initializes the time iteration counter n
636  it = 0
637  iloc = 0
638  ispec = 0
639 
640  ! 1.4.1 Creates the NetCDF file
641 
642  CALL w3crnc(nctype,fnamenc,ncid,dimid,dimln,varid)
643 
644  ! put start date in global attribute
645  CALL t2d(time,startdate,ierr)
646  WRITE(strstartdate,'(I4.4,A,4(I2.2,A),I2.2)') startdate(1),'-',startdate(2),'-', &
647  startdate(3),' ',startdate(5),':',startdate(6),':',startdate(7)
648  !
649  iret=nf90_put_att(ncid,nf90_global,'start_date',strstartdate)
650  CALL check_err(iret)
651 
652  ! End of define mode of NetCDF file
653  iret = nf90_enddef(ncid)
654  CALL check_err(iret)
655 
656  ! Process lower band and higher band frequencies
657  freq1(1:mk)=freq(1:mk)-0.5*(freq(1:mk)-(freq(1:mk)/xfr))
658  freq2(1:mk)=freq(1:mk)+0.5*(-freq(1:mk)+(freq(1:mk)*xfr))
659  freq1(1)=freq(1)
660  freq2(mk)=freq(mk)
661 
662  ! Converts direction unit in degree
663  dthd=360./mth
664  rth0=th1/dth
665  DO ith=1, mth
666  thd(ith)=dthd*(rth0+real(ith-1))
667  END DO
668  dir(1:mth)=mod(360-thd(1:mth),360.)
669 
670 
671  ! 1.4.2 Adds general variables to NetCDF file
672  iret=nf90_put_var(ncid,varid(2),freq)
673  CALL check_err(iret)
674 
675  iret=nf90_put_var(ncid,varid(3),freq1)
676  CALL check_err(iret)
677 
678  iret=nf90_put_var(ncid,varid(4),freq2)
679  CALL check_err(iret)
680 
681  iret=nf90_put_var(ncid,varid(5),dsip)
682  CALL check_err(iret)
683 
684  iret=nf90_put_var(ncid,varid(6),dir)
685  CALL check_err(iret)
686 
687  WRITE (ndso,973) fnamenc
688 
689  END IF ! IERR.EQ.0
690 
691 
692  ! 1.5 Defines the current time step and index
693 
694  CALL t2d(time,curdate,ierr)
695  outjulday=tsub(refdate,curdate)
696  WRITE(ndso,'(3A,I6,A,I4,A,I2.2,A,I2.2,A,I2.2,A,I2.2,A,I2.2,2A)') &
697  'Writing new record ', ename(2:) ,'number ',it, &
698  ' for ',curdate(1),':',curdate(2),':',curdate(3),'T',curdate(5), &
699  ':',curdate(6),':',curdate(7),' in file ',trim(fnamenc)
700 
701 
702  !
703  ! 1.6 Exit from W3EXNC if not sea point
704  !
705  IF ( tststr .NE. 'SEA' ) GOTO 888
706 
707 
708  !
709  ! 1.6.1 Process speed and direction components
710  !
711  wao = sqrt( wx**2 + wy**2 )
712  IF ( wao.GT.1.e-7 ) THEN
713  wdo = mod(270.-atan2(wy,wx)*rade,360.)
714  ELSE
715  wdo = 0.
716  END IF
717 
718  cao = sqrt( cx**2 + cy**2 )
719  IF ( cao.GT.1.e-7 ) THEN
720  cdo = mod(270.-atan2(cy,cx)*rade,360.)
721  ELSE
722  cdo = 0.
723  END IF
724 
725  !
726  ! 1.7.1 Puts dimensions variables in NetCDF file
727  !
728  it=it+1
729  IF ( ust .LT. 0. ) ust = -1.0
730 
731  ! time
732  iret=nf90_put_var(ncid,varid(1),outjulday,start=(/it/))
733  CALL check_err(iret)
734  ! longitude
735  iret=nf90_put_var(ncid,varid(7),m2km*x,start=(/it/))
736  CALL check_err(iret)
737  ! latitude
738  iret=nf90_put_var(ncid,varid(8),m2km*y,start=(/it/))
739  CALL check_err(iret)
740 
741 
742 
743  ! 1.7.2 Puts fields in NetCDF file
744 
745 
746  ! 1.7.2.a Write spectrum
747 
748  iret=nf90_put_var(ncid,varid(9), &
749  transpose(spec),start=(/1,1,it/), count=(/mth,mk,1/))
750  CALL check_err(iret)
751 
752  ! 1.7.2.b Write the basic stuff
753 
754  ! Write DW (depth)
755  iret=nf90_put_var(ncid, varid(10),dw ,start=(/it/))
756  CALL check_err(iret)
757  ! Write CAO (current - x direction)
758  iret=nf90_put_var(ncid, varid(11),cao ,start=(/it/))
759  CALL check_err(iret)
760  ! Write CDO (current - y direction)
761  iret=nf90_put_var(ncid,varid(12),cdo ,start=(/it/))
762  CALL check_err(iret)
763  ! Write WAO (wind velocity - x direction)
764  iret=nf90_put_var(ncid,varid(13),wao ,start=(/it/))
765  CALL check_err(iret)
766  ! Write WDO (wind velocity - y direction)
767  iret=nf90_put_var(ncid,varid(14),wdo ,start=(/it/))
768  CALL check_err(iret)
769  ! Write UST (friction velocity)
770  iret=nf90_put_var(ncid,varid(15),ust,start=(/it/))
771  CALL check_err(iret)
772  ! Write AS (air sea temperature difference)
773  iret=nf90_put_var(ncid,varid(16),as ,start=(/it/))
774  CALL check_err(iret)
775  ! Track name
776  iret=nf90_put_var(ncid,varid(18),trckid,start=(/1,it/),count=(/len_trim(trckid),1/))
777  CALL check_err(iret)
778 
779 
780  !
781 888 CONTINUE
782  !
783  RETURN
784 
785  !
786  ! Formats
787  !
788 973 FORMAT ( 'NEW NetCDF file was created ',a)
789 
790 
791  !/ End of W3EXNC ----------------------------------------------------- /
792  !/
793  END SUBROUTINE w3exnc
794 
795 
796 
797 
798  !--------------------------------------------------------------------------
810  SUBROUTINE w3crnc (NCTYPE,NCFILE,NCID,DIMID,DIMLN,VARID)
811 
812  USE netcdf
813 
814  IMPLICIT NONE
815 
816  INTEGER, INTENT(IN) :: NCTYPE
817  CHARACTER*(*), INTENT(IN) :: NCFILE
818  INTEGER, INTENT(IN) :: DIMLN(:)
819  INTEGER, INTENT(OUT) :: DIMID(:), VARID(:), NCID
820  INTEGER :: IRET
821  INTEGER :: DEFLATE=1
822 
823  !
824  ! Creation in netCDF3 or netCDF4
825  !
826  IF(nctype.EQ.3) iret = nf90_create(trim(ncfile), nf90_clobber, ncid)
827  IF(nctype.EQ.4) iret = nf90_create(trim(ncfile), nf90_netcdf4, ncid)
828  CALL check_err(iret)
829 
830  !
831  ! Define generals dimensions
832  !
833  iret = nf90_def_dim(ncid, 'time', dimln(1), dimid(1))
834  CALL check_err(iret)
835  iret = nf90_def_dim(ncid, 'frequency', dimln(2), dimid(2))
836  CALL check_err(iret)
837  iret = nf90_def_dim(ncid, 'direction', dimln(3), dimid(3))
838  CALL check_err(iret)
839  iret = nf90_def_dim(ncid, 'string32', dimln(4), dimid(4))
840  CALL check_err(iret)
841 
842  !
843  ! define generals variables
844  !
845 
846  ! time
847  iret=nf90_def_var(ncid, 'time', nf90_double, (/dimid(1)/), varid(1))
848  IF (nctype.EQ.4) iret=nf90_def_var_deflate(ncid, varid(1), 1, 1, deflate)
849  SELECT CASE (trim(caltype))
850  CASE ('360_day')
851  iret=nf90_put_att(ncid,varid(1),'long_name','time in 360 day calendar')
852  CASE ('365_day')
853  iret=nf90_put_att(ncid,varid(1),'long_name','time in 365 day calendar')
854  CASE ('standard')
855  iret=nf90_put_att(ncid,varid(1),'long_name','julian day (UT)')
856  END SELECT
857  iret=nf90_put_att(ncid,varid(1),'standard_name','time')
858  iret=nf90_put_att(ncid,varid(1),'units','days since 1990-01-01 00:00:00')
859  iret=nf90_put_att(ncid,varid(1),'conventions', &
860  'Relative julian days with decimal part (as parts of the day)')
861  iret=nf90_put_att(ncid,varid(1),'axis','T')
862  iret=nf90_put_att(ncid,varid(1),'calendar',trim(caltype))
863 
864  ! frequency
865  iret=nf90_def_var(ncid, 'frequency', nf90_float, (/dimid(2)/),varid(2))
866  IF (nctype.EQ.4) iret=nf90_def_var_deflate(ncid, varid(2), 1, 1, deflate)
867  iret=nf90_put_att(ncid,varid(2),'long_name','center frequencies for spectra')
868  iret=nf90_put_att(ncid,varid(2),'standard_name','frequency')
869  iret=nf90_put_att(ncid,varid(2),'units','s-1')
870  iret=nf90_put_att(ncid,varid(2),'scale_factor',1.)
871  iret=nf90_put_att(ncid,varid(2),'add_offset',0.)
872  iret=nf90_put_att(ncid,varid(2),'valid_min',0.)
873  iret=nf90_put_att(ncid,varid(2),'valid_max',10.)
874  iret=nf90_put_att(ncid,varid(2),'_FillValue',nf90_fill_float)
875  iret=nf90_put_att(ncid,varid(2),'axis','Y')
876 
877  !frequency1
878  iret=nf90_def_var(ncid, 'frequency1', nf90_float, (/dimid(2)/), varid(3))
879  CALL check_err(iret)
880  IF (nctype.EQ.4) iret=nf90_def_var_deflate(ncid, varid(3), 1, 1, deflate)
881  iret=nf90_put_att(ncid,varid(3),'long_name','frequency of lower band')
882  iret=nf90_put_att(ncid,varid(3),'standard_name','frequency_of_lower_band')
883  iret=nf90_put_att(ncid,varid(3),'globwave_name','frequency_lower_band')
884  iret=nf90_put_att(ncid,varid(3),'units','s-1')
885  iret=nf90_put_att(ncid,varid(3),'scale_factor',1.)
886  iret=nf90_put_att(ncid,varid(3),'add_offset',0.)
887  iret=nf90_put_att(ncid,varid(3),'valid_min',0.)
888  iret=nf90_put_att(ncid,varid(3),'valid_max',10.)
889  iret=nf90_put_att(ncid,varid(3),'_FillValue',nf90_fill_float)
890  iret=nf90_put_att(ncid,varid(3),'content','Y')
891  iret=nf90_put_att(ncid,varid(3),'associates','frequency')
892 
893  !frequency2
894  iret=nf90_def_var(ncid, 'frequency2', nf90_float, (/dimid(2)/), varid(4))
895  CALL check_err(iret)
896  IF (nctype.EQ.4) iret=nf90_def_var_deflate(ncid, varid(4), 1, 1, deflate)
897  iret=nf90_put_att(ncid,varid(4),'long_name','frequency of upper band')
898  iret=nf90_put_att(ncid,varid(4),'standard_name','frequency_of_upper_band')
899  iret=nf90_put_att(ncid,varid(4),'globwave_name','frequency_upper_band')
900  iret=nf90_put_att(ncid,varid(4),'units','s-1')
901  iret=nf90_put_att(ncid,varid(4),'scale_factor',1.)
902  iret=nf90_put_att(ncid,varid(4),'add_offset',0.)
903  iret=nf90_put_att(ncid,varid(4),'valid_min',0.)
904  iret=nf90_put_att(ncid,varid(4),'valid_max',10.)
905  iret=nf90_put_att(ncid,varid(4),'_FillValue',nf90_fill_float)
906  iret=nf90_put_att(ncid,varid(4),'content','Y')
907  iret=nf90_put_att(ncid,varid(4),'associates','frequency')
908 
909  ! frequency area
910  iret=nf90_def_var(ncid, 'frequency_area', nf90_float,(/dimid(2)/),varid(5))
911  IF (nctype.EQ.4) iret=nf90_def_var_deflate(ncid, varid(5), 1, 1, deflate)
912  iret=nf90_put_att(ncid,varid(5),'long_name','frequency spectral bin width')
913  iret=nf90_put_att(ncid,varid(5),'standard_name','frequency_area')
914  iret=nf90_put_att(ncid,varid(5),'units','s-2')
915  iret=nf90_put_att(ncid,varid(5),'scale_factor',1.)
916  iret=nf90_put_att(ncid,varid(5),'add_offset',0.)
917  iret=nf90_put_att(ncid,varid(5),'valid_min',0.)
918  iret=nf90_put_att(ncid,varid(5),'valid_max',10.)
919  iret=nf90_put_att(ncid,varid(5),'_FillValue',nf90_fill_float)
920  iret=nf90_put_att(ncid,varid(5),'content','Y')
921  iret=nf90_put_att(ncid,varid(5),'associates','frequency')
922 
923  ! direction
924  iret=nf90_def_var(ncid, 'direction', nf90_float, (/dimid(3)/),varid(6))
925  IF (nctype.EQ.4) iret=nf90_def_var_deflate(ncid, varid(6), 1, 1, deflate)
926  iret=nf90_put_att(ncid,varid(6),'long_name','sea surface wave to direction')
927  iret=nf90_put_att(ncid,varid(6),'standard_name','sea_surface_wave_to_direction')
928  iret=nf90_put_att(ncid,varid(6),'units','degree')
929  iret=nf90_put_att(ncid,varid(6),'scale_factor',1.)
930  iret=nf90_put_att(ncid,varid(6),'add_offset',0.)
931  iret=nf90_put_att(ncid,varid(6),'valid_min',0.)
932  iret=nf90_put_att(ncid,varid(6),'valid_max',360.)
933  iret=nf90_put_att(ncid,varid(6),'_FillValue',nf90_fill_float)
934  iret=nf90_put_att(ncid,varid(6),'axis','Z')
935 
936  IF (flagll) THEN
937  ! longitude
938  iret=nf90_def_var(ncid, 'longitude', nf90_float, (/dimid(1)/),varid(7))
939  IF (nctype.EQ.4) iret=nf90_def_var_deflate(ncid, varid(7), 1, 1, deflate)
940  iret=nf90_put_att(ncid,varid(7),'long_name','longitude')
941  iret=nf90_put_att(ncid,varid(7),'standard_name','longitude')
942  iret=nf90_put_att(ncid,varid(7),'units','degree_east')
943  iret=nf90_put_att(ncid,varid(7),'valid_min',-180.0)
944  iret=nf90_put_att(ncid,varid(7),'valid_max',360.)
945  iret=nf90_put_att(ncid,varid(7),'_FillValue',nf90_fill_float)
946  iret=nf90_put_att(ncid,varid(7),'content','T')
947  iret=nf90_put_att(ncid,varid(7),'associates','time')
948 
949 
950  ! latitude
951  iret=nf90_def_var(ncid, 'latitude', nf90_float, (/dimid(1)/),varid(8))
952  IF (nctype.EQ.4) iret=nf90_def_var_deflate(ncid, varid(8), 1, 1, deflate)
953  iret=nf90_put_att(ncid,varid(8),'long_name','latitude')
954  iret=nf90_put_att(ncid,varid(8),'standard_name','latitude')
955  iret=nf90_put_att(ncid,varid(8),'units','degree_north')
956  iret=nf90_put_att(ncid,varid(8),'valid_min',-90.0)
957  iret=nf90_put_att(ncid,varid(8),'valid_max',180.)
958  iret=nf90_put_att(ncid,varid(8),'_FillValue',nf90_fill_float)
959  iret=nf90_put_att(ncid,varid(8),'content','T')
960  iret=nf90_put_att(ncid,varid(8),'associates','time')
961  ELSE
962  ! longitude
963  iret=nf90_def_var(ncid, 'x', nf90_float, (/dimid(1)/),varid(7))
964  IF (nctype.EQ.4) iret=nf90_def_var_deflate(ncid, varid(7), 1, 1, deflate)
965  iret=nf90_put_att(ncid,varid(7),'long_name','x')
966  iret=nf90_put_att(ncid,varid(7),'standard_name','x')
967  iret=nf90_put_att(ncid,varid(7),'units','m')
968  iret=nf90_put_att(ncid,varid(7),'_FillValue',nf90_fill_float)
969  iret=nf90_put_att(ncid,varid(7),'content','T')
970  iret=nf90_put_att(ncid,varid(7),'associates','time')
971 
972  ! latitude
973  iret=nf90_def_var(ncid, 'y', nf90_float, (/dimid(1)/),varid(8))
974  IF (nctype.EQ.4) iret=nf90_def_var_deflate(ncid, varid(8), 1, 1, deflate)
975  iret=nf90_put_att(ncid,varid(8),'long_name','y')
976  iret=nf90_put_att(ncid,varid(8),'standard_name','y')
977  iret=nf90_put_att(ncid,varid(8),'units','m')
978  iret=nf90_put_att(ncid,varid(8),'_FillValue',nf90_fill_float)
979  iret=nf90_put_att(ncid,varid(8),'content','T')
980  iret=nf90_put_att(ncid,varid(8),'associates','time')
981 
982  END IF
983 
984 
985  ! Efth
986  iret=nf90_def_var(ncid,'efth',nf90_float,(/dimid(3),dimid(2),dimid(1)/),varid(9))
987  IF (nctype.EQ.4) iret=nf90_def_var_deflate(ncid, varid(9), 1, 1, deflate)
988  iret=nf90_put_att(ncid,varid(9),'long_name', &
989  'sea surface wave directional variance spectral density')
990  iret=nf90_put_att(ncid,varid(9),'standard_name', &
991  'sea_surface_wave_directional_variance_spectral_density')
992  iret=nf90_put_att(ncid,varid(9),'globwave_name', &
993  'directional_variance_spectral_density')
994  iret=nf90_put_att(ncid,varid(9),'units','m2 s rad-1')
995  iret=nf90_put_att(ncid,varid(9),'scale_factor',1.)
996  iret=nf90_put_att(ncid,varid(9),'add_offset',0.)
997  iret=nf90_put_att(ncid,varid(9),'valid_min',0.)
998  iret=nf90_put_att(ncid,varid(9),'valid_max',10.)
999  iret=nf90_put_att(ncid,varid(9),'_FillValue',nf90_fill_float)
1000  iret=nf90_put_att(ncid,varid(9),'content','TYZ')
1001  iret=nf90_put_att(ncid,varid(9),'associates','time frequency direction')
1002 
1003  ! DW - depth
1004  iret=nf90_def_var(ncid, 'dpt', nf90_float, (/dimid(1)/),varid(10))
1005  IF (nctype.EQ.4) iret=nf90_def_var_deflate(ncid, varid(10), 1, 1, deflate)
1006  iret=nf90_put_att(ncid,varid(10),'long_name','depth')
1007  iret=nf90_put_att(ncid,varid(10),'standard_name','depth')
1008  iret=nf90_put_att(ncid,varid(10),'globwave_name','depth')
1009  iret=nf90_put_att(ncid,varid(10),'units','m')
1010  iret=nf90_put_att(ncid,varid(10),'scale_factor',1.)
1011  iret=nf90_put_att(ncid,varid(10),'add_offset',0.)
1012  iret=nf90_put_att(ncid,varid(10),'_FillValue',nf90_fill_float)
1013  iret=nf90_put_att(ncid,varid(10),'content','T')
1014  iret=nf90_put_att(ncid,varid(10),'associates','time')
1015 
1016  ! CAO - current speed (m/s)
1017  iret=nf90_def_var(ncid, 'cur', nf90_float,(/dimid(1)/), varid(11))
1018  IF (nctype.EQ.4) iret=nf90_def_var_deflate(ncid, varid(11), 1, 1, deflate)
1019  iret=nf90_put_att(ncid,varid(11),'long_name','sea water speed')
1020  iret=nf90_put_att(ncid,varid(11),'standard_name','sea_water_speed')
1021  iret=nf90_put_att(ncid,varid(11),'globwave_name','sea_water_speed')
1022  iret=nf90_put_att(ncid,varid(11),'units','m s-1')
1023  iret=nf90_put_att(ncid,varid(11),'scale_factor',1.)
1024  iret=nf90_put_att(ncid,varid(11),'add_offset',0.)
1025  iret=nf90_put_att(ncid,varid(11),'_FillValue',nf90_fill_float)
1026  iret=nf90_put_att(ncid,varid(11),'content','T')
1027  iret=nf90_put_att(ncid,varid(11),'associates','time')
1028 
1029  ! CDO - current direction (degree)
1030  iret=nf90_def_var(ncid, 'curdir', nf90_float,(/dimid(1)/), varid(12))
1031  IF (nctype.EQ.4) iret=nf90_def_var_deflate(ncid, varid(12), 1, 1, deflate)
1032  iret=nf90_put_att(ncid,varid(12),'long_name','direction from of sea water velocity')
1033  iret=nf90_put_att(ncid,varid(12),'standard_name','direction_of_sea_water_velocity')
1034  iret=nf90_put_att(ncid,varid(12),'globwave_name','direction_of_sea_water_velocity')
1035  iret=nf90_put_att(ncid,varid(12),'units','degree')
1036  iret=nf90_put_att(ncid,varid(12),'scale_factor',1.)
1037  iret=nf90_put_att(ncid,varid(12),'add_offset',0.)
1038  iret=nf90_put_att(ncid,varid(12),'_FillValue',nf90_fill_float)
1039  iret=nf90_put_att(ncid,varid(12),'content','T')
1040  iret=nf90_put_att(ncid,varid(12),'associates','time')
1041 
1042  ! WAO - wind speed (m/s)
1043  iret=nf90_def_var(ncid, 'wnd', nf90_float,(/dimid(1)/), varid(13))
1044  IF (nctype.EQ.4) iret=nf90_def_var_deflate(ncid, varid(13), 1, 1, deflate)
1045  iret=nf90_put_att(ncid,varid(13),'long_name','wind speed at 10m')
1046  iret=nf90_put_att(ncid,varid(13),'standard_name','wind_speed')
1047  iret=nf90_put_att(ncid,varid(13),'globwave_name','wind_speed')
1048  iret=nf90_put_att(ncid,varid(13),'units','m s-1')
1049  iret=nf90_put_att(ncid,varid(13),'scale_factor',1.)
1050  iret=nf90_put_att(ncid,varid(13),'add_offset',0.)
1051  iret=nf90_put_att(ncid,varid(13),'_FillValue',nf90_fill_float)
1052  iret=nf90_put_att(ncid,varid(13),'content','T')
1053  iret=nf90_put_att(ncid,varid(13),'associates','time')
1054 
1055  ! WDO - wind direction (degree)
1056  iret=nf90_def_var(ncid, 'wnddir', nf90_float,(/dimid(1)/), varid(14))
1057  IF (nctype.EQ.4) iret=nf90_def_var_deflate(ncid, varid(14), 1, 1, deflate)
1058  iret=nf90_put_att(ncid,varid(14),'long_name','wind direction')
1059  iret=nf90_put_att(ncid,varid(14),'standard_name','wind_from_direction')
1060  iret=nf90_put_att(ncid,varid(14),'globwave_name','wind_from_direction')
1061  iret=nf90_put_att(ncid,varid(14),'units','m s-1')
1062  iret=nf90_put_att(ncid,varid(14),'scale_factor',1.)
1063  iret=nf90_put_att(ncid,varid(14),'add_offset',0.)
1064  iret=nf90_put_att(ncid,varid(14),'_FillValue',nf90_fill_float)
1065  iret=nf90_put_att(ncid,varid(14),'content','T')
1066  iret=nf90_put_att(ncid,varid(14),'associates','time')
1067 
1068  ! UST - friction velocity (m/s)
1069  iret=nf90_def_var(ncid, 'ust', nf90_float,(/dimid(1)/), varid(15))
1070  IF (nctype.EQ.4) iret=nf90_def_var_deflate(ncid, varid(15), 1, 1, deflate)
1071  iret=nf90_put_att(ncid,varid(15),'long_name','friction velocity')
1072  iret=nf90_put_att(ncid,varid(15),'standard_name','friction_velocity')
1073  iret=nf90_put_att(ncid,varid(15),'globwave_name','friction_velocity')
1074  iret=nf90_put_att(ncid,varid(15),'units','m s-1')
1075  iret=nf90_put_att(ncid,varid(15),'scale_factor',1.)
1076  iret=nf90_put_att(ncid,varid(15),'add_offset',0.)
1077  iret=nf90_put_att(ncid,varid(15),'_FillValue',nf90_fill_float)
1078  iret=nf90_put_att(ncid,varid(15),'content','T')
1079  iret=nf90_put_att(ncid,varid(15),'associates','time')
1080 
1081  ! AS - air-sea temperature difference (deg C)
1082  iret=nf90_def_var(ncid, 'ast',nf90_float,(/dimid(1)/), varid(16))
1083  IF (nctype.EQ.4) iret=nf90_def_var_deflate(ncid, varid(16), 1, 1, deflate)
1084  iret=nf90_put_att(ncid,varid(16),'long_name','air sea temperature difference')
1085  iret=nf90_put_att(ncid,varid(16),'standard_name','air_sea_temperature_difference')
1086  iret=nf90_put_att(ncid,varid(16),'globwave_name','air_sea_temperature_difference')
1087  iret=nf90_put_att(ncid,varid(16),'units','degree')
1088  iret=nf90_put_att(ncid,varid(16),'scale_factor',1.)
1089  iret=nf90_put_att(ncid,varid(16),'add_offset',0.)
1090  iret=nf90_put_att(ncid,varid(16),'_FillValue',nf90_fill_float)
1091  iret=nf90_put_att(ncid,varid(16),'content','T')
1092  iret=nf90_put_att(ncid,varid(16),'associates','time')
1093 
1094  ! string32
1095  iret=nf90_def_var(ncid, 'string32', nf90_int, (/dimid(4)/), varid(17))
1096  CALL check_err(iret)
1097  IF (nctype.EQ.4) iret=nf90_def_var_deflate(ncid, varid(17), 1, 1, deflate)
1098  iret=nf90_put_att(ncid,varid(17),'long_name','track_name number of characters')
1099  iret=nf90_put_att(ncid,varid(17),'_FillValue',nf90_fill_int)
1100  iret=nf90_put_att(ncid,varid(17),'axis','W')
1101 
1102  ! track_name
1103  iret=nf90_def_var(ncid, 'track_name', nf90_char, (/dimid(4),dimid(1)/), varid(18))
1104  CALL check_err(iret)
1105  IF (nctype.EQ.4) iret=nf90_def_var_deflate(ncid, varid(18), 1, 1, deflate)
1106  iret=nf90_put_att(ncid,varid(18),'long_name','track name')
1107  iret=nf90_put_att(ncid,varid(18),'content','TX')
1108  iret=nf90_put_att(ncid,varid(18),'associates','time string16')
1109 
1110  RETURN
1111 
1112  END SUBROUTINE w3crnc
1113 
1114  !==============================================================================
1121  SUBROUTINE check_err(IRET)
1123  USE netcdf
1124  USE w3odatmd, ONLY: ndse
1125  USE w3servmd, ONLY: extcde
1126 
1127  IMPLICIT NONE
1128 
1129  INTEGER IRET
1130 
1131  IF (iret .NE. nf90_noerr) THEN
1132  WRITE(ndse,*) ' *** WAVEWATCH III ERROR IN TRNC :'
1133  WRITE(ndse,*) ' NETCDF ERROR MESSAGE: '
1134  WRITE(ndse,*) nf90_strerror(iret)
1135  CALL extcde ( 59 )
1136  END IF
1137  RETURN
1138 
1139  END SUBROUTINE check_err
1140 
1141  !==============================================================================
1142 
1143 
1144  !/
1145  !/ End of W3TRNC ----------------------------------------------------- /
1146  !/
1147 END PROGRAM w3trnc
w3servmd::nextln
subroutine nextln(CHCKC, NDSI, NDSE)
Definition: w3servmd.F90:222
w3adatmd
Define data structures to set up wave model auxiliary data for several models simultaneously.
Definition: w3adatmd.F90:26
w3nmltrncmd::nml_file_t
Definition: w3nmltrncmd.F90:35
constants::rade
real, parameter rade
RADE Conversion factor from radians to degrees.
Definition: constants.F90:76
w3nmltrncmd
Definition: w3nmltrncmd.F90:3
w3gdatmd::gname
character(len=30), pointer gname
Definition: w3gdatmd.F90:1223
w3odatmd::fnmpre
character(len=80) fnmpre
Definition: w3odatmd.F90:330
w3gdatmd::w3setg
subroutine w3setg(IMOD, NDSE, NDST)
Definition: w3gdatmd.F90:2152
w3odatmd::ndse
integer, pointer ndse
Definition: w3odatmd.F90:456
w3adatmd::w3seta
subroutine w3seta(IMOD, NDSE, NDST)
Select one of the WAVEWATCH III grids / models.
Definition: w3adatmd.F90:2645
w3exnc
subroutine w3exnc(NX, NY, IX1, IXN, IY1, IYN, NSEA, FILEPREFIX, E3DF, P2MSF, US3DF, USSPF, NCTYPE, TOGETHER, NCVARTYPEI, FLG2D, NCIDS, S3, STRSTOPDATE)
Perform actual grid output in NetCDF file.
Definition: ww3_ounf.F90:863
w3servmd
Definition: w3servmd.F90:3
w3odatmd::w3seto
subroutine w3seto(IMOD, NDSERR, NDSTST)
Definition: w3odatmd.F90:1523
w3odatmd
Definition: w3odatmd.F90:3
w3adatmd::w3naux
subroutine w3naux(NDSE, NDST)
Set up the number of grids to be used.
Definition: w3adatmd.F90:704
w3iogrmd::w3iogr
subroutine w3iogr(INXOUT, NDSM, IMOD, FEXT ifdef W3_ASCII
Reading and writing of the model definition file.
Definition: w3iogrmd.F90:117
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
w3crnc
subroutine w3crnc(NCFILE, NCID, DIMID, DIMLN, VARID, EXTRADIM, NCTYPE, MAPSTAOUT)
Desc not available.
Definition: ww3_ounf.F90:3290
check_err
subroutine check_err(IRET)
Check input return status for error value.
Definition: ww3_bounc.F90:856
w3servmd::strace
subroutine strace(IENT, SNAME)
Definition: w3servmd.F90:148
w3odatmd::ndso
integer, pointer ndso
Definition: w3odatmd.F90:456
w3gdatmd::w3nmod
subroutine w3nmod(NUMBER, NDSE, NDST, NAUX)
Definition: w3gdatmd.F90:1433
w3gdatmd::xfr
real, pointer xfr
Definition: w3gdatmd.F90:1232
constants
Define some much-used constants for global use (all defined as PARAMETER).
Definition: constants.F90:20
w3gdatmd
Definition: w3gdatmd.F90:16
constants::file_endian
character(*), parameter file_endian
FILE_ENDIAN Filled by preprocessor with 'big_endian', 'little_endian', or 'native'.
Definition: constants.F90:86
w3servmd::extcde
subroutine extcde(IEXIT, UNIT, MSG, FILE, LINE, COMM)
Definition: w3servmd.F90:736
w3odatmd::w3nout
subroutine w3nout(NDSERR, NDSTST)
Definition: w3odatmd.F90:561
w3servmd::itrace
subroutine itrace(NDS, NMAX)
Definition: w3servmd.F90:91
w3trnc
program w3trnc
Convert direct access track output file to netCDF file.
Definition: ww3_trnc.F90:15
w3timemd
Definition: w3timemd.F90:3
w3nmltrncmd::w3nmltrnc
subroutine w3nmltrnc(NDSI, INFILE, NML_TRACK, NML_FILE, IERR)
Definition: w3nmltrncmd.F90:50
w3nmltrncmd::nml_track_t
Definition: w3nmltrncmd.F90:27
w3gdatmd::flagll
logical, pointer flagll
Definition: w3gdatmd.F90:1219