WAVEWATCH III  beta 0.0.1
wmupdtmd Module Reference

Update model input at the driver level of the multi-grid version of WAVEWATCH III. More...

Functions/Subroutines

subroutine wmupdt (IMOD, TDATA)
 Update inputs for selected wave model grid. More...
 
subroutine wmupd1 (IMOD, IDSTR, J, IERR)
 Update selected input using native input files. More...
 
subroutine wmupd2 (IMOD, J, JMOD, IERR)
 Update selected input using input grids. More...
 
subroutine wmupdv (IMOD, VX, VY, JMOD, VXI, VYI, UNDEF, CONSTP)
 Interpolate vector field from input grid to model grid. More...
 
subroutine wmupds (IMOD, FD, JMOD, FDI, UNDEF)
 Interpolate scalar field from input grid to model grid. More...
 
integer function xycurvisearch (LENGTH, ARRAY, VALUE, DELTA)
 Search the location of a point(XC,YC) in a regular grid. More...
 
real function interpolate (X_LEN, XARRAY, Y_LEN, YARRAY, FUNC, X, Y, DELTA)
 Perform interpolation from regular to curvilinear grid for a scalar field. More...
 
subroutine interpolate2d (X_LEN, XARRAY, Y_LEN, YARRAY, FUNC1, FUNC2, X, Y, DELTA, VAL1, VAL2)
 Perform interpolation from regular to curvilinear grid for a vector field. More...
 
real function averaging (X_LEN, XARRAY, Y_LEN, YARRAY, FUNC, X, Y, NPX, NPY)
 This function uses averaging to estimate the value of a function f at point (x,y). More...
 

Variables

integer, parameter swpmax = 5
 SWPMAX. More...
 

Detailed Description

Update model input at the driver level of the multi-grid version of WAVEWATCH III.

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

Function/Subroutine Documentation

◆ averaging()

real function wmupdtmd::averaging ( integer  X_LEN,
real, dimension(x_len)  XARRAY,
integer  Y_LEN,
real, dimension(y_len)  YARRAY,
real, dimension(x_len, y_len)  FUNC,
real  X,
real  Y,
integer  NPX,
integer  NPY 
)

This function uses averaging to estimate the value of a function f at point (x,y).

f is assumed to be on a regular grid, with the grid x values specified by xarray with dimension x_len and the grid y values specified by yarray with dimension y_len, the number of point to be taken into account in x and y.

Parameters
X_LEN
XARRAY
Y_LEN
YARRAY
FUNC
X
Y
NPX
NPY
Returns
AVERAGING
Author
H. L. Tolman
Date
25-Jul-2019

Definition at line 3028 of file wmupdtmd.F90.

3028  !/
3029  !/ +-----------------------------------+
3030  !/ | WAVEWATCH III NOAA/NCEP |
3031  !/ | H. L. Tolman |
3032  !/ | FORTRAN 90 |
3033  !/ | Last update : 25-July-2019 |
3034  !/ +-----------------------------------+
3035  !/
3036  !/ (R. Padilla-Hernandez, EMC/NOAA)
3037  !/
3038  !/ 29-July-2019 : ( version 7.13 )
3039  !/
3040 
3041  ! THIS FUNCTION USES AVERAGING TO ESTIMATE THE VALUE
3042  ! OF A FUNCTION F AT POINT (X,Y)
3043  ! F IS ASSUMED TO BE ON A REGULAR GRID, WITH THE GRID X VALUES SPECIFIED
3044  ! BY XARRAY WITH DIMENSION X_LEN
3045  ! AND THE GRID Y VALUES SPECIFIED BY YARRAY WITH DIMENSION Y_LEN
3046  ! ININI AND INEND, THE NUMBER OF POINT TO BE TAKEN INTO ACCOUNT
3047  ! IN X AND Y
3048 
3049  !IMPLICIT NONE
3050  INTEGER X_LEN, Y_LEN, INXEND, INYEND, NPX,NPY
3051  REAL, DIMENSION(X_LEN) :: XARRAY
3052  REAL, DIMENSION(Y_LEN) :: YARRAY
3053  REAL, DIMENSION(X_LEN, Y_LEN) :: FUNC
3054  REAL :: X,Y
3055 
3056  REAL :: X1, X2, Y1, Y2, SUM
3057  INTEGER :: INX,INY, INITIALX, INITIALY
3058  INTEGER :: INFINX, INFINY,ICOUNT,I,J
3059 
3060  inx = xycurvisearch(x_len, xarray, x)
3061  iny = xycurvisearch(y_len, yarray, y)
3062 
3063  x1 = xarray(inx)
3064  !X2 = XARRAY(INX+1)
3065 
3066  y1 = yarray(iny)
3067  !Y2 = YARRAY(INY+1)
3068 
3069  inxend=npx+1
3070  inyend=npy+1
3071  ! LETS FIX THE INITIAL INDEX =1 NEGATIVE INDEXES IN LONG
3072  IF (inx-npx .LT. 1) THEN
3073  initialx=1
3074  ELSE
3075  initialx=inx-npx
3076  END IF
3077  ! LETS FIX THE FINAL INDEX =NX IF LOOKING FOR INDEXES > NX
3078  IF (inx+inxend .GT. x_len) THEN
3079  infinx=x_len
3080  ELSE
3081  infinx=inx+inxend
3082  END IF
3083  ! LETS FIX THE INITIAL INDEX =1 FOR NEGATIVE INDEXES FOR LAT
3084  IF (iny-npy .LT. 1) THEN
3085  initialy=1
3086  ELSE
3087  initialy=iny-npy
3088  END IF
3089  ! LETS FIX THE FINAL INDEX =NX IF LOOKING FOR INDEXES > NX
3090  IF (iny+inyend .GT. y_len) THEN
3091  infiny=y_len
3092  ELSE
3093  infiny=iny+inyend
3094  END IF
3095 
3096 
3097  sum=0.0
3098  icount=0
3099  DO j=initialy,infiny
3100  DO i=initialx,infinx
3101  icount=icount+1
3102  sum=sum+func(i,j)
3103  END DO
3104  END DO
3105  averaging=sum/real(icount)
3106 

References xycurvisearch().

◆ interpolate()

real function wmupdtmd::interpolate ( integer, intent(in)  X_LEN,
real, dimension(x_len), intent(in)  XARRAY,
integer, intent(in)  Y_LEN,
real, dimension(y_len), intent(in)  YARRAY,
real, dimension(x_len, y_len), intent(in)  FUNC,
real, intent(in)  X,
real, intent(in)  Y,
real, intent(in), optional  DELTA 
)

Perform interpolation from regular to curvilinear grid for a scalar field.

This function uses bilinear interpolation to estimate the value of a function f at point (x,y). f is assumed to be on a regular grid, with the grid x values specified by xarray with dimension x_len and the grid y values specified by yarray with dimension y_len.

Parameters
X_LENDimension in X
XARRAY1D array for Longitudes
Y_LENDimension in Y
YARRAY1D array for Latitudes
FUNC1D Field
XLong for point in the curv grid
YLat for point in the curv grid
DELTAThreshold to determine if two values are equal
Returns
INTERPOLATE
Author
H. L. Tolman
Date
25-Jul-2019

Definition at line 2744 of file wmupdtmd.F90.

2744  !/
2745  !/ +-----------------------------------+
2746  !/ | WAVEWATCH III NOAA/NCEP |
2747  !/ | H. L. Tolman |
2748  !/ | FORTRAN 90 |
2749  !/ | Last update : 25-July-2019 |
2750  !/ +-----------------------------------+
2751  !/
2752  !/ (R. Padilla-Hernandez, EMC/NOAA)
2753  !/
2754  !/ 29-July-2019 : ( version 7.13 )
2755  !/
2756  ! 1. Purpose :
2757  !
2758  ! Perform interpolation from regular to curvilinear grid for a
2759  ! scalar field. THIS FUNCTION USES BILINEAR INTERPOLATION TO
2760  ! ESTIMATE THE VALUE OF A FUNCTION F AT POINT (X,Y) F IS ASSUMED
2761  ! TO BE ON A REGULAR GRID, WITH THE GRID X VALUES SPECIFIED BY
2762  ! XARRAY WITH DIMENSION X_LEN AND THE GRID Y VALUES SPECIFIED BY
2763  ! YARRAY WITH DIMENSION Y_LEN
2764  !
2765  ! 2. Parameters :
2766  !
2767  ! Parameter list
2768  ! ----------------------------------------------------------------
2769  ! X_LEN Int. Dimension in X
2770  ! XARRAY Int. 1D array for Longitudes
2771  ! Y_LEN Int. Dimension in Y
2772  ! YARRAY Int. 1D array for Latitudes
2773  ! FUNC Int. 1D Field
2774  ! X,Y Real Long-Lat for point in the curv grid
2775  ! DELTA Real Threshold to determine if two values are equal
2776  ! ----------------------------------------------------------------
2777  !
2778  ! Internal parameters
2779  ! ----------------------------------------------------------------
2780  ! INX Int. Index in X on the rectiliniear grid that is
2781  ! closest to, but less than, the given value for a
2782  ! point in the curvilinear grid.
2783  ! JNX Int. Idem INX for for Y.
2784  ! X1,Y1 Real (Long, Lat) left-bottom corner for the square in
2785  ! regular grid, where the given value for the point
2786  ! in the curvilinear grid lies
2787  ! X2,Y2 Real (Long, Lat) right-upper corner for the square in
2788  ! regular grid, where the given value for the point
2789  ! in the curvilinear grid lies
2790  ! ----------------------------------------------------------------
2791  !
2792  ! 3. Subroutines used :
2793  !
2794  ! Name Type Module Description
2795  ! ----------------------------------------------------------------
2796  ! XYCURVISEARCH Func. wmupdtmd Look for indexes in 1D array.
2797  ! ----------------------------------------------------------------
2798  !
2799  ! 4. Called by :
2800  !
2801  ! Main program in which it is contained.
2802  !
2803  ! 5. Error messages :
2804  !
2805  ! None.
2806  !
2807  ! 6. Remarks :
2808  !
2809  ! -
2810  !
2811  ! 7. Structure :
2812  !
2813  ! See source code.
2814  !
2815  ! 8. Switches :
2816  !
2817  ! -
2818  !
2819  ! 9. Source code :
2820  ! THIS FUNCTION USES BILINEAR INTERPOLATION TO ESTIMATE THE VALUE
2821  ! OF A FUNCTION F AT POINT (X,Y)
2822  ! F IS ASSUMED TO BE ON A REGULAR GRID, WITH THE GRID X VALUES SPECIFIED
2823  ! BY XARRAY WITH DIMENSION X_LEN
2824  ! AND THE GRID Y VALUES SPECIFIED BY YARRAY WITH DIMENSION Y_LEN
2825  !
2826  INTEGER, INTENT(IN) :: X_LEN, Y_LEN
2827  REAL, DIMENSION(X_LEN), INTENT(IN) :: XARRAY
2828  REAL, DIMENSION(Y_LEN), INTENT(IN) :: YARRAY
2829  REAL, DIMENSION(X_LEN, Y_LEN), INTENT(IN) :: FUNC
2830  REAL, INTENT(IN) :: X,Y
2831  REAL, INTENT(IN), OPTIONAL :: DELTA
2832  REAL :: DENOM, X1, X2, Y1, Y2
2833  INTEGER :: INX,JNX
2834 
2835  inx = xycurvisearch(x_len, xarray, x, delta)
2836  jnx = xycurvisearch(y_len, yarray, y, delta)
2837  !
2838  IF (inx .GE. x_len) THEN
2839  inx=inx-1
2840  END IF
2841  IF (jnx .GE. y_len) THEN
2842  jnx=jnx-1
2843  END IF
2844  !
2845  x1 = xarray(inx)
2846  x2 = xarray(inx+1)
2847  y1 = yarray(jnx)
2848  y2 = yarray(jnx+1)
2849  !
2850  denom = (x2 - x1)*(y2 - y1)
2851  !
2852  interpolate = (func(inx,jnx)*(x2-x)*(y2-y) + &
2853  func(inx+1,jnx)*(x-x1)*(y2-y) + &
2854  func(inx,jnx+1)*(x2-x)*(y-y1)+ &
2855  func(inx+1, jnx+1)*(x-x1)*(y-y1))/denom
2856  !

References xycurvisearch().

Referenced by wmupds().

◆ interpolate2d()

subroutine wmupdtmd::interpolate2d ( integer, intent(in)  X_LEN,
real, dimension(x_len), intent(in)  XARRAY,
integer, intent(in)  Y_LEN,
real, dimension(y_len), intent(in)  YARRAY,
real, dimension(x_len, y_len), intent(in)  FUNC1,
real, dimension(x_len, y_len), intent(in)  FUNC2,
real, intent(in)  X,
real, intent(in)  Y,
real, intent(in), optional  DELTA,
real, intent(out)  VAL1,
real, intent(out)  VAL2 
)

Perform interpolation from regular to curvilinear grid for a vector field.

This function uses bilinear interpolation to estimate the value of a function f at point (x,y). f is assumed to be on a regular grid, with the grid x values specified by xarray with dimension x_len and the grid y values specified by yarray with dimension y_len.

Parameters
[in]X_LENDimension in X
[in]XARRAY1D array for Longitudes
[in]Y_LENDimension in Y
[in]YARRAY1D array for Latitudes
[in]FUNC1First component of the 2D array
[in]FUNC2Second component of the 2D array
[in]XLong for point the curv grid
[in]YLat for point the curv grid
[in]DELTAThreshold to determine if two values are equal
[out]VAL1Interpolated value at a point in curv grid
[out]VAL2Interpolated value at a point in curv grid
Author
H. L. Tolman
Date
25-Jul-2019

Definition at line 2886 of file wmupdtmd.F90.

2886  !/
2887  !/ +-----------------------------------+
2888  !/ | WAVEWATCH III NOAA/NCEP |
2889  !/ | H. L. Tolman |
2890  !/ | FORTRAN 90 |
2891  !/ | Last update : 25-July-2019 |
2892  !/ +-----------------------------------+
2893  !/
2894  !/ (R. Padilla-Hernandez, EMC/NOAA)
2895  !/
2896  !/ 29-July-2019 : ( version 7.13 )
2897  !/
2898  ! 1. Purpose :
2899  !
2900  ! Perform interpolation from regular to curvilinear grid for a
2901  ! Vector field. THIS FUNCTION USES BILINEAR INTERPOLATION TO
2902  ! ESTIMATE THE VALUE OF A FUNCTION F AT POINT (X,Y) F IS ASSUMED
2903  ! TO BE ON A REGULAR GRID, WITH THE GRID X VALUES SPECIFIED BY
2904  ! XARRAY WITH DIMENSION X_LEN AND THE GRID Y VALUES SPECIFIED BY
2905  ! YARRAY WITH DIMENSION Y_LEN
2906  !
2907  ! 3. Parameters :
2908  !
2909  ! Parameter list
2910  ! ----------------------------------------------------------------
2911  ! X_LEN Int. Dimension in X
2912  ! XARRAY Int. 1D array for Longitudes
2913  ! Y_LEN Int. Dimension in Y
2914  ! YARRAY Int. 1D array for Latitudes
2915  ! FUNC1 Int. First componen of the 2D array
2916  ! FUNC2 Int. Second component of the 2D array
2917  ! X,Y Real Long-Lat for point in the curv grid
2918  ! DELTA Real Threshold to determine if two values are equal
2919  ! VAL1,VAL2 Real Interpolated values at a point in curvi grid
2920  ! ----------------------------------------------------------------
2921  !
2922  ! Internal parameters
2923  ! ----------------------------------------------------------------
2924  ! INX Int. Index in X on the rectiliniear grid that is
2925  ! closest to, but less than, the given value for a
2926  ! point in the curvilinear grid.
2927  ! JNX Int. Idem INX for for Y.
2928  ! X1,Y1 Real (Long, Lat) left-bottom corner for the square in
2929  ! regular grid, where the given value for the point
2930  ! in the curvilinear grid lies
2931  ! X2,Y2 Real (Long, Lat) right-upper corner for the square in
2932  ! regular grid, where the given value for the point
2933  ! in the curvilinear grid lies
2934  ! ----------------------------------------------------------------
2935  !
2936  ! 4. Subroutines used :
2937  !
2938  ! Name Type Module Description
2939  ! ----------------------------------------------------------------
2940  ! XYCURVISEARCH Func. wmupdtmd Look for indexes in 1D array.
2941  ! ----------------------------------------------------------------
2942  !
2943  ! 5. Called by :
2944  !
2945  ! Main program in which it is contained.
2946  !
2947  ! 6. Error messages :
2948  !
2949  ! None.
2950  !
2951  ! 7. Remarks :
2952  !
2953  ! -
2954  !
2955  ! 8. Structure :
2956  !
2957  ! See source code.
2958  !
2959  ! 9. Switches :
2960  !
2961  ! -
2962  !
2963  ! 10. Source code :
2964  !
2965  INTEGER, INTENT(IN) :: X_LEN, Y_LEN
2966  REAL, DIMENSION(X_LEN), INTENT(IN) :: XARRAY
2967  REAL, DIMENSION(Y_LEN), INTENT(IN) :: YARRAY
2968  REAL, DIMENSION(X_LEN, Y_LEN), INTENT(IN) :: FUNC1, FUNC2
2969  REAL, INTENT(IN) :: X,Y
2970  REAL, INTENT(IN), OPTIONAL :: DELTA
2971  REAL, INTENT(OUT) :: VAL1,VAL2
2972 
2973  REAL :: DENOM, X1, X2, Y1, Y2,C1,C2,C3,C4
2974  INTEGER :: INX,JNX
2975 
2976  inx = xycurvisearch(x_len, xarray, x, delta)
2977  jnx = xycurvisearch(y_len, yarray, y, delta)
2978  !
2979  IF (inx .GE. x_len) THEN
2980  inx=inx-1
2981  END IF
2982  IF (jnx .GE. y_len) THEN
2983  jnx=jnx-1
2984  END IF
2985  !
2986  x1 = xarray(inx)
2987  x2 = xarray(inx+1)
2988  y1 = yarray(jnx)
2989  y2 = yarray(jnx+1)
2990  !
2991  denom = (x2 - x1)*(y2 - y1)
2992  c1=(x2-x)*(y2-y)
2993  c2=(x-x1)*(y2-y)
2994  c3=(x2-x)*(y-y1)
2995  c4=(x-x1)*(y-y1)
2996  val1 = (func1(inx,jnx) *c1 + func1(inx+1,jnx) *c2 + &
2997  func1(inx,jnx+1)*c3 + func1(inx+1,jnx+1)*c4)/denom
2998 
2999  val2 = (func2(inx,jnx) *c1 + func2(inx+1,jnx) *c2 + &
3000  func2(inx,jnx+1)*c3 + func2(inx+1,jnx+1)*c4)/denom
3001  !

