1 |
C $Header: /u/gcmpack/MITgcm/pkg/obcs/OBCS_FIELDS.h,v 1.4 2012/11/15 15:55:42 dimitri Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#ifdef ALLOW_OBCS |
5 |
|
6 |
CBOP |
7 |
C !ROUTINE: OBCS_FIELDS.h |
8 |
C !INTERFACE: |
9 |
C #include "OBCS_FIELDS.h" |
10 |
|
11 |
C !DESCRIPTION: |
12 |
C *==========================================================* |
13 |
C | OBCS_FIELDS.h |
14 |
C | o Header file containing OB values of model fields |
15 |
C *==========================================================* |
16 |
CEOP |
17 |
|
18 |
#ifdef ALLOW_OBCS_PRESCRIBE |
19 |
C OBCS_ldRec :: time-record currently loaded (in temp arrays *[1]) |
20 |
COMMON /OBCS_LOAD_I/ OBCS_ldRec |
21 |
INTEGER OBCS_ldRec(nSx,nSy) |
22 |
#endif /* ALLOW_OBCS_PRESCRIBE */ |
23 |
|
24 |
C-- COMMON /OBCS_FIELDS/ Open boundary related stuff |
25 |
C OB[N,S,E,W][u,v,w,t,s,eta,am,ph] :: Fields with boundary conditions, |
26 |
C the letter combinations mean: |
27 |
C N/S/E/W :: northern/southern/eastern/western boundary |
28 |
C u/v/w/t/s :: ocean u/v/w velocities, temperature/salinity |
29 |
C eta :: sea surface height |
30 |
C am/ph :: tidal amplitude (m/s) / phase (s) |
31 |
C OBNu is the U value imposed at the Northern OB |
32 |
C OBNv is the V value imposed at the Northern OB |
33 |
C OBNt is the T value imposed at the Northern OB |
34 |
C OBNs is the S value imposed at the Northern OB |
35 |
C etc |
36 |
|
37 |
#ifdef ALLOW_OBCS_NORTH |
38 |
COMMON /OBCS_FIELDS_N/ |
39 |
& OBNu,OBNv,OBNt,OBNs |
40 |
_RL OBNu (1-Olx:sNx+Olx,Nr,nSx,nSy) |
41 |
_RL OBNv (1-Olx:sNx+Olx,Nr,nSx,nSy) |
42 |
_RL OBNt (1-Olx:sNx+Olx,Nr,nSx,nSy) |
43 |
_RL OBNs (1-Olx:sNx+Olx,Nr,nSx,nSy) |
44 |
# ifdef ALLOW_OBCS_PRESCRIBE |
45 |
COMMON /OBCS_FIELDS_AUX_N/ |
46 |
& OBNu0,OBNv0,OBNt0,OBNs0, |
47 |
& OBNu1,OBNv1,OBNt1,OBNs1 |
48 |
_RL OBNu0 (1-Olx:sNx+Olx,Nr,nSx,nSy) |
49 |
_RL OBNv0 (1-Olx:sNx+Olx,Nr,nSx,nSy) |
50 |
_RL OBNt0 (1-Olx:sNx+Olx,Nr,nSx,nSy) |
51 |
_RL OBNs0 (1-Olx:sNx+Olx,Nr,nSx,nSy) |
52 |
_RL OBNu1 (1-Olx:sNx+Olx,Nr,nSx,nSy) |
53 |
_RL OBNv1 (1-Olx:sNx+Olx,Nr,nSx,nSy) |
54 |
_RL OBNt1 (1-Olx:sNx+Olx,Nr,nSx,nSy) |
55 |
_RL OBNs1 (1-Olx:sNx+Olx,Nr,nSx,nSy) |
56 |
# endif /* ALLOW_OBCS_PRESCRIBE */ |
57 |
# ifdef ALLOW_OBCS_STEVENS |
58 |
COMMON /OBCS_FIELDS_STEVENS_N/ |
59 |
& OBNvStevens, OBNtStevens, OBNsStevens |
60 |
_RL OBNvStevens (1-Olx:sNx+Olx,Nr,nSx,nSy) |
61 |
_RL OBNtStevens (1-Olx:sNx+Olx,Nr,nSx,nSy) |
62 |
_RL OBNsStevens (1-Olx:sNx+Olx,Nr,nSx,nSy) |
63 |
# endif /* ALLOW_OBCS_STEVENS */ |
64 |
# ifdef ALLOW_OBCS_TIDES |
65 |
COMMON /OBCS_FIELDS_TIDES_N/ OBNam, OBNph |
66 |
_RL OBNam (1-Olx:sNx+Olx,tidalComponents,nSx,nSy) |
67 |
_RL OBNph (1-Olx:sNx+Olx,tidalComponents,nSx,nSy) |
68 |
# endif /* ALLOW_OBCS_TIDES */ |
69 |
#endif /* ALLOW_OBCS_NORTH */ |
70 |
|
71 |
#ifdef ALLOW_OBCS_SOUTH |
72 |
COMMON /OBCS_FIELDS_S/ |
73 |
& OBSu,OBSv,OBSt,OBSs |
74 |
_RL OBSu (1-Olx:sNx+Olx,Nr,nSx,nSy) |
75 |
_RL OBSv (1-Olx:sNx+Olx,Nr,nSx,nSy) |
76 |
_RL OBSt (1-Olx:sNx+Olx,Nr,nSx,nSy) |
77 |
_RL OBSs (1-Olx:sNx+Olx,Nr,nSx,nSy) |
78 |
# ifdef ALLOW_OBCS_PRESCRIBE |
79 |
COMMON /OBCS_FIELDS_AUX_S/ |
80 |
& OBSu0,OBSv0,OBSt0,OBSs0, |
81 |
& OBSu1,OBSv1,OBSt1,OBSs1 |
82 |
_RL OBSu0 (1-Olx:sNx+Olx,Nr,nSx,nSy) |
83 |
_RL OBSv0 (1-Olx:sNx+Olx,Nr,nSx,nSy) |
84 |
_RL OBSt0 (1-Olx:sNx+Olx,Nr,nSx,nSy) |
85 |
_RL OBSs0 (1-Olx:sNx+Olx,Nr,nSx,nSy) |
86 |
_RL OBSu1 (1-Olx:sNx+Olx,Nr,nSx,nSy) |
87 |
_RL OBSv1 (1-Olx:sNx+Olx,Nr,nSx,nSy) |
88 |
_RL OBSt1 (1-Olx:sNx+Olx,Nr,nSx,nSy) |
89 |
_RL OBSs1 (1-Olx:sNx+Olx,Nr,nSx,nSy) |
90 |
# endif /* ALLOW_OBCS_PRESCRIBE */ |
91 |
# ifdef ALLOW_OBCS_STEVENS |
92 |
COMMON /OBCS_FIELDS_STEVENS_S/ |
93 |
& OBSvStevens, OBStStevens, OBSsStevens |
94 |
_RL OBSvStevens (1-Olx:sNx+Olx,Nr,nSx,nSy) |
95 |
_RL OBStStevens (1-Olx:sNx+Olx,Nr,nSx,nSy) |
96 |
_RL OBSsStevens (1-Olx:sNx+Olx,Nr,nSx,nSy) |
97 |
# endif /* ALLOW_OBCS_STEVENS */ |
98 |
# ifdef ALLOW_OBCS_TIDES |
99 |
COMMON /OBCS_FIELDS_TIDES_S/ OBSam, OBSph |
100 |
_RL OBSam (1-Olx:sNx+Olx,tidalComponents,nSx,nSy) |
101 |
_RL OBSph (1-Olx:sNx+Olx,tidalComponents,nSx,nSy) |
102 |
# endif /* ALLOW_OBCS_TIDES */ |
103 |
#endif /* ALLOW_OBCS_SOUTH */ |
104 |
|
105 |
#ifdef ALLOW_OBCS_EAST |
106 |
COMMON /OBCS_FIELDS_E/ |
107 |
& OBEu,OBEv,OBEt,OBEs |
108 |
_RL OBEu (1-Oly:sNy+Oly,Nr,nSx,nSy) |
109 |
_RL OBEv (1-Oly:sNy+Oly,Nr,nSx,nSy) |
110 |
_RL OBEt (1-Oly:sNy+Oly,Nr,nSx,nSy) |
111 |
_RL OBEs (1-Oly:sNy+Oly,Nr,nSx,nSy) |
112 |
# ifdef ALLOW_OBCS_PRESCRIBE |
113 |
COMMON /OBCS_FIELDS_AUX_E/ |
114 |
& OBEu0,OBEv0,OBEt0,OBEs0, |
115 |
& OBEu1,OBEv1,OBEt1,OBEs1 |
116 |
_RL OBEu0 (1-Oly:sNy+Oly,Nr,nSx,nSy) |
117 |
_RL OBEv0 (1-Oly:sNy+Oly,Nr,nSx,nSy) |
118 |
_RL OBEt0 (1-Oly:sNy+Oly,Nr,nSx,nSy) |
119 |
_RL OBEs0 (1-Oly:sNy+Oly,Nr,nSx,nSy) |
120 |
_RL OBEu1 (1-Oly:sNy+Oly,Nr,nSx,nSy) |
121 |
_RL OBEv1 (1-Oly:sNy+Oly,Nr,nSx,nSy) |
122 |
_RL OBEt1 (1-Oly:sNy+Oly,Nr,nSx,nSy) |
123 |
_RL OBEs1 (1-Oly:sNy+Oly,Nr,nSx,nSy) |
124 |
# endif /* ALLOW_OBCS_PRESCRIBE */ |
125 |
# ifdef ALLOW_OBCS_STEVENS |
126 |
COMMON /OBCS_FIELDS_STEVENS_E/ |
127 |
& OBEuStevens, OBEtStevens, OBEsStevens |
128 |
_RL OBEuStevens (1-Oly:sNy+Oly,Nr,nSx,nSy) |
129 |
_RL OBEtStevens (1-Oly:sNy+Oly,Nr,nSx,nSy) |
130 |
_RL OBEsStevens (1-Oly:sNy+Oly,Nr,nSx,nSy) |
131 |
# endif /* ALLOW_OBCS_STEVENS */ |
132 |
# ifdef ALLOW_OBCS_TIDES |
133 |
COMMON /OBCS_FIELDS_TIDES_E/ OBEam, OBEph |
134 |
_RL OBEam (1-Oly:sNy+Oly,tidalComponents,nSx,nSy) |
135 |
_RL OBEph (1-Oly:sNy+Oly,tidalComponents,nSx,nSy) |
136 |
# endif /* ALLOW_OBCS_TIDES */ |
137 |
#endif /* ALLOW_OBCS_EAST */ |
138 |
|
139 |
#ifdef ALLOW_OBCS_WEST |
140 |
COMMON /OBCS_FIELDS_W/ |
141 |
& OBWu,OBWv,OBWt,OBWs |
142 |
_RL OBWu (1-Oly:sNy+Oly,Nr,nSx,nSy) |
143 |
_RL OBWv (1-Oly:sNy+Oly,Nr,nSx,nSy) |
144 |
_RL OBWt (1-Oly:sNy+Oly,Nr,nSx,nSy) |
145 |
_RL OBWs (1-Oly:sNy+Oly,Nr,nSx,nSy) |
146 |
# ifdef ALLOW_OBCS_PRESCRIBE |
147 |
COMMON /OBCS_FIELDS_AUX_W/ |
148 |
& OBWu0,OBWv0,OBWt0,OBWs0, |
149 |
& OBWu1,OBWv1,OBWt1,OBWs1 |
150 |
_RL OBWu0 (1-Oly:sNy+Oly,Nr,nSx,nSy) |
151 |
_RL OBWv0 (1-Oly:sNy+Oly,Nr,nSx,nSy) |
152 |
_RL OBWt0 (1-Oly:sNy+Oly,Nr,nSx,nSy) |
153 |
_RL OBWs0 (1-Oly:sNy+Oly,Nr,nSx,nSy) |
154 |
_RL OBWu1 (1-Oly:sNy+Oly,Nr,nSx,nSy) |
155 |
_RL OBWv1 (1-Oly:sNy+Oly,Nr,nSx,nSy) |
156 |
_RL OBWt1 (1-Oly:sNy+Oly,Nr,nSx,nSy) |
157 |
_RL OBWs1 (1-Oly:sNy+Oly,Nr,nSx,nSy) |
158 |
# endif /* ALLOW_OBCS_PRESCRIBE */ |
159 |
# ifdef ALLOW_OBCS_STEVENS |
160 |
COMMON /OBCS_FIELDS_STEVENS_W/ |
161 |
& OBWuStevens, OBWtStevens, OBWsStevens |
162 |
_RL OBWuStevens (1-Oly:sNy+Oly,Nr,nSx,nSy) |
163 |
_RL OBWtStevens (1-Oly:sNy+Oly,Nr,nSx,nSy) |
164 |
_RL OBWsStevens (1-Oly:sNy+Oly,Nr,nSx,nSy) |
165 |
# endif /* ALLOW_OBCS_STEVENS */ |
166 |
# ifdef ALLOW_OBCS_TIDES |
167 |
COMMON /OBCS_FIELDS_TIDES_W/ OBWam, OBWph |
168 |
_RL OBWam (1-Oly:sNy+Oly,tidalComponents,nSx,nSy) |
169 |
_RL OBWph (1-Oly:sNy+Oly,tidalComponents,nSx,nSy) |
170 |
# endif /* ALLOW_OBCS_TIDES */ |
171 |
#endif /* ALLOW_OBCS_WEST */ |
172 |
|
173 |
#ifdef ALLOW_NONHYDROSTATIC |
174 |
COMMON /OBCS_NH_FIELDS/ |
175 |
& OBNw, OBSw, OBEw, OBWw |
176 |
_RL OBNw (1-Olx:sNx+Olx,Nr,nSx,nSy) |
177 |
_RL OBSw (1-Olx:sNx+Olx,Nr,nSx,nSy) |
178 |
_RL OBEw (1-Oly:sNy+Oly,Nr,nSx,nSy) |
179 |
_RL OBWw (1-Oly:sNy+Oly,Nr,nSx,nSy) |
180 |
#ifdef ALLOW_OBCS_PRESCRIBE |
181 |
COMMON /OBCS_NH_FIELDS_AUX/ |
182 |
& OBNw0, OBSw0, OBEw0, OBWw0, |
183 |
& OBNw1, OBSw1, OBEw1, OBWw1 |
184 |
_RL OBNw0(1-Olx:sNx+Olx,Nr,nSx,nSy) |
185 |
_RL OBSw0(1-Olx:sNx+Olx,Nr,nSx,nSy) |
186 |
_RL OBEw0(1-Oly:sNy+Oly,Nr,nSx,nSy) |
187 |
_RL OBWw0(1-Oly:sNy+Oly,Nr,nSx,nSy) |
188 |
_RL OBNw1(1-Olx:sNx+Olx,Nr,nSx,nSy) |
189 |
_RL OBSw1(1-Olx:sNx+Olx,Nr,nSx,nSy) |
190 |
_RL OBEw1(1-Oly:sNy+Oly,Nr,nSx,nSy) |
191 |
_RL OBWw1(1-Oly:sNy+Oly,Nr,nSx,nSy) |
192 |
#endif /* ALLOW_OBCS_PRESCRIBE */ |
193 |
#endif /* ALLOW_NONHYDROSTATIC */ |
194 |
|
195 |
#ifdef NONLIN_FRSURF |
196 |
COMMON /OBCS_NLFS_FIELDS/ |
197 |
& OBNeta, OBSeta, OBEeta, OBWeta |
198 |
_RL OBNeta (1-Olx:sNx+Olx,nSx,nSy) |
199 |
_RL OBSeta (1-Olx:sNx+Olx,nSx,nSy) |
200 |
_RL OBEeta (1-Oly:sNy+Oly,nSx,nSy) |
201 |
_RL OBWeta (1-Oly:sNy+Oly,nSx,nSy) |
202 |
#ifdef ALLOW_OBCS_PRESCRIBE |
203 |
COMMON /OBCS_NLFS_FIELDS_AUX/ |
204 |
& OBNeta0,OBSeta0,OBEeta0,OBWeta0, |
205 |
& OBNeta1,OBSeta1,OBEeta1,OBWeta1 |
206 |
_RL OBNeta0(1-Olx:sNx+Olx,nSx,nSy) |
207 |
_RL OBSeta0(1-Olx:sNx+Olx,nSx,nSy) |
208 |
_RL OBEeta0(1-Oly:sNy+Oly,nSx,nSy) |
209 |
_RL OBWeta0(1-Oly:sNy+Oly,nSx,nSy) |
210 |
_RL OBNeta1(1-Olx:sNx+Olx,nSx,nSy) |
211 |
_RL OBSeta1(1-Olx:sNx+Olx,nSx,nSy) |
212 |
_RL OBEeta1(1-Oly:sNy+Oly,nSx,nSy) |
213 |
_RL OBWeta1(1-Oly:sNy+Oly,nSx,nSy) |
214 |
#endif /* ALLOW_OBCS_PRESCRIBE */ |
215 |
#endif /* NONLIN_FRSURF */ |
216 |
|
217 |
#endif /* ALLOW_OBCS */ |