WAVEWATCH III  beta 0.0.1
w3wdatmd Module Reference

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

Data Types

type  wdata
 

Functions/Subroutines

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

Variables

integer nwdata = -1
 
integer iwdata = -1
 
type(wdata), dimension(:), allocatable, target wdatas
 
integer, dimension(:), pointer time
 
integer, dimension(:), pointer tlev
 
integer, dimension(:), pointer tice
 
integer, dimension(:), pointer trho
 
integer, dimension(:), pointer tic1
 
integer, dimension(:), pointer tic5
 
integer, dimension(:), pointer time00
 
integer, dimension(:), pointer timeend
 
integer, dimension(:), pointer qi5tbeg
 
real, dimension(:), pointer qr5tim0
 
real, dimension(:, :), pointer qr5cvk0
 
real, dimension(:), pointer qr5tmix
 
complex, dimension(:, :), pointer qc5int0
 
real, dimension(:,:), pointer va
 
real, dimension(:), pointer wlv
 
real, dimension(:), pointer ice
 
real, dimension(:), pointer rhoair
 
real, dimension(:), pointer ust
 
real, dimension(:), pointer ustdir
 
real, dimension(:), pointer asf
 
real, dimension(:), pointer fpis
 
real, dimension(:), pointer berg
 
real, dimension(:), pointer iceh
 
real, dimension(:), pointer icef
 
real, dimension(:), pointer icedmax
 
real, dimension(:), pointer zeta_setup
 
real, dimension(:), pointer fx_zs
 
real, dimension(:), pointer fy_zs
 
real, dimension(:), pointer sxx_zs
 
real, dimension(:), pointer sxy_zs
 
real, dimension(:), pointer syy_zs
 
real, dimension(:,:), pointer vstot
 
real, dimension(:,:), pointer vdtot
 
real, dimension(:,:), pointer vaold
 
logical, dimension(:), pointer shavetot
 
logical, pointer dinit
 
logical, pointer fl_all
 

Detailed Description

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

The number of grids is taken from W3GDATMD, and needs to be set first with W3DIMG.

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

Function/Subroutine Documentation

◆ w3dimw()

subroutine w3wdatmd::w3dimw ( integer, intent(in)  IMOD,
integer, intent(in)  NDSE,
integer, intent(in)  NDST,
logical, intent(in), optional  F_ONLY 
)

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]F_ONLYFLag for initializing field arrays only.
Author
H. L. Tolman
Date
22-Mar-2021

Definition at line 343 of file w3wdatmd.F90.