References xycurvisearch().

Referenced by wmupdv().

◆ wmupd1()

subroutine wmupdtmd::wmupd1 ( integer, intent(in)  IMOD,
character(len=3), intent(in)  IDSTR,
integer, intent(in)  J,
integer, intent(out)  IERR 
)

Update selected input using native input files.

Reading from native grid files.

Parameters
[in]IMODModel number.
[in]IDSTRID string corresponding to J.
[in]JInput type.
[out]IERRError indicator.
Author
H. L. Tolman
Date
22-Mar-2021

Definition at line 516 of file wmupdtmd.F90.

516  !/
517  !/ +-----------------------------------+
518  !/ | WAVEWATCH III NOAA/NCEP |
519  !/ | H. L. Tolman |
520  !/ | FORTRAN 90 |
521  !/ | Last update : 22-Mar-2021 |
522  !/ +-----------------------------------+
523  !/
524  !/ 07-Oct-2006 : Origination. ( version 3.10 )
525  !/ 22-Mar-2021 : Add momentum and air density input ( version 7.13 )
526  !/
527  ! 1. Purpose :
528  !
529  ! Update selected input using native input files.
530  !
531  ! 2. Method :
532  !
533  ! Reading from native grid files.
534  !
535  ! 3. Parameters :
536  !
537  ! Parameter list
538  ! ----------------------------------------------------------------
539  ! IMOD Int. I Model number,
540  ! IDSTR C*3 I ID string corresponding to J.
541  ! J Int. I Input type.
542  ! IERR Int. O Error indicator.
543  ! ----------------------------------------------------------------
544  !
545  ! 4. Subroutines used :
546  !
547  ! Name Type Module Description
548  ! ----------------------------------------------------------------
549  ! WMDIMD Subr. WMMDATMD Set dimension of data grids.
550  ! W3FLDG Subr. W3FLDSMD Get input field.
551  ! W3FLDD Subr. Id. Get input data.
552  ! W3FLDM Subr. Id. Get grid speed data.
553  ! STRACE Subr. W3ERVMD Subroutine tracing.
554  ! ----------------------------------------------------------------
555  !
556  ! 5. Called by :
557  !
558  ! Name Type Module Description
559  ! ----------------------------------------------------------------
560  ! WMUPDT Subr. WMUPDTMD Master inpu update routine.
561  ! ----------------------------------------------------------------
562  !
563  ! 6. Error messages :
564  !
565  ! 7. Remarks :
566  !
567  ! - Pointers set in calling routine.
568  !
569  ! 8. Structure :
570  !
571  ! See source code.
572  !
573  ! 9. Switches :
574  !
575  ! !/S Enable subroutine tracing.
576  ! !/T Enable test output
577  !
578  ! 10. Source code :
579  !
580  !/ ------------------------------------------------------------------- /
581  !/
582  USE wmmdatmd, ONLY: wmdimd
583  USE w3fldsmd, ONLY: w3fldg, w3fldd, w3fldm
584 #ifdef W3_S
585  USE w3servmd, ONLY: strace
586 #endif
587  !/
588  USE w3gdatmd, ONLY: nx, ny
589 #ifdef W3_SMC
590  USE w3gdatmd, ONLY: fswnd, nsea
591 #endif
592  USE w3wdatmd, ONLY: time
593  USE w3idatmd, ONLY: tln, wlev, tc0, tcn, cx0, cxn, cy0, cyn, &
594  tw0, twn, tu0, tun, tr0, trn, wx0, wxn, &
595  wy0, wyn, dt0, dtn, tin, trn, icei, ux0, &
596  uxn, uy0, uyn, rh0, rhn, t0n, t1n, t2n, &
597  tdn, inflags1, tg0, tgn, ga0, gd0, gan, &
598  gdn, bergi, ttn, mudt, tvn, mudv, tzn, &
599  mudd, ti1, ti2, ti3, ti4, ti5, icep1, &
600  icep2, icep3, icep4, icep5
601  USE wmmdatmd, ONLY: improc, nmperr, mdst, mdse, mdsf, etime, &
602  fllstl, fllsti, fllstr, rcld, ndt, data0, &
603  data1, data2, nmv, nmvmax, tmv, amv, dmv
604  !/
605  IMPLICIT NONE
606  !/
607  !/ ------------------------------------------------------------------- /
608  !/ Parameter list
609  !/
610  INTEGER, INTENT(IN) :: IMOD, J
611  INTEGER, INTENT(OUT) :: IERR
612  CHARACTER(LEN=3), INTENT(IN) :: IDSTR
613  !/
614  !/ ------------------------------------------------------------------- /
615  !/ Local parameters
616  !/
617  INTEGER :: MDSEN, DTIME(2), NDTNEW
618  REAL :: XXX(NY,NX)
619 #ifdef W3_S
620  INTEGER, SAVE :: IENT = 0
621 #endif
622  !/
623  !/ ------------------------------------------------------------------- /
624  ! 0. Initialization
625  ! 0.a Subroutine tracing and echo of input
626  !
627 #ifdef W3_S
628  CALL strace (ient, 'WMUPD1')
629 #endif
630 #ifdef W3_T
631  WRITE (mdst,9000) imod, j
632 #endif
633  !
634  IF ( improc .EQ. nmperr ) THEN
635  mdsen = mdse
636  ELSE
637  mdsen = -1
638  END IF
639  !
640  ! 0.b Start case selection
641  !
642  SELECT CASE (j)
643  !
644  ! -7. Ice parameter 1 ---------------------------------------------- /
645  !
646  CASE (-7)
647  CALL w3fldg ('READ', idstr, mdsf(imod,j), mdst, mdsen, &
648  nx, ny, nx, ny, time, etime, dtime, &
649  xxx, xxx, xxx, ti1, xxx, xxx, icep1, ierr)
650  !
651  ! -6. Ice parameter 2 ---------------------------------------------- /
652  !
653  CASE (-6)
654  CALL w3fldg ('READ', idstr, mdsf(imod,j), mdst, mdsen, &
655  nx, ny, nx, ny, time, etime, dtime, &
656  xxx, xxx, xxx, ti2, xxx, xxx, icep2, ierr)
657  !
658  ! -5. Ice parameter 3 ---------------------------------------------- /
659  !
660  CASE (-5)
661  CALL w3fldg ('READ', idstr, mdsf(imod,j), mdst, mdsen, &
662  nx, ny, nx, ny, time, etime, dtime, &
663  xxx, xxx, xxx, ti3, xxx, xxx, icep3, ierr)
664  !
665  ! -4. Ice parameter 4 ---------------------------------------------- /
666  !
667  CASE (-4)
668  CALL w3fldg ('READ', idstr, mdsf(imod,j), mdst, mdsen, &
669  nx, ny, nx, ny, time, etime, dtime, &
670  xxx, xxx, xxx, ti4, xxx, xxx, icep4, ierr)
671  !
672  ! -3. Ice parameter 5 ---------------------------------------------- /
673  !
674  CASE (-3)
675  CALL w3fldg ('READ', idstr, mdsf(imod,j), mdst, mdsen, &
676  nx, ny, nx, ny, time, etime, dtime, &
677  xxx, xxx, xxx, ti5, xxx, xxx, icep5, ierr)
678  !
679  ! -2. Mud Density -------------------------------------------------- /
680  !
681  CASE (-2)
682  CALL w3fldg ('READ', idstr, mdsf(imod,j), mdst, mdsen, &
683  nx, ny, nx, ny, time, etime, dtime, &
684  xxx, xxx, xxx, tzn, xxx, xxx, mudd, ierr)
685  !
686  ! -1. Mud Thickness -------------------------------------------------- /
687  !
688  CASE (-1)
689  CALL w3fldg ('READ', idstr, mdsf(imod,j), mdst, mdsen, &
690  nx, ny, nx, ny, time, etime, dtime, &
691  xxx, xxx, xxx, ttn, xxx, xxx, mudt, ierr)
692  !
693  ! 0. Mud Viscosity -------------------------------------------------- /
694  !
695  CASE (0)
696  CALL w3fldg ('READ', idstr, mdsf(imod,j), mdst, mdsen, &
697  nx, ny, nx, ny, time, etime, dtime, &
698  xxx, xxx, xxx, tvn, xxx, xxx, mudv, ierr)
699  !
700  ! 1. Water levels --------------------------------------------------- /
701  !
702  CASE (1)
703  CALL w3fldg ('READ', idstr, mdsf(imod,j), mdst, mdsen, &
704  nx, ny, nx, ny, time, etime, dtime, &
705  xxx, xxx, xxx, tln, xxx, xxx, wlev, ierr)
706  IF ( ierr .LT. 0 ) fllstl = .true.
707  !
708  ! 2. Currents ------------------------------------------------------- /
709  !
710  CASE (2)
711 #ifdef W3_SMC
712  !!Li For sea point current option FSWND. JGLi08Feb2021
713  IF( fswnd ) THEN
714  CALL w3fldg ('READ', idstr, mdsf(imod,j), mdst, mdsen, &
715  nsea, 1, nsea, 1, time, etime, tc0, &
716  cx0, cy0, xxx, tcn, cxn, cyn, xxx, ierr)
717  ELSE
718 #endif
719  CALL w3fldg ('READ', idstr, mdsf(imod,j), mdst, mdsen, &
720  nx, ny, nx, ny, time, etime, tc0, &
721  cx0, cy0, xxx, tcn, cxn, cyn, xxx, ierr)
722 #ifdef W3_SMC
723  END IF
724 #endif
725  !
726  ! 3. Winds ---------------------------------------------------------- /
727  !
728  CASE (3)
729 #ifdef W3_SMC
730  !!Li For sea point wind option FSWND. JGLi08Feb2021
731  IF( fswnd ) THEN
732  CALL w3fldg ('READ', idstr, mdsf(imod,j), mdst, mdsen, &
733  nsea, 1, nsea, 1, time, etime, tw0, &
734  wx0, wy0, dt0, twn, wxn, wyn, dtn, ierr)
735  ELSE
736 #endif
737  CALL w3fldg ('READ', idstr, mdsf(imod,j), mdst, mdsen, &
738  nx, ny, nx, ny, time, etime, tw0, &
739  wx0, wy0, dt0, twn, wxn, wyn, dtn, ierr)
740 #ifdef W3_SMC
741  END IF
742 #endif
743  !
744  ! 4. Ice ------------------------------------------------------------ /
745  !
746  CASE (4)
747  CALL w3fldg ('READ', idstr, mdsf(imod,j), mdst, mdsen, &
748  nx, ny, nx, ny, time, etime, dtime, &
749  xxx, xxx, xxx, tin, xxx , bergi, icei, ierr)
750  IF ( ierr .LT. 0 ) fllsti = .true.
751  !
752  ! 5. Momentum ------------------------------------------------------- /
753  !
754  CASE (5)
755  CALL w3fldg ('READ', idstr, mdsf(imod,j), mdst, mdsen, &
756  nx, ny, nx, ny, time, etime, tu0, &
757  ux0, uy0, xxx, tun, uxn, uyn, xxx, ierr)
758  !
759  ! 6. Air density ---------------------------------------------------- /
760  !
761  CASE (6)
762  CALL w3fldg ('READ', idstr, mdsf(imod,j), mdst, mdsen, &
763  nx, ny, nx, ny, time, etime, tr0, &
764  xxx, xxx, rh0, trn, xxx, xxx, rhn, ierr)
765  IF ( ierr .LT. 0 ) fllstr = .true.
766  !
767  ! 7. Data type 0 ---------------------------------------------------- /
768  !
769  CASE (7)
770  CALL w3fldd ('SIZE', idstr, mdsf(imod,j), mdst, mdsen, &
771  time, t0n, rcld(1), ndt(1), ndtnew, &
772  data0, ierr )
773  IF ( ierr .LT. 0 ) THEN
774  inflags1(j) = .false.
775  rcld(1) = 1
776  ndt(1) = 1
777  CALL wmdimd ( imod, mdse, mdst, 1 )
778  ELSE
779  ndt(j) = ndtnew
780  CALL wmdimd ( imod, mdse, mdst, 1 )
781  CALL w3fldd ('SIZE', idstr, mdsf(imod,j), mdst, &
782  mdsen, time, t0n, rcld(1), ndt(1), &
783  ndtnew, data0, ierr )
784  END IF
785  !
786  ! 8. Data type 1 ---------------------------------------------------- /
787  !
788  CASE ( 8 )
789  CALL w3fldd ('SIZE', idstr, mdsf(imod,j), mdst, mdsen, &
790  time, t1n, rcld(2), ndt(2), ndtnew, &
791  data1, ierr )
792  IF ( ierr .LT. 0 ) THEN
793  inflags1(j) = .false.
794  rcld(2) = 1
795  ndt(2) = 1
796  CALL wmdimd ( imod, mdse, mdst, 2 )
797  ELSE
798  ndt(j) = ndtnew
799  CALL wmdimd ( imod, mdse, mdst, 2 )
800  CALL w3fldd ('SIZE', idstr, mdsf(imod,j), mdst, &
801  mdsen, time, t1n, rcld(2), ndt(2), &
802  ndtnew, data1, ierr )
803  END IF
804  !
805  ! 9. Data type 2 ---------------------------------------------------- /
806  !
807  CASE ( 9 )
808  CALL w3fldd ('SIZE', idstr, mdsf(imod,j), mdst, mdsen, &
809  time, t2n, rcld(3), ndt(3), ndtnew, &
810  data2, ierr )
811  IF ( ierr .LT. 0 ) THEN
812  inflags1(j) = .false.
813  rcld(3) = 1
814  ndt(3) = 1
815  CALL wmdimd ( imod, mdse, mdst, 3 )
816  ELSE
817  ndt(j) = ndtnew
818  CALL wmdimd ( imod, mdse, mdst, 3 )
819  CALL w3fldd ('SIZE', idstr, mdsf(imod,j), mdst, &
820  mdsen, time, t2n, rcld(3), ndt(3), &
821  ndtnew, data2, ierr )
822  END IF
823  !
824  ! 10. Moving grid data ---------------------------------------------- /
825  !
826  CASE ( 10 )
827  ! notes:
828  ! SUBROUTINE W3FLDM in w3fldsmd.ftn :
829  ! INTEGER, INTENT(INOUT) :: NH, THO(2,6,NHM), TF0(2), TFN(2)
830  ! INTEGER, INTENT(INOUT) :: NH, THO(2,-7:6,NHM), TF0(2), TFN(2)
831  ! REAL, INTENT(INOUT) :: HA(NHM,6), HD(NHM,6), A0, AN, D0, DN
832  ! REAL, INTENT(INOUT) :: HA(NHM,-7:6), HD(NHM,-7:6), A0, AN, D0, DN
833  ! Arguments #
834  ! THO 8
835  ! HA 9
836  ! HD 10
837  ! Here, that is TMV AMV DMV
838  CALL w3fldm ( 4, mdst, mdsen, time, etime, nmv, nmvmax, tmv,&
839  amv, dmv, tg0, ga0, gd0, tgn, gan, gdn, ierr )
840  !
841  END SELECT
842  !
843  ! 9. End of routine -------------------------------------------------- /
844  !
845  RETURN
846  !
847  ! Formats
848  !
849 #ifdef W3_T
850 9000 FORMAT ( ' TEST WMUPD1 : INPUT : ',2i4)
851 #endif
852  !/
853  !/ End of WMUPD1 ----------------------------------------------------- /
854  !/

References wmmdatmd::amv, wmmdatmd::data0, wmmdatmd::data1, wmmdatmd::data2, wmmdatmd::dmv, w3idatmd::dt0, w3idatmd::dtn, wmmdatmd::etime, wmmdatmd::fllsti, wmmdatmd::fllstl, wmmdatmd::fllstr, w3gdatmd::fswnd, w3idatmd::ga0, w3idatmd::gan, w3idatmd::gd0, w3idatmd::gdn, wmmdatmd::improc, w3idatmd::inflags1, wmmdatmd::mdse, wmmdatmd::mdsf, wmmdatmd::mdst, wmmdatmd::ndt, wmmdatmd::nmperr, wmmdatmd::nmv, wmmdatmd::nmvmax, w3gdatmd::nsea, w3gdatmd::nx, w3gdatmd::ny, wmmdatmd::rcld, w3servmd::strace(), w3idatmd::t0n, w3idatmd::t1n, w3idatmd::t2n, w3idatmd::tc0, w3idatmd::tcn, w3idatmd::tdn, w3idatmd::tg0, w3idatmd::tgn, w3idatmd::ti1, w3idatmd::ti2, w3idatmd::ti3, w3idatmd::ti4, w3idatmd::ti5, w3wdatmd::time, w3idatmd::tin, w3idatmd::tln, wmmdatmd::tmv, w3idatmd::tr0, w3idatmd::trn, w3idatmd::ttn, w3idatmd::tu0, w3idatmd::tun, w3idatmd::tvn, w3idatmd::tw0, w3idatmd::twn, w3idatmd::tzn, w3fldsmd::w3fldd(), w3fldsmd::w3fldg(), w3fldsmd::w3fldm(), wmmdatmd::wmdimd(), w3idatmd::wx0, w3idatmd::wxn, w3idatmd::wy0, and w3idatmd::wyn.

Referenced by wmupdt().

◆ wmupd2()

subroutine wmupdtmd::wmupd2 ( integer, intent(in)  IMOD,
integer, intent(in)  J,
integer, intent(in)  JMOD,
integer, intent(out)  IERR 
)

Update selected input using input grids.

Managing data, interpolation done in other routines.

Parameters
[in]IMODModel number.
[in]JInput type.
[in]JMODModel number source grid.
[out]IERRError indicator.
Author
H. L. Tolman
Date
22-Mar-2021

Definition at line 870 of file wmupdtmd.F90.

