WAVEWATCH III  beta 0.0.1
w3uqckmd Module Reference

Portable ULTIMATE QUICKEST schemes. More...

Functions/Subroutines

subroutine w3qck1 (MX, MY, NX, NY, CFLL, Q, CLOSE, INC, MAPACT, NACT, MAPBOU, NB0, NB1, NB2, NDSE, NDST)
 Preform one-dimensional propagation in a two-dimensional space with irregular boundaries and regular grid. More...
 
subroutine w3qck2 (MX, MY, NX, NY, VELO, DT, DX1, DX2, Q, CLOSE, INC, MAPACT, NACT, MAPBOU, NB0, NB1, NB2, NDSE, NDST)
 Like W3QCK1 with variable grid spacing. More...
 
subroutine w3qck3 (MX, MY, NX, NY, CFLL, TRANS, Q, CLOSE, INC, MAPACT, NACT, MAPBOU, NB0, NB1, NB2, NDSE, NDST)
 Like W3QCK1 with cell transparencies added. More...
 

Detailed Description

Portable ULTIMATE QUICKEST schemes.

Author
H. L. Tolman
Date
27-May-2014

Function/Subroutine Documentation

◆ w3qck1()

subroutine w3uqckmd::w3qck1 ( integer, intent(in)  MX,
integer, intent(in)  MY,
integer, intent(in)  NX,
integer, intent(in)  NY,
real, dimension(my*(mx+1)), intent(inout)  CFLL,
real, dimension(1-my:my*(mx+2)), intent(inout)  Q,
logical, intent(in)  CLOSE,
integer, intent(in)  INC,
integer, dimension(my*mx), intent(in)  MAPACT,
integer, intent(in)  NACT,
integer, dimension(my*mx), intent(in)  MAPBOU,
integer, intent(in)  NB0,
integer, intent(in)  NB1,
integer, intent(in)  NB2,
integer, intent(in)  NDSE,
integer, intent(in)  NDST 
)

Preform one-dimensional propagation in a two-dimensional space with irregular boundaries and regular grid.

ULTIMATE QUICKEST scheme (see manual).

Note that the check on monotonous behavior of QCN is performed using weights CFAC, to avoid the need for IF statements.

Called by: W3KTP2 Propagation in spectral space.

This routine can be used independently from WAVEWATCH III.

Parameters
[in]MXField dimensions, if grid is 'closed' or circular, MX is the closed dimension.
[in]MYField dimension (See MX)
[in]NXPart of field actually used
[in]NYPart of field actually used
[in]INCIncrement in 1-D array corresponding to increment in 2-D space.
[in]MAPACTList of active grid points.
[in]NACTSize of MAPACT.
[in]MAPBOUMap with boundary information (See W3MAP2).
[in]NB0Counter in MAPBOU
[in]NB1Counter in MAPBOU
[in]NB2Counter in MAPBOU
[in]NDSEError output unit number.
[in]NDSTTest output unit number.
[in,out]CFLLLocal Courant numbers (MY, MX+1)
[in,out]QPropagated quantity (MY,0:MX+2)
[in]CLOSEFlag for closed 'X' dimension.
Author
H. L. Tolman
Date
30-Oct-2009

Definition at line 120 of file w3uqckmd.F90.

