WAVEWATCH III  beta 0.0.1
w3idatmd Module Reference

Define data structures to set up wave model input data for several models simultaneously. More...

Data Types

type  input
 

Functions/Subroutines

subroutine w3ninp (NDSE, NDST)
 Set up the number of grids to be used. More...
 
subroutine w3dimi (IMOD, NDSE, NDST, FLAGSTIDEIN)
 Initialize an individual data grid at the proper dimensions. More...
 
subroutine w3seti (IMOD, NDSE, NDST)
 Select one of the WAVEWATCH III grids / models. More...
 

Variables

integer nidata = -1
 
integer iidata = -1
 
integer jfirst = 1
 
integer ntide
 
real, dimension(:), allocatable tidefreq
 
type(input), dimension(:), allocatable, target inputs
 
integer, dimension(:,:), pointer tfn
 
integer, dimension(:), pointer tln
 
integer, dimension(:), pointer tc0
 
integer, dimension(:), pointer tcn
 
integer, dimension(:), pointer tw0
 
integer, dimension(:), pointer twn
 
integer, dimension(:), pointer tu0
 
integer, dimension(:), pointer tun
 
integer, dimension(:), pointer tin
 
integer, dimension(:), pointer tr0
 
integer, dimension(:), pointer trn
 
integer, dimension(:), pointer t0n
 
integer, dimension(:), pointer t1n
 
integer, dimension(:), pointer t2n
 
integer, dimension(:), pointer tdn
 
integer, dimension(:), pointer tg0
 
integer, dimension(:), pointer tgn
 
integer, dimension(:), pointer ttn
 
integer, dimension(:), pointer tvn
 
integer, dimension(:), pointer tzn
 
integer, dimension(:), pointer ti1
 
integer, dimension(:), pointer ti2
 
integer, dimension(:), pointer ti3
 
integer, dimension(:), pointer ti4
 
integer, dimension(:), pointer ti5
 
real, pointer ga0
 
real, pointer gd0
 
real, pointer gan
 
real, pointer gdn
 
real, dimension(:,:), pointer wx0
 
real, dimension(:,:), pointer wy0
 
real, dimension(:,:), pointer dt0
 
real, dimension(:,:), pointer wxn
 
real, dimension(:,:), pointer wyn
 
real, dimension(:,:), pointer dtn
 
real, pointer ifdef
 
real, pointer w3_wrst
 
real, dimension(:,:,:,:), pointer cxtide
 
real, dimension(:,:,:,:), pointer cytide
 
real, dimension(:,:,:,:), pointer wltide
 
logical, pointer iinit
 
logical, dimension(:), pointer inflags1
 
logical, dimension(:), pointer inflags2
 
logical, dimension(:), pointer flagsc
 
logical, pointer fllev
 
logical, pointer flcur
 
logical, pointer flwind
 
logical, pointer flice
 
logical, pointer fltaua
 
logical, pointer flrhoa
 
logical, pointer flmth
 
logical, pointer flmvs
 
logical, pointer flmdn
 
logical, pointer flic1
 
logical, pointer flic2
 
logical, pointer flic3
 
logical, pointer flic4
 
logical, pointer flic5
 
logical, pointer fllevtide
 
logical, pointer flcurtide
 
logical, pointer fllevresi
 
logical, pointer flcurresi
 

Detailed Description

Define data structures to set up wave model input data for several models simultaneously.

Author
H. L. Tolman
Date
22-Mar-2021

Function/Subroutine Documentation

◆ w3dimi()

subroutine w3idatmd::w3dimi ( integer, intent(in)  IMOD,
integer, intent(in)  NDSE,
integer, intent(in)  NDST,
logical, dimension(4), intent(in), optional  FLAGSTIDEIN 
)

Initialize an individual data grid at the proper dimensions.

Allocate directly into the structure array. Note that this cannot be done through the pointer alias!

Parameters
[in]IMODModel number to point to.
[in]NDSEError output unit number.
[in]NDSTTest output unit number.
[in]FLAGSTIDEIN
Author
H. L. Tolman
Date
22-Mar-2021

Definition at line 435 of file w3idatmd.F90.

