WAVEWATCH III  beta 0.0.1
ww3_gspl.F90 File Reference

Contains the grid splitting program. More...

Go to the source code of this file.

Data Types

type  stats_grid
 
type  stats_mean
 
type  part_grid
 

Functions/Subroutines

program w3gspl
 Grid splitting program. More...
 
subroutine grinfo
 Compile statistical info on all sub grids (no halo). More...
 
subroutine grtrim
 Trim edges of all grids where they are next to another grid or next to unassigned grid points. More...
 
subroutine grfill (ND)
 Reassign unassigned grid points to grids, starting with the smallest grids. More...
 
subroutine grlost
 Reassign unassigned grid points to grids. More...
 
subroutine grsqrg
 Attempt to square-up grid. More...
 
subroutine grsngl (OK)
 Remove seapoints with only one adjacent point in same grid. More...
 
subroutine grsepa (OK, FRAC)
 Remove smaller grid parts. More...
 
subroutine grfsml
 Subroutine called when lowest grid size is stuck. More...
 
subroutine grflrg
 Like GRFSML for largest grid ... More...
 
subroutine gr1grd
 Extract single grid from master map. More...
 

Detailed Description

Contains the grid splitting program.

Author
H. L. Tolman
Date
18-Nov-2013

Definition in file ww3_gspl.F90.

Function/Subroutine Documentation

◆ gr1grd()

subroutine w3gspl::gr1grd

Extract single grid from master map.

Extract single grid from master map, including halo needed for grid overlap in ww3_multi.

Author
H. L. Tolman
Date
18-Nov-2012

Definition at line 3084 of file ww3_gspl.F90.

3084  !/
3085  !/ +-----------------------------------+
3086  !/ | WAVEWATCH III NOAA/NCEP |
3087  !/ | H. L. Tolman |
3088  !/ | FORTRAN 90 |
3089  !/ | Last update : 18-Nov-2012 |
3090  !/ +-----------------------------------+
3091  !/
3092  !/ 23-Sep-2012 : Origination. ( version 4.10 )
3093  !/ 24-Jan-2013 : Correct X0 to be in range. ( version 4.10 )
3094  !/ 04-Feb-2013 : Add corner point to halo. ( version 4.10 )
3095  !/ 18-Nov-2012 : Add user-defined halo extension. ( version 4.14 )
3096  !/
3097  ! 1. Purpose :
3098  !
3099  ! Extract single grid from master map, including halo needed for
3100  ! grid overlap in ww3_multi.
3101  !
3102  ! 10. Source code :
3103  !
3104  !/ ------------------------------------------------------------------- /
3105  !/
3106  IMPLICIT NONE
3107  !/
3108  !/ ------------------------------------------------------------------- /
3109  !/ Parameter list
3110  !/
3111  !/ ------------------------------------------------------------------- /
3112  !/ Local parameters
3113  !/
3114  INTEGER :: NIT, IIT, IXL, IXH, IYL, IYH, NOCNT,&
3115  NOCNTM, NOCNTL, JX, JY, ISEA, MX, MY
3116  INTEGER :: MTMP2(NY,NX)
3117 #ifdef W3_S
3118  INTEGER, SAVE :: IENT = 0
3119 #endif
3120  REAL :: XOFF
3121  LOGICAL :: MASK(NY,NX), LEFT, RIGHT, THERE
3122  !/
3123  !/ ------------------------------------------------------------------- /
3124  !/
3125 #ifdef W3_S
3126  CALL strace (ient, 'GR1GRD')
3127 #endif
3128  !
3129 #ifdef W3_T7
3130  WRITE (ndst,9000) ig
3131 #endif
3132  !
3133  ! 1. Set up MTEMP with MAPSTA 0,1,3 for grid ------------------------ *
3134  !
3135  DO ix=1, nx
3136  DO iy=1, ny
3137  IF ( msplit(iy,ix) .EQ. ig ) THEN
3138  mtemp(iy,ix) = 1
3139  ELSE IF ( msplit(iy,ix) .GT. 0 ) THEN
3140  mtemp(iy,ix) = 3
3141  ELSE
3142  mtemp(iy,ix) = 0
3143  END IF
3144  END DO
3145  END DO
3146  !
3147  ! 2. Add ALL MAPSTA = 2 points to grid ------------------------------ *
3148  !
3149  DO ix=1, nx
3150  DO iy=1, ny
3151  IF ( mapsta(iy,ix) .EQ. 2 ) THEN
3152  mtemp(iy,ix) = 2
3153  END IF
3154  END DO
3155  END DO
3156  !
3157  ! 3. Add halo ------------------------------------------------------- *
3158  ! 3.a Set up halo width depending on scheme and time steps
3159  ! NEEDED TO SET UP A LITTLE WIDER. NOT SURE WHY. NEED TO CHECK WITH
3160  ! WMEQL SUBROUTINE.
3161  !
3162 #ifdef W3_PR0
3163  nit = 0
3164 #endif
3165 #ifdef W3_PR1
3166  nit = 1 + nhext + ( 1 + int(dtmax/dtcfl-0.001) ) * 1
3167 #endif
3168 #ifdef W3_UQ
3169  nit = 1 + nhext + ( 1 + int(dtmax/dtcfl-0.001) ) * 3
3170 #endif
3171 #ifdef W3_UNO
3172  nit = 1 + nhext + ( 1 + int(dtmax/dtcfl-0.001) ) * 3
3173 #endif
3174  !
3175  ! 3.b Exand halo
3176  !
3177  DO iit=1, nit
3178  !
3179  mask = .false.
3180  !
3181  DO ix=1, nx
3182  ixl = 1 + mod(ix-2+nx,nx)
3183  ixh = 1 + mod(ix,nx)
3184  DO iy=2, ny-1
3185  IF ( mtemp(iy,ix) .EQ. 3 ) mask(iy,ix) = &
3186  ( ( mtemp(iy+1,ix ) .EQ. 1 ) .OR. &
3187  ( mtemp(iy-1,ix ) .EQ. 1 ) .OR. &
3188  ( mtemp(iy ,ixh) .EQ. 1 ) .OR. &
3189  ( mtemp(iy ,ixl) .EQ. 1 ) ) &
3190  .OR. ( ( mtemp(iy+1,ixl) .EQ. 1 ) .AND. &
3191  ( ( mtemp(iy ,ixl) .EQ. 1 ) .OR. &
3192  ( mtemp(iy+1,ix ) .EQ. 1 ) ) ) &
3193  .OR. ( ( mtemp(iy+1,ixh) .EQ. 1 ) .AND. &
3194  ( ( mtemp(iy ,ixh) .EQ. 1 ) .OR. &
3195  ( mtemp(iy+1,ix ) .EQ. 1 ) ) ) &
3196  .OR. ( ( mtemp(iy-1,ixh) .EQ. 1 ) .AND. &
3197  ( ( mtemp(iy ,ixh) .EQ. 1 ) .OR. &
3198  ( mtemp(iy-1,ix ) .EQ. 1 ) ) ) &
3199  .OR. ( ( mtemp(iy-1,ixl) .EQ. 1 ) .AND. &
3200  ( ( mtemp(iy ,ixl) .EQ. 1 ) .OR. &
3201  ( mtemp(iy-1,ix ) .EQ. 1 ) ) )
3202  END DO
3203  END DO
3204  !
3205  DO ix=1, nx
3206  DO iy=1, ny
3207  IF ( mask(iy,ix) ) mtemp(iy,ix) = 1
3208  END DO
3209  END DO
3210  !
3211  END DO
3212  !
3213  ! 3.c Contract halo
3214  !
3215  ! MTMP2 = MTEMP
3216  !
3217  ! DO IIT=1, NIT
3218  !
3219  ! MASK = .FALSE.
3220  !
3221  ! DO IX=1, NX
3222  ! IXL = 1 + MOD(IX-2+NX,NX)
3223  ! IXH = 1 + MOD(IX,NX)
3224  ! DO IY=2, NY-1
3225  ! IF ( MTMP2(IY,IX) .EQ. 1 ) MASK(IY,IX) = &
3226  ! ( ( MTMP2(IY+1,IX ) .EQ. 3 ) .OR. &
3227  ! ( MTMP2(IY-1,IX ) .EQ. 3 ) .OR. &
3228  ! ( MTMP2(IY ,IXH) .EQ. 3 ) .OR. &
3229  ! ( MTMP2(IY ,IXL) .EQ. 3 ) ) &
3230  ! .OR. ( ( MTMP2(IY+1,IXL) .EQ. 3 ) .AND. &
3231  ! ( ( MTMP2(IY ,IXL) .EQ. 3 ) .OR. &
3232  ! ( MTMP2(IY+1,IX ) .EQ. 3 ) ) ) &
3233  ! .OR. ( ( MTMP2(IY+1,IXH) .EQ. 3 ) .AND. &
3234  ! ( ( MTMP2(IY ,IXH) .EQ. 3 ) .OR. &
3235  ! ( MTMP2(IY+1,IX ) .EQ. 3 ) ) ) &
3236  ! .OR. ( ( MTMP2(IY-1,IXH) .EQ. 3 ) .AND. &
3237  ! ( ( MTMP2(IY ,IXH) .EQ. 3 ) .OR. &
3238  ! ( MTMP2(IY-1,IX ) .EQ. 3 ) ) ) &
3239  ! .OR. ( ( MTMP2(IY-1,IXL) .EQ. 3 ) .AND. &
3240  ! ( ( MTMP2(IY ,IXL) .EQ. 3 ) .OR. &
3241  ! ( MTMP2(IY-1,IX ) .EQ. 3 ) ) )
3242  ! END DO
3243  ! END DO
3244  !
3245  ! DO IX=1, NX
3246  ! DO IY=1, NY
3247  ! IF ( MASK(IY,IX) ) MTMP2(IY,IX) = 3
3248  ! END DO
3249  ! END DO
3250  !
3251  ! END DO
3252  !
3253  ! 3.d Check if consistent .....
3254  !
3255  ! DO IX=1, NX
3256  ! DO IY=1, NY
3257  ! IF ( MSPLIT(IY,IX).EQ.IG .OR. MTMP2(IY,IX).EQ.1 ) THEN
3258  ! IF ( MSPLIT(IY,IX).EQ.IG .AND. MTMP2(IY,IX).NE.1 ) THEN
3259  ! write (ndst,*) ix, iy, ' in grid, not in e-c grid'
3260  ! END IF
3261  ! IF ( MSPLIT(IY,IX).NE.IG .AND. MTMP2(IY,IX).EQ.1 ) THEN
3262  ! write (ndst,*) ix, iy, ' in e-c grid, not in grid'
3263  ! END IF
3264  ! END IF
3265  ! END DO
3266  ! END DO
3267  !
3268  ! 4. Remove extraeneous MAPSTA = 2 ---------------------------------- *
3269  !
3270  DO ix=1, nx
3271  !
3272  IF ( global ) THEN
3273  ixl = 1 + mod(ix-2+nx,nx)
3274  ixh = 1 + mod(ix,nx)
3275  ELSE
3276  ixl = max( 1 , ix-1 )
3277  ixh = min( nx , ix+1 )
3278  END IF
3279  !
3280  DO iy=1, ny
3281  IF ( mtemp(iy,ix) .EQ. 2 ) THEN
3282  iyl = max( 1 , iy-1 )
3283  iyh = min( ny , iy+1 )
3284  IF ( .NOT. ( ( mtemp(iyl,ix ) .EQ. 1 ) .OR. &
3285  ( mtemp(iyh,ix ) .EQ. 1 ) .OR. &
3286  ( mtemp(iy ,ixl) .EQ. 1 ) .OR. &
3287  ( mtemp(iy ,ixh) .EQ. 1 ) ) ) &
3288  mtemp(iy,ix) = 3
3289  END IF
3290  END DO
3291  !
3292  END DO
3293  !
3294 #ifdef W3_T7
3295  WRITE (ndst,9040)
3296 #endif
3297  !
3298  ! 5. Recompute grid range ------------------------------------------- *
3299  ! Using GSTOLD to store info for modified grid
3300  !
3301 #ifdef W3_T7
3302  WRITE (ndst,9050) gstats(ig)%STRADLE, gstats(ig)%NPTS, &
3303  gstats(ig)%NXL, gstats(ig)%NXH, &
3304  gstats(ig)%NYL, gstats(ig)%NYH
3305 #endif
3306  !
3307  gstold(ig)%STRADLE = .false.
3308  gstold(ig)%NPTS = 0
3309  gstold(ig)%NXL = nx
3310  gstold(ig)%NXH = 1
3311  gstold(ig)%NYL = ny
3312  gstold(ig)%NYH = 1
3313  !
3314  IF ( global ) THEN
3315  !
3316  left = .false.
3317  right = .false.
3318  !
3319  DO iy=1, ny
3320  IF ( mtemp(iy, 1).EQ.1 .OR. mtemp(iy, 1).EQ.2 ) left = .true.
3321  IF ( mtemp(iy,nx).EQ.1 .OR. mtemp(iy,nx).EQ.2 ) right = .true.
3322  END DO
3323  gstold(ig)%STRADLE = left .AND. right
3324  !
3325  END IF
3326  !
3327  DO iy=1, ny
3328  DO ix=1, nx
3329  IF ( mtemp(iy,ix).EQ.1 .OR. mtemp(iy,ix).EQ.2 ) THEN
3330  gstold(ig)%NPTS = gstold(ig)%NPTS + 1
3331  gstold(ig)%NXL = min( gstold(ig)%NXL , ix )
3332  gstold(ig)%NXH = max( gstold(ig)%NXH , ix )
3333  gstold(ig)%NYL = min( gstold(ig)%NYL , iy )
3334  gstold(ig)%NYH = max( gstold(ig)%NYH , iy )
3335  END IF
3336  END DO
3337  END DO
3338  !
3339  IF ( gstold(ig)%STRADLE ) THEN
3340  nocnt = 0
3341  nocntm = 0
3342  nocntl = 0
3343  DO ix=1, nx
3344  there = .false.
3345  DO iy=1, ny
3346  IF ( mtemp(iy,ix).EQ.1 .OR. mtemp(iy,ix).EQ.2 ) THEN
3347  there = .true.
3348  EXIT
3349  END IF
3350  END DO
3351  IF ( there ) THEN
3352  nocnt = 0
3353  ELSE
3354  nocnt = nocnt + 1
3355  IF ( nocnt .GT. nocntm ) THEN
3356  nocntm = nocnt
3357  nocntl = ix
3358  END IF
3359  END IF
3360  END DO
3361  gstold(ig)%NXL = nocntl + 1
3362  gstold(ig)%NXH = nocntl - nocntm
3363  END IF
3364  !
3365  ! ... Make sure outside of grid is 2 or 3
3366  !
3367 #ifdef W3_T7
3368  WRITE (ndst,9051) gstold(ig)%STRADLE, gstold(ig)%NPTS, &
3369  gstold(ig)%NXL, gstold(ig)%NXH, &
3370  gstold(ig)%NYL, gstold(ig)%NYH
3371 #endif
3372  left = .false.
3373  right = .false.
3374  !
3375  DO ix=1, nx
3376  left = left .OR. ( mtemp(gstold(ig)%NYL,ix) .EQ. 1 )
3377  right = right .OR. ( mtemp(gstold(ig)%NYH,ix) .EQ. 1 )
3378  END DO
3379  !
3380  IF ( left ) gstold(ig)%NYL = gstold(ig)%NYL - 1
3381  IF ( right ) gstold(ig)%NYH = gstold(ig)%NYH + 1
3382  !
3383  DO iy=1, ny
3384  left = left .OR. ( mtemp(iy,gstold(ig)%NXL) .EQ. 1 )
3385  right = right .OR. ( mtemp(iy,gstold(ig)%NXH) .EQ. 1 )
3386  END DO
3387  !
3388  IF ( left ) gstold(ig)%NXL = gstold(ig)%NXL - 1
3389  IF ( right ) gstold(ig)%NXH = gstold(ig)%NXH + 1
3390  !
3391  IF ( global .AND. gstold(ig)%NXL.EQ.0 ) THEN
3392  gstold(ig)%NXL = nx
3393  gstold(ig)%STRADLE = .true.
3394  END IF
3395  !
3396  IF ( global .AND. gstold(ig)%NXH.EQ.nx+1 ) THEN
3397  gstold(ig)%NXH = 1
3398  gstold(ig)%STRADLE = .true.
3399  END IF
3400  !
3401 #ifdef W3_T7
3402  WRITE (ndst,9052) gstold(ig)%STRADLE, gstold(ig)%NPTS, &
3403  gstold(ig)%NXL, gstold(ig)%NXH, &
3404  gstold(ig)%NYL, gstold(ig)%NYH
3405 #endif
3406  !
3407  ! 6. Extract reduced grid data -------------------------------------- *
3408  !
3409  my = 1 + gstold(ig)%NYH - gstold(ig)%NYL
3410  mx = 1 + gstold(ig)%NXH - gstold(ig)%NXL
3411  IF ( gstold(ig)%STRADLE ) mx = mx + nx
3412  pgrid(ig)%NY = my
3413  pgrid(ig)%NX = mx
3414  pgrid(ig)%NSEA = gstold(ig)%NPTS
3415  pgrid(ig)%X0 = x0 + real(gstold(ig)%NXL-1)*sx
3416  pgrid(ig)%Y0 = y0 + real(gstold(ig)%NYL-1)*sy
3417  pgrid(ig)%SX = sx
3418  pgrid(ig)%SY = sy
3419  !
3420  xoff = 360. * real( nint((pgrid(ig)%X0+0.5*real(mx-1)*sx)/360.) )
3421  pgrid(ig)%X0 = pgrid(ig)%X0 - xoff
3422  !
3423 #ifdef W3_T7
3424  WRITE (ndst,9060) pgrid(ig)%NX, pgrid(ig)%NY, pgrid(ig)%NSEA, &
3425  pgrid(ig)%X0, pgrid(ig)%Y0, pgrid(ig)%SX, pgrid(ig)%SY
3426 #endif
3427  !
3428  ALLOCATE ( pgrid(ig)%ZBIN(mx,my) , &
3429  pgrid(ig)%OBSX(mx,my) , &
3430  pgrid(ig)%OBSY(mx,my) , &
3431  pgrid(ig)%MASK(mx,my) )
3432  !
3433  pgrid(ig)%ZBIN = zbdum
3434  pgrid(ig)%OBSX = 0.
3435  pgrid(ig)%OBSY = 0.
3436  pgrid(ig)%MASK = 99
3437  !
3438  DO ix=1, pgrid(ig)%NX
3439  jx = 1 + mod( ix+gstold(ig)%NXL-2 , nx )
3440  DO iy=1, pgrid(ig)%NY
3441  jy = iy + gstold(ig)%NYL - 1
3442  isea = mapfs(jy,jx)
3443  IF ( mtemp(jy,jx) .NE. 0 ) THEN
3444  pgrid(ig)%ZBIN(ix,iy) = zb(isea)
3445  END IF
3446  IF ( trflag .NE. 0 ) THEN
3447  pgrid(ig)%OBSX(ix,iy) = 1. - trnx(jy,jx)
3448  pgrid(ig)%OBSY(ix,iy) = 1. - trny(jy,jx)
3449  END IF
3450  pgrid(ig)%MASK(ix,iy) = mtemp(jy,jx)
3451  END DO
3452  END DO
3453  !
3454  RETURN
3455  !
3456  ! Formats
3457  !
3458 #ifdef W3_T7
3459 9000 FORMAT ( 'TEST GR1GRD: EXTRACTING GRID:',i4)
3460 9040 FORMAT ( ' MASK ON FULL GRID COMPUTED')
3461 9050 FORMAT ( 'TEST GR1GRD: GRID STATS :'/ &
3462  ' GRID MAP :',l2,2x,i8,4i5)
3463 9051 FORMAT ( ' HALO ADDED :',l2,2x,i8,4i5)
3464 9052 FORMAT ( ' BORDER ADDED :',l2,2x,i8,4i5)
3465 9060 FORMAT ( 'TEST GR1GRD: EXTRACTED GRID :',2i5,i8/ &
3466  ' ',4e12.5)
3467 #endif
3468  !
3469  !/ End of GR1GRD ----------------------------------------------------- /
3470  !/

