WAVEWATCH III  beta 0.0.1
ww3_grib.F90 File Reference

Contains program W3GRIB for post-processing grid output. More...

Go to the source code of this file.

Functions/Subroutines

program w3grib
 Post-processing of grid output. More...
 
subroutine w3exgb (NX, NY, NSEA)
 Perform actual GRIB output. More...
 

Detailed Description

Contains program W3GRIB for post-processing grid output.

Author
H. L. Tolman
A. Chawla
J.-H. Alves
Date
22-Mar-2021

Definition in file ww3_grib.F90.

Function/Subroutine Documentation

◆ w3exgb()

subroutine w3grib::w3exgb ( integer, intent(in)  NX,
integer, intent(in)  NY,
integer, intent(in)  NSEA 
)

Perform actual GRIB output.

Parameters
[in]NXX array dimension
[in]NYY array dimension
[in]NSEASeapoint array dimension
Author
H. L. Tolman
A. Chawla
Date
22-Mar-2021

Definition at line 915 of file ww3_grib.F90.

915  !/
916  !/ +-----------------------------------+
917  !/ | WAVEWATCH III NOAA/NCEP |
918  !/ | H. L. Tolman |
919  !/ | A. Chawla |
920  !/ | FORTRAN 90 |
921  !/ | Last update : 22-Mar-2021 |
922  !/ +-----------------------------------+
923  !/
924  !/ 10-Jun-1999 : Final FORTRAN 77 ( version 1.18 )
925  !/ 24-Jan-2000 : Upgrade to FORTRAN 90 ( version 2.00 )
926  !/ Massive changes to logistics.
927  !/ 29-Apr-2002 : Adding output fields 17-18. ( version 2.20 )
928  !/ 24-Dec-2004 : Multiple grid version. ( version 3.06 )
929  !/ 18-May-2007 : Update GRIB1 for partitioning. ( version 3.11 )
930  !/ 16-Jul-2007 : Adding GRIB2 capability ( version 3.11 )
931  !/ (A. Chawla)
932  !/ 22-Mar-2021 : New coupling fields output ( version 7.13 )
933  !/
934  ! 1. Purpose :
935  !
936  ! Perform actual GRIB output.
937  !
938  ! 3. Parameters :
939  !
940  ! Parameter list
941  ! ----------------------------------------------------------------
942  ! NX, NY, NSEA
943  ! Int. I Array dimensions.
944  ! ----------------------------------------------------------------
945  !
946  ! Internal parameters
947  ! ----------------------------------------------------------------
948  ! X1, X2, XX, XY
949  ! R.A. Output fields
950  ! BITMAP L.A. Data / no data bitmap
951  ! ----------------------------------------------------------------
952  !
953  ! 4. Subroutines used :
954  !
955  ! Name Type Module Description
956  ! ----------------------------------------------------------------
957  ! STRACE Subr. W3SERVMD Subroutine tracing.
958  ! EXTCDE Subr. Id. Abort program as graceful as possible.
959  ! W3S2XY Subr. Id. Convert from storage to spatial grid.
960  ! PUTGB Subr. NCEP GRIB1 library routine.
961  ! GRIBCREATE Subr. NCEP GRIB2 library routine.
962  ! ADDGRID Subr. NCEP GRIB2 library routine.
963  ! ADDFIELD Subr. NCEP GRIB2 library routine.
964  ! GRIBEND Subr. NCEP GRIB2 library routine.
965  ! WRYTE Subr. NCEP GRIB2 library routine.
966  ! ----------------------------------------------------------------
967  !
968  ! 5. Called by :
969  !
970  ! Program in which it is contained.
971  !
972  ! 6. Error messages :
973  !
974  ! None.
975  !
976  ! 7. Remarks :
977  !
978  ! - Note that arrays CX and CY of the main program now contain
979  ! the absolute current speed and direction respectively.
980  !
981  ! 8. Structure :
982  !
983  ! See source code.
984  !
985  ! 9. Switches :
986  !
987  ! !/S Enable subroutine tracing.
988  ! !/T Enable test output.
989  ! !/NCEP2 NCEP IBM calls to GRIB2 packer (follows updated grib2
990  ! tables under verification as of 02/10/2012).
991  !
992  ! 10. Source code :
993  !
994  !/ ------------------------------------------------------------------- /
995  USE w3servmd, ONLY : w3s2xy
996 #ifdef W3_RTD
997  USE w3servmd, ONLY : w3thrtn, w3xyrtn
998 #endif
999  !/
1000  !/ ------------------------------------------------------------------- /
1001  !/ Parameter list
1002  !/
1003  INTEGER, INTENT(IN) :: NX, NY, NSEA
1004  !/
1005  !/ ------------------------------------------------------------------- /
1006  !/ Local parameters
1007  !/
1008  INTEGER :: J, IXY, NDATA
1009  INTEGER :: IO
1010 #ifdef W3_S
1011  INTEGER, SAVE :: IENT = 0
1012 #endif
1013  REAL :: X1(NX*NY), X2(NX*NY), XX(NX*NY), &
1014  XY(NX*NY), CABS, UABS, &
1015  YY(NX*NY,0:NOSWLL), KPDS5A, KPDS5B, &
1016  KPDS5A1(3)
1017  LOGICAL*1 :: BITMAP(NX*NY)
1018  LOGICAL :: FLONE, FLTWO, FLDIR, FLTRI, FLPRT
1019  !/
1020  !/ ------------------------------------------------------------------- /
1021  !/
1022 #ifdef W3_S
1023  CALL strace (ient, 'W3EXGB')
1024 #endif
1025  !
1026 #ifdef W3_T
1027  WRITE (ndst,9000) ((flreq(ifi,ifj),ifj=1,ngrpp), ifi=1,nogrp)
1028  WRITE (ndst,9001) ndsdat, kpds, kgds
1029 #endif
1030  !
1031  !--- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1032  ! 1. Preparations
1033  !
1034  x1 = undef
1035  x2 = undef
1036  xx = undef
1037  xy = undef
1038  yy = undef
1039  !
1040  !
1041  !--- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1042  ! 2. Loop over output fields.
1043  !
1044  DO ifi=1, nogrp
1045  DO ifj=1, ngrpp
1046  IF ( flreq(ifi,ifj) ) THEN
1047 
1048 
1049  !
1050  ! Initialize array dimension flags
1051  !
1052  flone = .false.
1053  fltwo = .false.
1054  fldir = .false.
1055  fltri = .false.
1056  flprt = .false.
1057  !
1058 #ifdef W3_T
1059  WRITE (ndst,9020) idout(ifi,ifj)
1060 #endif
1061  !
1062  ! 2.a Set output arrays and parameters
1063  !
1064  ! Water depth
1065  !
1066  IF ( ifi .EQ. 1 .AND. ifj .EQ. 1 ) THEN
1067  flone = .true.
1068 #ifdef W3_NCEP2
1069  kpds(2) = 14
1070  kpds(1) = 4
1071 #endif
1072  CALL w3s2xy ( nsea, nsea, nx, ny, dw(1:nsea) &
1073  , mapsf, x1 )
1074  !
1075  ! Current
1076  !
1077  ELSE IF ( ifi .EQ. 1 .AND. ifj .EQ. 2 ) THEN
1078  fltwo = .true.
1079 #ifdef W3_NCEP2
1080  kpds(2) = 1
1081  kpds(1) = 1
1082 #endif
1083 #ifdef W3_RTD
1084  ! Rotate x,y vector back to standard pole
1085  IF ( flagunr ) CALL w3xyrtn(nsea, cx, cy, angld)
1086 #endif
1087  CALL w3s2xy ( nsea, nsea, nx, ny, cx(1:nsea) &
1088  , mapsf, xx )
1089  CALL w3s2xy ( nsea, nsea, nx, ny, cy(1:nsea) &
1090  , mapsf, xy )
1091  DO isea=1, nsea
1092  IF (cx(isea) .NE. undef) THEN
1093  cabs = sqrt(cx(isea)**2+cy(isea)**2)
1094  IF ( cabs .GT. 0.001 ) THEN
1095  cy(isea) = mod( 630. - &
1096  rade*atan2(cy(isea),cx(isea)) , 360. )
1097  ELSE
1098  cy(isea) = 0.
1099  END IF
1100  ELSE
1101  cabs = undef
1102  cy(isea) = undef
1103  END IF
1104  cx(isea) = cabs
1105  END DO
1106  CALL w3s2xy ( nsea, nsea, nx, ny, cx(1:nsea) &
1107  , mapsf, x1 )
1108  CALL w3s2xy ( nsea, nsea, nx, ny, cy(1:nsea) &
1109  , mapsf, x2 )
1110  !
1111  ! Wind speed
1112  !
1113  ELSE IF ( ifi .EQ. 1 .AND. ifj .EQ. 3 ) THEN
1114  fltwo = .true.
1115 #ifdef W3_NCEP2
1116  kpds(2) = 1
1117  kpds(1) = 2
1118  listsec0(1) = 0
1119 #endif
1120 #ifdef W3_RTD
1121  ! Rotate x,y vector back to standard pole
1122  IF ( flagunr ) CALL w3xyrtn(nsea, ua, ud, angld)
1123 #endif
1124  CALL w3s2xy ( nsea, nsea, nx, ny, ua(1:nsea) &
1125  , mapsf, xx )
1126  CALL w3s2xy ( nsea, nsea, nx, ny, ud(1:nsea) &
1127  , mapsf, xy )
1128  DO isea=1, nsea
1129  IF (ua(isea) .NE. undef) THEN
1130  uabs = sqrt(ua(isea)**2+ud(isea)**2)
1131  IF ( uabs .GT. 0.001 ) THEN
1132  ud(isea) = mod( 630. - &
1133  rade*atan2(ud(isea),ua(isea)) , 360. )
1134  ELSE
1135  ud(isea) = 0.
1136  END IF
1137  ELSE
1138  uabs = undef
1139  ud(isea) = undef
1140  END IF
1141  ua(isea) = uabs
1142  END DO
1143  CALL w3s2xy ( nsea, nsea, nx, ny, ua(1:nsea) &
1144  , mapsf, x1 )
1145  CALL w3s2xy ( nsea, nsea, nx, ny, ud(1:nsea) &
1146  , mapsf, x2 )
1147  !
1148  ! Air-sea temp. dif.
1149  !
1150  ELSE IF ( ifi .EQ. 1 .AND. ifj .EQ. 4 ) THEN
1151  flone = .true.
1152 #ifdef W3_NCEP2
1153  kpds(2) = 255
1154  kpds(1) = 3
1155 #endif
1156  CALL w3s2xy ( nsea, nsea, nx, ny, as(1:nsea) &
1157  , mapsf, x1 )
1158  !
1159  ! Water level
1160  !
1161  ELSE IF ( ifi .EQ. 1 .AND. ifj .EQ. 5 ) THEN
1162  flone = .true.
1163 #ifdef W3_NCEP2
1164  kpds(2) = 1
1165  kpds(1) = 3
1166 #endif
1167  CALL w3s2xy ( nsea, nsea, nx, ny, wlv , mapsf, x1 )
1168  !
1169  ! Ice concentration
1170  !
1171  ELSE IF ( ifi .EQ. 1 .AND. ifj .EQ. 6 ) THEN
1172  flone = .true.
1173 #ifdef W3_NCEP2
1174  kpds(2) = 0
1175  kpds(1) = 2
1176 #endif
1177  CALL w3s2xy ( nsea, nsea, nx, ny, ice , mapsf, x1 )
1178  !
1179  ! Atmospheric momentum
1180  !
1181  ELSE IF ( ifi .EQ. 1 .AND. ifj .EQ. 8 ) THEN
1182  fltwo = .true.
1183 #ifdef W3_NCEP2
1184  kpds(2) = 1
1185  kpds(1) = 2
1186  listsec0(1) = 0
1187 #endif
1188 #ifdef W3_RTD
1189  ! Rotate x,y vector back to standard pole
1190  IF ( flagunr ) CALL w3xyrtn(nsea, taua, tauadir, angld)
1191 #endif
1192  CALL w3s2xy ( nsea, nsea, nx, ny, taua(1:nsea) &
1193  , mapsf, xx )
1194  CALL w3s2xy ( nsea, nsea, nx, ny, tauadir(1:nsea) &
1195  , mapsf, xy )
1196  DO isea=1, nsea
1197  IF (taua(isea) .NE. undef) THEN
1198  uabs = sqrt(taua(isea)**2+tauadir(isea)**2)
1199  IF ( uabs .GT. 0.001 ) THEN
1200  tauadir(isea) = mod( 630. - &
1201  rade*atan2(tauadir(isea),taua(isea)) , 360. )
1202  ELSE
1203  tauadir(isea) = 0.
1204  END IF
1205  ELSE
1206  uabs = undef
1207  tauadir(isea) = undef
1208  END IF
1209  taua(isea) = uabs
1210  END DO
1211  CALL w3s2xy ( nsea, nsea, nx, ny, taua(1:nsea) &
1212  , mapsf, x1 )
1213  CALL w3s2xy ( nsea, nsea, nx, ny, tauadir(1:nsea) &
1214  , mapsf, x2 )
1215  !
1216  ! Air density
1217  !
1218  ELSE IF ( ifi .EQ. 1 .AND. ifj .EQ. 9 ) THEN
1219  flone = .true.
1220 #ifdef W3_NCEP2
1221  kpds(2) = 0
1222  kpds(1) = 2
1223 #endif
1224  CALL w3s2xy ( nsea, nsea, nx, ny, rhoair, mapsf, x1 )
1225  !
1226  ! Significant wave height
1227  !
1228  ELSE IF ( ifi .EQ. 2 .AND. ifj .EQ. 1 ) THEN
1229  flone = .true.
1230 #ifdef W3_NCEP2
1231  kpds(2) = 3
1232 #endif
1233  CALL w3s2xy ( nsea, nsea, nx, ny, hs , mapsf, x1 )
1234  !
1235  ! Mean wave length
1236  !
1237  ELSE IF ( ifi .EQ. 2 .AND. ifj .EQ. 2 ) THEN
1238  flone = .true.
1239 #ifdef W3_NCEP2
1240  kpds(2) = 193
1241 #endif
1242  CALL w3s2xy ( nsea, nsea, nx, ny, wlm , mapsf, x1 )
1243  !
1244  ! Mean wave period (based on second moment)
1245  !
1246  ELSE IF ( ifi .EQ. 2 .AND. ifj .EQ. 3 ) THEN
1247  flone = .true.
1248 #ifdef W3_NCEP2
1249  if ((gen_pro.eq.1) .or. (gen_pro.eq.0)) then
1250  kpds(2) = 28
1251  else
1252  kpds(2) = 25
1253  endif
1254 #endif
1255 
1256  CALL w3s2xy ( nsea, nsea, nx, ny, t02 , mapsf, x1 )
1257  !
1258  ! Mean wave period (based on first moment)
1259  !
1260  ELSE IF ( ifi .EQ. 2 .AND. ifj .EQ. 4 ) THEN
1261  flone = .true.
1262 #ifdef W3_NCEP2
1263  kpds(2) = 15
1264 #endif
1265  CALL w3s2xy ( nsea, nsea, nx, ny, t0m1 , mapsf, x1 )
1266  !
1267  ! Mean wave period (based on first inverse moment)
1268  !
1269  ELSE IF ( ifi .EQ. 2 .AND. ifj .EQ. 5 ) THEN
1270  flone = .true.
1271 #ifdef W3_NCEP2
1272  if ((gen_pro.eq.1) .or. (gen_pro.eq.0)) then
1273  kpds(2) = 34
1274  else
1275  kpds(2) = 15
1276  endif
1277 #endif
1278  CALL w3s2xy ( nsea, nsea, nx, ny, t01 , mapsf, x1 )
1279  !
1280  ! Peak frequency
1281  !
1282  ELSE IF ( ifi .EQ. 2 .AND. ifj .EQ. 6 ) THEN
1283  flone = .true.
1284 #ifdef W3_NCEP2
1285  kpds(2) = 11
1286 #endif
1287  DO isea=1, nsea
1288  IF ( fp0(isea) .NE. undef .AND. fp0(isea) .NE. 0 ) THEN
1289  fp0(isea) = 1. / max(fr1,fp0(isea)) ! Limit FP to lowest discrete frequency
1290  END IF
1291  END DO
1292  CALL w3s2xy ( nsea, nsea, nx, ny, fp0 , mapsf, x1 )
1293  !
1294  !
1295  ! Mean wave direction
1296  !
1297  ELSE IF ( ifi .EQ. 2 .AND. ifj .EQ. 7 ) THEN
1298  flone = .true.
1299 #ifdef W3_NCEP2
1300  kpds(2) = 14
1301 #endif
1302 #ifdef W3_RTD
1303  ! Rotate direction back to standard pole
1304  IF ( flagunr ) CALL w3thrtn(nsea, thm, angld, .false.)
1305 #endif
1306  DO isea=1, nsea
1307  IF ( thm(isea) .NE. undef ) &
1308  thm(isea) = mod( 630. - rade*thm(isea) , 360. )
1309  END DO
1310  CALL w3s2xy ( nsea, nsea, nx, ny, thm , mapsf, x1 )
1311  !
1312  ! Directional spread
1313  !
1314  ELSE IF ( ifi .EQ. 2 .AND. ifj .EQ. 8 ) THEN
1315  flone = .true.
1316 #ifdef W3_NCEP2
1317  kpds(2) = 31
1318 #endif
1319  CALL w3s2xy ( nsea, nsea, nx, ny, ths , mapsf, x1 )
1320  !
1321  ! Peak direction
1322  !
1323  ELSE IF ( ifi .EQ. 2 .AND. ifj .EQ. 9 ) THEN
1324  flone = .true.
1325 #ifdef W3_NCEP2
1326  if ((gen_pro.eq.1) .or. (gen_pro.eq.0)) then
1327  kpds(2) = 46
1328  else
1329  kpds(2) = 10
1330  endif
1331 #endif
1332 #ifdef W3_RTD
1333  ! Rotate direction back to standard pole
1334  IF ( flagunr ) CALL w3thrtn(nsea, thp0, angld, .false.)
1335 #endif
1336  DO isea=1, nsea
1337  IF ( thp0(isea) .NE. undef ) THEN
1338  thp0(isea) = mod( 630-rade*thp0(isea) , 360. )
1339  END IF
1340  END DO
1341  CALL w3s2xy ( nsea, nsea, nx, ny, thp0 , mapsf, x1 )
1342  !
1343  ! Mean wave number
1344  !
1345  ELSE IF ( ifi .EQ. 2 .AND. ifj .EQ. 19 ) THEN
1346  flone = .true.
1347 #ifdef W3_NCEP2
1348  kpds(2) = 255
1349 #endif
1350  CALL w3s2xy ( nsea, nsea, nx, ny, wnmean, mapsf, x1 )
1351  !
1352  ! Partitioned wave height
1353  !
1354  ELSE IF ( ifi .EQ. 4 .AND. ifj .EQ. 1 ) THEN
1355  flprt = .true.
1356 #ifdef W3_NCEP2
1357  kpds5a = 5
1358  kpds5b = 8
1359  if ((gen_pro.eq.1) .or. (gen_pro.eq.0)) then
1360  kpds5a1(1) =47
1361  kpds5a1(2) =48
1362  kpds5a1(3) =49
1363  else
1364  kpds5b = 8
1365  endif
1366 #endif
1367  CALL w3s2xy &
1368  ( nsea, nsea, nx, ny, phs(:,0), mapsf, yy(:,0) )
1369  DO i=1, noswll
1370  CALL w3s2xy &
1371  ( nsea, nsea, nx, ny, phs(:,i), mapsf, yy(:,i) )
1372  END DO
1373  !
1374  ! Partitioned peak period
1375  !
1376  ELSE IF ( ifi .EQ. 4 .AND. ifj .EQ. 2 ) THEN
1377  flprt = .true.
1378 #ifdef W3_NCEP2
1379  kpds5a = 6
1380  kpds5b = 9
1381  if ((gen_pro.eq.1) .or. (gen_pro.eq.0)) then
1382  kpds5a1(1) = 50
1383  kpds5a1(2) = 51
1384  kpds5a1(3) = 52
1385  else
1386  kpds5b = 9
1387  endif
1388 #endif
1389  CALL w3s2xy &
1390  ( nsea, nsea, nx, ny, ptp(:,0), mapsf, yy(:,0) )
1391  DO i=1, noswll
1392  CALL w3s2xy &
1393  ( nsea, nsea, nx, ny, ptp(:,i), mapsf, yy(:,i) )
1394  END DO
1395  !
1396  ! Partitioned peak wave length
1397  !
1398  ELSE IF ( ifi .EQ. 4 .AND. ifj .EQ. 3 ) THEN
1399  flprt = .true.
1400 #ifdef W3_NCEP2
1401  kpds5a = 193
1402  kpds5b = 193
1403 #endif
1404  CALL w3s2xy &
1405  ( nsea, nsea, nx, ny, plp(:,0), mapsf, yy(:,0) )
1406  DO i=1, noswll
1407  CALL w3s2xy &
1408  ( nsea, nsea, nx, ny, plp(:,i), mapsf, yy(:,i) )
1409  END DO
1410  !
1411  ! Partitioned mean direction
1412  !
1413  ELSE IF ( ifi .EQ. 4 .AND. ifj .EQ. 4 ) THEN
1414  flprt = .true.
1415 #ifdef W3_NCEP2
1416  kpds5a = 4
1417  kpds5b = 7
1418  if ((gen_pro.eq.1) .or. (gen_pro.eq.0)) then
1419  kpds5a1(1) = 53
1420  kpds5a1(2) = 54
1421  kpds5a1(3) = 55
1422  else
1423  kpds5b = 7
1424  endif
1425 #endif
1426 #ifdef W3_RTD
1427  DO i = 0,noswll
1428  ! Rotate direction back to standard pole
1429  IF ( flagunr ) CALL w3thrtn(nsea, pdir(:,i), angld, .false.)
1430  END DO
1431 #endif
1432  DO isea = 1,nsea
1433  DO i = 0,noswll
1434  IF ( pdir(isea,i) .NE. undef ) THEN
1435  pdir(isea,i) = mod( 630 - rade*pdir(isea,i) , 360. )
1436  END IF
1437  END DO
1438  END DO
1439  CALL w3s2xy &
1440  ( nsea, nsea, nx, ny, pdir(:,0), mapsf, yy(:,0) )
1441  DO i=1, noswll
1442  CALL w3s2xy &
1443  ( nsea, nsea, nx, ny, pdir(:,i), mapsf, yy(:,i) )
1444  END DO
1445  !
1446  ! Partitioned Directional spread
1447  !
1448  ELSE IF ( ifi .EQ. 4 .AND. ifj .EQ. 5 ) THEN
1449  flprt = .true.
1450 #ifdef W3_NCEP2
1451  kpds5a = 32
1452  kpds5b = 33
1453 #endif
1454  CALL w3s2xy &
1455  ( nsea, nsea, nx, ny, psi(:,0), mapsf, yy(:,0) )
1456  DO i=1, noswll
1457  CALL w3s2xy &
1458  ( nsea, nsea, nx, ny, psi(:,i), mapsf, yy(:,i) )
1459  END DO
1460  !
1461  ! Partitioned wind sea fraction
1462  !
1463  ELSE IF ( ifi .EQ. 4 .AND. ifj .EQ. 6 ) THEN
1464  flprt = .true.
1465 #ifdef W3_NCEP2
1466  kpds5a = 255
1467  kpds5b = 255
1468 #endif
1469  CALL w3s2xy &
1470  ( nsea, nsea, nx, ny, pws(:,0), mapsf, yy(:,0) )
1471  DO i=1, noswll
1472  CALL w3s2xy &
1473  ( nsea, nsea, nx, ny, pws(:,i), mapsf, yy(:,i) )
1474  END DO
1475  !
1476  ! Total wind sea fraction
1477  !
1478  ELSE IF ( ifi .EQ. 4 .AND. ifj .EQ. 16 ) THEN
1479  flone = .true.
1480 #ifdef W3_NCEP2
1481  kpds(2) = 255
1482 #endif
1483  CALL w3s2xy ( nsea, nsea, nx, ny, pwst , mapsf, x1 )
1484  !
1485  ! Number of fields in partition
1486  !
1487  ELSE IF ( ifi .EQ. 4 .AND. ifj .EQ. 17 ) THEN
1488  flone = .true.
1489 #ifdef W3_NCEP2
1490  kpds(2) = 255
1491 #endif
1492  CALL w3s2xy ( nsea, nsea, nx, ny, pnr , mapsf, x1 )
1493  !
1494  ! Friction velocity
1495  !
1496  ELSE IF ( ifi .EQ. 5 .AND. ifj .EQ. 1 ) THEN
1497  fltwo = .true.
1498 #ifdef W3_NCEP2
1499  kpds(2) = 17
1500  kpds(1) = 1
1501 #endif
1502 #ifdef W3_RTD
1503  ! Rotate x,y vector back to standard pole
1504  IF ( flagunr ) CALL w3xyrtn(nsea, ust, ustdir, angld)
1505 #endif
1506  CALL w3s2xy ( nsea, nsea, nx, ny, ust(1:nsea) &
1507  , mapsf, x1 )
1508  CALL w3s2xy ( nsea, nsea, nx, ny, ustdir(1:nsea) &
1509  , mapsf, x2 )
1510  !
1511  ! Average source term time step
1512  !
1513  ELSE IF ( ifi .EQ. 9 .AND. ifj .EQ. 1 ) THEN
1514  flone = .true.
1515 #ifdef W3_NCEP2
1516  kpds(2) = 255
1517 #endif
1518  DO isea=1, nsea
1519  IF ( dtdyn(isea) .NE. undef ) &
1520  dtdyn(isea) = dtdyn(isea) / 60.
1521  END DO
1522  CALL w3s2xy ( nsea, nsea, nx, ny, dtdyn , mapsf, x1 )
1523  !
1524  ! Cut-off frequency
1525  !
1526  ELSE IF ( ifi .EQ. 9 .AND. ifj .EQ. 2 ) THEN
1527  flone = .true.
1528 #ifdef W3_NCEP2
1529  kpds(2) = 255
1530 #endif
1531  CALL w3s2xy ( nsea, nsea, nx, ny, fcut , mapsf, x1 )
1532  !
1533  ! CFL Maximum (in spatial space)
1534  !
1535  ELSE IF ( ifi .EQ. 9 .AND. ifj .EQ. 3 ) THEN
1536  flone = .true.
1537 #ifdef W3_NCEP2
1538  kpds(2) = 255
1539 #endif
1540  CALL w3s2xy ( nsea, nsea, nx, ny, cflxymax , mapsf, x1 )
1541  !
1542  ! CFL Maximum (in spectral space)
1543  !
1544  ELSE IF ( ifi .EQ. 9 .AND. ifj .EQ. 4 ) THEN
1545  flone = .true.
1546 #ifdef W3_NCEP2
1547  kpds(2) = 255
1548 #endif
1549  CALL w3s2xy ( nsea, nsea, nx, ny, cflthmax , mapsf, x1 )
1550  !
1551  ELSE
1552  WRITE (ndse,999)
1553  CALL extcde ( 1 )
1554  !
1555  END IF
1556  !
1557  ! 3 Perform output
1558  !
1559  ndata = nx*ny
1560  !
1561  ! 3.a Partitioned data
1562  !
1563  IF ( flprt ) THEN
1564  !
1565 #ifdef W3_NCEP2
1566  kpds(2) = kpds5a
1567 #endif
1568  DO ixy=1, nx*ny
1569  bitmap(ixy) = yy(ixy,0) .NE. undef
1570  END DO
1571 #ifdef W3_NCEP2
1572  CALL gribcreate (cgrib,lcgrib,listsec0,listsec1,io)
1573  IF (io .NE. 0) GOTO 810
1574  CALL addgrid (cgrib,lcgrib,igds,kgds,200,ideflist, &
1575  idefnum, io)
1576  IF (io .NE. 0) GOTO 820
1577  CALL addfield (cgrib,lcgrib,kpdsnum,kpds,200, &
1578  coordlist, numcoord, idrsnum, idrs, &
1579  200,yy(:,0), ndata, ibmp, bitmap, io)
1580  IF (io .NE. 0) GOTO 820
1581  CALL gribend (cgrib, lcgrib, lengrib, io)
1582  IF (io .NE. 0) GOTO 830
1583  CALL wryte (ndsdat, lengrib, cgrib)
1584 #endif
1585  !
1586 #ifdef W3_NCEP2
1587  if ((gen_pro.eq.0) .or. (gen_pro.eq.1)) then
1588  kpds(10) = 241
1589 #endif
1590  DO i=1, noswll
1591 #ifdef W3_NCEP2
1592  kpds(2) = kpds5a1(i)
1593  kpds(12) = i
1594 #endif
1595  DO ixy=1, nx*ny
1596  bitmap(ixy) = yy(ixy,i) .NE. undef
1597  END DO
1598 #ifdef W3_NCEP2
1599  CALL gribcreate (cgrib,lcgrib,listsec0,listsec1,io)
1600  IF (io .NE. 0) GOTO 810
1601  CALL addgrid (cgrib,lcgrib,igds,kgds,200,ideflist, &
1602  idefnum, io)
1603  IF (io .NE. 0) GOTO 820
1604  CALL addfield (cgrib,lcgrib,kpdsnum,kpds,200, &
1605  coordlist, numcoord, idrsnum, idrs, &
1606  200,yy(:,i), ndata, ibmp, bitmap, io)
1607  IF (io .NE. 0) GOTO 820
1608  CALL gribend (cgrib, lcgrib, lengrib, io)
1609  IF (io .NE. 0) GOTO 830
1610  CALL wryte (ndsdat, lengrib, cgrib)
1611 #endif
1612  END DO
1613 #ifdef W3_NCEP2
1614  ELSE
1615  kpds(2) = kpds5b
1616  kpds(10) = 241
1617 #endif
1618  DO i=1, noswll
1619 #ifdef W3_NCEP2
1620  kpds(12) = i
1621 #endif
1622  DO ixy=1, nx*ny
1623  bitmap(ixy) = yy(ixy,i) .NE. undef
1624  END DO
1625 #ifdef W3_NCEP2
1626  CALL gribcreate (cgrib,lcgrib,listsec0,listsec1,io)
1627  IF (io .NE. 0) GOTO 810
1628  CALL addgrid (cgrib,lcgrib,igds,kgds,200,ideflist, &
1629  idefnum, io)
1630  IF (io .NE. 0) GOTO 820
1631  CALL addfield (cgrib,lcgrib,kpdsnum,kpds,200, &
1632  coordlist, numcoord, idrsnum, idrs, &
1633  200,yy(:,i), ndata, ibmp, bitmap, io)
1634  IF (io .NE. 0) GOTO 820
1635  CALL gribend (cgrib, lcgrib, lengrib, io)
1636  IF (io .NE. 0) GOTO 830
1637  CALL wryte (ndsdat, lengrib, cgrib)
1638 #endif
1639  END DO
1640 #ifdef W3_NCEP2
1641  ENDIF
1642  kpds(10) = 1
1643  kpds(12) = 1
1644 #endif
1645  !
1646  ! 3.b Other data
1647  !
1648  ELSE IF (flone) THEN
1649  !
1650  DO ixy=1, nx*ny
1651  bitmap(ixy) = x1(ixy) .NE. undef
1652  END DO
1653  !
1654 #ifdef W3_NCEP2
1655  CALL gribcreate (cgrib,lcgrib,listsec0,listsec1,io)
1656  IF (io .NE. 0) GOTO 810
1657  CALL addgrid (cgrib,lcgrib,igds,kgds,200,ideflist, &
1658  idefnum, io)
1659  IF (io .NE. 0) GOTO 820
1660  CALL addfield (cgrib,lcgrib,kpdsnum,kpds,200, &
1661  coordlist, numcoord, idrsnum, idrs, &
1662  200,x1, ndata, ibmp, bitmap, io)
1663  IF (io .NE. 0) GOTO 820
1664  CALL gribend (cgrib, lcgrib, lengrib, io)
1665  IF (io .NE. 0) GOTO 830
1666  CALL wryte (ndsdat, lengrib, cgrib)
1667 #endif
1668  !
1669  ELSE IF ( fltwo ) THEN
1670  !
1671  DO ixy=1, nx*ny
1672  bitmap(ixy) = x1(ixy) .NE. undef
1673  END DO
1674 #ifdef W3_NCEP2
1675  CALL gribcreate (cgrib,lcgrib,listsec0,listsec1,io)
1676  IF (io .NE. 0) GOTO 810
1677  CALL addgrid (cgrib,lcgrib,igds,kgds,200,ideflist, &
1678  idefnum, io)
1679  IF (io .NE. 0) GOTO 820
1680  CALL addfield (cgrib,lcgrib,kpdsnum,kpds,200, &
1681  coordlist, numcoord, idrsnum, idrs, &
1682  200,x1, ndata, ibmp, bitmap, io)
1683  IF (io .NE. 0) GOTO 820
1684  CALL gribend (cgrib, lcgrib, lengrib, io)
1685  IF (io .NE. 0) GOTO 830
1686  CALL wryte (ndsdat, lengrib, cgrib)
1687 #endif
1688 
1689 #ifdef W3_NCEP2
1690  kpds(2) = 0
1691  CALL gribcreate (cgrib,lcgrib,listsec0,listsec1,io)
1692  IF (io .NE. 0) GOTO 810
1693  CALL addgrid (cgrib,lcgrib,igds,kgds,200,ideflist, &
1694  idefnum, io)
1695  IF (io .NE. 0) GOTO 820
1696  CALL addfield (cgrib,lcgrib,kpdsnum,kpds,200, &
1697  coordlist, numcoord, idrsnum, idrs, &
1698  200,x2, ndata, ibmp, bitmap, io)
1699  IF (io .NE. 0) GOTO 820
1700  CALL gribend (cgrib, lcgrib, lengrib, io)
1701  IF (io .NE. 0) GOTO 830
1702  CALL wryte (ndsdat, lengrib, cgrib)
1703  kpds(2) = 2
1704  CALL gribcreate (cgrib,lcgrib,listsec0,listsec1,io)
1705  IF (io .NE. 0) GOTO 810
1706  CALL addgrid (cgrib,lcgrib,igds,kgds,200,ideflist, &
1707  idefnum, io)
1708  IF (io .NE. 0) GOTO 820
1709  CALL addfield (cgrib,lcgrib,kpdsnum,kpds,200, &
1710  coordlist, numcoord, idrsnum, idrs, &
1711  200,xx, ndata, ibmp, bitmap, io)
1712  IF (io .NE. 0) GOTO 820
1713  CALL gribend (cgrib, lcgrib, lengrib, io)
1714  IF (io .NE. 0) GOTO 830
1715  CALL wryte (ndsdat, lengrib, cgrib)
1716  kpds(2) = 3
1717  CALL gribcreate (cgrib,lcgrib,listsec0,listsec1,io)
1718  IF (io .NE. 0) GOTO 810
1719  CALL addgrid (cgrib,lcgrib,igds,kgds,200,ideflist, &
1720  idefnum, io)
1721  IF (io .NE. 0) GOTO 820
1722  CALL addfield (cgrib,lcgrib,kpdsnum,kpds,200, &
1723  coordlist, numcoord, idrsnum, idrs, &
1724  200,xy, ndata, ibmp, bitmap, io)
1725  IF (io .NE. 0) GOTO 820
1726  CALL gribend (cgrib, lcgrib, lengrib, io)
1727  IF (io .NE. 0) GOTO 830
1728  CALL wryte (ndsdat, lengrib, cgrib)
1729 #endif
1730  !
1731  END IF
1732 #ifdef W3_NCEP2
1733  listsec0(1) = 10
1734  kpds(1) = 0
1735 #endif
1736  !
1737  ! ... End of fields loop
1738  !
1739  END IF
1740  END DO
1741  END DO
1742  !
1743  RETURN
1744  !
1745  ! Error escape locations
1746  !
1747 #ifdef W3_NCEP2
1748 810 CONTINUE
1749  WRITE (ndse,1010) io
1750  CALL extcde ( 20 )
1751 820 CONTINUE
1752  WRITE (ndse,1020) io
1753  CALL extcde ( 30 )
1754 830 CONTINUE
1755  WRITE (ndse,1030) io
1756  CALL extcde ( 40 )
1757 #endif
1758  !
1759  ! Formats
1760  !
1761 999 FORMAT (/' *** WAVEWATCH III ERROR IN W3EXGB :'/ &
1762  ' PLEASE UPDATE FIELDS !!! '/)
1763  !
1764 #ifdef W3_NCEP2
1765 1000 FORMAT (/' *** WAVEWATCH III ERROR IN W3EXGB : '/ &
1766  ' ERROR IN OPENING OUTPUT FILE'/ &
1767  ' IOSTAT =',i5/)
1768 #endif
1769  !
1770 #ifdef W3_NCEP2
1771 1010 FORMAT (/' *** WAVEWATCH III ERROR IN W3EXGB : '/ &
1772  ' ERROR CREATING NEW GRIB2 FIELD'/ &
1773  ' IOSTAT =',i5/)
1774 #endif
1775  !
1776 #ifdef W3_NCEP2
1777 1020 FORMAT (/' *** WAVEWATCH III ERROR IN W3EXGB : '/ &
1778  ' ERROR ADDING GRIB2 FIELD'/ &
1779  ' IOSTAT =',i5/)
1780 #endif
1781  !
1782 #ifdef W3_NCEP2
1783 1030 FORMAT (/' *** WAVEWATCH III ERROR IN W3EXGB : '/ &
1784  ' ERROR ENDING GRIB2 MESSAGE'/ &
1785  ' IOSTAT =',i5/)
1786 #endif
1787  !
1788 #ifdef W3_T
1789 9000 FORMAT (' TEST W3EXGB : FLAGS :',40l2)
1790 9001 FORMAT (' TEST W3EXGB : NDSDAT :',i4/ &
1791  ' KPDS :',13i4/ &
1792  ' ',12i4/ &
1793  ' KGDS :',8i6/ &
1794  ' ',8i6/ &
1795  ' ',6i6)
1796 #endif
1797  !
1798 #ifdef W3_T
1799 9012 FORMAT (' TEST W3EXGB : BLOK PARS : ',3i4)
1800 9014 FORMAT (' BASE NAME : ',a)
1801 #endif
1802  !
1803 #ifdef W3_T
1804 9020 FORMAT (' TEST W3EXGB : OUTPUT FIELD : ',a)
1805 #endif
1806  !/
1807  !/ End of W3EXGB ----------------------------------------------------- /
1808  !/

