WAVEWATCH III  beta 0.0.1
w3iobcmd Module Reference

Processing of boundary data output. More...

Functions/Subroutines

subroutine w3iobc (INXOUT, NDSB, TIME1, TIME2, IOTST, IMOD)
 Write/read boundary conditions file(s). More...
 

Variables

character(len=10), parameter verbptbc = '2018-03-01'
 
character(len=32), parameter idstrbc = 'WAVEWATCH III BOUNDARY DATA FILE'
 

Detailed Description

Processing of boundary data output.

Author
H. L. Tolman
Date
01-Mar-2018

Function/Subroutine Documentation

◆ w3iobc()

subroutine w3iobcmd::w3iobc ( character, dimension(*), intent(in)  INXOUT,
integer, intent(in)  NDSB,
integer, dimension(2), intent(inout)  TIME1,
integer, dimension(2), intent(out)  TIME2,
integer, intent(out)  IOTST,
integer, intent(in), optional  IMOD 
)

Write/read boundary conditions file(s).

The file(s) are opened within the routine, the names are pre-defined as nest.FILEXT for the input file and nest1.FILEXT through nest9.FILEXT for up to 9 output files.

Parameters
[in,out]INXOUTTest string for read/write.
[in,out]NDSBData set unit number.
[in,out]TIME1Present time (w), time of first field (r).
[in,out]TIME2Time of second field.
[in,out]IOTSTTest indictor for reading.
[in,out]IMODOptional grid number, defaults to 1.
Author
H. L. Tolman
Date
20-Jan-2017

Definition at line 99 of file w3iobcmd.F90.