◆ grfill()

subroutine w3gspl::grfill ( integer, intent(in)  ND)

Reassign unassigned grid points to grids, starting with the smallest grids.

Parameters
[in]NDDepth of halo for first sweep.
Author
H. L. Tolman
Date
01-Feb-2013

Definition at line 1652 of file ww3_gspl.F90.

1652  !/
1653  !/ +-----------------------------------+
1654  !/ | WAVEWATCH III NOAA/NCEP |
1655  !/ | H. L. Tolman |
1656  !/ | FORTRAN 90 |
1657  !/ | Last update : 01-Feb-2013 |
1658  !/ +-----------------------------------+
1659  !/
1660  !/ 07-Sep-2012 : Origination. ( version 4.10 )
1661  !/ 18-Sep-2012 : Include edge points of grid. ( version 4.10 )
1662  !/ Add convergence check.
1663  !/ 29-Jan-2013 : Add error code on stop. ( version 4.10 )
1664  !/ 29-Jan-2013 : Add error test output. ( version 4.10 )
1665  !/ 01-Feb-2013 : Loop over selected sea points only. ( version 4.10 )
1666  !/
1667  ! 1. Purpose :
1668  !
1669  ! Reassign unassigned grid points to grids, starting with the
1670  ! smallest grids.
1671  !
1672  ! 3. Parameters :
1673  !
1674  ! Parameter list
1675  ! ----------------------------------------------------------------
1676  ! ND Int. I Depth of halo for first sweep.
1677  ! ----------------------------------------------------------------
1678  !
1679  ! 10. Source code :
1680  !
1681  !/ ------------------------------------------------------------------- /
1682  !/
1683  IMPLICIT NONE
1684  !/
1685  !/ ------------------------------------------------------------------- /
1686  !/ Parameter list
1687  !/
1688  INTEGER, INTENT(IN) :: ND
1689  !/
1690  !/ ------------------------------------------------------------------- /
1691  !/ Local parameters
1692  !/
1693  INTEGER :: NMIN, I, NDEPTH, NITT, NADD, IXL, IXR,&
1694  NLEFT, NRIGHT, NXL, NXH, NYL, NYH
1695  INTEGER :: NXYOFF = 3
1696  INTEGER :: IIX(NSEA), IIY(NSEA), ISEA, NSEAL
1697 #ifdef W3_S
1698  INTEGER, SAVE :: IENT = 0
1699 #endif
1700  LOGICAL :: DONE(NG), MASK(NY,NX), FLOST(NG), &
1701  XFL(NX), YFL(NY)
1702  !/
1703  !/ ------------------------------------------------------------------- /
1704  !/
1705 #ifdef W3_S
1706  CALL strace (ient, 'GRFILL')
1707 #endif
1708  !
1709  ! 1. Loop to assure all reassigned ---------------------------------- *
1710  !
1711  ndepth = nd
1712  nitt = 0
1713  nleft = -1
1714  flost = .false.
1715  !
1716  nseal = 0
1717  DO ix=1, nx
1718  DO iy=1, ny
1719  IF ( msplit(iy,ix) .EQ. -1 ) THEN
1720  nseal = nseal + 1
1721  iix(nseal) = ix
1722  iiy(nseal) = iy
1723  END IF
1724  END DO
1725  END DO
1726  !
1727  DO
1728  nitt = nitt + 1
1729  !
1730  ! 2. Loop over all grids -------------------------------------------- *
1731  !
1732  done = .false.
1733  !
1734  DO j=1, ng
1735  !
1736  ! 3. Find smallest unprocessed grid --------------------------------- *
1737  !
1738  nmin = nsea + 1
1739  ig = 0
1740  !
1741  DO i=1, ng
1742  IF ( .NOT.done(i) .AND. gstats(i)%NPTS.LT.nmin ) THEN
1743  ig = i
1744  nmin = gstats(i)%NPTS
1745  END IF
1746  END DO
1747  !
1748  done(ig) = .true.
1749  !
1750 #ifdef W3_T2
1751  WRITE (ndst,9030) ig, j, nmin
1752 #endif
1753  !
1754  ! 4. Loop for halos per grid ---------------------------------------- *
1755  !
1756  DO, i=1, ndepth
1757  !
1758  mask = .false.
1759  !
1760  ! 5. Mark grid point for adding ------------------------------------- *
1761  !
1762  DO isea=1, nseal
1763  ix = iix(isea)
1764  iy = iiy(isea)
1765  ixl = 1 + mod(ix-2+nx,nx)
1766  ixr = 1 + mod(ix,nx)
1767  IF ( msplit(iy,ix) .EQ. -1 ) mask(iy,ix) = &
1768  ( msplit(iy+1,ix ) .EQ. ig ) &
1769  .OR. ( msplit(iy-1,ix ) .EQ. ig ) &
1770  .OR. ( msplit(iy ,ixr) .EQ. ig ) &
1771  .OR. ( msplit(iy ,ixl) .EQ. ig )
1772  END DO
1773  !
1774  ! 6. Add marked grid point ------------------------------------------ *
1775  !
1776  nadd = 0
1777  !
1778  DO isea=1, nseal
1779  ix = iix(isea)
1780  iy = iiy(isea)
1781  IF ( mask(iy,ix) ) THEN
1782  msplit(iy,ix) = ig
1783  nadd = nadd + 1
1784  END IF
1785  END DO
1786  !
1787  IF ( nadd .EQ. 0 ) EXIT
1788  !
1789  ! ... End loop started in 4.
1790  !
1791  END DO
1792  !
1793  ! ... End loop started in 2.
1794  !
1795  END DO
1796  !
1797  ndepth = 1
1798  !
1799  ! 7. Check convergence ---------------------------------------------- *
1800  ! 7.a Find number of points left
1801  !
1802  nright = nleft
1803  nleft = 0
1804  !
1805  DO isea=1, nseal
1806  ix = iix(isea)
1807  iy = iiy(isea)
1808  IF ( msplit(iy,ix) .EQ. -1 ) nleft = nleft + 1
1809  END DO
1810  !
1811 #ifdef W3_T2
1812  WRITE (ndst,9070) nitt, nleft
1813 #endif
1814  !
1815  ! 7.b No point left, exit loop
1816  !
1817  IF ( nleft .EQ. 0 ) EXIT
1818  !
1819  ! 7.c Stuck with points left
1820  !
1821  IF ( nright .GT. 0 ) THEN
1822  IF ( nleft .EQ. nright ) THEN
1823  !
1824  ! 7.d Do lost point correction once
1825  !
1826  IF ( .NOT. flost(ig) ) THEN
1827  CALL grlost
1828  flost(ig) = .true.
1829  ELSE
1830  !
1831  ! 7.e Got stuck for good, error message and ouput
1832  !
1833  WRITE (ndse,1000) ig, nitt, nleft
1834  !
1835  xfl = .false.
1836  yfl = .false.
1837  !
1838  DO isea=1, nseal
1839  ix = iix(isea)
1840  iy = iiy(isea)
1841  IF ( msplit(iy,ix) .EQ. -1 ) THEN
1842  xfl(max(1,ix-nxyoff):min(nx,ix+nxyoff)) = .true.
1843  yfl(max(1,iy-nxyoff):min(ny,iy+nxyoff)) = .true.
1844  END IF
1845  END DO
1846  !
1847  nxl = 0
1848  nxh = 0
1849  DO ix=1, nx
1850  IF ( xfl(ix) .AND. nxl.EQ. 0 ) nxl = ix
1851  IF ( xfl(ix) .AND. ix.EQ. nx ) nxh = ix
1852  IF ( .NOT. xfl(ix) .AND. nxl.NE. 0 ) nxh = ix-1
1853  IF ( nxh .NE. 0 ) THEN
1854  nyl = 0
1855  nyh = 0
1856  DO iy=1, ny
1857  IF ( yfl(iy) .AND. nyl.EQ. 0 ) nyl = iy
1858  IF ( yfl(iy) .AND. iy.EQ. ny ) nyh = iy
1859  IF ( .NOT. yfl(iy) .AND. nyl.NE. 0 ) &
1860  nyh = iy-1
1861  IF ( nyh .NE. 0 ) THEN
1862  WRITE (ndst,1001) nxl, nxh, nyh, nyl
1863  DO i=nyh, nyl, -1
1864  WRITE (ndst,1002) msplit(i,nxl:nxh)
1865  END DO
1866  nyl = 0
1867  nyh = 0
1868  END IF
1869  END DO
1870  nxl = 0
1871  nxh = 0
1872  END IF
1873  END DO
1874  !
1875  ! ... Stop program with error output ...
1876  !
1877  stop 01
1878  ENDIF
1879  !
1880  END IF
1881  END IF
1882  !
1883  ! ... End loop started in 1.
1884  !
1885  END DO
1886  !
1887  RETURN
1888  !
1889  ! Formats
1890  !
1891 1000 FORMAT (/' *** ERROR GRFILL : NO MORE CONVERGENCE, ', &
1892  'NITT, NLEFT:',2i8,' ***'/)
1893 1001 FORMAT ( ' MAP OUTPUT FOR GRID',i3,' AND X AND Y RANGE :',4i6/)
1894 1002 FORMAT ( ' ',60i2)
1895  !
1896 #ifdef W3_T2
1897 9030 FORMAT ( 'TEST GRFILL: PROCESSING GRID',i5,' (',i5,')',i8)
1898 9060 FORMAT ( 'TEST GRFILL: GRID, HALO, NADD :',i5,i2,i8)
1899 9070 FORMAT ( 'TEST GRFILL: NITT, NLEFT :',2i6)
1900 #endif
1901  !
1902  !/ End of GRFILL ----------------------------------------------------- /
1903  !/