343  !/
344  !/ +-----------------------------------+
345  !/ | WAVEWATCH III NOAA/NCEP |
346  !/ | H. L. Tolman |
347  !/ | FORTRAN 90 |
348  !/ | Last update : 22-Mar-2021 !
349  !/ +-----------------------------------+
350  !/
351  !/ 22-Oct-2004 : Origination. ( version 3.06 )
352  !/ 13-Jun-2006 : Allocate VA consistent with MPI ( version 3.09 )
353  !/ data types and initialize as needed.
354  !/ 05-Jul-2006 : Consolidate stress vector. ( version 3.09 )
355  !/ 04-Oct-2006 : Add filter to array pointers. ( version 3.10 )
356  !/ 14-Nov-2013 : Initialize UST and USTDIR. ( version 4.13 )
357  !/ 10-Dec-2014 : Add checks for allocate status ( version 5.04 )
358  !/ 22-Mar-2021 : Support for variable air density ( version 7.13 )
359  !/
360  ! 1. Purpose :
361  !
362  ! Initialize an individual data grid at the proper dimensions.
363  !
364  ! 2. Method :
365  !
366  ! Allocate directly into the structure array. Note that
367  ! this cannot be done through the pointer alias!
368  !
369  ! 3. Parameters :
370  !
371  ! Parameter list
372  ! ----------------------------------------------------------------
373  ! IMOD Int. I Model number to point to.
374  ! NDSE Int. I Error output unit number.
375  ! NDST Int. I Test output unit number.
376  ! F_ONLY L.O. I FLag for initializing field arrays only.
377  ! ----------------------------------------------------------------
378  !
379  ! 4. Subroutines used :
380  !
381  ! See module documentation.
382  !
383  ! 5. Called by :
384  !
385  ! Name Type Module Description
386  ! ----------------------------------------------------------------
387  ! W3IOGO Subr. W3IOGOMD Grid output IO routine.
388  ! W3IORS Subr. W3IORSMD Restart file IO routine.
389  ! WW3_SHEL Prog. N/A Main wave model driver.
390  ! WW3_STRT Prog. N/A Initial conditions program.
391  ! ----------------------------------------------------------------
392  !
393  ! 6. Error messages :
394  !
395  ! - Check on input parameters.
396  ! - Check on previous allocation.
397  !
398  ! 7. Remarks :
399  !
400  ! - W3SETW needs to be called after allocation to point to
401  ! proper allocated arrays.
402  !
403  ! 8. Structure :
404  !
405  ! See source code.
406  !
407  ! 9. Switches :
408  !
409  ! !/S Enable subroutine tracing.
410  ! !/T Enable test output
411  !
412  ! 10. Source code :
413  !
414  !/ ------------------------------------------------------------------- /
415  USE w3gdatmd, ONLY: ngrids, igrid, w3setg, nspec, nsea, nseal, grids
416  USE w3odatmd, ONLY: naproc, iaproc
417  USE w3servmd, ONLY: extcde
418  USE constants, ONLY : lpdlib, dair
420 #ifdef W3_NL5
421  USE w3gdatmd, ONLY: qi5nnz
422 #endif
423 #ifdef W3_PDLIB
424  use yownodepool, only: npa, np
425  use yowrankmodule, only : rank
426  USE w3gdatmd, ONLY: gtype, ungtype
427 #endif
428 #ifdef W3_S
429  USE w3servmd, ONLY: strace
430 #endif
431  !
432  IMPLICIT NONE
433  !
434  !/
435  !/ ------------------------------------------------------------------- /
436  !/ Parameter list
437  !/
438  INTEGER, INTENT(IN) :: IMOD, NDSE, NDST
439  LOGICAL, INTENT(IN), OPTIONAL :: F_ONLY
440  !/
441  !/ ------------------------------------------------------------------- /
442  !/ Local parameters
443  !/
444  INTEGER :: JGRID, NSEALM, NSEATM
445  INTEGER :: NSEAL_DUMMY, ISEA
446 #ifdef W3_PDLIB
447  INTEGER IRANK
448 #endif
449 #ifdef W3_S
450  INTEGER, SAVE :: IENT = 0
451 #endif
452  !/
453 #ifdef W3_S
454  CALL strace (ient, 'W3DIMW')
455 #endif
456 
457  !
458  ! -------------------------------------------------------------------- /
459  ! 1. Test input and module status
460  !
461  IF ( PRESENT(f_only) ) THEN
462  fl_all = .NOT. f_only
463  ELSE
464  fl_all = .true.
465  END IF
466  !
467  IF ( ngrids .EQ. -1 ) THEN
468  WRITE (ndse,1001)
469  CALL extcde (1)
470  END IF
471  !
472  IF ( imod.LT.1 .OR. imod.GT.nwdata ) THEN
473  WRITE (ndse,1002) imod, nwdata
474  CALL extcde (2)
475  END IF
476  !
477  IF ( wdatas(imod)%DINIT ) THEN
478  WRITE (ndse,1003)
479  CALL extcde (3)
480  END IF
481  !
482 #ifdef W3_T
483  WRITE (ndst,9000) imod
484 #endif
485  !
486  jgrid = igrid
487  IF ( jgrid .NE. imod ) CALL w3setg ( imod, ndse, ndst )
488  !
489  ! -------------------------------------------------------------------- /
490  ! 2. Allocate arrays
491  !
492  CALL set_up_nseal_nsealm(nseal_dummy, nsealm)
493  nseatm = nsealm * naproc
494  !
495  IF ( fl_all ) THEN
496  ALLOCATE ( wdatas(imod)%VA(nspec,0:nsealm), stat=istat ); wdatas(imod)%VA = 0.
497  check_alloc_status( istat )
498 #ifdef W3_PDLIB
499  ALLOCATE ( wdatas(imod)%SHAVETOT(nseal), stat=istat )
500 #endif
501 #ifdef W3_PDLIB
502  IF (.not. lsloc) THEN
503  ALLOCATE ( wdatas(imod)%VSTOT(nspec,nseal), stat=istat )
504 #endif
505 #ifdef W3_PDLIB
506  ALLOCATE ( wdatas(imod)%VDTOT(nspec,nseal), stat=istat )
507 #endif
508 #ifdef W3_PDLIB
509  ENDIF ! LSLOC
510  ALLOCATE ( wdatas(imod)%VAOLD(nspec,nseal), stat=istat )
511 #endif
512 #ifdef W3_PDLIB
513  IF (.not. lsloc) THEN
514  wdatas(imod)%VSTOT=0
515 #endif
516 #ifdef W3_PDLIB
517  wdatas(imod)%VDTOT=0
518 #endif
519 #ifdef W3_PDLIB
520  ENDIF ! LSLOC
521  wdatas(imod)%SHAVETOT=.false.
522 #endif
523  !
524  ! * Four arrays for NL5 (QL)
525  ! * AFAIK, the set up of QR5TIM0, QR5CVK0, QC5INT0 should be similar
526  ! * to VA, though I am not really clear about how FL_ALL works.
527  ! *
528 #ifdef W3_NL5
529  ALLOCATE ( wdatas(imod)%QR5TIM0(0:nsealm), &
530  wdatas(imod)%QR5CVK0(nspec, 0:nsealm), &
531  wdatas(imod)%QC5INT0(qi5nnz, 0:nsealm), &
532  wdatas(imod)%QR5TMIX(0:nsealm), stat=istat)
533  check_alloc_status( istat )
534 #endif
535  !
536  ! * Initialized NL5 arrays with zero (QL)
537 #ifdef W3_NL5
538  wdatas(imod)%QR5TIM0 = 0.0
539  wdatas(imod)%QR5CVK0 = 0.0
540  wdatas(imod)%QC5INT0 = (0.0, 0.0)
541  wdatas(imod)%QR5TMIX = 0.0
542 #endif
543  !
544 #ifdef W3_NL5
545  WRITE(*, *)
546  WRITE(*, '(A, I4, I12)') '⊚ → [WW3 WDAT]: IMOD & QI5NNZ: ', imod, qi5nnz
547  WRITE(*, *)
548 #endif
549  !
550  IF ( nseal .NE. nsealm ) THEN
551  DO isea=nseal+1,nsealm
552  wdatas(imod)%VA(:,isea) = 0.
553  !
554 #ifdef W3_NL5
555  wdatas(imod)%QR5TIM0(isea) = 0.0
556  wdatas(imod)%QR5CVK0(:,isea) = 0.0
557  wdatas(imod)%QC5INT0(:,isea) = (0.0, 0.0)
558  wdatas(imod)%QR5TMIX(isea) = 0.0
559 #endif
560  END DO
561  END IF
562  END IF
563  !
564  ! ICE, ICEH, ICEF must be defined from 0:NSEA
565  ALLOCATE ( wdatas(imod)%WLV(nsea), &
566  wdatas(imod)%ICE(0:nsea), &
567  wdatas(imod)%RHOAIR(nsea), &
568 #ifdef W3_SETUP
569  wdatas(imod)%ZETA_SETUP(nsea), &
570 #endif
571  wdatas(imod)%BERG(nsea), &
572  wdatas(imod)%ICEH(0:nsea), &
573  wdatas(imod)%ICEF(0:nsea), &
574  wdatas(imod)%ICEDMAX(nsea), &
575  wdatas(imod)%UST(0:nseatm), &
576  wdatas(imod)%USTDIR(0:nseatm), &
577  wdatas(imod)%ASF(nseatm), &
578  wdatas(imod)%FPIS(nseatm), stat=istat )
579  check_alloc_status( istat )
580 
581  wdatas(imod)%WLV (:) = 0.
582  wdatas(imod)%ICE (0:nsea) = 0.
583  wdatas(imod)%RHOAIR(:) = dair
584 #ifdef W3_SETUP
585  wdatas(imod)%ZETA_SETUP(:) = 0.
586 #endif
587  wdatas(imod)%BERG (:) = 0.
588  wdatas(imod)%ICEH (0:nsea) = grids(imod)%IICEHINIT
589  wdatas(imod)%ICEF (0:nsea) = 1000.
590  wdatas(imod)%ICEDMAX(:) = 1000.
591  wdatas(imod)%UST (0:nseatm) = 1.e-5
592  wdatas(imod)%USTDIR(0:nseatm) = 0.
593  wdatas(imod)%ASF (:) = 0.
594  wdatas(imod)%FPIS (:) = 0.
595  wdatas(imod)%DINIT = .true.
596  CALL w3setw ( imod, ndse, ndst )
597  !
598 #ifdef W3_T
599  WRITE (ndst,9003)
600 #endif
601  !
602  ! -------------------------------------------------------------------- /
603  ! 5. Restore previous grid setting if necessary
604  !
605  IF ( jgrid .NE. imod ) CALL w3setg ( jgrid, ndse, ndst )
606  !
607  RETURN
608  !
609  ! Formats
610  !
611 1001 FORMAT (/' *** ERROR W3DIMW : GRIDS NOT INITIALIZED *** '/ &
612  ' RUN W3NMOD FIRST '/)
613 1002 FORMAT (/' *** ERROR W3DIMW : ILLEGAL MODEL NUMBER *** '/ &
614  ' IMOD = ',i10/ &
615  ' NWDATA = ',i10/)
616 1003 FORMAT (/' *** ERROR W3DIMW : ARRAY(S) ALREADY ALLOCATED *** ')
617  !
618 #ifdef W3_T
619 9000 FORMAT (' TEST W3DIMW : MODEL ',i4,' DIM. AT ',2i5,i7)
620 #endif
621  !
622 #ifdef W3_T
623  WRITE (ndst,9001)
624 #endif
625  !
626  ! -------------------------------------------------------------------- /
627  ! 3. Point to allocated arrays
628  !
629  CALL w3setw ( imod, ndse, ndst )
630  !
631 #ifdef W3_T
632  WRITE (ndst,9002)
633 #endif
634  !
635  ! -------------------------------------------------------------------- /
636  ! 4. Update counters in grid
637 #ifdef W3_T
638 9001 FORMAT (' TEST W3DIMW : ARRAYS ALLOCATED')
639 9002 FORMAT (' TEST W3DIMW : POINTERS RESET')
640 9003 FORMAT (' TEST W3DIMW : DIMENSIONS STORED')
641 #endif
642  !/
643  !/ End of W3DIMW ----------------------------------------------------- /
644  !/

