WAVEWATCH III  beta 0.0.1
ww3_prtide.F90
Go to the documentation of this file.
1 
5 !
6 #include "w3macros.h"
7 !/ ------------------------------------------------------------------- /
12 PROGRAM w3prtide
13  !/
14  !/ +-----------------------------------+
15  !/ | WAVEWATCH III NOAA/NCEP |
16  !/ | F. Ardhuin |
17  !/ | FORTRAN 90 |
18  !/ | Last update : 21-Apr-2020 |
19  !/ +-----------------------------------+
20  !/
21  !/ 29-Mar-2013 : Creation ( version 4.11 )
22  !/ 17-Oct-2013 : Manages missing data for UNST grids ( version 4.12 )
23  !/ 06-Jun-2018 : COMPUTE VNEIGH: calculate the number of connected
24  !/ triangles for a given point ( version 6.04 )
25  !/ 21-Apr-2020 : MPI implementation ( version 7.13 )
26  !/ 21-Apr-2020 : bug fix for rectilinear grid ( version 7.13 )
27  !/ 1-Feb-2020 : Improve indexing, A.Roland ( version 7.14 )
28 
29  !/
30  !/ Copyright 2013 National Weather Service (NWS),
31  !/ National Oceanic and Atmospheric Administration. All rights
32  !/ reserved. WAVEWATCH III is a trademark of the NWS.
33  !/ No unauthorized use without permission.
34  !/
35  ! 1. Purpose :
36  !
37  ! Predicts tides (current or water level) to be used during
38  ! run by ww3_shel or ww3_multi (this takes much less memory).
39  !
40  ! 2. Method :
41  !
42  ! See documented input file.
43  !
44  ! 3. Parameters :
45  !
46  ! Local parameters.
47  ! ----------------------------------------------------------------
48  ! NDSI Int. Input unit number ("ww3_prtide.inp").
49  ! NDSLL Int. Unit number(s) of long-lat file(s)
50  ! NDSF I.A. Unit number(s) of input file(s).
51  ! NDSDAT Int. Unit number for output data file.
52  ! FLTIME Log. Time flag for input fields, if false, single
53  ! field, time read from NDSI.
54  ! IDLALL Int. Layout indicator used by INA2R. +
55  ! IDFMLL Int. Id. FORMAT indicator. |
56  ! FORMLL C*16 Id. FORMAT. | Long-lat
57  ! FROMLL C*4 'UNIT' / 'NAME' indicator | file(s)
58  ! NAMELL C*40 Name of long-lat file(s) +
59  ! IDLAF I.A. +
60  ! IDFMF I.A. |
61  ! FORMF C.A. | Idem. fields file(s)
62  ! FROMF C*4 |
63  ! NAMEF C*50 +
64  ! FORMT C.A. Format or time in field.
65  ! XC R.A. Components of input vector field or first
66  ! input scalar field
67  ! YC R.A. Components of input vector field or second
68  ! input scalar field
69  ! ----------------------------------------------------------------
70  !
71  ! 4. Subroutines used :
72  !
73  ! Name Type Module Description
74  ! ----------------------------------------------------------------
75  ! W3NMOD Subr. W3GDATMD Set number of model.
76  ! W3SETG Subr. Id. Point to selected model.
77  ! W3NOUT Subr. W3ODATMD Set number of model for output.
78  ! W3SETO Subr. Id. Point to selected model for output.
79  ! ITRACE Subr. W3SERVMD Subroutine tracing initialization.
80  ! STRACE Subr. Id. Subroutine tracing.
81  ! NEXTLN Subr. Id. Get next line from input filw
82  ! EXTCDE Subr. Id. Abort program as graceful as possible.
83  ! STME21 Subr. W3TIMEMD Convert time to string.
84  ! W3IOGR Subr. W3IOGRMD Reading/writing model definition file.
85  ! W3FLDO Subr. W3FLDSMD Opening of WAVEWATCH III generic shell
86  ! data file.
87  ! W3FLDG Subr. Id. Reading/writing shell input data.
88  ! ----------------------------------------------------------------
89  !
90  ! 5. Called by :
91  !
92  ! None, stand-alone program.
93  !
94  ! 6. Error messages :
95  !
96  ! - Checks on files and reading from file.
97  ! - Checks on validity of input parameters.
98  !
99  ! 7. Remarks :
100  !
101  ! - Input fields need to be continuous in longitude and latitude.
102  !
103  ! 8. Structure :
104  !
105  ! ----------------------------------------------------
106  ! To be updated ...
107  ! ----------------------------------------------------
108  !
109  ! 9. Switches :
110  !
111  ! !/S Enable subroutine tracing.
112  ! !/T Enable test output,
113  !
114  ! !/NCO NCEP NCO modifications for operational implementation.
115  !
116  ! 10. Source code :
117  !
118  !/ ------------------------------------------------------------------- /
119  USE constants
120  !/
121  ! USE W3GDATMD, ONLY: W3NMOD, W3SETG
122 #ifdef W3_NL1
123  USE w3adatmd,ONLY: w3naux, w3seta
124 #endif
125  USE w3odatmd, ONLY: iaproc, naproc, naperr, napout
126  USE w3odatmd, ONLY: w3nout, w3seto
127  USE w3servmd, ONLY : itrace, nextln, extcde, strsplit
128 #ifdef W3_S
129  USE w3servmd, ONLY : strace
130 #endif
131  USE w3timemd
132  USE w3arrymd, ONLY : ina2r, ina2i
133  USE w3iogrmd, ONLY: w3iogr
135  !/
136  USE w3gdatmd
137  USE w3gsrumd
138  USE w3odatmd, ONLY: ndse, ndst, ndso, fnmpre
139  USE w3tidemd
140  USE w3idatmd
141  !
142  IMPLICIT NONE
143  !
144 #ifdef W3_MPI
145  include "mpif.h"
146 #endif
147  !/
148  !/ ------------------------------------------------------------------- /
149  !/ Local parameters
150  !/
151  INTEGER :: ndsi, ndsf, ndsm, ndsdat, ndstrc, ntrace
152  INTEGER :: ierr, ifld, i, jj, j, ix, iy
153  INTEGER :: dttst, ndsen, prtide_dt
154  INTEGER :: tide_prmf, flagtide, tindex
155  INTEGER :: tideok, tide_max, tide_maxi
156  INTEGER :: k, icon, ix2, sumok, nbad, iter
157  INTEGER :: ie, ip, ip2, ii, ifound, alreadyfound
158  INTEGER :: tide_kd0, int24, intdys ! "Gregorian day constant"
159 #ifdef W3_MPI
160  INTEGER :: ierr_mpi, ind, rest, slice
161 #endif
162  INTEGER :: time(2), tide_start(2), tide_end(2)
163  INTEGER :: indmax(70), pr_inds(70)
164  !
165  INTEGER, ALLOCATABLE :: badpoints(:,:), vneigh(:,:), conn(:)
166 #ifdef W3_MPI
167  INTEGER, ALLOCATABLE :: nelem(:), cumul(:)
168 #endif
169  !
170  REAL :: wcurtidex, wcurtidey, tide_argx, tide_argy
171  REAL :: ampcos, ampsin
172  !
173  REAL :: tide_fx(44),ux(44),vx(44), maxvalcon(70)
174  !
175  REAL, ALLOCATABLE :: fx(:,:), fy(:,:), fa(:,:)
176 #ifdef W3_MPI
177  REAL, ALLOCATABLE :: fx1d(:), fy1d(:), fa1d(:)
178  REAL, ALLOCATABLE :: fx1dl(:), fy1dl(:), fa1dl(:)
179 #endif
180  !
181  DOUBLE PRECISION :: d1,h,tide_hour,hh,pp,s,p,enp,dh,dpp,ds,dp,dnp,tau
182  !
183  CHARACTER*256 :: filenamext
184  CHARACTER :: tideconstnames*1024
185  CHARACTER*23 :: idtime
186  CHARACTER :: comstr*1, idfld*3
187  !
188  CHARACTER(LEN=100) :: tidecon_prnames(70), tidecon_maxnames(70)
189  CHARACTER(LEN=100) :: tidecon_maxvals(70)
190  !
191  LOGICAL :: tidefill
192  !
193  !/
194  !/ ------------------------------------------------------------------- /
195  !/
196 
197  !==========================================================
198  !
199  ! Initialization
200  !
201  !==========================================================
202 
203  !
204  ! 1.a Set number of models
205  !
206  CALL w3nmod ( 1, 6, 6 )
207  CALL w3setg ( 1, 6, 6 )
208 #ifdef W3_NL1
209  CALL w3naux ( 6, 6 )
210  CALL w3seta ( 1, 6, 6 )
211 #endif
212  CALL w3nout ( 6, 6 )
213  CALL w3seto ( 1, 6, 6 )
214  !
215  ! 1.b IO set-up.
216  !
217  ndsi = 10
218  ndso = 6
219  ndse = 6
220  ndst = 6
221  ndsm = 11
222  ndsdat = 12
223  ndsf = 13
224  !
225  ndstrc = 6
226  ntrace = 10
227  tidefill =.true.
228  CALL itrace ( ndstrc, ntrace )
229  !
230 #ifdef W3_NCO
231  !
232  ! Redo according to NCO
233  !
234  ndsi = 11
235  ndso = 6
236  ndse = ndso
237  ndst = ndso
238  ndsm = 12
239  ndsdat = 51
240  ndstrc = ndso
241 #endif
242 
243 #ifdef W3_S
244  CALL strace (ient, 'W3PRTIDE')
245 #endif
246 
247  !
248  ! 1.c MPP initializations
249  !
250 #ifdef W3_SHRD
251  naproc = 1
252  iaproc = 1
253 #endif
254  !
255 #ifdef W3_MPI
256  CALL mpi_init ( ierr_mpi )
257  CALL mpi_comm_size ( mpi_comm_world, naproc, ierr_mpi )
258  CALL mpi_comm_rank ( mpi_comm_world, iaproc, ierr_mpi )
259  iaproc = iaproc + 1 ! this is to have IAPROC between 1 and NAPROC
260 #endif
261  !
262  IF ( iaproc .EQ. naperr ) THEN
263  ndsen = ndse
264  ELSE
265  ndsen = -1
266  END IF
267  !
268  IF ( iaproc .EQ. napout ) WRITE (ndso,900)
269  !
270  OPEN (ndsi,file=trim(fnmpre)//'ww3_prtide.inp',status='OLD', &
271  err=800,iostat=ierr)
272  rewind(ndsi)
273  READ (ndsi,'(A)',END=801,ERR=802,IOSTAT=IERR) comstr
274  IF (comstr.EQ.' ') comstr = '$'
275  IF ( iaproc .EQ. napout ) WRITE (ndso,901) comstr
276 
277  !==========================================================
278  !
279  ! Read model definition file.
280  !
281  !==========================================================
282 
283  CALL w3iogr ( 'READ', ndsm )
284  IF ( iaproc .EQ. napout ) WRITE (ndso,902) gname
285  ALLOCATE ( fx(nx,ny), fy(nx,ny), fa(nx,ny), badpoints(nx,ny) )
286 
287  !==========================================================
288  !
289  ! Read types from input file.
290  !
291  !==========================================================
292 
293  CALL nextln ( comstr , ndsi , ndsen )
294  READ (ndsi,*,END=801,ERR=802,IOSTAT=IERR) idfld
295  !
296  IF ( idfld.EQ.'LEV' ) THEN
297  ifld = 2
298  ELSE IF ( idfld.EQ.'CUR' ) THEN
299  ifld = 4
300  ELSE
301  WRITE (ndse,1030) idfld
302  CALL extcde ( 1 )
303  END IF
304  !
305  CALL nextln ( comstr , ndsi , ndsen )
306  READ (ndsi,'(A)',END=801,ERR=802,IOSTAT=IERR) tideconstnames
307  CALL nextln ( comstr , ndsi , ndse )
308  tidecon_prnames(:)=''
309  CALL strsplit(tideconstnames,tidecon_prnames)
310  !
311  CALL nextln ( comstr , ndsi , ndsen )
312  READ (ndsi,'(A)',END=801,ERR=802,IOSTAT=IERR) tideconstnames
313  tidecon_maxnames(:)=''
314  CALL strsplit(tideconstnames,tidecon_maxnames)
315  !
316  CALL nextln ( comstr , ndsi , ndsen )
317  tidecon_maxvals(:)=''
318  READ (ndsi,'(A)',END=801,ERR=802,IOSTAT=IERR) tideconstnames
319  CALL strsplit(tideconstnames,tidecon_maxvals)
320  !
321  CALL nextln ( comstr , ndsi , ndsen )
322  READ (ndsi,*,END=801,ERR=802,IOSTAT=IERR) TIDE_START,PRTIDE_DT,tide_end
323  CALL nextln ( comstr , ndsi , ndsen )
324  READ (ndsi,*,END=801,ERR=802,IOSTAT=IERR) filenamext
325  !
326  CALL w3fldo ('READ', idfld, ndsf, ndst, &
327  ndse, nx, ny, gtype, &
328  ierr, filenamext, '', tideflagin=flagtide )
329  !
330  IF (flagtide.NE.1) GOTO 803
331  !
332  CALL vuf_set_parameters
333 
334  !==========================================================
335  !
336  ! Read tidal amplitudes and phases
337  !
338  !==========================================================
339 
340  CALL w3fldtide1 ( 'READ', ndsf, ndst, ndse, nx, ny, idfld, ierr )
341  CALL w3fldtide2 ( 'READ', ndsf, ndst, ndse, nx, ny, idfld, 0, ierr )
342  CLOSE(ndsf)
343 
344  !
345 
346  IF (gtype.EQ.ungtype) THEN
347 
348  countri = maxval(ccon)
349  ALLOCATE(vneigh(nx,2*countri))
350  ALLOCATE(conn(nx))
351  vneigh(:,:) = 0
352  conn(:) = 0
353  !
354  j = 0
355  DO ip = 1, nx
356  ifound = 0
357  DO ii = 1, ccon(ip)
358  j = j + 1
359  ie = ie_cell(j)
360  IF (ip == trigp(1,ie)) THEN
361  DO ip2=2,3
362  alreadyfound = 0
363  DO i=1,ifound
364  IF (vneigh(ip,i).EQ.trigp(ip2,ie)) alreadyfound=alreadyfound+1
365  END DO
366  IF (alreadyfound.EQ.0) THEN
367  ifound=ifound+1
368  vneigh(ip,ifound)=trigp(ip2,ie)
369  END IF
370  END DO
371  END IF
372 
373  IF (ip == trigp(2,ie)) THEN
374  DO ip2=3,4
375  alreadyfound = 0
376  DO i=1,ifound
377  IF (vneigh(ip,i).EQ.trigp(mod(ip2-1,3)+1,ie)) alreadyfound=alreadyfound+1
378  END DO
379  IF (alreadyfound.EQ.0) THEN
380  ifound=ifound+1
381  vneigh(ip,ifound)=trigp(mod(ip2-1,3)+1,ie)
382  END IF
383  END DO
384  END IF
385 
386  IF (ip == trigp(3,ie)) THEN
387  DO ip2=1,2
388  alreadyfound = 0
389  DO i=1,ifound
390  IF (vneigh(ip,i).EQ.trigp(ip2,ie)) alreadyfound=alreadyfound+1
391  END DO
392  IF (alreadyfound.EQ.0) THEN
393  ifound=ifound+1
394  vneigh(ip,ifound)=trigp(ip2,ie)
395  END IF
396  END DO
397  END IF
398  END DO ! CCON
399  ! CONN is a counter on connected points. In comparison with the number of connected triangle
400  ! CCON, it will enable to spot whether a point belong to the contour
401  !
402  conn(ip)=ifound
403  DO i=2,ifound
404  DO jj=1,i-1
405  IF (vneigh(ip,jj).EQ. vneigh(ip,i)) THEN
406  countcon(ip)=countcon(ip)-1
407  vneigh(ip,i:ifound)=vneigh(ip,i+1:ifound+1) ! removes the double point
408  END IF
409  END DO
410  END DO
411  END DO !NX
412 
413  END IF ! UNGTYPE
414 
415  !==========================================================
416  !
417  ! Apply the maximum threshold value to tidal constituents
418  !
419  !==========================================================
420 
421  CALL tide_find_indices_prediction(tidecon_prnames,pr_inds,tide_prmf)
422  tide_max=0
423  tide_maxi=0
424  DO WHILE (len_trim(tidecon_maxnames(tide_maxi+1)).NE.0)
425  tide_maxi=tide_maxi+1
426  DO j=1,tide_mf
427  IF (trim(tidecon_name(j)).EQ.trim(tidecon_maxnames(tide_maxi))) THEN
428  tide_max=tide_max+1
429  indmax(tide_max)=j
430  READ(tidecon_maxvals(tide_maxi),*) maxvalcon(tide_max)
431  IF (iaproc.EQ.napout) THEN
432  WRITE(ndso,'(A,I8,A,F10.2)') &
433  'Maximum allowed value for amplitude:',&
434  j,trim(tidecon_name(j)),maxvalcon(tide_max)
435  END IF
436  END IF
437  END DO
438  END DO
439 
440  !==========================================================
441  !
442  ! Create the binary output file
443  !
444  !==========================================================
445 
446  flagtide = 0
447  IF (iaproc .EQ. napout) THEN
448  CALL w3fldo ('WRITE', idfld, ndsdat, ndst, ndse, nx, ny, &
449  gtype, ierr, 'ww3', tideflagin=flagtide)
450  END IF
451 
452  !==========================================================
453  !
454  ! Set arrays for MPI exchanges
455  !
456  !==========================================================
457 
458 #ifdef W3_MPI
459  slice=nx/naproc
460  rest=mod(nx,naproc)
461  IF(rest.GE.iaproc) slice=slice+1
462 #endif
463 
464 #ifdef W3_MPI
465  ! set total 1D array (nx)
466  ALLOCATE ( fx1d(nx), fy1d(nx), fa1d(nx))
467  fx1d(:)=0.
468  fy1d(:)=0.
469  fa1d(:)=0.
470 #endif
471 
472 #ifdef W3_MPI
473  ! set local 1D array (slice)
474  ALLOCATE(fx1dl(slice))
475  ALLOCATE(fy1dl(slice))
476  ALLOCATE(fa1dl(slice))
477  fx1dl(:)=0.
478  fy1dl(:)=0.
479  fa1dl(:)=0.
480 #endif
481 
482 
483 #ifdef W3_MPI
484  ! set arrays for number of elements per MPI proc
485  ALLOCATE(nelem(naproc))
486  ALLOCATE(cumul(naproc))
487  nelem(1) = nx / naproc
488  IF (rest .GT. 0) nelem(1) = nelem(1) + 1
489  cumul(1) = 0
490  DO i=2,naproc
491  cumul(i)=cumul(i-1)+nelem(i-1)
492  nelem(i) = nx / naproc
493  IF (rest .GT. i-1) nelem(i) = nelem(i) + 1
494  END DO
495 #endif
496 
497 #ifdef W3_MPIT
498  WRITE(100+iaproc,*) "Number of points for this processor ", iaproc, " : ", nelem(iaproc), ' / ', nx
499  WRITE(100+iaproc,*) "Cumul of points for this processor ", iaproc, " : ", cumul(iaproc), ' / ', nx
500 #endif
501 
502  !==========================================================
503  !
504  ! Loop on time steps
505  !
506  !==========================================================
507 
508  dttst = dsec21( tide_start , tide_end )
509  IF ( dttst .LE. 0. .OR. prtide_dt .LT. 1 ) GOTO 888
510  time = tide_start
511  tide_kd0= 2415020
512  !
513  tindex = 1
514  !
515  DO
516  dttst = dsec21( time, tide_end )
517  IF ( dttst .LT. 0. ) GOTO 888
518  !
519  CALL stme21 ( time , idtime )
520  IF ( iaproc .EQ. napout ) WRITE (ndso,973) idtime
521 
522  tide_hour = time2hours(time)
523  !
524  !* THE ASTRONOMICAL ARGUMENTS ARE CALCULATED
525  !
526  d1=tide_hour/24.d0
527  d1=d1-dfloat(tide_kd0)-0.5d0
528  CALL astr(d1,h,pp,s,p,enp,dh,dpp,ds,dp,dnp)
529  int24=24
530  intdys=int((tide_hour+0.00001)/int24)
531  hh=tide_hour-dfloat(intdys*int24)
532  tau=hh/24.d0+h-s
533 
534  !==========================================================
535  !
536  ! Treatment of 'bad points' at first time step
537  !
538  !==========================================================
539 
540  badpoints(:,:)=0
541  nbad =0
542 
543  IF (tindex.EQ.1) THEN
544  DO iy = 1, ny
545  DO ix=1, nx
546  tideok=1
547  DO i=1,tide_max
548  IF (abs(tidal_const(ix,iy,indmax(i),1,1)) .GT.maxvalcon(i) .OR. &
549  abs(tidal_const(ix,iy,indmax(i),2,1)) .GT.maxvalcon(i)) THEN
550  tideok = 0
551  WRITE(ndso,'(A,I8,F10.2,A,2F10.2)') &
552  '[BAD POINT] GREATER THAN THRESHOLD ', maxvalcon(i), &
553  ' AT INDEX ', indmax(i), &
554  ' WITH X-Y COMPONENTS : ', abs(tidal_const(ix,iy,indmax(i),1:2,1))
555  END IF
556  badpoints(ix,iy) = badpoints(ix,iy) + (1-tideok)
557  END DO
558 
559  IF (badpoints(ix,iy).GT.0) THEN
560  nbad = nbad +1
561  WRITE(ndse,*) 'BAD POINT:',ix,iy,nbad, &
562  tidal_const(ix,iy,:,1,1),'##',tidal_const(ix,iy,:,2,1)
563  END IF
564  END DO
565  END DO
566  !
567  DO iter=1,2
568  DO iy = 1, ny
569  DO ix= 1, nx
570  IF (badpoints(ix,iy).GT.0) THEN
571  tidal_const(ix,iy,:,1,1)=0
572  tidal_const(ix,iy,:,2,1)=0
573 
574 
575  IF (tidefill.AND.(gtype.EQ.ungtype)) THEN
576 
577  !
578  ! Performs a vector sum of tidal constituents over neighbor nodes
579  !
580  DO j=1, tide_mf
581  DO k=1, 2
582  ampcos = 0
583  ampsin = 0
584  sumok = 0
585  DO icon=1,countcon(ix)
586  ix2=vneigh(ix,icon)
587  IF (badpoints(ix2,iy).EQ.0) THEN
588  sumok = sumok + 1
589  ampcos = ampcos+tidal_const(ix2,iy,j,k,1)*cos(tidal_const(ix2,iy,j,k,2)*dera)
590  ampsin = ampsin+tidal_const(ix2,iy,j,k,1)*sin(tidal_const(ix2,iy,j,k,2)*dera)
591  END IF
592  END DO
593  IF (sumok.GT.1) THEN
594  !
595  ! Finalizes the amplitude and phase calculation from COS and SIN. Special case for mean value Z0.
596  !
597  IF (tidecon_name(j).NE.'Z0 ') THEN
598  tidal_const(ix,iy,j,k,1) = sqrt(ampcos**2+ampsin**2)/sumok
599  tidal_const(ix,iy,j,k,2) = atan2(ampsin,ampcos)/dera
600  ELSE
601  tidal_const(ix,iy,j,k,1) = ampcos/sumok
602  tidal_const(ix,iy,j,k,2) = 0.
603  END IF
604  IF(k.EQ.2.AND.j.EQ.tide_mf) THEN
605  nbad=nbad-1
606  badpoints(ix,iy) = 0
607  END IF
608  ENDIF
609  END DO
610  END DO
611  END IF
612  END IF
613  END DO
614  END DO
615  END DO
616  IF ( iaproc .EQ. napout ) WRITE(ndse,*) 'Number of remaining bad points:',nbad
617  END IF
618 
619  !==========================================================
620  !
621  ! For currents: 2 components
622  !
623  !==========================================================
624 
625  IF (ifld.EQ.4) THEN
626  DO iy = 1, ny
627 #ifdef W3_MPI
628  ind=0
629  DO ix=cumul(iaproc)+1,cumul(iaproc)+nelem(iaproc)
630 #endif
631 #ifdef W3_SHRD
632  DO ix=1,nx
633 #endif
634  CALL setvuf_fast(h,pp,s,p,enp,dh,dpp,ds,dp,dnp,tau,real(ygrd(iy,ix)),tide_fx,ux,vx)
635  wcurtidex = 0.
636  wcurtidey = 0.
637  DO i=1,tide_prmf
638  j=pr_inds(i)
639  IF (trim(tidecon_name(j)).EQ.'Z0') THEN
640  wcurtidex = wcurtidex+tidal_const(ix,iy,j,1,1)
641  wcurtidey = wcurtidey+tidal_const(ix,iy,j,2,1)
642  ELSE
643  tide_argx=(vx(j)+ux(j))*twpi-tidal_const(ix,iy,j,1,2)*dera
644  tide_argy=(vx(j)+ux(j))*twpi-tidal_const(ix,iy,j,2,2)*dera
645  wcurtidex = wcurtidex+tide_fx(j)*tidal_const(ix,iy,j,1,1)*cos(tide_argx)
646  wcurtidey = wcurtidey+tide_fx(j)*tidal_const(ix,iy,j,2,1)*cos(tide_argy)
647  END IF
648  END DO
649  IF (abs(wcurtidex).GT.10..OR.abs(wcurtidey).GT.10.) THEN
650  WRITE(ndse,*) &
651  'WARNING: VERY STRONG CURRENT... BAD CONSTITUENTS?', &
652  ix, wcurtidex, wcurtidey , tidal_const(ix,iy,:,1,1),'##',tidal_const(ix,iy,:,2,1)
653  stop
654  END IF
655 #ifdef W3_MPI
656  ind=ind+1
657  fx1dl(ind) = wcurtidex
658  fy1dl(ind) = wcurtidey
659  fa1dl(ind) = 0.
660  END DO ! NX
661 #endif
662 #ifdef W3_SHRD
663  fx(ix,iy) = wcurtidex
664  fy(ix,iy) = wcurtidey
665  fa(ix,iy) = 0.
666  END DO ! NX
667 #endif
668 
669  !
670  ! Gather from other MPI tasks
671  !
672 
673 #ifdef W3_MPI
674  IF (naproc.GT.1) THEN
675  CALL mpi_gatherv(fx1dl, slice, mpi_real, fx1d, nelem, &
676  cumul, mpi_real, napout-1, mpi_comm_world, ierr_mpi)
677  CALL mpi_gatherv(fy1dl, slice, mpi_real, fy1d, nelem, &
678  cumul, mpi_real, napout-1, mpi_comm_world, ierr_mpi)
679  CALL mpi_gatherv(fa1dl, slice, mpi_real, fa1d, nelem, &
680  cumul, mpi_real, napout-1, mpi_comm_world, ierr_mpi)
681  ELSE
682  fx1d = fx1dl
683  fy1d = fy1dl
684  fa1d = fa1dl
685  END IF
686 #endif
687 
688  !
689  ! Convert from 1D to 2D array
690  !
691 #ifdef W3_MPI
692  IF (iaproc .EQ. napout) THEN
693  ind=0
694  DO ix=1,nx
695  ind=ind+1
696  fx(ix,iy)=fx1d(ind)
697  fy(ix,iy)=fy1d(ind)
698  fa(ix,iy)=fa1d(ind)
699  END DO
700  END IF
701 #endif
702 
703  END DO ! NY
704  END IF ! IFLD.EQ.4
705 
706 
707  !==========================================================
708  !
709  ! For water levels: only 1 component
710  !
711  !==========================================================
712 
713  IF (ifld.EQ.2) THEN
714  DO iy = 1, ny
715 #ifdef W3_MPI
716  ind=0
717  DO ix=cumul(iaproc)+1,cumul(iaproc)+nelem(iaproc)
718 #endif
719 #ifdef W3_SHRD
720  DO ix=1,nx
721 #endif
722  CALL setvuf_fast(h,pp,s,p,enp,dh,dpp,ds,dp,dnp,tau,real(ygrd(iy,ix)),tide_fx,ux,vx)
723  !
724  ! Removes unlikely values ...
725  !
726  IF (tindex.EQ.1) THEN
727  tideok=1
728  DO i=1,tide_max
729  IF (abs(tidal_const(ix,iy,indmax(i),1,1)) .GT.maxvalcon(i)) &
730  tideok = 0
731  END DO
732  IF (tideok.EQ.0) THEN
733  WRITE(ndse,*) 'BAD POINT:',ix,iy, tidal_const(ix,iy,:,1,1)
734  tidal_const(ix,iy,:,1,1)=0
735  END IF
736  END IF
737 
738  wcurtidex = 0.
739  DO i=1,tide_prmf
740  j=pr_inds(i)
741  IF (trim(tidecon_name(j)).EQ.'Z0') THEN
742  wcurtidex = wcurtidex+tidal_const(ix,iy,j,1,1)
743  ELSE
744  tide_argx=(vx(j)+ux(j))*twpi-tidal_const(ix,iy,j,1,2)*dera
745  wcurtidex = wcurtidex+tide_fx(j)*tidal_const(ix,iy,j,1,1)*cos(tide_argx)
746  END IF
747  END DO
748 #ifdef W3_MPI
749  ind=ind+1
750  fx1dl(ind) = 0.
751  fy1dl(ind) = 0.
752  fa1dl(ind) = wcurtidex
753  END DO ! NX
754 #endif
755 #ifdef W3_SHRD
756  fx(ix,iy) = 0.
757  fy(ix,iy) = 0.
758  fa(ix,iy) = wcurtidex
759  END DO ! NX
760 #endif
761 
762 
763 
764  !
765  ! Gather from other MPI tasks
766  !
767 
768 #ifdef W3_MPI
769  IF (naproc.GT.1) THEN
770  CALL mpi_gatherv(fx1dl, slice, mpi_real, fx1d, nelem,&
771  cumul, mpi_real, napout-1, mpi_comm_world, ierr_mpi)
772  CALL mpi_gatherv(fy1dl, slice, mpi_real, fy1d, nelem,&
773  cumul, mpi_real, napout-1, mpi_comm_world, ierr_mpi)
774  CALL mpi_gatherv(fa1dl, slice, mpi_real, fa1d, nelem,&
775  cumul, mpi_real, napout-1, mpi_comm_world, ierr_mpi)
776  ELSE
777  fx1d = fx1dl
778  fy1d = fy1dl
779  fa1d = fa1dl
780  END IF
781 #endif
782 
783  !
784  ! Convert from 1D to 2D array
785  !
786 #ifdef W3_MPI
787  IF (iaproc .EQ. napout) THEN
788  ind=0
789  DO ix=1,nx
790  ind=ind+1
791  fx(ix,iy)=fx1d(ind)
792  fy(ix,iy)=fy1d(ind)
793  fa(ix,iy)=fa1d(ind)
794  END DO
795  END IF
796 #endif
797 
798  END DO ! NY
799  END IF ! IFLD.EQ.2
800 
801 
802  !==========================================================
803  !
804  ! Write into binary output file
805  !
806  !==========================================================
807 
808  IF (iaproc .EQ. napout) THEN
809 
810  ! WHERE(FX.NE.FX) FX = 0.
811  ! WHERE(FY.NE.FY) FY = 0.
812  ! WHERE(FA.NE.FA) FA = 0.
813 
814  CALL w3fldg ('WRITE', idfld, ndsdat, ndst, ndse, nx, ny, &
815  nx, ny, time, time, time, fx, fy, fa, time, &
816  fx, fy, fa, ierr)
817  END IF
818 
819  !==========================================================
820  !
821  ! Increment the clock
822  !
823  !==========================================================
824 
825  CALL tick21 ( time, float(prtide_dt) )
826  tindex = tindex +1
827 
828  END DO
829  !
830  GOTO 888
831  !
832  ! Error escape locations
833  !
834 800 CONTINUE
835  WRITE (ndse,1000) ierr
836  CALL extcde ( 40 )
837  !
838 801 CONTINUE
839  WRITE (ndse,1001)
840  CALL extcde ( 41 )
841  !
842 802 CONTINUE
843  WRITE (ndse,1002) ierr
844  CALL extcde ( 42 )
845  !
846 803 CONTINUE
847  WRITE (ndse,1003)
848  CALL extcde ( 43 )
849  !
850 888 CONTINUE
851  IF ( iaproc .EQ. napout ) WRITE (ndso,999)
852 #ifdef W3_MPI
853  CALL mpi_finalize ( ierr_mpi )
854 #endif
855 
856  !
857  ! Formats
858  !
859 900 FORMAT (/15x,' *** WAVEWATCH III tide prediction *** '/ &
860  15x,'==============================================='/)
861 901 FORMAT ( ' Comment character is ''',a,''''/)
862 902 FORMAT ( ' Grid name : ',a/)
863 973 FORMAT ( ' Time : ',a)
864  !
865 999 FORMAT(/' End of program '/ &
866  ' ========================================='/ &
867  ' WAVEWATCH III Input preprocessing '/)
868  !
869 1000 FORMAT (/' *** WAVEWATCH III ERROR IN W3PRTIDE : '/ &
870  ' ERROR IN OPENING INPUT FILE'/ &
871  ' IOSTAT =',i5/)
872  !
873 1001 FORMAT (/' *** WAVEWATCH III ERROR IN W3PRTIDE : '/ &
874  ' PREMATURE END OF INPUT FILE'/)
875  !
876 1002 FORMAT (/' *** WAVEWATCH III ERROR IN W3PRTIDE : '/ &
877  ' ERROR IN READING FROM INPUT FILE'/ &
878  ' IOSTAT =',i5/)
879  !
880 1003 FORMAT (/' *** WAVEWATCH III ERROR IN W3PRTIDE : '/ &
881  ' THE INPUT FILE DOES NOT CONTAIN TIDAL DATA'/)
882  !
883 1030 FORMAT (/' *** WAVEWATCH III ERROR IN W3PRTIDE : '/ &
884  ' ILLEGAL FIELD ID -->',a,'<--'/)
885  !
886  !/
887  !/ End of W3PRTIDE ----------------------------------------------------- /
888  !/
889 END PROGRAM w3prtide
w3servmd::nextln
subroutine nextln(CHCKC, NDSI, NDSE)
Definition: w3servmd.F90:222
w3fldsmd::w3fldd
subroutine w3fldd(INXOUT, IDFLD, NDS, NDST, NDSE, TIME, TD, NR, ND, NDOUT, DATA, IERR)
Definition: w3fldsmd.F90:1474
include
cmake src_list cmake include(${CMAKE_CURRENT_SOURCE_DIR}/cmake/check_switches.cmake) check_switches("$
Definition: CMakeLists.txt:15
w3adatmd
Define data structures to set up wave model auxiliary data for several models simultaneously.
Definition: w3adatmd.F90:26
w3fldsmd::w3fldtide1
subroutine w3fldtide1(INXOUT, NDS, NDST, NDSE, NX, NY, IDFLD, IERR)
Definition: w3fldsmd.F90:531
w3tidemd
Definition: w3tidemd.F90:3
constants::dera
real, parameter dera
DERA Conversion factor from degrees to radians.
Definition: constants.F90:77
w3gsrumd
Definition: w3gsrumd.F90:17
w3odatmd::iaproc
integer, pointer iaproc
Definition: w3odatmd.F90:457
w3prtide
program w3prtide
Predicts tides (current or water level) to be used during run by ww3_shel or ww3_multi (this takes mu...
Definition: ww3_prtide.F90:12
w3odatmd::fnmpre
character(len=80) fnmpre
Definition: w3odatmd.F90:330
w3arrymd::ina2i
subroutine ina2i(ARRAY, MX, MY, LX, HX, LY, HY, NDS, NDST, NDSE, IDFM, RFORM, IDLA, VSC, VOF)
Definition: w3arrymd.F90:295
w3servmd::strsplit
subroutine strsplit(STRING, TAB)
Definition: w3servmd.F90:1440
w3arrymd::ina2r
subroutine ina2r(ARRAY, MX, MY, LX, HX, LY, HY, NDS, NDST, NDSE, IDFM, RFORM, IDLA, VSC, VOF)
Definition: w3arrymd.F90:78
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
w3odatmd::naperr
integer, pointer naperr
Definition: w3odatmd.F90:457
w3servmd
Definition: w3servmd.F90:3
w3odatmd::w3seto
subroutine w3seto(IMOD, NDSERR, NDSTST)
Definition: w3odatmd.F90:1523
w3odatmd
Definition: w3odatmd.F90:3
w3adatmd::w3naux
subroutine w3naux(NDSE, NDST)
Set up the number of grids to be used.
Definition: w3adatmd.F90:704
w3odatmd::naproc
integer, pointer naproc
Definition: w3odatmd.F90:457
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
w3idatmd
Define data structures to set up wave model input data for several models simultaneously.
Definition: w3idatmd.F90:16
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
w3odatmd::napout
integer, pointer napout
Definition: w3odatmd.F90:457
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
w3fldsmd::w3fldg
subroutine w3fldg(INXOUT, IDFLD, NDS, NDST, NDSE, MX, MY, NX, NY, T0, TN, TF0, FX0, FY0, FA0, TFN, FXN, FYN, FAN, IERR, FLAGSC ifdef W3_OASIS
Definition: w3fldsmd.F90:958
w3servmd::extcde
subroutine extcde(IEXIT, UNIT, MSG, FILE, LINE, COMM)
Definition: w3servmd.F90:736
w3odatmd::w3nout
subroutine w3nout(NDSERR, NDSTST)
Definition: w3odatmd.F90:561
w3fldsmd::w3fldtide2
subroutine w3fldtide2(INXOUT, NDS, NDST, NDSE, NX, NY, IDFLD, IDAT, IERR)
Definition: w3fldsmd.F90:722
w3servmd::itrace
subroutine itrace(NDS, NMAX)
Definition: w3servmd.F90:91
w3fldsmd::w3fldo
subroutine w3fldo(INXOUT, IDFLD, NDS, NDST, NDSE, NX, NY, GTYPE, IERR, FEXT, FPRE, FHDR, TIDEFLAGIN)
Definition: w3fldsmd.F90:90
w3timemd
Definition: w3timemd.F90:3
w3fldsmd
Definition: w3fldsmd.F90:3