References grlost().

◆ grflrg()

subroutine w3gspl::grflrg

Like GRFSML for largest grid ...

Author
H. L. Tolman
Date
29-Jan-2013

Definition at line 2925 of file ww3_gspl.F90.

2925  !/
2926  !/ +-----------------------------------+
2927  !/ | WAVEWATCH III NOAA/NCEP |
2928  !/ | H. L. Tolman |
2929  !/ | FORTRAN 90 |
2930  !/ | Last update : 29-Jan-2013 |
2931  !/ +-----------------------------------+
2932  !/
2933  !/ 19-Sep-2012 : Origination. ( version 4.10 )
2934  !/ 29-Jan-2013 : Add error code on stop. ( version 4.10 )
2935  !/
2936  ! 1. Purpose :
2937  !
2938  ! Like GRFSML for largest grid ...
2939  !
2940  ! 10. Source code :
2941  !
2942  !/ ------------------------------------------------------------------- /
2943  !/
2944  IMPLICIT NONE
2945  !/
2946  !/ ------------------------------------------------------------------- /
2947  !/ Parameter list
2948  !/
2949  !/ ------------------------------------------------------------------- /
2950  !/ Local parameters
2951  !/
2952  INTEGER :: NBIG, IGMAX(NG), NNEXT, JG
2953 !!! INTEGER :: NSMALL, IGMIN(NG), NNEXT, JG, IGADD, &
2954 !!! IGTEST, FREE(NG), NFREE, NBIG, IGB, &
2955 !!! MX, MY, NX0, NXN, NY0, NYN
2956 #ifdef W3_S
2957  INTEGER, SAVE :: IENT = 0
2958 #endif
2959  CHARACTER(LEN=1) :: NEXTTO(0:NG,0:NG), TEMP(NG)
2960  !/
2961  !/ ------------------------------------------------------------------- /
2962  !/
2963 #ifdef W3_S
2964  CALL strace (ient, 'GRFLRG')
2965 #endif
2966  !
2967  ! 1. Find big(s) ---------------------------------------------------- *
2968  !
2969  nbig = 0
2970  igmax = 0
2971  !
2972  DO ig=1,ng
2973  IF ( gstats(ig)%INSTAT .AND. &
2974  gstats(ig)%NPTS .EQ. mstats%NMAX ) THEN
2975  nbig = nbig + 1
2976  igmax(nbig) = ig
2977  END IF
2978  END DO
2979  !
2980 #ifdef W3_T6
2981  WRITE (ndst,9010) nbig, igmax(:nbig)
2982 #endif
2983  !
2984  ! 2. Find neighbours ------------------------------------------------ *
2985  !
2986  nextto = '.'
2987  !
2988  DO ix=1, nx-1
2989  DO iy=1, ny-1
2990  nextto(msplit(iy ,ix ),msplit(iy+1,ix )) = 'X'
2991  nextto(msplit(iy+1,ix ),msplit(iy ,ix )) = 'X'
2992  nextto(msplit(iy ,ix+1),msplit(iy ,ix )) = 'X'
2993  nextto(msplit(iy ,ix ),msplit(iy ,ix+1)) = 'X'
2994  END DO
2995  END DO
2996  !
2997  IF ( global ) THEN
2998  DO iy=1, ny-1
2999  nextto(msplit(iy ,nx),msplit(iy+1,nx)) = 'X'
3000  nextto(msplit(iy+1,nx),msplit(iy ,nx)) = 'X'
3001  nextto(msplit(iy , 1),msplit(iy ,nx)) = 'X'
3002  nextto(msplit(iy ,nx),msplit(iy , 1)) = 'X'
3003  END DO
3004  END IF
3005  !
3006  DO ig=0,ng
3007  nextto(ig,ig) = '-'
3008  END DO
3009  !
3010 #ifdef W3_T6
3011  WRITE (ndst,9020)
3012  DO ig=1, ng
3013  temp = nextto(ig,1:)
3014  WRITE (ndst,9021) ig, temp
3015  END DO
3016 #endif
3017  !
3018  ! 3. Loop over big grids -------------------------------------------- *
3019  !
3020  DO j=1, nbig
3021  !
3022 #ifdef W3_T6
3023  WRITE (ndst,9030) igmax(j)
3024 #endif
3025  !
3026  ! 3.a Find neighbours
3027  !
3028  ig = igmax(j)
3029  nnext = 0
3030  DO jg=1, ng
3031  IF ( nextto(ig,jg) .EQ. 'X' ) nnext = nnext + 1
3032  END DO
3033  !
3034 #ifdef W3_T6
3035  WRITE (ndst,9031) nnext
3036 #endif
3037  !
3038  ! 3.b Enough neighbours found, mark as 'not to be processed further'
3039  !
3040  IF ( nnext .GE. 1 ) THEN
3041  gstats(ig)%INSTAT = .false.
3042 #ifdef W3_T6
3043  WRITE (ndst,9032)
3044 #endif
3045  ELSE
3046  !
3047  ! 3.c Biggest grid is isolated, should split
3048  !
3049  WRITE (ndse,930)
3050  stop 11
3051  !
3052  END IF
3053  !
3054  END DO
3055  !
3056  RETURN
3057  !
3058  ! Formats
3059  !
3060 930 FORMAT ( ' *** ERROR GRFLRG: LARGEST GRID IS ISOLATED ***' &
3061  ' SPLITTING NOT YET IMPLEMENTED '/)
3062  !
3063 #ifdef W3_T6
3064 9010 FORMAT ( 'TEST GRFLRG:',i2,' BIG GRIDS:',10i4)
3065 9020 FORMAT ( 'TEST GRFLRG: NEIGHBOUR MAP PER GRID')
3066 9021 FORMAT (2x,i3,2x,120a1)
3067 9030 FORMAT ( 'TEST GRFLRG: PROCESSING BIG GRID',i4)
3068 9031 FORMAT ( ' GRID HAS',i3,' NEIGHBOURS')
3069 9032 FORMAT ( ' NO ACTION')
3070 #endif
3071  !
3072  !/ End of GRFLRG ----------------------------------------------------- /
3073  !/