435  !/
436  !/ +-----------------------------------+
437  !/ | WAVEWATCH III NOAA/NCEP |
438  !/ | H. L. Tolman |
439  !/ | FORTRAN 90 |
440  !/ | Last update : 22-Mar-2021 !
441  !/ +-----------------------------------+
442  !/
443  !/ 02-Apr-2004 : Origination. ( version 3.06 )
444  !/ 19-Jul-2006 : Adding auxiliary grids. ( version 3.10 )
445  !/ 04-Oct-2006 : Add filter to array pointers. ( version 3.10 )
446  !/ 10-Dec-2014 : Add checks for allocate status ( version 5.04 )
447  !/ 21-Jun-2018 : Add FSWND input for SMC grid. JGLi ( version 6.04 )
448  !/ 22-Mar-2021 : Momentum and air density support ( version 7.13 )
449  !/
450  ! 1. Purpose :
451  !
452  ! Initialize an individual data grid at the proper dimensions.
453  !
454  ! 2. Method :
455  !
456  ! Allocate directly into the structure array. Note that
457  ! this cannot be done through the pointer alias!
458  !
459  ! 3. Parameters :
460  !
461  ! Parameter list
462  ! ----------------------------------------------------------------
463  ! IMOD Int. I Model number to point to.
464  ! NDSE Int. I Error output unit number.
465  ! NDST Int. I Test output unit number.
466  ! ----------------------------------------------------------------
467  !
468  ! 4. Subroutines used :
469  !
470  ! See module documentation.
471  !
472  ! 5. Called by :
473  !
474  ! Main wave model drivers.
475  !
476  ! 6. Error messages :
477  !
478  ! - Check on input parameters.
479  ! - Check on previous allocation.
480  !
481  ! 7. Remarks :
482  !
483  ! - W3SETI needs to be called after allocation to point to
484  ! proper allocated arrays.
485  !
486  ! 8. Structure :
487  !
488  ! See source code.
489  !
490  ! 9. Switches :
491  !
492  ! !/S Enable subroutine tracing.
493  ! !/T Enable test output
494  !
495  ! 10. Source code :
496  !
497  !/ ------------------------------------------------------------------- /
498  USE w3gdatmd, ONLY: ngrids, nauxgr, igrid, w3setg, nx, ny
499 #ifdef W3_SMC
500  USE w3gdatmd, ONLY: fswnd, nsea
501 #endif
502  USE w3servmd, ONLY: extcde
503 #ifdef W3_S
504  USE w3servmd, ONLY: strace
505 #endif
506  !
507  IMPLICIT NONE
508  !/
509  !/ ------------------------------------------------------------------- /
510  !/ Parameter list
511  !/
512  INTEGER, INTENT(IN) :: IMOD, NDSE, NDST
513  LOGICAL, INTENT(IN), OPTIONAL :: FLAGSTIDEIN(4)
514  !/
515  !/ ------------------------------------------------------------------- /
516  !/ Local parameters
517  !/
518  INTEGER :: JGRID
519  LOGICAL :: FLAGSTIDE(4)=.false.
520 #ifdef W3_S
521  INTEGER, SAVE :: IENT = 0
522  CALL strace (ient, 'W3DIMI')
523 #endif
524  !
525  ! -------------------------------------------------------------------- /
526  ! 1. Test input and module status
527  !
528  IF ( ngrids .EQ. -1 ) THEN
529  WRITE (ndse,1001)
530  CALL extcde (1)
531  END IF
532  !
533  IF ( imod.LT.-nauxgr .OR. imod.GT.nidata ) THEN
534  WRITE (ndse,1002) imod, -nauxgr, nidata
535  CALL extcde (2)
536  END IF
537  !
538  IF ( inputs(imod)%IINIT ) THEN
539  WRITE (ndse,1003)
540  CALL extcde (3)
541  END IF
542  !
543 #ifdef W3_T
544  WRITE (ndst,9000) imod
545 #endif
546  !
547  jgrid = igrid
548  IF ( jgrid .NE. imod ) CALL w3setg ( imod, ndse, ndst )
549  !
550  ! -------------------------------------------------------------------- /
551  ! 2. Allocate arrays
552  !
553 #ifdef W3_TIDE
554  IF ( PRESENT(flagstidein) ) THEN
555  flagstide(:) = flagstidein(:)
556  END IF
557 #endif
558 
559  flic1 => inputs(imod)%INFLAGS1(-7)
560  flic2 => inputs(imod)%INFLAGS1(-6)
561  flic3 => inputs(imod)%INFLAGS1(-5)
562  flic4 => inputs(imod)%INFLAGS1(-4)
563  flic5 => inputs(imod)%INFLAGS1(-3)
564  !
565  flmdn => inputs(imod)%INFLAGS1(-2)
566  flmth => inputs(imod)%INFLAGS1(-1)
567  flmvs => inputs(imod)%INFLAGS1(0)
568  !
569  fllev => inputs(imod)%INFLAGS1(1)
570  flcur => inputs(imod)%INFLAGS1(2)
571 #ifdef W3_TIDE
572  fllevtide => inputs(imod)%INFLAGS1(11)
573  flcurtide => inputs(imod)%INFLAGS1(12)
574  fllevresi => inputs(imod)%INFLAGS1(13)
575  flcurresi => inputs(imod)%INFLAGS1(14)
576  !
577  fllevtide = flagstide(1)
578  flcurtide = flagstide(2)
579  fllevresi = flagstide(3)
580  flcurresi = flagstide(4)
581 #endif
582 
583  flwind => inputs(imod)%INFLAGS1(3)
584  flice => inputs(imod)%INFLAGS1(4)
585  fltaua => inputs(imod)%INFLAGS1(5)
586  flrhoa => inputs(imod)%INFLAGS1(6)
587  !
588  ! notes: future improvement: flags for ICEPx should be
589  ! "all or nothing" rather than 5 individual flags
590 
591  IF ( flic1 ) THEN
592  ALLOCATE ( inputs(imod)%ICEP1(nx,ny), stat=istat )
593  check_alloc_status( istat )
594  END IF
595  IF ( flic2 ) THEN
596  ALLOCATE ( inputs(imod)%ICEP2(nx,ny), stat=istat )
597  check_alloc_status( istat )
598  END IF
599  IF ( flic3 ) THEN
600  ALLOCATE ( inputs(imod)%ICEP3(nx,ny), stat=istat )
601  check_alloc_status( istat )
602  END IF
603  IF ( flic4 ) THEN
604  ALLOCATE ( inputs(imod)%ICEP4(nx,ny), stat=istat )
605  check_alloc_status( istat )
606  END IF
607  IF ( flic5 ) THEN
608  ALLOCATE ( inputs(imod)%ICEP5(nx,ny), stat=istat )
609  check_alloc_status( istat )
610  END IF
611  !
612  IF ( flmdn ) THEN
613  ALLOCATE ( inputs(imod)%MUDD(nx,ny), stat=istat )
614  check_alloc_status( istat )
615  END IF
616  IF ( flmth ) THEN
617  ALLOCATE ( inputs(imod)%MUDT(nx,ny), stat=istat )
618  check_alloc_status( istat )
619  END IF
620  IF ( flmvs ) THEN
621  ALLOCATE ( inputs(imod)%MUDV(nx,ny), stat=istat )
622  check_alloc_status( istat )
623  END IF
624  !
625  IF ( fllev ) THEN
626  ALLOCATE ( inputs(imod)%WLEV(nx,ny), stat=istat )
627  check_alloc_status( istat )
628  END IF
629  !
630  IF ( flcur ) THEN
631 #ifdef W3_SMC
632  IF( fswnd ) THEN
633  ALLOCATE ( inputs(imod)%CX0(nsea,1) , &
634  inputs(imod)%CY0(nsea,1) , &
635  inputs(imod)%CXN(nsea,1) , &
636  inputs(imod)%CYN(nsea,1) , stat=istat )
637  ELSE
638 #endif
639  ALLOCATE ( inputs(imod)%CX0(nx,ny) , &
640  inputs(imod)%CY0(nx,ny) , &
641  inputs(imod)%CXN(nx,ny) , &
642  inputs(imod)%CYN(nx,ny) , stat=istat )
643 #ifdef W3_SMC
644  ENDIF
645 #endif
646  check_alloc_status( istat )
647  END IF
648  !
649 #ifdef W3_TIDE
650  IF ( fllevtide ) THEN
651  ALLOCATE ( inputs(imod)%WLTIDE(nx,ny,ntide,2), stat=istat )
652  check_alloc_status( istat )
653  END IF
654  !
655  IF ( flcurtide ) THEN
656  ALLOCATE ( inputs(imod)%CXTIDE(nx,ny,ntide,2), &
657  inputs(imod)%CYTIDE(nx,ny,ntide,2), stat=istat )
658  check_alloc_status( istat )
659  END IF
660 #endif
661  !
662 
663 #ifdef W3_WRST
664  IF(.NOT.(inputs(imod)%WRSTIINIT)) THEN
665  ALLOCATE ( inputs(imod)%WXNwrst(nx,ny) , &
666  inputs(imod)%WYNwrst(nx,ny) , stat=istat )
667  inputs(imod)%WRSTIINIT=.true.
668  ENDIF
669 #endif
670 
671  IF ( flwind ) THEN
672 #ifdef W3_SMC
673  IF( fswnd ) THEN
674  ALLOCATE ( inputs(imod)%WX0(nsea,1) , &
675  inputs(imod)%WY0(nsea,1) , &
676  inputs(imod)%DT0(nsea,1) , &
677  inputs(imod)%WXN(nsea,1) , &
678  inputs(imod)%WYN(nsea,1) , &
679  inputs(imod)%DTN(nsea,1) , stat=istat )
680  ELSE
681 #endif
682  ALLOCATE ( inputs(imod)%WX0(nx,ny) , &
683  inputs(imod)%WY0(nx,ny) , &
684  inputs(imod)%DT0(nx,ny) , &
685  inputs(imod)%WXN(nx,ny) , &
686  inputs(imod)%WYN(nx,ny) , &
687  inputs(imod)%DTN(nx,ny) , stat=istat )
688 #ifdef W3_SMC
689  ENDIF
690 #endif
691  check_alloc_status( istat )
692  inputs(imod)%DT0 = 0.
693  inputs(imod)%DTN = 0.
694  END IF
695  !
696  IF ( flice ) THEN
697  ALLOCATE ( inputs(imod)%ICEI(nx,ny), &
698  inputs(imod)%BERGI(nx,ny), stat=istat )
699  check_alloc_status( istat )
700  inputs(imod)%BERGI = 0.
701  END IF
702  !
703  IF ( fltaua ) THEN
704 #ifdef W3_SMC
705  IF( fswnd ) THEN
706  ALLOCATE ( inputs(imod)%UX0(nsea,1) , &
707  inputs(imod)%UY0(nsea,1) , &
708  inputs(imod)%UXN(nsea,1) , &
709  inputs(imod)%UYN(nsea,1) , stat=istat )
710  ELSE
711 #endif
712  ALLOCATE ( inputs(imod)%UX0(nx,ny) , &
713  inputs(imod)%UY0(nx,ny) , &
714  inputs(imod)%UXN(nx,ny) , &
715  inputs(imod)%UYN(nx,ny) , stat=istat )
716 #ifdef W3_SMC
717  ENDIF
718 #endif
719  check_alloc_status( istat )
720  END IF
721  !
722  IF ( flrhoa ) THEN
723 #ifdef W3_SMC
724  IF( fswnd ) THEN
725  ALLOCATE ( inputs(imod)%RH0(nsea,1) , &
726  inputs(imod)%RHN(nsea,1) , stat=istat )
727  ELSE
728 #endif
729  ALLOCATE ( inputs(imod)%RH0(nx,ny) , &
730  inputs(imod)%RHN(nx,ny) , stat=istat )
731 #ifdef W3_SMC
732  ENDIF
733 #endif
734  check_alloc_status( istat )
735  END IF
736  !
737  inputs(imod)%IINIT = .true.
738  !
739 #ifdef W3_T
740  WRITE (ndst,9001)
741 #endif
742  !
743  ! -------------------------------------------------------------------- /
744  ! 3. Point to allocated arrays
745  !
746  CALL w3seti ( imod, ndse, ndst )
747  !
748 #ifdef W3_T
749  WRITE (ndst,9002)
750 #endif
751  !
752  ! -------------------------------------------------------------------- /
753  ! 4. Update counters in grid
754  !
755 #ifdef W3_T
756  WRITE (ndst,9003)
757 #endif
758  !
759  ! -------------------------------------------------------------------- /
760  ! 5. Restore previous grid setting if necessary
761  !
762  IF ( jgrid .NE. imod ) CALL w3setg ( jgrid, ndse, ndst )
763  !
764  RETURN
765  !
766  ! Check inputs for stresses
767  IF(fltaua) THEN
768 #ifdef W3_FLX0
769  WRITE (ndse,*) " *** WARNING W3DIMI : TAUA NOT USED *** "
770 #endif
771 #ifdef W3_FLX1
772  WRITE (ndse,*) " *** WARNING W3DIMI : TAUA NOT USED *** "
773 #endif
774 #ifdef W3_FLX2
775  WRITE (ndse,*) " *** WARNING W3DIMI : TAUA NOT USED *** "
776 #endif
777 #ifdef W3_FLX3
778  WRITE (ndse,*) " *** WARNING W3DIMI : TAUA NOT USED *** "
779 #endif
780 #ifdef W3_FLX4
781  WRITE (ndse,*) " *** WARNING W3DIMI : TAUA NOT USED *** "
782 #endif
783  END IF
784  !
785  ! Formats
786  !
787 1001 FORMAT (/' *** ERROR W3DIMI : GRIDS NOT INITIALIZED *** '/ &
788  ' RUN W3NMOD FIRST '/)
789 1002 FORMAT (/' *** ERROR W3DIMI : ILLEGAL MODEL NUMBER *** '/ &
790  ' IMOD = ',i10/ &
791  ' NAUXGR = ',i10/ &
792  ' NIDATA = ',i10/)
793 1003 FORMAT (/' *** ERROR W3DIMI : ARRAY(S) ALREADY ALLOCATED *** ')
794  !
795 #ifdef W3_T
796 9000 FORMAT (' TEST W3DIMI : MODEL ',i4,' DIM. AT ',2i5,i7)
797 9001 FORMAT (' TEST W3DIMI : ARRAYS ALLOCATED')
798 9002 FORMAT (' TEST W3DIMI : POINTERS RESET')
799 9003 FORMAT (' TEST W3DIMI : DIMENSIONS STORED')
800 #endif
801  !/
802  !/ End of W3DIMI ----------------------------------------------------- /
803  !/