99  !/
100  !/ +-----------------------------------+
101  !/ | WAVEWATCH III NOAA/NCEP |
102  !/ | H. L. Tolman |
103  !/ | FORTRAN 90 |
104  !/ | Last update : 20-Jan-2017 |
105  !/ +-----------------------------------+
106  !/
107  !/ 12-Jan-1999 : Distributed FORTRAN 77 version. ( version 1.18 )
108  !/ 20-May-1999 : Remove read bug for IPBP and RDBP ( see web page )
109  !/ 30-Dec-1999 : Upgrade to FORTRAN 90 ( version 2.00 )
110  !/ Major changes to logistics.
111  !/ 13-Dec-2004 : Multiple grid version. ( version 3.06 )
112  !/ 19-Sep-2005 : Allow for change of spec. res. ( version 3.08 )
113  !/ (on read only).
114  !/ 30-Sep-2005 : Add 'DUMP' option. ( version 3.08 )
115  !/ 27-Jun-2006 : Adding file name preamble. ( version 3.09 )
116  !/ 29-May-2009 : Preparing distribution version. ( version 3.14 )
117  !/ 30-Oct-2009 : Implement curvilinear grid type. ( version 3.14 )
118  !/ (W. E. Rogers & T. J. Campbell, NRL)
119  !/ 28-Jul-2010 : Moving NKI, NTHI, XFRI, FR1I and
120  !/ TH1I to W3ODATMD. ( version 3.14.3 )
121  !/ 31-Oct-2010 : Implementing unstructured grid ( version 3.14.3 )
122  !/ (A. Roland and F. Ardhuin)
123  !/ 05-Apr-2011 : Moved the W3CSPC call into loop ( version 3.14.3 )
124  !/ 12-Jun-2012 : Add /RTD option or rotated grid option.
125  !/ (Jian-Guo Li) ( version 4.06 )
126  !/ 03-Jul-2013 : Corrected ABPIN indices ( version 4.11 )
127  !/ 14-Jan-2014 : Corrected ABPIN indices for W3CSPC ( version 4.18 )
128  !/ 20-Jan-2017 : Allow input boundary points to lie outside the grid
129  !/ within a distance of 0.1 times the grid cell size.
130  !/ (T.J. Campbell, NRL) ( version 6.02 )
131  !/ 01-Mar-2018 : Rotate boundary points and directions
132  !/ of input spectra for rotated grids ( version 6.02 )
133  !/ 07-Oct-2019 : RTD option with standard lat-lon
134  !/ grid when nesting to rotated grid ( version 7.11 )
135  !/
136  ! 1. Purpose :
137  !
138  ! Write/read boundary conditions file(s).
139  !
140  ! 2. Method :
141  !
142  ! The file(s) are opened within the routine, the names are
143  ! pre-defined as nest.FILEXT for the input file and nest1.FILEXT
144  ! through nest9.FILEXT for up to 9 output files.
145  !
146  ! 3. Parameters :
147  !
148  ! Parameter list
149  ! ----------------------------------------------------------------
150  ! INXOUT C*(*) I Test string for read/write, valid are:
151  ! 'READ', 'WRITE' or 'DUMP'.
152  ! NDSB Int. I Data set unit number.
153  ! TIME1 I.A. I/O Present time. (w)
154  ! Time of first field. (r)
155  ! TIME2 I.A. O Time of second field. (r)
156  ! IOTST Int. O Test indictor for reading.
157  ! 1 : File not found.
158  ! 0 : Fields read.
159  ! -1 : Past end of file.
160  ! IMOD Int. I Optional grid number, defaults to 1.
161  ! ----------------------------------------------------------------
162  ! (w) used for write only
163  ! (r) used for write only
164  !
165  ! 4. Subroutines used :
166  !
167  ! See module documentation.
168  !
169  ! 5. Called by :
170  !
171  ! Name Type Module Description
172  ! ----------------------------------------------------------------
173  ! W3WAVE Subr. W3WAVEMD Actual wave model routine.
174  ! ----------------------------------------------------------------
175  !
176  ! 6. Error messages :
177  !
178  ! Tests on INXOUT, file status and data present in file.
179  !
180  ! 7. Remarks :
181  !
182  ! - Array dimensions are tested in W3IOGR.
183  ! - Spectra are stored as frequency (sigma) spectra to guarantee
184  ! conservation under grid transformation.
185  ! - At the moment it is mplicitly assumed that the number of
186  ! spectral components is larger that the number of spectra
187  ! per time step per file.
188  ! - Dump option used in multi-grid model.
189  !
190  ! 8. Structure :
191  !
192  ! See source code.
193  !
194  ! 9. Switches :
195  !
196  ! !/SHRD Switch for shared / distributed memory architecture.
197  ! !/DIST Id.
198  !
199  ! !/S Enable subroutine tracing.
200  ! !/T General test output.
201  ! !/T0 Point info test output.
202  ! !/T1 Wave heights at input/output points.
203  !
204  ! 10. Source code :
205  !
206  !/ ------------------------------------------------------------------- /
207  USE constants
208  !
209  USE w3gdatmd, ONLY: w3setg
210  USE w3wdatmd, ONLY: w3setw
211  USE w3adatmd, ONLY: w3seta
212  USE w3odatmd, ONLY: w3seto, w3dmo5
213  USE w3cspcmd, ONLY: w3cspc
214  USE w3triamd, ONLY: w3nestug
215  !
216  USE w3gdatmd, ONLY: nk, nth, nspec, nsea, nseal, nx, ny, &
217  x0, y0, sx, sy, gsu, mapsta, mapfs, mapsf, &
218  xfr, fr1, sig2, th, dth, filext, fachfe, &
220  USE w3gdatmd, ONLY: dxymax
221 #ifdef W3_T1
222  USE w3gdatmd, ONLY: sig
223 #endif
224 #ifdef W3_RTD
225  !! Use rotated N-Pole lat/lon and conversion sub. JGLi12Jun2012
226  USE w3gdatmd, ONLY: polat, polon, angld
227  USE w3servmd, ONLY: w3lltoeq, w3eqtoll, w3acturn
228 #endif
229  USE w3wdatmd, ONLY: va
230  USE w3adatmd, ONLY: cg
231  USE w3odatmd, ONLY: ndse, ndst, iaproc, naproc, naperr, napbpt, &
232  nbi, nbi2, nfbpo, nbo, nbo2, ndsl, &
233  nki, nthi, xfri, fr1i, th1i, &
234  ipbpi, isbpi, xbpi, ybpi, rdbpi, &
235  ipbpo, isbpo, xbpo, ybpo, rdbpo, &
236  abpi0, abpin, abpos, flbpi, filer, filew, &
238  USE w3gsrumd
239  !
240  USE w3servmd, ONLY: extcde
241 #ifdef W3_S
242  USE w3servmd, ONLY: strace
243 #endif
244  !
245 #ifdef W3_SMC
246  USE w3psmcmd, ONLY: w3smcgmp
247 #endif
248  !
249  IMPLICIT NONE
250  !/
251  !/ ------------------------------------------------------------------- /
252  !/ Parameter list
253  !/
254  INTEGER, INTENT(IN) :: NDSB
255  INTEGER, INTENT(INOUT) :: TIME1(2)
256  INTEGER, INTENT(OUT) :: TIME2(2), IOTST
257  INTEGER, INTENT(IN), OPTIONAL :: IMOD
258  CHARACTER, INTENT(IN) :: INXOUT*(*)
259  !/
260  !/
261  !/ ------------------------------------------------------------------- /
262  !/ Local parameters
263  !/
264  INTEGER :: IFILE, IERR, I, J, IX, IY, ISEA, &
265  IP, ISP, NPTS, ISOUT, IS, IGRD
266 #ifdef W3_T1
267  INTEGER :: IK, ITH
268 #endif
269 #ifdef W3_S
270  INTEGER, SAVE :: IENT = 0
271 #endif
272 #ifdef W3_T1
273  REAL :: HS, HS0
274 #endif
275 #ifdef W3_RTD
276  !! Declare rotation angle and rotated lat/lon variables for
277  !! boundary points. JGLi12Jun2012
278  REAL, ALLOCATABLE :: Anglbdy(:), ELatbdy(:), ELonbdy(:)
279  REAL :: Spectr(NK*NTH)
280  REAL :: XRLIM, YRLIM
281 #endif
282  REAL, ALLOCATABLE :: TMPSPC(:,:)
283  LOGICAL :: FLOK
284  CHARACTER(LEN=18) :: FILEN
285  CHARACTER(LEN=10) :: VERTST
286  CHARACTER(LEN=32) :: IDTST
287  !/
288  !/ ------------------------------------------------------------------- /
289  !/
290 #ifdef W3_S
291  CALL strace (ient, 'W3IOBC')
292 #endif
293  !
294  iotst = 0
295  !
296  ! test parameter list input ------------------------------------------ *
297  !
298  IF ( PRESENT(imod) ) THEN
299  igrd = imod
300  ELSE
301  igrd = 1
302  END IF
303  !
304  CALL w3seto ( igrd, ndse, ndst )
305  CALL w3setg ( igrd, ndse, ndst )
306  CALL w3setw ( igrd, ndse, ndst )
307  CALL w3seta ( igrd, ndse, ndst )
308  !
309  IF (inxout.NE.'READ' .AND. inxout.NE.'WRITE' .AND. &
310  inxout.NE.'DUMP' ) THEN
311  IF ( iaproc .EQ. naperr ) WRITE (ndse,900) inxout
312  CALL extcde ( 1 )
313  END IF
314  !
315 #ifdef W3_T
316  WRITE (ndst,9000) inxout, filer, filew, filed, ndsb
317 #endif
318  !
319  ! open file ---------------------------------------------------------- *
320  !
321  i = len_trim(filext)
322  j = len_trim(fnmpre)
323  !
324  IF ( inxout.EQ.'READ' .AND. filer ) THEN
325  WRITE (filen,'(A5,A)') 'nest.', filext(:i)
326 #ifdef W3_T
327  WRITE (ndst,9001) filen(:5+i), ndsb
328 #endif
329  OPEN (ndsb,file=fnmpre(:j)//filen(:5+i),form='UNFORMATTED', convert=file_endian, &
330  err=801,iostat=ierr,status='OLD')
331  END IF
332  !
333  IF ( inxout.EQ.'WRITE' .AND. filew ) THEN
334  DO ifile=1, nfbpo
335  ndsl(ifile) = ndsb + ifile - 1
336  WRITE (filen,'(A4,I1,A1,A)') 'nest', ifile, '.', &
337  filext(:i)
338 #ifdef W3_T
339  WRITE (ndst,9001) filen(:6+i), ndsl(ifile)
340 #endif
341  OPEN (ndsl(ifile),file=fnmpre(:j)//filen(:6+i), &
342  form='UNFORMATTED', convert=file_endian,err=800,iostat=ierr)
343  END DO
344  END IF
345  !
346  IF ( inxout.EQ.'DUMP' .AND. filed ) THEN
347  WRITE (filen,'(A5,A)') 'nest.', filext(:i)
348 #ifdef W3_T
349  WRITE (ndst,9001) filen(:5+i), ndsb
350 #endif
351  OPEN (ndsb,file=fnmpre(:j)//filen(:5+i),form='UNFORMATTED', convert=file_endian, &
352  err=800,iostat=ierr)
353  END IF
354  !
355  ! test info ---------------------------------------------------------- *
356  ! ( new files only )
357  ! ... writing
358  !
359  IF ( inxout.EQ.'WRITE' .AND. filew ) THEN
360  IF ( iaproc .EQ. napbpt ) THEN
361  DO ifile=1, nfbpo
362  WRITE (ndsl(ifile)) &
363  idstrbc, verbptbc, nk, nth, xfr, fr1, th(1), &
364  nbo(ifile)-nbo(ifile-1)
365  !
366 #ifdef W3_T
367  WRITE (ndst,9002) ifile, ndsl(ifile), idstrbc, &
368  verbptbc, nbo(ifile)-nbo(ifile-1)
369 #endif
370  !
371 #ifdef W3_RTD
372  ! By running the ww3_grid program the arrays XBPO, YBPO have been
373  ! remapped to standard lat-lon and stored in mod_def.*
374  !
375 #endif
376  WRITE (ndsl(ifile)) &
377  (xbpo(i),i=nbo(ifile-1)+1,nbo(ifile)), &
378  (ybpo(i),i=nbo(ifile-1)+1,nbo(ifile)), &
379  ((ipbpo(i,j),i=nbo(ifile-1)+1,nbo(ifile)),j=1,4),&
380  ((rdbpo(i,j),i=nbo(ifile-1)+1,nbo(ifile)),j=1,4)
381  !
382 #ifdef W3_T0
383  WRITE (ndst,9003)
384  DO i=nbo(ifile-1)+1, nbo(ifile)
385  WRITE (ndst,9004) i-nbo(ifile-1), xbpo(i), &
386  ybpo(i), (ipbpo(i,j),j=1,4), &
387  (rdbpo(i,j),j=1,4)
388  END DO
389 #endif
390  !
391  END DO
392  END IF
393  END IF
394  !
395  ! ... dumping
396  !
397  IF ( inxout.EQ.'DUMP' .AND. filed ) THEN
398  IF ( iaproc .EQ. napbpt ) THEN
399  WRITE (ndsb) idstrbc, verbptbc, nk, nth, xfr, fr1, th(1), nbi
400  !
401 #ifdef W3_T
402  WRITE (ndst,9002) 1, ndsb, idstrbc, verbptbc, nbi
403 #endif
404  !
405  WRITE (ndsb) (xbpi(i),i=1,nbi), (ybpi(i),i=1,nbi), &
406  ((ipbpi(i,j),i=1,nbi),j=1,4), &
407  ((rdbpi(i,j),i=1,nbi),j=1,4)
408  !
409 #ifdef W3_T0
410  WRITE (ndst,9003)
411  DO i=1, nbi
412  WRITE (ndst,9004) i, xbpi(i), ybpi(i), &
413  (ipbpi(i,j),j=1,4), (rdbpi(i,j),j=1,4)
414  END DO
415 #endif
416  !
417  END IF
418  END IF
419  !
420  ! ... reading
421  !
422  IF ( inxout.EQ.'READ' .AND. filer ) THEN
423  !
424  READ (ndsb,err=803,iostat=ierr) &
425  idtst, vertst, nki, nthi, xfri, fr1i, th1i, nbi
426  !
427 #ifdef W3_T
428  WRITE (ndst,9002) 1, ndsb, idtst, vertst, nbi
429 #endif
430  !
431  IF ( idtst .NE. idstrbc ) THEN
432  IF ( iaproc .EQ. naperr ) &
433  WRITE (ndse,901) idtst, idstrbc
434  CALL extcde ( 10 )
435  END IF
436  IF ( vertst .NE. verbptbc ) THEN
437  IF ( iaproc .EQ. naperr ) &
438  WRITE (ndse,902) vertst, verbptbc
439  CALL extcde ( 11 )
440  END IF
441  !
442  ! Determines if the spectrum in nest file needs to be converted
443  !
444  spconv = nki.NE.nk .OR. nthi.NE.nth .OR. &
445  abs(xfri/xfr-1.).GT.0.01 .OR. &
446  abs(fr1i/fr1-1.).GT.0.01 .OR. &
447  abs(th1i-th(1)).GT.0.01*dth
448  !
449  CALL w3dmo5 ( igrd, ndse, ndst, 1 )
450  !
451  READ (ndsb,err=803,iostat=ierr) &
452  (xbpi(i),i=1,nbi), (ybpi(i),i=1,nbi), &
453  ((ipbpi(i,j),i=1,nbi),j=1,4), &
454  ((rdbpi(i,j),i=1,nbi),j=1,4)
455  !
456 #ifdef W3_RTD
457  ! All boundary conditions position arrays XBPI, YBPI are defined
458  ! in standard lat/lon coordinates. If Polat = 90. (and Polon = -180.),
459  ! the b.c. positions don't need to be remapped
460  IF ( polat < 90. ) THEN
461  !! Convert standard into rotated lat/lon. JGLi12Jun2012
462  ALLOCATE ( anglbdy(nbi), elatbdy(nbi), elonbdy(nbi) )
463 
464  CALL w3lltoeq ( ybpi, xbpi, elatbdy, elonbdy, &
465  & anglbdy, polat, polon, nbi )
466 
467  xbpi = elonbdy
468  ybpi = elatbdy
469  !! W3LLTOEQ outputs longitudes on 0->360 degree grid
470  !! Next section will revise to -180->180 convention if required
471  !! by nested model rotated grid; determined by X0 lon value
472  IF ( x0 .LT. 0.0 ) THEN
473  DO i=1, nbi
474  IF ( xbpi(i) .GT. 180.0) xbpi(i) = xbpi(i) - 360.0
475  ENDDO
476  END IF
477  !! The old (4.18) W3GFPT was very strict so this loop reassigns RTD
478  !! values to within a tolerance of the boundary - possibly this is
479  !! no longer required after the 20-Jan-2017 change?
480  xrlim = x0 + (nx-1) * sx
481  yrlim = y0 + (ny-1) * sy
482  DO i=1, nbi
483  IF ( abs(xbpi(i) - x0) .LT. sx/4.0 ) xbpi(i) = x0
484  IF ( abs(ybpi(i) - y0) .LT. sy/4.0 ) ybpi(i) = y0
485  IF ( abs(xbpi(i) - xrlim) .LT. sx/4.0 ) xbpi(i) = xrlim
486  IF ( abs(ybpi(i) - yrlim) .LT. sy/4.0 ) ybpi(i) = yrlim
487  ENDDO
488 
489  DEALLOCATE ( anglbdy, elatbdy, elonbdy )
490 
491  END IF ! ( Polat < 90. )
492 
493 #endif
494  flok = .true.
495  IF (gtype .EQ. ungtype) THEN
496  CALL w3nestug(dxymax,flok)
497 #ifdef W3_SMC
498  !Li For SMC grid check whether boundary points are within cell area.
499  ELSE IF( gtype .EQ. smctype ) THEN
500  CALL w3smcgmp( igrd, nbi, xbpi, ybpi, isbpi )
501  IF ( iaproc .EQ. naperr ) WRITE (ndse,920) &
502  ( isbpi(i), xbpi(i), ybpi(i), i=1,nbi )
503 #endif
504  ELSE
505  DO i=1, nbi
506  ! W3GFTP: find the nearest grid point to the input boundary point
507  ! DCIN=0.1 is the distance outside of source grid in units of
508  ! cell width to treat target point as inside the source grid.
509  IF ( w3gfpt( gsu, xbpi(i), ybpi(i), ix, iy, dcin=0.1 ) ) THEN
510  IF ( abs(mapsta(iy,ix)) .NE. 2 ) THEN
511  IF ( iaproc .EQ. naperr ) &
512  WRITE (ndse,909) ix, iy, abs(mapsta(iy,ix))
513  flok = .false.
514  END IF
515  ELSE
516  IF ( iaproc .EQ. naperr ) &
517  WRITE (ndse,910) i, xbpi(i), ybpi(i)
518  CALL extcde ( 12 )
519  END IF
520  isbpi(i) = mapfs(iy,ix)
521  END DO
522  END IF
523  !
524 #ifdef W3_T0
525  WRITE (ndst,9003)
526  DO i=1, nbi
527  WRITE (ndst,9005) i, isbpi(i), xbpi(i), ybpi(i), &
528  (ipbpi(i,j),j=1,4), (rdbpi(i,j),j=1,4)
529  END DO
530 #endif
531  !
532  IF ( .NOT.flok ) CALL extcde ( 20 )
533  !
534  DO isea=1, nsea
535  ix = mapsf(isea,1)
536  iy = mapsf(isea,2)
537  IF ( abs(mapsta(iy,ix)) .EQ. 2 ) THEN
538  flok = .false.
539  DO i=1, nbi
540  IF ( isea .EQ. isbpi(i) ) flok = .true.
541  END DO
542  IF ( .NOT.flok .AND. iaproc.EQ.naperr ) &
543  WRITE (ndse,911) ix, iy
544  END IF
545  END DO
546  !
547  ! Read first time and allocate ABPI0/N
548  !
549  READ (ndsb,END=810,ERR=810) TIME2, nbi2
550  backspace(ndsb)
551 #ifdef W3_T
552  WRITE (ndst,9012) ndsb, time2, nbi2
553 #endif
554  CALL w3dmo5 ( igrd, ndse, ndst, 3 )
555  !
556  END IF
557  !
558  ! Save previous spectra on read -------------------------------------- *
559  !
560  IF ( inxout.EQ.'READ' .AND. .NOT.filer ) THEN
561 #ifdef W3_T
562  WRITE (ndst,9020)
563 #endif
564  time1 = time2
565  abpi0(:,1:nbi2) = abpin(:,1:nbi2)
566  END IF
567  !
568  ! TIME --------------------------------------------------------------- *
569  !
570  IF ( inxout .EQ. 'WRITE' ) THEN
571  DO ifile=1, nfbpo
572  npts = nbo2(ifile) - nbo2(ifile-1)
573  WRITE (ndsl(ifile)) time1, npts
574 #ifdef W3_T
575  WRITE (ndst,9010) ifile, ndsl(ifile), time1, npts
576 #endif
577  END DO
578  END IF
579  !
580  IF ( inxout .EQ. 'DUMP' ) THEN
581  WRITE (ndsb) time1, nbi2
582 #ifdef W3_T
583  WRITE (ndst,9011) ndsb, time1, nbi2
584 #endif
585  END IF
586  !
587  IF ( inxout .EQ. 'READ' ) THEN
588  READ (ndsb,err=810,END=810) TIME2, nbi2
589 #ifdef W3_T
590  WRITE (ndst,9011) ndsb, time2, nbi2
591 #endif
592  END IF
593  !
594  ! Spectra ------------------------------------------------------------ *
595  !
596  IF ( inxout .EQ. 'WRITE' ) THEN
597  !
598 #ifdef W3_T1
599  WRITE (ndst,9040)
600 #endif
601  !
602  DO ifile=1, nfbpo
603  DO isout=nbo2(ifile-1)+1, nbo2(ifile)
604  !
605  isea = isbpo(isout)
606  !
607  ! ... Shared memory version data gather
608  !
609 #ifdef W3_SHRD
610  DO is=1, nspec
611  abpos(is,isout) = va(is,isea) * sig2(is) / &
612  cg(1+(is-1)/nth,isea)
613  END DO
614 #endif
615  !
616  ! ... Distributed memory version data gather
617  ! ( Array pre-filled in W3WAVE )
618  !
619 #ifdef W3_DIST
620  DO is=1, nspec
621  abpos(is,isout) = abpos(is,isout) * sig2(is) / &
622  cg(1+(is-1)/nth,isea)
623  END DO
624 #endif
625  !
626 #ifdef W3_RTD
627  ! Polat == 90. means the grid is standard lat-lon, and the spectra
628  ! need not be rotated back
629  IF ( polat < 90. ) THEN
630  ! Added spectral turning for rotated grid
631  ! (rotate back to standard pole)
632  spectr = abpos(:,isout)
633  CALL w3acturn( nth, nk, -angld(isea), spectr )
634  abpos(:,isout) = spectr
635  END IF
636 #endif
637  !
638  WRITE (ndsl(ifile)) (abpos(is,isout),is=1,nspec)
639  !
640 #ifdef W3_T1
641  hs = 0.
642  DO ik=1, nk
643  DO ith=1, nth
644  is = ith + (ik-1)*nth
645  hs = hs + abpos(is,isout)*sig(ik)
646  END DO
647  END DO
648  hs = 4. * sqrt( hs * dth * 0.5 * (xfr-1./xfr) )
649  WRITE (ndst,9041) ndsl(ifile), isout, isea, hs
650 #endif
651  !
652  END DO
653  END DO
654  !
655  END IF
656  !
657  IF ( inxout .EQ. 'DUMP' ) THEN
658  DO i=1, nbi2
659  WRITE (ndsb) abpin(:,i)
660  END DO
661  END IF
662  !
663  IF ( inxout .EQ. 'READ' ) THEN
664  !
665  IF ( .NOT. spconv ) THEN
666  DO ip=1, nbi2
667  READ (ndsb,err=803,iostat=ierr) abpin(:,ip)
668  END DO
669  ELSE
670  !
671  ! In this case the spectral resolution is not compatible and
672  ! the spectrum TMPSPC in nest file must be re-gridded into ABPIN to fit the model run
673  ! spectral conversion is done by W3CSPC in w3cspcmd.ftn
674  !
675  ALLOCATE ( tmpspc(nki*nthi,nbi2) )
676  DO ip=1, nbi2
677  READ (ndsb,err=803,iostat=ierr) tmpspc(:,ip)
678  END DO
679  CALL w3cspc ( tmpspc , nki, nthi, xfri, fr1i, th1i, &
680  abpin(:,1:nbi2),nk, nth, xfr, fr1, th(1),&
681  nbi2, ndst, ndse, fachfe )
682  DEALLOCATE ( tmpspc )
683  END IF
684  !
685 #ifdef W3_T1
686  WRITE (ndst,9042)
687  DO ip=1, nbi2
688  hs = 0.
689  hs0 = 0.
690  DO isp=1, nspec
691  hs = hs + abpin(isp,ip)*sig2(isp)
692  IF ( .NOT.filer ) hs0 = hs0 + abpi0(isp,ip)*sig2(isp)
693  END DO
694  hs = 4. * sqrt( hs * dth * 0.5 * (xfr-1./xfr) )
695  hs0 = 4. * sqrt( hs0 * dth * 0.5 * (xfr-1./xfr) )
696  WRITE (ndst,9043) ip, hs0, hs
697  END DO
698 #endif
699  !
700  END IF
701  !
702  ! Set first spectra on first read ------------------------------------ *
703  !
704  IF ( inxout.EQ.'READ' .AND. filer ) THEN
705 #ifdef W3_T
706  WRITE (ndst,9021)
707 #endif
708  time1 = time2
709  DO ip=1, nbi2
710  abpi0(:,ip) = abpin(:,ip)
711  END DO
712  abpi0(:,0) = 0.
713  abpin(:,0) = 0.
714  END IF
715  !
716  ! Reset flags -------------------------------------------------------- *
717  !
718  IF ( inxout .EQ. 'WRITE' ) filew = .false.
719  IF ( inxout .EQ. 'DUMP' ) filed = .false.
720  IF ( inxout .EQ. 'READ' ) filer = .false.
721  !
722  RETURN
723  !
724  ! Escape locations IO errors
725  !
726 800 CONTINUE
727  IF ( iaproc .EQ. naperr ) WRITE (ndse,1000) filen, ierr
728  CALL extcde ( 40 )
729  !
730 801 CONTINUE
731  IF ( iaproc .EQ. naperr ) WRITE (ndse,1001) imod
732  iotst = 1
733  flbpi = .false.
734  RETURN
735  !
736 802 CONTINUE
737  IF ( iaproc .EQ. naperr ) WRITE (ndse,1002)
738  CALL extcde ( 41 )
739  !
740 803 CONTINUE
741  IF ( iaproc .EQ. naperr ) WRITE (ndse,1003) ierr
742  CALL extcde ( 42 )
743  !
744 810 CONTINUE
745  IF ( filer ) THEN
746  IF ( iaproc .EQ. naperr ) WRITE (ndse,1010)
747  CALL extcde ( 43 )
748  END IF
749  !
750 #ifdef W3_T
751  WRITE (ndst,9022)
752 #endif
753  time1(1) = time2(1)
754  time1(2) = time2(2)
755  DO ip=0, nbi2
756  DO isp=1, nspec
757  abpi0(isp,ip) = abpin(isp,ip)
758  END DO
759  END DO
760  !
761  iotst = -1
762  flbpi = .false.
763  RETURN
764  !
765  ! Formats
766  !
767 900 FORMAT (/' *** WAVEWATCH III ERROR IN W3IOBC :'/ &
768  ' ILLEGAL INXOUT VALUE: ',a/)
769 901 FORMAT (/' *** WAVEWATCH III ERROR IN W3IOBC :'/ &
770  ' ILLEGAL IDSTRBC, READ : ',a/ &
771  ' CHECK : ',a/)
772 902 FORMAT (/' *** WAVEWATCH III ERROR IN W3IOBC :'/ &
773  ' ILLEGAL VEROGR, READ : ',a/ &
774  ' CHECK : ',a/)
775  !
776 909 FORMAT (/' *** WAVEWATCH III ERROR IN W3IOBC :'/ &
777  ' POINT',2i4,' NOT ACTIVE SEA POINT (',i1,')')
778 910 FORMAT (/' *** WAVEWATCH III ERROR IN W3IOBC :'/ &
779  ' POINT',i4,2e14.6,' NOT LOCATED IN GRID')
780 911 FORMAT ( ' *** WAVEWATCH III WARNING : POINT',2i7, &
781  ' WILL NOT BE UPDATED')
782 920 FORMAT (/' *** SMCTYPE mapped boundary cells:'/ ((i8,2f9.3)) )
783  !
784 1000 FORMAT (/' *** WAVEWATCH III ERROR IN W3IOBC : '/ &
785  ' ERROR IN OPENING FILE ',a/ &
786  ' IOSTAT =',i5/)
787  !
788  ! Note: This 1001 error can occur when multi-grid time steps are not
789  ! compatible.
790 1001 FORMAT (/' *** WAVEWATCH III WARNING IN W3IOBC : '/ &
791  ' INPUT FILE WITH BOUNDARY CONDITIONS NOT FOUND'/ &
792  ' BOUNDARY CONDITIONS WILL NOT BE UPDATED ',i5/)
793 1002 FORMAT (/' *** WAVEWATCH III ERROR IN W3IOBC : '/ &
794  ' PREMATURE END OF FILE'/)
795 1003 FORMAT (/' *** WAVEWATCH III ERROR IN W3IOBC : '/ &
796  ' ERROR IN READING FROM FILE'/ &
797  ' IOSTAT =',i5/)
798  !
799 1010 FORMAT (/' *** WAVEWATCH III ERROR IN W3IOBC : '/ &
800  ' NO DATA IN INPUT FILE'/)
801  !
802 #ifdef W3_T
803 9000 FORMAT (' TEST W3IOBC : INXOUT : ',a5/ &
804  ' FLAGS : ',3l2/ &
805  ' UNIT : ',i4)
806 9001 FORMAT (' TEST W3IOBC : OPENING FILE ',a,' (',i2,')')
807 9002 FORMAT (' TEST W3IOBC : FILE # : ',i4/ &
808  ' UNIT : ',i4/ &
809  ' ID : ',a/ &
810  ' VERSION : ',a/ &
811  ' POINTS : ',i4)
812 #endif
813  !
814 #ifdef W3_T0
815 9003 FORMAT (' TEST W3IOBC : POINT DATA ')
816 9004 FORMAT (' ',i3,2e10.3,2x,4i4,2x,4f5.2)
817 9005 FORMAT (' ',i3,i4,2e10.3,2x,4i4,2x,4f5.2)
818 #endif
819  !
820 #ifdef W3_T
821 9010 FORMAT (' TEST W3IOBC : OUTPUT FILE ',i1,' UNIT',i3,' TIME', &
822  i9.8,i7.6,',',i5,' SPECTRA')
823 9011 FORMAT (' TEST W3IOBC : INPUT FILE UNIT',i3,' TIME', &
824  i9.8,i7.6,',',i5,' SPECTRA')
825 9012 FORMAT (' TEST W3IOBC : INPUT FILE UNIT',i3,' TIME', &
826  i9.8,i7.6,',',i5,' SPECTRA (TEST READ)')
827  !
828 9020 FORMAT (' TEST W3IOBC : SAVING OLD DATA')
829 9021 FORMAT (' TEST W3IOBC : SAVING FIRST DATA')
830 9022 FORMAT (' TEST W3IOBC : EOF REACHED')
831 #endif
832  !
833 #ifdef W3_T1
834 9040 FORMAT (' TEST W3IOBC : UNIT, ISOUT, ISEA, HS(NO TAIL) ')
835 9041 FORMAT ( ' ',i3,2i6,f8.2)
836 9042 FORMAT (' TEST W3IOBC : IP, HS(NO TAIL) ')
837 9043 FORMAT ( ' ',i6,2f8.2)
838 #endif
839  !/
840  !/ End of W3IOBC ----------------------------------------------------- /
841  !/