◆ grfsml()

subroutine w3gspl::grfsml

Subroutine called when lowest grid size is stuck.

Attempting to join to neighbor grid, otherwise mark as a accepted small grid. note that a small grid does not influence parallel scaling like a big grid does .....

Author
H. L. Tolman
Date
04-Feb-2013

Definition at line 2639 of file ww3_gspl.F90.

2639  !/
2640  !/ +-----------------------------------+
2641  !/ | WAVEWATCH III NOAA/NCEP |
2642  !/ | H. L. Tolman |
2643  !/ | FORTRAN 90 |
2644  !/ | Last update : 04-Feb-2013 |
2645  !/ +-----------------------------------+
2646  !/
2647  !/ 13-Sep-2012 : Origination. ( version 4.10 )
2648  !/ 04-Feb-2013 : Bug fix grid splitting. ( version 4.10 )
2649  !/
2650  ! 1. Purpose :
2651  !
2652  ! Subroutine called when lowest grid size is stuck. Attempting to
2653  ! joint to neighbor grid, otherwise mark as accepted small grid.
2654  ! note that small grid does not influence parallel scaling like a
2655  ! big grid does .....
2656  !
2657  ! 1-Feb-2013: Also used for early small-grid merging.
2658  !
2659  ! 10. Source code :
2660  !
2661  !/ ------------------------------------------------------------------- /
2662  !/
2663  IMPLICIT NONE
2664  !/
2665  !/ ------------------------------------------------------------------- /
2666  !/ Parameter list
2667  !/
2668  !/ ------------------------------------------------------------------- /
2669  !/ Local parameters
2670  !/
2671  INTEGER :: NSMALL, IGMIN(NG), NNEXT, JG, IGADD, &
2672  IGTEST, FREE(NG), NFREE, NBIG, IGB, &
2673  MX, MY, NX0, NXN, NY0, NYN, JX
2674 #ifdef W3_T5
2675  INTEGER :: NXNT
2676 #endif
2677 #ifdef W3_S
2678  INTEGER, SAVE :: IENT = 0
2679 #endif
2680  CHARACTER(LEN=1) :: NEXTTO(0:NG,0:NG), TEMP(NG)
2681  !/
2682  !/ ------------------------------------------------------------------- /
2683  !/
2684 #ifdef W3_S
2685  CALL strace (ient, 'GRFSML')
2686 #endif
2687  !
2688  ! 1. Find small(s) -------------------------------------------------- *
2689  !
2690  nsmall = 0
2691  igmin = 0
2692  !
2693  DO ig=1,ng
2694  IF ( gstats(ig)%INSTAT .AND. &
2695  gstats(ig)%NPTS .EQ. mstats%NMIN ) THEN
2696  nsmall = nsmall + 1
2697  igmin(nsmall) = ig
2698  END IF
2699  END DO
2700  !
2701 #ifdef W3_T5
2702  WRITE (ndst,9010) nsmall, igmin(:nsmall)
2703 #endif
2704  !
2705  ! 2. Find neighbours ------------------------------------------------ *
2706  !
2707  nextto = '.'
2708  !
2709  DO ix=1, nx-1
2710  DO iy=1, ny-1
2711  nextto(msplit(iy ,ix ),msplit(iy+1,ix )) = 'X'
2712  nextto(msplit(iy+1,ix ),msplit(iy ,ix )) = 'X'
2713  nextto(msplit(iy ,ix+1),msplit(iy ,ix )) = 'X'
2714  nextto(msplit(iy ,ix ),msplit(iy ,ix+1)) = 'X'
2715  END DO
2716  END DO
2717  !
2718  IF ( global ) THEN
2719  DO iy=1, ny-1
2720  nextto(msplit(iy ,nx),msplit(iy+1,nx)) = 'X'
2721  nextto(msplit(iy+1,nx),msplit(iy ,nx)) = 'X'
2722  nextto(msplit(iy , 1),msplit(iy ,nx)) = 'X'
2723  nextto(msplit(iy ,nx),msplit(iy , 1)) = 'X'
2724  END DO
2725  END IF
2726  !
2727  DO ig=0,ng
2728  nextto(ig,ig) = '-'
2729  END DO
2730  !
2731 #ifdef W3_T5
2732  WRITE (ndst,9020)
2733  DO ig=1, ng
2734  temp = nextto(ig,1:)
2735  WRITE (ndst,9021) ig, temp
2736  END DO
2737 #endif
2738  !
2739  ! 3. Loop over small grids ------------------------------------------ *
2740  !
2741  free = 0
2742  nfree = 0
2743  !
2744  DO j=1, nsmall
2745  !
2746 #ifdef W3_T5
2747  WRITE (ndst,9030) igmin(j)
2748 #endif
2749  !
2750  ! 3.a Find neighbours
2751  !
2752  ig = igmin(j)
2753  igadd = 0
2754  igtest = nsea + 1
2755  nnext = 0
2756  DO jg=1, ng
2757  IF ( nextto(ig,jg) .EQ. 'X' ) THEN
2758  nnext = nnext + 1
2759  IF ( gstats(jg)%NPTS .LT. igtest ) THEN
2760  igtest = gstats(jg)%NPTS
2761  igadd = jg
2762  END IF
2763  END IF
2764  END DO
2765  !
2766 #ifdef W3_T5
2767  WRITE (ndst,9031) nnext
2768 #endif
2769  !
2770  ! 3.b No neighbours found, mark as 'not to be processed further'
2771  !
2772  IF ( nnext .EQ. 0 ) THEN
2773  gstats(ig)%INSTAT = .false.
2774 #ifdef W3_T5
2775  WRITE (ndst,9032) ig
2776 #endif
2777  ELSE
2778  !
2779  ! 3.c Check smallest neighbor
2780  !
2781 #ifdef W3_T5
2782  WRITE (ndst,9033) igadd, igtest, igtest+ingmin, nint(xmean)
2783 #endif
2784  !
2785  IF ( igtest + ingmin .LT. nint(xmean) ) THEN
2786  !
2787  ! ... Merge grids
2788  !
2789  DO ix=1, nx
2790  DO iy=1, ny
2791  IF ( msplit(iy,ix) .EQ. ig ) msplit(iy,ix) = igadd
2792  END DO
2793  END DO
2794  !
2795  nfree = nfree + 1
2796  free(nfree) = ig
2797  !
2798  ELSE
2799  !
2800  ! ... Remove grid(s) from stats
2801  !
2802 #ifdef W3_T5
2803  WRITE (ndst,9034)
2804 #endif
2805  !
2806  gstats(ig)%INSTAT = .false.
2807 #ifdef W3_T5
2808  WRITE (ndst,9032) ig
2809 #endif
2810  nnext = 0
2811  DO jg=1, ng
2812  IF ( nextto(igadd,jg) .EQ. 'X' ) nnext = nnext + 1
2813  END DO
2814  IF ( nnext .EQ. 1 ) THEN
2815  gstats(igadd)%INSTAT = .false.
2816 #ifdef W3_T5
2817  WRITE (ndst,9032) igadd
2818 #endif
2819  END IF
2820  !
2821  END IF
2822  !
2823  END IF
2824  !
2825  END DO
2826  !
2827  ! 4. Make new grids as needed --------------------------------------- *
2828  !
2829 #ifdef W3_T5
2830  WRITE (ndst,9040) nfree
2831 #endif
2832  !
2833  DO j=1, nfree
2834  !
2835 #ifdef W3_T5
2836  WRITE (ndst,9041) free(j)
2837 #endif
2838  !
2839  ! 4.a Find biggest grid
2840  !
2841  nbig = 0
2842  igb = 0
2843  !
2844  DO ig=1, ng
2845  IF ( gstats(ig)%NPTS .GT. nbig ) THEN
2846  nbig = gstats(ig)%NPTS
2847  igb = ig
2848  END IF
2849  END DO
2850  !
2851  ! 4.a Split biggest grid
2852  !
2853  nx0 = gstats(igb)%NXL
2854  nxn = gstats(igb)%NXH
2855  ny0 = gstats(igb)%NYL
2856  nyn = gstats(igb)%NYH
2857  !
2858  my = 1 + gstats(igb)%NYH - gstats(igb)%NYL
2859  mx = 1 + gstats(igb)%NXH - gstats(igb)%NXL
2860  IF ( gstats(igb)%STRADLE ) mx = mx + nx
2861  !
2862  IF ( my .GE. mx ) THEN
2863 #ifdef W3_T5
2864  WRITE (ndst,9042) igb, 'VERTICAL', mx, my
2865 #endif
2866  nyn = ny0 + my/2
2867  ELSE
2868 #ifdef W3_T5
2869  WRITE (ndst,9042) igb, 'HORIZONTAL', mx, my
2870 #endif
2871  nxn = nx0 + mx/2
2872 #ifdef W3_T5
2873  nxnt = 1 + mod(nxn-1,nx)
2874 #endif
2875  END IF
2876 #ifdef W3_T5
2877  WRITE (ndst,9043) gstats(igb)%NXL, gstats(igb)%NXH, &
2878  gstats(igb)%NYL, gstats(igb)%NYH, &
2879  gstats(igb)%STRADLE, nx0, nxn, ny0, nyn
2880 #endif
2881  !
2882  DO ix=nx0, nxn
2883  jx = 1 + mod(ix-1,nx)
2884  DO iy=ny0, nyn
2885  IF ( msplit(iy,jx) .EQ. igb ) msplit(iy,jx) = free(j)
2886  END DO
2887  END DO
2888  !
2889  gstats(igb)%NPTS = 0
2890  gstats(free(j))%NPTS = 0
2891  !
2892  END DO
2893  !
2894  RETURN
2895  !
2896  ! Formats
2897  !
2898 #ifdef W3_T5
2899 9010 FORMAT ( 'TEST GRFSML:',i2,' SMALL GRIDS:',10i4)
2900 9020 FORMAT ( 'TEST GRFSML: NEIGHBOUR MAP PER GRID')
2901 9021 FORMAT (2x,i3,2x,120a1)
2902 9030 FORMAT ( 'TEST GRFSML: PROCESSING SMALL GRID',i4)
2903 9031 FORMAT ( ' GRID HAS',i3,' NEIGHBOURS')
2904 9032 FORMAT ( ' REMOVED GRID',i4,' FROM STATS')
2905 9033 FORMAT ( ' SMALLEST NEIGHBOUR AND SIZE',i4,i6/ &
2906  ' SIZE OF COMBINED GRIDS',i8,' (',i8,')')
2907 9034 FORMAT ( ' GRIDS TOO LARGE TO MERGE')
2908 9040 FORMAT ( 'TEST GRFSML: GENERATING',i3,' NEW GRIDS')
2909 9041 FORMAT ( ' MAKING GRID NR.:',i4)
2910 9042 FORMAT ( ' SPLITTING GRID',i3,' ',a,', MX,MY:',2i6)
2911 9043 FORMAT ( ' OLD RANGE :',4i6,l4/ &
2912  ' NEW RANGE :',4i6)
2913 #endif
2914  !
2915  !/ End of GRFSML ----------------------------------------------------- /
2916  !/