References w3servmd::extcde(), flcur, flcurresi, flcurtide, flic1, flic2, flic3, flic4, flic5, flice, fllev, fllevresi, fllevtide, flmdn, flmth, flmvs, flrhoa, fltaua, flwind, w3gdatmd::fswnd, w3gdatmd::igrid, inputs, w3gdatmd::nauxgr, w3gdatmd::ngrids, nidata, w3gdatmd::nsea, ntide, w3gdatmd::nx, w3gdatmd::ny, w3servmd::strace(), w3gdatmd::w3setg(), and w3seti().

Referenced by w3initmd::w3init(), wminitmd::wminit(), and wminitmd::wminitnml().

◆ w3ninp()

subroutine w3idatmd::w3ninp ( integer, intent(in)  NDSE,
integer, intent(in)  NDST 
)

Set up the number of grids to be used.

Use data stored in NGRIDS in W3GDATMD.

Parameters
[in]NDSEError output unit number.
[in]NDSTTest output unit number.
Author
H. L. Tolman
Date
22-Mar-2021

Definition at line 283 of file w3idatmd.F90.

283  !/
284  !/ +-----------------------------------+
285  !/ | WAVEWATCH III NOAA/NCEP |
286  !/ | H. L. Tolman |
287  !/ | FORTRAN 90 |
288  !/ | Last update : 22-Mar-2021 !
289  !/ +-----------------------------------+
290  !/
291  !/ 02-Apr-2004 : Origination. ( version 3.06 )
292  !/ 19-Jul-2006 : Adding auxiliary grids. ( version 3.10 )
293  !/ 10-Dec-2014 : Add checks for allocate status ( version 5.04 )
294  !/ 22-Mar-2021 : Momentum and air density support ( version 7.13 )
295  !/
296  ! 1. Purpose :
297  !
298  ! Set up the number of grids to be used.
299  !
300  ! 2. Method :
301  !
302  ! Use data stored in NGRIDS in W3GDATMD.
303  !
304  ! 3. Parameters :
305  !
306  ! Parameter list
307  ! ----------------------------------------------------------------
308  ! NDSE Int. I Error output unit number.
309  ! NDST Int. I Test output unit number.
310  ! ----------------------------------------------------------------
311  !
312  ! 4. Subroutines used :
313  !
314  ! See module documentation.
315  !
316  ! 5. Called by :
317  !
318  ! Any program that uses this grid structure.
319  !
320  ! 6. Error messages :
321  !
322  ! - Error checks on previous setting of variable NGRIDS.
323  !
324  ! 7. Remarks :
325  !
326  ! 8. Structure :
327  !
328  ! 9. Switches :
329  !
330  ! !/S Enable subroutine tracing.
331  ! !/T Enable test output
332  !
333  ! 10. Source code :
334  !
335  !/ ------------------------------------------------------------------- /
336  USE w3gdatmd, ONLY: ngrids, nauxgr
337  USE w3servmd, ONLY: extcde
338 #ifdef W3_S
339  USE w3servmd, ONLY: strace
340 #endif
341  !
342  IMPLICIT NONE
343  !/
344  !/ ------------------------------------------------------------------- /
345  !/ Parameter list
346  !/
347  INTEGER, INTENT(IN) :: NDSE, NDST
348  !/
349  !/ ------------------------------------------------------------------- /
350  !/ Local parameters
351  !/
352  INTEGER :: I
353 #ifdef W3_S
354  INTEGER, SAVE :: IENT = 0
355  CALL strace (ient, 'W3NINP')
356 #endif
357  !
358  ! -------------------------------------------------------------------- /
359  ! 1. Test input and module status
360  !
361  IF ( ngrids .EQ. -1 ) THEN
362  WRITE (ndse,1001) ngrids
363  CALL extcde (1)
364  END IF
365  !
366  ! -------------------------------------------------------------------- /
367  ! 2. Set variable and allocate arrays
368  !
369  ALLOCATE ( inputs(-nauxgr:ngrids), stat=istat )
370  check_alloc_status( istat )
371  nidata = ngrids
372  !
373  ! -------------------------------------------------------------------- /
374  ! 3. Initialize parameters
375  !
376  DO i=-nauxgr, ngrids
377  inputs(i)%TFN(1,:) = -1
378  inputs(i)%TFN(2,:) = 0
379  inputs(i)%TC0(1) = -1
380  inputs(i)%TC0(2) = 0
381  inputs(i)%TW0(1) = -1
382  inputs(i)%TW0(2) = 0
383  inputs(i)%TU0(1) = -1
384  inputs(i)%TU0(2) = 0
385  inputs(i)%TR0(1) = -1
386  inputs(i)%TR0(2) = 0
387  inputs(i)%TDN(1) = -1
388  inputs(i)%TDN(2) = 0
389  inputs(i)%TG0(1) = -1
390  inputs(i)%TG0(2) = 0
391  inputs(i)%GA0 = 0.
392  inputs(i)%GD0 = 0.
393  inputs(i)%GAN = 0.
394  inputs(i)%GDN = 0.
395  inputs(i)%IINIT = .false.
396  inputs(i)%INFLAGS1 = .false.
397  inputs(i)%INFLAGS2 = .false.
398  inputs(i)%FLAGSC = .false.
399  END DO
400  !
401 #ifdef W3_T
402  WRITE (ndst,9000) -nauxgr, ngrids
403 #endif
404  !
405  RETURN
406  !
407  ! Formats
408  !
409 1001 FORMAT (/' *** ERROR W3NINP : NGRIDS NOT YET SET *** '/ &
410  ' NGRIDS = ',i10/ &
411  ' RUN W3NMOD FIRST'/)
412  !
413 #ifdef W3_T
414 9000 FORMAT (' TEST W3NINP : SETTING UP FOR ',i2,' -',i3,' GRIDS')
415 #endif
416  !/
417  !/ End of W3NINP ----------------------------------------------------- /
418  !/