120  !/
121  !/ +-----------------------------------+
122  !/ | WAVEWATCH III NOAA/NCEP |
123  !/ | H. L. Tolman |
124  !/ | FORTRAN 90 |
125  !/ | Last update : 30-Oct-2009 |
126  !/ +-----------------------------------+
127  !/
128  !/ 11-Mar-1997 : Final FORTRAN 77 ( version 1.18 )
129  !/ 15-Dec-1999 : Upgrade to FORTRAN 90 ( version 2.00 )
130  !/ 15-Feb-2001 : Unit numbers added to par list. ( version 2.08 )
131  !/ 05-Mar-2008 : Added NEC sxf90 compiler directives.
132  !/ (Chris Bunney, UK Met Office) ( version 3.13 )
133  !/ 30-Oct-2009 : Fixed "Called by" doc line. ( version 3.14 )
134  !/ (T. J. Campbell, NRL)
135  !/
136  ! 1. Purpose :
137  !
138  ! Preform one-dimensional propagation in a two-dimensional space
139  ! with irregular boundaries and regular grid.
140  !
141  ! 2. Method :
142  !
143  ! ULTIMATE QUICKEST scheme (see manual).
144  !
145  ! Note that the check on monotonous behavior of QCN is performed
146  ! using weights CFAC, to avoid the need for IF statements.
147  !
148  ! 3. Parameters :
149  !
150  ! Parameter list
151  ! ----------------------------------------------------------------
152  ! MX,MY Int. I Field dimensions, if grid is 'closed' or
153  ! circular, MX is the closed dimension.
154  ! NX,NY Int. I Part of field actually used.
155  ! CFLL R.A. I Local Courant numbers. (MY, MX+1)
156  ! Q R.A. I/O Propagated quantity. (MY,0:MX+2)
157  ! CLOSE Log. I Flag for closed 'X' dimension'
158  ! INC Int. I Increment in 1-D array corresponding to
159  ! increment in 2-D space.
160  ! MAPACT I.A. I List of active grid points.
161  ! NACT Int. I Size of MAPACT.
162  ! MAPBOU I.A. I Map with boundary information (see W3MAP2).
163  ! NBn Int. I Counters in MAPBOU.
164  ! NDSE Int. I Error output unit number.
165  ! NDST Int. I Test output unit number.
166  ! ----------------------------------------------------------------
167  ! - CFLL amd Q need only bee filled in the (MY,MX) range,
168  ! extension is used internally for closure.
169  ! - CFLL and Q are defined as 1-D arrays internally.
170  !
171  ! 4. Subroutines used :
172  !
173  ! STRACE Service routine.
174  !
175  ! 5. Called by :
176  !
177  ! W3KTP2 Propagation in spectral space
178  !
179  ! 6. Error messages :
180  !
181  ! None.
182  !
183  ! 7. Remarks :
184  !
185  ! - This routine can be used independently from WAVEWATCH III.
186  !
187  ! 8. Structure :
188  !
189  ! ------------------------------------------------------
190  ! 1. Initialize aux. array FLA.
191  ! 2. Fluxes for central points (3rd order + limiter).
192  ! 3. Fluxes boundary point above (1st order).
193  ! 4. Fluxes boundary point below (1st order).
194  ! 5. Closure of 'X' if required
195  ! 6. Propagate.
196  ! ------------------------------------------------------
197  !
198  ! 9. Switches :
199  !
200  ! !/S Enable subroutine tracing.
201  ! !/T Enable test output.
202  ! !/T0 Test output input/output fields.
203  ! !/T1 Test output fluxes.
204  ! !/T2 Test output integration.
205  !
206  ! 10. Source code :
207  !
208  !/ ------------------------------------------------------------------- /
209  IMPLICIT NONE
210  !/
211  !/ ------------------------------------------------------------------- /
212  !/ Parameter list
213  !/
214  INTEGER, INTENT(IN) :: MX, MY, NX, NY, INC, MAPACT(MY*MX), &
215  NACT, MAPBOU(MY*MX), NB0, NB1, NB2, &
216  NDSE, NDST
217  REAL, INTENT(INOUT) :: CFLL(MY*(MX+1)), Q(1-MY:MY*(MX+2))
218  LOGICAL, INTENT(IN) :: CLOSE
219  !/
220  !/ ------------------------------------------------------------------- /
221  !/ Local parameters
222  !/
223  INTEGER :: IXY, IP, IXYC, IXYU, IXYD, IY, IX, &
224  IAD00, IAD02, IADN0, IADN1, IADN2
225 #ifdef W3_S
226  INTEGER, SAVE :: IENT = 0
227 #endif
228 #ifdef W3_T1
229  INTEGER :: IX2, IY2
230 #endif
231  REAL :: CFL, QB, DQ, DQNZ, QCN, QBN, QBR, CFAC
232  REAL :: FLA(1-MY:MY*MX)
233 #ifdef W3_T0
234  REAL :: QMAX
235 #endif
236 #ifdef W3_T1
237  REAL :: QBO, QN
238 #endif
239 #ifdef W3_T2
240  REAL :: QOLD
241 #endif
242  !/
243  !/ ------------------------------------------------------------------- /
244  !/
245 #ifdef W3_S
246  CALL strace (ient, 'W3QCK1')
247 #endif
248  !
249 #ifdef W3_T
250  WRITE (ndst,9000) mx, my, nx, ny, CLOSE, inc, nb0, nb1, nb2
251 #endif
252  !
253 #ifdef W3_T0
254  qmax = 0.
255  DO iy=1, ny
256  DO ix=1, nx
257  qmax = max( qmax , q(iy+(ix-1)*my) )
258  END DO
259  END DO
260  qmax = max( 0.01*qmax , 1.e-10 )
261 #endif
262  !
263 #ifdef W3_T0
264  WRITE (ndst,9001) 'CFLL'
265  DO iy=ny,1,-1
266  WRITE (ndst,9002) (nint(100.*cfll(iy+(ix-1)*my)),ix=1,nx)
267  END DO
268  WRITE (ndst,9001) 'Q'
269  DO iy=ny,1,-1
270  WRITE (ndst,9002) (nint(q(iy+(ix-1)*my)/qmax),ix=1,nx)
271  END DO
272  WRITE (ndst,9001) 'MAPACT'
273  WRITE (ndst,9003) (mapact(ixy),ixy=1,nact)
274 #endif
275  !
276  ! 1. Initialize aux. array FLA and closure ------------------------- *
277  !
278  fla = 0.
279  !
280  IF ( CLOSE ) THEN
281 #ifdef W3_T
282  WRITE (ndst,9005)
283 #endif
284  iad00 = -my
285  iad02 = my
286  iadn0 = iad00 + my*nx
287  iadn1 = my*nx
288  iadn2 = iad02 + my*nx
289  DO iy=1, ny
290  q(iy+iad00) = q(iy+iadn0)
291  q(iy+iadn1) = q( iy )
292  q(iy+iadn2) = q(iy+iad02)
293  cfll(iy+iadn1) = cfll( iy )
294  END DO
295  END IF
296  !
297  ! 2. Fluxes for central points ------------------------------------- *
298  ! ( 3rd order + limiter )
299  !
300 #ifdef W3_T1
301  WRITE (ndst,9010)
302  WRITE (ndst,9011) nb0, 'CENTRAL'
303 #endif
304  !
305  DO ip=1, nb0
306  !
307  ixy = mapbou(ip)
308  cfl = 0.5 * ( cfll(ixy) + cfll(ixy+inc) )
309  ixyc = ixy - inc * int( min( 0. , sign(1.1,cfl) ) )
310  qb = 0.5 * ( (1.-cfl)*q(ixy+inc) + (1.+cfl)*q(ixy) ) &
311  - (1.-cfl**2)/6. * (q(ixyc-inc)-2.*q(ixyc)+q(ixyc+inc))
312  !
313  ixyu = ixyc - inc * int( sign(1.1,cfl) )
314  ixyd = 2*ixyc - ixyu
315  dq = q(ixyd) - q(ixyu)
316  dqnz = sign( max(1.e-15,abs(dq)) , dq )
317  qcn = ( q(ixyc) - q(ixyu) ) / dqnz
318  qcn = min( 1.1, max( -0.1 , qcn ) )
319  !
320 #ifdef W3_T1
321  qbo = qb
322 #endif
323  qbn = max( (qb-q(ixyu))/dqnz , qcn )
324  qbn = min( qbn , 1. , qcn/max(1.e-10,abs(cfl)) )
325  qbr = q(ixyu) + qbn*dq
326  cfac = real( int( 2. * abs(qcn-0.5) ) )
327  qb = (1.-cfac)*qbr + cfac*q(ixyc)
328  !
329  fla(ixy) = cfl * qb
330  !
331 #ifdef W3_T1
332  iy = mod( ixy , my )
333  ix = 1 + ixy/my
334  iy2 = mod( ixy+inc , my )
335  ix2 = 1 + (ixy+inc)/my
336  qn = max( qb, qbo, q(ixy-inc), q( ixy ), &
337  q(ixy+inc), q(ixy+2*inc) )
338  IF ( qn .GT. 1.e-10 ) THEN
339  qn = 1. /qn
340  WRITE (ndst,9012) ip, ix, iy, ix2, iy2, &
341  cfl, cfll(ixy), cfll(ixy+inc), &
342  qbo*qn, qb*qn, q(ixy-inc)*qn, q( ixy )*qn, &
343  q(ixy+inc)*qn, q(ixy+2*inc)*qn
344  END IF
345 #endif
346  !
347  END DO
348  !
349  ! 3. Fluxes for points with boundary above ------------------------- *
350  ! ( 1st order without limiter )
351  !
352 #ifdef W3_T1
353  WRITE (ndst,9011) nb1-nb0, 'BOUNDARY ABOVE'
354 #endif
355  !
356  DO ip=nb0+1, nb1
357  ixy = mapbou(ip)
358  cfl = cfll(ixy)
359  ixyc = ixy - inc * int( min( 0. , sign(1.1,cfl) ) )
360  fla(ixy) = cfl * q(ixyc)
361 #ifdef W3_T1
362  iy = mod( ixy , my )
363  ix = 1 + ixy/my
364  iy2 = mod( ixy+inc , my )
365  ix2 = 1 + (ixy+inc)/my
366  qn = max( q(ixy+inc), q(ixy) )
367  IF ( qn .GT. 1.e-10 ) THEN
368  qn = 1. /qn
369  WRITE (ndst,9013) ip, ix, iy, ix2, iy2, cfl, &
370  cfll(ixy), q(ixyc)*qn, q(ixy)*qn, q(ixy+inc)*qn
371  END IF
372 #endif
373  END DO
374  !
375  ! 4. Fluxes for points with boundary below ------------------------- *
376  ! ( 1st order without limiter )
377  !
378 #ifdef W3_T1
379  WRITE (ndst,9011) nb2-nb1, 'BOUNDARY BELOW'
380 #endif
381  !
382  DO ip=nb1+1, nb2
383  ixy = mapbou(ip)
384  cfl = cfll(ixy+inc)
385  ixyc = ixy - inc * int( min( 0. , sign(1.1,cfl) ) )
386  fla(ixy) = cfl * q(ixyc)
387 #ifdef W3_T1
388  iy = mod( ixy , my )
389  ix = 1 + ixy/my
390  iy2 = mod( ixy+inc , my )
391  ix2 = 1 + (ixy+inc)/my
392  qn = max( q(ixy+inc), q(ixy) )
393  IF ( qn .GT. 1.e-10 ) THEN
394  qn = 1. /qn
395  WRITE (ndst,9014) ip, ix, iy, ix2, iy2, cfl, &
396  cfll(ixy+inc), q(ixyc)*qn, q(ixy)*qn, q(ixy+inc)*qn
397  END IF
398 #endif
399  END DO
400  !
401  ! 5. Global closure ----------------------------------------------- *
402  !
403  IF ( CLOSE ) THEN
404 #ifdef W3_T
405  WRITE (ndst,9015)
406 #endif
407  DO iy=1, ny
408  fla(iy+iad00) = fla(iy+iadn0)
409  END DO
410  END IF
411  !
412  ! 6. Propagation -------------------------------------------------- *
413  !
414 #ifdef W3_T2
415  WRITE (ndst,9020)
416 #endif
417  DO ip=1, nact
418  ixy = mapact(ip)
419 #ifdef W3_T2
420  qold = q(ixy)
421 #endif
422  q(ixy) = max( 0. , q(ixy) + fla(ixy-inc) - fla(ixy) )
423 #ifdef W3_T2
424  IF ( qold + q(ixy) .GT. 1.e-10 ) &
425  WRITE (ndst,9021) ip, ixy, qold, q(ixy), &
426  fla(ixy-inc), fla(ixy)
427 #endif
428  END DO
429  !
430 #ifdef W3_T0
431  WRITE (ndst,9001) 'Q'
432  DO iy=ny,1,-1
433  WRITE (ndst,9002) (nint(q(iy+(ix-1)*my)/qmax),ix=1,nx)
434  END DO
435 #endif
436  !
437  RETURN
438  !
439  ! Formats
440  !
441 #ifdef W3_T
442 9000 FORMAT ( ' TEST W3QCK1 : ARRAY DIMENSIONS :',2i6/ &
443  ' USED :',2i6/ &
444  ' CLOSE, INC :',l6,i6/ &
445  ' NB0, NB1, NB2 :',3i6)
446 #endif
447 #ifdef W3_T0
448 9001 FORMAT ( ' TEST W3QCK1 : DUMP ARRAY ',a,' :')
449 9002 FORMAT ( 1x,43i3)
450 9003 FORMAT ( 1x,21i6)
451 #endif
452 #ifdef W3_T
453 9005 FORMAT (' TEST W3QCK1 : GLOBAL CLOSURE (1)')
454 #endif
455  !
456 #ifdef W3_T1
457 9010 FORMAT (' TEST W3QCK1 : IP, 2x(IX,IY), CFL (b,i,i+1), ', &
458  ' Q (b,b,i-1,i,i+1,i+2)')
459 9011 FORMAT (' TEST W3QCK1 :',i6,' POINTS OF TYPE ',a)
460 9012 FORMAT (10x,i6,4i4,1x,3f6.2,1x,f7.2,f6.2,1x,4f6.2)
461 9013 FORMAT (10x,i6,4i4,1x,f6.2,f6.2,' --- ',1x,f7.2,1x,' --- ',&
462  2f6.2,' --- ')
463 9014 FORMAT (10x,i6,4i4,1x,f6.2,' --- ',f6.2,1x,f7.2,1x,' --- ',&
464  2f6.2,' --- ')
465 #endif
466 #ifdef W3_T
467 9015 FORMAT (' TEST W3QCK1 : GLOBAL CLOSURE (2)')
468 #endif
469  !
470 #ifdef W3_T2
471 9020 FORMAT (' TEST W3QCK1 : IP, IXY, 2Q, 2FL')
472 9021 FORMAT (' ',2i6,2(1x,2e11.3))
473 #endif
474  !/
475  !/ End of W3QCK1 ----------------------------------------------------- /
476  !/