◆ grinfo()

subroutine w3gspl::grinfo

Compile statistical info on all sub grids (no halo).

Author
H. L. Tolman
Date
13-Sep-2012

Definition at line 1302 of file ww3_gspl.F90.

1302  !/
1303  !/ +-----------------------------------+
1304  !/ | WAVEWATCH III NOAA/NCEP |
1305  !/ | H. L. Tolman |
1306  !/ | FORTRAN 90 |
1307  !/ | Last update : 13-Sep-2012 |
1308  !/ +-----------------------------------+
1309  !/
1310  !/ 06-Sep-2012 : Origination. ( version 4.10 )
1311  !/ 13-Sep-2012 : Option to exclude grids from stats. ( version 4.10 )
1312  !/
1313  ! 1. Purpose :
1314  !
1315  ! Compile statistical info on all sub grids (no halo).
1316  !
1317  ! 10. Source code :
1318  !
1319  !/ ------------------------------------------------------------------- /
1320  !/
1321  IMPLICIT NONE
1322  !/
1323  !/ ------------------------------------------------------------------- /
1324  !/ Parameter list
1325  !/
1326  !/ ------------------------------------------------------------------- /
1327  !/ Local parameters
1328  !/
1329  INTEGER :: NOCNT, NOCNTM, NOCNTL, NGC, NSEAC
1330 #ifdef W3_S
1331  INTEGER, SAVE :: IENT = 0
1332 #endif
1333  REAL :: SUMSQR
1334  LOGICAL :: LEFT, RIGHT, THERE
1335  !/
1336  !/ ------------------------------------------------------------------- /
1337  !/
1338 #ifdef W3_S
1339  CALL strace (ient, 'GRINFO')
1340 #endif
1341  !
1342  ! 1. Initialization ------------------------------------------------- *
1343  !
1344  gstats(:)%STRADLE = .false.
1345  gstats(:)%NPTS = 0
1346  gstats(:)%NXL = nx
1347  gstats(:)%NXH = 1
1348  gstats(:)%NYL = ny
1349  gstats(:)%NYH = 1
1350  !
1351  ! 2. Get STRADLE, NGC ----------------------------------------------- *
1352  !
1353  ngc = 0
1354  !
1355  DO ig=1, ng
1356  left = .false.
1357  right = .false.
1358  IF ( gstats(ig)%INSTAT ) ngc = ngc + 1
1359  DO iy=1, ny
1360  IF ( msplit(iy, 1) .EQ. ig ) left = .true.
1361  IF ( msplit(iy,nx) .EQ. ig ) right = .true.
1362  END DO
1363  gstats(ig)%STRADLE = left .AND. right
1364  END DO
1365  !
1366  IF ( ngc .EQ. 0 ) THEN
1367  ngc = 1
1368  done = .true.
1369  END IF
1370  !
1371  ! 3. Run grid stats ------------------------------------------------- *
1372  ! 3.a General
1373  !
1374  DO iy=1, ny
1375  DO ix=1, nx
1376  ig = msplit(iy,ix)
1377  IF ( msplit(iy,ix) .GT. 0 ) THEN
1378  gstats(ig)%NPTS = gstats(ig)%NPTS + 1
1379  gstats(ig)%NXL = min( gstats(ig)%NXL , ix )
1380  gstats(ig)%NXH = max( gstats(ig)%NXH , ix )
1381  gstats(ig)%NYL = min( gstats(ig)%NYL , iy )
1382  gstats(ig)%NYH = max( gstats(ig)%NYH , iy )
1383  END IF
1384  END DO
1385  END DO
1386  !
1387  ! 3.b Stradled grids
1388  !
1389  IF ( ng .GT. 1) THEN
1390  DO ig=1, ng
1391  IF ( gstats(ig)%STRADLE ) THEN
1392  nocnt = 0
1393  nocntm = 0
1394  nocntl = 0
1395  DO ix=1, nx
1396  there = .false.
1397  DO iy=1, ny
1398  IF ( msplit(iy,ix) .EQ. ig ) THEN
1399  there = .true.
1400  EXIT
1401  END IF
1402  END DO
1403  IF ( there ) THEN
1404  nocnt = 0
1405  ELSE
1406  nocnt = nocnt + 1
1407  IF ( nocnt .GT. nocntm ) THEN
1408  nocntm = nocnt
1409  nocntl = ix
1410  END IF
1411  END IF
1412  END DO
1413  gstats(ig)%NXL = nocntl + 1
1414  gstats(ig)%NXH = nocntl - nocntm
1415  END IF
1416  END DO
1417  ELSE
1418  gstats(1)%STRADLE = .false.
1419  END IF
1420  !
1421  ! 3.c Corrected NSEA
1422  !
1423  nseac = 0
1424  !
1425  DO ig=1, ng
1426  IF ( gstats(ig)%INSTAT ) nseac = nseac + gstats(ig)%NPTS
1427  END DO
1428  !
1429  ! 4. Run overall stats ---------------------------------------------- *
1430  !
1431  mstats%NMIN = nsea + 1
1432  mstats%NMAX = 0
1433  xmean = real(nseac) / real(ngc)
1434  sumsqr = 0.
1435  !
1436  DO ig=1, ng
1437  IF ( .NOT. gstats(ig)%INSTAT ) cycle
1438  mstats%NMIN = min( mstats%NMIN , gstats(ig)%NPTS )
1439  mstats%NMAX = max( mstats%NMAX , gstats(ig)%NPTS )
1440  sumsqr = sumsqr + ( real(gstats(ig)%NPTS) - xmean )**2
1441  END DO
1442  !
1443  mstats%RSTD = sqrt( sumsqr / real(ngc) )
1444  !
1445  ! 5. Test output ---------------------------------------------------- *
1446  !
1447 #ifdef W3_T1
1448  WRITE (ndst,9000)
1449  DO ig=1, ng
1450  WRITE (ndst,9001) ig, gstats(ig)%STRADLE, gstats(ig)%NPTS, &
1451  gstats(ig)%NXL, gstats(ig)%NXH, &
1452  gstats(ig)%NYL, gstats(ig)%NYH
1453  END DO
1454  WRITE (ndst,9010) mstats%NMIN, mstats%NMAX, mstats%RSTD
1455 #endif
1456  !
1457  RETURN
1458  !
1459  ! Formats
1460  !
1461 #ifdef W3_T1
1462 9000 FORMAT ( 'TEST GRINFO: J, STRADLE, NPTS,NXL-H, NYL-H')
1463 9001 FORMAT ( ' ',i4,2x,l1,2x,i7,4i5)
1464 9010 FORMAT ( 'TEST GRINFO: MIN, MAX, STD:',2i8,f10.2)
1465 #endif
1466  !
1467  !/ End of GRINFO ----------------------------------------------------- /
1468  !/

◆ grlost()

subroutine w3gspl::grlost

Reassign unassigned grid points to grids.

Dealing with lost point by finding closest grids.

Author
H. L. Tolman
Date
09-Jan-2013

Definition at line 1913 of file ww3_gspl.F90.

1913  !/
1914  !/ +-----------------------------------+
1915  !/ | WAVEWATCH III NOAA/NCEP |
1916  !/ | H. L. Tolman |
1917  !/ | FORTRAN 90 |
1918  !/ | Last update : .9-Jan-2013 |
1919  !/ +-----------------------------------+
1920  !/
1921  !/ 31-Jan-2013 : Origination. ( version 4.10 )
1922  !/
1923  ! 1. Purpose :
1924  !
1925  ! Reassign unassigned grid points to grids. Dealing with lost
1926  ! point by finding closest grids.
1927  !
1928  ! 10. Source code :
1929  !
1930  !/ ------------------------------------------------------------------- /
1931  !/
1932  IMPLICIT NONE
1933  !/
1934  !/ ------------------------------------------------------------------- /
1935  !/ Local parameters
1936  !/
1937  INTEGER :: IX, IY, IOFF, JJX, JX, JY, IG, I
1938 #ifdef W3_S
1939  INTEGER, SAVE :: IENT = 0
1940 #endif
1941  INTEGER :: IFOUND(-1:NG)
1942  !/
1943  !/ ------------------------------------------------------------------- /
1944  !/
1945 #ifdef W3_S
1946  CALL strace (ient, 'GRLOST')
1947 #endif
1948  !
1949  ! 1. Loop over all grid points -------------------------------------- *
1950  !
1951  DO ix=1, nx
1952  DO iy=1, ny
1953  !
1954  IF ( msplit(iy,ix) .EQ. -1 ) THEN
1955  !
1956  ! 2. Find nearest grid(s) ------------------------------------------- *
1957  !
1958  ioff = 1
1959  !
1960  DO
1961  !
1962  ifound = 0
1963  DO jjx=ix-ioff, ix+ioff
1964  IF ( global ) THEN
1965  jx = 1 + mod(jjx-1+2*nx,nx)
1966  ELSE
1967  jx = jjx
1968  END IF
1969  IF ( jx.LT.1 .OR. jx.GT.nx ) cycle
1970  DO jy=iy-ioff, iy+ioff
1971  IF ( jy.LT.1 .OR. jy.GT.ny ) cycle
1972  ifound(msplit(jy,jx)) = ifound(msplit(jy,jx)) + 1
1973  END DO
1974  END DO
1975  !
1976  ig = 0
1977  DO i=1, ng
1978  IF ( ifound(i) .GT. 0 ) THEN
1979  ig = i
1980  EXIT
1981  END IF
1982  END DO
1983  !
1984  IF ( ig .NE. 0 ) THEN
1985  msplit(iy,ix) = ig
1986  EXIT
1987  END IF
1988  !
1989  ioff = ioff + 1
1990  IF ( ioff .GT. nx .AND. ioff.GT.ny ) EXIT
1991  END DO
1992  !
1993  ! ... End of loops and logic started in 1.
1994  !
1995  END IF
1996  !
1997  END DO
1998  END DO
1999  !
2000  RETURN
2001  !
2002  ! Formats
2003  !
2004  !/ End of GRLOST ----------------------------------------------------- /
2005  !/