References w3servmd::extcde(), inputs, w3gdatmd::nauxgr, w3gdatmd::ngrids, nidata, and w3servmd::strace().

Referenced by w3shel(), w3strt(), w3uprstr(), wminitmd::wminit(), and wminitmd::wminitnml().

◆ w3seti()

subroutine w3idatmd::w3seti ( integer, intent(in)  IMOD,
integer, intent(in)  NDSE,
integer, intent(in)  NDST 
)

Select one of the WAVEWATCH III grids / models.

Point pointers to the proper variables in the proper element of the GRIDS array.

Parameters
[in]IMODModel number to point to.
[in]NDSEError output unit number.
[in]NDSTTest output unit number.
Author
H. L. Tolman
Date
22-Mar-2021

Definition at line 819 of file w3idatmd.F90.

819  !/
820  !/ +-----------------------------------+
821  !/ | WAVEWATCH III NOAA/NCEP |
822  !/ | H. L. Tolman |
823  !/ | FORTRAN 90 |
824  !/ | Last update : 22-Mar-2021 !
825  !/ +-----------------------------------+
826  !/
827  !/ 02-Apr-2004 : Origination. ( version 3.06 )
828  !/ 19-Jul-2006 : Adding auxiliary grids. ( version 3.10 )
829  !/ 04-Oct-2006 : Add filter to array pointers. ( version 3.10 )
830  !/ 22-Mar-2021 : Momentum and air density support ( version 7.13 )
831  !/
832  ! 1. Purpose :
833  !
834  ! Select one of the WAVEWATCH III grids / models.
835  !
836  ! 2. Method :
837  !
838  ! Point pointers to the proper variables in the proper element of
839  ! the GRIDS array.
840  !
841  ! 3. Parameters :
842  !
843  ! Parameter list
844  ! ----------------------------------------------------------------
845  ! IMOD Int. I Model number to point to.
846  ! NDSE Int. I Error output unit number.
847  ! NDST Int. I Test output unit number.
848  ! ----------------------------------------------------------------
849  !
850  ! 4. Subroutines used :
851  !
852  ! See module documentation.
853  !
854  ! 5. Called by :
855  !
856  ! Any subroutine.
857  !
858  ! 6. Error messages :
859  !
860  ! Many subroutines in the WAVEWATCH system.
861  !
862  ! 7. Remarks :
863  !
864  ! 8. Structure :
865  !
866  ! 9. Switches :
867  !
868  ! !/S Enable subroutine tracing.
869  ! !/T Enable test output
870  !
871  ! 10. Source code :
872  !
873  !/ ------------------------------------------------------------------- /
874  USE w3gdatmd, ONLY: nauxgr
875  !
876  USE w3servmd, ONLY: extcde
877 #ifdef W3_S
878  USE w3servmd, ONLY: strace
879 #endif
880  !
881  IMPLICIT NONE
882  !/
883  !/ ------------------------------------------------------------------- /
884  !/ Parameter list
885  !/
886  INTEGER, INTENT(IN) :: IMOD, NDSE, NDST
887  !/
888  !/ ------------------------------------------------------------------- /
889  !/ Local parameters
890  !/
891 #ifdef W3_S
892  INTEGER, SAVE :: IENT = 0
893  CALL strace (ient, 'W3SETI')
894 #endif
895  !
896  ! -------------------------------------------------------------------- /
897  ! 1. Test input and module status
898  !
899  IF ( nidata .EQ. -1 ) THEN
900  WRITE (ndse,1001)
901  CALL extcde (1)
902  END IF
903  !
904  IF ( imod.LT.-nauxgr .OR. imod.GT.nidata ) THEN
905  WRITE (ndse,1002) imod, -nauxgr, nidata
906  CALL extcde (2)
907  END IF
908  !
909 #ifdef W3_T
910  WRITE (ndst,9000) imod
911 #endif
912  !
913  ! -------------------------------------------------------------------- /
914  ! 2. Set model numbers
915  !
916  iidata = imod
917  !
918  ! -------------------------------------------------------------------- /
919  ! 3. Set pointers
920  !
921  tfn => inputs(imod)%TFN
922  tc0 => inputs(imod)%TC0
923  tw0 => inputs(imod)%TW0
924  tu0 => inputs(imod)%TU0
925  tr0 => inputs(imod)%TR0
926  tg0 => inputs(imod)%TG0
927  tdn => inputs(imod)%TDN
928  !
929  ti1 => inputs(imod)%TFN(:,-7)
930  ti2 => inputs(imod)%TFN(:,-6)
931  ti3 => inputs(imod)%TFN(:,-5)
932  ti4 => inputs(imod)%TFN(:,-4)
933  ti5 => inputs(imod)%TFN(:,-3)
934  !
935  tzn => inputs(imod)%TFN(:,-2)
936  ttn => inputs(imod)%TFN(:,-1)
937  tvn => inputs(imod)%TFN(:,0)
938  !
939  tln => inputs(imod)%TFN(:,1)
940  tcn => inputs(imod)%TFN(:,2)
941  twn => inputs(imod)%TFN(:,3)
942  tin => inputs(imod)%TFN(:,4)
943  tun => inputs(imod)%TFN(:,5)
944  trn => inputs(imod)%TFN(:,6)
945  t0n => inputs(imod)%TFN(:,7)
946  t1n => inputs(imod)%TFN(:,8)
947  t2n => inputs(imod)%TFN(:,9)
948  tgn => inputs(imod)%TFN(:,10)
949  !
950  ga0 => inputs(imod)%GA0
951  gd0 => inputs(imod)%GD0
952  gan => inputs(imod)%GAN
953  gdn => inputs(imod)%GDN
954  !
955  iinit => inputs(imod)%IINIT
956  inflags1 => inputs(imod)%INFLAGS1
957  inflags2 => inputs(imod)%INFLAGS2
958  flagsc => inputs(imod)%FLAGSC
959  !
960  flic1 => inputs(imod)%INFLAGS1(-7)
961  flic2 => inputs(imod)%INFLAGS1(-6)
962  flic3 => inputs(imod)%INFLAGS1(-5)
963  flic4 => inputs(imod)%INFLAGS1(-4)
964  flic5 => inputs(imod)%INFLAGS1(-3)
965  !
966  flmdn => inputs(imod)%INFLAGS1(-2)
967  flmth => inputs(imod)%INFLAGS1(-1)
968  flmvs => inputs(imod)%INFLAGS1(0)
969  !
970  fllev => inputs(imod)%INFLAGS1(1)
971  flcur => inputs(imod)%INFLAGS1(2)
972 #ifdef W3_TIDE
973  fllevtide => inputs(imod)%INFLAGS1(11)
974  flcurtide => inputs(imod)%INFLAGS1(12)
975  fllevresi => inputs(imod)%INFLAGS1(13)
976  flcurresi => inputs(imod)%INFLAGS1(14)
977 #endif
978 
979  flwind => inputs(imod)%INFLAGS1(3)
980  flice => inputs(imod)%INFLAGS1(4)
981  fltaua => inputs(imod)%INFLAGS1(5)
982  flrhoa => inputs(imod)%INFLAGS1(6)
983  !
984  IF ( iinit ) THEN
985  !
986  IF ( flic1 ) THEN
987  icep1 => inputs(imod)%ICEP1
988  END IF
989  IF ( flic2 ) THEN
990  icep2 => inputs(imod)%ICEP2
991  END IF
992  IF ( flic3 ) THEN
993  icep3 => inputs(imod)%ICEP3
994  END IF
995  IF ( flic4 ) THEN
996  icep4 => inputs(imod)%ICEP4
997  END IF
998  IF ( flic5 ) THEN
999  icep5 => inputs(imod)%ICEP5
1000  END IF
1001  !
1002  IF ( flmdn ) THEN
1003  mudd => inputs(imod)%MUDD
1004  END IF
1005  IF ( flmth ) THEN
1006  mudt => inputs(imod)%MUDT
1007  END IF
1008  IF ( flmvs ) THEN
1009  mudv => inputs(imod)%MUDV
1010  END IF
1011  !
1012  IF ( fllev ) THEN
1013  wlev => inputs(imod)%WLEV
1014  END IF
1015  !
1016  IF ( flcur ) THEN
1017  cx0 => inputs(imod)%CX0
1018  cy0 => inputs(imod)%CY0
1019  cxn => inputs(imod)%CXN
1020  cyn => inputs(imod)%CYN
1021  END IF
1022 #ifdef W3_TIDE
1023  IF ( fllevtide ) THEN
1024  wltide => inputs(imod)%WLTIDE
1025  END IF
1026  IF ( flcurtide ) THEN
1027  cxtide => inputs(imod)%CXTIDE
1028  cytide => inputs(imod)%CYTIDE
1029  END IF
1030 #endif
1031  !
1032 #ifdef W3_WRST
1033  wxnwrst => inputs(imod)%WXNwrst
1034  wynwrst => inputs(imod)%WYNwrst
1035 #endif
1036 
1037  IF ( flwind ) THEN
1038  wx0 => inputs(imod)%WX0
1039  wy0 => inputs(imod)%WY0
1040  dt0 => inputs(imod)%DT0
1041  wxn => inputs(imod)%WXN
1042  wyn => inputs(imod)%WYN
1043  dtn => inputs(imod)%DTN
1044  END IF
1045  !
1046  IF ( flice ) THEN
1047  icei => inputs(imod)%ICEI
1048  bergi => inputs(imod)%BERGI
1049  END IF
1050  !
1051  IF ( fltaua ) THEN
1052  ux0 => inputs(imod)%UX0
1053  uy0 => inputs(imod)%UY0
1054  uxn => inputs(imod)%UXN
1055  uyn => inputs(imod)%UYN
1056  END IF
1057  !
1058  IF ( flrhoa ) THEN
1059  rh0 => inputs(imod)%RH0
1060  rhn => inputs(imod)%RHN
1061  END IF
1062  !
1063  END IF
1064  !
1065  RETURN
1066  !
1067  ! Formats
1068  !
1069 1001 FORMAT (/' *** ERROR W3SETI : GRIDS NOT INITIALIZED *** '/ &
1070  ' RUN W3NMOD FIRST '/)
1071 1002 FORMAT (/' *** ERROR W3SETI : ILLEGAL MODEL NUMBER *** '/ &
1072  ' IMOD = ',i10/ &
1073  ' NAUXGR = ',i10/ &
1074  ' NIDATA = ',i10/)
1075  !
1076 #ifdef W3_T
1077 9000 FORMAT (' TEST W3SETI : MODEL ',i4,' SELECTED')
1078 #endif
1079  !/
1080  !/ End of W3SETI ----------------------------------------------------- /
1081  !/

