1 |
C $Header: /u/gcmpack/MITgcm/pkg/timeave/TIMEAVE_STATV.h,v 1.16 2004/12/04 01:54:04 dimitri Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "TIMEAVE_OPTIONS.h" |
5 |
|
6 |
#ifdef ALLOW_TIMEAVE |
7 |
|
8 |
CBOP |
9 |
C !ROUTINE: TIMEAVE_STATV.h |
10 |
C !INTERFACE: |
11 |
C include "TIMEAVE_STATV.h" |
12 |
C !DESCRIPTION: \bw |
13 |
C *================================================================* |
14 |
C | TIMEAVE_STATV.h |
15 |
C | o Time averages of model state-variables |
16 |
C | (common block TAVE_STATEVARS) |
17 |
C *================================================================* |
18 |
C | Time average of state variables is (generally) centered on the |
19 |
C | middle of the time step (time average interval = TimeAve_half) |
20 |
C | Time average of intermediate and tendency variables is centered |
21 |
C | on the time step (time average interval=TimeAve_full) |
22 |
C *================================================================* |
23 |
C \ev |
24 |
CEOP |
25 |
|
26 |
C TimeAve_* :: time of temporal integration (s) *** for each thread *** |
27 |
C TimeAve_half :: half time_step multiple (used for state variables) |
28 |
C TimeAve_full :: full time_step multiple (used for for intermediate var.) |
29 |
C uFluxtave :: zonal surface wind stress (N/m^2, |
30 |
C >0 for increase in uVel, i=1 held at western face) |
31 |
C vFluxtave :: meridional surface wind stress (N/m^2, |
32 |
C >0 for increase in vVel, j=1 held at southern face) |
33 |
C tFluxtave :: net surface heat flux (W/m^2, >0 for increase in theta) |
34 |
C sFluxtave :: net surface salt flux (g/m^2/s, >0 for increase in salt) |
35 |
C etatave :: surface displacement (r unit, i.e. ocean:z, atmos:p) |
36 |
C uVeltave :: zonal velocity (m/s, i=1 held at western face) |
37 |
C vVeltave :: meridional velocity (m/s, j=1 held at southern face) |
38 |
C wVeltave :: vertical velocity ([r]/s, i.e.: ocean:m/s atmos:Pa/s) |
39 |
C thetatave :: potential temperature (oC, held at pressure/tracer point) |
40 |
C salttave :: salinity (ppt, held at pressure/tracer point) |
41 |
C Eta2tave :: eta * eta |
42 |
C TTtave :: theta * theta |
43 |
C UUtave :: uVel * uVel (used to compute the averaged KE) |
44 |
C VVtave :: vVel * vVel (used to compute the averaged KE) |
45 |
C UVtave :: uVel * vVel (at vorticity point, i.e. grid-corner) |
46 |
C KEtave :: Kinetic Energy |
47 |
C UTtave :: uVel * theta (* hFacW) |
48 |
C VTtave :: vVel * theta (* hFacS) |
49 |
C WTtave :: wVel * theta |
50 |
C UStave :: uVel * salt (* hFacW) |
51 |
C VStave :: vVel * salt (* hFacS) |
52 |
C WStave :: wVel * salt |
53 |
C tDiffRtave :: vertical diffusion flux of Temperature (theta) |
54 |
C uZetatave :: uVel*Relativ_Vorticity_3 (computed at v point) |
55 |
C vZetatave :: vVel*Relativ_Vorticity_3 (computed at u point) |
56 |
C phiHydtave :: Hydrostatic (ocean) pressure / (atmos) geo- Potential |
57 |
C phiHydLowtave:: Hydrostatic (ocean) pressure / (atmos) geo- Potential |
58 |
C at the fixed boundary: (ocean) bottom pressure |
59 |
C (atmos) geo- Potential |
60 |
C ConvectCountTave :: Average number of convective adjustment event |
61 |
|
62 |
COMMON /TAVE_TIME/ TimeAve_half,TimeAve_full |
63 |
_RL TimeAve_half(Nr,nSx,nSy) |
64 |
_RL TimeAve_full(Nr,nSx,nSy) |
65 |
|
66 |
COMMON /TAVE_STATEVARS/ |
67 |
& uFluxtave,vFluxtave,tFluxtave,sFluxtave |
68 |
& ,etatave,uVeltave,vVeltave,wVeltave |
69 |
& ,thetatave,salttave,phiHydLowtave |
70 |
& ,UTtave,VTtave,WTtave,UStave,VStave,WStave |
71 |
& ,Eta2tave,TTtave,UUtave,VVtave,UVtave |
72 |
& ,TdiffRtave |
73 |
#ifdef ALLOW_MOM_VECINV |
74 |
& ,uZetatave, vZetatave |
75 |
#endif /* ALLOW_MOM_VECINV */ |
76 |
& ,phiHydtave |
77 |
& ,phiHydLow2Tave |
78 |
& ,ConvectCountTave |
79 |
_RL uFluxtave(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
80 |
_RL vFluxtave(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
81 |
_RL tFluxtave(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
82 |
_RL sFluxtave(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
83 |
_RL etatave (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
84 |
_RL uVeltave (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
85 |
_RL vVeltave (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
86 |
_RL wVeltave (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
87 |
_RL thetatave(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
88 |
_RL salttave (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
89 |
_RL phiHydLowtave (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
90 |
_RL UTtave (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
91 |
_RL VTtave (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
92 |
_RL WTtave (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
93 |
_RL UStave (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
94 |
_RL VStave (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
95 |
_RL WStave (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
96 |
_RL eta2Tave (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
97 |
_RL TTtave (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
98 |
_RL UUtave (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
99 |
_RL VVtave (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
100 |
_RL UVtave (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
101 |
_RL TdiffRtave(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
102 |
#ifdef ALLOW_MOM_VECINV |
103 |
_RL uZetatave (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
104 |
_RL vZetatave (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
105 |
#endif /* ALLOW_MOM_VECINV */ |
106 |
_RL phiHydtave(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
107 |
_RL phiHydLow2Tave(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
108 |
_RL ConvectCountTave(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
109 |
|
110 |
#ifdef NONLIN_FRSURF |
111 |
C hUtave :: average zonal flow (=hFacW*uVel) (still in m/s !) |
112 |
C hVtave :: average merid.flow (=hFacS*vVel) (still in m/s !) |
113 |
C hFacCtave :: average thickness fraction of open water, Center |
114 |
C hFacWtave :: average thickness fraction of open water, West side |
115 |
C hFacStave :: average thickness fraction of open water, South side |
116 |
|
117 |
COMMON /TAVE_THICKNESS/ |
118 |
& hUtave, hVtave |
119 |
c & , hFacCtave, hFacWtave, hFacStave |
120 |
_RL hUtave (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
121 |
_RL hVtave (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
122 |
c _RL hFacCtave(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
123 |
c _RL hFacWtave(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
124 |
c _RL hFacStave(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
125 |
#endif /* NONLIN_FRSURF */ |
126 |
|
127 |
#endif /* ALLOW_TIMEAVE */ |