Referenced by grfill().

◆ grsepa()

subroutine w3gspl::grsepa ( logical, intent(inout)  OK,
real, intent(in)  FRAC 
)

Remove smaller grid parts.

Remove smaller parts of a grid that are separated from the main body, and that can be attached to other grids.

Parameters
[in,out]OK
[in]FRACFraction of average size used to remove grid part.
Author
H. L. Tolman
Date
01-Feb-2013

Definition at line 2294 of file ww3_gspl.F90.

2294  !/
2295  !/ +-----------------------------------+
2296  !/ | WAVEWATCH III NOAA/NCEP |
2297  !/ | H. L. Tolman |
2298  !/ | FORTRAN 90 |
2299  !/ | Last update : 01-Feb-2013 |
2300  !/ +-----------------------------------+
2301  !/
2302  !/ 10-Sep-2012 : Origination. ( version 4.10 )
2303  !/ 18-Sep-2012 : Include edge points of grid. ( version 4.10 )
2304  !/ 01-Feb-2013 : Much faster algorithms. ( version 4.10 )
2305  !/
2306  ! 1. Purpose :
2307  !
2308  ! Remove smller parts of a grid that are separated from the main
2309  ! body, and that can be attached to other grids.
2310  !
2311  ! 3. Parameters :
2312  !
2313  ! Parameter list
2314  ! ----------------------------------------------------------------
2315  ! OK Log. I/O Flag for grid status, .F. if values of
2316  ! -1 are left in MSPLIT.
2317  ! FRAC Real I Fraction of average size used to remove grid
2318  ! part.
2319  ! ----------------------------------------------------------------
2320  !
2321  ! 10. Source code :
2322  !
2323  !/ ------------------------------------------------------------------- /
2324  !/
2325  IMPLICIT NONE
2326  !/
2327  !/ ------------------------------------------------------------------- /
2328  !/ Parameter list
2329  !/
2330  REAL, INTENT(IN) :: FRAC
2331  LOGICAL, INTENT(INOUT) :: OK
2332  !/
2333  !/ ------------------------------------------------------------------- /
2334  !/ Local parameters
2335  !/
2336  INTEGER :: IPAVG, IPCHCK, ID, IPTOT, IX, IY, &
2337  IXL, IYL, IDL, JX, JY, KY, IPT, &
2338  IXH, IYH, I, J, K, L, IMIN, LMIN
2339 #ifdef W3_S
2340  INTEGER, SAVE :: IENT = 0
2341 #endif
2342  INTEGER :: GMASK(NY,NX), IIX(NSEA), IIY(NSEA)
2343  INTEGER, ALLOCATABLE :: PMAP(:), INGRD(:)
2344  LOGICAL :: PREV
2345  LOGICAL,ALLOCATABLE :: FLNEXT(:), NEXTTO(:,:)
2346  !/
2347  !/ ------------------------------------------------------------------- /
2348  !/
2349 #ifdef W3_S
2350  CALL strace (ient, 'GRSEPA')
2351 #endif
2352  !
2353  ipavg = nint( real(nsea) / real(ng) )
2354  ipchck = nint( frac * real(nsea) / real(ng) )
2355  !
2356 #ifdef W3_T4
2357  WRITE (ndst,9000) ipavg, ipchck
2358 #endif
2359  !
2360  ! 1. Loop over grids ------------------------------------------------ *
2361  !
2362  DO ig=1, ng
2363  !
2364  gmask = 0
2365  id = 0
2366  !
2367 #ifdef W3_T4
2368  WRITE (ndst,9010) ig
2369 #endif
2370  !
2371  ! 2. Find all parts ------------------------------------------------- *
2372  ! 2.a First loop, partial parts
2373  !
2374  iptot = 0
2375  !
2376  DO ix=1, nx
2377  !
2378  ixl = 1 + mod(ix-2+nx,nx)
2379  prev = .false.
2380  !
2381  DO iy=1, ny
2382  IF (msplit(iy,ix) .EQ. ig ) THEN
2383  iptot = iptot + 1
2384  iix(iptot) = ix
2385  iiy(iptot) = iy
2386  IF ( .NOT. prev) THEN
2387  id = id + 1
2388  prev = .true.
2389  END IF
2390  gmask(iy,ix) = id
2391  ELSE IF ( prev ) THEN
2392  prev = .false.
2393  idl = 0
2394  DO jy=iy-1, 1, -1
2395  IF ( gmask(jy,ix) .EQ. 0 ) EXIT
2396  IF ( gmask(jy,ixl).NE.0 .AND. idl.EQ.0 ) &
2397  idl = gmask(jy,ixl)
2398  END DO
2399  IF ( idl .NE. 0 ) THEN
2400  DO ky=jy+1, iy-1
2401  IF ( gmask(ky,ix).EQ.id ) gmask(ky,ix) = idl
2402  END DO
2403  id = id - 1
2404  END IF
2405  !
2406  END IF
2407  END DO
2408  END DO
2409  !
2410  ! 2.b Grid too small, do not cut
2411  !
2412  IF ( iptot .LE. ipavg ) THEN
2413 #ifdef W3_T4
2414  WRITE (ndst,9020) iptot, ipavg
2415 #endif
2416  cycle
2417  END IF
2418  !
2419  ! 2.c Neighbouring grid parts
2420  ! Raw data
2421  !
2422  ALLOCATE ( nextto(0:id,0:id), pmap(0:id) )
2423  nextto = .false.
2424  !
2425  DO ipt=1, iptot
2426  ix = iix(ipt)
2427  iy = iiy(ipt)
2428  ixl = 1 + mod(ix-2+nx,nx)
2429  iyl = iy - 1
2430  ixh = 1 + mod(ix,nx)
2431  iyh = iy + 1
2432  nextto( gmask(iy,ix) , gmask(iy ,ixl) ) = .true.
2433  nextto( gmask(iy,ix) , gmask(iy ,ixh) ) = .true.
2434  nextto( gmask(iy,ix) , gmask(iyl,ix ) ) = .true.
2435  nextto( gmask(iy,ix) , gmask(iyh,ix ) ) = .true.
2436  END DO
2437  !
2438  ! Make symmetric
2439  !
2440  DO i=1, id
2441  DO j=1, id
2442  nextto(i,j) = nextto(i,j) .OR. nextto(j,i)
2443  END DO
2444  END DO
2445  !
2446  ! Connect accross neighbours
2447  !
2448  DO i=1, id
2449  DO j=1, id
2450  IF ( nextto(i,j) ) THEN
2451  DO k=1, id
2452  IF ( nextto(k,j) ) THEN
2453  nextto(k,i) = .true.
2454  nextto(i,k) = .true.
2455  END IF
2456  END DO
2457  END IF
2458  END DO
2459  END DO
2460  !
2461  ! Map the parts
2462  !
2463  idl = id
2464  pmap = 0
2465  id = 0
2466  !
2467  DO i=1, idl
2468  IF ( pmap(i) .EQ. 0 ) THEN
2469  id = id + 1
2470  DO j=1, idl
2471  IF ( nextto(j,i) ) EXIT
2472  END DO
2473  IF ( j .GT. idl ) THEN
2474  pmap(i) = id
2475  ELSE
2476  DO k=i, idl
2477  IF ( pmap(k).EQ.0 .AND. nextto(j,k) ) pmap(k) = id
2478  END DO
2479  END IF
2480  END IF
2481  END DO
2482  !
2483  DEALLOCATE ( nextto )
2484  !
2485  ! 3. Grid is contiguous --------------------------------------------- *
2486  !
2487  IF ( id .EQ. 1 ) THEN
2488 #ifdef W3_T4
2489  WRITE (ndst,9030) ig
2490 #endif
2491  DEALLOCATE ( pmap )
2492  cycle
2493  END IF
2494  !
2495  ! 4. Grid is split, get stats --------------------------------------- *
2496  !
2497 #ifdef W3_T4
2498  WRITE (ndst,9040) ig
2499 #endif
2500  !
2501  ! 4.a Construct final map for grid
2502  !
2503  DO ipt=1, iptot
2504  ix = iix(ipt)
2505  iy = iiy(ipt)
2506  gmask(iy,ix) = pmap(gmask(iy,ix))
2507  END DO
2508  !
2509  DEALLOCATE ( pmap )
2510  !
2511  ! 4.b Run stats
2512  !
2513  ALLOCATE ( ingrd(id), flnext(id) )
2514  ingrd = 0
2515  flnext = .false.
2516  iptot = 0
2517  !
2518  DO jx=1, nx
2519  DO jy=1, ny
2520  IF ( gmask(jy,jx) .GT. 0 ) THEN
2521  ingrd(gmask(jy,jx)) = ingrd(gmask(jy,jx)) + 1
2522  iptot = iptot + 1
2523  END IF
2524  END DO
2525  END DO
2526  !
2527  DO jx=1, nx
2528  DO jy=1, ny-1
2529  IF ( ( gmask(jy ,jx) .GT. 0 ) .AND. &
2530  ( sea(jy+1,jx) .AND. msplit(jy+1,jx).NE.ig ) ) &
2531  flnext(gmask(jy ,jx)) = .true.
2532  IF ( ( gmask(jy+1,jx) .GT. 0 ) .AND. &
2533  ( sea(jy ,jx) .AND. msplit(jy ,jx).NE.ig ) ) &
2534  flnext(gmask(jy+1,jx)) = .true.
2535  END DO
2536  END DO
2537  !
2538  DO jy=1, ny
2539  DO jx=1, nx-1
2540  IF ( ( gmask(jy,jx ) .GT. 0 ) .AND. &
2541  ( sea(jy,jx+1) .AND. msplit(jy,jx+1).NE.ig ) ) &
2542  flnext(gmask(jy,jx )) = .true.
2543  IF ( ( gmask(jy,jx+1) .GT. 0 ) .AND. &
2544  ( sea(jy,jx ) .AND. msplit(jy,jx ).NE.ig ) ) &
2545  flnext(gmask(jy,jx+1)) = .true.
2546  END DO
2547  IF ( global ) THEN
2548  IF ( ( gmask(jy,nx) .GT. 0 ) .AND. &
2549  ( sea(jy, 1) .AND. msplit(jy, 1).NE.ig ) ) &
2550  flnext(gmask(jy,nx)) = .true.
2551  IF ( ( gmask(jy, 1) .GT. 0 ) .AND. &
2552  ( sea(jy,nx) .AND. msplit(jy,nx).NE.ig ) ) &
2553  flnext(gmask(jy, 1)) = .true.
2554  END IF
2555  END DO
2556  !
2557 #ifdef W3_T4
2558  DO j=1, id
2559  WRITE (ndst,9041) j, ingrd(j), flnext(j)
2560  END DO
2561 #endif
2562  !
2563  ! 5. Grid large enough, find smallest part -------------------------- *
2564  !
2565  imin = nsea
2566  lmin = 0
2567  !
2568  DO j=1, id
2569  IF ( flnext(j) .AND. ingrd(j).LT.imin ) THEN
2570  imin = ingrd(j)
2571  lmin = j
2572  END IF
2573  END DO
2574  !
2575  IF ( lmin .EQ. 0 ) THEN
2576 #ifdef W3_T4
2577  WRITE (ndst,9050)
2578 #endif
2579  DEALLOCATE ( ingrd, flnext )
2580  cycle
2581  END IF
2582  !
2583  IF ( imin .GT. ipchck ) THEN
2584 #ifdef W3_T4
2585  WRITE (ndst,9051)
2586 #endif
2587  DEALLOCATE ( ingrd, flnext )
2588  cycle
2589  END IF
2590  !
2591  ! 6. Part to cut has been identified -------------------------------- *
2592  !
2593 #ifdef W3_T4
2594  WRITE (ndst,9060) lmin
2595 #endif
2596  !
2597  DO jx=1, nx
2598  DO jy=1, ny
2599  IF ( gmask(jy,jx) .EQ. lmin ) msplit(jy,jx) = -1
2600  END DO
2601  END DO
2602  !
2603  DEALLOCATE ( ingrd, flnext )
2604  ok = .false.
2605  !
2606  ! ... End loops started in 1.
2607  !
2608  END DO
2609  !
2610  RETURN
2611  !
2612  ! Formats
2613  !
2614 #ifdef W3_T4
2615 9000 FORMAT ( 'TEST GRSEPA: IPAVG/CHCK:',2i8)
2616 9010 FORMAT ( 'TEST GRSEPA: WORKING ON GRID'i4)
2617 9020 FORMAT ( ' GRID TOO SMALL TO CUT',2i8)
2618 9030 FORMAT ( 'TEST GRSEPA: GRID',i4,' IS CONTIGUOUS')
2619 9040 FORMAT ( 'TEST GRSEPA: GRID',i4,' CONTAINS PARTS')
2620 9041 FORMAT ( ' PART, SIZE, NEIGHBOUR:',i4,i8,l4)
2621 9050 FORMAT ( ' NO PART NEXT TO OTHER')
2622 9051 FORMAT ( ' NO PART SMALL ENOUGH')
2623 9060 FORMAT ( ' CUTTING PART',i4)
2624 #endif
2625  !
2626  !/ End of GRSEPA ----------------------------------------------------- /
2627  !/