References cxtide, cytide, dt0, dtn, w3servmd::extcde(), flagsc, flcur, flcurresi, flcurtide, flic1, flic2, flic3, flic4, flic5, flice, fllev, fllevresi, fllevtide, flmdn, flmth, flmvs, flrhoa, fltaua, flwind, ga0, gan, gd0, gdn, iidata, iinit, inflags1, inflags2, inputs, w3gdatmd::nauxgr, nidata, w3servmd::strace(), t0n, t1n, t2n, tc0, tcn, tdn, tfn, tg0, tgn, ti1, ti2, ti3, ti4, ti5, tin, tln, tr0, trn, ttn, tu0, tun, tvn, tw0, twn, tzn, wltide, wx0, wxn, wy0, and wyn.

Referenced by wmesmfmd::createexpmesh(), w3dimi(), w3initmd::w3init(), w3iorsmd::w3iors(), w3sbs1(), w3shel(), w3uprstr(), w3wavemd::w3wave(), wminitmd::wminit(), wminitmd::wminitnml(), and wmupdtmd::wmupdt().

Variable Documentation

◆ cxtide

real, dimension(:,:,:,:), pointer w3idatmd::cxtide

Definition at line 256 of file w3idatmd.F90.

256  REAL, POINTER :: CXTIDE(:,:,:,:), &
257  CYTIDE(:,:,:,:), WLTIDE(:,:,:,:)