References constants::rade, constants::undef, w3servmd::w3s2xy(), w3servmd::w3thrtn(), and w3servmd::w3xyrtn().

Referenced by w3grib().

◆ w3grib()

program w3grib

Post-processing of grid output.

Data is read from the grid output file out_grd.ww3 (raw data) and from the file ww3_grib.inp ( NDSI, output requests ). Model definition and raw data files are read using WAVEWATCH III subroutines. GRIB packing is performed using NCEP's W3 library (not supplied).

When adding new parameters to GRIB packing, keep in mind that packing is done differently for scalar and vector quantities

Author
H. L. Tolman
A. Chawla
J.-H. Alves
Date
22-Mar-2021

Definition at line 27 of file ww3_grib.F90.

References w3timemd::dsec21(), w3servmd::extcde(), file(), w3odatmd::flogd, w3odatmd::flogrd, w3odatmd::fnmpre, w3wdatmd::ice, w3odatmd::idout, w3servmd::itrace(), w3odatmd::ndse, w3odatmd::ndso, w3odatmd::ndst, w3servmd::nextln(), w3odatmd::ngrpp, w3odatmd::noge, w3odatmd::nogrp, w3odatmd::noswll, w3wdatmd::rhoair, w3timemd::stme21(), w3servmd::strace(), w3timemd::tick21(), w3wdatmd::time, constants::undef, w3wdatmd::ust, w3wdatmd::ustdir, w3exgb(), w3iogomd::w3iogo(), w3iogrmd::w3iogr(), w3adatmd::w3naux(), w3wdatmd::w3ndat(), w3gdatmd::w3nmod(), w3odatmd::w3nout(), w3iogomd::w3readflgrd(), w3adatmd::w3seta(), w3gdatmd::w3setg(), w3odatmd::w3seto(), w3wdatmd::w3setw(), and w3wdatmd::wlv.