◆ grsngl()

subroutine w3gspl::grsngl ( logical, intent(inout)  OK)

Remove seapoints with only one adjacent point in same grid.

Remove points from a grid that are in the middle of the sea, but that have only one adjacent point in the same grid. Directly select a new grid for this point rather than deactivate and use GRFILL.

Parameters
[in,out]OKFlag for grid status
Author
H. L. Tolman
Date
09-Sep-2012

Definition at line 2123 of file ww3_gspl.F90.

2123  !/
2124  !/ +-----------------------------------+
2125  !/ | WAVEWATCH III NOAA/NCEP |
2126  !/ | H. L. Tolman |
2127  !/ | FORTRAN 90 |
2128  !/ | Last update : 09-Sep-2012 |
2129  !/ +-----------------------------------+
2130  !/
2131  !/ 09-Sep-2012 : Origination. ( version 4.10 )
2132  !/
2133  ! 1. Purpose :
2134  !
2135  ! Remove points from a grid that are in the middle of the sea, but
2136  ! that have omly one adjacent point in the same grid. Directly
2137  ! select a new grid for this point rather than deactivate and use
2138  ! GRFILL.
2139  !
2140  ! 3. Parameters :
2141  !
2142  ! Parameter list
2143  ! ----------------------------------------------------------------
2144  ! OK Log. I/O Flag for grid status, .F. if values of
2145  ! -1 are left in MSPLIT.
2146  ! ----------------------------------------------------------------
2147  !
2148  ! 10. Source code :
2149  !
2150  !/ ------------------------------------------------------------------- /
2151  !/
2152  IMPLICIT NONE
2153  !/
2154  !/ ------------------------------------------------------------------- /
2155  !/ Parameter list
2156  !/
2157  LOGICAL, INTENT(INOUT) :: OK
2158  !/
2159  !/ ------------------------------------------------------------------- /
2160  !/ Local parameters
2161  !/
2162  INTEGER :: NX0, NXN, IXL, IXH, COUNT(-1:NG), &
2163  INEW1, INEW2, INEW
2164 #ifdef W3_S
2165  INTEGER, SAVE :: IENT = 0
2166 #endif
2167  !/
2168  !/ ------------------------------------------------------------------- /
2169  !/
2170 #ifdef W3_S
2171  CALL strace (ient, 'GRSNGL')
2172 #endif
2173  !
2174  ! 1. Set up looping ------------------------------------------------- *
2175  !
2176  IF ( global ) THEN
2177  nx0 = 1
2178  nxn = nx
2179  ELSE
2180  nx0 = 2
2181  nxn = nx-1
2182  END IF
2183  !
2184  ! 2. Loops over 2D grid --------------------------------------------- *
2185  !
2186  DO ix=nx0, nxn
2187  !
2188  ixl = ix - 1
2189  ixh = ix + 1
2190  IF ( ix .EQ. 1 ) ixl = nx
2191  IF ( ix .EQ. nx ) ixh = 1
2192  !
2193  DO iy=2, ny-1
2194  !
2195  ! 3. Central sea points only ---------------------------------------- *
2196  !
2197  IF ( sea(iy,ix) .AND. sea(iy-1,ix ) .AND. sea(iy+1,ix ) &
2198  .AND. sea(iy ,ixl) .AND. sea(iy ,ixh) ) THEN
2199  !
2200  ! 4. Check for 'lost points' ---------------------------------------- *
2201  !
2202  count = 0
2203  ig = msplit(iy,ix)
2204  !
2205  count(ig) = 1
2206  count(msplit(iy-1,ix )) = count(msplit(iy-1,ix )) + 1
2207  count(msplit(iy+1,ix )) = count(msplit(iy+1,ix )) + 1
2208  count(msplit(iy ,ixl)) = count(msplit(iy ,ixl)) + 1
2209  count(msplit(iy ,ixh)) = count(msplit(iy ,ixh)) + 1
2210  !
2211  IF ( count(ig) .LE. 2 ) THEN
2212  !
2213 #ifdef W3_T3
2214  WRITE (ndst,9040) ix, iy, ig
2215 #endif
2216  !
2217  inew1 = -1
2218  inew2 = -1
2219  !
2220  DO j=1, ng
2221  IF ( count(j) .GE. 2 ) THEN
2222 #ifdef W3_T3
2223  WRITE (ndst,9041) j
2224 #endif
2225  IF ( inew1 .EQ. -1 ) THEN
2226  inew1 = j
2227  ELSE
2228  inew2 = j
2229  EXIT
2230  END IF
2231  END IF
2232  END DO
2233  !
2234  IF ( inew1 .EQ. -1 ) THEN
2235  inew = -1
2236  ok = .false.
2237 #ifdef W3_T3
2238  WRITE (ndst,9043)
2239 #endif
2240  ELSE IF ( inew2 .EQ. -1 ) THEN
2241  inew = inew1
2242 #ifdef W3_T3
2243  WRITE (ndst,9042) inew
2244 #endif
2245  ELSE
2246  IF ( gstats(inew1)%NPTS .GT. &
2247  gstats(inew2)%NPTS ) THEN
2248  inew = inew2
2249  ELSE
2250  inew = inew1
2251  END IF
2252 #ifdef W3_T3
2253  WRITE (ndst,9042) inew
2254 #endif
2255  END IF
2256  !
2257  msplit(iy,ix) = inew
2258  !
2259  END IF
2260  !
2261  END IF
2262  !
2263  ! ... End loops started in 2.
2264  !
2265  END DO
2266  !
2267  END DO
2268  !
2269  RETURN
2270  !
2271  ! Formats
2272  !
2273 #ifdef W3_T3
2274 9040 FORMAT ( 'TEST GRSNGL: POINT FOUND, IX, IY, IG:',2i5,i4)
2275 9041 FORMAT ( ' CANDIDATE GRID :',10x,i4)
2276 9042 FORMAT ( ' GRID USED :',10x,i4)
2277 9043 FORMAT ( ' GRID LEFT UNDIFINED')
2278 #endif
2279  !
2280  !/ End of GRSNGL ----------------------------------------------------- /
2281  !/

◆ grsqrg()

subroutine w3gspl::grsqrg

Attempt to square-up grid.

Attempt to square-up grid by taking off grid point in outermost grid point in X and Y only, after which GRFILL is to be run to re-assign grid points.

Author
H. L. Tolman
Date
07-Sep-2012

Definition at line 2017 of file ww3_gspl.F90.

2017  !/
2018  !/ +-----------------------------------+
2019  !/ | WAVEWATCH III NOAA/NCEP |
2020  !/ | H. L. Tolman |
2021  !/ | FORTRAN 90 |
2022  !/ | Last update : 07-Sep-2012 |
2023  !/ +-----------------------------------+
2024  !/
2025  !/ 07-Sep-2012 : Origination. ( version 4.10 )
2026  !/
2027  ! 1. Purpose :
2028  !
2029  ! Attemp to square-up grid by taking off grid point in outermost
2030  ! grid point in X and Y only, after which GRFILL is to be run to
2031  ! re-assign grid points,
2032  !
2033  ! 10. Source code :
2034  !
2035  !/ ------------------------------------------------------------------- /
2036  !/
2037  IMPLICIT NONE
2038  !/
2039  !/ ------------------------------------------------------------------- /
2040  !/ Parameter list
2041  !/
2042  !/ ------------------------------------------------------------------- /
2043  !/ Local parameters
2044  !/
2045  INTEGER :: MX, MY
2046 #ifdef W3_S
2047  INTEGER, SAVE :: IENT = 0
2048 #endif
2049  !/
2050  !/ ------------------------------------------------------------------- /
2051  !/
2052 #ifdef W3_S
2053  CALL strace (ient, 'GRSQRG')
2054 #endif
2055  !
2056  ! 1. Loop over grids ------------------------------------------------ *
2057  !
2058  DO ig=1, ng
2059  !
2060  my = 1 + gstats(ig)%NYH - gstats(ig)%NYL
2061  mx = 1 + gstats(ig)%NXH - gstats(ig)%NXL
2062  IF ( gstats(ig)%STRADLE ) mx = mx + nx
2063  !
2064  ! 2. Top ------------------------------------------------------------ *
2065  !
2066  IF ( my .GE. 5 ) THEN
2067  !
2068  DO ix=1, nx
2069  IF (msplit(gstats(ig)%NYH,ix) .EQ. ig ) &
2070  msplit(gstats(ig)%NYH,ix) = -1
2071  END DO
2072  !
2073  ! 3. Bottom --------------------------------------------------------- *
2074  !
2075  DO ix=1, nx
2076  IF (msplit(gstats(ig)%NYL,ix) .EQ. ig ) &
2077  msplit(gstats(ig)%NYL,ix) = -1
2078  END DO
2079  !
2080  END IF
2081  !
2082  ! 4. Left ----------------------------------------------------------- *
2083  !
2084  IF ( mx .GE. 5 ) THEN
2085  !
2086  DO iy=gstats(ig)%NYL, gstats(ig)%NYH
2087  IF (msplit(iy,gstats(ig)%NXL) .EQ. ig ) &
2088  msplit(iy,gstats(ig)%NXL) = -1
2089  END DO
2090  !
2091  ! 5. Right ---------------------------------------------------------- *
2092  !
2093  DO iy=gstats(ig)%NYH, gstats(ig)%NYH
2094  IF (msplit(iy,gstats(ig)%NXH) .EQ. ig ) &
2095  msplit(iy,gstats(ig)%NXH) = -1
2096  END DO
2097  !
2098  END IF
2099  !
2100  ! ... End loop started in 1.
2101  !
2102  END DO
2103  !
2104  RETURN
2105  !
2106  ! Formats
2107  !
2108  !/ End of GRSQRG ----------------------------------------------------- /
2109  !/