870  !/
871  !/ +-----------------------------------+
872  !/ | WAVEWATCH III NOAA/NCEP |
873  !/ | H. L. Tolman |
874  !/ | FORTRAN 90 |
875  !/ | Last update : 22-Mar-2021 |
876  !/ +-----------------------------------+
877  !/
878  !/ 14-Oct-2006 : Origination. ( version 3.10 )
879  !/ 10-Dec-2006 : Bug fix WMUPD2 initial fields. ( version 3.10 )
880  !/ 22-Mar-2021 : Add momentum and air density input ( version 7.13 )
881  !/
882  ! 1. Purpose :
883  !
884  ! Update selected input using input grids.
885  !
886  ! 2. Method :
887  !
888  ! Managing data, interpolation done inother routines.
889  !
890  ! 3. Parameters :
891  !
892  ! Parameter list
893  ! ----------------------------------------------------------------
894  ! IMOD Int. I Model number,
895  ! J Int. I Input type.
896  ! JMOD Int. I Model number source grid.
897  ! IERR Int. O Error indicator.
898  ! ----------------------------------------------------------------
899  !
900  ! 4. Subroutines used :
901  !
902  ! Name Type Module Description
903  ! ----------------------------------------------------------------
904  ! STRACE Subr. W3ERVMD Subroutine tracing.
905  ! EXTCDE Subr. Id. Program abort.
906  ! WMUPDV Subr. local Interpolation of vector fields.
907  ! WMUPDS Subr. local Interpolation of scalar fields.
908  ! ----------------------------------------------------------------
909  !
910  ! 5. Called by :
911  !
912  ! Name Type Module Description
913  ! ----------------------------------------------------------------
914  ! WMUPDT Subr. WMUPDTMD Master input update routine.
915  ! ----------------------------------------------------------------
916  !
917  ! 6. Error messages :
918  !
919  ! 7. Remarks :
920  !
921  ! 8. Structure :
922  !
923  ! See source code.
924  !
925  ! 9. Switches :
926  !
927  ! !/CRX0 Current vector component conservation.
928  ! !/CRX1 Current speed conservation.
929  ! !/CRX2 Current exenrgy conservation.
930  !
931  ! !/WNX0 Wind vector component conservation.
932  ! !/WNX1 Wind speed conservation.
933  ! !/WNX2 Wind exenrgy conservation.
934  !
935  ! !/S Enable subroutine tracing.
936  ! !/T Enable test output
937  !
938  ! 10. Source code :
939  !
940  !/ ------------------------------------------------------------------- /
941  !/
942  USE w3servmd, ONLY: extcde
943 #ifdef W3_S
944  USE w3servmd, ONLY: strace
945 #endif
946  !/
947  USE w3wdatmd, ONLY: time
948  USE w3idatmd, ONLY: inputs
949  USE wmmdatmd, ONLY: improc, nmperr, nmpscr, mdst, mdse, mdss, &
950  mdso, etime, idinp
951  USE constants, ONLY: dair
952  !/
953  IMPLICIT NONE
954  !/
955  !/ ------------------------------------------------------------------- /
956  !/ Parameter list
957  !/
958  INTEGER, INTENT(IN) :: IMOD, J, JMOD
959  INTEGER, INTENT(OUT) :: IERR
960  !/
961  !/ ------------------------------------------------------------------- /
962  !/ Local parameters
963  !/
964  INTEGER :: ICONSC, ICONSW, ICONSU
965 #ifdef W3_S
966  INTEGER, SAVE :: IENT = 0
967 #endif
968  !/
969  !/ ------------------------------------------------------------------- /
970  ! 0. Initialization
971  ! 0.a Subroutine tracing and echo of input
972  !
973 #ifdef W3_S
974  CALL strace (ient, 'WMUPD2')
975 #endif
976  !
977 #ifdef W3_T
978  WRITE (mdst,9000) imod, j, jmod
979  WRITE (mdst,9001) inputs(imod)%TFN(:,j), &
980  inputs(jmod)%TFN(:,j), etime
981 #endif
982  !
983  ierr = 0
984 #ifdef W3_CRX0
985  iconsc = 0
986 #endif
987 #ifdef W3_CRX1
988  iconsc = 1
989 #endif
990 #ifdef W3_CRX2
991  iconsc = 2
992 #endif
993 #ifdef W3_WNX0
994  iconsw = 0
995 #endif
996 #ifdef W3_WNX1
997  iconsw = 1
998 #endif
999 #ifdef W3_WNX2
1000  iconsw = 2
1001 #endif
1002 #ifdef W3_WNX0
1003  iconsu = 0
1004 #endif
1005 #ifdef W3_WNX1
1006  iconsu = 1
1007 #endif
1008 #ifdef W3_WNX2
1009  iconsu = 2
1010 #endif
1011  !
1012  ! 1. Shift fields ( currents and winds only ) ------------------------ /
1013  !
1014  SELECT CASE (j)
1015  !
1016  ! 1.a Currents
1017  !
1018  CASE (2)
1019  IF ( inputs(imod)%TFN(1,j) .GT. 0 ) THEN
1020  inputs(imod)%TC0(:) = inputs(imod)%TFN(:,j)
1021  inputs(imod)%CX0 = inputs(imod)%CXN
1022  inputs(imod)%CY0 = inputs(imod)%CYN
1023 #ifdef W3_T
1024  WRITE (mdst,9010) j, inputs(imod)%TFN(:,j)
1025  ELSE
1026  WRITE (mdst,9011) j
1027 #endif
1028  END IF
1029  !
1030  ! 1.b Winds
1031  !
1032  CASE (3)
1033  IF ( inputs(imod)%TFN(1,j) .GT. 0 ) THEN
1034  inputs(imod)%TW0(:) = inputs(imod)%TFN(:,j)
1035  inputs(imod)%WX0 = inputs(imod)%WXN
1036  inputs(imod)%WY0 = inputs(imod)%WYN
1037  inputs(imod)%DT0 = inputs(imod)%DTN
1038 #ifdef W3_T
1039  WRITE (mdst,9010) j, inputs(imod)%TFN(:,j)
1040  ELSE
1041  WRITE (mdst,9011) j
1042 #endif
1043  END IF
1044  !
1045  ! 1.c Momentum
1046  !
1047  CASE (5)
1048  IF ( inputs(imod)%TFN(1,j) .GT. 0 ) THEN
1049  inputs(imod)%TU0(:) = inputs(imod)%TFN(:,j)
1050  inputs(imod)%UX0 = inputs(imod)%UXN
1051  inputs(imod)%UY0 = inputs(imod)%UYN
1052 #ifdef W3_T
1053  WRITE (mdst,9010) j, inputs(imod)%TFN(:,j)
1054  ELSE
1055  WRITE (mdst,9011) j
1056 #endif
1057  END IF
1058  !
1059  END SELECT
1060  !
1061  ! 2. Process fields at ending time ----------------------------------- /
1062  !
1063 #ifdef W3_T
1064  WRITE (mdst,9020) j, inputs(jmod)%TFN(:,j)
1065 #endif
1066  inputs(imod)%TFN(:,j) = inputs(jmod)%TFN(:,j)
1067  !
1068  SELECT CASE (j)
1069  !
1070  ! 2.a-3 Ice parameter 1
1071  !
1072  CASE (-7)
1073  CALL wmupds ( imod, inputs(imod)%ICEP1, &
1074  jmod, inputs(jmod)%ICEP1, 0. )
1075  !
1076  ! 2.a-3 Ice parameter 2
1077  !
1078  CASE (-6)
1079  CALL wmupds ( imod, inputs(imod)%ICEP2, &
1080  jmod, inputs(jmod)%ICEP2, 0. )
1081  !
1082  ! 2.a-3 Ice parameter 3
1083  !
1084  CASE (-5)
1085  CALL wmupds ( imod, inputs(imod)%ICEP3, &
1086  jmod, inputs(jmod)%ICEP3, 0. )
1087  !
1088  ! 2.a-3 Ice parameter 4
1089  !
1090  CASE (-4)
1091  CALL wmupds ( imod, inputs(imod)%ICEP4, &
1092  jmod, inputs(jmod)%ICEP4, 0. )
1093  !
1094  ! 2.a-3 Ice parameter 5
1095  !
1096  CASE (-3)
1097 
1098  CALL wmupds ( imod, inputs(imod)%ICEP5, &
1099  jmod, inputs(jmod)%ICEP5, 0. )
1100  !
1101  ! 2.a-2 Mud densities
1102  !
1103  CASE (-2)
1104  CALL wmupds ( imod, inputs(imod)%MUDD, &
1105  jmod, inputs(jmod)%MUDD, 0. )
1106  !
1107  ! 2.a-1 Mud viscosities
1108  !
1109  CASE (-1)
1110  CALL wmupds ( imod, inputs(imod)%MUDT, &
1111  jmod, inputs(jmod)%MUDT, 0. )
1112  !
1113  ! 2.a-0 Mud thicknesses
1114  !
1115  CASE (0)
1116  CALL wmupds ( imod, inputs(imod)%MUDV, &
1117  jmod, inputs(jmod)%MUDV, 0. )
1118  !
1119  ! 2.a Water levels
1120  !
1121  CASE (1)
1122  CALL wmupds ( imod, inputs(imod)%WLEV, &
1123  jmod, inputs(jmod)%WLEV, 0. )
1124  !
1125  ! 2.b Curents
1126  !
1127  CASE (2)
1128  CALL wmupdv ( imod, inputs(imod)%CXN, inputs(imod)%CYN, &
1129  jmod, inputs(jmod)%CXN, inputs(jmod)%CYN, &
1130  0., iconsc )
1131  !
1132  ! 2.c Wind speeds
1133  !
1134  CASE (3)
1135  CALL wmupdv ( imod, inputs(imod)%WXN, inputs(imod)%WYN, &
1136  jmod, inputs(jmod)%WXN, inputs(jmod)%WYN, &
1137  0., iconsw )
1138  !
1139  IF ( idinp(imod,j) .EQ. 'WNS' ) CALL wmupds &
1140  ( imod, inputs(imod)%DTN, &
1141  jmod, inputs(jmod)%DTN, 0. )
1142  !
1143  ! 2.d Ice concentrations
1144  !
1145  CASE (4)
1146  CALL wmupds ( imod, inputs(imod)%ICEI, &
1147  jmod, inputs(jmod)%ICEI, 0. )
1148  IF ( idinp(imod,j) .EQ. 'ISI' ) CALL wmupds &
1149  ( imod, inputs(imod)%BERGI, &
1150  jmod, inputs(jmod)%BERGI, 0. )
1151  !
1152  ! 2.e Momentum
1153  !
1154  CASE (5)
1155  CALL wmupdv ( imod, inputs(imod)%UXN, inputs(imod)%UYN, &
1156  jmod, inputs(jmod)%UXN, inputs(jmod)%UYN, &
1157  0., iconsu )
1158  !
1159  ! 2.f Air density
1160  !
1161  CASE (6)
1162  CALL wmupds ( imod, inputs(imod)%RHN, &
1163  jmod, inputs(jmod)%RHN, dair )
1164  !
1165  ! 2.g Assimilation data 0
1166  !
1167  CASE (7)
1168  GOTO 2999
1169  !
1170  ! 2.h Assimilation data 1
1171  !
1172  CASE (8)
1173  GOTO 2999
1174  !
1175  ! 2.i Assimilation data 2
1176  !
1177  CASE (9)
1178  GOTO 2999
1179  !
1180  END SELECT
1181  !
1182  ! 3. Check and update first fields ( currents and winds only ) ------- /
1183  !
1184  SELECT CASE (j)
1185  !
1186  ! 3.a Currents
1187  !
1188  CASE (2)
1189  IF ( inputs(imod)%TC0(1) .LT. 0 ) THEN
1190  inputs(imod)%TC0(:) = inputs(jmod)%TC0(:)
1191 #ifdef W3_T
1192  WRITE (mdst,9030) j, inputs(imod)%TC0(:)
1193 #endif
1194 #ifdef W3_CRX0
1195  iconsc = 0
1196 #endif
1197 #ifdef W3_CRX1
1198  iconsc = 1
1199 #endif
1200 #ifdef W3_CRX2
1201  iconsc = 2
1202 #endif
1203  CALL wmupdv ( imod, inputs(imod)%CX0, inputs(imod)%CY0, &
1204  jmod, inputs(jmod)%CX0, inputs(jmod)%CY0, &
1205  0., iconsc )
1206  END IF
1207  !
1208  ! 3.b Winds
1209  !
1210  CASE (3)
1211  IF ( inputs(imod)%TW0(1) .LT. 0 ) THEN
1212  inputs(imod)%TW0(:) = inputs(jmod)%TW0(:)
1213 #ifdef W3_T
1214  WRITE (mdst,9030) j, inputs(imod)%TW0(:)
1215 #endif
1216 #ifdef W3_WNX0
1217  iconsw = 0
1218 #endif
1219 #ifdef W3_WNX1
1220  iconsw = 1
1221 #endif
1222 #ifdef W3_WNX2
1223  iconsw = 2
1224 #endif
1225  CALL wmupdv ( imod, inputs(imod)%WX0, inputs(imod)%WY0, &
1226  jmod, inputs(jmod)%WX0, inputs(jmod)%WY0, &
1227  0., iconsw )
1228  IF ( idinp(imod,j) .EQ. 'WNS' ) CALL wmupds &
1229  ( imod, inputs(imod)%DT0, &
1230  jmod, inputs(jmod)%DT0, 0. )
1231  END IF
1232  !
1233  ! 3.c Momentum
1234  !
1235  CASE (5)
1236  IF ( inputs(imod)%TU0(1) .LT. 0 ) THEN
1237  inputs(imod)%TU0(:) = inputs(jmod)%TU0(:)
1238 #ifdef W3_T
1239  WRITE (mdst,9030) j, inputs(imod)%TU0(:)
1240 #endif
1241 #ifdef W3_WNX0
1242  iconsu = 0
1243 #endif
1244 #ifdef W3_WNX1
1245  iconsu = 1
1246 #endif
1247 #ifdef W3_WNX2
1248  iconsu = 2
1249 #endif
1250  CALL wmupdv ( imod, inputs(imod)%UX0, inputs(imod)%UY0, &
1251  jmod, inputs(jmod)%UX0, inputs(jmod)%UY0, &
1252  0., iconsu )
1253  END IF
1254  !
1255  END SELECT
1256  !
1257  ! 4. End of routine -------------------------------------------------- /
1258  !
1259  RETURN
1260  !
1261  ! Error escape locations
1262  !
1263 2999 CONTINUE
1264  IF ( mdss.NE.mdso .AND. nmpscr.EQ.improc ) WRITE (mdse,1999)
1265  CALL extcde ( 2999 )
1266  RETURN
1267  !
1268  ! Formats
1269  !
1270 1999 FORMAT (/' *** ERROR WMUPD2: OPTION NOT YET IMPLEMENTED ***'/)
1271  !
1272 #ifdef W3_T
1273 9000 FORMAT ( ' TEST WMUPD2 : INPUT : ',3i4)
1274 9001 FORMAT ( ' TEST WMUPD2 : TIME OF IMOD : ',i9.8,1x,i6.6/ &
1275  ' TIME OF JMOD : ',i9.8,1x,i6.6/ &
1276  ' ENDING TIME : ',i9.8,1x,i6.6)
1277 9010 FORMAT ( ' TEST WMUPD2 : SHIFTING ',i1,' TIME = ',i8.8,i7.6)
1278 9011 FORMAT ( ' TEST WMUPD2 : NO DATA FOR ',i1,' TO SHIFT')
1279 9020 FORMAT ( ' TEST WMUPD2 : PROCESSING ',i1,' TIME = ',i8.8,i7.6)
1280 9030 FORMAT ( ' TEST WMUPD2 : INITIAL FIELD FOR ',i1, &
1281  ' TIME = ',i8.8,i7.6)
1282 #endif
1283  !/
1284  !/ End of WMUPD2 ----------------------------------------------------- /
1285  !/

References constants::dair, wmmdatmd::etime, w3servmd::extcde(), wmmdatmd::idinp, wmmdatmd::improc, w3idatmd::inputs, wmmdatmd::mdse, wmmdatmd::mdso, wmmdatmd::mdss, wmmdatmd::mdst, wmmdatmd::nmperr, wmmdatmd::nmpscr, w3servmd::strace(), w3wdatmd::time, wmupds(), and wmupdv().

Referenced by wmupdt().

◆ wmupds()

subroutine wmupdtmd::wmupds ( integer, intent(in)  IMOD,
real, dimension(nx,ny), intent(out)  FD,
integer, intent(in)  JMOD,
real, dimension(grids(jmod)%nx,grids(jmod)%ny), intent(in)  FDI,
real, intent(in)  UNDEF 
)

Interpolate scalar field from input grid to model grid.

Interpolating or averaging from input grid.

Parameters
[in]IMODOutput model number
[out]FDOutput scalar field
[in]JMODInput model number
[in]FDIInput scalar field
[in]UNDEFValue for mapped out point and points not covered.
Author
H. L. Tolman
Date
06-Dec-2010

Definition at line 1986 of file wmupdtmd.F90.

