Processing of boundary data output.
More...
|
| subroutine | w3iobc (INXOUT, NDSB, TIME1, TIME2, IOTST, IMOD) |
| | Write/read boundary conditions file(s). More...
|
| |
|
| character(len=10), parameter | verbptbc = '2018-03-01' |
| |
| character(len=32), parameter | idstrbc = 'WAVEWATCH III BOUNDARY DATA FILE' |
| |
Processing of boundary data output.
- Author
- H. L. Tolman
- Date
- 01-Mar-2018
◆ 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] | INXOUT | Test string for read/write. |
| [in,out] | NDSB | Data set unit number. |
| [in,out] | TIME1 | Present time (w), time of first field (r). |
| [in,out] | TIME2 | Time of second field. |
| [in,out] | IOTST | Test indictor for reading. |
| [in,out] | IMOD | Optional grid number, defaults to 1. |
- Author
- H. L. Tolman
- Date
- 20-Jan-2017
Definition at line 99 of file w3iobcmd.F90.
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*(*)
264 INTEGER :: IFILE, IERR, I, J, IX, IY, ISEA, &
265 IP, ISP, NPTS, ISOUT, IS, IGRD
270 INTEGER,
SAVE :: IENT = 0
278 REAL,
ALLOCATABLE :: Anglbdy(:), ELatbdy(:), ELonbdy(:)
279 REAL :: Spectr(NK*NTH)
282 REAL,
ALLOCATABLE :: TMPSPC(:,:)
284 CHARACTER(LEN=18) :: FILEN
285 CHARACTER(LEN=10) :: VERTST
286 CHARACTER(LEN=32) :: IDTST
291 CALL strace (ient,
'W3IOBC')
298 IF (
PRESENT(imod) )
THEN
309 IF (inxout.NE.
'READ' .AND. inxout.NE.
'WRITE' .AND. &
310 inxout.NE.
'DUMP' )
THEN
324 IF ( inxout.EQ.
'READ' .AND.
filer )
THEN
325 WRITE (filen,
'(A5,A)')
'nest.', filext(:i)
327 WRITE (
ndst,9001) filen(:5+i), ndsb
330 err=801,iostat=ierr,
status=
'OLD')
333 IF ( inxout.EQ.
'WRITE' .AND.
filew )
THEN
335 ndsl(ifile) = ndsb + ifile - 1
336 WRITE (filen,
'(A4,I1,A1,A)')
'nest', ifile,
'.', &
339 WRITE (
ndst,9001) filen(:6+i),
ndsl(ifile)
342 form=
'UNFORMATTED', convert=
file_endian,err=800,iostat=ierr)
346 IF ( inxout.EQ.
'DUMP' .AND.
filed )
THEN
347 WRITE (filen,
'(A5,A)')
'nest.', filext(:i)
349 WRITE (
ndst,9001) filen(:5+i), ndsb
359 IF ( inxout.EQ.
'WRITE' .AND.
filew )
THEN
362 WRITE (
ndsl(ifile)) &
363 idstrbc, verbptbc, nk, nth, xfr, fr1, th(1), &
367 WRITE (
ndst,9002) ifile,
ndsl(ifile), idstrbc, &
368 verbptbc,
nbo(ifile)-
nbo(ifile-1)
376 WRITE (
ndsl(ifile)) &
384 DO i=
nbo(ifile-1)+1,
nbo(ifile)
397 IF ( inxout.EQ.
'DUMP' .AND.
filed )
THEN
399 WRITE (ndsb) idstrbc, verbptbc, nk, nth, xfr, fr1, th(1),
nbi
402 WRITE (
ndst,9002) 1, ndsb, idstrbc, verbptbc,
nbi
422 IF ( inxout.EQ.
'READ' .AND.
filer )
THEN
424 READ (ndsb,err=803,iostat=ierr) &
428 WRITE (
ndst,9002) 1, ndsb, idtst, vertst,
nbi
431 IF ( idtst .NE. idstrbc )
THEN
433 WRITE (
ndse,901) idtst, idstrbc
436 IF ( vertst .NE. verbptbc )
THEN
438 WRITE (
ndse,902) vertst, verbptbc
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
451 READ (ndsb,err=803,iostat=ierr) &
460 IF (
polat < 90. )
THEN
462 ALLOCATE ( anglbdy(
nbi), elatbdy(
nbi), elonbdy(
nbi) )
464 CALL w3lltoeq (
ybpi,
xbpi, elatbdy, elonbdy, &
472 IF ( x0 .LT. 0.0 )
THEN
480 xrlim = x0 + (nx-1) * sx
481 yrlim = y0 + (ny-1) * sy
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
489 DEALLOCATE ( anglbdy, elatbdy, elonbdy )
495 IF (gtype .EQ. ungtype)
THEN
499 ELSE IF( gtype .EQ. smctype )
THEN
509 IF ( w3gfpt( gsu,
xbpi(i),
ybpi(i), ix, iy, dcin=0.1 ) )
THEN
510 IF ( abs(mapsta(iy,ix)) .NE. 2 )
THEN
512 WRITE (
ndse,909) ix, iy, abs(mapsta(iy,ix))
520 isbpi(i) = mapfs(iy,ix)
532 IF ( .NOT.flok )
CALL extcde ( 20 )
537 IF ( abs(mapsta(iy,ix)) .EQ. 2 )
THEN
540 IF ( isea .EQ.
isbpi(i) ) flok = .true.
543 WRITE (
ndse,911) ix, iy
549 READ (ndsb,
END=810,ERR=810) TIME2,
nbi2
560 IF ( inxout.EQ.
'READ' .AND. .NOT.
filer )
THEN
570 IF ( inxout .EQ.
'WRITE' )
THEN
580 IF ( inxout .EQ.
'DUMP' )
THEN
581 WRITE (ndsb) time1,
nbi2
587 IF ( inxout .EQ.
'READ' )
THEN
588 READ (ndsb,err=810,
END=810) TIME2,
nbi2
596 IF ( inxout .EQ.
'WRITE' )
THEN
603 DO isout=
nbo2(ifile-1)+1,
nbo2(ifile)
611 abpos(is,isout) =
va(is,isea) * sig2(is) / &
612 cg(1+(is-1)/nth,isea)
621 abpos(is,isout) =
abpos(is,isout) * sig2(is) / &
622 cg(1+(is-1)/nth,isea)
629 IF (
polat < 90. )
THEN
632 spectr =
abpos(:,isout)
633 CALL w3acturn( nth, nk, -
angld(isea), spectr )
634 abpos(:,isout) = spectr
638 WRITE (
ndsl(ifile)) (
abpos(is,isout),is=1,nspec)
644 is = ith + (ik-1)*nth
645 hs = hs +
abpos(is,isout)*sig(ik)
648 hs = 4. * sqrt( hs * dth * 0.5 * (xfr-1./xfr) )
649 WRITE (
ndst,9041)
ndsl(ifile), isout, isea, hs
657 IF ( inxout .EQ.
'DUMP' )
THEN
659 WRITE (ndsb)
abpin(:,i)
663 IF ( inxout .EQ.
'READ' )
THEN
667 READ (ndsb,err=803,iostat=ierr)
abpin(:,ip)
677 READ (ndsb,err=803,iostat=ierr) tmpspc(:,ip)
680 abpin(:,1:
nbi2),nk, nth, xfr, fr1, th(1),&
682 DEALLOCATE ( tmpspc )
691 hs = hs +
abpin(isp,ip)*sig2(isp)
692 IF ( .NOT.
filer ) hs0 = hs0 +
abpi0(isp,ip)*sig2(isp)
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
704 IF ( inxout.EQ.
'READ' .AND.
filer )
THEN
718 IF ( inxout .EQ.
'WRITE' )
filew = .false.
719 IF ( inxout .EQ.
'DUMP' )
filed = .false.
720 IF ( inxout .EQ.
'READ' )
filer = .false.
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/ &
772 902
FORMAT (/
' *** WAVEWATCH III ERROR IN W3IOBC :'/ &
773 ' ILLEGAL VEROGR, READ : ',a/ &
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)) )
784 1000
FORMAT (/
' *** WAVEWATCH III ERROR IN W3IOBC : '/ &
785 ' ERROR IN OPENING FILE ',a/ &
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'/ &
799 1010
FORMAT (/
' *** WAVEWATCH III ERROR IN W3IOBC : '/ &
800 ' NO DATA IN INPUT FILE'/)
803 9000
FORMAT (
' TEST W3IOBC : INXOUT : ',a5/ &
806 9001
FORMAT (
' TEST W3IOBC : OPENING FILE ',a,
' (',i2,
')')
807 9002
FORMAT (
' TEST W3IOBC : FILE # : ',i4/ &
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)
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)')
828 9020
FORMAT (
' TEST W3IOBC : SAVING OLD DATA')
829 9021
FORMAT (
' TEST W3IOBC : SAVING FIRST DATA')
830 9022
FORMAT (
' TEST W3IOBC : EOF REACHED')
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)
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().
◆ idstrbc
| character(len=32), parameter w3iobcmd::idstrbc = 'WAVEWATCH III BOUNDARY DATA FILE' |
◆ verbptbc
| character(len=10), parameter w3iobcmd::verbptbc = '2018-03-01' |
integer, dimension(:), pointer nbo
Reads triangle and unstructured grid information.
Define data structures to set up wave model auxiliary data for several models simultaneously.
integer, parameter ungtype
Define data structures to set up wave model dynamic data for several models simultaneously.
integer, dimension(:,:), pointer ipbpo
real, dimension(:,:), pointer cg
real, dimension(:), pointer sig
subroutine w3nestug(DISTMIN, FLOK)
UGTYPE nesting initialization.
real, dimension(:), pointer ybpi
real, dimension(:,:), pointer abpin
Convert spectra to new discrete spectral grid.
integer, dimension(:), pointer isbpo
subroutine w3eqtoll(PHI_EQ, LAMBDA_EQ, PHI, LAMBDA, ANGLED, PHI_POLE, LAMBDA_POLE, POINTS)
real, dimension(:,:), pointer va
real, dimension(:), pointer th
subroutine w3setg(IMOD, NDSE, NDST)
subroutine w3acturn(NDirc, NFreq, Alpha, Spectr)
character(len=8), dimension(max_timers), save status
subroutine w3seta(IMOD, NDSE, NDST)
Select one of the WAVEWATCH III grids / models.
real, dimension(:,:), pointer abpos
integer, dimension(:,:), pointer mapfs
subroutine w3dmo5(IMOD, NDSE, NDST, IBLOCK)
subroutine w3lltoeq(PHI, LAMBDA, PHI_EQ, LAMBDA_EQ, ANGLED, PHI_POLE, LAMBDA_POLE, POINTS)
integer, dimension(:), pointer nbo2
subroutine w3setw(IMOD, NDSE, NDST)
Select one of the WAVEWATCH III grids / models.
integer, dimension(:,:), pointer ipbpi
subroutine w3seto(IMOD, NDSERR, NDSTST)
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.
integer, dimension(:,:), pointer mapsf
real, dimension(:), pointer xbpi
integer, parameter smctype
file(STRINGS ${CMAKE_BINARY_DIR}/switch switch_strings) separate_arguments(switches UNIX_COMMAND $
real, dimension(:), pointer sig2
real, dimension(:), pointer xbpo
subroutine strace(IENT, SNAME)
real, dimension(:,:), pointer rdbpi
real, dimension(:,:), pointer abpi0
Spherical Multiple-Cell (SMC) grid routines.
Define some much-used constants for global use (all defined as PARAMETER).
real, dimension(:), pointer ybpo
integer, dimension(:), pointer ndsl
character(*), parameter file_endian
FILE_ENDIAN Filled by preprocessor with 'big_endian', 'little_endian', or 'native'.
real, dimension(:), pointer angld
integer, dimension(:), pointer isbpi
subroutine w3smcgmp(IMOD, NC, XLon, YLat, IDCl)
Map lat-lon points to SMC grid cells.
integer, dimension(:,:), pointer mapsta
real, dimension(:,:), pointer rdbpo
character(len=13), pointer filext