◆ grtrim()

subroutine w3gspl::grtrim

Trim edges of all grids where they are next to another grid or next to unassigned grid points.

This trimming is done in preparation for reassigning edges of grids to smaller adjacent grids.

Author
H. L. Tolman
Date
01-Feb-2013

Definition at line 1480 of file ww3_gspl.F90.

1480  !/
1481  !/ +-----------------------------------+
1482  !/ | WAVEWATCH III NOAA/NCEP |
1483  !/ | H. L. Tolman |
1484  !/ | FORTRAN 90 |
1485  !/ | Last update : 01-Feb-2013 |
1486  !/ +-----------------------------------+
1487  !/
1488  !/ 07-Sep-2012 : Origination. ( version 4.10 )
1489  !/ 18-Sep-2012 : Include edge points of grid. ( version 4.10 )
1490  !/ 01-Feb-2013 : Add dynamic trim range. ( version 4.10 )
1491  !/
1492  ! 1. Purpose :
1493  !
1494  ! Trim edges of all grids where they are next to another grid or next
1495  ! to unassigned grid points. This is done in preparation for
1496  ! reassigning edges of grids to smaller adjacent grids.
1497  !
1498  ! 10. Source code :
1499  !
1500  !/ ------------------------------------------------------------------- /
1501  !/
1502  IMPLICIT NONE
1503  !/
1504  !/ ------------------------------------------------------------------- /
1505  !/ Parameter list
1506  !/
1507  !/ ------------------------------------------------------------------- /
1508  !/ Local parameters
1509  !/
1510  INTEGER :: ITARG, ITL, IPTS, MX, MY, ICIRC, NWDTH
1511 #ifdef W3_S
1512  INTEGER, SAVE :: IENT = 0
1513 #endif
1514  LOGICAL :: MASK(NY,NX)
1515  !/
1516  !/ ------------------------------------------------------------------- /
1517  !/
1518 #ifdef W3_S
1519  CALL strace (ient, 'GRTRIM')
1520 #endif
1521  !
1522  itarg = nsea / ng
1523  !
1524  ! 1. Loop over grids ------------------------------------------------ *
1525  !
1526  DO ig=1, ng
1527  !
1528  ipts = gstats(ig)%NPTS
1529  my = 1 + gstats(ig)%NYH - gstats(ig)%NYL
1530  mx = 1 + gstats(ig)%NXH - gstats(ig)%NXL
1531  IF ( gstats(ig)%STRADLE ) mx = mx + nx
1532  icirc = 2 * ( mx + my )
1533  !
1534  nwdth = 1
1535  !
1536  itl = min( itarg , max( itarg-2*icirc , 3*icirc ) )
1537  IF ( ipts .LT. itl ) nwdth = 0
1538  !
1539  IF ( ipts.GT.itarg ) THEN
1540  nwdth = 1 + &
1541  max(0,+nint((real((ipts-itarg))/real(icirc)-1.)/3.))
1542  ENDIF
1543  !
1544  DO j=1, nwdth
1545  !
1546  mask = .false.
1547  !
1548  ! 2. Mark points to be removed -------------------------------------- *
1549  !
1550  DO ix=2, nx-1
1551  IF ( msplit( 1,ix) .EQ. ig ) mask( 1,ix) = &
1552  (sea( 2,ix ).AND.(msplit( 2,ix ).NE.ig)) &
1553  .OR. (sea( 1,ix+1).AND.(msplit( 1,ix+1).NE.ig)) &
1554  .OR. (sea( 1,ix-1).AND.(msplit( 1,ix-1).NE.ig))
1555  DO iy=2, ny-1
1556  IF ( msplit(iy,ix) .EQ. ig ) mask(iy,ix) = &
1557  (sea(iy+1,ix ).AND.(msplit(iy+1,ix ).NE.ig)) &
1558  .OR. (sea(iy-1,ix ).AND.(msplit(iy-1,ix ).NE.ig)) &
1559  .OR. (sea(iy ,ix+1).AND.(msplit(iy ,ix+1).NE.ig)) &
1560  .OR. (sea(iy ,ix-1).AND.(msplit(iy ,ix-1).NE.ig))
1561  END DO
1562  IF ( msplit(ny,ix) .EQ. ig ) mask(ny,ix) = &
1563  (sea(ny-1,ix ).AND.(msplit(ny-1,ix ).NE.ig)) &
1564  .OR. (sea(ny ,ix+1).AND.(msplit(ny ,ix+1).NE.ig)) &
1565  .OR. (sea(ny ,ix-1).AND.(msplit(ny ,ix-1).NE.ig))
1566  END DO
1567  !
1568  IF ( global ) THEN
1569  IF ( msplit( 1, 1) .EQ. ig ) mask( 1, 1) = &
1570  (sea( 2, 1).AND.(msplit( 2, 1).NE.ig)) &
1571  .OR. (sea( 1, 2).AND.(msplit( 1, 2).NE.ig)) &
1572  .OR. (sea( 1,nx).AND.(msplit( 1,nx).NE.ig))
1573  IF ( msplit( 1,nx) .EQ. ig ) mask( 1,nx) = &
1574  (sea( 2,nx ).AND.(msplit( 2,nx ).NE.ig)) &
1575  .OR. (sea( 1, 1 ).AND.(msplit( 1, 1 ).NE.ig)) &
1576  .OR. (sea( 1,nx-1).AND.(msplit( 1,nx-1).NE.ig))
1577  DO iy=2, ny-1
1578  IF ( msplit(iy, 1) .EQ. ig ) mask(iy, 1) = &
1579  (sea(iy+1, 1).AND.(msplit(iy+1, 1).NE.ig)) &
1580  .OR. (sea(iy-1, 1).AND.(msplit(iy-1, 1).NE.ig)) &
1581  .OR. (sea(iy , 2).AND.(msplit(iy , 2).NE.ig)) &
1582  .OR. (sea(iy ,nx).AND.(msplit(iy ,nx).NE.ig))
1583  IF ( msplit(iy,nx) .EQ. ig ) mask(iy,nx) = &
1584  (sea(iy+1,nx).AND.(msplit(iy+1,nx).NE.ig)) &
1585  .OR. (sea(iy-1,nx).AND.(msplit(iy-1,nx).NE.ig)) &
1586  .OR. (sea(iy , 1).AND.(msplit(iy , 1).NE.ig)) &
1587  .OR. (sea(iy,nx-1).AND.(msplit(iy,nx-1).NE.ig))
1588  END DO
1589  IF ( msplit(ny, 1) .EQ. ig ) mask(ny, 1) = &
1590  (sea(ny-1, 1).AND.(msplit(ny-1, 1).NE.ig)) &
1591  .OR. (sea(ny , 2).AND.(msplit(ny , 2).NE.ig)) &
1592  .OR. (sea(ny ,nx).AND.(msplit(ny ,nx).NE.ig))
1593  IF ( msplit(ny,nx) .EQ. ig ) mask(ny,nx) = &
1594  (sea(ny-1,nx).AND.(msplit(ny-1,nx).NE.ig)) &
1595  .OR. (sea(ny , 1).AND.(msplit(ny , 1).NE.ig)) &
1596  .OR. (sea(ny,nx-1).AND.(msplit(ny,nx-1).NE.ig))
1597  ELSE
1598  IF ( msplit( 1, 1) .EQ. ig ) mask( 1, 1) = &
1599  (sea( 2, 1).AND.(msplit( 2, 1).NE.ig)) &
1600  .OR. (sea( 1, 2).AND.(msplit( 1, 2).NE.ig))
1601  IF ( msplit( 1,nx) .EQ. ig ) mask( 1,nx) = &
1602  (sea( 2,nx ).AND.(msplit( 2,nx ).NE.ig)) &
1603  .OR. (sea( 1,nx-1).AND.(msplit( 1,nx-1).NE.ig))
1604  DO iy=2, ny-1
1605  IF ( msplit(iy, 1) .EQ. ig ) mask(iy, 1) = &
1606  (sea(iy+1, 1).AND.(msplit(iy+1, 1).NE.ig)) &
1607  .OR. (sea(iy-1, 1).AND.(msplit(iy-1, 1).NE.ig)) &
1608  .OR. (sea(iy , 2).AND.(msplit(iy , 2).NE.ig))
1609  IF ( msplit(iy,nx) .EQ. ig ) mask(iy,nx) = &
1610  (sea(iy+1,nx).AND.(msplit(iy+1,nx).NE.ig)) &
1611  .OR. (sea(iy-1,nx).AND.(msplit(iy-1,nx).NE.ig)) &
1612  .OR. (sea(iy,nx-1).AND.(msplit(iy,nx-1).NE.ig))
1613  END DO
1614  IF ( msplit(ny, 1) .EQ. ig ) mask(ny, 1) = &
1615  (sea(ny-1, 1).AND.(msplit(ny-1, 1).NE.ig)) &
1616  .OR. (sea(ny , 2).AND.(msplit(ny , 2).NE.ig))
1617  IF ( msplit(ny,nx) .EQ. ig ) mask(ny,nx) = &
1618  (sea(ny-1,nx).AND.(msplit(ny-1,nx).NE.ig)) &
1619  .OR. (sea(ny,nx-1).AND.(msplit(ny,nx-1).NE.ig))
1620  END IF
1621  !
1622  ! 3. Remove marked points ------------------------------------------- *
1623  !
1624  DO ix=1, nx
1625  DO iy=1, ny
1626  IF ( mask(iy,ix) ) THEN
1627  msplit(iy,ix) = -1
1628  END IF
1629  END DO
1630  END DO
1631  !
1632  ! ... End loops started in 1.
1633  !
1634  END DO
1635  END DO
1636  !
1637  RETURN
1638  !
1639  ! Formats
1640  !
1641  !/ End of GRTRIM ----------------------------------------------------- /
1642  !/

◆ w3gspl()

program w3gspl

Grid splitting program.

Take an existing grid and create from this the grid data for a set of overlapping grids to be used in the ww3_multi code for hybid paralellization.

Author
H. L. Tolman
Date
18-Nov-2013

Definition at line 26 of file ww3_gspl.F90.

References w3servmd::extcde(), w3odatmd::fnmpre, w3servmd::itrace(), w3odatmd::ndse, w3odatmd::ndso, w3odatmd::ndst, w3servmd::nextln(), w3arrymd::outa2i(), w3arrymd::outa2r(), w3servmd::strace(), w3iogrmd::w3iogr(), w3adatmd::w3naux(), w3odatmd::w3nout(), w3adatmd::w3seta(), and w3odatmd::w3seto().

grlost
subroutine grlost
Reassign unassigned grid points to grids.
Definition: ww3_gspl.F90:1913