Referenced by w3seti(), and w3updtmd::w3ucur().

◆ cytide

real, dimension(:,:,:,:), pointer w3idatmd::cytide

Definition at line 256 of file w3idatmd.F90.

Referenced by w3seti(), and w3updtmd::w3ucur().

◆ dt0

real, dimension(:,:), pointer w3idatmd::dt0

Definition at line 243 of file w3idatmd.F90.

Referenced by w3seti(), w3updtmd::w3uwnd(), and wmupdtmd::wmupd1().

◆ dtn

real, dimension(:,:), pointer w3idatmd::dtn

Definition at line 243 of file w3idatmd.F90.

Referenced by w3seti(), w3updtmd::w3uwnd(), and wmupdtmd::wmupd1().

◆ flagsc

logical, dimension(:), pointer w3idatmd::flagsc

Definition at line 260 of file w3idatmd.F90.

Referenced by w3seti(), and w3shel().

◆ flcur

◆ flcurresi

logical, pointer w3idatmd::flcurresi

Definition at line 266 of file w3idatmd.F90.

Referenced by w3dimi(), and w3seti().

◆ flcurtide

logical, pointer w3idatmd::flcurtide

Definition at line 266 of file w3idatmd.F90.

Referenced by w3dimi(), w3seti(), and w3updtmd::w3ucur().

◆ flic1

logical, pointer w3idatmd::flic1

Definition at line 264 of file w3idatmd.F90.

264  LOGICAL, POINTER :: FLIC1, FLIC2, FLIC3, FLIC4, FLIC5

Referenced by w3dimi(), w3initmd::w3init(), w3seti(), w3updtmd::w3uic1(), and w3wavemd::w3wave().

◆ flic2

logical, pointer w3idatmd::flic2

Definition at line 264 of file w3idatmd.F90.

Referenced by w3dimi(), w3initmd::w3init(), w3seti(), and w3wavemd::w3wave().

◆ flic3

logical, pointer w3idatmd::flic3

Definition at line 264 of file w3idatmd.F90.

Referenced by w3dimi(), w3initmd::w3init(), w3seti(), and w3wavemd::w3wave().

◆ flic4

logical, pointer w3idatmd::flic4

Definition at line 264 of file w3idatmd.F90.

Referenced by w3dimi(), w3initmd::w3init(), w3seti(), and w3wavemd::w3wave().

◆ flic5

logical, pointer w3idatmd::flic5

Definition at line 264 of file w3idatmd.F90.

Referenced by w3dimi(), w3initmd::w3init(), w3seti(), and w3wavemd::w3wave().

◆ flice

logical, pointer w3idatmd::flice

Definition at line 261 of file w3idatmd.F90.

Referenced by w3dimi(), w3initmd::w3init(), w3seti(), and w3wavemd::w3wave().

◆ fllev

◆ fllevresi

logical, pointer w3idatmd::fllevresi

Definition at line 266 of file w3idatmd.F90.

Referenced by w3dimi(), and w3seti().

◆ fllevtide

logical, pointer w3idatmd::fllevtide

Definition at line 266 of file w3idatmd.F90.

266  LOGICAL, POINTER :: FLLEVTIDE, FLCURTIDE, &
267  FLLEVRESI, FLCURRESI

Referenced by w3dimi(), w3seti(), and w3updtmd::w3ulev().

◆ flmdn

logical, pointer w3idatmd::flmdn

Definition at line 263 of file w3idatmd.F90.

Referenced by w3dimi(), w3initmd::w3init(), and w3seti().

◆ flmth

logical, pointer w3idatmd::flmth

Definition at line 263 of file w3idatmd.F90.

263  LOGICAL, POINTER :: FLMTH, FLMVS, FLMDN

Referenced by w3dimi(), w3initmd::w3init(), and w3seti().

◆ flmvs

logical, pointer w3idatmd::flmvs

Definition at line 263 of file w3idatmd.F90.

Referenced by w3dimi(), w3initmd::w3init(), and w3seti().

◆ flrhoa

logical, pointer w3idatmd::flrhoa

Definition at line 261 of file w3idatmd.F90.

Referenced by w3adatmd::w3dima(), w3dimi(), w3initmd::w3init(), w3seti(), and w3wavemd::w3wave().

◆ fltaua

logical, pointer w3idatmd::fltaua

Definition at line 261 of file w3idatmd.F90.

Referenced by w3adatmd::w3dima(), w3dimi(), w3initmd::w3init(), w3seti(), and w3wavemd::w3wave().

◆ flwind

logical, pointer w3idatmd::flwind

◆ ga0

real, pointer w3idatmd::ga0

Definition at line 242 of file w3idatmd.F90.

242  REAL, POINTER :: GA0, GD0, GAN, GDN

Referenced by w3seti(), w3wavemd::w3wave(), and wmupdtmd::wmupd1().

◆ gan

real, pointer w3idatmd::gan

Definition at line 242 of file w3idatmd.F90.

Referenced by w3seti(), w3wavemd::w3wave(), and wmupdtmd::wmupd1().

◆ gd0

real, pointer w3idatmd::gd0

Definition at line 242 of file w3idatmd.F90.

Referenced by w3seti(), w3wavemd::w3wave(), and wmupdtmd::wmupd1().

◆ gdn

real, pointer w3idatmd::gdn

Definition at line 242 of file w3idatmd.F90.

Referenced by w3seti(), w3wavemd::w3wave(), and wmupdtmd::wmupd1().

◆ ifdef

real, pointer w3idatmd::ifdef

Definition at line 243 of file w3idatmd.F90.

◆ iidata

integer w3idatmd::iidata = -1

Definition at line 160 of file w3idatmd.F90.

Referenced by w3seti(), and w3wavemd::w3wave().

◆ iinit

logical, pointer w3idatmd::iinit

Definition at line 259 of file w3idatmd.F90.

259  LOGICAL, POINTER :: IINIT

Referenced by w3seti(), wminitmd::wminit(), and wminitmd::wminitnml().

◆ inflags1

logical, dimension(:), pointer w3idatmd::inflags1

Definition at line 260 of file w3idatmd.F90.

260  LOGICAL, POINTER :: INFLAGS1(:), INFLAGS2(:), FLAGSC(:)

Referenced by w3sbt8md::w3sbt8(), w3sbt9md::w3sbt9(), w3seti(), w3shel(), w3wavemd::w3wave(), wminitmd::wminit(), wminitmd::wminitnml(), wmupdtmd::wmupd1(), wmupdtmd::wmupdt(), and wmwavemd::wmwave().

◆ inflags2

◆ inputs

type(input), dimension(:), allocatable, target w3idatmd::inputs

Definition at line 232 of file w3idatmd.F90.

232  TYPE(INPUT), TARGET, ALLOCATABLE :: INPUTS(:)

Referenced by w3dimi(), w3ninp(), w3adatmd::w3seta(), w3seti(), w3adatmd::w3xeta(), wminitmd::wminit(), wminitmd::wminitnml(), and wmupdtmd::wmupd2().

◆ jfirst

integer w3idatmd::jfirst = 1

