Contains wave model subroutine, w3wave.
More...
|
| subroutine | w3wave (IMOD, ODAT, TEND, STAMP, NO_OUT ifdef W3_OASIS |
| | Run WAVEWATCH III for a given time interval. More...
|
| |
| subroutine | w3gath (ISPEC, FIELD) |
| | Gather spectral bin information into a propagation field array. More...
|
| |
| subroutine | w3scat (ISPEC, MAPSTA, FIELD) |
| | Scatter data back to spectral storage after propagation. More...
|
| |
| subroutine | w3nmin (MAPSTA, FLAG0) |
| | Check minimum number of active sea points at given processor to evaluate the need for a MPI_BARRIER call. More...
|
| |
Contains wave model subroutine, w3wave.
- Author
- H. L. Tolman
- Date
- 22-Mar-2021
◆ w3gath()
| subroutine w3wavemd::w3gath |
( |
integer, intent(in) |
ISPEC, |
|
|
real, dimension(1-ny:ny*(nx+2)), intent(out) |
FIELD |
|
) |
| |
Gather spectral bin information into a propagation field array.
Direct copy or communication calls (MPP version). The field is extracted but not converted.
MPI version requires posing of send and receive calls in W3WAVE to match local calls.
MPI version does not require an MPI_TESTALL call for the posted gather operation as MPI_WAITALL is mandatory to reset persistent communication for next time step.
MPI version allows only two new pre-fetch postings per call to minimize chances to be slowed down by gathers that are not yet needed, while maximizing the pre-loading during the early (low-frequency) calls to the routine where the amount of calculation needed for proagation is the largest.
- Parameters
-
| [in] | ISPEC | Spectral bin considered. |
| [out] | FIELD | Full field to be propagated. |
- Author
- H. L. Tolman
- Date
- 26-Dec-2012
Definition at line 2938 of file w3wavemd.F90.
3041 INTEGER,
INTENT(IN) :: ISPEC
3042 REAL,
INTENT(OUT) :: FIELD(1-NY:NY*(NX+2))
3048 INTEGER :: ISEA, IXY
3051 INTEGER :: STATUS(MPI_STATUS_SIZE,NSPEC), &
3052 IOFF, IERR_MPI, JSEA, ISEA, &
3053 IXY, IS0, IB0, NPST, J
3056 INTEGER,
SAVE :: IENT
3059 CHARACTER(LEN=15) :: STR(MPIBUF), STRT
3065 CALL strace (ient,
'W3GATH')
3075 field(ixy) = a(ispec,isea)
3093 IF (
isploc .EQ. 1 )
THEN
3094 str =
'--------------+'
3095 WRITE (
ndst,9000) str
3176 IF ( is0 .GT.
nsploc )
EXIT
3177 ib0 = 1 + mod(ib0,
mpibuf)
3178 IF (
bstat(ib0) .EQ. 0 )
THEN
3181 ioff = 1 + (is0-1)*
nrqsg2
3188 WRITE (strt(1:7),
'(I2,I5)')
bstat(ib0),
bispl(ib0)
3193 IF ( npst .GE. 2 )
EXIT
3202 IF ( strt(2:2) .EQ.
' ' )
THEN
3203 IF (
bstat(ib0) .EQ. 0 )
THEN
3204 WRITE (strt(1:2),
'(I2)')
bstat(ib0)
3206 WRITE (strt(1:7),
'(I2,I5)')
bstat(ib0),
bispl(ib0)
3221 9000
FORMAT (
' TEST OF BUFFER MANAGEMENT MPI :'/ &
3222 ' -------------------------------'/ &
3223 ' RECORDS ALTERNATELY WRITTEN BY W3GATH AND W3SCAT'/ &
3224 ' FRIST COLLUMN : LOCAL ISPEC'/ &
3225 ' OTHER COLLUMNS : BUFFER STATUS INDICATOR '/ &
3229 ' LOCAL ISPEC FOR BUFFER'/ &
3230 ' A : ACTIVE BUFFER'/ &
3231 ' g/G: START/FINISH RECIEVE'/ &
3232 ' s/S: START/FINISH SEND'/ &
3234 9010
FORMAT (
' |',i4,
' |',8a15)
References w3adatmd::bispl, w3adatmd::bstat, w3gdatmd::dmin, w3adatmd::gstore, w3odatmd::iaproc, w3adatmd::ibfloc, include(), w3parall::init_get_isea(), w3adatmd::irqsg2, w3adatmd::isploc, w3gdatmd::mapsf, w3adatmd::mpibuf, w3odatmd::naproc, w3odatmd::ndst, w3odatmd::notype, w3adatmd::nrqsg2, w3gdatmd::nsea, w3gdatmd::nseal, w3gdatmd::nspec, w3adatmd::nsploc, w3gdatmd::nx, w3gdatmd::ny, w3servmd::strace(), and w3wdatmd::va.
Referenced by w3wave().
◆ w3nmin()
| subroutine w3wavemd::w3nmin |
( |
integer, dimension(ny*nx), intent(in) |
MAPSTA, |
|
|
logical, intent(out) |
FLAG0 |
|
) |
| |
Check minimum number of active sea points at given processor to evaluate the need for a MPI_BARRIER call.
- Parameters
-
| [in] | MAPSTA | Status map for spatial grid. |
| [out] | FLAG0 | Flag to identify 0 as minimum. |
- Author
- H. L. Tolman
- Date
- 28-Dec-2004
Definition at line 3568 of file w3wavemd.F90.
3640 INTEGER,
INTENT(IN) :: MAPSTA(NY*NX)
3641 LOGICAL,
INTENT(OUT) :: FLAG0
3646 INTEGER :: NMIN, IPROC, NLOC, ISEA, IXY
3647 INTEGER :: JSEA, ISPROC
3649 INTEGER,
SAVE :: IENT
3655 CALL strace (ient,
'W3NMIN')
3668 IF (isproc .eq. iproc)
THEN
3670 IF ( mapsta(ixy) .EQ. 1 ) nloc = nloc + 1
3680 WRITE (
ndst,9000) iproc, nloc
3682 nmin = min( nmin , nloc )
3691 WRITE (
ndst,9001) nmin, flag0
3699 9000
FORMAT (
' TEST W3NMIN : IPROC =',i3,
' NLOC =',i5)
3700 9001
FORMAT (
' TEST W3NMIN : NMIN =',i5,
' FLAG0 =',l2)
References w3parall::init_get_jsea_isproc(), w3gdatmd::mapsf, w3odatmd::naproc, w3odatmd::ndst, w3gdatmd::nsea, w3gdatmd::nx, w3gdatmd::ny, and w3servmd::strace().
Referenced by w3wave().
◆ w3scat()
| subroutine w3wavemd::w3scat |
( |
integer, intent(in) |
ISPEC, |
|
|
integer, dimension(ny*nx), intent(in) |
MAPSTA, |
|
|
real, dimension(1-ny:ny*(nx+2)), intent(in) |
FIELD |
|
) |
| |
Scatter data back to spectral storage after propagation.
Direct copy or communication calls (MPP version). See also W3GATH. The field is put back but not converted! MPI persistent communication calls initialize in W3MPII. See W3GATH and W3MPII for additional comments on data buffering.
- Parameters
-
| [in,out] | ISPEC | Spectral bin considered. |
| [in,out] | MAPSTA | Status map for spatial grid. |
| [in,out] | FIELD | Full field to be propagated. |
- Author
- H. L. Tolman
- Date
- 13-Jun-2006
Definition at line 3256 of file w3wavemd.F90.
3357 INTEGER,
INTENT(IN) :: ISPEC, MAPSTA(NY*NX)
3358 REAL,
INTENT(IN) :: FIELD(1-NY:NY*(NX+2))
3364 INTEGER :: ISEA, IXY
3367 INTEGER :: ISEA, IXY, IOFF, IERR_MPI, J, &
3368 STATUS(MPI_STATUS_SIZE,NSPEC), &
3372 INTEGER,
SAVE :: IENT
3375 CHARACTER(LEN=15) :: STR(MPIBUF), STRT
3384 CALL strace (ient,
'W3SCAT')
3392 IF ( mapsta(ixy) .NE. 0 ) a(ispec,isea) = field(ixy)
3419 IF ( mapsta(ixy) .NE. 0 )
sstore(isea,
ibfloc) = field(ixy)
3442 IF (mapsta(ixy) .NE. 0) a(ispec,jsea) =
sstore(isea,
ibfloc)
3454 ib0 = 1 + mod(ib0,
mpibuf)
3455 IF (
bstat(ib0) .EQ. 2 )
THEN
3457 IF (
nrqsg2 .GT. 0 )
THEN
3462 IF ( done .AND.
nrqsg2.GT.0 )
THEN
3470 WRITE (strt(1:7),
'(I2,I5)')
bstat(ib0),
bispl(ib0)
3488 IF (
bstat(ib0) .EQ. 2 )
THEN
3495 WRITE (strt(1:7),
'(I2,I5)')
bstat(ib0),
bispl(ib0)
3518 IF ( strt(2:2) .EQ.
' ' )
THEN
3519 IF (
bstat(ib0) .EQ. 0 )
THEN
3520 WRITE (strt(1:2),
'(I2)')
bstat(ib0)
3522 WRITE (strt(1:7),
'(I2,I5)')
bstat(ib0),
bispl(ib0)
3530 WRITE (ndst,9000) str
3534 IF (
isploc .EQ. 0 )
THEN
3536 str(ib0) =
'--------------+'
3538 WRITE (ndst,9010) str
3550 9000
FORMAT (
' | |',8a15)
3551 9010
FORMAT (
' +-----+',8a15)
References w3adatmd::bispl, w3adatmd::bstat, w3odatmd::iaproc, w3adatmd::ibfloc, include(), w3parall::init_get_isea(), w3adatmd::irqsg2, w3adatmd::isploc, constants::lpdlib, w3gdatmd::mapsf, w3adatmd::mpibuf, w3odatmd::naproc, w3odatmd::ndst, w3adatmd::nrqsg2, w3gdatmd::nsea, w3gdatmd::nseal, w3gdatmd::nspec, w3adatmd::nsploc, w3gdatmd::nx, w3gdatmd::ny, w3adatmd::sstore, w3servmd::strace(), and w3wdatmd::va.
Referenced by w3wave().
◆ w3wave()
| subroutine w3wavemd::w3wave |
( |
integer, intent(in) |
IMOD, |
|
|
integer, dimension(35), intent(in) |
ODAT, |
|
|
integer, dimension(2), intent(in) |
TEND, |
|
|
logical, intent(in), optional |
STAMP, |
|
|
logical, intent(in), optional |
NO_OUT, |
|
|
|
ifdef, |
|
|
|
W3_OASIS |
|
) |
| |
Run WAVEWATCH III for a given time interval.
Currents are updated before winds as currents are used in wind and USTAR processing.
Ice and water levels can be updated only once per call.
If ice or water level time are undefined, the update takes place asap, otherwise around the "half-way point" between the old and new times.
To increase accuracy, the calculation of the intra-spectral propagation is performed in two parts around the spatial propagation.
- Parameters
-
| [in] | IMOD | Model number. |
| [in] | TEND | Ending time of integration. |
| [in] | STAMP | Write time stamp (optional, defaults to T). |
| [in] | NO_OUT | Skip output (optional, defaults to F). |
| [in] | ODAT | |
| [in] | ID_LCOMM | Present only when using W3_OASIS. |
| [in] | TIMEN | Present only when using W3_OASIS. |
- Author
- H. L. Tolman
- Date
- 22-Mar-2021
Definition at line 230 of file w3wavemd.F90.
502 INTEGER,
INTENT(IN) :: IMOD, TEND(2),ODAT(35)
503 LOGICAL,
INTENT(IN),
OPTIONAL :: STAMP, NO_OUT
505 INTEGER,
INTENT(IN),
OPTIONAL :: ID_LCOMM
506 INTEGER,
INTENT(IN),
OPTIONAL :: TIMEN(2)
516 INTEGER,
SAVE :: IENT = 0
519 INTEGER :: TCALC(2), IT, IT0, NT, ITEST, &
520 ITLOC, ITLOCH, NTLOC, ISEA, JSEA, &
521 IX, IY, ISPEC, J, TOUT(2), TLST(2), &
522 REFLED(6), IK, ITH, IS, NKCFL
523 INTEGER :: ISP, IP_glob
524 INTEGER :: TTEST(2),DTTEST
534 INTEGER :: JJ, NDSOFLG
537 INTEGER :: IERR_MPI, NRQMAX
538 INTEGER,
ALLOCATABLE :: STATCO(:,:), STATIO(:,:)
541 REAL :: DTTST, DTTST1, DTTST2, DTTST3, &
542 DTL0, DTI0, DTR0, DTI10, DTI50, &
543 DTGA, DTG, DTGpre, DTRES, &
544 FAC, VGX, VGY, FACK, FACTH, &
545 FACX, XXX, REFLEC(4), &
546 DELX, DELY, DELA, DEPTH, D50, PSIC
547 REAL :: VSioDummy(NSPEC), VDioDummy(NSPEC), VAoldDummy(NSPEC)
548 LOGICAL :: SHAVETOTioDummy
553 REAL,
ALLOCATABLE :: FIELD(:)
554 REAL :: TMP1(4), TMP2(3), TMP3(2), TMP4(2)
556 REAL,
ALLOCATABLE :: WN_I(:)
559 REAL,
ALLOCATABLE :: CIK(:)
564 REAL,
ALLOCATABLE :: TAUWX(:), TAUWY(:)
566 LOGICAL :: FLACT, FLZERO, FLFRST, FLMAP, TSTAMP,&
567 SKIP_O, FLAG_O, FLDDIR, READBC, &
568 FLAG0 = .false., floutg, flpfld, &
569 flpart, local, floutg2
572 LOGICAL :: FLGMPI(0:8)
575 REAL :: FIXEDVISC,FIXEDDENS,FIXEDELAS
576 REAL :: USE_CHENG, USE_CGICE, HICE
578 LOGICAL :: UGDTUPDATE
579 CHARACTER(LEN=8) :: STTIME
580 CHARACTER(LEN=21) :: IDACT
581 CHARACTER(LEN=16) :: OUTID
582 CHARACTER(LEN=23) :: IDTIME
586 REAL :: VS_SPEC(NSPEC)
587 REAL :: VD_SPEC(NSPEC)
591 CHARACTER(LEN=30) :: FOUTNAME
595 REAL :: INDSORT(NSEA), DTCFL1(NSEA)
600 REAL,
ALLOCATABLE :: BACSPEC(:)
625 CALL all_va_integral_print(imod,
"W3WAVEMD, step 1", 1)
634 IF (
PRESENT(stamp) )
THEN
640 IF (
PRESENT(no_out) )
THEN
646 CALL all_va_integral_print(imod,
"W3WAVEMD, step 2", 1)
652 CALL strace (ient,
'W3WAVE')
681 ALLOCATE ( field(
ncel) )
684 ALLOCATE ( field(1-
ny:
ny*(
nx+2)) )
719 flzero = dttst .EQ. 0.
721 WRITE (
ndst,9010) dttst, flzero
723 IF ( dttst .LT. 0. )
THEN
731 IF (
tlev(1) .GE. 0. )
THEN
737 WRITE (
ndst,9011) dtl0
739 IF ( dtl0 .LT. 0. )
THEN
747 CALL all_va_integral_print(imod,
"W3WAVEMD, step 4", 1)
757 WRITE (
ndst,9012) dttst1, dttst2, dttst3
759 IF ( dttst1.LT.0. .OR. dttst2.LT.0. .OR. dttst3.LT.0. )
THEN
763 IF ( dttst2.EQ.0..AND.
itime.EQ.0 )
THEN
776 WRITE (
ndst,9013) dttst1, dttst2, dttst3
778 IF ( dttst1.LT.0. .OR. dttst2.LT.0. .OR. dttst3.LT.0. )
THEN
782 IF ( dttst2.EQ.0..AND.
itime.EQ.0 )
THEN
788 CALL all_va_integral_print(imod,
"W3WAVEMD, step 5", 1)
794 IF (
tice(1) .GE. 0 )
THEN
800 WRITE (
ndst,9014) dti0
802 IF ( dti0 .LT. 0. )
THEN
810 CALL all_va_integral_print(imod,
"W3WAVEMD, step 6", 1)
820 WRITE (
ndst,9017) dttst1, dttst2, dttst3
822 IF ( dttst1.LT.0. .OR. dttst2.LT.0. .OR. dttst3.LT.0. )
THEN
826 IF ( dttst2.EQ.0..AND.
itime.EQ.0 )
THEN
839 WRITE (
ndst,9018) dttst1, dttst2, dttst3
841 IF ( dttst1.LT.0. .OR. dttst2.LT.0. .OR. dttst3.LT.0. )
THEN
845 IF ( dttst2.EQ.0..AND.
itime.EQ.0 )
THEN
854 IF (
tic1(1) .GE. 0 )
THEN
860 WRITE (
ndst,9015) dti10
862 IF ( dti10 .LT. 0. )
THEN
874 IF (
tic5(1) .GE. 0 )
THEN
880 WRITE (
ndst,9016) dti50
882 IF ( dti50 .LT. 0. )
THEN
901 CALL all_va_integral_print(imod,
"W3WAVEMD, step 6.1", 1)
908 IF( use_cheng==1.0 )
THEN
912 IF ( (fixedvisc.LT.0.0).OR.(fixeddens.LT.0.0) .OR. (fixedelas.LT.0.0) )
THEN
914 WRITE(
ndse,*)
'Cheng method requires stationary', &
915 ' and uniform rheology from namelist.'
933 IF ( use_cgice==1.0 )
THEN
947 ALLOCATE(wn_i(
SIZE(
wn(:,isea))))
949 depth = max(
dmin ,
dw(isea) )
958 IF ( use_cheng==1.0 )
THEN
961 ELSEIF (
ic3pars(13).GE.0.0)
THEN
965 WRITE(
ndse,*)
'ICE THICKNESS NOT AVAILABLE FOR CG CALC'
971 fixeddens, fixedelas, depth)
983 IF (icep1(ix,iy)>0.0)
THEN
985 icep2(ix,iy), icep3(ix,iy),icep4(ix,iy),depth)
989 WRITE(
ndse,*)
'ICE PARAMETERS NOT AVAILABLE FOR CG CALC'
1002 IF (
tofrst(1) .GT. 0 )
THEN
1008 IF ( dttst.GE.0. )
THEN
1015 nt = 1 + int( dttst /
dtmax - 0.001 )
1016 dtga = dttst / real(nt)
1017 IF ( dttst .EQ. 0. )
THEN
1027 WRITE (
ndst,9020) it0, nt, dtga
1054 CALL all_va_integral_print(imod,
"Beginning time loop", 1)
1060 call print_memcheck(memunit,
'memcheck_____:'//
' WW3_WAVE TIME LOOP 0')
1064 dtg = real(nint(dtga+dtres+0.0001))
1065 dtres = dtres + dtga - dtg
1066 IF ( abs(dtres) .LT. 0.001 ) dtres = 0.
1069 call print_memcheck(memunit,
'memcheck_____:'//
' WW3_WAVE TIME LOOP 1')
1074 WRITE (
screen,950) idtime, sttime
1077 call print_memcheck(memunit,
'memcheck_____:'//
' WW3_WAVE TIME LOOP 2')
1084 fac = dttst1 / max( 1. , dttst2 )
1093 WRITE (
ndst,9021)
itime, it,
time, flmap, flddir, vgx, vgy, dtg, dtres
1099 #ifdef W3_DEBUGDCXDX
1102 call print_memcheck(memunit,
'memcheck_____:'//
' WW3_WAVE TIME LOOP 3a')
1106 CALL all_va_integral_print(imod,
"Before UCUR", 1)
1113 call print_memcheck(memunit,
'memcheck_____:'//
' WW3_WAVE TIME LOOP 3b')
1122 #ifdef W3_DEBUGDCXDX
1123 WRITE(740+
iaproc,*)
'Before call to UG_GRADIENT for assigning DCXDX/DCXDY array'
1134 call print_memcheck(memunit,
'memcheck_____:'//
' WW3_WAVE TIME LOOP 4')
1136 ELSE IF ( flfrst )
THEN
1146 call print_memcheck(memunit,
'memcheck_____:'//
' WW3_WAVE TIME LOOP 5')
1149 IF ( flfrst )
asf = 1.
1150 CALL w3uwnd ( flfrst, vgx, vgy )
1151 ELSE IF ( flfrst )
THEN
1161 IF (
va(is, jsea) .LT. 0.)
THEN
1162 WRITE(740+
iaproc,*)
'TEST W3WAVE 5',
va(is,jsea)
1167 IF (sum(
va) .NE. sum(
va))
THEN
1168 WRITE(740+
iaproc,*)
'NAN in ACTION 5', ix, iy, sum(
va)
1173 call print_memcheck(memunit,
'memcheck_____:'//
' WW3_WAVE TIME LOOP 6')
1180 CALL all_va_integral_print(imod,
"Before call to W3UINI", 1)
1189 ELSE IF ( flfrst )
THEN
1196 ELSE IF ( flfrst )
THEN
1203 CALL all_va_integral_print(imod,
"Before boundary update", 1)
1208 call print_memcheck(memunit,
'memcheck_____:'//
' WW3_WAVE TIME LOOP 7')
1210 IF (
flbpi .AND. local )
THEN
1213 IF (
tbpin(1) .EQ. -1 )
THEN
1218 IF (readbc.AND.idact(1:1).EQ.
' ') idact(1:1) =
'X'
1220 flact = readbc .OR. flact
1224 IF ( itest .NE. 1 )
CALL w3ubpt
1228 IF ( itest .LT. 0 ) idact(1:1) =
'L'
1229 IF ( itest .GT. 0 ) idact(1:1) =
' '
1230 IF ( .NOT. (readbc.AND.
flbpi) )
EXIT
1234 call print_memcheck(memunit,
'memcheck_____:'//
' WW3_WAVE TIME LOOP 7')
1237 CALL apply_boundary_condition_va
1239 CALL all_va_integral_print(imod,
"After FLBPI and LOCAL", 1)
1245 call print_memcheck(memunit,
'memcheck_____:'//
' WW3_WAVE TIME LOOP 8')
1250 IF (
flice .AND. dti0.NE.0. )
THEN
1252 IF (
tice(1).GE.0 )
THEN
1253 IF ( dti0 .LT. 0. )
THEN
1257 IF ( dttst .LE. 0.5*dti0 ) idact(13:13) =
'U'
1263 IF ( idact(13:13).NE.
' ' )
THEN
1271 CALL all_va_integral_print(imod,
"After FLICE and DTI0", 1)
1276 call print_memcheck(memunit,
'memcheck_____:'//
' WW3_WAVE TIME LOOP 9')
1280 IF (
flic1 .AND. dti10.NE.0. )
THEN
1282 IF (
tic1(1).GE.0 )
THEN
1283 IF ( dti10 .LT. 0. )
THEN
1287 IF ( dttst .LE. 0.5*dti10 ) idact(15:15) =
'U'
1294 IF ( idact(15:15).NE.
' ' )
THEN
1302 call print_memcheck(memunit,
'memcheck_____:'//
' WW3_WAVE TIME LOOP 10')
1307 IF (
flic5 .AND. dti50.NE.0. )
THEN
1309 IF (
tic5(1).GE.0 )
THEN
1310 IF ( dti50 .LT. 0. )
THEN
1314 IF ( dttst .LE. 0.5*dti50 ) idact(18:18) =
'U'
1320 IF ( idact(18:18).NE.
' ' )
THEN
1329 call print_memcheck(memunit,
'memcheck_____:'//
' WW3_WAVE TIME LOOP 11a')
1333 IF (
fllev .AND. dtl0 .NE.0. )
THEN
1335 IF (
tlev(1) .GE. 0 )
THEN
1336 IF ( dtl0 .LT. 0. )
THEN
1340 IF ( dttst .LE. 0.5*dtl0 ) idact(5:5) =
'U'
1346 IF ( idact(5:5).NE.
' ' )
THEN
1359 CALL all_va_integral_print(imod,
"After FFLEV and DTL0", 1)
1364 call print_memcheck(memunit,
'memcheck_____:'//
' WW3_WAVE TIME LOOP 11b')
1388 CALL w3nmin (
mapsta, flag0 )
1410 call print_memcheck(memunit,
'memcheck_____:'//
' WW3_WAVE TIME LOOP 12')
1437 CALL all_va_integral_print(imod,
"VA before W3SRCE_IMP_PRE", 1)
1438 CALL all_field_integral_print(
vstot,
"VSTOT before W3SRCE_IMP_PRE")
1439 CALL all_field_integral_print(
vdtot,
"VDTOT before W3SRCE_IMP_PRE")
1453 call print_memcheck(memunit,
'memcheck_____:'//
' WW3_WAVE TIME LOOP 13')
1485 CALL init_get_isea(isea, jsea)
1501 delx=
hpfac(iy,ix)/ facx
1502 dely=
hqfac(iy,ix)/ facx
1505 reflec=
reflc(:,isea)
1506 reflec(4)=
berg(isea)*reflec(4)
1507 refled=
refld(:,isea)
1517 WRITE(740+
iaproc,*)
'NODE_SRCE_IMP_PRE : IX=', ix,
' JSEA=', jsea
1519 WRITE(740+
iaproc,*)
'IT/IX/IY/IMOD=', it, ix, iy, imod
1520 WRITE(740+
iaproc,*)
'ISEA/JSEA=', isea, jsea
1521 WRITE(740+
iaproc,*)
'Before sum(VA)=', sum(
va(:,jsea))
1526 vsiodummy, vdiodummy,
shavetot(jsea), &
1533 as(isea),
ust(isea), &
1537 reflec, refled, delx, dely, dela, &
1540 fcut(jsea), dtgpre, tauwx(jsea), tauwy(jsea), &
1544 tws(jsea),
phioc(jsea), tmp1, d50, psic, tmp2, &
1548 IF (.not. lsloc)
THEN
1549 vstot(:,jsea) = vsiodummy
1550 vdtot(:,jsea) = vdiodummy
1553 WRITE(740+
iaproc,*)
'After sum(VA)=', sum(
va(:,jsea))
1554 WRITE(740+
iaproc,*)
' sum(VSTOT)=', sum(
vstot(:,jsea))
1555 WRITE(740+
iaproc,*)
' sum(VDTOT)=', sum(
vdtot(:,jsea))
1567 CALL all_va_integral_print(imod,
"VA after W3SRCE_IMP_PRE", 1)
1568 CALL all_field_integral_print(
vstot,
"VSTOT after W3SRCE_IMP_PRE")
1569 CALL all_field_integral_print(
vdtot,
"VDTOT after W3SRCE_IMP_PRE")
1578 call print_memcheck(memunit,
'memcheck_____:'//
' WW3_WAVE TIME LOOP 14')
1601 WRITE(
ndse,*)
'Computing CFLs .... ',
nseal
1603 IF (
flogrd(9,3).AND. ugdtupdate )
THEN
1615 CALL init_get_isea(isea, jsea)
1621 IF (mod(isea,100).EQ.0)
WRITE(
ndse,*)
'COMPUTING CFL FOR NODE:',isea
1647 call print_memcheck(memunit,
'memcheck_____:'//
' WW3_WAVE TIME LOOP 15')
1656 indsort(jsea)=float(jsea)
1660 IF (
iaproc .EQ.
naperr )
WRITE(
ndse,*)
'Nodes requesting smallest timesteps:'
1661 IF (
iaproc .EQ.
naperr )
WRITE(
ndse,
'(A,10I10)')
'Nodes ',nint(indsort(1:10))
1662 IF (
iaproc .EQ.
naperr )
WRITE(
ndse,
'(A,10F10.2)')
'time steps ',dtcfl1(1:10)
1663 DO jsea = 1, min(
nseal,200)
1664 isea = nint(indsort(jsea))
1667 WRITE(995,*)
' IP dtmax_exp(ip) x-coord y-coord z-coord'
1669 WRITE(995,
'(I10,F10.2,3F10.4)') ix, dtcfl1(jsea),
xgrd(1,ix),
ygrd(2,ix),
zb(ix)
1685 ntloc = 1 + int( dtg/
dtcfli - 0.001 )
1688 WRITE(
ndse,
'(A,I4,A,I4)')
' SUBSECOND STEP:',isec1,
' out of ',
nitersec1
1692 facth = dtg / (
dth*real(ntloc))
1693 fack = dtg / real(ntloc)
1698 itloch = ( ntloc + 1 - mod(nint(dttest/dtg),2) ) / 2
1703 CALL all_va_integral_print(imod,
"Before intraspectral part 1", 1)
1717 CALL init_get_isea(isea, jsea)
1728 IF (
iobp(isea) .NE. 1) cycle
1732 IF (
mapsta(iy,ix) .EQ. 1 )
THEN
1733 depth = max(
dmin ,
dw(isea) )
1745 cg(:,isea),
wn(:,isea), depth, &
1757 cg(:,isea),
wn(:,isea), depth, &
1761 dcdx(:,iy,ixrel),
dcdy(:,iy,ixrel),
va(:,jsea))
1765 cg(:,isea),
wn(:,isea), depth, &
1769 dcdx(:,iy,ixrel),
dcdy(:,iy,ixrel),
va(:,jsea))
1773 cg(:,isea),
wn(:,isea), depth, &
1777 dcdx(:,iy,ixrel),
dcdy(:,iy,ixrel),
va(:,jsea), &
1794 call print_memcheck(memunit,
'memcheck_____:'//
' WW3_WAVE TIME LOOP 16')
1797 CALL all_va_integral_print(imod,
"Before spatial advection", 1)
1820 CALL pdlib_w3xypug ( ispec, facx, facx, dtg, vgx, vgy, ugdtupdate )
1830 CALL all_va_integral_print(imod,
"Before Block implicit", 1)
1833 CALL pdlib_w3xypug_block_implicit(imod, facx, facx, dtg, vgx, vgy, ugdtupdate )
1839 CALL pdlib_w3xypug_block_explicit(imod, facx, facx, dtg, vgx, vgy, ugdtupdate )
1848 IF (
nrqsg1 .GT. 0 )
THEN
1867 ELSE IF (.NOT.
lpdlib )
THEN
1868 CALL w3gath ( ispec, field )
1875 CALL w3psmc ( ispec, dtg, field )
1884 CALL w3xypug ( ispec, facx, facx, dtg, field, vgx, vgy, ugdtupdate )
1887 CALL w3xypug ( ispec, facx, facx, dtg, field, vgx, vgy, ugdtupdate )
1890 CALL w3xypug ( ispec, facx, facx, dtg, field, vgx, vgy, ugdtupdate )
1916 ELSE IF (.NOT.
lpdlib )
THEN
1917 CALL w3scat ( ispec,
mapsta, field )
1924 IF (
nrqsg1 .GT. 0 )
THEN
1925 ALLOCATE ( statco(mpi_status_size,
nrqsg1) )
1926 CALL mpi_waitall (
nrqsg1,
irqsg1(1,1), statco, ierr_mpi)
1927 CALL mpi_waitall (
nrqsg1,
irqsg1(1,2), statco, ierr_mpi)
1928 DEALLOCATE ( statco )
1931 call print_memcheck(memunit,
'memcheck_____:'//
' WW3_WAVE TIME LOOP 17')
1950 ispec = mod( iy-1,
naproc )
1951 jsea = 1 + (iy - ispec - 1)/
naproc
1962 IF(
iaproc .EQ. ispec+1 )
THEN
1981 ALLOCATE ( bacspec(
nspec) )
1993 ispec = mod( ix-1,
naproc )
1994 jsea = 1 + (ix - ispec - 1)/
naproc
2003 IF(
iaproc .EQ. ispec+1 )
THEN
2008 va(:,jsea) = bacspec
2013 DEALLOCATE ( bacspec )
2024 CALL all_va_integral_print(imod,
"After spatial advection", 1)
2033 DO itloc=itloch+1, ntloc
2042 CALL init_get_isea(isea, jsea)
2045 depth = max(
dmin ,
dw(isea) )
2053 IF (
iobp(isea) .NE. 1) cycle
2057 IF (
mapsta(iy,ix) .EQ. 1 )
THEN
2069 cg(:,isea),
wn(:,isea), depth, &
2080 cg(:,isea),
wn(:,isea), depth, &
2084 dcdx(:,iy,ixrel),
dcdy(:,iy,ixrel),
va(:,jsea))
2088 cg(:,isea),
wn(:,isea), depth, &
2092 dcdx(:,iy,ixrel),
dcdy(:,iy,ixrel),
va(:,jsea))
2096 cg(:,isea),
wn(:,isea), depth, &
2100 dcdx(:,iy,ixrel),
dcdy(:,iy,ixrel),
va(:,jsea), &
2117 CALL all_va_integral_print(imod,
"After intraspectral adv.", 1)
2124 ugdtupdate = .false.
2140 CALL all_vaold_integral_print(
"VAOLD before W3SRCE_IMP_POST", 1)
2141 CALL all_va_integral_print(imod,
"VA before W3SRCE_IMP_POST", 1)
2160 CALL init_get_isea(isea, jsea)
2174 delx=
hpfac(iy,ix)/ facx
2175 dely=
hqfac(iy,ix)/ facx
2181 reflec=
reflc(:,isea)
2182 reflec(4)=
berg(isea)*reflec(4)
2183 refled=
refld(:,isea)
2200 vsiodummy,vdiodummy,
shavetot(jsea), &
2207 as(isea),
ust(isea), &
2211 reflec, refled, delx, dely, dela, &
2214 fcut(jsea), dtg, tauwx(jsea), tauwy(jsea), &
2218 tws(jsea),
phioc(jsea), tmp1, d50, psic, tmp2, &
2225 vaolddummy,
va(:,jsea), &
2226 vsiodummy, vdiodummy, shavetotiodummy, &
2233 as(isea),
ust(isea), &
2237 reflec, refled, delx, dely, dela, &
2240 fcut(jsea), dtg, tauwx(jsea), tauwy(jsea), &
2244 tws(jsea),
phioc(jsea), tmp1, d50, psic,tmp2, &
2273 CALL all_vaold_integral_print(
"VAOLD after W3SRCE_IMP_PRE_POST", 1)
2274 CALL all_va_integral_print(imod,
"VA after W3SRCE_IMP_PRE_POST", 1)
2284 CALL all_va_integral_print(imod,
"After source terms", 1)
2295 IF (it.GT.0) dtg=dtgtemp
2310 dtg = dttst / real(nt-it)
2313 IF ( flact .AND. it.NE.nt .AND.
iaproc.EQ.
naplog )
THEN
2326 CALL all_va_integral_print(imod,
"end of time loop", 1)
2342 call print_memcheck(memunit,
'memcheck_____:'//
' WW3_WAVE END TIME LOOP')
2354 IF (
tofrst(1) .EQ. -1 )
THEN
2360 IF (
tdn(1) .EQ. -1 )
THEN
2367 flag_o = .NOT.skip_o .OR. ( skip_o .AND. dttst2.NE.0. )
2373 IF ( dttst.LE.0. .AND. dttst1.NE.0. .AND. flag_o )
THEN
2381 IF (
flout(1) )
THEN
2387 IF (
flout(7) )
THEN
2398 WRITE (
ndst,9042) local, flpart, floutg
2401 IF ( local .AND. flpart )
CALL w3cprt ( imod )
2402 IF ( local .AND. (floutg .OR. floutg2) )
then
2403 CALL w3outg (
va, flpfld, floutg, floutg2 )
2412 IF ( (floutg) .OR. (floutg2 .AND. sbsed) )
THEN
2414 IF (
nrqgo.NE.0 )
THEN
2417 CALL mpi_startall (
nrqgo,
irqgo , ierr_mpi )
2422 nrqmax = max( nrqmax ,
nrqgo )
2439 nrqmax = max( nrqmax ,
nrqgo2 )
2455 call print_memcheck(memunit,
'memcheck_____:'//
' WW3_WAVE AFTER TIME LOOP 1')
2462 nrqmax = max( nrqmax ,
nrqpo )
2475 CALL mpi_startall (
nrqrs,
irqrs , ierr_mpi )
2477 nrqmax = max( nrqmax ,
nrqrs )
2490 CALL mpi_startall (
nrqrs,
irqrs , ierr_mpi )
2492 nrqmax = max( nrqmax ,
nrqrs )
2507 nrqmax = max( nrqmax ,
nrqbp )
2521 nrqmax = max( nrqmax ,
nrqbp2 )
2532 IF ( nrqmax .NE. 0 )
ALLOCATE ( statio(mpi_status_size,nrqmax) )
2534 call print_memcheck(memunit,
'memcheck_____:'//
' WW3_WAVE AFTER TIME LOOP 2')
2544 IF (
flout(j) )
THEN
2554 IF ( dttst .EQ. 0. )
THEN
2562 IF ( flgmpi(1) )
CALL mpi_waitall (
nrqgo2,
irqgo2, statio, ierr_mpi )
2567 IF ( j .EQ. 1 )
THEN
2569 CALL w3iogo(
'WRITE',
nds(7), itest, imod &
2584 foutname =
'Field_done.' // idtime(1:4) &
2585 // idtime(6:7) // idtime(9:10) &
2586 // idtime(12:13) //
'.' //
filext(1:jj)
2590 OPEN( unit=ndsoflg,
file=foutname)
2595 ELSE IF ( j .EQ. 2 )
THEN
2605 CALL w3iopon (
'WRITE',
nds(8), itest, imod )
2607 CALL w3iopo (
'WRITE',
nds(8), itest, imod &
2615 ELSE IF ( j .EQ. 3 )
THEN
2620 ELSE IF ( j .EQ. 4 )
THEN
2623 ELSE IF ( j .EQ. 5 )
THEN
2631 ELSE IF ( j .EQ. 6 )
THEN
2634 ELSE IF ( j .EQ. 7 )
THEN
2638 IF (
dtout(7).NE.0)
THEN
2666 dttst =
dsec21( tout , tlst )
2667 flout(j) = dttst.GE.0.
2668 IF (
flout(j) )
THEN
2669 outid(2*j-1:2*j-1) =
'X'
2671 IF ( (
dtout(7).NE.0) .AND. &
2676 outid(2*j-1:2*j-1) =
'L'
2682 IF (
flout(j) )
THEN
2683 IF (
tofrst(1).EQ.-1 )
THEN
2687 IF ( dttst.GT.0.)
THEN
2697 call print_memcheck(memunit,
'memcheck_____:'//
' WW3_WAVE AFTER TIME LOOP 3')
2701 IF (
flout(j) )
THEN
2707 IF ( dttst .EQ. 0. )
THEN
2713 dttst =
dsec21( tout , tlst )
2714 flout(j) = dttst.GE.0.
2715 IF (
flout(j) )
THEN
2716 outid(2*j-1:2*j-1) =
'X'
2718 IF ( (
dtout(7).NE.0) .AND. &
2723 outid(2*j-1:2*j-1) =
'L'
2729 IF (
flout(j) )
THEN
2730 IF (
tofrst(1).EQ.-1 )
THEN
2734 IF ( dttst.GT.0.)
THEN
2742 call print_memcheck(memunit,
'memcheck_____:'//
' WW3_WAVE AFTER TIME LOOP 3')
2745 IF ( flgmpi(0) )
CALL mpi_waitall (
nrqgo,
irqgo , statio, ierr_mpi )
2746 IF ( flgmpi(2) )
CALL mpi_waitall (
nrqpo,
irqpo1, statio, ierr_mpi )
2747 IF ( flgmpi(4) )
CALL mpi_waitall (
nrqrs,
irqrs , statio, ierr_mpi )
2748 IF ( flgmpi(8) )
CALL mpi_waitall (
nrqrs,
irqrs , statio, ierr_mpi )
2749 IF ( flgmpi(5) )
CALL mpi_waitall (
nrqbp,
irqbp1, statio, ierr_mpi )
2750 IF ( nrqmax .NE. 0 )
DEALLOCATE ( statio )
2760 call print_memcheck(memunit,
'memcheck_____:'//
' WW3_WAVE AFTER TIME LOOP 4')
2769 IF ( dttst .EQ. 0. ) idact(7:7) =
'X'
2773 IF ( dttst .EQ. 0. ) idact(3:3) =
'X'
2777 IF ( dttst .EQ. 0. ) idact(9:9) =
'X'
2781 IF ( dttst .EQ. 0. ) idact(11:11) =
'X'
2783 IF (
tdn(1) .GT. 0 )
THEN
2785 IF ( dttst .EQ. 0. ) idact(21:21) =
'X'
2804 IF ( dttst .EQ. 0. )
EXIT
2809 call print_memcheck(memunit,
'memcheck_____:'//
' WW3_WAVE AFTER TIME LOOP 5')
2814 WRITE (
screen,951) sttime
2820 DEALLOCATE(tauwx, tauwy)
2822 call print_memcheck(memunit,
'memcheck_____:'//
' WW3_WAVE END W3WAVE')
2828 900
FORMAT (4x,i6,
'|',i6,
'| ', a19 ,
' | ',a,
' | ',a,
' |')
2829 901
FORMAT (4x,i6,
'|',i6,
'| ',11x,a8,
' | ',a,
' | ',a,
' |')
2830 902
FORMAT (2x,
'--------+------+---------------------+' &
2831 ,
'-----------------------+------------------+')
2834 920
FORMAT (
' Updating k and Cg from ice param. 1,2,3,4.'/)
2836 950
FORMAT (
' WAVEWATCH III calculating for ',a,
' at ',a)
2837 951
FORMAT (
' WAVEWATCH III reached the end of a computation', &
2839 1000
FORMAT (/
' *** WAVEWATCH III ERROR IN W3WAVE :'/ &
2840 ' ENDING TIME BEFORE STARTING TIME '/)
2841 1001
FORMAT (/
' *** WAVEWATCH III ERROR IN W3WAVE :'/ &
2842 ' NEW WATER LEVEL BEFORE OLD WATER LEVEL '/)
2843 1002
FORMAT (/
' *** WAVEWATCH III ERROR IN W3WAVE :'/ &
2844 ' ILLEGAL CURRENT INTERVAL '/)
2845 1003
FORMAT (/
' *** WAVEWATCH III ERROR IN W3WAVE :'/ &
2846 ' ILLEGAL WIND INTERVAL '/)
2847 1004
FORMAT (/
' *** WAVEWATCH III ERROR IN W3WAVE :'/ &
2848 ' NEW ICE FIELD BEFORE OLD ICE FIELD '/)
2849 1005
FORMAT (/
' *** WAVEWATCH III ERROR IN W3WAVE :'/ &
2850 ' NEW IC1 FIELD BEFORE OLD IC1 FIELD '/)
2851 1007
FORMAT (/
' *** WAVEWATCH III ERROR IN W3WAVE :'/ &
2852 ' NEW ATM MOMENTUM BEFORE OLD ATM MOMENTUM '/)
2853 1008
FORMAT (/
' *** WAVEWATCH III ERROR IN W3WAVE :'/ &
2854 ' NEW AIR DENSITY BEFORE OLD AIR DENSITY '/)
2856 1006
FORMAT (/
' *** WAVEWATCH III ERROR IN W3WAVE :'/ &
2857 ' NEW IC5 FIELD BEFORE OLD IC5 FIELD '/)
2859 1030
FORMAT (/
' *** WAVEWATCH III WARING IN W3WAVE :'/ &
2860 ' AT LEAST ONE PROCESSOR HAS 0 ACTIVE POINTS', &
2863 1040
FORMAT (/
' *** WAVEWATCH III ERROR IN W3WAVE :'/ &
2864 ' EXPERIMENTAL FEATURE !/REFRX NOT FULLY IMPLEMENTED.'/)
2869 '============================================================', &
2870 '===================='/ &
2871 ' TEST W3WAVE : RUN MODEL',i3,
' FILEXT [',a, &
2872 '] UP TO ',i8.8,i7.6 / &
2873 '====================', &
2874 '============================================================')
2875 9010
FORMAT (
' TEST W3WAVE : DT INT. =',f12.1,
' FLZERO = ',l1)
2876 9011
FORMAT (
' TEST W3WAVE : DT LEV. =',f12.1)
2877 9012
FORMAT (
' TEST W3WAVE : DT CUR. =',f12.1/ &
2880 9013
FORMAT (
' TEST W3WAVE : DT WIND =',f12.1/ &
2883 9014
FORMAT (
' TEST W3WAVE : DT ICE =',f12.1)
2884 9015
FORMAT (
' TEST W3WAVE : DT IC1 =',f12.1)
2885 9016
FORMAT (
' TEST W3WAVE : DT IC5 =',f12.1)
2886 9017
FORMAT (
' TEST W3WAVE : DT TAU =',f12.1)
2887 9018
FORMAT (
' TEST W3WAVE : DT RHO =',f12.1)
2888 9020
FORMAT (
' TEST W3WAVE : IT0, NT, DTG :',2i4,f8.1)
2889 9021
FORMAT (
' TEST W3WAVE : ITIME etc',i6,i4,i10.8,i7.6,1x,2l1, &
2891 9022
FORMAT (
' TEST W3WAVE : SKIP TO 400 IN 3.5')
2892 9023
FORMAT (
' TEST W3WAVE : SKIP TO 380 IN 3.5')
2893 9030
FORMAT (
' TEST W3WAVE : END OF COMPUTATION LOOP')
2894 9040
FORMAT (
' TEST W3WAVE : CHECKING FOR OUTPUT'/ &
2895 ' TOFRST :',i9.8,i7.6/ &
2896 ' TND :',i9.8,i7.6/ &
2897 ' DTTST[1], FLAG_O :',2f8.1,l4)
2898 9041
FORMAT (
' TEST W3WAVE : PERFORMING OUTPUT')
2899 9042
FORMAT (
' TEST W3WAVE : OUTPUT COMPUTATION FLAGS: ',3l2)
2902 9043
FORMAT (
' TEST W3WAVE : TYPE, NRQ, NRQMAX, NA : ',a2,3i6)
2905 9044
FORMAT (
' TEST W3WAVE : END OF OUTPUT')
References pdlib_w3profsmd::all_field_integral_print(), pdlib_w3profsmd::all_va_integral_print(), pdlib_w3profsmd::all_vaold_integral_print(), w3adatmd::alpha, w3gdatmd::angarc, pdlib_w3profsmd::apply_boundary_condition_va(), w3gdatmd::arctc, w3adatmd::as, w3wdatmd::asf, pdlib_w3profsmd::aspar_diag_all, pdlib_w3profsmd::aspar_jac, pdlib_w3profsmd::b_jac, w3adatmd::bedforms, w3wdatmd::berg, w3sic3md::calledic3table, w3adatmd::cflkmax, w3adatmd::cflthmax, w3adatmd::cflxymax, w3adatmd::cg, w3adatmd::charn, w3gdatmd::clats, w3gdatmd::clgtype, w3oacpmd::cplt0, w3gdatmd::cthg0s, w3adatmd::cx, w3adatmd::cy, constants::dair, w3adatmd::dcdx, w3adatmd::dcdy, w3adatmd::dcxdx, w3adatmd::dcxdy, w3adatmd::dcydx, w3adatmd::dcydy, w3adatmd::dddx, w3adatmd::dddy, constants::debug_node, constants::dera, w3adatmd::dhdx, w3adatmd::dhdy, w3adatmd::dhlmt, w3gdatmd::dmin, pdlib_field_vec::do_output_exchanges(), w3timemd::dsec21(), w3gdatmd::dtcfli, w3adatmd::dtdyn, w3gdatmd::dth, w3gdatmd::dtmax, w3odatmd::dtout, w3adatmd::dw, w3servmd::extcde(), w3adatmd::fcut, file(), w3gdatmd::filext, w3gdatmd::flagll, w3gdatmd::flagst, w3odatmd::flbpi, w3gdatmd::flck, w3adatmd::flcold, w3gdatmd::flcth, w3idatmd::flcur, w3gdatmd::flcx, w3gdatmd::flcy, w3gdatmd::fldry, w3idatmd::flic1, w3idatmd::flic2, w3idatmd::flic3, w3idatmd::flic4, w3idatmd::flic5, w3idatmd::flice, w3adatmd::fliwnd, w3idatmd::fllev, w3odatmd::flogr2, w3odatmd::flogrd, w3odatmd::flout, w3idatmd::flrhoa, w3gdatmd::flsou, w3idatmd::fltaua, w3idatmd::flwind, w3wdatmd::fpis, w3gdatmd::fsfreqshift, w3gdatmd::fsrefraction, w3gdatmd::fssource, w3gdatmd::fstotalexp, w3gdatmd::fstotalimp, w3idatmd::ga0, w3idatmd::gan, w3idatmd::gd0, w3idatmd::gdn, w3gdatmd::gtype, w3gdatmd::hpfac, w3gdatmd::hqfac, w3adatmd::iadata, w3adatmd::iappro, w3odatmd::iaproc, w3gdatmd::ic3pars, w3sic3md::ic3table_cheng(), w3wdatmd::ice, w3wdatmd::icedmax, w3wdatmd::icef, w3wdatmd::iceh, w3gdatmd::iclbac, w3oacpmd::id_oasis_time, w3adatmd::idlast, w3gdatmd::igrid, w3idatmd::iidata, include(), w3idatmd::inflags1, w3parall::init_get_isea(), w3gdatmd::iobp, w3gdatmd::iobp_loc, w3odatmd::ioutp, w3adatmd::ipass, yownodepool::iplg, w3odatmd::irqbp1, w3odatmd::irqbp2, w3odatmd::irqgo, w3odatmd::irqgo2, w3odatmd::irqpo1, w3odatmd::irqrs, w3adatmd::irqsg1, w3adatmd::itime, w3wdatmd::iwdata, constants::lpdlib, w3parall::lsloc, w3gdatmd::mapfs, w3gdatmd::mapsf, w3gdatmd::mapsta, w3adatmd::mpi_comm_wave, w3odatmd::napbpt, w3odatmd::naperr, w3odatmd::napfld, w3odatmd::naplog, w3odatmd::napout, w3odatmd::nappnt, w3odatmd::naproc, w3odatmd::naprst, w3gdatmd::nbac, w3gdatmd::nbgl, w3gdatmd::ncel, w3odatmd::nds, w3odatmd::ndse, w3odatmd::ndso, w3odatmd::ndst, w3gdatmd::nglo, w3gdatmd::nitersec1, w3gdatmd::nk, w3odatmd::noge, w3odatmd::notype, yownodepool::np, yownodepool::npa, w3odatmd::nrqbp, w3odatmd::nrqbp2, w3odatmd::nrqgo, w3odatmd::nrqgo2, w3odatmd::nrqpo, w3odatmd::nrqrs, w3adatmd::nrqsg1, w3gdatmd::nseal, w3gdatmd::nth, w3gdatmd::nx, w3gdatmd::ny, w3parall::pdlib_nseal, w3parall::pdlib_nsealm, pdlib_w3profsmd::pdlib_w3xypug(), pdlib_w3profsmd::pdlib_w3xypug_block_explicit(), pdlib_w3profsmd::pdlib_w3xypug_block_implicit(), w3adatmd::phiaw, w3adatmd::phibbl, w3adatmd::phice, w3adatmd::phioc, w3servmd::print_memcheck(), w3parall::print_my_time(), constants::radius, w3gdatmd::reflc, w3gdatmd::refld, w3wdatmd::rhoair, w3gdatmd::rlgtype, w3gdatmd::rstype, w3odatmd::screen, w3gdatmd::sed_d50, w3gdatmd::sed_psic, w3wdatmd::shavetot, w3gdatmd::sig, w3psmcmd::smcdcxy(), w3psmcmd::smcdhxy(), w3gdatmd::smctype, w3agcmmd::snd_fields_to_atmos(), w3igcmmd::snd_fields_to_ice(), w3ogcmmd::snd_fields_to_ocean(), w3gdatmd::spcbac, constants::srce_direct, constants::srce_imp_post, constants::srce_imp_pre, w3servmd::ssort1(), w3timemd::stme21(), w3odatmd::stop, w3servmd::strace(), w3gdatmd::sx, w3gdatmd::sy, w3adatmd::taua, w3adatmd::tauadir, w3adatmd::taubbl, w3adatmd::tauice, w3adatmd::tauocx, w3adatmd::tauocy, w3adatmd::tauox, w3adatmd::tauoy, w3adatmd::tauwix, w3adatmd::tauwiy, w3adatmd::tauwnx, w3adatmd::tauwny, w3odatmd::tbpi0, w3odatmd::tbpin, w3idatmd::tc0, w3idatmd::tcn, w3idatmd::tdn, w3idatmd::tg0, w3idatmd::tgn, w3idatmd::ti1, w3idatmd::ti5, w3wdatmd::tic1, w3wdatmd::tic5, w3wdatmd::tice, w3timemd::tick21(), w3wdatmd::time, w3wdatmd::time00, w3wdatmd::timeend, w3idatmd::tin, w3wdatmd::tlev, w3idatmd::tln, w3odatmd::tofrst, w3odatmd::tolast, w3odatmd::tonext, w3odatmd::tosnl5, constants::tpiinv, w3idatmd::tr0, w3idatmd::trn, w3gdatmd::trnx, w3gdatmd::trny, w3idatmd::tu0, w3idatmd::tun, w3idatmd::tw0, w3idatmd::twn, w3adatmd::tws, w3adatmd::u10, w3adatmd::u10d, w3triamd::ug_gradients(), constants::undef, w3gdatmd::ungtype, w3uostmd::uost_setgrid(), w3wdatmd::ust, w3wdatmd::ustdir, w3wdatmd::va, w3wdatmd::vaold, w3wdatmd::vdtot, w3wdatmd::vstot, w3servmd::w3acturn(), w3profsmd::w3cflug(), w3pro3md::w3cflxy(), w3iosfmd::w3cprt(), w3updtmd::w3dzxy(), w3gath(), w3psmcmd::w3gathsmc(), w3sic3md::w3ic3wncg_cheng(), w3sic3md::w3ic3wncg_v1(), w3iobcmd::w3iobc(), w3iogomd::w3iogo(), w3iopomd::w3iope(), w3iopomd::w3iopo(), w3iopomd::w3iopon(), w3iorsmd::w3iors(), w3iosfmd::w3iosf(), w3iotrmd::w3iotr(), w3psmcmd::w3krtn(), w3pro1md::w3ktp1(), w3pro2md::w3ktp2(), w3pro3md::w3ktp3(), w3pro1md::w3map1(), w3pro2md::w3map2(), w3pro3md::w3map3(), w3pro3md::w3mapt(), w3nmin(), w3iogomd::w3outg(), w3psmcmd::w3psmc(), w3scat(), w3psmcmd::w3scatsmc(), w3adatmd::w3seta(), w3gdatmd::w3setg(), w3idatmd::w3seti(), w3odatmd::w3seto(), w3wdatmd::w3setw(), w3srcemd::w3srce(), w3updtmd::w3ubpt(), w3updtmd::w3ucur(), w3updtmd::w3uic1(), w3updtmd::w3uic5(), w3updtmd::w3uice(), w3updtmd::w3uini(), w3updtmd::w3ulev(), w3updtmd::w3urho(), w3updtmd::w3utau(), w3updtmd::w3utrn(), w3updtmd::w3uwnd(), w3pro1md::w3xyp1(), w3pro2md::w3xyp2(), w3pro3md::w3xyp3(), w3profsmd::w3xypug(), w3wavset::wave_setup_computation(), w3adatmd::whitecap, w3adatmd::wn, w3adatmd::wnmean, w3servmd::wwtime(), w3gdatmd::xgrd, w3gdatmd::ygrd, and w3gdatmd::zb.
Referenced by wmwavemd::wmwave().
real, dimension(:,:), allocatable aspar_diag_all
integer, dimension(:), pointer tbpi0
integer, parameter srce_imp_pre
srce_imp_pre
real function dsec21(TIME1, TIME2)
subroutine ug_gradients(PARAM, DIFFX, DIFFY)
Calculate gradients at a point via its connection.
real, dimension(:), pointer phice
cmake src_list cmake include(${CMAKE_CURRENT_SOURCE_DIR}/cmake/check_switches.cmake) check_switches("$
integer, pointer nitersec1
subroutine all_field_integral_print(FIELD, string)
real, dimension(:), pointer iceh
integer, dimension(:), pointer twn
Module used for coupling applications between atmospheric model and WW3 with OASIS3-MCT.
logical, dimension(:), pointer inflags1
real, dimension(:), pointer charn
real, dimension(:), pointer dtdyn
subroutine w3scatsmc(ISPEC, MAPSTA, FIELD)
SMC version of W3GATH.
real, dimension(:), pointer fpis
Reads triangle and unstructured grid information.
real, dimension(:), pointer as
logical, dimension(:), pointer shavetot
integer, dimension(:), pointer tosnl5
double precision, dimension(:,:), pointer ygrd
Parmeterization of the unresolved obstacles.
Define data structures to set up wave model auxiliary data for several models simultaneously.
subroutine w3cprt(IMOD)
Partitioning of spectra into fields for all grid points that are locally stored.
real, dimension(:,:), pointer trnx
real, dimension(:,:,:), pointer dcdx
real, parameter dair
DAIR Density of air (kg/m3).
real, dimension(:), pointer zb
real, parameter dera
DERA Conversion factor from degrees to radians.
subroutine w3ktp3(ISEA, FACTH, FACK, CTHG0, CG, WN, DW, DDDX, DDDY, CX, CY, DCXDX, DCXDY, DCYDX, DCYDY, DCDX, DCDY, VA, CFLTHMAX, CFLKMAX)
Propagation in spectral space.
logical, dimension(:), pointer flagst
integer, parameter ungtype
subroutine w3ucur(FLFRST)
Interpolate the current field to the present time.
integer, parameter srce_direct
srce_direct
Implicit solution of wave setup problem following Dingemans for structured and unstructured grids.
Define data structures to set up wave model dynamic data for several models simultaneously.
subroutine w3urho(FLFRST)
Interpolate air density field to the given time.
real, dimension(:), pointer dhdy
real, dimension(:,:), pointer dhlmt
real, dimension(:,:), pointer tauice
subroutine w3uini(A)
Initialize the wave field with fetch-limited spectra before the actual calculation start.
real, dimension(:), pointer sed_d50
real, dimension(:), pointer cflxymax
real, dimension(:,:), pointer cg
real, dimension(:), pointer tws
subroutine smcdhxy
Calculates water-depth gradient for refraction.
real, dimension(:), pointer fcut
integer, parameter rlgtype
logical, dimension(:,:), pointer flogr2
real, dimension(:,:), allocatable aspar_jac
subroutine w3ktp2(ISEA, FACTH, FACK, CTHG0, CG, WN, DEPTH, DDDX, DDDY, CX, CY, DCXDX, DCXDY, DCYDX, DCYDY, DCDX, DCDY, VA)
Propagation in spectral space.
integer, dimension(:), pointer bispl
real, dimension(:), pointer dw
subroutine w3iors(INXOUT, NDSR, DUMFPI, IMOD, FLRSTRT)
Reads/writes restart files.
real, dimension(:), pointer u10d
real, dimension(:,:), pointer reflc
real, dimension(:), pointer dtout
real, dimension(:), pointer sig
real, dimension(:), pointer icef
double precision, dimension(:,:), pointer xgrd
subroutine w3xyp3(ISP, DTG, MAPSTA, MAPFS, VQ, VGX, VGY)
Propagation in phyiscal space for a given spectral component.
integer, dimension(:), pointer ti5
logical, pointer fsrefraction
Bundles all input updating routines for WAVEWATCH III.
integer, dimension(:), pointer time
real, dimension(:), pointer tauocy
integer, dimension(:), pointer irqpo1
integer, parameter srce_imp_post
srce_imp_post
subroutine w3ubpt
Update spectra at the active boundary points.
integer, dimension(:), allocatable, public iplg
Node local to global mapping.
subroutine w3srce(srce_call, IT, ISEA, JSEA, IX, IY, IMOD, SPECOLD, SPEC, VSIO, VDIO, SHAVEIO, ALPHA, WN1, CG1, CLATSL, D_INP, U10ABS, U10DIR, ifdef W3_FLX5
Calculate and integrate source terms for a single grid point.
integer, dimension(:), pointer tbpin
integer, dimension(:), pointer tg0
integer, dimension(:), pointer iclbac
integer, dimension(:), pointer iappro
logical, pointer fssource
integer *2, dimension(:), pointer iobp_loc
integer, public npa
number of ghost + resident nodes this partition holds
integer, save calledic3table
real, dimension(:,:), pointer va
subroutine w3iobc(INXOUT, NDSB, TIME1, TIME2, IOTST, IMOD)
Write/read boundary conditions file(s).
real, dimension(:), pointer tauwix
integer, dimension(:), pointer tu0
real, dimension(:,:), pointer spcbac
subroutine w3ulev(A, VA)
Update the water level.
real, dimension(:,:), pointer dcydx
Bundles routines for third order propagation scheme in single module.
real, dimension(:,:), pointer hqfac
integer, dimension(:), pointer tlev
subroutine, public snd_fields_to_ice()
Send coupling fields to ice model.
subroutine w3setg(IMOD, NDSE, NDST)
Bundles routines for first order propagation scheme in single module.
Bundles routines for third order porpagation scheme in single module.
subroutine w3acturn(NDirc, NFreq, Alpha, Spectr)
real, dimension(:), pointer tauwiy
real, dimension(:,:), pointer taubbl
character(len=8), dimension(max_timers), save status
real, dimension(:), pointer berg
subroutine w3seta(IMOD, NDSE, NDST)
Select one of the WAVEWATCH III grids / models.
real, dimension(:), pointer tauwny
integer, dimension(:,:), pointer mapfs
integer, dimension(:,:), pointer tonext
real, dimension(:,:,:), pointer dcdy
subroutine, public ic3table_cheng(ICE2, ICE3, ICE4)
real, dimension(:), pointer phiaw
subroutine w3iotr(NDSINP, NDSOUT, A, IMOD)
Perform output of spectral information along provided tracks.
real, dimension(:,:), pointer vdtot
real, dimension(:), pointer phibbl
integer, parameter mpibuf
integer, dimension(:,:), pointer refld
real, dimension(:), pointer sed_psic
subroutine, public snd_fields_to_atmos()
Send coupling fields to atmospheric model.
Has data that belong to nodes.
logical lpdlib
LPDLIB Logical for using the PDLIB library.
subroutine do_output_exchanges(IMOD)
real, dimension(:), pointer cflthmax
subroutine w3cflxy(ISEA, DTG, MAPSTA, MAPFS, CFLXYMAX, VGX, VGY)
Computes the maximum CFL number for spatial advection.
integer, parameter clgtype
real, dimension(:,:), pointer vstot
real, dimension(:,:), pointer dcxdx
logical, dimension(:,:), pointer flogrd
subroutine all_vaold_integral_print(string, choice)
real, dimension(:,:), pointer bedforms
subroutine tick21(TIME, DTIME)
integer, dimension(:), pointer tic1
subroutine pdlib_w3xypug(ISP, FACX, FACY, DTG, VGX, VGY, LCALC)
subroutine w3setw(IMOD, NDSE, NDST)
Select one of the WAVEWATCH III grids / models.
integer, dimension(nogrp) noge
integer, dimension(:,:), pointer irqsg1
real, parameter tpiinv
TPIINV Inverse of 2*Pi.
integer, dimension(:), pointer irqbp1
subroutine w3seto(IMOD, NDSERR, NDSTST)
subroutine stme21(TIME, DTME21)
subroutine print_memcheck(iun, msg)
Write memory statistics if requested.
real, dimension(:,:), pointer gstore
real, dimension(:), pointer tauwnx
real, dimension(:,:), pointer alpha
subroutine ssort1(X, Y, N, KFLAG)
subroutine w3xypug(ISP, FACX, FACY, DTG, VQ, VGX, VGY, LCALC)
real, dimension(:,:), pointer dddy
real, dimension(:), pointer clats
integer, dimension(:), pointer irqgo2
integer, dimension(:), pointer nds
subroutine w3uice(VA)
Update ice map in the wave model.
subroutine pdlib_w3xypug_block_explicit(IMOD, FACX, FACY, DTG, VGX, VGY, LCALC)
real, dimension(:), pointer cy
Source term integration routine.
subroutine w3map1(MAPSTA)
Generate 'map' arrays for the first order upstream scheme.
Processing of boundary data output.
real, dimension(:), pointer taua
real, dimension(:), pointer angarc
integer, dimension(:,:), pointer tolast
subroutine, public w3ic3wncg_cheng(WN_R, WN_I, CG, ICE1, ICE2, ICE3, ICE4, DPT)
subroutine pdlib_w3xypug_block_implicit(IMOD, FACX, FACY, DTG, VGX, VGY, LCALC)
integer, dimension(:,:), pointer mapsf
integer, dimension(:), pointer ti1
subroutine print_my_time(string)
Print timings.
integer, public np
number of nodes, local
real, dimension(:), pointer wnmean
real, dimension(:), pointer cflkmax
integer, dimension(:), pointer irqbp2
real, dimension(:,:), pointer sstore
subroutine w3iogo(INXOUT, NDSOG, IOTST, IMOD ifdef W3_ASCII
Read/write gridded output.
subroutine w3seti(IMOD, NDSE, NDST)
Select one of the WAVEWATCH III grids / models.
subroutine w3uic5(FLFRST)
Update ice floe mean and max diameters in the wave model.
integer, parameter smctype
real, dimension(:), pointer cthg0s
real, dimension(:), pointer phioc
file(STRINGS ${CMAKE_BINARY_DIR}/switch switch_strings) separate_arguments(switches UNIX_COMMAND $
Reading/writing of model definition file.
real, parameter radius
RADIUS Radius of the earth (m).
subroutine w3iopon(INXOUT, NDSOP, IOTST, IMOD)
Read or write the netCDF point output file, depending on the value of the first parameter.
subroutine wave_setup_computation
General driver.
real, dimension(:,:), pointer wn
subroutine w3dzxy(ZZ, ZUNIT, DZZDX, DZZDY)
Calculate derivatives of a field.
real, dimension(:), pointer u10
logical, pointer fstotalexp
Module used for coupling applications between ice model and WW3 with OASIS3-MCT.
integer, dimension(:), pointer irqrs
integer, dimension(:), pointer timeend
real, dimension(:,:), pointer dcydy
integer, dimension(:), pointer tdn
logical, dimension(:), pointer flout
Gridded output of mean wave parameters.
subroutine strace(IENT, SNAME)
real, dimension(:), pointer icedmax
integer, public id_oasis_time
Define data structures to set up wave model input data for several models simultaneously.
real, dimension(:), pointer ice
subroutine w3utau(FLFRST)
Interpolate atmosphere momentum fields to the given time.
Floe-size dependant scattering of waves in the marginal ice zone.
real, dimension(:,:), pointer whitecap
subroutine, public uost_setgrid(IGRID)
Sets the current grid in the sourceterm object.
subroutine w3xyp1(ISP, DTG, MAPSTA, FIELD, VGX, VGY)
Propagation in physical space for a given spectral component.
subroutine w3map3
Generate 'map' arrays for the ULTIMATE QUICKEST scheme.
integer, dimension(:), pointer tc0
subroutine, public snd_fields_to_ocean()
real, dimension(:), pointer tauox
integer, dimension(:), pointer tw0
integer, dimension(:), pointer tin
Spherical Multiple-Cell (SMC) grid routines.
real, dimension(:,:), pointer hpfac
subroutine init_get_jsea_isproc(ISEA, JSEA, ISPROC)
Set JSEA for all schemes.
subroutine w3xyp2(ISP, DTG, MAPSTA, MAPFS, VQ, VGX, VGY)
Propagation in physical space for a given spectral component.
integer, dimension(:), pointer time00
integer *2, dimension(:), pointer iobp
integer, dimension(:), pointer tofrst
real, dimension(:), pointer ust
integer, pointer mpi_comm_wave
integer, parameter debug_node
DEBUG_NODE Node number used for debugging.
Define some much-used constants for global use (all defined as PARAMETER).
real, dimension(:), pointer tauoy
subroutine w3cflug(ISEA, NKCFL, FACX, FACY, DT, MAPFS, CFLXYMAX, VGX, VGY)
real, dimension(:), pointer ic3pars
subroutine w3gathsmc(ISPEC, FIELD)
SMC version of W3GATH.
real, dimension(:,:), pointer dcxdy
subroutine w3uic1(FLFRST)
Update ice thickness in the wave model.
real, dimension(:), pointer tauadir
subroutine extcde(IEXIT, UNIT, MSG, FILE, LINE, COMM)
real, dimension(:), pointer rhoair
integer, dimension(:), pointer bstat
real, dimension(:), pointer ustdir
subroutine w3outg(A, FLPART, FLOUTG, FLOUTG2)
Fill necessary arrays with gridded data for output.
subroutine all_va_integral_print(IMOD, string, choice)
subroutine w3iosf(NDSPT, IMOD)
Write partitioned spectral data to file.
real, dimension(:,:), pointer dddx
subroutine, public w3ic3wncg_v1(WN_R, WN_I, CG, ICE1, ICE2, ICE3, ICE4, DPT)
integer, dimension(:,:), pointer irqsg2
integer, dimension(:), pointer tln
real, dimension(:), pointer cx
integer, dimension(:), pointer tic5
integer, dimension(:), pointer trn
subroutine w3psmc(ISP, DTG, VQ)
Propagation in phyiscal space for a given spectral component.
logical, pointer fsfreqshift
real undef
UNDEF Value for undefined variable in output.
double precision, parameter fac
Parallel routines for implicit solver.
subroutine w3iopo(INXOUT, NDSOP, IOTST, IMOD ifdef W3_ASCII
Read or write point output.
real, dimension(:,:), allocatable b_jac
integer, dimension(:), pointer irqgo
logical, pointer fstotalimp
integer, dimension(:), pointer tun
integer, dimension(:), pointer tr0
subroutine w3mapt
Generate 'map' arrays for the ULTIMATE QUICKEST scheme to combine GSE alleviation with obstructions.
subroutine w3iope(A)
Extract point output data and store in output COMMONs.
subroutine w3utrn(TRNX, TRNY)
Update cell boundary transparencies for general use in propagation routines.
integer, dimension(:,:), pointer mapsta
subroutine w3uwnd(FLFRST, VGX, VGY)
Interpolate wind fields to the given time.
subroutine w3map2
Generate 'map' arrays for the ULTIMATE QUICKEST scheme.
real, dimension(:,:), pointer vaold
Read/write restart files.
integer, dimension(:), pointer tice
real, dimension(:,:), pointer trny
I/O and computational routines for the wave-field separation output.
real, dimension(:), pointer tauocx
integer, dimension(:), pointer tgn
subroutine apply_boundary_condition_va
subroutine init_get_isea(ISEA, JSEA)
Set ISEA for all schemes.
integer, dimension(:), pointer tcn
real, dimension(:), pointer dhdx
subroutine w3krtn(ISEA, FACTH, FACK, CTHG0, CG, WN, DEPTH, DDDX, DDDY, ALFLMT, CX, CY, DCXDX, DCXDY, DCYDX, DCYDY, DCDX, DCDY, VA)
Refraction and great-circle turning by spectral rotation.
subroutine w3ktp1(ISEA, FACTH, FACK, CTHG0, CG, WN, DEPTH, DDDX, DDDY, CX, CY, DCXDX, DCXDY, DCYDX, DCYDY, DCDX, DCDY, VA)
Propagation in spectral space.
real, dimension(:), pointer asf
character(len=13), pointer filext
subroutine smcdcxy
Calculates current velocity gradient for refraction.