1986  !/
1987  !/ +-----------------------------------+
1988  !/ | WAVEWATCH III NOAA/NCEP |
1989  !/ | H. L. Tolman |
1990  !/ | FORTRAN 90 |
1991  !/ | Last update : 06-Dec-2010 |
1992  !/ +-----------------------------------+
1993  !/
1994  !/ 14-Oct-2006 : Origination. ( version 3.10 )
1995  !/ 12-Jan-2007 : General clean-up and bug fixes. ( version 3.10 )
1996  !/ 30-Oct-2009 : Implement run-time grid selection. ( version 3.14 )
1997  !/ (W. E. Rogers & T. J. Campbell, NRL)
1998  !/ 06-Dec-2010 : Change from GLOBAL (logical) to ICLOSE (integer) to
1999  !/ specify index closure for a grid. ( version 3.14 )
2000  !/ (T. J. Campbell, NRL)
2001  !/ 11-May-2015 : Updates to 2-ways nestings for UG ( version 5.08 )
2002  !/ 01-Jul-2019 : Generalize output to curv grids ( version 7.13 )
2003  !/ (R. Padilla-Hernandez, J.H. Alves, EMC/NOAA)
2004  !/
2005  ! 1. Purpose :
2006  !
2007  ! Interpolate scalar field from input grid to model grid.
2008  !
2009  ! 2. Method :
2010  !
2011  ! Interpolating or averaging from input grid.
2012  !
2013  ! 3. Parameters :
2014  !
2015  ! Parameter list
2016  ! ----------------------------------------------------------------
2017  ! IMOD Int. I Output model number,
2018  ! FD Int. O Output scaler field.
2019  ! JMOD Int. I Input model number,
2020  ! FDI Int. I Input scaler field.
2021  ! UNDEF Int. I Value for mapped out point and points not
2022  ! covered.
2023  ! ----------------------------------------------------------------
2024  !
2025  ! 4. Subroutines used :
2026  !
2027  ! Name Type Module Description
2028  ! ----------------------------------------------------------------
2029  ! STRACE Subr. W3ERVMD Subroutine tracing.
2030  ! EXTCDE Subr. Id. Program abort.
2031  ! ----------------------------------------------------------------
2032  !
2033  ! 5. Called by :
2034  !
2035  ! Name Type Module Description
2036  ! ----------------------------------------------------------------
2037  ! WMUPD2 Subr. WMUPDTMD Input update routine.
2038  ! ----------------------------------------------------------------
2039  !
2040  ! 6. Error messages :
2041  !
2042  ! 7. Remarks :
2043  !
2044  ! - Grid pointers for output grid need to be set externally.
2045  ! - If input grid does not cover point of target grid, target grid
2046  ! values are set to UNDEF.
2047  !
2048  ! 8. Structure :
2049  !
2050  ! See source code.
2051  !
2052  ! 9. Switches :
2053  !
2054  ! !/S Enable subroutine tracing.
2055  ! !/T Enable test output
2056  ! !/T1 Test output interpolation data.
2057  !
2058  ! 10. Source code :
2059  !
2060  !/ ------------------------------------------------------------------- /
2061  !/
2062  USE w3servmd, ONLY: extcde
2063 #ifdef W3_S
2064  USE w3servmd, ONLY: strace
2065 #endif
2066  !/
2067  USE w3gdatmd, ONLY: nx, ny, x0, y0, sx, sy, grids, flagll, &
2070  hpfac, hqfac, xgrd, ygrd
2071  USE wmmdatmd, ONLY: improc, nmperr, nmpscr, mdst, mdse, mdso, &
2072  mdss
2073  !/
2074  IMPLICIT NONE
2075  !/
2076  !/ ------------------------------------------------------------------- /
2077  !/ Parameter list
2078  !/
2079  INTEGER, INTENT(IN) :: IMOD, JMOD
2080  REAL, INTENT(OUT) :: FD(NX,NY)
2081  REAL, INTENT(IN) :: FDI(GRIDS(JMOD)%NX,GRIDS(JMOD)%NY), &
2082  UNDEF
2083  !/
2084  !/ ------------------------------------------------------------------- /
2085  !/ Local parameters
2086  !/
2087  INTEGER :: IXO, IYO, IX, IY, IXF0, IXFN, IYF0, &
2088  IYFN, IXS0, IXSN, IYS0, IYSN, IXS, &
2089  MXA, MYA, J, J1, J2, IXC, IYC, JJ, &
2090  JX, JY
2091 
2092  INTEGER :: NPOIX, NPOIY, I, CURVI !RP
2093 
2094 #ifdef W3_S
2095  INTEGER, SAVE :: IENT = 0
2096 #endif
2097  INTEGER, ALLOCATABLE :: NXA(:,:), NYA(:,:)
2098  REAL :: XR, YR, R1, R2, RT, XFL, XFR, XSL, &
2099  XSR, YFL, YFR, YSL, YSR
2100  REAL :: FDL, WTOT, WL
2101 
2102  REAL :: LONC, LATC, SXYC, &
2103  XDI, DTOLER, VALUEINTER
2104 
2105  REAL, ALLOCATABLE :: RXA(:,:), RYA(:,:)
2106 
2107 
2108  LOGICAL :: MAP1(NX,NY), MAP2(NX,NY), &
2109  MAP3(NX,NY), FLAGUP
2110  !
2111  INTEGER, POINTER :: NXI, NYI, MAP(:,:), MAPI(:,:)
2112 
2113  DOUBLE PRECISION, POINTER :: XGRDI(:,:), YGRDI(:,:), XGRDC(:,:), YGRDC(:,:)
2114  REAL, POINTER :: HPFACI(:,:), HQFACI(:,:) !RP
2115 
2116  REAL, POINTER :: X0I, Y0I, SXI, SYI !RPXXX , HPFACI, HQFACI
2117  INTEGER, POINTER :: ICLOSE
2118 #ifdef W3_T1
2119  CHARACTER(LEN=17) :: FORMAT1
2120 #endif
2121  !/
2122  !/ ------------------------------------------------------------------- /
2123  ! 0. Initialization
2124  ! 0.a Subroutine tracing and test output
2125  !
2126 #ifdef W3_S
2127  CALL strace (ient, 'WMUPDS')
2128 #endif
2129  !
2130  nxi => grids(jmod)%NX
2131  nyi => grids(jmod)%NY
2132  x0i => grids(jmod)%X0
2133  y0i => grids(jmod)%Y0
2134  sxi => grids(jmod)%SX
2135  syi => grids(jmod)%SY
2136  hpfaci => grids(jmod)%HPFAC
2137  hqfaci => grids(jmod)%HQFAC
2138  map => grids(imod)%MAPSTA
2139  mapi => grids(jmod)%MAPSTA
2140  iclose => grids(jmod)%ICLOSE
2141  !
2142  IF ( iclose .EQ. iclose_trpl ) THEN
2143  IF ( improc.EQ.nmperr ) WRITE(mdse,*)'SUBROUTINE WMUPDS IS'// &
2144  ' NOT YET ADAPTED FOR TRIPOLE GRIDS. STOPPING NOW.'
2145  CALL extcde ( 1 )
2146  END IF
2147  !
2148 #ifdef W3_T
2149  WRITE (mdst,9000) imod, nx, ny, x0, y0, sx, sy, &
2150  jmod, nxi, nyi, x0i, y0i, sxi, syi, undef
2151 #endif
2152  !
2153  ! 0.b Initialize fields
2154  !
2155  fd = undef
2156  !
2157  curvi=0
2158  IF ( grids(imod)%GTYPE .EQ. clgtype .OR. &
2159  grids(jmod)%GTYPE .EQ. clgtype ) THEN
2160  curvi=1
2161  END IF
2162 
2163  ! 1. Case of identical resolution and coinciding grids --------------- /
2164  !
2165  IF(curvi .EQ. 0) THEN
2166  IF ( abs(sx/sxi-1.) .LT. 1.e-3 .AND. &
2167  abs(sy/syi-1.) .LT. 1.e-3 .AND. &
2168  abs(mod((abs(x0-x0i))/sx+0.5,1.)-0.5) .LT. 1.e-2 .AND. &
2169  abs(mod((abs(y0-y0i))/sy+0.5,1.)-0.5) .LT. 1.e-2 ) THEN
2170  !
2171  ! 1.a Offsets
2172  !
2173  ixo = nint((x0-x0i)/sx)
2174  !
2175  IF ( flagll ) THEN
2176  ixf0 = 1
2177  ixfn = nx
2178  ixs0 = -999
2179  ixsn = -999
2180  ELSE
2181  ixf0 = max( 1 , 1-ixo )
2182  ixfn = min( nx , nxi-ixo )
2183  ixs0 = max( 1 , 1+ixo )
2184  ixsn = ixs0 + ixfn - ixf0
2185  END IF
2186  !
2187  iyo = nint((y0-y0i)/sy)
2188  !
2189  iyf0 = max( 1 , 1-iyo )
2190  iyfn = min( ny , nyi-iyo )
2191  iys0 = max( 1 , 1+iyo )
2192  iysn = iys0 + iyfn - iyf0
2193  !
2194 #ifdef W3_T
2195  WRITE (mdst,9010) ixo, iyo, ixf0, ixfn, iyf0, iyfn, &
2196  ixs0, ixsn, iys0, iysn
2197 #endif
2198  !
2199  ! 1.b Fill arrays for sea points only
2200  !
2201  IF ( flagll ) THEN
2202  DO ix=ixf0, ixfn
2203  ixs = 1 + nint( mod( &
2204  1080.+x0+(real(ix)-0.5)*sx-x0i , 360. ) / sx - 0.5 )
2205  IF ( ixs .GT. nxi ) cycle
2206  fd(ix,iyf0:iyfn) = fdi(ixs,iys0:iysn)
2207  END DO
2208  ELSE
2209  DO ix=ixf0, ixfn
2210  ixs = ix + ixo
2211  fd(ix,iyf0:iyfn) = fdi(ixs,iys0:iysn)
2212  END DO
2213  END IF
2214  !
2215  ! 1.c Return to calling routine
2216  !
2217  RETURN
2218  !
2219  END IF
2220  END IF !CURVI
2221  !
2222  ! 2. General case --------------------------------------------------- /
2223  !
2224  !
2225  ! 2.a Curvilinear grids
2226  !
2227  IF ( grids(imod)%GTYPE .EQ. clgtype .OR. &
2228  grids(jmod)%GTYPE .EQ. clgtype ) THEN
2229 
2230  ! 2.a.1 Getting the info for reg and curvi grids
2231  xgrdi => grids(jmod)%XGRD !LONS FOR INPUT FIELD
2232  ygrdi => grids(jmod)%YGRD !LATS FOR INPUT FIELD
2233 
2234  ! GETTING THE INFO FOR THE CURVILINEAR GRID
2235  xgrdc => grids(imod)%XGRD !LONS FOR CURVI GRID
2236  ygrdc => grids(imod)%YGRD !LATS FOR CURVI GRID
2237  !HPFAC => GRIDS(IMOD)%HPFAC !DELTAYGRDC(:,:)YGRDC(:,:)YGRDC(:,:)YGRDC(:,:)S IN LON FOR CURVI GRID
2238  !HQFAC => GRIDS(IMOD)%HQFAC !DELTAS IN LAT FOR CURVI GRID
2239 
2240  ! FOR NOW ONLY INTERPOLATION NOT AVERAGING THEN MXA=2
2241  mxa=2
2242  mya=2
2243  ALLOCATE ( nxa(nx,0:mxa) , rxa(nx,mxa) )
2244  nxa = 0
2245  rxa = 0.
2246  ALLOCATE ( nya(ny,0:mya) , rya(ny,mya) )
2247  nya = 0
2248  rya = 0.
2249  !
2250  !IS THE TOLERANCE USED TO DETERMINE IF TWO VALUES ARE EQUAL IN LOCATION
2251  dtoler = 1e-5
2252  ! 2.a.2 running over the curvilinear grid
2253  DO j=1,ny
2254  DO i=1,nx
2255  lonc=real(xgrdc(j,i)) !LON FOR EVERY CURVL GRID POINT
2256  latc=real(ygrdc(j,i)) !LAT FOR EVERY CURVL GRID POINT
2257  !SXC =HPFAC(J,I) !DELTA IN LON FOR CURVI GRID
2258  !SYC =HQFAC(J,I) !DELTA IN LAT FOR CURVI GRID
2259 
2260  valueinter=interpolate(nxi,real(xgrdi(1,:)),nyi,real(ygrdi(:,1)), &
2261  fdi,lonc,latc,dtoler)
2262  fd(i,j)=valueinter
2263  END DO !END I
2264  END DO !END J
2265 
2266  ELSE
2267  !
2268  ! 2.b Rectilinear grids
2269  !
2270  ! 2.b.1 Interpolation / averaging data for X axis
2271  !
2272  IF ( sx/sxi .LT. 1.0001 ) THEN
2273  mxa = 2
2274  ELSE
2275  mxa = 2 + int(sx/sxi)
2276  END IF
2277  !
2278 #ifdef W3_T
2279  WRITE (mdst,9020) 'X'
2280 #endif
2281 #ifdef W3_T1
2282  format1 = '(10X, I5, F6.2)'
2283  WRITE (format1,'(A,I2,A,I2,A)') "'(10X,",mxa+1,'I5,',mxa+1,"F6.2)'"
2284  WRITE (mdst,9021) nx, mxa
2285 #endif
2286  !
2287  ALLOCATE ( nxa(nx,0:mxa) , rxa(nx,mxa) )
2288  nxa = 0
2289  rxa = 0.
2290  !
2291  !
2292  IF ( mxa .EQ. 2 ) THEN
2293  !
2294  DO ix=1, nx
2295  IF ( flagll ) THEN
2296  xr = 1. + mod &
2297  ( 1080.+x0+real(ix-1)*sx-x0i , 360. ) / sxi
2298  ELSE
2299  xr = 1. + ( x0+real(ix-1)*sx - x0i ) / sxi
2300  END IF
2301  IF ( xr.GT.0. ) THEN
2302  j1 = int(xr)
2303  j2 = j1 + 1
2304  r2 = max( 0. , xr-real(j1) )
2305  r1 = 1. - r2
2306  IF ( flagll .AND. iclose.NE.iclose_none ) THEN
2307  j1 = 1 + mod(j1-1,nxi)
2308  j2 = 1 + mod(j2-1,nxi)
2309  END IF
2310  IF ( j1.GE.1 .AND. j1.LE.nxi .AND. r1.GT.0.05 ) THEN
2311  nxa(ix,0) = nxa(ix,0) + 1
2312  nxa(ix,nxa(ix,0)) = j1
2313  rxa(ix,nxa(ix,0)) = r1
2314  END IF
2315  IF ( j2.GE.1 .AND. j2.LE.nxi .AND. r2.GT.0.05 ) THEN
2316  nxa(ix,0) = nxa(ix,0) + 1
2317  nxa(ix,nxa(ix,0)) = j2
2318  rxa(ix,nxa(ix,0)) = r2
2319  END IF
2320  IF ( nxa(ix,0) .GT. 0 ) THEN
2321  rt = sum( rxa(ix,:) )
2322  IF ( rt .LT. 0.7 ) THEN
2323  nxa(ix,:) = 0
2324  rxa(ix,:) = 0.
2325  END IF
2326  END IF
2327  END IF
2328  END DO
2329  !
2330  ELSE
2331  !
2332  DO ix=1, nx
2333  !
2334  xfl = x0 + real(ix-1)*sx - 0.5*sx
2335  xfr = x0 + real(ix-1)*sx + 0.5*sx
2336  IF ( flagll ) THEN
2337  ixc = 1 + nint( mod( &
2338  1080.+x0+real(ix-1)*sx-x0i , 360. ) / sxi )
2339  ixs0 = ixc - 1 - mxa/2
2340  ixsn = ixc + 1 + mxa/2
2341  ELSE
2342  ixc = nint( 1. + ( x0+real(ix-1)*sx - x0i ) / sxi )
2343  ixs0 = max( 1 , ixc - 1 - mxa/2 )
2344  ixsn = min( nxi , ixc + 1 + mxa/2 )
2345  END IF
2346  DO j=ixs0, ixsn
2347  IF ( flagll ) THEN
2348  IF ( iclose.NE.iclose_none ) jj = 1 + mod(j-1+nxi,nxi)
2349  IF ( jj.LT.1 .OR. jj.GT. nxi ) cycle
2350  ixc = nint((0.5*(xfl+xfr)-x0i-real(jj-1)*sxi)/360.)
2351  IF ( ixc .NE. 0 ) THEN
2352  xfl = xfl - real(ixc) * 360.
2353  xfr = xfr - real(ixc) * 360.
2354  END IF
2355  ELSE
2356  jj = j
2357  END IF
2358  xsl = max( xfl , x0i + real(jj-1)*sxi - 0.5*sxi )
2359  xsr = min( xfr , x0i + real(jj-1)*sxi + 0.5*sxi )
2360  r1 = max( 0. , xsr - xsl ) / sx
2361  IF ( r1 .GT. 0 ) THEN
2362  nxa(ix,0) = nxa(ix,0) + 1
2363  nxa(ix,nxa(ix,0)) = jj
2364  rxa(ix,nxa(ix,0)) = r1
2365  END IF
2366  END DO
2367  IF ( nxa(ix,0) .GT. 0 ) THEN
2368  rt = sum( rxa(ix,:) )
2369  IF ( rt .LT. 0.7 ) THEN
2370  nxa(ix,:) = 0
2371  rxa(ix,:) = 0.
2372  END IF
2373  END IF
2374  END DO
2375  !
2376  END IF
2377  !
2378 #ifdef W3_T1
2379  DO, ix=1, nx
2380  IF ( nxa(ix,0) .GT. 0 ) WRITE (mdst,format1) &
2381  ix, nxa(ix,1:mxa), rxa(ix,1:mxa), sum(rxa(ix,1:mxa))
2382  END DO
2383 #endif
2384  !
2385  ! 2.b.2 Interpolation / averaging data for Y axis
2386  !
2387  IF ( sy/syi .LT. 1.0001 ) THEN
2388  mya = 2
2389  ELSE
2390  mya = 2 + int(sy/syi)
2391  END IF
2392  !
2393 #ifdef W3_T
2394  WRITE (mdst,9020) 'Y'
2395 #endif
2396 #ifdef W3_T1
2397  WRITE (format1,'(A,I2,A,I2,A)') "'(10X,",mya+1,'I5,',mya+1,"F6.2)'"
2398  WRITE (mdst,9021) ny, mya
2399 #endif
2400  !
2401  ALLOCATE ( nya(ny,0:mya) , rya(ny,mya) )
2402  nya = 0
2403  rya = 0.
2404  !
2405  !
2406  IF ( mya .EQ. 2 ) THEN
2407  !
2408  DO iy=1, ny
2409  yr = 1. + ( y0+real(iy-1)*sy - y0i ) / syi
2410  IF ( yr.GT.0. ) THEN
2411  j1 = int(yr)
2412  j2 = j1 + 1
2413  r2 = max( 0. , yr-real(j1) )
2414  r1 = 1. - r2
2415  IF ( j1.GE.1 .AND. j1.LE.nyi .AND. r1.GT.0.05 ) THEN
2416  nya(iy,0) = nya(iy,0) + 1
2417  nya(iy,nya(iy,0)) = j1
2418  rya(iy,nya(iy,0)) = r1
2419  END IF
2420  IF ( j2.GE.1 .AND. j2.LE.nyi .AND. r2.GT.0.05 ) THEN
2421  nya(iy,0) = nya(iy,0) + 1
2422  nya(iy,nya(iy,0)) = j2
2423  rya(iy,nya(iy,0)) = r2
2424  END IF
2425  IF ( nya(iy,0) .GT. 0 ) THEN
2426  rt = sum( rya(iy,:) )
2427  IF ( rt .LT. 0.7 ) THEN
2428  nya(iy,:) = 0
2429  rya(iy,:) = 0.
2430  END IF
2431  END IF
2432  END IF
2433  END DO
2434  !
2435  ELSE
2436  !
2437  DO iy=1, ny
2438  yfl = y0 + real(iy-1)*sy - 0.5*sy
2439  yfr = y0 + real(iy-1)*sy + 0.5*sy
2440  iyc = nint( 1. + ( y0+real(iy-1)*sy - y0i ) / syi )
2441  iys0 = max( 1 , iyc - 1 - mya/2 )
2442  iysn = min( nyi , iyc + 1 + mya/2 )
2443  DO j=iys0, iysn
2444  ysl = max( yfl , y0i + real(j-1)*syi - 0.5*syi )
2445  ysr = min( yfr , y0i + real(j-1)*syi + 0.5*syi )
2446  r1 = max( 0. , ysr - ysl ) / sy
2447  IF ( r1 .GT. 0 ) THEN
2448  nya(iy,0) = nya(iy,0) + 1
2449  nya(iy,nya(iy,0)) = j
2450  rya(iy,nya(iy,0)) = r1
2451  END IF
2452  END DO
2453  IF ( nya(iy,0) .GT. 0 ) THEN
2454  rt = sum( rya(iy,:) )
2455  IF ( rt .LT. 0.7 ) THEN
2456  nya(iy,:) = 0
2457  rya(iy,:) = 0.
2458  END IF
2459  END IF
2460  END DO
2461  !
2462  END IF
2463  !
2464  END IF
2465  !
2466 #ifdef W3_T1
2467  DO, iy=1, ny
2468  IF ( nya(iy,0) .GT. 0 ) WRITE (mdst,format1) &
2469  iy, nya(iy,1:mya), rya(iy,1:mya), sum(rya(iy,1:mya))
2470  END DO
2471 #endif
2472  !
2473  ! 2.c Process grid
2474  !
2475  map1 = .false.
2476  map2 = .false.
2477  !
2478  DO ix=1, nx
2479  IF ( nxa(ix,0) .EQ. 0 ) cycle
2480  DO iy=1, ny
2481  IF ( nya(iy,0) .EQ. 0 ) cycle
2482  IF ( map(iy,ix).NE.0 ) THEN
2483  fdl = 0.
2484  wtot = 0.
2485  DO j1=1, nxa(ix,0)
2486  jx = nxa(ix,j1)
2487  DO j2=1, nya(iy,0)
2488  jy = nya(iy,j2)
2489  IF ( mapi(jy,jx) .NE. 0 ) THEN
2490  wl = rxa(ix,j1) * rya(iy,j2)
2491  wtot = wtot + wl
2492  fdl = fdl + wl * fdi(jx,jy)
2493  END IF
2494  END DO
2495  END DO
2496  IF ( wtot .LT. 0.05 ) THEN
2497  map1(ix,iy) = .true.
2498  ELSE
2499  map2(ix,iy) = .true.
2500  fdl = fdl / wtot
2501  fd(ix,iy) = fdl
2502  END IF
2503  END IF
2504  END DO
2505  END DO
2506  !
2507  ! 2.d Reconcile mask differences
2508  !
2509 #ifdef W3_T
2510  WRITE (mdst,9022)
2511 #endif
2512  !
2513  jj = 0
2514  iclose => grids(imod)%ICLOSE
2515  !
2516  DO
2517  IF ( jj .GT. swpmax ) EXIT
2518  flagup = .false.
2519  map3 = .false.
2520  jj = jj + 1
2521 #ifdef W3_T
2522  WRITE (mdst,9023) jj
2523 #endif
2524  DO ix=1, nx
2525  DO iy=1, ny
2526  IF ( map1(ix,iy) ) THEN
2527  fdl = 0.
2528  j1 = 0
2529  IF ( flagll ) THEN
2530  DO j2=ix-1, ix+1
2531  IF ( (j2.GT.1 .AND. j2.LE.nx) .OR. iclose.NE.iclose_none ) THEN
2532  jx = 1 + mod(nx+j2-1,nx)
2533  DO jy=iy-1, iy+1
2534  IF ( jy.GT.1 .AND. jy.LE.ny ) THEN
2535  IF ( map2(jx,jy) ) THEN
2536  fdl = fdl + fd(jx,jy)
2537  j1 = j1 + 1
2538  END IF
2539  END IF
2540  END DO
2541  END IF
2542  END DO
2543  ELSE
2544  DO jx=ix-1, ix+1
2545  IF ( jx.GT.1 .AND. jx.LE.nx ) THEN
2546  DO jy=iy-1, iy+1
2547  IF ( jy.GT.1 .AND. jy.LE.ny ) THEN
2548  IF ( map2(jx,jy) ) THEN
2549  fdl = fdl + fd(jx,jy)
2550  j1 = j1 + 1
2551  END IF
2552  END IF
2553  END DO
2554  END IF
2555  END DO
2556  END IF !FLAGLL
2557  IF ( j1 .GT. 0 ) THEN
2558  fd(ix,iy) = fdl / real(j1)
2559  map1(ix,iy) = .false.
2560  map3(ix,iy) = .true.
2561  flagup = .true.
2562  END IF
2563  END IF
2564  END DO
2565  END DO
2566  IF ( flagup ) THEN
2567  map2 = map2 .OR. map3
2568  ELSE
2569  EXIT
2570  END IF
2571  END DO
2572  !
2573  ! 3. End of routine -------------------------------------------------- /
2574  !
2575  DEALLOCATE ( nxa, nya, rxa, rya )
2576  !
2577  RETURN
2578  !
2579  ! Formats
2580  !
2581 #ifdef W3_T
2582 9000 FORMAT ( ' TEST WMUPDS : GRID INFORMATION : '/ &
2583  ' ',3i5,4e11.3/ &
2584  ' ',3i5,4e11.3/ &
2585  ' UNDEFINED = ',e10.3)
2586 9010 FORMAT ( ' TEST WMUPDS : COINCIDING GRIDS, OFFSETS :',2i6/ &
2587  ' TARGET GRID RANGES :',4i6/ &
2588  ' SOURCE GRID RANGES :',4i6)
2589 9020 FORMAT ( ' TEST WMUPDS : WEIGHTS FOR ',a,' INTERPOATION')
2590 #endif
2591 #ifdef W3_T1
2592 9021 FORMAT ( ' TEST WMUPDS : ARAY DIMENSIONED AS : ',2i6)
2593 #endif
2594 #ifdef W3_T
2595 9022 FORMAT ( ' TEST WMUPDS : RECONCILING MASKS')
2596 9023 FORMAT ( ' SWEEP NR ',i4)
2597 #endif
2598  !/
2599  !/ End of WMUPDS ----------------------------------------------------- /
2600  !/