Definition at line 162 of file w3idatmd.F90.

162  INTEGER :: JFIRST = 1

Referenced by w3shel(), wminitmd::wminit(), wminitmd::wminitnml(), and wmupdtmd::wmupdt().

◆ nidata

integer w3idatmd::nidata = -1

Definition at line 160 of file w3idatmd.F90.

160  INTEGER :: NIDATA = -1, iidata = -1

Referenced by w3dimi(), w3ninp(), and w3seti().

◆ ntide

integer w3idatmd::ntide

Definition at line 165 of file w3idatmd.F90.

165  INTEGER :: NTIDE ! number of tidal constituents

Referenced by w3dimi(), w3updtmd::w3ucur(), and w3updtmd::w3ulev().

◆ t0n

integer, dimension(:), pointer w3idatmd::t0n

Definition at line 236 of file w3idatmd.F90.

Referenced by w3seti(), wmupdtmd::wmupd1(), and wmupdtmd::wmupdt().

◆ t1n

integer, dimension(:), pointer w3idatmd::t1n

Definition at line 236 of file w3idatmd.F90.

Referenced by w3seti(), wmupdtmd::wmupd1(), and wmupdtmd::wmupdt().

◆ t2n

integer, dimension(:), pointer w3idatmd::t2n

Definition at line 236 of file w3idatmd.F90.

Referenced by w3seti(), wmupdtmd::wmupd1(), and wmupdtmd::wmupdt().

◆ tc0

integer, dimension(:), pointer w3idatmd::tc0

◆ tcn

integer, dimension(:), pointer w3idatmd::tcn

◆ tdn

integer, dimension(:), pointer w3idatmd::tdn

Definition at line 236 of file w3idatmd.F90.

Referenced by w3seti(), w3wavemd::w3wave(), wmupdtmd::wmupd1(), and wmupdtmd::wmupdt().

◆ tfn

integer, dimension(:,:), pointer w3idatmd::tfn

Definition at line 236 of file w3idatmd.F90.

236  INTEGER, POINTER :: TFN(:,:), TLN(:), TC0(:), TCN(:), &
237  TW0(:), TWN(:), TU0(:), TUN(:), &
238  TIN(:), TR0(:), TRN(:), T0N(:), &
239  T1N(:), T2N(:), TDN(:), TG0(:), &
240  TGN(:), TTN(:), TVN(:), TZN(:), &
241  TI1(:), TI2(:), TI3(:), TI4(:), TI5(:)

Referenced by w3seti(), and wmupdtmd::wmupdt().

◆ tg0

integer, dimension(:), pointer w3idatmd::tg0

Definition at line 236 of file w3idatmd.F90.

Referenced by w3seti(), w3wavemd::w3wave(), wmupdtmd::wmupd1(), and wmupdtmd::wmupdt().

◆ tgn

integer, dimension(:), pointer w3idatmd::tgn

Definition at line 236 of file w3idatmd.F90.

Referenced by w3seti(), w3wavemd::w3wave(), wmupdtmd::wmupd1(), and wmupdtmd::wmupdt().

◆ ti1

integer, dimension(:), pointer w3idatmd::ti1

◆ ti2

integer, dimension(:), pointer w3idatmd::ti2

Definition at line 236 of file w3idatmd.F90.

Referenced by w3seti(), wmupdtmd::wmupd1(), and wmupdtmd::wmupdt().

◆ ti3

integer, dimension(:), pointer w3idatmd::ti3

Definition at line 236 of file w3idatmd.F90.

Referenced by w3seti(), wmupdtmd::wmupd1(), and wmupdtmd::wmupdt().

◆ ti4

integer, dimension(:), pointer w3idatmd::ti4

Definition at line 236 of file w3idatmd.F90.

Referenced by w3seti(), wmupdtmd::wmupd1(), and wmupdtmd::wmupdt().

◆ ti5

integer, dimension(:), pointer w3idatmd::ti5

◆ tidefreq

real, dimension(:), allocatable w3idatmd::tidefreq

Definition at line 166 of file w3idatmd.F90.

166  REAL, ALLOCATABLE :: TIDEFREQ(:)

◆ tin

integer, dimension(:), pointer w3idatmd::tin

◆ tln

integer, dimension(:), pointer w3idatmd::tln

◆ tr0

integer, dimension(:), pointer w3idatmd::tr0

◆ trn

integer, dimension(:), pointer w3idatmd::trn

◆ ttn

integer, dimension(:), pointer w3idatmd::ttn

Definition at line 236 of file w3idatmd.F90.

Referenced by w3seti(), wmupdtmd::wmupd1(), and wmupdtmd::wmupdt().

◆ tu0

integer, dimension(:), pointer w3idatmd::tu0

◆ tun

integer, dimension(:), pointer w3idatmd::tun

◆ tvn

integer, dimension(:), pointer w3idatmd::tvn

Definition at line 236 of file w3idatmd.F90.

Referenced by w3seti(), wmupdtmd::wmupd1(), and wmupdtmd::wmupdt().

◆ tw0

integer, dimension(:), pointer w3idatmd::tw0

◆ twn

integer, dimension(:), pointer w3idatmd::twn

◆ tzn

integer, dimension(:), pointer w3idatmd::tzn

Definition at line 236 of file w3idatmd.F90.

Referenced by w3seti(), wmupdtmd::wmupd1(), and wmupdtmd::wmupdt().

◆ w3_wrst

real, pointer w3idatmd::w3_wrst

Definition at line 243 of file w3idatmd.F90.

◆ wltide

real, dimension(:,:,:,:), pointer w3idatmd::wltide

Definition at line 256 of file w3idatmd.F90.

Referenced by w3seti(), and w3updtmd::w3ulev().

◆ wx0

real, dimension(:,:), pointer w3idatmd::wx0

Definition at line 243 of file w3idatmd.F90.

243  REAL, POINTER :: WX0(:,:), WY0(:,:), DT0(:,:), &
244  WXN(:,:), WYN(:,:), DTN(:,:), &
245 #ifdef W3_WRST

Referenced by w3seti(), w3updtmd::w3uwnd(), and wmupdtmd::wmupd1().

◆ wxn

real, dimension(:,:), pointer w3idatmd::wxn

Definition at line 243 of file w3idatmd.F90.

Referenced by w3iorsmd::w3iors(), w3seti(), w3updtmd::w3uwnd(), and wmupdtmd::wmupd1().

◆ wy0

real, dimension(:,:), pointer w3idatmd::wy0

Definition at line 243 of file w3idatmd.F90.

Referenced by w3seti(), w3updtmd::w3uwnd(), and wmupdtmd::wmupd1().

◆ wyn

real, dimension(:,:), pointer w3idatmd::wyn

Definition at line 243 of file w3idatmd.F90.

Referenced by w3iorsmd::w3iors(), w3seti(), w3updtmd::w3uwnd(), and wmupdtmd::wmupd1().