References w3servmd::strace().

Referenced by w3pro2md::w3ktp2(), and w3pro3md::w3ktp3().

◆ w3qck2()

subroutine w3uqckmd::w3qck2 ( integer, intent(in)  MX,
integer, intent(in)  MY,
integer, intent(in)  NX,
integer, intent(in)  NY,
real, dimension(my*(mx+1)), intent(inout)  VELO,
real, intent(in)  DT,
real, dimension(my*(mx+1)), intent(inout)  DX1,
real, dimension(1-my:my*(mx+1)), intent(inout)  DX2,
real, dimension(1-my:my*(mx+2)), intent(inout)  Q,
logical, intent(in)  CLOSE,
integer, intent(in)  INC,
integer, dimension(my*mx), intent(in)  MAPACT,
integer, intent(in)  NACT,
integer, dimension(my*mx), intent(in)  MAPBOU,
integer, intent(in)  NB0,
integer, intent(in)  NB1,
integer, intent(in)  NB2,
integer, intent(in)  NDSE,
integer, intent(in)  NDST 
)

Like W3QCK1 with variable grid spacing.

VELO amd Q need only bee filled in the (MY,MX) range, extension is used internally for closure. VELO and Q are defined as 1-D arrays internally.

Called by: W3KTP2 Propagation in spectral space.

This routine can be used independently from WAVEWATCH III.

Parameters
[in]MXField dimensions, if grid is 'closed' or circular, MX is the closed dimension.
[in]MYField dimension (See MX).
[in]NXPart of field actually used.
[in]NYPart of field actually used.
[in]MAPACTList of active grid points.
[in]NACTSize of MAPACT.
[in]MAPBOUMap with boundary information (See W3MAP2).
[in]NB0Counter in MAPBOU.
[in]NB1Counter in MAPBOU.
[in]NB2Counter in MAPBOU.
[in,out]VELOLocal velocities (MY, MX+1).
[in]DTTime step.
[in,out]DX1Band width at points (MY, MX+1).
[in,out]DX2Band width between points (MY,0:MX+1)
[in]NDSEError output unit number.
[in]NDSTTest output unit number.
[in,out]QPropagated quantity (MY,0:MX+2).
[in]CLOSEFlag for closed 'X' dimension.
[in]INCIncrement in 1-D array corresponding to increment in 2-D space.
Author
H. L. Tolman
Date
30-Oct-2009

Definition at line 517 of file w3uqckmd.F90.