References w3gdatmd::clgtype, w3servmd::extcde(), w3gdatmd::flagll, w3gdatmd::grids, w3gdatmd::gtype, w3gdatmd::hpfac, w3gdatmd::hqfac, w3gdatmd::iclose_none, w3gdatmd::iclose_smpl, w3gdatmd::iclose_trpl, wmmdatmd::improc, interpolate(), wmmdatmd::mdse, wmmdatmd::mdso, wmmdatmd::mdss, wmmdatmd::mdst, wmmdatmd::nmperr, wmmdatmd::nmpscr, w3gdatmd::nx, w3gdatmd::ny, w3gdatmd::rlgtype, w3servmd::strace(), swpmax, w3gdatmd::sx, w3gdatmd::sy, w3gdatmd::ungtype, w3gdatmd::x0, w3gdatmd::xgrd, w3gdatmd::y0, and w3gdatmd::ygrd.

Referenced by wmupd2().

◆ wmupdt()

subroutine wmupdtmd::wmupdt ( integer, intent(in)  IMOD,
integer, dimension(2), intent(inout)  TDATA 
)

Update inputs for selected wave model grid.

Reading from native grid files if update is needed based on time of data.

Parameters
[in]IMODModel number
[in,out]TDATATime for which all is data available.
Author
H. L. Tolman
Date
22-Mar-2021

Definition at line 103 of file wmupdtmd.F90.

103  !/
104  !/ +-----------------------------------+
105  !/ | WAVEWATCH III NOAA/NCEP |
106  !/ | H. L. Tolman |
107  !/ | FORTRAN 90 |
108  !/ | Last update : 22-Mar-2021 |
109  !/ +-----------------------------------+
110  !/
111  !/ 22-Feb-2005 : Origination. ( version 3.07 )
112  !/ 14-Oct-2006 : Adding separate input grids. ( version 3.10 )
113  !/ 12-Jan-2007 : General clean-up and bug fixes. ( version 3.10 )
114  !/ 20-Jan-2017 : Enable using input from coupler ( version 6.02 )
115  !/ (T. J. Campbell, NRL)
116  !/ 22-Mar-2021 : Add momentum and air density input ( version 7.13 )
117  !/
118  ! 1. Purpose :
119  !
120  ! Update inputs for seleceted wave model grid.
121  !
122  ! 2. Method :
123  !
124  ! Reading from native grid files if update is needed based on
125  ! time of data.
126  !
127  ! 3. Parameters :
128  !
129  ! Parameter list
130  ! ----------------------------------------------------------------
131  ! IMOD Int. I Model number,
132  ! TDATA I.A. I Time for which all is data available.
133  ! ----------------------------------------------------------------
134  !
135  ! 4. Subroutines used :
136  !
137  ! Name Type Module Description
138  ! ----------------------------------------------------------------
139  ! W3SETG Subr. W3GDATMD Point to grid/model.
140  ! W3SETW Subr. W3WDATMD Point to grid/model.
141  ! W3SETI Subr. W3IDATMD Point to grid/model.
142  ! WMSETM Subr. WMMDATMD Point to grid/model.
143  ! STRACE Subr. W3ERVMD Subroutine tracing.
144  ! EXTCDE Subr. Id. Program abort.
145  ! DSEC21 Func. W3TIMEMD Time difference.
146  ! STME21 Subr. Id. Write time string.
147  ! TICK21 Subr. Id. Advancing time.
148  ! ----------------------------------------------------------------
149  !
150  ! 5. Called by :
151  !
152  ! Name Type Module Description
153  ! ----------------------------------------------------------------
154  ! WMWAVE Subr. WMWAVEMD Multi-grid model main routine.
155  ! ----------------------------------------------------------------
156  !
157  ! 6. Error messages :
158  !
159  ! 7. Remarks :
160  !
161  ! - IDFLDS dimensioning is hardwired as IDFLDS(-7:9) where
162  ! lowest possible value of JFIRST is JFIRST=-7
163  !
164  ! 8. Structure :
165  !
166  ! See source code.
167  !
168  ! 9. Switches :
169  !
170  ! !/S Enable subroutine tracing.
171  ! !/T Enable test output
172  !
173  ! 10. Source code :
174  !
175  !/ ------------------------------------------------------------------- /
176  !/
177  USE w3gdatmd, ONLY: w3setg
178  USE w3wdatmd, ONLY: w3setw
179  USE w3idatmd, ONLY: w3seti
180  USE wmmdatmd, ONLY: wmsetm
181  USE w3servmd, ONLY: extcde
182 #ifdef W3_S
183  USE w3servmd, ONLY: strace
184 #endif
185  USE w3timemd, ONLY: dsec21, stme21, tick21
186  !/
187  USE w3gdatmd, ONLY: nx, ny, filext
188  USE w3wdatmd, ONLY: time
189 
190  USE w3idatmd, ONLY: inflags1, tln, tc0, tcn, tw0, twn, tu0, &
191  tun, tr0, trn, tin, t0n, t1n, t2n, tg0, &
192  tgn, tfn, tdn, ttn, tvn, tzn, ti1, ti2, &
193  ti3, ti4, ti5, jfirst
194 
195  USE wmmdatmd, ONLY: improc, mdso, mdss, mdst, mdse, nmpscr, &
198  !/
199  IMPLICIT NONE
200  !/
201  !/ ------------------------------------------------------------------- /
202  !/ Parameter list
203  !/
204  INTEGER, INTENT(IN) :: IMOD
205  INTEGER, INTENT(INOUT) :: TDATA(2)
206  !/
207  !/ ------------------------------------------------------------------- /
208  !/ Local parameters
209  !/
210  INTEGER :: MDSEN, J, DTIME(2), IERR, NDTNEW, JJ
211 #ifdef W3_S
212  INTEGER, SAVE :: IENT = 0
213 #endif
214  REAL :: DTTST
215  LOGICAL :: FIRST
216  CHARACTER(LEN=13) :: IDFLDS(-7:10)
217  CHARACTER(LEN=23) :: DTME21
218  !
219  DATA idflds / 'ice param. 1 ' , 'ice param. 2 ' , &
220  'ice param. 3 ' , 'ice param. 4 ' , &
221  'ice param. 5 ' , &
222  'mud density ' , 'mud thkness ' , &
223  'mud viscos. ' , &
224  'water levels ' , 'currents ' , &
225  'winds ' , 'ice fields ' , &
226  'momentum ' , 'air density ' , &
227  'mean param. ' , '1D spectra ' , &
228  '2D spectra ' , 'grid speed ' /
229  !/
230  !/ ------------------------------------------------------------------- /
231  ! 0. Initialization
232  ! 0.a Subroutine tracing and echo of input
233  !
234 #ifdef W3_S
235  CALL strace (ient, 'WMUPDT')
236 #endif
237 #ifdef W3_T
238  WRITE (mdst,9000) imod, tdata
239 #endif
240  !
241  IF ( improc .EQ. nmperr ) THEN
242  mdsen = mdse
243  ELSE
244  mdsen = -1
245  END IF
246  !
247  ! 0.b Point to proper grids and initialize
248  !
249  CALL w3setg ( imod, mdse, mdst )
250  CALL w3setw ( imod, mdse, mdst )
251  CALL w3seti ( imod, mdse, mdst )
252  CALL wmsetm ( imod, mdse, mdst )
253  !
254  fllstl = .false.
255  fllsti = .false.
256  fllstr = .false.
257  ierr = 0
258  !
259  ! 0.c Output
260  !
261  IF ( mdss.NE.mdso .AND. nmpscr.EQ.improc ) THEN
262  CALL stme21 ( time , dtme21 )
263  WRITE (mdss,900) imod, dtme21
264  END IF
265  !
266 #ifdef W3_T
267  WRITE (mdst,9001) ' J', '0-N', time, etime
268  IF (lbound(idinp,2).LE.-7) WRITE (mdst,9002) -7, idinp(imod,-7), inflags1(-7), ti1
269  IF (lbound(idinp,2).LE.-6) WRITE (mdst,9002) -6, idinp(imod,-6), inflags1(-6), ti2
270  IF (lbound(idinp,2).LE.-5) WRITE (mdst,9002) -5, idinp(imod,-5), inflags1(-5), ti3
271  IF (lbound(idinp,2).LE.-4) WRITE (mdst,9002) -4, idinp(imod,-4), inflags1(-4), ti4
272  IF (lbound(idinp,2).LE.-3) WRITE (mdst,9002) -3, idinp(imod,-3), inflags1(-3), ti5
273  IF (lbound(idinp,2).LE.-2) WRITE (mdst,9002) -2, idinp(imod,-2), inflags1(-2), tzn
274  IF (lbound(idinp,2).LE.-1) WRITE (mdst,9002) -1, idinp(imod,-1), inflags1(-1), ttn
275  IF (lbound(idinp,2).LE. 0) WRITE (mdst,9002) 0, idinp(imod, 0), inflags1( 0), tvn
276  WRITE (mdst,9002) 1, idinp(imod,1), inflags1(1), tln
277  WRITE (mdst,9003) 2, idinp(imod,2), inflags1(2), tc0, tcn
278  WRITE (mdst,9003) 3, idinp(imod,3), inflags1(3), tw0, twn
279  WRITE (mdst,9002) 4, idinp(imod,4), inflags1(4), tin
280  WRITE (mdst,9003) 5, idinp(imod,5), inflags1(5), tu0, tun
281  WRITE (mdst,9003) 6, idinp(imod,6), inflags1(6), tr0, trn
282  WRITE (mdst,9002) 7, idinp(imod,7), inflags1(7), t0n
283  WRITE (mdst,9002) 8, idinp(imod,8), inflags1(8), t1n
284  WRITE (mdst,9002) 9, idinp(imod,9), inflags1(9), t2n
285  WRITE (mdst,9003) 10, 'MOV' , inflags1(10), tg0, tgn
286  WRITE (mdst,9004) 'GRD', nx, ny
287 #endif
288  !
289  ! 1. Loop over input types ------------------------------------------ /
290  !
291  DO j=jfirst, 10
292  !
293  ! 1.a Check if update needed
294  !
295  IF ( .NOT. inflags1(j) ) cycle
296  !
297 #ifdef W3_T
298  WRITE (mdst,9010) j, inflags1(j), inpmap(imod,j)
299 #endif
300  !
301  ! 1.b Test time
302  !
303  IF ( tfn(1,j) .EQ. -1 ) THEN
304  first = .true.
305  dttst = 0.
306  ELSE
307  first = .false.
308  dttst = dsec21( time , tfn(:,j) )
309  END IF
310  !
311 #ifdef W3_T
312  WRITE (mdst,9011) idinp(imod,j), dttst, tfn(:,j)
313 #endif
314  !
315  IF ( dttst .GT. 0. ) cycle
316  IF ( mdss.NE.mdso .AND. nmpscr.EQ.improc ) &
317  WRITE (mdss,901) idflds(j)
318  !
319  ! 2. Forcing input from file & defined on the native grid ----------- /
320  !
321  IF ( inpmap(imod,j) .EQ. 0 ) THEN
322  !
323 #ifdef W3_T
324  WRITE (mdst,9020)
325 #endif
326  !
327  CALL wmupd1 ( imod, idinp(imod,j), j, ierr )
328  !
329  ! 3. Forcing input from file & defined on an input grid ------------- /
330  !
331  ELSE IF ( inpmap(imod,j) .GT. 0 ) THEN
332  !
333 #ifdef W3_T
334  WRITE (mdst,9030) inpmap(imod,j)
335 #endif
336  !
337  ! 3.a Check if input grid is available
338  !
339  jj = -inpmap(imod,j)
340  CALL w3setg ( jj, mdse, mdst )
341  CALL w3seti ( jj, mdse, mdst )
342  !
343  IF ( tfn(1,j) .EQ. -1 ) THEN
344  dttst = 0.
345  ELSE
346  IF ( first .OR. ( j.EQ.1 .AND. iflstl(-jj) ) &
347  .OR. ( j.EQ.4 .AND. iflsti(-jj) ) &
348  .OR. ( j.EQ.6 .AND. iflstr(-jj) ) ) THEN
349  dttst = 1.
350  ELSE
351  dttst = dsec21( time , tfn(:,j) )
352  END IF
353  END IF
354  !
355  IF ( j .EQ. 1 ) fllstl = iflstl(-jj)
356  IF ( j .EQ. 4 ) fllsti = iflsti(-jj)
357  IF ( j .EQ. 6 ) fllstr = iflstr(-jj)
358  !
359 #ifdef W3_T
360  WRITE (mdst,9031) j, idinp(jj,j), dttst, tfn(:,j)
361 #endif
362  !
363  ! 3.b If needed, update input grid
364  ! Note: flags in WMMDATMD set for grid IMOD !
365  !
366  IF ( dttst .LE. 0. ) THEN
367  !
368  IF ( mdss.NE.mdso .AND. nmpscr.EQ.improc ) &
369  WRITE (mdss,930) filext
370  !
371  CALL wmupd1 ( jj, idinp(jj,j), j, ierr )
372  !
373  IF ( j .EQ. 1 ) iflstl(-jj) = fllstl
374  IF ( j .EQ. 4 ) iflsti(-jj) = fllsti
375  IF ( j .EQ. 6 ) iflstr(-jj) = fllstr
376  !
377  END IF
378  !
379  ! 3.c Set up for update, and call updating routine
380  !
381  CALL w3setg ( imod, mdse, mdst )
382  CALL w3seti ( imod, mdse, mdst )
383  !
384  CALL wmupd2 ( imod, j, jj, ierr )
385  !
386  ! 4. Forcing input from CPL ----------------------------------------- /
387  !
388  ELSE ! INPMAP(IMOD,J) .LT. 0
389  ! Data input and time stamp settings for forcing input from
390  ! CPL are handled in wmesmfmd.ftn:GetImport
391  !
392 #ifdef W3_T
393  IF ( inpmap(imod,j) .EQ. -999 ) THEN
394  ! *** Forcing input from CPL & defined on native grid ***
395  WRITE (mdst,9040)
396  ELSE
397  ! *** Forcing input from CPL & defined on an input grid ***
398  WRITE (mdst,9050) -inpmap(imod,j)
399  END IF
400 #endif
401  !
402  END IF
403  !
404  ! 5. Finalize for each type ----------------------------------------- /
405  ! 5.a Process IERR output
406  !
407  IF ( ierr.GT.0 ) GOTO 2000
408  IF ( ierr.LT.0 .AND. mdss.NE.mdso .AND. nmpscr.EQ.improc ) &
409  WRITE (mdss,950) idflds(j)
410  !
411  ! 5.b End of master loop
412  !
413  END DO
414  !
415  ! 6. Compute TDATA -------------------------------------------------- /
416  !
417  tdata = etime
418  !
419  DO j=jfirst, 10
420  IF ( .NOT. inflags1(j) ) cycle
421  dttst = dsec21( tfn(:,j) , tdata )
422  IF ( dttst.GT.0. .AND. .NOT. ( (fllstl .AND. j.EQ.1) .OR. &
423  (fllsti .AND. j.EQ.4) .OR. &
424  (fllstr .AND. j.EQ.6) ) ) THEN
425  tdata = tfn(:,j)
426  END IF
427  END DO
428  !
429  ! 6. Compute TDN ---------------------------------------------------- /
430  !
431  tdn = tdata
432  CALL tick21 ( tdn, 1. )
433  DO j=7, 9
434  IF ( inflags1(j) ) THEN
435  dttst = dsec21( tfn(:,j) , tdn )
436  IF ( dttst.GT.0. ) tdn = tfn(:,j)
437  END IF
438  END DO
439  !
440  ! 7. Final test output ---------------------------------------------- /
441  !
442 #ifdef W3_T
443  WRITE (mdst,9070) ' J', '0-N', time, etime, tdata
444  IF (lbound(idinp,2).LE.-7) WRITE (mdst,9071) -7, idinp(imod,-7), inflags1(-7), ti1
445  IF (lbound(idinp,2).LE.-6) WRITE (mdst,9071) -6, idinp(imod,-6), inflags1(-6), ti2
446  IF (lbound(idinp,2).LE.-5) WRITE (mdst,9071) -5, idinp(imod,-5), inflags1(-5), ti3
447  IF (lbound(idinp,2).LE.-4) WRITE (mdst,9071) -4, idinp(imod,-4), inflags1(-4), ti4
448  IF (lbound(idinp,2).LE.-3) WRITE (mdst,9071) -3, idinp(imod,-3), inflags1(-3), ti5
449  IF (lbound(idinp,2).LE.-2) WRITE (mdst,9071) -2, idinp(imod,-2), inflags1(-2), tzn
450  IF (lbound(idinp,2).LE.-1) WRITE (mdst,9071) -1, idinp(imod,-1), inflags1(-1), ttn
451  IF (lbound(idinp,2).LE. 0) WRITE (mdst,9071) 0, idinp(imod, 0), inflags1( 0), tvn
452  WRITE (mdst,9071) 1, idinp(imod,1), inflags1(1), tln
453  WRITE (mdst,9072) 2, idinp(imod,2), inflags1(2), tc0, tcn
454  WRITE (mdst,9072) 3, idinp(imod,3), inflags1(3), tw0, twn
455  WRITE (mdst,9071) 4, idinp(imod,4), inflags1(4), tin
456  WRITE (mdst,9072) 5, idinp(imod,5), inflags1(5), tu0, tun
457  WRITE (mdst,9072) 6, idinp(imod,6), inflags1(6), tr0, trn
458  WRITE (mdst,9071) 7, idinp(imod,7), inflags1(7), t0n
459  WRITE (mdst,9071) 8, idinp(imod,8), inflags1(8), t1n
460  WRITE (mdst,9073) 9, idinp(imod,9), inflags1(9), t2n, tdn
461  WRITE (mdst,9072) 10, 'MOV' , inflags1(10), tg0, tgn
462 #endif
463  !
464  RETURN
465  !
466  ! Error escape locations
467  !
468 2000 CONTINUE
469  CALL extcde ( 2000 )
470  RETURN
471  !
472  ! Formats
473  !
474 900 FORMAT ( ' Updating input for grid',i3,' at ',a)
475 901 FORMAT ( ' Updating ',a)
476 930 FORMAT ( ' First updating ',a)
477 950 FORMAT ( ' Past last ',a)
478  !
479 #ifdef W3_T
480 9000 FORMAT ( ' TEST WMUPDT : INPUT : ',i4,i10.8,i7.6, &
481  ' <============================')
482 9001 FORMAT ( ' TEST WMUPDT : ',a2,1x,a3,3x, 2(i10.8,i7.6))
483 9002 FORMAT ( ' ',i2,1x,a3,l3,17x,1(i10.8,i7.6))
484 9003 FORMAT ( ' ',i2,1x,a3,l3, 2(i10.8,i7.6))
485 9004 FORMAT ( ' ',2x,1x,a3,3x,2i10 )
486 9010 FORMAT ( ' TEST WMUPDT : J, FLAG, INPMAP : ',i2,l2,i4)
487 9011 FORMAT ( ' TEST WMUPDT : ',a,', DTTST = ',e10.3,2x,i9.8,i7.6)
488 9020 FORMAT ( ' TEST WMUPDT : FORCING INPUT FROM FILE & DEFINED ON THE NATIVE GRID')
489 9030 FORMAT ( ' TEST WMUPDT : FORCING INPUT FROM FILE & DEFINED ON INPUT GRID',i4)
490 9031 FORMAT ( ' TEST WMUPDT : J =',i4,3xa,', DTTST = ', &
491  e10.3,2x,i9.8,i7.6)
492 9040 FORMAT ( ' TEST WMUPDT : FORCING INPUT FROM CPL & DEFINED ON THE NATIVE GRID')
493 9050 FORMAT ( ' TEST WMUPDT : FORCING INPUT FROM CPL & DEFINED ON INPUT GRID',i4)
494 9070 FORMAT ( ' TEST WMUPDT : ',a2,1x,a3,3x, 3(i10.8,i7.6))
495 9071 FORMAT ( ' ',i2,1x,a3,l3,17x,1(i10.8,i7.6))
496 9072 FORMAT ( ' ',i2,1x,a3,l3, 2(i10.8,i7.6))
497 9073 FORMAT ( ' ',i2,1x,a3,l3,17x,2(i10.8,i7.6))
498 #endif
499  !/
500  !/ End of WMUPDT ----------------------------------------------------- /
501  !/

