WAVEWATCH III  beta 0.0.1
ww3_gspl.F90
Go to the documentation of this file.
1 
7 
8 #include "w3macros.h"
9 !/ ------------------------------------------------------------------- /
25 
26 PROGRAM w3gspl
27  !/
28  !/ +-----------------------------------+
29  !/ | WAVEWATCH III NOAA/NCEP |
30  !/ | H. L. Tolman |
31  !/ | FORTRAN 90 |
32  !/ | Last update : 18-Nov-2013 |
33  !/ +-----------------------------------+
34  !/
35  !/ 24-Sep-2012 : Origination. ( version 4.10 )
36  !/ 16-Jan-2013 : Add output of mask file (no halo). ( version 4.10 )
37  !/ 19-Jan-2013 : Tweaking the template file. ( version 4.10 )
38  !/ 24-Jan-2013 : Set up for minimum of 2 grids. ( version 4.10 )
39  !/ Add XOFF to grid origin in X.
40  !/ Fix IDCLSE for partial grids.
41  !/ Add FRFLAG option to disable side-by-side
42  !/ running of grids in ww3_multi.
43  !/ 29-Jan-2013 : Add error code on stop. ( version 4.10 )
44  !/ 31-Jan-2013 : Add routine GRLOST. ( version 4.10 )
45  !/ 01-Feb-2013 : Speed up GRSEPA. ( version 4.10 )
46  !/ Add dynamic trim range in GRTRIM.
47  !/ Speed up GRFILL.
48  !/ Add small grid merge (GRFSML) early in loop.
49  !/ 04-Feb-2013 : Testing on zero grid size added. ( version 4.10 )
50  !/ Corner point in halo for GR1GRD.
51  !/ 04-Mar-2013 : Adding GrADS output. ( version 4.10 )
52  !/ 05-Aug-2013 : Add UQ/UNO for distances. ( version 4.12 )
53  !/ 18-Nov-2013 : Add user-defined halo extension. ( version 4.14 )
54  !/
55  !/ Copyright 2012-2013 National Weather Service (NWS),
56  !/ National Oceanic and Atmospheric Administration. All rights
57  !/ reserved. WAVEWATCH III is a trademark of the NWS.
58  !/ No unauthorized use without permission.
59  !/
60  ! 1. Purpose :
61  !
62  ! Take an existing grid and create from this the grid data for a set
63  ! of overlapping grids to be used in the ww3_multi code for hybid
64  ! paralellization.
65  !
66  ! 2. Method :
67  !
68  ! See Section 8.
69  !
70  ! 3. Parameters :
71  !
72  ! Local parameters.
73  ! ----------------------------------------------------------------
74  ! NDSI Int. Input unit number ("ww3_prep.inp").
75  ! NDSO Int. Output unit number.
76  ! NDSE Int. Error unit number.
77  ! NDST Int. Test output unit number.
78  ! NDSM Int. Unit number for mod_def file.
79  ! NG Int. Number of grids to be generated.
80  ! NITMAX Int. Maximum number of iterations on grid ref.
81  ! STARG Real std target in percent.
82  ! GLOBAL Log. Closure flag.
83  ! SEA L.A. Sea point map.
84  ! ----------------------------------------------------------------
85  !
86  ! 4. Subroutines used :
87  !
88  ! Name Type Module Description
89  ! ----------------------------------------------------------------
90  ! W3NMOD Subr. W3GDATMD Set number of model.
91  ! W3SETG Subr. Id. Point to selected model.
92  ! W3NDAT Subr. W3WDATMD Set number of model for wave data.
93  ! W3SETW Subr. Id. Point to selected model for wave data.
94  ! W3NOUT Subr. W3ODATMD Set number of model for output.
95  ! W3SETO Subr. Id. Point to selected model for output.
96  ! ITRACE Subr. W3SERVMD Subroutine tracing initialization.
97  ! STRACE Subr. Id. Subroutine tracing.
98  ! NEXTLN Subr. Id. Get next line from input filw
99  ! EXTCDE Subr. Id. Abort program as graceful as possible.
100  ! W3IOGR Subr. W3IOGRMD Reading/writing model definition file.
101  !
102  ! GRINFO Subr. Internal Compile info on all grids.
103  ! GRTRIM Subr. Internal Trim edges of grids.
104  ! GRFILL Subr. Internal Fill unassigned space in grid.
105  ! GRLOST Subr. Internal Assign "lost points".
106  ! GRSQRG Subr. Internal Attempt to square-up grid.
107  ! GRSNGL Subr. Internal Remove grid points that stick out.
108  ! GRSEPA Subr. Internal Remove separated grid pieces.
109  ! GRFSML Subr. Internal Deal with fixed minimum size.
110  ! GRFLRG Subr. Internal Deal with fixed maximum size.
111  ! GR1GRD Subr. Internal Extract single grid from map.
112  ! ----------------------------------------------------------------
113  !
114  ! 5. Called by :
115  !
116  ! None, stand-alone program.
117  !
118  ! 6. Error messages :
119  !
120  ! 7. Remarks :
121  !
122  ! 8. Structure :
123  !
124  ! ----------------------------------------------------
125  ! 1.a Number of models.
126  ! ( W3NMOD , W3NOUT , W3SETG , W3SETO )
127  ! b I-O setup.
128  ! c Print heading(s).
129  ! 2. Read model definition file. ( W3IOGR )
130  ! 3. Read options from file.
131  ! 4. Generate first-guess map of sub-grids
132  ! a Set up array
133  ! b First cut with regular grid set up
134  ! 1 Set up 'checkerboard'
135  ! 2 Fill checkerboard
136  ! 3 Remove smallest grids as necessary
137  ! 4 Store first guess in MSPLIT
138  ! 5. Refine map of sub-grids (no halo).
139  ! a Set up loop. ( GRINFO )
140  ! b Remove small grids. ( GRFSML )
141  ! c Trim edges of grids ( GRTRIM )
142  ! ( GRFILL, GRLOST )
143  ! d Attempt to square-up grid ( GRINFO )
144  ! ( GRSQRG )
145  ! ( GRFILL, GRLOST )
146  ! e Remove mid-sea points sticking out of grid
147  ! ( GRSNGL )
148  ! f Remove detached grid parts. ( GRSEPA )
149  ! g Recompute stats ( GRINFO )
150  ! h Optional GrADS output.
151  ! i Test convergence
152  ! Check if stuck on min or max. ( GRFSML )
153  ! ( GRFLRG )
154  ! j Test output
155  ! 6. Output info for all sub grids.
156  ! a Set up loop. ( GRINFO )
157  ! b Extract grid including halo. ( GR1GRD )
158  ! 7. End of program.
159  ! ----------------------------------------------------
160  !
161  ! 9. Switches :
162  !
163  ! !/PRn Select propgation scheme.
164  !
165  ! !/O16 Generate GrADS output of grid partitioning.
166  !
167  ! !/S Enable subroutine tracing.
168  ! !/T Enable test output (main).
169  ! !/T1 Enable test output (GRINFO).
170  ! !/T2 Enable test output (GRFILL).
171  ! !/T3 Enable test output (GRSNGL).
172  ! !/T4 Enable test output (GRSEPA).
173  ! !/T5 Enable test output (GRFSML).
174  ! !/T6 Enable test output (GRFLRG).
175  ! !/T7 Enable test output (GR1GRD).
176  !
177  ! 10. Source code :
178  !
179  !/ ------------------------------------------------------------------- /
180  USE constants
181  !/
182  ! USE W3GDATMD, ONLY: W3NMOD, W3SETG
183  USE w3adatmd, ONLY: w3naux, w3seta
184  USE w3odatmd, ONLY: w3nout, w3seto
185  USE w3servmd, ONLY : itrace, nextln, extcde
186  USE w3arrymd, ONLY : outa2i, outa2r
187 #ifdef W3_S
188  USE w3servmd, ONLY : strace
189 #endif
190  USE w3iogrmd, ONLY: w3iogr
191  !/
192  USE w3gdatmd
193  USE w3odatmd, ONLY: ndse, ndst, ndso, fnmpre
194  !
195  IMPLICIT NONE
196  !/
197  !/ ------------------------------------------------------------------- /
198  !/ Local parameters
199  !/
200  INTEGER :: ndsi, ndsm, ndstrc, ntrace, j, ierr, &
201  ng, ix, iy, ngb, ngx, ngy, ig, igg, &
202  igx, igy, igy0, igyn, igx0, igxn, &
203  mingrd, minnr, minnxt, minnnr, &
204  nitmax, iit, ingmin, ingmax, &
205  ingmnc, ingmxc, inglag, jj, &
206  nstdlg, mstdlg = 5, nseat, j1, j2, &
207  j3, j4, j5, idfm1, idfm2, idfm3, &
208  idla1, idla2, idla3, vsc3, nhext
209 #ifdef W3_S
210  INTEGER, SAVE :: ient = 0
211 #endif
212 #ifdef W3_O16
213  INTEGER :: ndsg = 35, ntgrds = 0
214 #endif
215  INTEGER, ALLOCATABLE :: msplit(:,:), mtemp(:,:), ingrd(:)
216  REAL :: ratio1, xmean, starg, stdmin, &
217  zbdum, zbmin, vsc1, vsc2, fracl, frach
218  LOGICAL :: global, ok, done, frflag
219  LOGICAL, ALLOCATABLE :: isnext(:), sea(:,:)
220  CHARACTER(LEN=1) :: comstr
221  CHARACTER(LEN=3) :: g0id
222  CHARACTER(LEN=4) :: idgrid, idclse, ptclse
223  CHARACTER(LEN=6) :: nrfmt
224  CHARACTER(LEN=11) :: fext, aext
225  CHARACTER(LEN=16) :: rform1, rform2, rform3
226  CHARACTER(LEN=20) :: fname, iname
227  !
229  LOGICAL :: stradle, instat
230  INTEGER :: npts, nyl, nyh, nxl, nxh
231  END TYPE stats_grid
232  !
234  INTEGER :: nmin, nmax
235  REAL :: rstd
236  END TYPE stats_mean
237  !
239  INTEGER :: nx, ny, nsea
240  INTEGER, POINTER :: mask(:,:)
241  REAL :: x0, y0, sx, sy
242  REAL, POINTER :: zbin(:,:), obsx(:,:), obsy(:,:)
243  LOGICAL :: global
244  END TYPE part_grid
245  !
246  TYPE(stats_grid), POINTER :: gstats(:), gstold(:)
247  TYPE(stats_mean) :: mstats , mstold
248  TYPE(part_grid), POINTER :: pgrid(:)
249  !/
250  !/ ------------------------------------------------------------------- /
251  !/
252  ! 1.a Set number of models
253  !
254  CALL w3nmod ( 1, 6, 6 )
255  CALL w3setg ( 1, 6, 6 )
256  CALL w3naux ( 6, 6 )
257  CALL w3seta ( 1, 6, 6 )
258  CALL w3nout ( 6, 6 )
259  CALL w3seto ( 1, 6, 6 )
260  !
261  ! 1.b IO set-up.
262  !
263  ndsi = 10
264  ndso = 6
265  ndse = 6
266  ndst = 6
267  ndsm = 11
268  !
269  ndstrc = 6
270  ntrace = 100
271  CALL itrace ( ndstrc, ntrace )
272  !
273 #ifdef W3_O16
274  OPEN ( ndsg, file='./ww3.ww3_gspl', form='UNFORMATTED', convert=file_endian)
275 #endif
276  !
277  ! 1.c Print header
278  !
279  WRITE (ndso,900)
280 #ifdef W3_S
281  CALL strace (ient, 'W3GSPL')
282 #endif
283  !
284  j = len_trim(fnmpre)
285  OPEN (ndsi,file=fnmpre(:j)//'ww3_gspl.inp',status='OLD', &
286  err=800,iostat=ierr)
287  rewind(ndsi)
288  READ (ndsi,'(A)',END=801,ERR=802,IOSTAT=IERR) comstr
289  IF (comstr.EQ.' ') comstr = '$'
290  WRITE (ndso,901) comstr
291  !
292  !--- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
293  ! 2. Read model definition file.
294  !
295  CALL nextln ( comstr , ndsi , ndse )
296  READ (ndsi,*,END=801,ERR=802,IOSTAT=IERR) fext
297  !
298  CALL w3iogr ( 'READ', ndsm, 1, fext )
299  CLOSE (ndsm)
300  !
301  WRITE (ndso,902) fext, gname
302  !
303  SELECT CASE (gtype)
304  CASE (rlgtype)
305  WRITE ( ndso,903) 'rectilinear'
306  idgrid = 'RECT'
307  CASE (clgtype)
308  WRITE ( ndso,903) 'curvictilinear'
309  idgrid = 'CURV'
310  CASE (ungtype)
311  WRITE ( ndso,903) 'unstructured'
312  idgrid = 'UNST'
313  GOTO 820
314  CASE DEFAULT
315  WRITE ( ndso,903) 'not recognized'
316  GOTO 821
317  END SELECT
318  !
319  SELECT CASE (iclose)
320  CASE (iclose_none)
321  WRITE ( ndso,904) 'none'
322  idclse = 'NONE'
323  global = .false.
324  CASE (iclose_smpl)
325  WRITE ( ndso,904) 'global (simple)'
326  idclse = 'SMPL'
327  global = .true.
328  CASE (iclose_trpl)
329  WRITE ( ndso,904) 'global (tripolar)'
330  idclse = 'TRPL'
331  global = .true.
332  GOTO 822
333  CASE DEFAULT
334  WRITE ( ndso,904) 'not recognized'
335  GOTO 823
336  END SELECT
337  !
338  WRITE (ndso,905) nx, ny, nsea
339  IF ( nsea .EQ. 0 ) GOTO 824
340  !
341  !--- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
342  ! 3. Read options from input file.
343  !
344  CALL nextln ( comstr , ndsi , ndse )
345  READ (ndsi,*,END=801,ERR=802,IOSTAT=IERR) NG, NITMAX, STARG, nhext
346  ng = max( 2, ng )
347  nitmax = max( 1, nitmax )
348  starg = max( 0. , starg )
349  nhext = max( 0, nhext )
350  WRITE (ndso,930) ng, nitmax, starg, nhext
351  !
352  CALL nextln ( comstr , ndsi , ndse )
353  READ (ndsi,*,END=801,ERR=802,IOSTAT=IERR) IDLA1, IDFM1, &
354  vsc1, rform1
355  IF (idla1.LT.1 .OR. idla1.GT.4) idla1 = 1
356  IF (idfm1.LT.1 .OR. idfm1.GT.3) idfm1 = 1
357  IF ( abs(vsc1) .LT. 1.e-15 ) vsc1 = 1.
358 
359  WRITE (ndso,931) idla1, idfm1, vsc1, rform1
360  !
361  CALL nextln ( comstr , ndsi , ndse )
362  READ (ndsi,*,END=801,ERR=802,IOSTAT=IERR) IDLA2, IDFM2, &
363  vsc2, rform2
364  IF (idla2.LT.1 .OR. idla2.GT.4) idla2 = 1
365  IF (idfm2.LT.1 .OR. idfm2.GT.3) idfm2 = 1
366  IF ( abs(vsc2) .LT. 1.e-15 ) vsc2 = 1.
367  IF ( trflag .EQ. 0 ) THEN
368  WRITE (ndso,932)
369  ELSE
370  WRITE (ndso,933) idla2, idfm2, vsc2, rform2
371  END IF
372  !
373  CALL nextln ( comstr , ndsi , ndse )
374  READ (ndsi,*,END=801,ERR=802,IOSTAT=IERR) IDLA3, IDFM3, &
375  vsc3, rform3
376  IF (idla3.LT.1 .OR. idla3.GT.4) idla3 = 1
377  IF (idfm3.LT.1 .OR. idfm3.GT.3) idfm3 = 1
378  IF ( vsc3 .EQ. 0 ) vsc3 = 1
379  WRITE (ndso,934) idla3, idfm3, vsc3, rform3
380  !
381  CALL nextln ( comstr , ndsi , ndse )
382  READ (ndsi,*,END=801,ERR=802,IOSTAT=IERR) FRACL, FRACH, frflag
383  fracl = max( 0. , fracl )
384  frach = min( 1. , frach )
385  WRITE (ndso,935) fracl, frach
386  IF ( fracl .GT. frach ) GOTO 830
387  IF ( .NOT. frflag ) WRITE (ndso,936)
388  !
389  !--- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
390  ! 4. Generate map of sub-grids (no halo)
391  ! 4.a Set up array
392  !
393  ALLOCATE ( msplit(ny,nx) , mtemp(ny,nx), sea(ny,nx) )
394  !
395  DO iy=1, ny
396  DO ix=1, nx
397  IF ( mapsta(iy,ix) .EQ. 0 ) THEN
398  msplit(iy,ix) = 0
399  sea(iy,ix) = .false.
400  ELSE
401  msplit(iy,ix) = -1
402  sea(iy,ix) = .true.
403  END IF
404  END DO
405  END DO
406 
407  !
408  ! 4.b First cut with regular grid set up
409  ! 4.b.1 Set up 'checkerboard'
410  !
411  ratio1 = real(nx) / real(ny)
412  !
413  ngx = 1
414  ngy = 1
415  !
416  DO
417  IF ( ngx*ngy .GE. ng ) EXIT
418  IF ( real(ngx)/real(ngy) .GT. ratio1 ) THEN
419  ngy = ngy + 1
420  ELSE
421  ngx = ngx + 1
422  END IF
423  END DO
424  !
425  IF ( ngx .GT. ngy ) THEN
426  IF ( (ngy-1)*ngx .GE. ng ) ngy = ngy - 1
427  IF ( (ngx-1)*ngy .GE. ng ) ngx = ngx - 1
428  ELSE
429  IF ( (ngy-1)*ngx .GE. ng ) ngy = ngy - 1
430  IF ( (ngx-1)*ngy .GE. ng ) ngx = ngx - 1
431  END IF
432  !
433 #ifdef W3_T
434  WRITE (ndst,9040) ngx, ngy
435 #endif
436  !
437  ! 4.b.2 Fill checkerboard
438  !
439  j = 0
440  DO
441  !
442  mtemp = msplit
443  ig = 1
444  igyn = 0
445  j = j + 1
446  ALLOCATE ( ingrd(ngx*ngy) )
447  ingrd = 0
448  !
449 #ifdef W3_T
450  WRITE (ndst,9041) j
451 #endif
452  !
453  DO igy=1, ngy
454  !
455  igy0 = igyn + 1
456  IF ( igy .EQ. ngy ) THEN
457  igyn = ny
458  ELSE
459  igyn = nint( real(ny*igy) / real(ngy) )
460  END IF
461  igxn = 0
462  !
463  DO igx=1, ngx
464  !
465  igx0 = igxn + 1
466  IF ( igx .EQ. ngx ) THEN
467  igxn = nx
468  ELSE
469  igxn = nint( real(nx*igx) / real(ngx) )
470  END IF
471  !
472  DO ix=igx0, igxn
473  DO iy=igy0, igyn
474  IF ( mtemp(iy,ix) .EQ. -1 ) THEN
475  mtemp(iy,ix) = ig
476  ingrd(ig) = ingrd(ig) + 1
477  END IF
478  END DO
479  END DO
480  !
481  IF ( ingrd(ig) .GT. 0 ) THEN
482 #ifdef W3_T
483  WRITE (ndst,9042) ig, igx0, igxn, igy0, igyn, &
484  ingrd(ig), 'OK'
485 #endif
486  ig = ig + 1
487 #ifdef W3_T
488  ELSE
489  WRITE (ndst,9042) ig, igx0, igxn, igy0, igyn, &
490  ingrd(ig), 'EMPTY (SKIPPED)'
491 #endif
492  END IF
493  !
494  END DO
495  !
496  END DO
497  !
498  ig = ig - 1
499  IF ( ig .LT. ng ) THEN
500  IF ( ngx .LT. ngy ) THEN
501  ngy = ngy + 1
502  ELSE
503  ngx = ngx + 1
504  END IF
505  DEALLOCATE ( ingrd )
506 #ifdef W3_T
507  WRITE (ndst,9040) ngx, ngy
508 #endif
509  ELSE
510  EXIT
511  END IF
512  !
513  END DO
514  !
515  mingrd = 0
516  DO j=1, ig
517  mingrd = mingrd + ingrd(j)
518  END DO
519  IF ( mingrd .NE. nsea ) GOTO 825
520  !
521 #ifdef W3_T
522  WRITE (ndst,9043) ig, ng
523 #endif
524  !
525  ! 4.b.3 Merge smallest grids as necessary
526  !
527  igg = ig
528  !
529  DO
530  !
531  IF ( igg .EQ. ng ) EXIT
532  !
533  mingrd = nsea
534  minnr = 0
535  DO j=1, ig
536  IF ( ingrd(j) .LT. mingrd ) THEN
537  mingrd = ingrd(j)
538  minnr = j
539  END IF
540  END DO
541  ingrd(minnr) = nsea + 1
542  !
543 #ifdef W3_T
544  WRITE (ndst,9044) mingrd, minnr
545 #endif
546  !
547  ALLOCATE ( isnext(0:ig) )
548  isnext = .false.
549  !
550  DO iy=1, ny-1
551  DO ix=1, nx-1
552  IF ( ( mtemp(iy ,ix ) - minnr ) * &
553  ( mtemp(iy+1,ix ) - minnr ) * &
554  ( mtemp(iy ,ix+1) - minnr ) * &
555  ( mtemp(iy+1,ix+1) - minnr ) .EQ. 0 ) THEN
556  isnext(mtemp(iy ,ix )) = .true.
557  isnext(mtemp(iy+1,ix )) = .true.
558  isnext(mtemp(iy ,ix+1)) = .true.
559  isnext(mtemp(iy+1,ix+1)) = .true.
560  END IF
561  END DO
562  END DO
563  !
564  IF ( global ) THEN
565  DO iy=1, ny-1
566  IF ( ( mtemp(iy ,nx) - minnr ) * &
567  ( mtemp(iy+1,nx) - minnr ) * &
568  ( mtemp(iy , 1) - minnr ) * &
569  ( mtemp(iy+1, 1) - minnr ) .EQ. 0 ) THEN
570  isnext(mtemp(iy ,nx)) = .true.
571  isnext(mtemp(iy+1,nx)) = .true.
572  isnext(mtemp(iy , 1)) = .true.
573  isnext(mtemp(iy+1, 1)) = .true.
574  END IF
575  END DO
576  END IF
577  !
578  minnxt = nsea
579  minnnr = 0
580  DO j=1, ig
581  IF ( isnext(j) .AND. ( ingrd(j) .LT. minnxt ) ) THEN
582  minnxt = ingrd(j)
583  minnnr = j
584  END IF
585  END DO
586  !
587 #ifdef W3_T
588  WRITE (ndst,9045) minnxt, minnnr
589 #endif
590  !
591  IF ( minnnr .GT. 0 ) THEN
592  DO iy=1, ny
593  DO ix=1, nx
594  IF ( mtemp(iy,ix) .EQ. minnr ) THEN
595  mtemp(iy,ix) = minnnr
596  ingrd(minnnr) = ingrd(minnnr) + 1
597  END IF
598  END DO
599  END DO
600  igg = igg - 1
601 #ifdef W3_T
602  WRITE (ndst,9046) minnr, minnnr
603  DO j=1, ig
604  WRITE (ndst,9047) j, ingrd(j)
605  END DO
606  ELSE
607  WRITE (ndst,9048) minnr
608 #endif
609  END IF
610  !
611  DEALLOCATE ( isnext)
612 #ifdef W3_T
613  WRITE (ndst,9043) igg, ng
614 #endif
615  !
616  END DO
617  !
618 #ifdef W3_T
619  WRITE (ndst,9049) ng
620 #endif
621  !
622  DO j=1, ig
623  IF ( ingrd(j) .GT. nsea ) ingrd(j) = 0
624 #ifdef W3_T
625  WRITE (ndso,9047) j, ingrd(j)
626 #endif
627  END DO
628  !
629  ! 4.b.4 Store first guess in MSPLT
630  !
631  igg = 0
632  DO j=1, ig
633  IF ( ingrd(j) .NE. 0 ) THEN
634  igg = igg + 1
635  DO iy=1, ny
636  DO ix=1, nx
637  IF ( mtemp(iy,ix) .EQ. j ) msplit(iy,ix) = igg
638  END DO
639  END DO
640  END IF
641  END DO
642  !
643  ! 5.b.5 Optional GrADS output
644  !
645 #ifdef W3_O16
646  WRITE ( ndsg ) ((real(msplit(iy,ix)),ix=1,nx),iy=1,ny)
647  ntgrds = ntgrds + 1
648 #endif
649  !
650  !--- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
651  ! 5. Refine grids
652  ! 5.a Set up loop
653  !
654  ALLOCATE ( gstats(ng), gstold(ng), pgrid(ng) )
655  gstats(:)%INSTAT = .true.
656  WRITE (ndso,950)
657  done = .false.
658  !
659  CALL grinfo
660  WRITE (ndso,951) 0, mstats%NMIN, mstats%NMAX, &
661  100.*mstats%RSTD/xmean
662  g0id = '5.a'
663  IF ( mstats%NMIN .EQ. 0 ) GOTO 850
664  ingmin = mstats%NMIN
665  ingmax = mstats%NMAX
666  ingmnc = 0
667  ingmxc = 0
668  inglag = 3
669  stdmin = 100.*mstats%RSTD/xmean
670  nstdlg = 0
671  !
672  DO iit=1, nitmax
673  !
674  IF ( ng .EQ. 1 ) EXIT
675  !
676  mstold = mstats
677  gstold = gstats
678  !
679 #ifdef W3_T
680  WRITE (ndst,9050) 'a', mstats%NMIN, mstats%NMAX, mstats%RSTD
681 #endif
682  !
683  ! 5.b Small grid attempt to merge
684  !
685  IF ( mstats%NMIN .LT. nint(0.45*xmean) ) THEN
686  !
687  CALL grfsml
688  CALL grinfo
689  !
690  g0id = '5.b'
691  IF ( mstold%NMIN .NE. mstats%NMIN ) THEN
692  WRITE (ndso,951) iit, mstats%NMIN, mstats%NMAX, &
693  100.*mstats%RSTD/xmean
694  IF ( mstats%NMIN .EQ. 0 ) GOTO 850
695 #ifdef W3_O16
696  WRITE ( ndsg ) ((real(msplit(iy,ix)),ix=1,nx),iy=1,ny)
697  ntgrds = ntgrds + 1
698 #endif
699 
700  cycle
701  ELSE
702  WRITE (ndso,952) mstats%NMIN, mstats%NMAX, &
703  100.*mstats%RSTD/xmean
704  IF ( mstats%NMIN .EQ. 0 ) GOTO 850
705  END IF
706  !
707  END IF
708  !
709  ! 5.c Trim edges of grids and reassign
710  !
711  CALL grtrim
712  CALL grfill ( 2 )
713  !
714  ! 5.d Attempt to quare-up grid
715  !
716  CALL grinfo ! call needed as GRSQRG uses grid ranges
717 #ifdef W3_T
718  WRITE (ndst,9051) 'd', mstats%NMIN, mstats%NMAX, mstats%RSTD
719 #endif
720  CALL grsqrg
721  CALL grfill ( 1 )
722  !
723  ! 5.e Remove mid-sea points sticking out of grid
724  ! Call more than once to remove most .....
725  !
726  ok = .true.
727  !
728  DO jj=1, 4
729  CALL grsngl ( ok )
730  END DO
731  !
732  ! 5.f Remove parts of grid separated from main body, and attachable to
733  ! other grids.
734  !
735  CALL grsepa ( ok , 0.10 )
736  IF ( .NOT. ok ) THEN
737  CALL grfill ( 1 )
738  ok = .true.
739  END IF
740  !
741  ! 5.g Re-compute grid stats
742  !
743  CALL grinfo
744  WRITE (ndso,951) iit, mstats%NMIN, mstats%NMAX, &
745  100.*mstats%RSTD/xmean
746 #ifdef W3_T
747  WRITE (ndst,9051) 'g', mstats%NMIN, mstats%NMAX, mstats%RSTD
748 #endif
749  !
750  g0id = '5.g'
751  IF ( mstats%NMIN .EQ. 0 ) GOTO 850
752  !
753  ! 5.h Optional GrADS output
754  !
755 #ifdef W3_O16
756  WRITE ( ndsg ) ((real(msplit(iy,ix)),ix=1,nx),iy=1,ny)
757  ntgrds = ntgrds + 1
758 #endif
759  !
760  ! 5.i Convergence tests
761  ! ... The quick one
762  !
763  IF ( 100.*mstats%RSTD/xmean .LE. starg ) THEN
764  WRITE (ndso,959)
765  EXIT
766  END IF
767  !
768  ! ... Monitoring convergence ....
769  !
770  IF ( 100.*mstats%RSTD/xmean .LT. 1.0001*stdmin ) THEN
771  IF ( nstdlg .LT. mstdlg ) THEN
772  nstdlg = 0
773  ELSE
774  WRITE (ndso,959)
775  EXIT
776  END IF
777  stdmin = 100.*mstats%RSTD/xmean
778  ELSE
779  nstdlg = nstdlg + 1
780  IF ( nstdlg .GT. mstdlg ) stdmin = 1.01*stdmin
781  END IF
782  !
783  ! ... Check if stuck on min or max
784  !
785  IF ( mstats%NMAX .LT. ingmax ) THEN
786  ingmax = mstats%NMAX
787  ingmxc = 0
788  ELSE
789  ingmxc = ingmxc + 1
790  END IF
791  !
792  IF ( mstats%NMIN .GT. ingmin ) THEN
793  ingmin = mstats%NMIN
794  ingmnc = 0
795  ELSE
796  ingmnc = ingmnc + 1
797  END IF
798  !
799  ! ... Stuck in min ...
800  !
801  IF ( ingmnc .GE. inglag ) THEN
802  !
803 #ifdef W3_T
804  WRITE (ndst,9052) 'MINIMUM'
805 #endif
806  !
807  IF ( real(ingmin) .LT. 0.85*xmean ) THEN
808  !
809 #ifdef W3_T
810  WRITE (ndst,9053) 0.85*xmean / real(ingmin)
811 #endif
812  CALL grfsml
813  CALL grinfo
814  WRITE (ndso,952) mstats%NMIN, mstats%NMAX, &
815  100.*mstats%RSTD/xmean
816  ingmin = mstats%NMIN
817  ingmax = mstats%NMAX
818  ingmnc = 0
819  ingmxc = 0
820  IF ( done ) EXIT
821  !
822 #ifdef W3_T
823  ELSE
824  WRITE (ndst,9054)
825 
826 #endif
827  END IF
828  !
829  END IF
830  !
831  ! ... Stuck in max ...
832  !
833  IF ( ingmxc .GE. inglag ) THEN
834  !
835 #ifdef W3_T
836  WRITE (ndst,9052) 'MAXIMUM'
837 #endif
838  !
839  IF ( real(ingmax) .GT. 1.075*xmean ) THEN
840  !
841 #ifdef W3_T
842  WRITE (ndst,9053) real(ingmax) / ( 1.075*xmean )
843 #endif
844  CALL grinfo
845  WRITE (ndso,952) mstats%NMIN, mstats%NMAX, &
846  100.*mstats%RSTD/xmean
847  ingmin = mstats%NMIN
848  ingmax = mstats%NMAX
849  ingmnc = 0
850  ingmxc = 0
851  IF ( done ) EXIT
852  !
853 #ifdef W3_T
854  ELSE
855  WRITE (ndst,9054)
856 
857 #endif
858  END IF
859  !
860  END IF
861  !
862  END DO
863  !
864  ! 5.j Test output
865  !
866  WRITE (ndso,955)
867  ALLOCATE ( isnext(ng) )
868  isnext = .true.
869  CALL grinfo
870  !
871  DO jj=1, ng
872  minnr = nsea + 1
873  DO j=1, ng
874  IF ( isnext(j) .AND. gstats(j)%NPTS.LT.minnr ) THEN
875  minnr = gstats(j)%NPTS
876  ig = j
877  END IF
878  END DO
879  isnext(ig) = .false.
880  WRITE (ndst,956) ig, gstats(ig)%STRADLE, gstats(ig)%NPTS, &
881  gstats(ig)%NXL, gstats(ig)%NXH, &
882  gstats(ig)%NYL, gstats(ig)%NYH
883  END DO
884  !
885  DEALLOCATE ( isnext )
886  !
887  !--- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
888  ! 6. Generate output to make separate grids
889  ! 6.a Set up loop
890  !
891  WRITE (ndso,960)
892  !
893  zbdum = 999.
894  IF ( maxval(zb) .LT. -0.11 ) THEN
895  zbmin = -0.1
896  ELSE
897  zbmin = maxval(zb) + 1.
898  zbdum = max( zbdum , zbmin+1 )
899  END IF
900  !
901  j1 = len_trim(fext)
902  j2 = 1 + int(log10(real(ng)+0.5))
903  WRITE (nrfmt,'(A2,I1,A1,I1,A1)') '(I', j2, '.', j2, ')'
904  !
905  IF ( j1 + j2 + 2 .LE. 10 ) THEN
906  fname = fext(:j1) // '_p'
907  j3 = j1 + 3
908  ELSE
909  fname = 'part_'
910  j3 = 6
911  END IF
912  j4 = j3 + j2 - 1
913  !
914  nseat = 0
915  !
916  DO ig=1, ng
917  !
918  !
919  ! 6.b Extract grid including halo
920  !
921  WRITE (ndso,961) ig
922  CALL gr1grd
923  nseat = nseat + pgrid(ig)%NSEA
924  !
925  WRITE (aext,nrfmt) ig
926  fname(j3:j4) = aext(:j2)
927  j = len_trim(fnmpre)
928  !
929  ! 6.c Writing bottom file
930  !
931  j5 = j4 + 4
932  fname(j4+1:j5) = '.bot'
933  WRITE (ndso,962) fname(:j5)
934  !
935  IF ( idfm1 .EQ. 3 ) THEN
936  OPEN (ndsm,file=fnmpre(:j)//fname(:j5), &
937  form='UNFORMATTED', convert=file_endian,err=860,iostat=ierr)
938  ELSE
939  OPEN (ndsm,file=fnmpre(:j)//fname(:j5), err=860,iostat=ierr)
940  END IF
941  rewind(ndsm)
942  CALL outa2r ( pgrid(ig)%ZBIN, pgrid(ig)%NX, pgrid(ig)%NY, &
943  1, pgrid(ig)%NX, 1, pgrid(ig)%NY, ndsm, ndst, &
944  ndse, idfm1, rform1, idla1, vsc1, 0.0 )
945  CLOSE (ndsm)
946  !
947  ! 6.d Writing obstruction file
948  !
949  j5 = j4 + 5
950  fname(j4+1:j5) = '.obst'
951  !
952  IF ( trflag .EQ. 0 ) THEN
953  WRITE (ndso,963) fname(:j5)
954  ELSE
955  WRITE (ndso,962) fname(:j5)
956  !
957  IF ( idfm2 .EQ. 3 ) THEN
958  OPEN (ndsm,file=fnmpre(:j)//fname(:j5), &
959  form='UNFORMATTED', convert=file_endian,err=860,iostat=ierr)
960  ELSE
961  OPEN (ndsm,file=fnmpre(:j)//fname(:j5), &
962  err=860,iostat=ierr)
963  END IF
964  rewind(ndsm)
965  CALL outa2r ( pgrid(ig)%OBSX, pgrid(ig)%NX, pgrid(ig)%NY, &
966  1, pgrid(ig)%NX, 1, pgrid(ig)%NY, ndsm, &
967  ndst, ndse, idfm2, rform2, idla2, vsc2, 0.0 )
968  CALL outa2r ( pgrid(ig)%OBSY, pgrid(ig)%NX, pgrid(ig)%NY, &
969  1, pgrid(ig)%NX, 1, pgrid(ig)%NY, ndsm, &
970  ndst, ndse, idfm2, rform2, idla2, vsc2, 0.0 )
971  CLOSE (ndsm)
972  !
973  END IF
974  !
975  ! 6.e Writing mask file
976  !
977  j5 = j4 + 5
978  fname(j4+1:j5) = '.mask'
979  WRITE (ndso,962) fname(:j5)
980  !
981  IF ( idfm3 .EQ. 3 ) THEN
982  OPEN (ndsm,file=fnmpre(:j)//fname(:j5), &
983  form='UNFORMATTED', convert=file_endian,err=860,iostat=ierr)
984  ELSE
985  OPEN (ndsm,file=fnmpre(:j)//fname(:j5), err=860,iostat=ierr)
986  END IF
987  rewind(ndsm)
988  CALL outa2i ( pgrid(ig)%MASK, pgrid(ig)%NX, pgrid(ig)%NY, &
989  1, pgrid(ig)%NX, 1, pgrid(ig)%NY, ndsm, ndst, &
990  ndse, idfm3, rform3, idla3, vsc3, 0 )
991  CLOSE (ndsm)
992  !
993  ! 6.f Writing input file
994  !
995  j5 = j4 + 5
996  fname(j4+1:j5) = '.tmpl'
997  WRITE (ndso,962) fname(:j5)
998  !
999  OPEN (ndsm,file=fnmpre(:j)//fname(:j5), err=860,iostat=ierr)
1000  !
1001  gname(31-j2:30) = aext
1002  gname(30-j2:30-j2) = 'p'
1003  WRITE (ndsm,965) gname, sig(2)/sig(1), tpiinv*sig(1), nk, &
1004  nth, th(1)/dth, fldry, flcx, flcy, flcth, &
1005  flck, flsou, dtmax, dtcfl, dtcfli, dtmin
1006  j5 = len_trim(rform1)
1007  IF ( real(pgrid(ig)%NX) * pgrid(ig)%SX .LT. 359.9 ) THEN
1008  ptclse = 'NONE'
1009  ELSE
1010  ptclse = idclse
1011  END IF
1012  WRITE (ndsm,966) idgrid, flagll, ptclse, &
1013  pgrid(ig)%NX, pgrid(ig)%NY, &
1014  pgrid(ig)%SX, pgrid(ig)%SY, &
1015  pgrid(ig)%X0, pgrid(ig)%Y0, &
1016  zbmin, dmin, vsc1, idla1, idfm1, &
1017  rform1(:j5), fname(:j4)//'.bot'
1018  IF ( trflag .NE. 0 ) THEN
1019  j5 = len_trim(rform2)
1020  WRITE (ndsm,967) vsc2,idla2, idfm2, rform2(:j5), &
1021  fname(:j4)//'.obst'
1022  END IF
1023  j5 = len_trim(rform3)
1024  WRITE (ndsm,968) idla3, idfm3, rform3(:j5), fname(:j4)//'.mask'
1025  CLOSE (ndsm)
1026  !
1027  END DO
1028  !
1029  WRITE (ndso,969) 100. * (real(nseat)/real(nsea)-1.)
1030  !
1031  !--- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1032  ! 7. Write part of ww3_multi.inp
1033  !
1034  j5 = 11+j1+j2
1035  iname(:j5) = 'ww3_multi.'//fext(:j1)//'.'//aext(:j2)
1036  OPEN (ndsm,file=fnmpre(:j)//iname(:j5), err=870,iostat=ierr)
1037  !
1038  DO ig=1, ng
1039  WRITE (aext,nrfmt) ig
1040  fname(j3:j4) = aext(:j2)
1041  IF ( frflag ) THEN
1042  WRITE (ndsm,970) fname(:j4), &
1043  fracl + real(ig-1)*(frach-fracl)/real(ng), &
1044  fracl + real( ig )*(frach-fracl)/real(ng)
1045  ELSE
1046  WRITE (ndsm,970) fname(:j4), fracl, frach
1047  END IF
1048  END DO
1049  !
1050  CLOSE (ndsm)
1051  !
1052  !--- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1053  ! 8. Write mask file (no halo)
1054  !
1055  j5 = 10+j1+j2
1056  iname(:j5) = 'ww3_mask.'//fext(:j1)//'.'//aext(:j2)
1057  OPEN (ndsm,file=fnmpre(:j)//iname(:j5), err=870,iostat=ierr)
1058  !
1059  DO iy=1, ny
1060  WRITE (ndsm,980) msplit(iy,:)
1061  END DO
1062  !
1063 #ifdef W3_O16
1064  CLOSE ( ndsg )
1065 #endif
1066  !
1067 #ifdef W3_O16
1068  OPEN ( ndsg,file='ww3.ctl')
1069  WRITE (ndsg,985) nx, x0, sx, ny, y0, sy, ntgrds
1070  CLOSE ( ndsg )
1071 #endif
1072  !
1073  !--- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1074  ! 9. End of program
1075  !
1076  GOTO 888
1077  !
1078  ! Error escape locations
1079  !
1080 800 CONTINUE
1081  WRITE (ndse,1000) ierr
1082  CALL extcde ( 40 )
1083  !
1084 801 CONTINUE
1085  WRITE (ndse,1001)
1086  CALL extcde ( 41 )
1087  !
1088 802 CONTINUE
1089  WRITE (ndse,1002) ierr
1090  CALL extcde ( 42 )
1091  !
1092 820 CONTINUE
1093  WRITE (ndse,1020) gtype
1094  CALL extcde ( 20 )
1095  !
1096 821 CONTINUE
1097  WRITE (ndse,1021) gtype
1098  CALL extcde ( 21 )
1099  !
1100 822 CONTINUE
1101  WRITE (ndse,1022) iclose
1102  CALL extcde ( 22 )
1103  !
1104 823 CONTINUE
1105  WRITE (ndse,1023) iclose
1106  CALL extcde ( 23 )
1107  !
1108 824 CONTINUE
1109  WRITE (ndse,1024)
1110  CALL extcde ( 24 )
1111  !
1112 825 CONTINUE
1113  WRITE (ndse,1025) mingrd, nsea
1114  CALL extcde ( 25 )
1115  !
1116 830 CONTINUE
1117  WRITE (ndse,1030)
1118  CALL extcde ( 30 )
1119  !
1120 850 CONTINUE
1121  WRITE (ndse,1050) g0id
1122  CALL extcde ( 50 )
1123  !
1124 860 CONTINUE
1125  WRITE (ndse,1060) fnmpre(:j)//fname(:j5), ierr
1126  CALL extcde ( 60 )
1127  !
1128 870 CONTINUE
1129  WRITE (ndse,1070) fnmpre(:j)//iname(:j5), ierr
1130  CALL extcde ( 70 )
1131  !
1132 888 CONTINUE
1133  WRITE (ndso,999)
1134  !
1135  ! Formats
1136  !
1137 900 FORMAT (/15x,' *** WAVEWATCH III Grid splitting *** '/ &
1138  15x,'=========================================='/)
1139 901 FORMAT ( ' Comment character is ''',a,''''/)
1140 902 FORMAT ( ' Grid ID : ',a/ &
1141  ' Grid name : ',a)
1142 903 FORMAT ( ' Grid type : ',a)
1143 904 FORMAT ( ' Closure : ',a)
1144 905 FORMAT ( ' Grid size : ',i4,' x',i4,' (',i8,')'/)
1145  !
1146 930 FORMAT ( ' Generating ',i3,' grids'/ &
1147  ' No more than',i4,' refinement iterations'/ &
1148  ' Grid point count std target (%) :',f6.2/ &
1149  ' Halo per sub grid extended by',i3,' grid point.')
1150 931 FORMAT ( ' Format info for bottom file :',2i2,f12.4,2x,a)
1151 932 FORMAT ( ' Format info for obstruction file not used')
1152 933 FORMAT ( ' Format info for obstruction file :',2i2,f12.4,2x,a)
1153 934 FORMAT ( ' Format info for mask file :',2i2,i7,7x,a)
1154 935 FORMAT ( ' Part of cummunicator to be used :',2f7.4)
1155 936 FORMAT ( ' Not running grids side-by-side'/ &
1156  ' *** NON CONVENTIONAL OPERATION ***'/)
1157  !
1158 950 FORMAT (/' Iterations:'/ &
1159  ' nr min max std (%) '/ &
1160  ' ---------------------------------')
1161 951 FORMAT (2x,i5,2i8,2f10.2)
1162 952 FORMAT (2x,5x,2i8,2f10.2)
1163 955 FORMAT (/' Resulting grids:'/ &
1164  ' grid stradle points range X range Y '/ &
1165  ' ---------------------------------------------')
1166 956 FORMAT ( ' ',i4,5x,l1,2x,i7,4i5)
1167 959 FORMAT ( ' Convergence reached')
1168  !
1169 960 FORMAT (/' Generating grid data:'/ &
1170  ' ---------------------------------------------')
1171 961 FORMAT ( ' Extracting data for grid',i4)
1172 962 FORMAT ( ' Writing file ',a)
1173 963 FORMAT ( ' File ',a,' not requested')
1174  !
1175 970 FORMAT ( ' ''',a,''' ''LEV'' ''CUR'' ''WND'' ''ICE''', &
1176  ' ''D1'' ''D2'' ''D3'' RANK GROUP',2f10.7,' BFLAG')
1177  !
1178 980 FORMAT (1x,360i2)
1179  !
1180 #ifdef W3_O16
1181 985 FORMAT ( 'DSET ww3.ww3_gspl'/ &
1182  'TITLE WAVEWATCH III grid splitting data'/ &
1183  'OPTIONS sequential'/ &
1184  'UNDEF -999.9'/ &
1185  'XDEF ',i6,' LINEAR ',2f12.5/ &
1186  'YDEF ',i6,' LINEAR ',2f12.5/ &
1187  'ZDEF 1 LINEAR 1000.00000 1.00000'/ &
1188  'TDEF ',i6,' LINEAR 00:00 06JUN1968 1HR'/ &
1189  'VARS 1'/ &
1190  'MAP 0 99 grid use map '/ &
1191  'ENDVARS')
1192 #endif
1193  !
1194 965 FORMAT ( '$ -------------------------------------', &
1195  '------------------------------- $'/ &
1196  '$ WAVEWATCH III Grid preprocessor input', &
1197  ' file $'/ &
1198  '$ -------------------------------------', &
1199  '------------------------------- $'/ &
1200  ' ''',a,''''/'$'/ &
1201  ' ',f8.4,f10.6,2i6,f8.4/' ',6l2/' ',4f12.4/ &
1202  '$ NAMELISTS'/'$')
1203 966 FORMAT ( ' ''',a4,''' ',l1,' ''',a4,''''/1x,i8,i12/ &
1204  4x,2f12.6,' 1.0'/4x,2f12.6,' 1.0'/2f8.2,' 20', &
1205  f12.6,2i2,' ''',a,''' ''NAME'' ''',a,'''')
1206 967 FORMAT ( 18x,'30',f12.6,2i2,' ''',a,''' ''NAME'' ''',a,'''' )
1207 968 FORMAT ( 18x,'40',12x,2i2,' ''',a,''' ''NAME'' ''',a,''''/'$'/ &
1208  '$ Note: cannot make output boundary points here'/'$'/ &
1209  ' 0. 0. 0. 0. 0'/ &
1210  '$ -------------------------------------', &
1211  '------------------------------- $'/ &
1212  '$ End of input file ', &
1213  ' $'/ &
1214  '$ -------------------------------------', &
1215  '------------------------------- $')
1216  !
1217 969 FORMAT (/' Grid point inflation',f7.2,'%')
1218  !
1219 999 FORMAT(//' End of program '/ &
1220  ' ========================================='/ &
1221  ' WAVEWATCH III Grid splitting '/)
1222  !
1223 1000 FORMAT (/' *** WAVEWATCH III ERROR IN W3GSPL : '/ &
1224  ' ERROR IN OPENING INPUT FILE'/ &
1225  ' IOSTAT =',i5/)
1226  !
1227 1001 FORMAT (/' *** WAVEWATCH III ERROR IN W3GSPL : '/ &
1228  ' PREMATURE END OF INPUT FILE'/)
1229  !
1230 1002 FORMAT (/' *** WAVEWATCH III ERROR IN W3GSPL : '/ &
1231  ' ERROR IN READING FROM INPUT FILE'/ &
1232  ' IOSTAT =',i5/)
1233  !
1234 1020 FORMAT (/' *** WAVEWATCH III ERROR IN W3GSPL : '/ &
1235  ' SPLITTING NOT AVAILABLE FOR GRID TYPE'/ &
1236  ' GTYPE =',i5/)
1237  !
1238 1021 FORMAT (/' *** WAVEWATCH III ERROR IN W3GSPL : '/ &
1239  ' GRID TYPE NOT RECOGNIZED'/ &
1240  ' GTYPE =',i5/)
1241  !
1242 1022 FORMAT (/' *** WAVEWATCH III ERROR IN W3GSPL : '/ &
1243  ' SPLITTING NOT AVAILABLE FOR CLOSURE TYPE'/ &
1244  ' ICLOSE =',i5/)
1245  !
1246 1023 FORMAT (/' *** WAVEWATCH III ERROR IN W3GSPL : '/ &
1247  ' CLOSURE TYPE NOT RECOGNIZED'/ &
1248  ' ICLOSE =',i5/)
1249  !
1250 1024 FORMAT (/' *** WAVEWATCH III ERROR IN W3GSPL : '/ &
1251  ' NO ACTIVE SEA POINT IN GRID'/)
1252  !
1253 1025 FORMAT (/' *** WAVEWATCH III ERROR IN W3GSPL : '/ &
1254  ' WRONG NUMBER OF SEA POINTS'/ &
1255  ' MINGRD, NSEA =',2i7/)
1256  !
1257 1030 FORMAT (/' *** WAVEWATCH III ERROR IN W3GSPL : '/ &
1258  ' ILLEGAL PART OF COMMUNICATOR REQUESTED'/)
1259  !
1260 1050 FORMAT (/' *** WAVEWATCH III ERROR IN W3GSPL : '/ &
1261  ' SHOULD NOT HAVE ZERO GRID SIZE (',a,') ...'/)
1262  !
1263 1060 FORMAT (/' *** WAVEWATCH III ERROR IN W3GSPL : '/ &
1264  ' ERROR IN OPENING FILE ',a/ &
1265  ' IOSTAT =',i5/)
1266  !
1267 1070 FORMAT (/' *** WAVEWATCH III ERROR IN W3GSPL : '/ &
1268  ' ERROR IN OPENING FILE ',a/ &
1269  ' IOSTAT =',i5/)
1270  !
1271 #ifdef W3_T
1272 9040 FORMAT ( 'TEST W3GSPL: CHECKERBOARD X-Y:',2i8)
1273 9041 FORMAT ( 'TEST W3GSPL: FILLING CHECKERBOARD TRY:',i3/ &
1274  ' GRID, IGX0, IGXN, IGY0, IGYN, POINTS ')
1275 9042 FORMAT ( ' ',i6,2(2i8), i8,2x,a)
1276 9043 FORMAT ( 'TEST W3GSPL: CHECKERBOARD GRIDS:',i4,' (',i4,')')
1277 9044 FORMAT ( ' SMALLEST SIZE/GRID:',i8,i4)
1278 9045 FORMAT ( ' SMALLEST NEIGHBOR :',i8,i4)
1279 9046 FORMAT ( ' GRID',i4', MERGED WITH GRID',i4)
1280 9047 FORMAT ( ' ',i6,i8)
1281 9048 FORMAT ( ' GRID',i4', IS ISOLATED, LEFT UNCHANGED')
1282 9049 FORMAT ( 'TEST W3GSPL: CHECKERBOARD CONSOLIDATED ON',i4,' GRIDS')
1283 #endif
1284  !
1285 #ifdef W3_T
1286 9050 FORMAT ( 'TEST W3GSPL',a,': MIN, MAX, STD:',2i8,f10.2)
1287 9051 FORMAT ( ' ',a,': MIN, MAX, STD:',2i8,f10.2)
1288 #endif
1289 9052 FORMAT ( 'TEST W3GSPL: STUCK ON ',a,' GRID SIZE')
1290 9053 FORMAT ( ' OUT OF RANGE, PROCESSING (',f6.3,')')
1291 9054 FORMAT ( ' IN RANGE, NO ACTION')
1292  !/
1293  !/ Embedded subroutines ---------------------------------------------- /
1294  !/
1295 CONTAINS
1296  !/ ------------------------------------------------------------------- /
1297 
1301  SUBROUTINE grinfo
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  !/
1469  END SUBROUTINE grinfo
1470  !/ ------------------------------------------------------------------- /
1471 
1479  SUBROUTINE grtrim
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  !/
1643  END SUBROUTINE grtrim
1644  !/ ------------------------------------------------------------------- /
1645 
1651  SUBROUTINE grfill ( ND )
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  !/
1904  END SUBROUTINE grfill
1905  !/ ------------------------------------------------------------------- /
1906 
1912  SUBROUTINE grlost
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  !/
2006  END SUBROUTINE grlost
2007  !/ ------------------------------------------------------------------- /
2008 
2016  SUBROUTINE grsqrg
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  !/
2110  END SUBROUTINE grsqrg
2111  !/ ------------------------------------------------------------------- /
2112 
2122  SUBROUTINE grsngl ( OK )
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  !/
2282  END SUBROUTINE grsngl
2283  !/ ------------------------------------------------------------------- /
2284 
2293  SUBROUTINE grsepa ( OK, FRAC )
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  !/
2628  END SUBROUTINE grsepa
2629  !/ ------------------------------------------------------------------- /
2630 
2638  SUBROUTINE grfsml
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  !/
2917  END SUBROUTINE grfsml
2918 
2919 
2923  !/ ------------------------------------------------------------------- /
2924  SUBROUTINE grflrg
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  !/
3074  END SUBROUTINE grflrg
3075  !/ ------------------------------------------------------------------- /
3076 
3083  SUBROUTINE gr1grd
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  !/
3471  END SUBROUTINE gr1grd
3472  !/
3473  !/ End of W3GSPL ----------------------------------------------------- /
3474  !/
3475 END PROGRAM w3gspl
w3servmd::nextln
subroutine nextln(CHCKC, NDSI, NDSE)
Definition: w3servmd.F90:222
w3adatmd
Define data structures to set up wave model auxiliary data for several models simultaneously.
Definition: w3adatmd.F90:26
w3arrymd::outa2i
subroutine outa2i(ARRAY, MX, MY, LX, HX, LY, HY, NDS, NDST, NDSE, IDFM, RFORM, IDLA, VSC, VOF)
Definition: w3arrymd.F90:627
grflrg
subroutine grflrg
Like GRFSML for largest grid ...
Definition: ww3_gspl.F90:2925
w3arrymd::outa2r
subroutine outa2r(ARRAY, MX, MY, LX, HX, LY, HY, NDS, NDST, NDSE, IDFM, RFORM, IDLA, VSC, VOF)
Definition: w3arrymd.F90:465
grfsml
subroutine grfsml
Subroutine called when lowest grid size is stuck.
Definition: ww3_gspl.F90:2639
grsepa
subroutine grsepa(OK, FRAC)
Remove smaller grid parts.
Definition: ww3_gspl.F90:2294
grsqrg
subroutine grsqrg
Attempt to square-up grid.
Definition: ww3_gspl.F90:2017
w3gdatmd::sy
real, pointer sy
Definition: w3gdatmd.F90:1183
w3gdatmd::ny
integer, pointer ny
Definition: w3gdatmd.F90:1097
w3odatmd::fnmpre
character(len=80) fnmpre
Definition: w3odatmd.F90:330
grsngl
subroutine grsngl(OK)
Remove seapoints with only one adjacent point in same grid.
Definition: ww3_gspl.F90:2123
grinfo
subroutine grinfo
Compile statistical info on all sub grids (no halo).
Definition: ww3_gspl.F90:1302
part_grid
Definition: ww3_gspl.F90:238
w3gdatmd::w3setg
subroutine w3setg(IMOD, NDSE, NDST)
Definition: w3gdatmd.F90:2152
w3odatmd::ndse
integer, pointer ndse
Definition: w3odatmd.F90:456
w3adatmd::w3seta
subroutine w3seta(IMOD, NDSE, NDST)
Select one of the WAVEWATCH III grids / models.
Definition: w3adatmd.F90:2645
w3gdatmd::x0
real, pointer x0
Definition: w3gdatmd.F90:1183
w3gdatmd::nsea
integer, pointer nsea
Definition: w3gdatmd.F90:1097
w3servmd
Definition: w3servmd.F90:3
constants::tpiinv
real, parameter tpiinv
TPIINV Inverse of 2*Pi.
Definition: constants.F90:74
w3odatmd::w3seto
subroutine w3seto(IMOD, NDSERR, NDSTST)
Definition: w3odatmd.F90:1523
w3odatmd
Definition: w3odatmd.F90:3
grlost
subroutine grlost
Reassign unassigned grid points to grids.
Definition: ww3_gspl.F90:1913
w3adatmd::w3naux
subroutine w3naux(NDSE, NDST)
Set up the number of grids to be used.
Definition: w3adatmd.F90:704
w3iogrmd::w3iogr
subroutine w3iogr(INXOUT, NDSM, IMOD, FEXT ifdef W3_ASCII
Reading and writing of the model definition file.
Definition: w3iogrmd.F90:117
file
file(STRINGS ${CMAKE_BINARY_DIR}/switch switch_strings) separate_arguments(switches UNIX_COMMAND $
Definition: CMakeLists.txt:3
w3iogrmd
Reading/writing of model definition file.
Definition: w3iogrmd.F90:20
w3servmd::strace
subroutine strace(IENT, SNAME)
Definition: w3servmd.F90:148
w3odatmd::ndso
integer, pointer ndso
Definition: w3odatmd.F90:456
w3gdatmd::w3nmod
subroutine w3nmod(NUMBER, NDSE, NDST, NAUX)
Definition: w3gdatmd.F90:1433
w3arrymd
Definition: w3arrymd.F90:3
w3gdatmd::y0
real, pointer y0
Definition: w3gdatmd.F90:1183
gr1grd
subroutine gr1grd
Extract single grid from master map.
Definition: ww3_gspl.F90:3084
w3gdatmd::sx
real, pointer sx
Definition: w3gdatmd.F90:1183
w3odatmd::ndst
integer, pointer ndst
Definition: w3odatmd.F90:456
constants
Define some much-used constants for global use (all defined as PARAMETER).
Definition: constants.F90:20
w3gdatmd
Definition: w3gdatmd.F90:16
grtrim
subroutine grtrim
Trim edges of all grids where they are next to another grid or next to unassigned grid points.
Definition: ww3_gspl.F90:1480
grfill
subroutine grfill(ND)
Reassign unassigned grid points to grids, starting with the smallest grids.
Definition: ww3_gspl.F90:1652
stats_grid
Definition: ww3_gspl.F90:228
constants::file_endian
character(*), parameter file_endian
FILE_ENDIAN Filled by preprocessor with 'big_endian', 'little_endian', or 'native'.
Definition: constants.F90:86
w3servmd::extcde
subroutine extcde(IEXIT, UNIT, MSG, FILE, LINE, COMM)
Definition: w3servmd.F90:736
w3odatmd::w3nout
subroutine w3nout(NDSERR, NDSTST)
Definition: w3odatmd.F90:561
w3servmd::itrace
subroutine itrace(NDS, NMAX)
Definition: w3servmd.F90:91
w3gspl
program w3gspl
Grid splitting program.
Definition: ww3_gspl.F90:26
w3gdatmd::nx
integer, pointer nx
Definition: w3gdatmd.F90:1097
stats_mean
Definition: ww3_gspl.F90:233