References constants::dair, w3servmd::extcde(), fl_all, w3gdatmd::grids, w3gdatmd::gtype, w3odatmd::iaproc, w3gdatmd::igrid, constants::lpdlib, w3parall::lsloc, w3odatmd::naproc, w3gdatmd::ngrids, yownodepool::np, yownodepool::npa, w3gdatmd::nsea, w3gdatmd::nseal, w3gdatmd::nspec, nwdata, w3gdatmd::qi5nnz, yowrankmodule::rank, w3parall::set_up_nseal_nsealm(), w3servmd::strace(), w3gdatmd::ungtype, w3gdatmd::w3setg(), w3setw(), and wdatas.

Referenced by w3grid_interp(), w3initmd::w3init(), w3iogomd::w3iogo(), w3iorsmd::w3iors(), and w3shel().

◆ w3ndat()

subroutine w3wdatmd::w3ndat ( 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
10-Dec-2014

Definition at line 210 of file w3wdatmd.F90.

210  !/
211  !/ +-----------------------------------+
212  !/ | WAVEWATCH III NOAA/NCEP |
213  !/ | H. L. Tolman |
214  !/ | FORTRAN 90 |
215  !/ | Last update : 10-Dec-2014 !
216  !/ +-----------------------------------+
217  !/
218  !/ 31-Mar-2004 : Origination. ( version 3.06 )
219  !/ 10-Dec-2014 : Add checks for allocate status ( version 5.04 )
220  !/
221  ! 1. Purpose :
222  !
223  ! Set up the number of grids to be used.
224  !
225  ! 2. Method :
226  !
227  ! Use data stored in NGRIDS in W3GDATMD.
228  !
229  ! 3. Parameters :
230  !
231  ! Parameter list
232  ! ----------------------------------------------------------------
233  ! NDSE Int. I Error output unit number.
234  ! NDST Int. I Test output unit number.
235  ! ----------------------------------------------------------------
236  !
237  ! 4. Subroutines used :
238  !
239  ! See module documentation.
240  !
241  ! 5. Called by :
242  !
243  ! Any program that uses this grid structure.
244  !
245  ! 6. Error messages :
246  !
247  ! - Error checks on previous setting of variable NGRIDS.
248  !
249  ! 7. Remarks :
250  !
251  ! 8. Structure :
252  !
253  ! 9. Switches :
254  !
255  ! !/S Enable subroutine tracing.
256  ! !/T Enable test output
257  !
258  ! 10. Source code :
259  !
260  !/ ------------------------------------------------------------------- /
261  USE w3gdatmd, ONLY: ngrids
262  USE w3servmd, ONLY: extcde
263 #ifdef W3_S
264  USE w3servmd, ONLY: strace
265 #endif
266  !
267  IMPLICIT NONE
268  !/
269  !/ ------------------------------------------------------------------- /
270  !/ Parameter list
271  !/
272  INTEGER, INTENT(IN) :: NDSE, NDST
273  !/
274  !/ ------------------------------------------------------------------- /
275  !/ Local parameters
276  !/
277  INTEGER :: I
278 #ifdef W3_S
279  INTEGER, SAVE :: IENT = 0
280 #endif
281  !/
282 #ifdef W3_S
283  CALL strace (ient, 'W3NDAT')
284 #endif
285  !
286  ! -------------------------------------------------------------------- /
287  ! 1. Test input and module status
288  !
289  IF ( ngrids .EQ. -1 ) THEN
290  WRITE (ndse,1001) ngrids
291  CALL extcde (1)
292  END IF
293  !
294  ! -------------------------------------------------------------------- /
295  ! 2. Set variable and allocate arrays
296  !
297  ALLOCATE ( wdatas(0:ngrids), stat=istat )
298  check_alloc_status( istat )
299  nwdata = ngrids
300  !
301  ! -------------------------------------------------------------------- /
302  ! 3. Initialize parameters
303  !
304  DO i=0, ngrids
305  wdatas(i)%DINIT = .false.
306  wdatas(i)%FL_ALL = .false.
307  END DO
308  !
309 #ifdef W3_T
310  WRITE (ndst,9000) ngrids
311 #endif
312  !
313  RETURN
314  !
315  ! Formats
316  !
317 1001 FORMAT (/' *** ERROR W3NDAT : NGRIDS NOT YET SET *** '/ &
318  ' NGRIDS = ',i10/ &
319  ' RUN W3NMOD FIRST'/)
320  !
321 #ifdef W3_T
322 9000 FORMAT (' TEST W3NDAT : SETTING UP FOR ',i4,' GRIDS')
323 #endif
324  !/
325  !/ End of W3NDAT ----------------------------------------------------- /
326  !/

References w3servmd::extcde(), w3gdatmd::ngrids, nwdata, w3servmd::strace(), and wdatas.

Referenced by gxoutf(), gxoutp(), w3bounc(), w3bound(), w3grib(), w3grid_interp(), w3ounf(), w3ounp(), w3outf(), w3outp(), w3shel(), w3strt(), w3uprstr(), wminitmd::wminit(), and wminitmd::wminitnml().

◆ w3setw()

subroutine w3wdatmd::w3setw ( 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 660 of file w3wdatmd.F90.

660  !/
661  !/ +-----------------------------------+
662  !/ | WAVEWATCH III NOAA/NCEP |
663  !/ | H. L. Tolman |
664  !/ | FORTRAN 90 |
665  !/ | Last update : 22-Mar-2021 !
666  !/ +-----------------------------------+
667  !/
668  !/ 31-Mar-2004 : Origination. ( version 3.06 )
669  !/ 05-Jul-2006 : Consolidate stress vector. ( version 3.09 )
670  !/ 04-Oct-2006 : Add filter to array pointers. ( version 3.10 )
671  !/ 22-Mar-2021 : Support for variable air density ( version 7.13 )
672  !/
673  ! 1. Purpose :
674  !
675  ! Select one of the WAVEWATCH III grids / models.
676  !
677  ! 2. Method :
678  !
679  ! Point pointers to the proper variables in the proper element of
680  ! the GRIDS array.
681  !
682  ! 3. Parameters :
683  !
684  ! Parameter list
685  ! ----------------------------------------------------------------
686  ! IMOD Int. I Model number to point to.
687  ! NDSE Int. I Error output unit number.
688  ! NDST Int. I Test output unit number.
689  ! ----------------------------------------------------------------
690  !
691  ! 4. Subroutines used :
692  !
693  ! See module documentation.
694  !
695  ! 5. Called by :
696  !
697  ! Many subroutines in the WAVEWATCH system.
698  !
699  ! 6. Error messages :
700  !
701  ! Checks on parameter list IMOD.
702  !
703  ! 7. Remarks :
704  !
705  ! 8. Structure :
706  !
707  ! 9. Switches :
708  !
709  ! !/S Enable subroutine tracing.
710  ! !/T Enable test output
711  !
712  ! 10. Source code :
713  !
714  !/ ------------------------------------------------------------------- /
715  USE w3servmd, ONLY: extcde
716 #ifdef W3_S
717  USE w3servmd, ONLY: strace
718 #endif
719  !
720  IMPLICIT NONE
721  !/
722  !/ ------------------------------------------------------------------- /
723  !/ Parameter list
724  !/
725  INTEGER, INTENT(IN) :: IMOD, NDSE, NDST
726  !/
727  !/ ------------------------------------------------------------------- /
728  !/ Local parameters
729  !/
730 #ifdef W3_S
731  INTEGER, SAVE :: IENT = 0
732 #endif
733  !/
734 #ifdef W3_S
735  CALL strace (ient, 'W3SETW')
736 #endif
737  !
738  ! -------------------------------------------------------------------- /
739  ! 1. Test input and module status
740  !
741  IF ( nwdata .EQ. -1 ) THEN
742  WRITE (ndse,1001)
743  CALL extcde (1)
744  END IF
745  !
746  IF ( imod.LT.0 .OR. imod.GT.nwdata ) THEN
747  WRITE (ndse,1002) imod, nwdata
748  CALL extcde (2)
749  END IF
750  !
751 #ifdef W3_T
752  WRITE (ndst,9000) imod
753 #endif
754  !
755  ! -------------------------------------------------------------------- /
756  ! 2. Set model numbers
757  !
758  iwdata = imod
759  !
760  ! -------------------------------------------------------------------- /
761  ! 3. Set pointers
762  !
763  time => wdatas(imod)%TIME
764 #ifdef W3_OASIS
765  time00 => wdatas(imod)%TIME00
766  timeend => wdatas(imod)%TIMEEND
767 #endif
768 #ifdef W3_NL5
769  qi5tbeg => wdatas(imod)%QI5TBEG
770 #endif
771  tlev => wdatas(imod)%TLEV
772  tice => wdatas(imod)%TICE
773  trho => wdatas(imod)%TRHO
774  tic1 => wdatas(imod)%TIC1
775  tic5 => wdatas(imod)%TIC5
776  dinit => wdatas(imod)%DINIT
777  fl_all => wdatas(imod)%FL_ALL
778  !
779  IF ( dinit ) THEN
780  IF ( fl_all ) THEN
781  va => wdatas(imod)%VA
782 #ifdef W3_NL5
783  qr5tim0 => wdatas(imod)%QR5TIM0
784  qr5cvk0 => wdatas(imod)%QR5CVK0
785  qc5int0 => wdatas(imod)%QC5INT0
786  qr5tmix => wdatas(imod)%QR5TMIX
787 #endif
788 #ifdef W3_PDLIB
789  shavetot => wdatas(imod)%SHAVETOT
790  vstot => wdatas(imod)%VSTOT
791  vdtot => wdatas(imod)%VDTOT
792  vaold => wdatas(imod)%VAOLD
793 #endif
794  END IF
795  wlv => wdatas(imod)%WLV
796  ice => wdatas(imod)%ICE
797  rhoair => wdatas(imod)%RHOAIR
798 #ifdef W3_SETUP
799  zeta_setup => wdatas(imod)%ZETA_SETUP
800  fx_zs => wdatas(imod)%FX_zs
801  fy_zs => wdatas(imod)%FY_zs
802  sxx_zs => wdatas(imod)%SXX_zs
803  sxy_zs => wdatas(imod)%SXY_zs
804  syy_zs => wdatas(imod)%SYY_zs
805 #endif
806  berg => wdatas(imod)%BERG
807  iceh => wdatas(imod)%ICEH
808  icef => wdatas(imod)%ICEF
809  icedmax=> wdatas(imod)%ICEDMAX
810  ust => wdatas(imod)%UST
811  ustdir => wdatas(imod)%USTDIR
812  asf => wdatas(imod)%ASF
813  fpis => wdatas(imod)%FPIS
814  END IF
815  !
816  RETURN
817  !
818  ! Formats
819  !
820 1001 FORMAT (/' *** ERROR W3SETW : GRIDS NOT INITIALIZED *** '/ &
821  ' RUN W3NMOD FIRST '/)
822 1002 FORMAT (/' *** ERROR W3SETW : ILLEGAL MODEL NUMBER *** '/ &
823  ' IMOD = ',i10/ &
824  ' NWDATA = ',i10/)
825  !
826 #ifdef W3_T
827 9000 FORMAT (' TEST W3SETW : MODEL ',i4,' SELECTED')
828 #endif
829  !/
830  !/ End of W3SETW ----------------------------------------------------- /
831  !/

References asf, berg, dinit, w3servmd::extcde(), fl_all, fpis, fx_zs, fy_zs, ice, icedmax, icef, iceh, iwdata, nwdata, qc5int0, qi5tbeg, qr5cvk0, qr5tim0, qr5tmix, rhoair, shavetot, w3servmd::strace(), sxx_zs, sxy_zs, syy_zs, tic1, tic5, tice, time, time00, timeend, tlev, trho, ust, ustdir, va, vaold, vdtot, vstot, wdatas, wlv, and zeta_setup.

Referenced by wmesmfmd::createexpmesh(), gxoutf(), gxoutp(), w3bounc(), w3bound(), w3dimw(), w3grib(), w3grid_interp(), w3initmd::w3init(), w3iobcmd::w3iobc(), w3iogomd::w3iogo(), w3iopomd::w3iopo(), w3iopomd::w3iopon(), w3iorsmd::w3iors(), w3iotrmd::w3iotr(), w3ounf(), w3ounp(), w3outf(), w3outp(), w3shel(), w3strt(), w3uprstr(), w3wavemd::w3wave(), wminitmd::wminit(), wminitmd::wminitnml(), wminiomd::wmiobg(), wminiomd::wmiobs(), wminiomd::wmioeg(), wminiomd::wmioes(), wminiomd::wmiohg(), wminiomd::wmiohs(), wmiopomd::wmiopo(), wmiopomd::wmiopp(), wmupdtmd::wmupdt(), and wmwavemd::wmwave().

Variable Documentation

◆ asf

◆ berg

real, dimension(:), pointer w3wdatmd::berg

◆ dinit

logical, pointer w3wdatmd::dinit

Definition at line 195 of file w3wdatmd.F90.

195  LOGICAL, POINTER :: DINIT, FL_ALL

Referenced by w3iogomd::w3iogo(), w3iorsmd::w3iors(), and w3setw().

◆ fl_all

logical, pointer w3wdatmd::fl_all

Definition at line 195 of file w3wdatmd.F90.

Referenced by w3dimw(), and w3setw().

◆ fpis

real, dimension(:), pointer w3wdatmd::fpis

◆ fx_zs

real, dimension(:), pointer w3wdatmd::fx_zs

Definition at line 187 of file w3wdatmd.F90.

Referenced by w3setw().

◆ fy_zs

real, dimension(:), pointer w3wdatmd::fy_zs

Definition at line 187 of file w3wdatmd.F90.

Referenced by w3setw().

◆ ice

◆ icedmax

real, dimension(:), pointer w3wdatmd::icedmax

Definition at line 183 of file w3wdatmd.F90.

Referenced by w3setw(), w3updtmd::w3uic5(), and w3wavemd::w3wave().

◆ icef

◆ iceh

real, dimension(:), pointer w3wdatmd::iceh

◆ iwdata

integer w3wdatmd::iwdata = -1

Definition at line 134 of file w3wdatmd.F90.

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

◆ nwdata

integer w3wdatmd::nwdata = -1

Definition at line 134 of file w3wdatmd.F90.

134  INTEGER :: NWDATA = -1, iwdata = -1

Referenced by w3dimw(), w3ndat(), and w3setw().

◆ qc5int0

complex, dimension(:, :), pointer w3wdatmd::qc5int0

Definition at line 181 of file w3wdatmd.F90.

181  COMPLEX, POINTER :: QC5INT0(:, :)

Referenced by w3setw(), and w3snl5md::w3snl5().

◆ qi5tbeg

integer, dimension(:), pointer w3wdatmd::qi5tbeg

Definition at line 179 of file w3wdatmd.F90.

179  INTEGER, POINTER :: QI5TBEG(:)

Referenced by w3setw(), w3shel(), and w3snl5md::w3snl5().

◆ qr5cvk0

real, dimension(:, :), pointer w3wdatmd::qr5cvk0

Definition at line 180 of file w3wdatmd.F90.

Referenced by w3setw(), and w3snl5md::w3snl5().

◆ qr5tim0

real, dimension(:), pointer w3wdatmd::qr5tim0

Definition at line 180 of file w3wdatmd.F90.

180  REAL, POINTER :: QR5TIM0(:), QR5CVK0(:, :), QR5TMIX(:)

Referenced by w3setw(), and w3snl5md::w3snl5().

◆ qr5tmix

real, dimension(:), pointer w3wdatmd::qr5tmix

Definition at line 180 of file w3wdatmd.F90.

Referenced by w3setw(), and w3snl5md::w3snl5().

◆ rhoair

real, dimension(:), pointer w3wdatmd::rhoair

◆ shavetot

logical, dimension(:), pointer w3wdatmd::shavetot

◆ sxx_zs

real, dimension(:), pointer w3wdatmd::sxx_zs

Definition at line 188 of file w3wdatmd.F90.

188  REAL, POINTER :: SXX_zs(:), SXY_zs(:), SYY_zs(:)

Referenced by w3setw().

◆ sxy_zs

real, dimension(:), pointer w3wdatmd::sxy_zs

Definition at line 188 of file w3wdatmd.F90.

Referenced by w3setw().

◆ syy_zs

real, dimension(:), pointer w3wdatmd::syy_zs

Definition at line 188 of file w3wdatmd.F90.

Referenced by w3setw().

◆ tic1

integer, dimension(:), pointer w3wdatmd::tic1

Definition at line 172 of file w3wdatmd.F90.

Referenced by w3iorsmd::w3iors(), w3setw(), w3updtmd::w3uic1(), and w3wavemd::w3wave().

◆ tic5

integer, dimension(:), pointer w3wdatmd::tic5

Definition at line 172 of file w3wdatmd.F90.

Referenced by w3iorsmd::w3iors(), w3setw(), w3updtmd::w3uic5(), and w3wavemd::w3wave().

◆ tice

integer, dimension(:), pointer w3wdatmd::tice

◆ time

integer, dimension(:), pointer w3wdatmd::time

Definition at line 172 of file w3wdatmd.F90.

172  INTEGER, POINTER :: TIME(:), TLEV(:), TICE(:), TRHO(:), &
173  TIC1(:), TIC5(:)

Referenced by pdlib_w3profsmd::apply_boundary_condition(), pdlib_w3profsmd::apply_boundary_condition_va(), gxoutf(), gxoutp(), pdlib_w3profsmd::pdlib_jacobi_gauss_seidel_block(), pdlib_w3profsmd::pdlib_w3xypfsfct2(), pdlib_w3profsmd::pdlib_w3xypfsn2(), pdlib_w3profsmd::pdlib_w3xypfspsi2(), pdlib_w3profsmd::pdlib_w3xypug(), wmesmfmd::setservices(), w3bullmd::w3bull(), w3pro3md::w3cflxy(), w3grib(), w3grid_interp(), w3initmd::w3init(), w3iogomd::w3iogo(), w3iopomd::w3iopo(), w3iopomd::w3iopon(), w3iopomd::w3iopon_read(), w3iopomd::w3iopon_write(), w3iorsmd::w3iors(), w3iosfmd::w3iosf(), w3iotrmd::w3iotr(), w3ounf(), w3ounp(), w3outf(), w3outp(), w3psmcmd::w3psmc(), w3src6md::w3sds6(), w3setw(), w3shel(), w3srcemd::w3srce(), w3updtmd::w3ucur(), w3updtmd::w3uic1(), w3updtmd::w3uic5(), w3updtmd::w3uice(), w3updtmd::w3ulev(), w3uprstr(), w3updtmd::w3urho(), w3updtmd::w3utau(), w3updtmd::w3uwnd(), w3wavemd::w3wave(), w3wdasmd::w3wdas(), w3pro1md::w3xyp1(), w3pro2md::w3xyp2(), w3pro3md::w3xyp3(), w3profsmd::w3xypfsfct2(), w3profsmd::w3xypfsn2(), w3profsmd::w3xypfsnimp(), w3profsmd::w3xypfspsi2(), w3profsmd::w3xypug(), wminitmd::wminit(), wminitmd::wminitnml(), wminiomd::wmiobg(), wminiomd::wmiobs(), wminiomd::wmioeg(), wminiomd::wmioes(), wminiomd::wmiohg(), wminiomd::wmiohs(), wmiopomd::wmiopo(), wmupdtmd::wmupd1(), wmupdtmd::wmupd2(), wmupdtmd::wmupdt(), and wmwavemd::wmwave().

◆ time00

integer, dimension(:), pointer w3wdatmd::time00

Definition at line 175 of file w3wdatmd.F90.

175  INTEGER, POINTER :: TIME00(:)

Referenced by w3setw(), w3shel(), and w3wavemd::w3wave().

◆ timeend

integer, dimension(:), pointer w3wdatmd::timeend

Definition at line 176 of file w3wdatmd.F90.

176  INTEGER, POINTER :: TIMEEND(:)

Referenced by w3setw(), w3shel(), and w3wavemd::w3wave().

◆ tlev

integer, dimension(:), pointer w3wdatmd::tlev

◆ trho

integer, dimension(:), pointer w3wdatmd::trho

Definition at line 172 of file w3wdatmd.F90.

Referenced by w3initmd::w3init(), w3iorsmd::w3iors(), w3setw(), and w3updtmd::w3urho().

◆ ust

◆ ustdir

◆ va

real, dimension(:,:), pointer w3wdatmd::va

◆ vaold

◆ vdtot

◆ vstot

real, dimension(:,:), pointer w3wdatmd::vstot

◆ wdatas

type(wdata), dimension(:), allocatable, target w3wdatmd::wdatas

Definition at line 168 of file w3wdatmd.F90.

168  TYPE(WDATA), TARGET, ALLOCATABLE :: WDATAS(:)

Referenced by w3dimw(), w3exgi(), w3grid_interp(), w3ndat(), w3setw(), and wminiomd::wmiohg().

◆ wlv

real, dimension(:), pointer w3wdatmd::wlv

◆ zeta_setup

real, dimension(:), pointer w3wdatmd::zeta_setup

Definition at line 187 of file w3wdatmd.F90.

187  REAL, POINTER :: ZETA_SETUP(:), FX_zs(:), FY_zs(:)

Referenced by w3wavset::fd_wave_setup_computation(), w3wavset::trig_wave_setup_computation(), w3initmd::w3init(), w3iogomd::w3iogo(), w3iopomd::w3iope(), w3ounf(), w3setw(), and w3updtmd::w3ulev().

w3wdatmd::qr5tim0
real, dimension(:), pointer qr5tim0
Definition: w3wdatmd.F90:180
w3wdatmd::iwdata
integer iwdata
Definition: w3wdatmd.F90:134
w3gdatmd::nseal
integer, pointer nseal
Definition: w3gdatmd.F90:1097
w3wdatmd::iceh
real, dimension(:), pointer iceh
Definition: w3wdatmd.F90:183
w3wdatmd::fpis
real, dimension(:), pointer fpis
Definition: w3wdatmd.F90:183
w3adatmd::nsealm
integer, pointer nsealm
Definition: w3adatmd.F90:686
w3wdatmd::shavetot
logical, dimension(:), pointer shavetot
Definition: w3wdatmd.F90:193
w3gdatmd::nspec
integer, pointer nspec
Definition: w3gdatmd.F90:1230
yowrankmodule::rank
type(t_rank), dimension(:), allocatable, public rank
Provides access to some information of all threads e.g.
Definition: yowrankModule.F90:68
constants::dair
real, parameter dair
DAIR Density of air (kg/m3).
Definition: constants.F90:63
w3gdatmd::ungtype
integer, parameter ungtype
Definition: w3gdatmd.F90:626
w3parall::set_up_nseal_nsealm
subroutine set_up_nseal_nsealm(NSEALout, NSEALMout)
Setup NSEAL, NSEALM in context of PDLIB.
Definition: w3parall.F90:1040
w3parall::lsloc
logical, parameter lsloc
Definition: w3parall.F90:89
w3wdatmd::icef
real, dimension(:), pointer icef
Definition: w3wdatmd.F90:183
w3wdatmd::qc5int0
complex, dimension(:, :), pointer qc5int0
Definition: w3wdatmd.F90:181
w3wdatmd::wlv
real, dimension(:), pointer wlv
Definition: w3wdatmd.F90:183
w3odatmd::iaproc
integer, pointer iaproc
Definition: w3odatmd.F90:457
w3wdatmd::wdatas
type(wdata), dimension(:), allocatable, target wdatas
Definition: w3wdatmd.F90:168
w3wdatmd::time
integer, dimension(:), pointer time
Definition: w3wdatmd.F90:172
w3gdatmd::grids
type(grid), dimension(:), allocatable, target grids
Definition: w3gdatmd.F90:1088
w3adatmd::fl_all
logical, pointer fl_all
Definition: w3adatmd.F90:688
yownodepool::npa
integer, public npa
number of ghost + resident nodes this partition holds
Definition: yownodepool.F90:99
w3wdatmd::va
real, dimension(:,:), pointer va
Definition: w3wdatmd.F90:183
w3wdatmd::trho
integer, dimension(:), pointer trho
Definition: w3wdatmd.F90:172
w3wdatmd::tlev
integer, dimension(:), pointer tlev
Definition: w3wdatmd.F90:172
w3gdatmd::w3setg
subroutine w3setg(IMOD, NDSE, NDST)
Definition: w3gdatmd.F90:2152
w3wdatmd::fx_zs
real, dimension(:), pointer fx_zs
Definition: w3wdatmd.F90:187
w3wdatmd::berg
real, dimension(:), pointer berg
Definition: w3wdatmd.F90:183
w3wdatmd::qr5tmix
real, dimension(:), pointer qr5tmix
Definition: w3wdatmd.F90:180
w3wdatmd::vdtot
real, dimension(:,:), pointer vdtot
Definition: w3wdatmd.F90:191
yownodepool
Has data that belong to nodes.
Definition: yownodepool.F90:39
constants::lpdlib
logical lpdlib
LPDLIB Logical for using the PDLIB library.
Definition: constants.F90:101
w3gdatmd::nsea
integer, pointer nsea
Definition: w3gdatmd.F90:1097
w3servmd
Definition: w3servmd.F90:3
w3wdatmd::vstot
real, dimension(:,:), pointer vstot
Definition: w3wdatmd.F90:191
w3wdatmd::tic1
integer, dimension(:), pointer tic1
Definition: w3wdatmd.F90:172
yowrankmodule
Provides access to some information of all threads e.g.
Definition: yowrankModule.F90:44
w3wdatmd::w3setw
subroutine w3setw(IMOD, NDSE, NDST)
Select one of the WAVEWATCH III grids / models.
Definition: w3wdatmd.F90:660
w3wdatmd::qr5cvk0
real, dimension(:, :), pointer qr5cvk0
Definition: w3wdatmd.F90:180
w3gdatmd::qi5nnz
integer(kind=8), pointer qi5nnz
Definition: w3gdatmd.F90:1373
w3odatmd
Definition: w3odatmd.F90:3
w3odatmd::naproc
integer, pointer naproc
Definition: w3odatmd.F90:457
w3wdatmd::sxy_zs
real, dimension(:), pointer sxy_zs
Definition: w3wdatmd.F90:188
yownodepool::np
integer, public np
number of nodes, local
Definition: yownodepool.F90:93
w3gdatmd::igrid
integer igrid
Definition: w3gdatmd.F90:618
w3wdatmd::timeend
integer, dimension(:), pointer timeend
Definition: w3wdatmd.F90:176
w3wdatmd::qi5tbeg
integer, dimension(:), pointer qi5tbeg
Definition: w3wdatmd.F90:179
w3servmd::strace
subroutine strace(IENT, SNAME)
Definition: w3servmd.F90:148
w3wdatmd::icedmax
real, dimension(:), pointer icedmax
Definition: w3wdatmd.F90:183
w3gdatmd::gtype
integer, pointer gtype
Definition: w3gdatmd.F90:1094
w3wdatmd::zeta_setup
real, dimension(:), pointer zeta_setup
Definition: w3wdatmd.F90:187
w3wdatmd::ice
real, dimension(:), pointer ice
Definition: w3wdatmd.F90:183
w3wdatmd::time00
integer, dimension(:), pointer time00
Definition: w3wdatmd.F90:175
w3wdatmd::ust
real, dimension(:), pointer ust
Definition: w3wdatmd.F90:183
w3wdatmd::nwdata
integer nwdata
Definition: w3wdatmd.F90:134
constants
Define some much-used constants for global use (all defined as PARAMETER).
Definition: constants.F90:20
w3gdatmd
Definition: w3gdatmd.F90:16
w3wdatmd::fy_zs
real, dimension(:), pointer fy_zs
Definition: w3wdatmd.F90:187
w3servmd::extcde
subroutine extcde(IEXIT, UNIT, MSG, FILE, LINE, COMM)
Definition: w3servmd.F90:736
w3wdatmd::rhoair
real, dimension(:), pointer rhoair
Definition: w3wdatmd.F90:183
w3wdatmd::ustdir
real, dimension(:), pointer ustdir
Definition: w3wdatmd.F90:183
w3wdatmd::syy_zs
real, dimension(:), pointer syy_zs
Definition: w3wdatmd.F90:188
w3wdatmd::tic5
integer, dimension(:), pointer tic5
Definition: w3wdatmd.F90:172
w3gdatmd::ngrids
integer ngrids
Definition: w3gdatmd.F90:618
w3parall
Parallel routines for implicit solver.
Definition: w3parall.F90:22
w3wdatmd::vaold
real, dimension(:,:), pointer vaold
Definition: w3wdatmd.F90:192
w3wdatmd::tice
integer, dimension(:), pointer tice
Definition: w3wdatmd.F90:172
w3wdatmd::sxx_zs
real, dimension(:), pointer sxx_zs
Definition: w3wdatmd.F90:188
w3wdatmd::dinit
logical, pointer dinit
Definition: w3wdatmd.F90:195
w3wdatmd::asf
real, dimension(:), pointer asf
Definition: w3wdatmd.F90:183