517  !/
518  !/ +-----------------------------------+
519  !/ | WAVEWATCH III NOAA/NCEP |
520  !/ | H. L. Tolman |
521  !/ | FORTRAN 90 |
522  !/ | Last update : 30-Oct-2009 |
523  !/ +-----------------------------------+
524  !/
525  !/ 07-Sep-1997 : Final FORTRAN 77 ( version 1.18 )
526  !/ 16-Dec-1999 : Upgrade to FORTRAN 90 ( version 2.00 )
527  !/ 14-Feb-2001 : Unit numbers added to par list. ( version 2.08 )
528  !/ 05-Mar-2008 : Added NEC sxf90 compiler directives.
529  !/ (Chris Bunney, UK Met Office) ( version 3.13 )
530  !/ 30-Oct-2009 : Fixed "Called by" doc line. ( version 3.14 )
531  !/ (T. J. Campbell, NRL)
532  !/
533  ! 1. Purpose :
534  !
535  ! Like W3QCK1 with variable grid spacing.
536  !
537  ! 3. Parameters :
538  !
539  ! Parameter list
540  ! ----------------------------------------------------------------
541  ! MX,MY Int. I Field dimensions, if grid is 'closed' or
542  ! circular, MX is the closed dimension.
543  ! NX,NY Int. I Part of field actually used.
544  ! VELO R.A. I Local velocities. (MY, MX+1)
545  ! DT Real I Time step.
546  ! DX1 R.A. I/O Band width at points. (MY, MX+1)
547  ! DX2 R.A. I/O Band width between points. (MY,0:MX+1)
548  ! (local counter and counter+INC)
549  ! Q R.A. I/O Propagated quantity. (MY,0:MX+2)
550  ! CLOSE Log. I Flag for closed 'X' dimension'
551  ! INC Int. I Increment in 1-D array corresponding to
552  ! increment in 2-D space.
553  ! MAPACT I.A. I List of active grid points.
554  ! NACT Int. I Size of MAPACT.
555  ! MAPBOU I.A. I Map with boundary information (see W3MAP2).
556  ! NBn Int. I Counters in MAPBOU.
557  ! NDSE Int. I Error output unit number.
558  ! NDST Int. I Test output unit number.
559  ! ----------------------------------------------------------------
560  ! - VELO amd Q need only bee filled in the (MY,MX) range,
561  ! extension is used internally for closure.
562  ! - VELO and Q are defined as 1-D arrays internally.
563  !
564  ! 4. Subroutines used :
565  !
566  ! STRACE Service routine.
567  !
568  ! 5. Called by :
569  !
570  ! W3KTP2 Propagation in spectral space
571  !
572  ! 6. Error messages :
573  !
574  ! None.
575  !
576  ! 7. Remarks :
577  !
578  ! - This routine can be used independently from WAVEWATCH III.
579  !
580  ! 8. Structure :
581  !
582  ! ------------------------------------------------------
583  ! 1. Initialize aux. array FLA.
584  ! 2. Fluxes for central points (3rd order + limiter).
585  ! 3. Fluxes boundary point above (1st order).
586  ! 4. Fluxes boundary point below (1st order).
587  ! 5. Closure of 'X' if required
588  ! 6. Propagate.
589  ! ------------------------------------------------------
590  !
591  ! 9. Switches :
592  !
593  ! !/S Enable subroutine tracing.
594  ! !/T Enable test output.
595  ! !/T0 Test output input/output fields.
596  ! !/T1 Test output fluxes.
597  ! !/T2 Test output integration.
598  !
599  ! 10. Source code :
600  !
601  !/ ------------------------------------------------------------------- /
602  IMPLICIT NONE
603  !/
604  !/ ------------------------------------------------------------------- /
605  !/ Parameter list
606  !/
607  INTEGER, INTENT(IN) :: MX, MY, NX, NY, INC, MAPACT(MY*MX), &
608  NACT, MAPBOU(MY*MX), NB0, NB1, NB2, &
609  NDSE, NDST
610  REAL, INTENT(IN) :: DT
611  REAL, INTENT(INOUT) :: VELO(MY*(MX+1)), DX1(MY*(MX+1)), &
612  DX2(1-MY:MY*(MX+1)), Q(1-MY:MY*(MX+2))
613  LOGICAL, INTENT(IN) :: CLOSE
614  !/
615  !/ ------------------------------------------------------------------- /
616  !/ Local parameters
617  !/
618  INTEGER :: IXY, IP, IXYC, IXYU, IXYD, IY, IX, &
619  IAD00, IAD02, IADN0, IADN1, IADN2
620 #ifdef W3_S
621  INTEGER, SAVE :: IENT
622 #endif
623 #ifdef W3_T1
624  INTEGER :: IX2, IY2
625 #endif
626  REAL :: CFL, VEL, QB, DQ, DQNZ, QCN, QBN, &
627  QBR, CFAC, FLA(1-MY:MY*MX)
628 #ifdef W3_T0
629  REAL :: QMAX
630 #endif
631 #ifdef W3_T1
632  REAL :: QBO, QN, XCFL
633 #endif
634 #ifdef W3_T2
635  REAL :: QOLD
636 #endif
637  !/
638  !/ ------------------------------------------------------------------- /
639  !/
640 #ifdef W3_S
641  CALL strace (ient, 'W3QCK2')
642 #endif
643  !
644 #ifdef W3_T
645  WRITE (ndst,9000) mx, my, nx, ny, dt, CLOSE, inc, nb0, nb1, nb2
646 #endif
647  !
648 #ifdef W3_T0
649  qmax = 0.
650  DO iy=1, ny
651  DO ix=1, nx
652  qmax = max( qmax , q(iy+(ix-1)*my) )
653  END DO
654  END DO
655  qmax = max( 0.01*qmax , 1.e-10 )
656 #endif
657  !
658 #ifdef W3_T0
659  WRITE (ndst,9001) 'VELO'
660  DO iy=ny,1,-1
661  WRITE (ndst,9002) (nint(100.*velo(iy+(ix-1)*my) &
662  *dt/dx1(iy+(ix-1)*my)),ix=1,nx)
663  END DO
664  WRITE (ndst,9001) 'Q'
665  DO iy=ny,1,-1
666  WRITE (ndst,9002) (nint(q(iy+(ix-1)*my)/qmax),ix=1,nx)
667  END DO
668  WRITE (ndst,9001) 'MAPACT'
669  WRITE (ndst,9003) (mapact(ixy),ixy=1,nact)
670 #endif
671  !
672  ! 1. Initialize aux. array FLA and closure ------------------------- *
673  !
674  fla = 0.
675  !
676  IF ( CLOSE ) THEN
677 #ifdef W3_T
678  WRITE (ndst,9005)
679 #endif
680  iad00 = -my
681  iad02 = my
682  iadn0 = iad00 + my*nx
683  iadn1 = my*nx
684  iadn2 = iad02 + my*nx
685  DO iy=1, ny
686  q(iy+iad00) = q(iy+iadn0)
687  q(iy+iadn1) = q( iy )
688  q(iy+iadn2) = q(iy+iad02)
689  velo(iy+iadn1) = velo( iy )
690  dx1(iy+iadn1) = dx1( iy )
691  dx2(iy+iad00) = dx1(iy+iadn0)
692  dx2(iy+iadn1) = dx1( iy )
693  END DO
694  END IF
695  !
696  ! 2. Fluxes for central points ------------------------------------- *
697  ! ( 3rd order + limiter )
698  !
699 #ifdef W3_T1
700  WRITE (ndst,9010)
701  WRITE (ndst,9011) nb0, 'CENTRAL'
702 #endif
703  !
704  DO ip=1, nb0
705  !
706  ixy = mapbou(ip)
707  vel = 0.5 * ( velo(ixy) + velo(ixy+inc) )
708  cfl = dt * vel / dx2(ixy)
709  ixyc = ixy - inc * int( min( 0. , sign(1.1,cfl) ) )
710  qb = 0.5 * ( (1.-cfl)*q(ixy+inc) + (1.+cfl)*q(ixy) ) &
711  - dx2(ixy)**2 / dx1(ixyc) * (1.-cfl**2) / 6. &
712  * ( (q(ixyc+inc)-q(ixyc))/dx2(ixyc) &
713  - (q(ixyc)-q(ixyc-inc))/dx2(ixyc-inc) )
714  !
715  ixyu = ixyc - inc * int( sign(1.1,cfl) )
716  ixyd = 2*ixyc - ixyu
717  dq = q(ixyd) - q(ixyu)
718  dqnz = sign( max(1.e-15,abs(dq)) , dq )
719  qcn = ( q(ixyc) - q(ixyu) ) / dqnz
720  qcn = min( 1.1, max( -0.1 , qcn ) )
721  !
722 #ifdef W3_T1
723  qbo = qb
724 #endif
725  qbn = max( (qb-q(ixyu))/dqnz , qcn )
726  qbn = min( qbn , 1. , qcn/max(1.e-10,abs(cfl)) )
727  qbr = q(ixyu) + qbn*dq
728  cfac = real( int( 2. * abs(qcn-0.5) ) )
729  qb = (1.-cfac)*qbr + cfac*q(ixyc)
730  !
731  fla(ixy) = vel * qb
732  !
733 #ifdef W3_T1
734  iy = mod( ixy , my )
735  ix = 1 + ixy/my
736  iy2 = mod( ixy+inc , my )
737  ix2 = 1 + (ixy+inc)/my
738  qn = max( qb, qbo, q(ixy-inc), q( ixy ), &
739  q(ixy+inc), q(ixy+2*inc) )
740  IF ( qn .GT. 1.e-10 ) THEN
741  qn = 1. /qn
742  WRITE (ndst,9012) ip, ix, iy, ix2, iy2, &
743  cfl, dt*velo(ixy)/dx1(ixy), &
744  dt*velo(ixy+inc)/dx1(ixy+inc), &
745  qbo*qn, qb*qn, q(ixy-inc)*qn, q( ixy )*qn, &
746  q(ixy+inc)*qn, q(ixy+2*inc)*qn
747  END IF
748 #endif
749  !
750  END DO
751  !
752  ! 3. Fluxes for points with boundary above ------------------------- *
753  ! ( 1st order without limiter )
754  !
755 #ifdef W3_T1
756  WRITE (ndst,9011) nb1-nb0, 'BOUNDARY ABOVE'
757 #endif
758  !
759  DO ip=nb0+1, nb1
760  ixy = mapbou(ip)
761  vel = velo(ixy)
762  ixyc = ixy - inc * int( min( 0. , sign(1.1,vel) ) )
763  fla(ixy) = vel * q(ixyc)
764 #ifdef W3_T1
765  iy = mod( ixy , my )
766  ix = 1 + ixy/my
767  iy2 = mod( ixy+inc , my )
768  ix2 = 1 + (ixy+inc)/my
769  qn = max( q(ixy+inc), q(ixy) )
770  IF ( qn .GT. 1.e-10 ) THEN
771  qn = 1. /qn
772  WRITE (ndst,9013) ip, ix, iy, ix2, iy2, xcfl, &
773  dt*velo(ixy)/dx2(ixy), &
774  q(ixyc)*qn, q(ixy)*qn, q(ixy+inc)*qn
775  END IF
776 #endif
777  END DO
778  !
779  ! 4. Fluxes for points with boundary below ------------------------- *
780  ! ( 1st order without limiter )
781  !
782 #ifdef W3_T1
783  WRITE (ndst,9011) nb2-nb1, 'BOUNDARY BELOW'
784 #endif
785  !
786  DO ip=nb1+1, nb2
787  ixy = mapbou(ip)
788  vel = velo(ixy+inc)
789  ixyc = ixy - inc * int( min( 0. , sign(1.1,vel) ) )
790  fla(ixy) = vel * q(ixyc)
791 #ifdef W3_T1
792  iy = mod( ixy , my )
793  ix = 1 + ixy/my
794  iy2 = mod( ixy+inc , my )
795  ix2 = 1 + (ixy+inc)/my
796  qn = max( q(ixy+inc), q(ixy) )
797  IF ( qn .GT. 1.e-10 ) THEN
798  qn = 1. /qn
799  WRITE (ndst,9014) ip, ix, iy, ix2, iy2, xcfl, &
800  dt*velo(ixy+inc)/dx2(ixy), &
801  q(ixyc)*qn, q(ixy)*qn, q(ixy+inc)*qn
802  END IF
803 #endif
804  END DO
805  !
806  ! 5. Global closure ----------------------------------------------- *
807  !
808  IF ( CLOSE ) THEN
809 #ifdef W3_T
810  WRITE (ndst,9015)
811 #endif
812  DO iy=1, ny
813  fla(iy+iad00) = fla(iy+iadn0)
814  END DO
815  END IF
816  !
817  ! 6. Propagation -------------------------------------------------- *
818  !
819 #ifdef W3_T2
820  WRITE (ndst,9020)
821 #endif
822  DO ip=1, nact
823  ixy = mapact(ip)
824 #ifdef W3_T2
825  qold = q(ixy)
826 #endif
827  q(ixy) = max( 0. , q(ixy) + dt/dx1(ixy) * &
828  (fla(ixy-inc)-fla(ixy)) )
829 #ifdef W3_T2
830  IF ( qold + q(ixy) .GT. 1.e-10 ) &
831  WRITE (ndst,9021) ip, ixy, qold, q(ixy), &
832  dt*fla(ixy-inc)/dx1(ixy), &
833  dt*fla(ixy)/dx1(ixy)
834 #endif
835  END DO
836  !
837 #ifdef W3_T0
838  WRITE (ndst,9001) 'Q'
839  DO iy=ny,1,-1
840  WRITE (ndst,9002) (nint(q(iy+(ix-1)*my)/qmax),ix=1,nx)
841  END DO
842 #endif
843  !
844  RETURN
845  !
846  ! Formats
847  !
848 #ifdef W3_T
849 9000 FORMAT ( ' TEST W3QCK2 : ARRAY DIMENSIONS :',2i6/ &
850  ' USED :',2i6/ &
851  ' TIME STEP :',f8.1/ &
852  ' CLOSE, INC :',l6,i6/ &
853  ' NB0, NB1, NB2 :',3i6)
854 #endif
855 #ifdef W3_T0
856 9001 FORMAT ( ' TEST W3QCK2 : DUMP ARRAY ',a,' :')
857 9002 FORMAT ( 1x,43i3)
858 9003 FORMAT ( 1x,21i6)
859 #endif
860 #ifdef W3_T
861 9005 FORMAT (' TEST W3QCK2 : GLOBAL CLOSURE (1)')
862 #endif
863  !
864 #ifdef W3_T1
865 9010 FORMAT (' TEST W3QCK2 : IP, 2x(IX,IY), CFL (b,i,i+1), ', &
866  ' Q (b,b,i-1,i,i+1,i+2)')
867 9011 FORMAT (' TEST W3QCK2 :',i6,' POINTS OF TYPE ',a)
868 9012 FORMAT (10x,i6,4i4,1x,3f6.2,1x,f7.2,f6.2,1x,4f6.2)
869 9013 FORMAT (10x,i6,4i4,1x,f6.2,f6.2,' --- ',1x,f7.2,1x,' --- ',&
870  2f6.2,' --- ')
871 9014 FORMAT (10x,i6,4i4,1x,f6.2,' --- ',f6.2,1x,f7.2,1x,' --- ',&
872  2f6.2,' --- ')
873 #endif
874 #ifdef W3_T
875 9015 FORMAT (' TEST W3QCK2 : GLOBAL CLOSURE (2)')
876 #endif
877  !
878 #ifdef W3_T2
879 9020 FORMAT (' TEST W3QCK2 : IP, IXY, 2Q, 2FL')
880 9021 FORMAT (' ',2i6,2(1x,2e11.3))
881 #endif
882  !/
883  !/ End of W3QCK2 ----------------------------------------------------- /
884  !/

