1 |
C $Header: /u/gcmpack/MITgcm/model/inc/DYNVARS.h,v 1.40 2010/02/17 00:21:32 gforget Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
CBOP |
5 |
C !ROUTINE: DYNVARS.h |
6 |
C !INTERFACE: |
7 |
C include "DYNVARS.h" |
8 |
C !DESCRIPTION: |
9 |
C \bv |
10 |
C *==========================================================* |
11 |
C | DYNVARS.h |
12 |
C | o Dynamical model variables (common block DYNVARS_R) |
13 |
C *==========================================================* |
14 |
C | The value and two levels of time tendency are held for |
15 |
C | each prognostic variable. |
16 |
C *==========================================================* |
17 |
C \ev |
18 |
CEOP |
19 |
|
20 |
C State Variables: |
21 |
C etaN :: free-surface r-anomaly (r unit) at current time level |
22 |
C uVel :: zonal velocity (m/s, i=1 held at western face) |
23 |
C vVel :: meridional velocity (m/s, j=1 held at southern face) |
24 |
C theta :: potential temperature (oC, held at pressure/tracer point) |
25 |
C salt :: salinity (ppt, held at pressure/tracer point) |
26 |
C gX, gxNm1 :: Time tendencies at current and previous time levels. |
27 |
C etaH :: surface r-anomaly, advanced in time consistently |
28 |
C with 2.D flow divergence (Exact-Conservation): |
29 |
C etaH^n+1 = etaH^n - delta_t*Div.(H^n U^n) |
30 |
C note: a) used with "exactConserv" but strictly necessary for NonLinFreeSurf |
31 |
C b) same as etaN but not necessarely at the same time, e.g.: |
32 |
C implicDiv2DFlow=0 => etaH=etaN ; =1 => etaH=etaNm1 ; |
33 |
|
34 |
#ifdef ALLOW_ADAMSBASHFORTH_3 |
35 |
COMMON /DYNVARS_R/ |
36 |
& etaN, |
37 |
& uVel,vVel,wVel,theta,salt, |
38 |
& gU, gV, gT, gS, |
39 |
& guNm, gvNm, gtNm, gsNm |
40 |
#else /* ALLOW_ADAMSBASHFORTH_3 */ |
41 |
COMMON /DYNVARS_R/ |
42 |
& etaN, |
43 |
& uVel,vVel,wVel,theta,salt, |
44 |
& gU, gV, gT, gS, |
45 |
& guNm1,gvNm1,gtNm1,gsNm1 |
46 |
#endif /* ALLOW_ADAMSBASHFORTH_3 */ |
47 |
_RL etaN (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
48 |
_RL uVel (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
49 |
_RL vVel (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
50 |
_RL wVel (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
51 |
_RL theta(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
52 |
_RL salt (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
53 |
_RL gU(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
54 |
_RL gV(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
55 |
_RL gT(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
56 |
_RL gS(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
57 |
#ifdef ALLOW_ADAMSBASHFORTH_3 |
58 |
_RL guNm(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy,2) |
59 |
_RL gvNm(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy,2) |
60 |
_RL gtNm(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy,2) |
61 |
_RL gsNm(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy,2) |
62 |
#else /* ALLOW_ADAMSBASHFORTH_3 */ |
63 |
_RL guNm1(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
64 |
_RL gvNm1(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
65 |
_RL gtNm1(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
66 |
_RL gsNm1(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
67 |
#endif /* ALLOW_ADAMSBASHFORTH_3 */ |
68 |
|
69 |
COMMON /DYNVARS_R_2/ |
70 |
& etaH |
71 |
_RL etaH (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
72 |
|
73 |
#ifdef ALLOW_3D_DIFFKR |
74 |
C diffKr :: full 3D specification of Laplacian diffusion coeff. |
75 |
C for mixing of tracers vertically ( units of r^2/s ) |
76 |
COMMON /DYNVARS_DIFFKR/ |
77 |
& diffKr |
78 |
_RL diffKr (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
79 |
#endif |
80 |
|
81 |
cph( |
82 |
cph the following block will eventually move to a separate |
83 |
cph header file containing requires anomaly fields of control vars. |
84 |
cph |
85 |
#if (defined (ALLOW_AUTODIFF) && defined (ALLOW_KAPGM_CONTROL)) |
86 |
COMMON /DYNVARS_KAPGM/ |
87 |
& kapgm |
88 |
_RL kapgm (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
89 |
#endif |
90 |
#if (defined (ALLOW_AUTODIFF) && defined (ALLOW_KAPREDI_CONTROL)) |
91 |
COMMON /DYNVARS_KAPREDI/ |
92 |
& kapredi |
93 |
_RL kapredi (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
94 |
#endif |
95 |
#if (defined (ALLOW_AUTODIFF) && defined (ALLOW_BOTTOMDRAG_CONTROL)) |
96 |
COMMON /DYNVARS_BOTTOMDRAG/ |
97 |
& bottomdragfld |
98 |
_RL bottomdragfld (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
99 |
#endif |
100 |
cph |
101 |
cph) |
102 |
|
103 |
#ifdef ALLOW_BL79_LAT_VARY |
104 |
C BL79LatArray :: is used for latitudinal dependence of |
105 |
C BryanLewis79 vertical diffusivity |
106 |
COMMON /DYNVARS_BL79LatArray/ BL79LatArray |
107 |
_RL BL79LatArray (1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
108 |
#endif |
109 |
|
110 |
#ifdef ALLOW_ADDFLUID |
111 |
C addMass :: source (<0: sink) of fluid in the domain interior |
112 |
C in mass per unit of time [kg/s] |
113 |
C (generalisation of oceanic real fresh-water flux) |
114 |
COMMON /DYNVARS_ADD_FLUID/ addMass |
115 |
_RL addMass(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
116 |
#endif |
117 |
|
118 |
C Diagnostic Variables: |
119 |
C phiHydLow :: Phi-Hydrostatic at r-lower boundary |
120 |
C (bottom in z-coordinates, top in p-coordinates) |
121 |
C totPhiHyd :: total hydrostatic Potential (anomaly, for now), |
122 |
C at cell center level ; includes surface contribution. |
123 |
C (for diagnostic + used in Z-coord with EOS_funct_P) |
124 |
C rhoInSitu :: In-Situ density anomaly [kg/m^3] at cell center level. |
125 |
C hMixLayer :: Mixed layer depth [m] |
126 |
C (for diagnostic + used GMRedi "fm07") |
127 |
C IVDConvCount :: Impl.Vert.Diffusion convection counter: |
128 |
C = 0 (not convecting) or 1 (convecting) |
129 |
COMMON /DYNVARS_DIAG/ |
130 |
& phiHydLow, totPhiHyd, |
131 |
& rhoInSitu, |
132 |
& hMixLayer, IVDConvCount |
133 |
_RL phiHydLow(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
134 |
_RL totPhiHyd(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
135 |
_RL rhoInSitu(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
136 |
_RL hMixLayer(1-OLx:sNx+OLx,1-OLy:sNy+OLy,nSx,nSy) |
137 |
_RL IVDConvCount(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy) |
138 |
|