References w3timemd::dsec21(), wmmdatmd::etime, w3servmd::extcde(), w3gdatmd::filext, wmmdatmd::fllsti, wmmdatmd::fllstl, wmmdatmd::fllstr, wmmdatmd::idinp, wmmdatmd::iflsti, wmmdatmd::iflstl, wmmdatmd::iflstr, wmmdatmd::improc, w3idatmd::inflags1, wmmdatmd::inpmap, w3idatmd::jfirst, wmmdatmd::mdse, wmmdatmd::mdso, wmmdatmd::mdss, wmmdatmd::mdst, wmmdatmd::nmperr, wmmdatmd::nmpscr, w3gdatmd::nx, w3gdatmd::ny, w3timemd::stme21(), w3servmd::strace(), w3idatmd::t0n, w3idatmd::t1n, w3idatmd::t2n, w3idatmd::tc0, w3idatmd::tcn, w3idatmd::tdn, w3idatmd::tfn, w3idatmd::tg0, w3idatmd::tgn, w3idatmd::ti1, w3idatmd::ti2, w3idatmd::ti3, w3idatmd::ti4, w3idatmd::ti5, w3timemd::tick21(), w3wdatmd::time, w3idatmd::tin, w3idatmd::tln, w3idatmd::tr0, w3idatmd::trn, w3idatmd::ttn, w3idatmd::tu0, w3idatmd::tun, w3idatmd::tvn, w3idatmd::tw0, w3idatmd::twn, w3idatmd::tzn, w3gdatmd::w3setg(), w3idatmd::w3seti(), w3wdatmd::w3setw(), wmmdatmd::wmsetm(), wmupd1(), and wmupd2().

Referenced by wmwavemd::wmwave().

◆ wmupdv()

subroutine wmupdtmd::wmupdv ( integer, intent(in)  IMOD,
real, dimension(nx,ny), intent(out)  VX,
real, dimension(nx,ny), intent(out)  VY,
integer, intent(in)  JMOD,
real, dimension(grids(jmod)%nx,grids(jmod)%ny), intent(in)  VXI,
real, dimension(grids(jmod)%nx,grids(jmod)%ny), intent(in)  VYI,
real, intent(in)  UNDEF,
integer, intent(in)  CONSTP 
)

Interpolate vector field from input grid to model grid.

Interpolating or averaging from input grid.

Parameters
[in]IMODOutput model number
[out]VXOutput vector field
[out]VYOutput vector field
[in]JMODInput model number
[in]VXIInput vector field
[in]VYIInput vector field
[in]UNDEFValue for mapped out point and points not covered.
[in]CONSTPConvervation type
Author
H. L. Tolman
Date
06-Dec-2010

Definition at line 1305 of file wmupdtmd.F90.