References w3odatmd::abpi0, w3odatmd::abpin, w3odatmd::abpos, w3gdatmd::angld, w3adatmd::cg, w3gdatmd::dth, w3gdatmd::dxymax, w3gdatmd::fachfe, file(), constants::file_endian, w3odatmd::filed, w3odatmd::filer, w3odatmd::filew, w3gdatmd::filext, w3odatmd::flbpi, w3odatmd::fnmpre, w3gdatmd::fr1, w3odatmd::fr1i, w3gdatmd::gsu, w3gdatmd::gtype, w3odatmd::iaproc, idstrbc, w3odatmd::ipbpi, w3odatmd::ipbpo, w3odatmd::isbpi, w3odatmd::isbpo, w3gdatmd::mapfs, w3gdatmd::mapsf, w3gdatmd::mapsta, w3odatmd::napbpt, w3odatmd::naperr, w3odatmd::naproc, w3odatmd::nbi, w3odatmd::nbi2, w3odatmd::nbo, w3odatmd::nbo2, w3odatmd::ndse, w3odatmd::ndsl, w3odatmd::ndst, w3odatmd::nfbpo, w3gdatmd::nk, w3odatmd::nki, w3gdatmd::nsea, w3gdatmd::nseal, w3gdatmd::nspec, w3gdatmd::nth, w3odatmd::nthi, w3gdatmd::nx, w3gdatmd::ny, w3gdatmd::polat, w3gdatmd::polon, w3odatmd::rdbpi, w3odatmd::rdbpo, w3gdatmd::sig, w3gdatmd::sig2, w3gdatmd::smctype, w3odatmd::spconv, w3servmd::strace(), w3gdatmd::sx, w3gdatmd::sy, w3gdatmd::th, w3odatmd::th1i, w3gdatmd::ungtype, w3wdatmd::va, verbptbc, w3servmd::w3acturn(), w3cspcmd::w3cspc(), w3odatmd::w3dmo5(), w3servmd::w3eqtoll(), w3servmd::w3lltoeq(), w3triamd::w3nestug(), w3adatmd::w3seta(), w3gdatmd::w3setg(), w3odatmd::w3seto(), w3wdatmd::w3setw(), w3psmcmd::w3smcgmp(), w3gdatmd::x0, w3odatmd::xbpi, w3odatmd::xbpo, w3gdatmd::xfr, w3odatmd::xfri, w3gdatmd::y0, w3odatmd::ybpi, and w3odatmd::ybpo.