References w3servmd::strace().

Referenced by w3pro2md::w3ktp2(), and w3pro3md::w3ktp3().

◆ w3qck3()

subroutine w3uqckmd::w3qck3 ( integer, intent(in)  MX,
integer, intent(in)  MY,
integer, intent(in)  NX,
integer, intent(in)  NY,
real, dimension(my*(mx+1)), intent(inout)  CFLL,
real, dimension(my*mx,-1:1), intent(in)  TRANS,
real, dimension(1-my:my*(mx+2)), intent(inout)  Q,
logical, intent(in)  CLOSE,
integer, intent(in)  INC,
integer, dimension(my*mx), intent(in)  MAPACT,
integer, intent(in)  NACT,
integer, dimension(my*mx), intent(in)  MAPBOU,
integer, intent(in)  NB0,
integer, intent(in)  NB1,
integer, intent(in)  NB2,
integer, intent(in)  NDSE,
integer, intent(in)  NDST 
)

Like W3QCK1 with cell transparencies added.

CFLL amd Q need only bee filled in the (MY,MX) range, extension is used internally for closure. CFLL and Q are defined as 1-D arrays internally.

Called by: W3XYP2 Propagation in physical space.

This routine can be used independently from WAVEWATCH III.

Parameters
[in]MXField dimensions, if grid is 'closed' or circular, MX is the closed dimension.
[in]MYField dimension (See MX)
[in]NXPart of field actually used
[in]NYPart of field actually used
[in,out]CFLLLocal Courant numbers (MY, MX+1).
[in]TRANS
[in]INCIncrement in 1-D array corresponding to increment in 2-D space.
[in]MAPACTList of active grid points.
[in]NACTSize of MAPACT.
[in]MAPBOUMap with boundary information (See W3MAP2).
[in]NB0Counter in MAPBOU
[in]NB1Counter in MAPBOU
[in]NB2Counter in MAPBOU
[in]NDSEError output unit number.
[in]NDSTTest output unit number.
[in,out]QPropagated quantity (MY,0:MX+2)
[in]CLOSEFlag for closed 'X' dimension.
Author
H. L. Tolman
Date
30-Oct-2009