1305  !/
1306  !/ +-----------------------------------+
1307  !/ | WAVEWATCH III NOAA/NCEP |
1308  !/ | H. L. Tolman |
1309  !/ | FORTRAN 90 |
1310  !/ | Last update : 06-Dec-2010 |
1311  !/ +-----------------------------------+
1312  !/
1313  !/ 14-Oct-2006 : Origination. ( version 3.10 )
1314  !/ 12-Jan-2007 : General clean-up and bug fixes. ( version 3.10 )
1315  !/ 30-Oct-2009 : Implement run-time grid selection. ( version 3.14 )
1316  !/ (W. E. Rogers & T. J. Campbell, NRL)
1317  !/ 06-Dec-2010 : Change from GLOBAL (logical) to ICLOSE (integer) to
1318  !/ specify index closure for a grid. ( version 3.14 )
1319  !/ (T. J. Campbell, NRL)
1320  !/ 01-Jul-2019 : Generalize output to curv grids ( version 7.13 )
1321  !/ (R. Padilla-Hernandez, J.H. Alves, EMC/NOAA)
1322  !/
1323  ! 1. Purpose :
1324  !
1325  ! Interpolate vector field from input grid to model grid.
1326  !
1327  ! 2. Method :
1328  !
1329  ! Interpolating or averaging from input grid.
1330  !
1331  ! 3. Parameters :
1332  !
1333  ! Parameter list
1334  ! ----------------------------------------------------------------
1335  ! IMOD Int. I Output model number,
1336  ! VX/Y Int. O Output vector field.
1337  ! JMOD Int. I Input model number,
1338  ! VX/YI Int. I Input vector field.
1339  ! UNDEF Int. I Value for mapped out point and points not
1340  ! covered.
1341  ! CONSTP Int. I Conservation type :
1342  ! 1: Vector speed.
1343  ! 2: Vector speed squared.
1344  ! *: Vector components.
1345  ! ----------------------------------------------------------------
1346  !
1347  ! 4. Subroutines used :
1348  !
1349  ! Name Type Module Description
1350  ! ----------------------------------------------------------------
1351  ! STRACE Subr. W3ERVMD Subroutine tracing.
1352  ! EXTCDE Subr. Id. Program abort.
1353  ! ----------------------------------------------------------------
1354  !
1355  ! 5. Called by :
1356  !
1357  ! Name Type Module Description
1358  ! ----------------------------------------------------------------
1359  ! WMUPD2 Subr. WMUPDTMD Input update routine.
1360  ! ----------------------------------------------------------------
1361  !
1362  ! 6. Error messages :
1363  !
1364  ! 7. Remarks :
1365  !
1366  ! - Grid pointers for output grid need to be set externally.
1367  ! - If input grid does not cover point of target grid, target grid
1368  ! values are set to UNDEF.
1369  !
1370  ! 8. Structure :
1371  !
1372  ! See source code.
1373  !
1374  ! 9. Switches :
1375  !
1376  ! !/S Enable subroutine tracing.
1377  ! !/T Enable test output
1378  ! !/T1 Test output interpolation data.
1379  !
1380  ! 10. Source code :
1381  !
1382  !/ ------------------------------------------------------------------- /
1383  !/
1384  USE w3servmd, ONLY: extcde
1385 #ifdef W3_S
1386  USE w3servmd, ONLY: strace
1387 #endif
1388  !/
1389  USE w3gdatmd, ONLY: nx, ny, x0, y0, sx, sy, grids, flagll, &
1392  hpfac, hqfac, xgrd, ygrd
1393  USE wmmdatmd, ONLY: improc, nmperr, nmpscr, mdst, mdse, mdso, &
1394  mdss
1395  !/
1396  IMPLICIT NONE
1397  !/
1398  !/ ------------------------------------------------------------------- /
1399  !/ Parameter list
1400  !/
1401  INTEGER, INTENT(IN) :: IMOD, JMOD, CONSTP
1402  REAL, INTENT(OUT) :: VX(NX,NY), VY(NX,NY)
1403  REAL, INTENT(IN) :: VXI(GRIDS(JMOD)%NX,GRIDS(JMOD)%NY), &
1404  VYI(GRIDS(JMOD)%NX,GRIDS(JMOD)%NY), &
1405  UNDEF
1406  !/
1407  !/ ------------------------------------------------------------------- /
1408  !/ Local parameters
1409  !/
1410  INTEGER :: IXO, IYO, IX, IY, IXF0, IXFN, IYF0, &
1411  IYFN, IXS0, IXSN, IYS0, IYSN, IXS, &
1412  MXA, MYA, J, J1, J2, IXC, IYC, JJ, &
1413  JX, JY
1414  INTEGER :: NPOIX, NPOIY, I, IFIELDS,CURVI !RP
1415 
1416 #ifdef W3_S
1417  INTEGER, SAVE :: IENT = 0
1418 #endif
1419  INTEGER, ALLOCATABLE :: NXA(:,:), NYA(:,:)
1420  REAL :: XR, YR, R1, R2, RT, XFL, XFR, XSL, &
1421  XSR, YFL, YFR, YSL, YSR
1422  REAL :: VXL, VYL, VA0, VA, VA2, FACTOR, &
1423  WTOT, WL
1424 
1425  REAL :: LONC, LATC, SXYC, &
1426  XDI, DTOLER, VALUEX, VALUEY
1427 
1428  REAL, ALLOCATABLE :: RXA(:,:), RYA(:,:)
1429 
1430  ! REAL, ALLOCATABLE :: VARIN(:,:) !RP
1431 
1432  LOGICAL :: MAP1(NX,NY), MAP2(NX,NY), &
1433  MAP3(NX,NY), FLAGUP
1434  !
1435  INTEGER, POINTER :: NXI, NYI, MAP(:,:), MAPI(:,:)
1436  REAL, POINTER :: X0I, Y0I, SXI, SYI !RP , HPFACI, HQFACI
1437 
1438  REAL, POINTER :: HPFACI(:,:), HQFACI(:,:)
1439  DOUBLE PRECISION, POINTER :: XGRDI(:,:), YGRDI(:,:), XGRDC(:,:), YGRDC(:,:)
1440 
1441  INTEGER, POINTER :: ICLOSE
1442  REAL, ALLOCATABLE :: XGRTMP(:),YGRTMP(:)
1443 #ifdef W3_T1
1444  CHARACTER(LEN=17) :: FORMAT1
1445 #endif
1446  !/
1447  !/ ------------------------------------------------------------------- /
1448  ! 0. Initialization
1449  ! 0.a Subroutine tracing and test output
1450  !
1451 #ifdef W3_S
1452  CALL strace (ient, 'WMUPDV')
1453 #endif
1454  !
1455  IF ( grids(imod)%GTYPE .EQ. ungtype .OR. &
1456  grids(jmod)%GTYPE .EQ. ungtype ) THEN
1457  WRITE (mdse,'(/2A)') ' *** ERROR WMUPDV: ', &
1458  'UNSTRUCTURED GRID SUPPORT NOT YET IMPLEMENTED ***'
1459  CALL extcde ( 999 )
1460  END IF
1461  !
1462  nxi => grids(jmod)%NX
1463  nyi => grids(jmod)%NY
1464  x0i => grids(jmod)%X0
1465  y0i => grids(jmod)%Y0
1466  sxi => grids(jmod)%SX
1467  syi => grids(jmod)%SY
1468  hpfaci => grids(jmod)%HPFAC
1469  hqfaci => grids(jmod)%HQFAC
1470  map => grids(imod)%MAPSTA
1471  mapi => grids(jmod)%MAPSTA
1472  iclose => grids(jmod)%ICLOSE
1473  !
1474  IF ( iclose .EQ. iclose_trpl ) THEN
1475  IF ( improc.EQ.nmperr ) WRITE(mdse,*)'SUBROUTINE WMUPDV IS'// &
1476  ' NOT YET ADAPTED FOR TRIPOLE GRIDS. STOPPING NOW.'
1477  CALL extcde ( 1 )
1478  END IF
1479 
1480  !
1481 #ifdef W3_T
1482  WRITE (mdst,9000) imod, nx, ny, x0, y0, sx, sy, &
1483  jmod, nxi, nyi, x0i, y0i, sxi, syi, undef
1484 #endif
1485  !
1486  ! 0.b Initialize fields
1487  !
1488  vx = undef
1489  vy = undef
1490  !
1491  curvi=0
1492  IF ( grids(imod)%GTYPE .EQ. clgtype .OR. &
1493  grids(jmod)%GTYPE .EQ. clgtype ) THEN
1494  curvi=1
1495  END IF
1496 
1497  ! 1. Case of identical resolution and coinciding grids --------------- /
1498  IF(curvi .EQ. 0) THEN
1499  !
1500  IF ( abs(sx/sxi-1.) .LT. 1.e-3 .AND. &
1501  abs(sy/syi-1.) .LT. 1.e-3 .AND. &
1502  abs(mod((abs(x0-x0i))/sx+0.5,1.)-0.5) .LT. 1.e-2 .AND. &
1503  abs(mod((abs(y0-y0i))/sy+0.5,1.)-0.5) .LT. 1.e-2 ) THEN
1504  !
1505  ! 1.a Offsets
1506  !
1507 
1508  ixo = nint((x0-x0i)/sx)
1509  !
1510  IF ( flagll ) THEN
1511  ixf0 = 1
1512  ixfn = nx
1513  ixs0 = -999
1514  ixsn = -999
1515  ELSE
1516  ixf0 = max( 1 , 1-ixo )
1517  ixfn = min( nx , nxi-ixo )
1518  ixs0 = max( 1 , 1+ixo )
1519  ixsn = ixs0 + ixfn - ixf0
1520  END IF
1521  !
1522  iyo = nint((y0-y0i)/sy)
1523  !
1524  iyf0 = max( 1 , 1-iyo )
1525  iyfn = min( ny , nyi-iyo )
1526  iys0 = max( 1 , 1+iyo )
1527  iysn = iys0 + iyfn - iyf0
1528  !
1529 #ifdef W3_T
1530  WRITE (mdst,9010) ixo, iyo, ixf0, ixfn, iyf0, iyfn, &
1531  ixs0, ixsn, iys0, iysn
1532 #endif
1533  !
1534  ! 1.b Fill arrays for sea points only
1535  !
1536 
1537  DO ix=ixf0, ixfn
1538  IF ( flagll ) THEN
1539  ixs = 1 + nint( mod( &
1540  1080.+x0+(real(ix)-0.5)*sx-x0i , 360. ) / sx - 0.5 )
1541  IF ( ixs .GT. nxi ) cycle
1542  ELSE
1543  ixs = ix + ixo
1544  END IF
1545  vx(ix,iyf0:iyfn) = vxi(ixs,iys0:iysn)
1546  vy(ix,iyf0:iyfn) = vyi(ixs,iys0:iysn)
1547  END DO
1548  !
1549  ! 1.c Return to calling routine
1550  !
1551  RETURN
1552  !
1553  END IF
1554  END IF !CURVI
1555  !
1556  ! 2. General case --------------------------------------------------- /
1557  !
1558  ! 2.a Curvilinear grids
1559  !
1560  IF ( grids(imod)%GTYPE .EQ. clgtype .OR. &
1561  grids(jmod)%GTYPE .EQ. clgtype ) THEN
1562 
1563  xgrdi => grids(jmod)%XGRD !LONS FOR INPUT FIELD
1564  ygrdi => grids(jmod)%YGRD !LATS FOR INPUT FIELD
1565 
1566  ! GETTING THE INFO FOR THE CURVILINEAR GRID
1567  xgrdc => grids(imod)%XGRD !LONS FOR CURVI GRID
1568  ygrdc => grids(imod)%YGRD !LATS FOR CURVI GRID
1569  !HPFAC => GRIDS(IMOD)%HPFAC !DELTAS IN LON FOR CURVI GRID
1570  !HQFAC => GRIDS(IMOD)%HQFAC !DELTAS IN LAT FOR CURVI GRID
1571  !
1572  !
1573  ! FOR NOW ONLY INTERPOLATION NOT AVERAGING THEN MXA=2
1574  mxa=2
1575  mya=2
1576  ALLOCATE ( nxa(nx,0:mxa) , rxa(nx,mxa) )
1577  nxa = 0
1578  rxa = 0.
1579  ALLOCATE ( nya(ny,0:mya) , rya(ny,mya) )
1580  nya = 0
1581  rya = 0.
1582 
1583  !IS THE TOLERANCE USED TO DETERMINE IF TWO VALUES ARE EQUAL IN LOCATION
1584  dtoler = 1e-5
1585  ! 2.a.1 running over the curvilinear grid
1586  ALLOCATE (xgrtmp(nxi),ygrtmp(nyi))
1587  xgrtmp=real(xgrdi(1,:))
1588  ygrtmp=real(ygrdi(:,1))
1589 #ifdef W3_OMPH
1590  !$OMP PARALLEL DO PRIVATE(J,I,LONC,LATC,VALUEX,VALUEY)
1591 #endif
1592  DO j=1,ny
1593  DO i=1,nx
1594  lonc=real(xgrdc(j,i)) !LON FOR EVERY CURVL GRID POINT
1595  latc=real(ygrdc(j,i)) !LAT FOR EVERY CURVL GRID POINT
1596 
1597  CALL interpolate2d(nxi,real(xgrtmp),nyi,real(ygrtmp), &
1598  vxi,vyi,lonc,latc,dtoler,valuex,valuey)
1599  vx(i,j)=valuex
1600  vy(i,j)=valuey
1601 
1602  END DO !END I
1603  END DO !END J
1604 #ifdef W3_OMPH
1605  !$OMP END PARALLEL DO
1606 #endif
1607  DEALLOCATE (xgrtmp, ygrtmp)
1608 
1609  ELSE
1610  !
1611  ! 2.b Rectilinear grids
1612  !
1613  ! 2.b.1 Interpolation / averaging data for X axis
1614  !
1615  IF ( sx/sxi .LT. 1.0001 ) THEN
1616  mxa = 2
1617  ELSE
1618  mxa = 2 + int(sx/sxi)
1619  END IF
1620  !
1621 #ifdef W3_T
1622  WRITE (mdst,9020) 'X'
1623 #endif
1624 #ifdef W3_T1
1625  WRITE (format1,'(A,I2,A,I2,A)') "'(10X,",mxa+1,'I5,',mxa+1,"F6.2)'"
1626  WRITE (mdst,9021) nx, mxa
1627 #endif
1628  !
1629  ALLOCATE ( nxa(nx,0:mxa) , rxa(nx,mxa) )
1630  nxa = 0
1631  rxa = 0.
1632  !
1633  IF ( mxa .EQ. 2 ) THEN
1634  !
1635  DO ix=1, nx
1636  IF ( flagll ) THEN
1637  xr = 1. + mod &
1638  ( 1080.+x0+real(ix-1)*sx-x0i , 360. ) / sxi
1639  ELSE
1640  xr = 1. + ( x0+real(ix-1)*sx - x0i ) / sxi
1641  END IF
1642  IF ( xr.GT.0. ) THEN
1643  j1 = int(xr)
1644  j2 = j1 + 1
1645  r2 = max( 0. , xr-real(j1) )
1646  r1 = 1. - r2
1647  IF ( flagll .AND. iclose.NE.iclose_none ) THEN
1648  j1 = 1 + mod(j1-1,nxi)
1649  j2 = 1 + mod(j2-1,nxi)
1650  END IF
1651  IF ( j1.GE.1 .AND. j1.LE.nxi .AND. r1.GT.0.05 ) THEN
1652  nxa(ix,0) = nxa(ix,0) + 1
1653  nxa(ix,nxa(ix,0)) = j1
1654  rxa(ix,nxa(ix,0)) = r1
1655  END IF
1656  IF ( j2.GE.1 .AND. j2.LE.nxi .AND. r2.GT.0.05 ) THEN
1657  nxa(ix,0) = nxa(ix,0) + 1
1658  nxa(ix,nxa(ix,0)) = j2
1659  rxa(ix,nxa(ix,0)) = r2
1660  END IF
1661  IF ( nxa(ix,0) .GT. 0 ) THEN
1662  rt = sum( rxa(ix,:) )
1663  IF ( rt .LT. 0.7 ) THEN
1664  nxa(ix,:) = 0
1665  rxa(ix,:) = 0.
1666  END IF
1667  END IF
1668  END IF
1669  END DO
1670  !
1671  ELSE
1672  !
1673  DO ix=1, nx
1674  !
1675  xfl = x0 + real(ix-1)*sx - 0.5*sx
1676  xfr = x0 + real(ix-1)*sx + 0.5*sx
1677  IF ( flagll ) THEN
1678  ixc = 1 + nint( mod( &
1679  1080.+x0+real(ix-1)*sx-x0i , 360. ) / sxi )
1680  ixs0 = ixc - 1 - mxa/2
1681  ixsn = ixc + 1 + mxa/2
1682  ELSE
1683  ixc = nint( 1. + ( x0+real(ix-1)*sx - x0i ) / sxi )
1684  ixs0 = max( 1 , ixc - 1 - mxa/2 )
1685  ixsn = min( nxi , ixc + 1 + mxa/2 )
1686  END IF
1687  DO j=ixs0, ixsn
1688  jj=j
1689  IF ( flagll ) THEN
1690  IF ( iclose.NE.iclose_none ) jj = 1 + mod(j-1+nxi,nxi)
1691  IF ( jj.LT.1 .OR. jj.GT. nxi ) cycle
1692  ixc = nint((0.5*(xfl+xfr)-x0i-real(jj-1)*sxi)/360.)
1693  IF ( ixc .NE. 0 ) THEN
1694  xfl = xfl - real(ixc) * 360.
1695  xfr = xfr - real(ixc) * 360.
1696  END IF
1697  ELSE
1698  jj = j
1699  END IF
1700  xsl = max( xfl , x0i + real(jj-1)*sxi - 0.5*sxi )
1701  xsr = min( xfr , x0i + real(jj-1)*sxi + 0.5*sxi )
1702  r1 = max( 0. , xsr - xsl ) / sx
1703  IF ( r1 .GT. 0 ) THEN
1704  nxa(ix,0) = nxa(ix,0) + 1
1705  nxa(ix,nxa(ix,0)) = jj
1706  rxa(ix,nxa(ix,0)) = r1
1707  END IF
1708  END DO
1709  IF ( nxa(ix,0) .GT. 0 ) THEN
1710  rt = sum( rxa(ix,:) )
1711  IF ( rt .LT. 0.7 ) THEN
1712  nxa(ix,:) = 0
1713  rxa(ix,:) = 0.
1714  END IF
1715  END IF
1716  END DO
1717  !
1718  END IF
1719  !
1720 #ifdef W3_T1
1721  DO, ix=1, nx
1722  IF ( nxa(ix,0) .GT. 0 ) WRITE (mdst,format1) &
1723  ix, nxa(ix,1:mxa), rxa(ix,1:mxa), sum(rxa(ix,1:mxa))
1724  END DO
1725 #endif
1726  !
1727  ! 2.b.2 Interpolation / averaging data for Y axis
1728  !
1729  IF ( sy/syi .LT. 1.0001 ) THEN
1730  mya = 2
1731  ELSE
1732  mya = 2 + int(sy/syi)
1733  END IF
1734  !
1735 #ifdef W3_T
1736  WRITE (mdst,9020) 'Y'
1737 #endif
1738 #ifdef W3_T1
1739  format1 = '(10X, I5, F6.2)'
1740  WRITE (format1,'(A,I2,A,I2,A)') "'(10X,",mya+1,'I5,',mya+1,"F6.2)'"
1741  WRITE (mdst,9021) ny, mya
1742 #endif
1743  !
1744  ALLOCATE ( nya(ny,0:mya) , rya(ny,mya) )
1745  nya = 0
1746  rya = 0.
1747  !
1748  IF ( mya .EQ. 2 ) THEN
1749  !
1750  DO iy=1, ny
1751  yr = 1. + ( y0+real(iy-1)*sy - y0i ) / syi
1752  IF ( yr.GT.0. ) THEN
1753  j1 = int(yr)
1754  j2 = j1 + 1
1755  r2 = max( 0. , yr-real(j1) )
1756  r1 = 1. - r2
1757  IF ( j1.GE.1 .AND. j1.LE.nyi .AND. r1.GT.0.05 ) THEN
1758  nya(iy,0) = nya(iy,0) + 1
1759  nya(iy,nya(iy,0)) = j1
1760  rya(iy,nya(iy,0)) = r1
1761  END IF
1762  IF ( j2.GE.1 .AND. j2.LE.nyi .AND. r2.GT.0.05 ) THEN
1763  nya(iy,0) = nya(iy,0) + 1
1764  nya(iy,nya(iy,0)) = j2
1765  rya(iy,nya(iy,0)) = r2
1766  END IF
1767  IF ( nya(iy,0) .GT. 0 ) THEN
1768  rt = sum( rya(iy,:) )
1769  IF ( rt .LT. 0.7 ) THEN
1770  nya(iy,:) = 0
1771  rya(iy,:) = 0.
1772  END IF
1773  END IF
1774  END IF
1775  END DO
1776  !
1777  ELSE
1778  !
1779  DO iy=1, ny
1780  yfl = y0 + real(iy-1)*sy - 0.5*sy
1781  yfr = y0 + real(iy-1)*sy + 0.5*sy
1782  iyc = nint( 1. + ( y0+real(iy-1)*sy - y0i ) / syi )
1783  iys0 = max( 1 , iyc - 1 - mya/2 )
1784  iysn = min( nyi , iyc + 1 + mya/2 )
1785  DO j=iys0, iysn
1786  ysl = max( yfl , y0i + real(j-1)*syi - 0.5*syi )
1787  ysr = min( yfr , y0i + real(j-1)*syi + 0.5*syi )
1788  r1 = max( 0. , ysr - ysl ) / sy
1789  IF ( r1 .GT. 0 ) THEN
1790  nya(iy,0) = nya(iy,0) + 1
1791  nya(iy,nya(iy,0)) = j
1792  rya(iy,nya(iy,0)) = r1
1793  END IF
1794  END DO
1795  IF ( nya(iy,0) .GT. 0 ) THEN
1796  rt = sum( rya(iy,:) )
1797  IF ( rt .LT. 0.7 ) THEN
1798  nya(iy,:) = 0
1799  rya(iy,:) = 0.
1800  END IF
1801  END IF
1802  END DO
1803  !
1804  END IF
1805  !
1806  END IF
1807 
1808  !
1809 #ifdef W3_T1
1810  DO, iy=1, ny
1811  IF ( nya(iy,0) .GT. 0 ) WRITE (mdst,format1) &
1812  iy, nya(iy,1:mya), rya(iy,1:mya), sum(rya(iy,1:mya))
1813  END DO
1814 #endif
1815  !
1816  ! 2.c Process grid
1817  !
1818  map1 = .false.
1819  map2 = .false.
1820  factor = 1.
1821  !
1822 
1823  DO ix=1, nx
1824  IF ( nxa(ix,0) .EQ. 0 ) cycle
1825  DO iy=1, ny
1826  IF ( nya(iy,0) .EQ. 0 ) cycle
1827  IF ( map(iy,ix).NE.0 ) THEN
1828  vxl = 0.
1829  vyl = 0.
1830  va = 0.
1831  va2 = 0.
1832  wtot = 0.
1833  DO j1=1, nxa(ix,0)
1834  jx = nxa(ix,j1)
1835  DO j2=1, nya(iy,0)
1836  jy = nya(iy,j2)
1837  IF ( mapi(jy,jx) .NE. 0 ) THEN
1838  wl = rxa(ix,j1) * rya(iy,j2)
1839  wtot = wtot + wl
1840  vxl = vxl + wl * vxi(jx,jy)
1841  vyl = vyl + wl * vyi(jx,jy)
1842  va = va + wl * sqrt &
1843  ( vxi(jx,jy)**2 + vyi(jx,jy)**2 )
1844  va2 = va2 + wl * &
1845  ( vxi(jx,jy)**2 + vyi(jx,jy)**2 )
1846  END IF
1847  END DO
1848  END DO
1849  IF ( wtot .LT. 0.05 ) THEN
1850  map1(ix,iy) = .true.
1851  ELSE
1852  map2(ix,iy) = .true.
1853  vxl = vxl / wtot
1854  vyl = vyl / wtot
1855  va = va / wtot
1856  va2 = sqrt( va2 / wtot )
1857  va0 = sqrt( vxl**2 + vyl**2 )
1858  IF ( constp .EQ. 1 ) THEN
1859  factor = min( 1.25 , va/max(1.e-7,va0) )
1860  ELSE IF ( constp .EQ. 2 ) THEN
1861  factor = min( 1.25 , va2/max(1.e-7,va0) )
1862  END IF
1863  vx(ix,iy) = factor * vxl
1864  vy(ix,iy) = factor * vyl
1865  END IF
1866  END IF
1867  END DO
1868  END DO
1869 
1870  !
1871  ! 2.d Reconcile mask differences
1872  !
1873 #ifdef W3_T
1874  WRITE (mdst,9022)
1875 #endif
1876  !
1877  jj = 0
1878  iclose => grids(imod)%ICLOSE
1879  !
1880  DO
1881  IF ( jj .GT. swpmax ) EXIT
1882  flagup = .false.
1883  map3 = .false.
1884  jj = jj + 1
1885 #ifdef W3_T
1886  WRITE (mdst,9023) jj
1887 #endif
1888  DO ix=1, nx
1889  DO iy=1, ny
1890  IF ( map1(ix,iy) ) THEN
1891  vxl = 0.
1892  vyl = 0.
1893  j1 = 0
1894  IF ( flagll ) THEN
1895  DO j2=ix-1, ix+1
1896  IF ( (j2.GT.1 .AND. j2.LE.nx) .OR. iclose.NE.iclose_none ) THEN
1897  jx = 1 + mod(nx+j2-1,nx)
1898  DO jy=iy-1, iy+1
1899  IF ( jy.GT.1 .AND. jy.LE.ny ) THEN
1900  IF ( map2(jx,jy) ) THEN
1901  vxl = vxl + vx(jx,jy)
1902  vyl = vyl + vy(jx,jy)
1903  j1 = j1 + 1
1904  END IF
1905  END IF
1906  END DO
1907  END IF
1908  END DO
1909  ELSE
1910  DO jx=ix-1, ix+1
1911  IF ( jx.GT.1 .AND. jx.LE.nx ) THEN
1912  DO jy=iy-1, iy+1
1913  IF ( jy.GT.1 .AND. jy.LE.ny ) THEN
1914  IF ( map2(jx,jy) ) THEN
1915  vxl = vxl + vx(jx,jy)
1916  vyl = vyl + vy(jx,jy)
1917  j1 = j1 + 1
1918  END IF
1919  END IF
1920  END DO
1921  END IF
1922  END DO
1923  END IF !FLAGLL
1924  IF ( j1 .GT. 0 ) THEN
1925  vx(ix,iy) = vxl / real(j1)
1926  vy(ix,iy) = vyl / real(j1)
1927  map1(ix,iy) = .false.
1928  map3(ix,iy) = .true.
1929  flagup = .true.
1930  END IF
1931  END IF
1932  END DO
1933  END DO
1934  IF ( flagup ) THEN
1935  map2 = map2 .OR. map3
1936  ELSE
1937  EXIT
1938  END IF
1939  END DO
1940 
1941  !
1942  ! 3. End of routine -------------------------------------------------- /
1943  !
1944  DEALLOCATE ( nxa, nya, rxa, rya )
1945  !
1946  RETURN
1947  !
1948  ! Formats
1949  !
1950 #ifdef W3_T
1951 9000 FORMAT ( ' TEST WMUPDV : GRID INFORMATION : '/ &
1952  ' ',3i5,4e11.3/ &
1953  ' ',3i5,4e11.3/ &
1954  ' UNDEFINED = ',e10.3)
1955 9010 FORMAT ( ' TEST WMUPDV : COINCIDING GRIDS, OFFSETS :',2i6/ &
1956  ' TARGET GRID RANGES :',4i6/ &
1957  ' SOURCE GRID RANGES :',4i6)
1958 9020 FORMAT ( ' TEST WMUPDV : WEIGHTS FOR ',a,' INTERPOATION')
1959 #endif
1960 #ifdef W3_T1
1961 9021 FORMAT ( ' TEST WMUPDV : ARAY DIMENSIONED AS : ',2i6)
1962 #endif
1963 #ifdef W3_T
1964 9022 FORMAT ( ' TEST WMUPDV : RECONCILING MASKS')
1965 9023 FORMAT ( ' SWEEP NR ',i4)
1966 #endif
1967  !/
1968  !/ End of WMUPDV ----------------------------------------------------- /
1969  !/