w3adatmd::dtdyn
real, dimension(:), pointer dtdyn
Definition: w3adatmd.F90:620
w3adatmd::as
real, dimension(:), pointer as
Definition: w3adatmd.F90:584
w3adatmd::t02
real, dimension(:), pointer t02
Definition: w3adatmd.F90:587
w3adatmd::cflxymax
real, dimension(:), pointer cflxymax
Definition: w3adatmd.F90:620
w3adatmd::fcut
real, dimension(:), pointer fcut
Definition: w3adatmd.F90:620
w3adatmd::ptp
real, dimension(:,:), pointer ptp
Definition: w3adatmd.F90:597
w3adatmd::dw
real, dimension(:), pointer dw
Definition: w3adatmd.F90:584
w3adatmd::t01
real, dimension(:), pointer t01
Definition: w3adatmd.F90:587
w3wdatmd::wlv
real, dimension(:), pointer wlv
Definition: w3wdatmd.F90:183
w3adatmd::pdir
real, dimension(:,:), pointer pdir
Definition: w3adatmd.F90:597
w3adatmd::thp0
real, dimension(:), pointer thp0
Definition: w3adatmd.F90:587
constants::rade
real, parameter rade
RADE Conversion factor from radians to degrees.
Definition: constants.F90:76
w3adatmd::phs
real, dimension(:,:), pointer phs
Definition: w3adatmd.F90:597
w3adatmd::hs
real, dimension(:), pointer hs
Definition: w3adatmd.F90:587
w3servmd::w3thrtn
subroutine w3thrtn(NSEA, THETA, AnglD, Degrees)
Definition: w3servmd.F90:1333
w3adatmd::pws
real, dimension(:,:), pointer pws
Definition: w3adatmd.F90:597
w3adatmd::plp
real, dimension(:,:), pointer plp
Definition: w3adatmd.F90:597
w3adatmd::cflthmax
real, dimension(:), pointer cflthmax
Definition: w3adatmd.F90:620
w3adatmd::psi
real, dimension(:,:), pointer psi
Definition: w3adatmd.F90:597
w3servmd
Definition: w3servmd.F90:3
w3adatmd::ths
real, dimension(:), pointer ths
Definition: w3adatmd.F90:587
w3adatmd::ud
real, dimension(:), pointer ud
Definition: w3adatmd.F90:584
w3adatmd::pwst
real, dimension(:), pointer pwst
Definition: w3adatmd.F90:597
w3adatmd::cy
real, dimension(:), pointer cy
Definition: w3adatmd.F90:584
w3adatmd::pnr
real, dimension(:), pointer pnr
Definition: w3adatmd.F90:597
w3adatmd::taua
real, dimension(:), pointer taua
Definition: w3adatmd.F90:584
w3adatmd::wlm
real, dimension(:), pointer wlm
Definition: w3adatmd.F90:587
w3adatmd::wnmean
real, dimension(:), pointer wnmean
Definition: w3adatmd.F90:587
w3wdatmd::ice
real, dimension(:), pointer ice
Definition: w3wdatmd.F90:183
w3wdatmd::ust
real, dimension(:), pointer ust
Definition: w3wdatmd.F90:183
w3adatmd::ua
real, dimension(:), pointer ua
Definition: w3adatmd.F90:584
w3adatmd::fp0
real, dimension(:), pointer fp0
Definition: w3adatmd.F90:587
w3servmd::w3xyrtn
subroutine w3xyrtn(NSEA, XVEC, YVEC, AnglD)
Definition: w3servmd.F90:1387
w3adatmd::tauadir
real, dimension(:), pointer tauadir
Definition: w3adatmd.F90:584
w3wdatmd::rhoair
real, dimension(:), pointer rhoair
Definition: w3wdatmd.F90:183
w3wdatmd::ustdir
real, dimension(:), pointer ustdir
Definition: w3wdatmd.F90:183
w3adatmd::thm
real, dimension(:), pointer thm
Definition: w3adatmd.F90:587
w3adatmd::cx
real, dimension(:), pointer cx
Definition: w3adatmd.F90:584
constants::undef
real undef
UNDEF Value for undefined variable in output.
Definition: constants.F90:84
w3servmd::w3s2xy
subroutine w3s2xy(NSEA, MSEA, MX, MY, S, MAPSF, XY)
Definition: w3servmd.F90:337
w3adatmd::t0m1
real, dimension(:), pointer t0m1
Definition: w3adatmd.F90:587