Definition at line 922 of file w3uqckmd.F90.

922  !/
923  !/ +-----------------------------------+
924  !/ | WAVEWATCH III NOAA/NCEP |
925  !/ | H. L. Tolman |
926  !/ | FORTRAN 90 |
927  !/ | Last update : 27-May-2014 |
928  !/ +-----------------------------------+
929  !/
930  !/ 13_nov-2001 : Origination. ( version 2.14 )
931  !/ 16-Oct-2002 : Fix INTENT for TRANS. ( version 3.00 )
932  !/ 05-Mar-2008 : Added NEC sxf90 compiler directives.
933  !/ (Chris Bunney, UK Met Office) ( version 3.13 )
934  !/ 27-May-2014 : Added OMPH switches in W3QCK3. ( version 5.02 )
935  !/
936  ! 1. Purpose :
937  !
938  ! Like W3QCK1 with cell transparencies added.
939  !
940  ! 2. Method :
941  !
942  ! 3. Parameters :
943  !
944  ! Parameter list
945  ! ----------------------------------------------------------------
946  ! MX,MY Int. I Field dimensions, if grid is 'closed' or
947  ! circular, MX is the closed dimension.
948  ! NX,NY Int. I Part of field actually used.
949  ! CFLL R.A. I Local Courant numbers. (MY, MX+1)
950  ! Q R.A. I/O Propagated quantity. (MY,0:MX+2)
951  ! CLOSE Log. I Flag for closed 'X' dimension'
952  ! INC Int. I Increment in 1-D array corresponding to
953  ! increment in 2-D space.
954  ! MAPACT I.A. I List of active grid points.
955  ! NACT Int. I Size of MAPACT.
956  ! MAPBOU I.A. I Map with boundary information (see W3MAP2).
957  ! NBn Int. I Counters in MAPBOU.
958  ! NDSE Int. I Error output unit number.
959  ! NDST Int. I Test output unit number.
960  ! ----------------------------------------------------------------
961  ! - CFLL amd Q need only bee filled in the (MY,MX) range,
962  ! extension is used internally for closure.
963  ! - CFLL and Q are defined as 1-D arrays internally.
964  !
965  ! 4. Subroutines used :
966  !
967  ! STRACE Service routine.
968  !
969  ! 5. Called by :
970  !
971  ! W3XYP2 Propagation in physical space
972  !
973  ! 6. Error messages :
974  !
975  ! None.
976  !
977  ! 7. Remarks :
978  !
979  ! - This routine can be used independently from WAVEWATCH III.
980  !
981  ! 8. Structure :
982  !
983  ! ------------------------------------------------------
984  ! 1. Initialize aux. array FLA.
985  ! 2. Fluxes for central points (3rd order + limiter).
986  ! 3. Fluxes boundary point above (1st order).
987  ! 4. Fluxes boundary point below (1st order).
988  ! 5. Closure of 'X' if required
989  ! 6. Propagate.
990  ! ------------------------------------------------------
991  !
992  ! 9. Switches :
993  !
994  ! !/OMPH Ading OMP directves for hybrid paralellization.
995  !
996  ! !/S Enable subroutine tracing.
997  ! !/T Enable test output.
998  ! !/T0 Test output input/output fields.
999  ! !/T1 Test output fluxes.
1000  ! !/T2 Test output integration.
1001  !
1002  ! 10. Source code :
1003  !
1004  !/ ------------------------------------------------------------------- /
1005  IMPLICIT NONE
1006  !/
1007  !/ ------------------------------------------------------------------- /
1008  !/ Parameter list
1009  !/
1010  INTEGER, INTENT(IN) :: MX, MY, NX, NY, INC, MAPACT(MY*MX), &
1011  NACT, MAPBOU(MY*MX), NB0, NB1, NB2, &
1012  NDSE, NDST
1013  REAL, INTENT(IN) :: TRANS(MY*MX,-1:1)
1014  REAL, INTENT(INOUT) :: CFLL(MY*(MX+1)), Q(1-MY:MY*(MX+2))
1015  LOGICAL, INTENT(IN) :: CLOSE
1016  !/
1017  !/ ------------------------------------------------------------------- /
1018  !/ Local parameters
1019  !/
1020  INTEGER :: IXY, IP, IXYC, IXYU, IXYD, IY, IX, &
1021  IAD00, IAD02, IADN0, IADN1, IADN2, &
1022  JN, JP
1023 #ifdef W3_S
1024  INTEGER, SAVE :: IENT = 0
1025 #endif
1026 #ifdef W3_T1
1027  INTEGER :: IX2, IY2
1028 #endif
1029  REAL :: CFL, QB, DQ, DQNZ, QCN, QBN, QBR, CFAC
1030  REAL :: FLA(1-MY:MY*MX)
1031 #ifdef W3_T0
1032  REAL :: QMAX
1033 #endif
1034 #ifdef W3_T1
1035  REAL :: QBO, QN
1036 #endif
1037 #ifdef W3_T2
1038  REAL :: QOLD
1039 #endif
1040  !/
1041  !/ ------------------------------------------------------------------- /
1042  !/
1043 #ifdef W3_S
1044  CALL strace (ient, 'W3QCK3')
1045 #endif
1046  !
1047 #ifdef W3_T
1048  WRITE (ndst,9000) mx, my, nx, ny, CLOSE, inc, nb0, nb1, nb2
1049 #endif
1050  !
1051 #ifdef W3_T0
1052  qmax = 0.
1053  DO iy=1, ny
1054  DO ix=1, nx
1055  qmax = max( qmax , q(iy+(ix-1)*my) )
1056  END DO
1057  END DO
1058  qmax = max( 0.01*qmax , 1.e-10 )
1059 #endif
1060  !
1061 #ifdef W3_T0
1062  WRITE (ndst,9001) 'CFLL'
1063  DO iy=ny,1,-1
1064  WRITE (ndst,9002) (nint(100.*cfll(iy+(ix-1)*my)),ix=1,nx)
1065  END DO
1066  WRITE (ndst,9001) 'Q'
1067  DO iy=ny,1,-1
1068  WRITE (ndst,9002) (nint(q(iy+(ix-1)*my)/qmax),ix=1,nx)
1069  END DO
1070  WRITE (ndst,9001) 'MAPACT'
1071  WRITE (ndst,9003) (mapact(ixy),ixy=1,nact)
1072 #endif
1073  !
1074  ! 1. Initialize aux. array FLA and closure ------------------------- *
1075  !
1076  fla = 0.
1077  !
1078  IF ( CLOSE ) THEN
1079 #ifdef W3_T
1080  WRITE (ndst,9005)
1081 #endif
1082  iad00 = -my
1083  iad02 = my
1084  iadn0 = iad00 + my*nx
1085  iadn1 = my*nx
1086  iadn2 = iad02 + my*nx
1087  !
1088 #ifdef W3_OMPH
1089  !$OMP PARALLEL DO PRIVATE (IY)
1090 #endif
1091  !
1092  DO iy=1, ny
1093  q(iy+iad00) = q(iy+iadn0) ! 1 ghost column to left
1094  q(iy+iadn1) = q( iy ) ! 1st ghost column to right
1095  q(iy+iadn2) = q(iy+iad02) ! 2nd ghost column to right
1096  cfll(iy+iadn1) = cfll( iy ) ! as for Q above, 1st to rt
1097  END DO
1098  !
1099 #ifdef W3_OMPH
1100  !$OMP END PARALLEL DO
1101 #endif
1102  !
1103  END IF
1104  !
1105  ! 2. Fluxes for central points ------------------------------------- *
1106  ! ( 3rd order + limiter )
1107  !
1108 #ifdef W3_T1
1109  WRITE (ndst,9010)
1110  WRITE (ndst,9011) nb0, 'CENTRAL'
1111 #endif
1112  !
1113 #ifdef W3_OMPH
1114  !$OMP PARALLEL DO PRIVATE (IP, IXY, CFL, IXYC, QB, IXYU, IXYD, &
1115 #ifdef W3_T1
1116  !$OMP QBO, QN, IX, IY, IX2, IY2, &
1117 #endif
1118  !$OMP& DQ, DQNZ, QCN, QBN, QBR, CFAC )
1119 #endif
1120  !
1121  DO ip=1, nb0
1122  !
1123  ixy = mapbou(ip)
1124  cfl = 0.5 * ( cfll(ixy) + cfll(ixy+inc) )
1125  ixyc = ixy - inc * int( min( 0. , sign(1.1,cfl) ) )
1126  qb = 0.5 * ( (1.-cfl)*q(ixy+inc) + (1.+cfl)*q(ixy) ) &
1127  - (1.-cfl**2)/6. * (q(ixyc-inc)-2.*q(ixyc)+q(ixyc+inc))
1128  !
1129  ixyu = ixyc - inc * int( sign(1.1,cfl) )
1130  ixyd = 2*ixyc - ixyu
1131  dq = q(ixyd) - q(ixyu)
1132  dqnz = sign( max(1.e-15,abs(dq)) , dq )
1133  qcn = ( q(ixyc) - q(ixyu) ) / dqnz
1134  qcn = min( 1.1, max( -0.1 , qcn ) )
1135  !
1136 #ifdef W3_T1
1137  qbo = qb
1138 #endif
1139  qbn = max( (qb-q(ixyu))/dqnz , qcn )
1140  qbn = min( qbn , 1. , qcn/max(1.e-10,abs(cfl)) )
1141  qbr = q(ixyu) + qbn*dq
1142  cfac = real( int( 2. * abs(qcn-0.5) ) )
1143  qb = (1.-cfac)*qbr + cfac*q(ixyc)
1144  !
1145  fla(ixy) = cfl * qb
1146  !
1147 #ifdef W3_T1
1148  iy = mod( ixy , my )
1149  ix = 1 + ixy/my
1150  iy2 = mod( ixy+inc , my )
1151  ix2 = 1 + (ixy+inc)/my
1152  qn = max( qb, qbo, q(ixy-inc), q( ixy ), &
1153  q(ixy+inc), q(ixy+2*inc) )
1154  IF ( qn .GT. 1.e-10 ) THEN
1155  qn = 1. /qn
1156  WRITE (ndst,9012) ip, ix, iy, ix2, iy2, &
1157  cfl, cfll(ixy), cfll(ixy+inc), &
1158  qbo*qn, qb*qn, q(ixy-inc)*qn, q( ixy )*qn, &
1159  q(ixy+inc)*qn, q(ixy+2*inc)*qn
1160  END IF
1161 #endif
1162  !
1163  END DO
1164  !
1165 #ifdef W3_OMPH
1166  !$OMP END PARALLEL DO
1167 #endif
1168  !
1169  ! 3. Fluxes for points with boundary above ------------------------- *
1170  ! ( 1st order without limiter )
1171  !
1172 #ifdef W3_T1
1173  WRITE (ndst,9011) nb1-nb0, 'BOUNDARY ABOVE'
1174 #endif
1175  !
1176 !!!/OMPH/!$OMP PARALLEL DO PRIVATE (IP, IXY, CFL, IXYC)
1177 !!!
1178  DO ip=nb0+1, nb1
1179  ixy = mapbou(ip)
1180  cfl = cfll(ixy)
1181  ixyc = ixy - inc * int( min( 0. , sign(1.1,cfl) ) )
1182  fla(ixy) = cfl * q(ixyc)
1183 #ifdef W3_T1
1184  iy = mod( ixy , my )
1185  ix = 1 + ixy/my
1186  iy2 = mod( ixy+inc , my )
1187  ix2 = 1 + (ixy+inc)/my
1188  qn = max( q(ixy+inc), q(ixy) )
1189  IF ( qn .GT. 1.e-10 ) THEN
1190  qn = 1. /qn
1191  WRITE (ndst,9013) ip, ix, iy, ix2, iy2, cfl, &
1192  cfll(ixy), q(ixyc)*qn, q(ixy)*qn, q(ixy+inc)*qn
1193  END IF
1194 #endif
1195  END DO
1196 !!!
1197 !!!/OMPH/!$OMP END PARALLEL DO
1198  !
1199  ! 4. Fluxes for points with boundary below ------------------------- *
1200  ! ( 1st order without limiter )
1201  !
1202 #ifdef W3_T1
1203  WRITE (ndst,9011) nb2-nb1, 'BOUNDARY BELOW'
1204 #endif
1205  !
1206 !!!/OMPH/!$OMP PARALLEL DO PRIVATE (IP, IXY, CFL, IXYC)
1207 !!!
1208  DO ip=nb1+1, nb2
1209  ixy = mapbou(ip)
1210  cfl = cfll(ixy+inc)
1211  ixyc = ixy - inc * int( min( 0. , sign(1.1,cfl) ) )
1212  fla(ixy) = cfl * q(ixyc)
1213 #ifdef W3_T1
1214  iy = mod( ixy , my )
1215  ix = 1 + ixy/my
1216  iy2 = mod( ixy+inc , my )
1217  ix2 = 1 + (ixy+inc)/my
1218  qn = max( q(ixy+inc), q(ixy) )
1219  IF ( qn .GT. 1.e-10 ) THEN
1220  qn = 1. /qn
1221  WRITE (ndst,9014) ip, ix, iy, ix2, iy2, cfl, cfll(ixy+inc), &
1222  q(ixyc)*qn, q(ixy)*qn, q(ixy+inc)*qn
1223  END IF
1224 #endif
1225  END DO
1226  !
1227 !!!/OMPH/!$OMP END PARALLEL DO
1228  !
1229  ! 5. Global closure ----------------------------------------------- *
1230  !
1231  IF ( CLOSE ) THEN
1232 #ifdef W3_T
1233  WRITE (ndst,9015)
1234 #endif
1235  DO iy=1, ny
1236  fla(iy+iad00) = fla(iy+iadn0)
1237  END DO
1238  END IF
1239  !
1240  ! 6. Propagation -------------------------------------------------- *
1241  !
1242 #ifdef W3_T2
1243  WRITE (ndst,9020)
1244 #endif
1245 #ifdef W3_OMPH
1246  !$OMP PARALLEL DO PRIVATE (IP, IXY, JN, JP )
1247 #endif
1248  !
1249  DO ip=1, nact
1250  !
1251  ixy = mapact(ip)
1252  IF ( fla(ixy-inc) .GT. 0. ) THEN
1253  jn = -1
1254  ELSE
1255  jn = 0
1256  END IF
1257  IF ( fla(ixy ) .LT. 0. ) THEN
1258  jp = 1
1259  ELSE
1260  jp = 0
1261  END IF
1262  !
1263 #ifdef W3_T2
1264  qold = q(ixy)
1265 #endif
1266  q(ixy) = max( 0. , q(ixy) + trans(ixy,jn) * fla(ixy-inc) &
1267  - trans(ixy,jp) * fla(ixy) )
1268 
1269 #ifdef W3_T2
1270  IF ( qold + q(ixy) .GT. 1.e-10 ) &
1271  WRITE (ndst,9021) ip, ixy, qold, q(ixy), &
1272  fla(ixy-inc), fla(ixy)
1273 #endif
1274  END DO
1275  !
1276 #ifdef W3_OMPH
1277  !$OMP END PARALLEL DO
1278 #endif
1279  !
1280 #ifdef W3_T0
1281  WRITE (ndst,9001) 'Q'
1282  DO iy=ny,1,-1
1283  WRITE (ndst,9002) (nint(q(iy+(ix-1)*my)/qmax),ix=1,nx)
1284  END DO
1285 #endif
1286  !
1287  RETURN
1288  !
1289  ! Formats
1290  !
1291 #ifdef W3_T
1292 9000 FORMAT ( ' TEST W3QCK3 : ARRAY DIMENSIONS :',2i6/ &
1293  ' USED :',2i6/ &
1294  ' CLOSE, INC :',l6,i6/ &
1295  ' NB0, NB1, NB2 :',3i6)
1296 #endif
1297 #ifdef W3_T0
1298 9001 FORMAT ( ' TEST W3QCK3 : DUMP ARRAY ',a,' :')
1299 9002 FORMAT ( 1x,43i3)
1300 9003 FORMAT ( 1x,21i6)
1301 #endif
1302 #ifdef W3_T
1303 9005 FORMAT (' TEST W3QCK3 : GLOBAL CLOSURE (1)')
1304 #endif
1305  !
1306 #ifdef W3_T1
1307 9010 FORMAT (' TEST W3QCK3 : IP, 2x(IX,IY), CFL (b,i,i+1), ', &
1308  ' Q (b,b,i-1,i,i+1,i+2)')
1309 9011 FORMAT (' TEST W3QCK3 :',i6,' POINTS OF TYPE ',a)
1310 9012 FORMAT (10x,i6,4i4,1x,3f6.2,1x,f7.2,f6.2,1x,4f6.2)
1311 9013 FORMAT (10x,i6,4i4,1x,f6.2,f6.2,' --- ',1x,f7.2,1x,' --- ',&
1312  2f6.2,' --- ')
1313 9014 FORMAT (10x,i6,4i4,1x,f6.2,' --- ',f6.2,1x,f7.2,1x,' --- ',&
1314  2f6.2,' --- ')
1315 #endif
1316 #ifdef W3_T
1317 9015 FORMAT (' TEST W3QCK3 : GLOBAL CLOSURE (2)')
1318 #endif
1319  !
1320 #ifdef W3_T2
1321 9020 FORMAT (' TEST W3QCK3 : IP, IXY, 2Q, 2FL')
1322 9021 FORMAT (' ',2i6,2(1x,2e11.3))
1323 #endif
1324  !/
1325  !/ End of W3QCK3 ----------------------------------------------------- /
1326  !/

References w3servmd::strace().

Referenced by w3pro2md::w3xyp2(), and w3pro3md::w3xyp3().