Referenced by w3wavemd::w3wave(), and wminiomd::wmiobg().

Variable Documentation

◆ idstrbc

character(len=32), parameter w3iobcmd::idstrbc = 'WAVEWATCH III BOUNDARY DATA FILE'

Definition at line 77 of file w3iobcmd.F90.

77  CHARACTER(LEN=32), PARAMETER :: &
78  IDSTRBC = 'WAVEWATCH III BOUNDARY DATA FILE'

Referenced by w3bounc(), w3bound(), and w3iobc().

◆ verbptbc

character(len=10), parameter w3iobcmd::verbptbc = '2018-03-01'

Definition at line 76 of file w3iobcmd.F90.

76  CHARACTER(LEN=10), PARAMETER :: VERBPTBC = '2018-03-01'

Referenced by w3bounc(), w3bound(), and w3iobc().

w3gdatmd::nk
integer, pointer nk
Definition: w3gdatmd.F90:1230
w3odatmd::nbo
integer, dimension(:), pointer nbo
Definition: w3odatmd.F90:531
w3gdatmd::nseal
integer, pointer nseal
Definition: w3gdatmd.F90:1097
w3gdatmd::dth
real, pointer dth
Definition: w3gdatmd.F90:1232
w3triamd
Reads triangle and unstructured grid information.
Definition: w3triamd.F90:21
w3odatmd::nki
integer, pointer nki
Definition: w3odatmd.F90:530
w3gdatmd::gsu
type(t_gsu), pointer gsu
Definition: w3gdatmd.F90:1226
w3odatmd::spconv
logical, pointer spconv
Definition: w3odatmd.F90:546
w3adatmd
Define data structures to set up wave model auxiliary data for several models simultaneously.
Definition: w3adatmd.F90:26
w3gdatmd::nspec
integer, pointer nspec
Definition: w3gdatmd.F90:1230
w3odatmd::th1i
real, pointer th1i
Definition: w3odatmd.F90:540
w3gdatmd::ungtype
integer, parameter ungtype
Definition: w3gdatmd.F90:626
w3wdatmd
Define data structures to set up wave model dynamic data for several models simultaneously.
Definition: w3wdatmd.F90:18
w3odatmd::ipbpo
integer, dimension(:,:), pointer ipbpo
Definition: w3odatmd.F90:535
w3adatmd::cg
real, dimension(:,:), pointer cg
Definition: w3adatmd.F90:575
w3gsrumd
Definition: w3gsrumd.F90:17
w3gdatmd::sig
real, dimension(:), pointer sig
Definition: w3gdatmd.F90:1234
w3gdatmd::sy
real, pointer sy
Definition: w3gdatmd.F90:1183
w3odatmd::iaproc
integer, pointer iaproc
Definition: w3odatmd.F90:457
w3gdatmd::fachfe
real, pointer fachfe
Definition: w3gdatmd.F90:1232
w3triamd::w3nestug
subroutine w3nestug(DISTMIN, FLOK)
UGTYPE nesting initialization.
Definition: w3triamd.F90:2184
w3odatmd::ybpi
real, dimension(:), pointer ybpi
Definition: w3odatmd.F90:541
w3gdatmd::ny
integer, pointer ny
Definition: w3gdatmd.F90:1097
w3odatmd::abpin
real, dimension(:,:), pointer abpin
Definition: w3odatmd.F90:541
w3cspcmd
Convert spectra to new discrete spectral grid.
Definition: w3cspcmd.F90:21
w3odatmd::fnmpre
character(len=80) fnmpre
Definition: w3odatmd.F90:330
w3odatmd::isbpo
integer, dimension(:), pointer isbpo
Definition: w3odatmd.F90:535
w3odatmd::nbi
integer, pointer nbi
Definition: w3odatmd.F90:530
w3odatmd::flbpi
logical, pointer flbpi
Definition: w3odatmd.F90:546
w3servmd::w3eqtoll
subroutine w3eqtoll(PHI_EQ, LAMBDA_EQ, PHI, LAMBDA, ANGLED, PHI_POLE, LAMBDA_POLE, POINTS)
Definition: w3servmd.F90:1224
w3wdatmd::va
real, dimension(:,:), pointer va
Definition: w3wdatmd.F90:183
w3odatmd::napbpt
integer, pointer napbpt
Definition: w3odatmd.F90:457
w3gdatmd::th
real, dimension(:), pointer th
Definition: w3gdatmd.F90:1234
w3snl4md::npts
integer, parameter npts
Definition: w3snl4md.F90:129
w3gdatmd::w3setg
subroutine w3setg(IMOD, NDSE, NDST)
Definition: w3gdatmd.F90:2152
w3servmd::w3acturn
subroutine w3acturn(NDirc, NFreq, Alpha, Spectr)
Definition: w3servmd.F90:977
scrip_timers::status
character(len=8), dimension(max_timers), save status
Definition: scrip_timers.f:63
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
w3odatmd::abpos
real, dimension(:,:), pointer abpos
Definition: w3odatmd.F90:541
w3odatmd::nbi2
integer, pointer nbi2
Definition: w3odatmd.F90:530
w3gdatmd::mapfs
integer, dimension(:,:), pointer mapfs
Definition: w3gdatmd.F90:1163
w3odatmd::naperr
integer, pointer naperr
Definition: w3odatmd.F90:457
w3odatmd::filew
logical, pointer filew
Definition: w3odatmd.F90:546
w3odatmd::w3dmo5
subroutine w3dmo5(IMOD, NDSE, NDST, IBLOCK)
Definition: w3odatmd.F90:1321
w3gdatmd::polat
real, pointer polat
Definition: w3gdatmd.F90:1191
w3gdatmd::x0
real, pointer x0
Definition: w3gdatmd.F90:1183
w3gdatmd::nsea
integer, pointer nsea
Definition: w3gdatmd.F90:1097
w3servmd
Definition: w3servmd.F90:3
w3servmd::w3lltoeq
subroutine w3lltoeq(PHI, LAMBDA, PHI_EQ, LAMBDA_EQ, ANGLED, PHI_POLE, LAMBDA_POLE, POINTS)
Definition: w3servmd.F90:1084
w3odatmd::nbo2
integer, dimension(:), pointer nbo2
Definition: w3odatmd.F90:531
w3wdatmd::w3setw
subroutine w3setw(IMOD, NDSE, NDST)
Select one of the WAVEWATCH III grids / models.
Definition: w3wdatmd.F90:660
w3odatmd::ipbpi
integer, dimension(:,:), pointer ipbpi
Definition: w3odatmd.F90:535
w3odatmd::w3seto
subroutine w3seto(IMOD, NDSERR, NDSTST)
Definition: w3odatmd.F90:1523
w3gdatmd::nth
integer, pointer nth
Definition: w3gdatmd.F90:1230
w3odatmd
Definition: w3odatmd.F90:3
w3odatmd::fr1i
real, pointer fr1i
Definition: w3odatmd.F90:540
w3cspcmd::w3cspc
subroutine w3cspc(SP1, NFR1, NTH1, XF1, FR1, TH1, SP2, NFR2, NTH2, XF2, FR2, TH2, NSP, NDST, NDSE, FTL)
Convert a set of spectra to a new spectral grid.
Definition: w3cspcmd.F90:146
w3odatmd::nfbpo
integer, pointer nfbpo
Definition: w3odatmd.F90:530
w3gdatmd::mapsf
integer, dimension(:,:), pointer mapsf
Definition: w3gdatmd.F90:1163
w3odatmd::naproc
integer, pointer naproc
Definition: w3odatmd.F90:457
w3odatmd::xbpi
real, dimension(:), pointer xbpi
Definition: w3odatmd.F90:541
w3gdatmd::smctype
integer, parameter smctype
Definition: w3gdatmd.F90:627
file
file(STRINGS ${CMAKE_BINARY_DIR}/switch switch_strings) separate_arguments(switches UNIX_COMMAND $
Definition: CMakeLists.txt:3
w3odatmd::xfri
real, pointer xfri
Definition: w3odatmd.F90:540
w3gdatmd::sig2
real, dimension(:), pointer sig2
Definition: w3gdatmd.F90:1234
w3gdatmd::fr1
real, pointer fr1
Definition: w3gdatmd.F90:1232
w3odatmd::xbpo
real, dimension(:), pointer xbpo
Definition: w3odatmd.F90:541
w3servmd::strace
subroutine strace(IENT, SNAME)
Definition: w3servmd.F90:148
w3odatmd::rdbpi
real, dimension(:,:), pointer rdbpi
Definition: w3odatmd.F90:541
w3gdatmd::gtype
integer, pointer gtype
Definition: w3gdatmd.F90:1094
w3odatmd::abpi0
real, dimension(:,:), pointer abpi0
Definition: w3odatmd.F90:541
w3odatmd::filer
logical, pointer filer
Definition: w3odatmd.F90:546
w3gdatmd::xfr
real, pointer xfr
Definition: w3gdatmd.F90:1232
w3gdatmd::y0
real, pointer y0
Definition: w3gdatmd.F90:1183
w3gdatmd::dxymax
real, pointer dxymax
Definition: w3gdatmd.F90:1133
w3psmcmd
Spherical Multiple-Cell (SMC) grid routines.
Definition: w3psmcmd.F90:18
w3gdatmd::sx
real, pointer sx
Definition: w3gdatmd.F90:1183
w3odatmd::ndst
integer, pointer ndst
Definition: w3odatmd.F90:456
constants
Define some much-used constants for global use (all defined as PARAMETER).
Definition: constants.F90:20
w3gdatmd
Definition: w3gdatmd.F90:16
w3odatmd::ybpo
real, dimension(:), pointer ybpo
Definition: w3odatmd.F90:541
w3odatmd::ndsl
integer, dimension(:), pointer ndsl
Definition: w3odatmd.F90:531
constants::file_endian
character(*), parameter file_endian
FILE_ENDIAN Filled by preprocessor with 'big_endian', 'little_endian', or 'native'.
Definition: constants.F90:86
w3gdatmd::angld
real, dimension(:), pointer angld
Definition: w3gdatmd.F90:1192
w3odatmd::filed
logical, pointer filed
Definition: w3odatmd.F90:546
w3odatmd::isbpi
integer, dimension(:), pointer isbpi
Definition: w3odatmd.F90:535
w3gdatmd::nx
integer, pointer nx
Definition: w3gdatmd.F90:1097
w3odatmd::nthi
integer, pointer nthi
Definition: w3odatmd.F90:530
w3psmcmd::w3smcgmp
subroutine w3smcgmp(IMOD, NC, XLon, YLat, IDCl)
Map lat-lon points to SMC grid cells.
Definition: w3psmcmd.F90:3804
w3gdatmd::mapsta
integer, dimension(:,:), pointer mapsta
Definition: w3gdatmd.F90:1163
w3gdatmd::polon
real, pointer polon
Definition: w3gdatmd.F90:1191
w3odatmd::rdbpo
real, dimension(:,:), pointer rdbpo
Definition: w3odatmd.F90:541
w3gdatmd::filext
character(len=13), pointer filext
Definition: w3gdatmd.F90:1224