References w3gdatmd::clgtype, w3servmd::extcde(), w3gdatmd::flagll, w3gdatmd::grids, w3gdatmd::gtype, w3gdatmd::hpfac, w3gdatmd::hqfac, w3gdatmd::iclose_none, w3gdatmd::iclose_smpl, w3gdatmd::iclose_trpl, wmmdatmd::improc, interpolate2d(), wmmdatmd::mdse, wmmdatmd::mdso, wmmdatmd::mdss, wmmdatmd::mdst, wmmdatmd::nmperr, wmmdatmd::nmpscr, w3gdatmd::nx, w3gdatmd::ny, w3gdatmd::rlgtype, w3servmd::strace(), swpmax, w3gdatmd::sx, w3gdatmd::sy, w3gdatmd::ungtype, w3gdatmd::x0, w3gdatmd::xgrd, w3gdatmd::y0, and w3gdatmd::ygrd.

Referenced by wmupd2().

◆ xycurvisearch()

integer function wmupdtmd::xycurvisearch ( integer, intent(in)  LENGTH,
real, dimension(length), intent(in)  ARRAY,
real, intent(in)  VALUE,
real, intent(in), optional  DELTA 
)

Search the location of a point(XC,YC) in a regular grid.

Given an array and a value to search, it returns the index of the element on the rectilinear grid that is closest to, but less than, the given value. "delta" is the threshold used to determine if two values are equal.

       if ( abs(x1 - x2) <= delta) then
          assume x1 = x2
       endif
Parameters
LENGTHDimension of input array
ARRAY1D array for lats or longs
VALUEValue to located in ARRAY
DELTAThreshold to determine if two values are equal
Returns
XYCURVISEARCH
Author
H. L. Tolman
Date
20-Jan-2017

Definition at line 2625 of file wmupdtmd.F90.

2625  !/ +-----------------------------------+
2626  !/ | WAVEWATCH III NOAA/NCEP |
2627  !/ | H. L. Tolman |
2628  !/ | FORTRAN 90 |
2629  !/ | Last update : 20-Jan-2017 |
2630  !/ +-----------------------------------+
2631  !/ (R. Padilla-Hernandez, EMC/NOAA)
2632  !/
2633  !/ 01-Jul-2019 : Origination. ( version 7.13 )
2634  !/ 01-Jul-2019 : Generalize output to curv grids ( version 7.13 )
2635  !/
2636  !/ Copyright 2009 National Weather Service (NWS),
2637  !/ National Oceanic and Atmospheric Administration. All rights
2638  !/ reserved. WAVEWATCH III is a trademark of the NWS.
2639  !/ No unauthorized use without permission.
2640  !/
2641  ! 1. Purpose :
2642  !
2643  ! Search the location of a point(XC,YC) in a regular grid
2644  !
2645  !
2646  ! 2. Parameters :
2647  !
2648  ! Parameter list
2649  ! ----------------------------------------------------------------
2650  ! LENGTH Int. Input Dimension of input array
2651  ! ARRAY Int. Input 1D array for lats or longs
2652  ! VALUE Real Input Value to be located in ARRAY
2653  ! DELTA Real Input Threshold to determine if two values
2654  ! are equal
2655  ! ----------------------------------------------------------------
2656  !
2657  ! Internal parameters
2658  ! ----------------------------------------------------------------
2659  !
2660  ! ----------------------------------------------------------------
2661  !
2662  ! 3. Subroutines and functions :
2663  !
2664  ! Name Type Scope Description
2665  ! ----------------------------------------------------------------
2666  ! XYCURVISEARCH Function Find indexes See bellow
2667  ! ----------------------------------------------------------------
2668  !
2669  ! 4. Subroutines and functions used :
2670  !
2671  ! -
2672  !
2673  ! 5. Remarks :
2674  ! GIVEN AN ARRAY AND A VALUE TO SEARCH, IT RETURNS THE INDEX OF
2675  ! THE ELEMENT ON THE RECTILINEAR GRID THAT IS CLOSEST TO, BUT
2676  ! LESS THAN, THE GIVEN VALUE.
2677  ! "DELTA" IS THE THERSHOLD USED TO DETERMINE IF TWO VALUES ARE EQUAL
2678  ! IF ( ABS(X1 - X2) <= DELTA) THEN
2679  ! ASSUME X1 = X2
2680  ! ENDIF
2681  !
2682  ! 6. Switches :
2683  ! -
2684  !
2685  ! 7. Source code :
2686 
2687  INTEGER, INTENT(IN) :: LENGTH
2688  REAL, DIMENSION(LENGTH), INTENT(IN) :: ARRAY
2689  REAL, INTENT(IN) :: VALUE
2690  REAL, INTENT(IN), OPTIONAL :: DELTA
2691 
2692  INTEGER :: XYCURVISEARCH
2693 
2694  INTEGER :: LEFT, MIDDLE, RIGHT
2695 
2696 
2697  left = 1
2698  right = length
2699  DO
2700  IF (left > right) THEN
2701  EXIT
2702  ENDIF
2703  middle = nint((left+right) / 2.0)
2704  IF ( abs(array(middle) - VALUE) <= delta) THEN
2705  xycurvisearch = middle
2706  RETURN
2707  ELSE IF (array(middle) > VALUE) THEN
2708  right = middle - 1
2709  ELSE
2710  left = middle + 1
2711  END IF
2712  END DO
2713  xycurvisearch = right
2714 

Referenced by averaging(), interpolate(), and interpolate2d().

Variable Documentation

◆ swpmax

integer, parameter wmupdtmd::swpmax = 5

SWPMAX.

Definition at line 88 of file wmupdtmd.F90.

88  INTEGER, PARAMETER :: SWPMAX = 5

Referenced by wmupds(), and wmupdv().

w3idatmd::gdn
real, pointer gdn
Definition: w3idatmd.F90:242
w3timemd::dsec21
real function dsec21(TIME1, TIME2)
Definition: w3timemd.F90:333
w3fldsmd::w3fldd
subroutine w3fldd(INXOUT, IDFLD, NDS, NDST, NDSE, TIME, TD, NR, ND, NDOUT, DATA, IERR)
Definition: w3fldsmd.F90:1474
wmmdatmd::mdse
integer mdse
MDSE.
Definition: wmmdatmd.F90:316
w3idatmd::twn
integer, dimension(:), pointer twn
Definition: w3idatmd.F90:236
w3idatmd::inflags1
logical, dimension(:), pointer inflags1
Definition: w3idatmd.F90:260
w3gdatmd::ygrd
double precision, dimension(:,:), pointer ygrd
Definition: w3gdatmd.F90:1205
w3idatmd::t1n
integer, dimension(:), pointer t1n
Definition: w3idatmd.F90:236
wmmdatmd::iflstl
logical, dimension(:), allocatable iflstl
IFLSTL.
Definition: wmmdatmd.F90:384
w3idatmd::dt0
real, dimension(:,:), pointer dt0
Definition: w3idatmd.F90:243
wmmdatmd::iflstr
logical, dimension(:), allocatable iflstr
IFLSTR.
Definition: wmmdatmd.F90:385
wmmdatmd::fllsti
logical, pointer fllsti
FLLSTI.
Definition: wmmdatmd.F90:562
wmmdatmd::dmv
real, dimension(:,:), pointer dmv
DMV.
Definition: wmmdatmd.F90:551
constants::dair
real, parameter dair
DAIR Density of air (kg/m3).
Definition: constants.F90:63
w3gdatmd::fswnd
logical, pointer fswnd
Definition: w3gdatmd.F90:1264
wmmdatmd::rcld
integer, dimension(:), pointer rcld
RCLD.
Definition: wmmdatmd.F90:535
w3gdatmd::ungtype
integer, parameter ungtype
Definition: w3gdatmd.F90:626
w3wdatmd
Define data structures to set up wave model dynamic data for several models simultaneously.
Definition: w3wdatmd.F90:18
wmmdatmd::nmpscr
integer nmpscr
NMPSCR.
Definition: wmmdatmd.F90:324
wmmdatmd::data2
real, dimension(:,:), pointer data2
DATA2.
Definition: wmmdatmd.F90:549
wmmdatmd::mdso
integer mdso
MDSO.
Definition: wmmdatmd.F90:313
w3gdatmd::rlgtype
integer, parameter rlgtype
Definition: w3gdatmd.F90:624
w3idatmd::wy0
real, dimension(:,:), pointer wy0
Definition: w3idatmd.F90:243
w3idatmd::inputs
type(input), dimension(:), allocatable, target inputs
Definition: w3idatmd.F90:232
wmmdatmd::fllstl
logical, pointer fllstl
FLLSTL.
Definition: wmmdatmd.F90:560
w3gdatmd::xgrd
double precision, dimension(:,:), pointer xgrd
Definition: w3gdatmd.F90:1205
w3gdatmd::sy
real, pointer sy
Definition: w3gdatmd.F90:1183
w3idatmd::t0n
integer, dimension(:), pointer t0n
Definition: w3idatmd.F90:236
wmmdatmd::nmv
integer, pointer nmv
NMV.
Definition: wmmdatmd.F90:537
w3idatmd::ti5
integer, dimension(:), pointer ti5
Definition: w3idatmd.F90:236
w3idatmd::gan
real, pointer gan
Definition: w3idatmd.F90:242
wmmdatmd::mdss
integer mdss
MDSS.
Definition: wmmdatmd.F90:314
w3wdatmd::time
integer, dimension(:), pointer time
Definition: w3wdatmd.F90:172
w3gdatmd::grids
type(grid), dimension(:), allocatable, target grids
Definition: w3gdatmd.F90:1088
w3gdatmd::ny
integer, pointer ny
Definition: w3gdatmd.F90:1097
w3idatmd::tg0
integer, dimension(:), pointer tg0
Definition: w3idatmd.F90:236
w3fldsmd::w3fldm
subroutine w3fldm(J, NDST, NDSE, T0, TN, NH, NHM, THO, HA, HD, TF0, A0, D0, TFN, AN, DN, IERR)
Definition: w3fldsmd.F90:2491
wmmdatmd::improc
integer improc
IMPROC.
Definition: wmmdatmd.F90:322
w3idatmd::ti3
integer, dimension(:), pointer ti3
Definition: w3idatmd.F90:236
w3idatmd::tu0
integer, dimension(:), pointer tu0
Definition: w3idatmd.F90:236
w3gdatmd::iclose_none
integer, parameter iclose_none
Definition: w3gdatmd.F90:629
w3gdatmd::hqfac
real, dimension(:,:), pointer hqfac
Definition: w3gdatmd.F90:1212
w3gdatmd::w3setg
subroutine w3setg(IMOD, NDSE, NDST)
Definition: w3gdatmd.F90:2152
wmmdatmd::data1
real, dimension(:,:), pointer data1
DATA1.
Definition: wmmdatmd.F90:548
w3idatmd::tzn
integer, dimension(:), pointer tzn
Definition: w3idatmd.F90:236
w3idatmd::jfirst
integer jfirst
Definition: w3idatmd.F90:162
w3idatmd::wxn
real, dimension(:,:), pointer wxn
Definition: w3idatmd.F90:243
w3gdatmd::x0
real, pointer x0
Definition: w3gdatmd.F90:1183
w3gdatmd::nsea
integer, pointer nsea
Definition: w3gdatmd.F90:1097
wmmdatmd::amv
real, dimension(:,:), pointer amv
AMV.
Definition: wmmdatmd.F90:550
wmmdatmd::fllstr
logical, pointer fllstr
FLLSTR.
Definition: wmmdatmd.F90:561
w3gdatmd::clgtype
integer, parameter clgtype
Definition: w3gdatmd.F90:625
w3servmd
Definition: w3servmd.F90:3
wmmdatmd::idinp
character(len=3), dimension(:,:), allocatable idinp
IDINP.
Definition: wmmdatmd.F90:386
w3timemd::tick21
subroutine tick21(TIME, DTIME)
Definition: w3timemd.F90:84
w3wdatmd::w3setw
subroutine w3setw(IMOD, NDSE, NDST)
Select one of the WAVEWATCH III grids / models.
Definition: w3wdatmd.F90:660
wmmdatmd::tmv
integer, dimension(:,:,:), pointer tmv
TMV.
Definition: wmmdatmd.F90:538
w3timemd::stme21
subroutine stme21(TIME, DTME21)
Definition: w3timemd.F90:682
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
wmmdatmd::nmperr
integer nmperr
NMPERR.
Definition: wmmdatmd.F90:326
w3idatmd::w3seti
subroutine w3seti(IMOD, NDSE, NDST)
Select one of the WAVEWATCH III grids / models.
Definition: w3idatmd.F90:819
wmmdatmd::ndt
integer, dimension(:), pointer ndt
NDT.
Definition: wmmdatmd.F90:536
wmmdatmd::mdst
integer mdst
MDST.
Definition: wmmdatmd.F90:315
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::ga0
real, pointer ga0
Definition: w3idatmd.F90:242
w3servmd::strace
subroutine strace(IENT, SNAME)
Definition: w3servmd.F90:148
wmmdatmd::wmsetm
subroutine wmsetm(IMOD, NDSE, NDST)
Select one of the WAVEWATCH III grids / models.
Definition: wmmdatmd.F90:1169
w3gdatmd::gtype
integer, pointer gtype
Definition: w3gdatmd.F90:1094
w3idatmd
Define data structures to set up wave model input data for several models simultaneously.
Definition: w3idatmd.F90:16
w3idatmd::t2n
integer, dimension(:), pointer t2n
Definition: w3idatmd.F90:236
wmmdatmd::wmdimd
subroutine wmdimd(IMOD, NDSE, NDST, J)
Initialize an individual data grid at the proper dimensions.
Definition: wmmdatmd.F90:787
wmmdatmd::data0
real, dimension(:,:), pointer data0
DATA0.
Definition: wmmdatmd.F90:547
w3idatmd::wx0
real, dimension(:,:), pointer wx0
Definition: w3idatmd.F90:243
w3gdatmd::y0
real, pointer y0
Definition: w3gdatmd.F90:1183
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
w3gdatmd::hpfac
real, dimension(:,:), pointer hpfac
Definition: w3gdatmd.F90:1211
w3gdatmd::sx
real, pointer sx
Definition: w3gdatmd.F90:1183
constants
Define some much-used constants for global use (all defined as PARAMETER).
Definition: constants.F90:20
w3gdatmd
Definition: w3gdatmd.F90:16
wmmdatmd
Define data structures to set up wave model dynamic data for several models simultaneously.
Definition: wmmdatmd.F90:16
w3gdatmd::iclose_trpl
integer, parameter iclose_trpl
Definition: w3gdatmd.F90:631
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
w3idatmd::ti2
integer, dimension(:), pointer ti2
Definition: w3idatmd.F90:236
w3idatmd::tln
integer, dimension(:), pointer tln
Definition: w3idatmd.F90:236
wmmdatmd::nmvmax
integer nmvmax
NMVMAX.
Definition: wmmdatmd.F90:333
wmmdatmd::inpmap
integer, dimension(:,:), allocatable inpmap
INPMAP.
Definition: wmmdatmd.F90:368
w3gdatmd::nx
integer, pointer nx
Definition: w3gdatmd.F90:1097
w3idatmd::trn
integer, dimension(:), pointer trn
Definition: w3idatmd.F90:236
w3timemd
Definition: w3timemd.F90:3
w3gdatmd::iclose_smpl
integer, parameter iclose_smpl
Definition: w3gdatmd.F90:630
w3idatmd::tun
integer, dimension(:), pointer tun
Definition: w3idatmd.F90:236
w3idatmd::tr0
integer, dimension(:), pointer tr0
Definition: w3idatmd.F90:236
wmmdatmd::mdsf
integer, dimension(:,:), allocatable mdsf
MDSF.
Definition: wmmdatmd.F90:352
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
wmmdatmd::etime
integer, dimension(2) etime
ETIME.
Definition: wmmdatmd.F90:329
w3fldsmd
Definition: w3fldsmd.F90:3
w3gdatmd::flagll
logical, pointer flagll
Definition: w3gdatmd.F90:1219
wmmdatmd::iflsti
logical, dimension(:), allocatable iflsti
IFLSTI.
Definition: wmmdatmd.F90:383
w3gdatmd::filext
character(len=13), pointer filext
Definition: w3gdatmd.F90:1224