w3idatmd::gdn
real, pointer gdn
Definition: w3idatmd.F90:242
w3idatmd::flmdn
logical, pointer flmdn
Definition: w3idatmd.F90:263
w3idatmd::twn
integer, dimension(:), pointer twn
Definition: w3idatmd.F90:236
w3idatmd::inflags1
logical, dimension(:), pointer inflags1
Definition: w3idatmd.F90:260
w3idatmd::t1n
integer, dimension(:), pointer t1n
Definition: w3idatmd.F90:236
w3idatmd::dt0
real, dimension(:,:), pointer dt0
Definition: w3idatmd.F90:243
w3gdatmd::fswnd
logical, pointer fswnd
Definition: w3gdatmd.F90:1264
w3idatmd::flmth
logical, pointer flmth
Definition: w3idatmd.F90:263
w3idatmd::inflags2
logical, dimension(:), pointer inflags2
Definition: w3idatmd.F90:260
w3idatmd::wy0
real, dimension(:,:), pointer wy0
Definition: w3idatmd.F90:243
w3idatmd::inputs
type(input), dimension(:), allocatable, target inputs
Definition: w3idatmd.F90:232
w3idatmd::nidata
integer nidata
Definition: w3idatmd.F90:160
w3idatmd::t0n
integer, dimension(:), pointer t0n
Definition: w3idatmd.F90:236
w3idatmd::ti5
integer, dimension(:), pointer ti5
Definition: w3idatmd.F90:236
w3idatmd::gan
real, pointer gan
Definition: w3idatmd.F90:242
w3idatmd::flic4
logical, pointer flic4
Definition: w3idatmd.F90:264
w3idatmd::fllevresi
logical, pointer fllevresi
Definition: w3idatmd.F90:266
w3gdatmd::ny
integer, pointer ny
Definition: w3gdatmd.F90:1097
w3idatmd::tg0
integer, dimension(:), pointer tg0
Definition: w3idatmd.F90:236
w3idatmd::flcurtide
logical, pointer flcurtide
Definition: w3idatmd.F90:266
w3idatmd::wltide
real, dimension(:,:,:,:), pointer wltide
Definition: w3idatmd.F90:256
w3idatmd::flcur
logical, pointer flcur
Definition: w3idatmd.F90:261
w3idatmd::ti3
integer, dimension(:), pointer ti3
Definition: w3idatmd.F90:236
w3idatmd::tu0
integer, dimension(:), pointer tu0
Definition: w3idatmd.F90:236
w3idatmd::flcurresi
logical, pointer flcurresi
Definition: w3idatmd.F90:266
w3idatmd::flic5
logical, pointer flic5
Definition: w3idatmd.F90:264
w3gdatmd::w3setg
subroutine w3setg(IMOD, NDSE, NDST)
Definition: w3gdatmd.F90:2152
w3odatmd::ndse
integer, pointer ndse
Definition: w3odatmd.F90:456
w3idatmd::tzn
integer, dimension(:), pointer tzn
Definition: w3idatmd.F90:236
w3idatmd::fllev
logical, pointer fllev
Definition: w3idatmd.F90:261
w3idatmd::wxn
real, dimension(:,:), pointer wxn
Definition: w3idatmd.F90:243
w3idatmd::ntide
integer ntide
Definition: w3idatmd.F90:165
w3gdatmd::nsea
integer, pointer nsea
Definition: w3gdatmd.F90:1097
w3servmd
Definition: w3servmd.F90:3
w3idatmd::fllevtide
logical, pointer fllevtide
Definition: w3idatmd.F90:266
w3idatmd::cxtide
real, dimension(:,:,:,:), pointer cxtide
Definition: w3idatmd.F90:256
w3idatmd::flic2
logical, pointer flic2
Definition: w3idatmd.F90:264
w3idatmd::flagsc
logical, dimension(:), pointer flagsc
Definition: w3idatmd.F90:260
w3idatmd::flic3
logical, pointer flic3
Definition: w3idatmd.F90:264
w3idatmd::iidata
integer iidata
Definition: w3idatmd.F90:160
w3idatmd::dtn
real, dimension(:,:), pointer dtn
Definition: w3idatmd.F90:243
w3idatmd::gd0
real, pointer gd0
Definition: w3idatmd.F90:242
w3idatmd::ti4
integer, dimension(:), pointer ti4
Definition: w3idatmd.F90:236
w3idatmd::ti1
integer, dimension(:), pointer ti1
Definition: w3idatmd.F90:236
w3gdatmd::igrid
integer igrid
Definition: w3gdatmd.F90:618
w3idatmd::w3seti
subroutine w3seti(IMOD, NDSE, NDST)
Select one of the WAVEWATCH III grids / models.
Definition: w3idatmd.F90:819
w3idatmd::flmvs
logical, pointer flmvs
Definition: w3idatmd.F90:263
w3idatmd::ttn
integer, dimension(:), pointer ttn
Definition: w3idatmd.F90:236
w3idatmd::tdn
integer, dimension(:), pointer tdn
Definition: w3idatmd.F90:236
w3idatmd::tvn
integer, dimension(:), pointer tvn
Definition: w3idatmd.F90:236
w3idatmd::flwind
logical, pointer flwind
Definition: w3idatmd.F90:261
w3idatmd::ga0
real, pointer ga0
Definition: w3idatmd.F90:242
w3servmd::strace
subroutine strace(IENT, SNAME)
Definition: w3servmd.F90:148
w3idatmd::t2n
integer, dimension(:), pointer t2n
Definition: w3idatmd.F90:236
w3idatmd::wx0
real, dimension(:,:), pointer wx0
Definition: w3idatmd.F90:243
w3idatmd::tc0
integer, dimension(:), pointer tc0
Definition: w3idatmd.F90:236
w3idatmd::tw0
integer, dimension(:), pointer tw0
Definition: w3idatmd.F90:236
w3idatmd::tin
integer, dimension(:), pointer tin
Definition: w3idatmd.F90:236
w3idatmd::iinit
logical, pointer iinit
Definition: w3idatmd.F90:259
w3idatmd::fltaua
logical, pointer fltaua
Definition: w3idatmd.F90:261
w3idatmd::cytide
real, dimension(:,:,:,:), pointer cytide
Definition: w3idatmd.F90:256
w3gdatmd
Definition: w3gdatmd.F90:16
w3idatmd::flrhoa
logical, pointer flrhoa
Definition: w3idatmd.F90:261
w3idatmd::flice
logical, pointer flice
Definition: w3idatmd.F90:261
w3idatmd::ti2
integer, dimension(:), pointer ti2
Definition: w3idatmd.F90:236
w3idatmd::tln
integer, dimension(:), pointer tln
Definition: w3idatmd.F90:236
w3gdatmd::nx
integer, pointer nx
Definition: w3gdatmd.F90:1097
w3gdatmd::ngrids
integer ngrids
Definition: w3gdatmd.F90:618
w3idatmd::trn
integer, dimension(:), pointer trn
Definition: w3idatmd.F90:236
w3idatmd::tun
integer, dimension(:), pointer tun
Definition: w3idatmd.F90:236
w3idatmd::tr0
integer, dimension(:), pointer tr0
Definition: w3idatmd.F90:236
w3gdatmd::nauxgr
integer nauxgr
Definition: w3gdatmd.F90:618
w3idatmd::wyn
real, dimension(:,:), pointer wyn
Definition: w3idatmd.F90:243
w3idatmd::tfn
integer, dimension(:,:), pointer tfn
Definition: w3idatmd.F90:236
w3idatmd::tgn
integer, dimension(:), pointer tgn
Definition: w3idatmd.F90:236
w3idatmd::tcn
integer, dimension(:), pointer tcn
Definition: w3idatmd.F90:236
w3idatmd::flic1
logical, pointer flic1
Definition: w3idatmd.F90:264