1 | MODULE obc_oce |
---|
2 | !!============================================================================== |
---|
3 | !! *** MODULE obc_oce *** |
---|
4 | !! Open Boundary Cond. : define related variables |
---|
5 | !!============================================================================== |
---|
6 | !! OPA 9.0 , LOCEAN-IPSL (2005) |
---|
7 | !! $Id$ |
---|
8 | !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt) |
---|
9 | !!---------------------------------------------------------------------- |
---|
10 | !!---------------------------------------------------------------------- |
---|
11 | !! 'key_obc' : Open Boundary Condition |
---|
12 | !!---------------------------------------------------------------------- |
---|
13 | !! history : |
---|
14 | !! 8.0 01/91 (CLIPPER) Original code |
---|
15 | !! 8.5 06/02 (C. Talandier) modules |
---|
16 | !! 06/04 (F. Durand) ORCA2_ZIND config |
---|
17 | !! |
---|
18 | !!---------------------------------------------------------------------- |
---|
19 | !! * Modules used |
---|
20 | USE par_oce ! ocean parameters |
---|
21 | USE obc_par ! open boundary condition parameters |
---|
22 | |
---|
23 | #if defined key_obc |
---|
24 | |
---|
25 | IMPLICIT NONE |
---|
26 | PUBLIC |
---|
27 | |
---|
28 | !!---------------------------------------------------------------------- |
---|
29 | !! open boundary variables |
---|
30 | !!---------------------------------------------------------------------- |
---|
31 | !! |
---|
32 | !!General variables for open boundaries: |
---|
33 | !!-------------------------------------- |
---|
34 | INTEGER :: & !: * namelist ??? * |
---|
35 | nbobc = 2 , & !: number of open boundaries ( 1=< nbobc =< 4 ) |
---|
36 | nobc_dta = 0 !: = 0 use the initial state as obc data |
---|
37 | ! ! = 1 read obc data in obcxxx.dta files |
---|
38 | |
---|
39 | LOGICAL :: ln_obc_clim = .true. !: obc data files are climatological |
---|
40 | LOGICAL :: ln_obc_fla = .false. !: Flather open boundary condition not used |
---|
41 | LOGICAL :: ln_vol_cst = .true. !: Conservation of the whole volume |
---|
42 | |
---|
43 | REAL(wp) :: & !!: open boundary namelist (namobc) |
---|
44 | rdpein = 1. , & !: damping time scale for inflow at East open boundary |
---|
45 | rdpwin = 1. , & !: " " at West open boundary |
---|
46 | rdpsin = 1. , & !: " " at South open boundary |
---|
47 | rdpnin = 1. , & !: " " at North open boundary |
---|
48 | rdpeob = 15. , & !: damping time scale for the climatology at East open boundary |
---|
49 | rdpwob = 15. , & !: " " at West open boundary |
---|
50 | rdpsob = 15. , & !: " " at South open boundary |
---|
51 | rdpnob = 15. , & !: " " at North open boundary |
---|
52 | volemp = 1. !: = 0 the total volume will have the variability of the |
---|
53 | ! surface Flux E-P else (volemp = 1) the volume will be constant |
---|
54 | ! = 1 the volume will be constant during all the integration. |
---|
55 | |
---|
56 | LOGICAL :: & !: |
---|
57 | lfbceast, lfbcwest, & !: logical flag for a fixed East and West open boundaries |
---|
58 | lfbcnorth, lfbcsouth !: logical flag for a fixed North and South open boundaries |
---|
59 | ! ! These logical flags are set to 'true' if damping time |
---|
60 | ! ! scale are set to 0 in the namelist, for both inflow and outflow). |
---|
61 | |
---|
62 | REAL(wp), PUBLIC :: & !: |
---|
63 | obcsurftot !: Total lateral surface of open boundaries |
---|
64 | |
---|
65 | REAL(wp), PUBLIC, DIMENSION(jpi,jpj) :: & !: |
---|
66 | obctmsk, & !: mask array identical to tmask, execpt along OBC where it is set to 0 |
---|
67 | ! ! it used to calculate the cumulate flux E-P in the obcvol.F90 routine |
---|
68 | obcumask, obcvmask !: u-, v- Force filtering mask for the open |
---|
69 | ! ! boundary condition on grad D |
---|
70 | |
---|
71 | !!-------------------- |
---|
72 | !! East open boundary: |
---|
73 | !!-------------------- |
---|
74 | INTEGER :: nie0 , nie1 !: do loop index in mpp case for jpieob |
---|
75 | INTEGER :: nie0p1, nie1p1 !: do loop index in mpp case for jpieob+1 |
---|
76 | INTEGER :: nie0m1, nie1m1 !: do loop index in mpp case for jpieob-1 |
---|
77 | INTEGER :: nje0 , nje1 !: do loop index in mpp case for jpjed, jpjef |
---|
78 | INTEGER :: nje0p1, nje1m1 !: do loop index in mpp case for jpjedp1,jpjefm1 |
---|
79 | INTEGER :: nje1m2, nje0m1 !: do loop index in mpp case for jpjefm1-1,jpjed |
---|
80 | |
---|
81 | REAL(wp), DIMENSION(jpj) :: & !: |
---|
82 | bsfeob !: now barotropic stream fuction computed at the OBC. The corres- |
---|
83 | ! ! ponding bsfn will be computed by the forward time step in dynspg. |
---|
84 | |
---|
85 | REAL(wp), DIMENSION(jpj,3,3) :: & !: |
---|
86 | bebnd !: east boundary barotropic streamfunction over 3 rows |
---|
87 | ! ! and 3 time step (now, before, and before before) |
---|
88 | |
---|
89 | REAL(wp), DIMENSION(jpjed:jpjef) :: & !: |
---|
90 | bfoe, & !: now climatology of the east boundary barotropic stream function |
---|
91 | sshfoe, & !: now climatology of the east boundary sea surface height |
---|
92 | ubtfoe,vbtfoe !: now climatology of the east boundary barotropic transport |
---|
93 | |
---|
94 | REAL(wp), DIMENSION(jpj,jpk) :: & !: |
---|
95 | ufoe, vfoe, & !: now climatology of the east boundary velocities |
---|
96 | tfoe, sfoe, & !: now climatology of the east boundary temperature and salinity |
---|
97 | uclie !: baroclinic componant of the zonal velocity after radiation |
---|
98 | ! ! in the obcdyn.F90 routine |
---|
99 | |
---|
100 | REAL(wp), DIMENSION(jpjed:jpjef,jpj) :: & !: |
---|
101 | sshfoe_b !: east boundary ssh correction averaged over the barotropic loop |
---|
102 | !: (if Flather's algoritm applied at open boundary) |
---|
103 | |
---|
104 | !!------------------------------- |
---|
105 | !! Arrays for radiative East OBC: |
---|
106 | !!------------------------------- |
---|
107 | REAL(wp), DIMENSION(jpj,jpk,3,3) :: & !: |
---|
108 | uebnd, vebnd !: baroclinic u & v component of the velocity over 3 rows |
---|
109 | ! and 3 time step (now, before, and before before) |
---|
110 | |
---|
111 | REAL(wp), DIMENSION(jpj,jpk,2,2) :: & !: |
---|
112 | tebnd, sebnd !: East boundary temperature and salinity over 2 rows |
---|
113 | ! and 2 time step (now and before) |
---|
114 | |
---|
115 | REAL(wp), DIMENSION(jpj,jpk) :: & !: |
---|
116 | u_cxebnd, v_cxebnd !: Zonal component of the phase speed ratio computed with |
---|
117 | ! radiation of u and v velocity (respectively) at the |
---|
118 | ! east open boundary (u_cxebnd = cx rdt ) |
---|
119 | |
---|
120 | REAL(wp), DIMENSION(jpj,jpk) :: & !: |
---|
121 | uemsk, vemsk, temsk !: 2D mask for the East OB |
---|
122 | |
---|
123 | ! Note that those arrays are optimized for mpp case |
---|
124 | ! (hence the dimension jpj is the size of one processor subdomain) |
---|
125 | |
---|
126 | !!-------------------- |
---|
127 | !! West open boundary |
---|
128 | !!-------------------- |
---|
129 | INTEGER :: niw0 , niw1 !: do loop index in mpp case for jpiwob |
---|
130 | INTEGER :: niw0p1, niw1p1 !: do loop index in mpp case for jpiwob+1 |
---|
131 | INTEGER :: njw0 , njw1 !: do loop index in mpp case for jpjwd, jpjwf |
---|
132 | INTEGER :: njw0p1, njw1m1 !: do loop index in mpp case for jpjwdp1,jpjwfm1 |
---|
133 | INTEGER :: njw1m2, njw0m1 !: do loop index in mpp case for jpjwfm2,jpjwd |
---|
134 | |
---|
135 | REAL(wp), DIMENSION(jpj) :: & !: |
---|
136 | bsfwob !: now barotropic stream fuction computed at the OBC. The corres- |
---|
137 | ! ! ponding bsfn will be computed by the forward time step in dynspg. |
---|
138 | |
---|
139 | REAL(wp), DIMENSION(jpj,3,3) :: & !: |
---|
140 | bwbnd !: West boundary barotropic streamfunction over |
---|
141 | ! ! 3 rows and 3 time step (now, before, and before before) |
---|
142 | |
---|
143 | REAL(wp), DIMENSION(jpjwd:jpjwf) :: & !: |
---|
144 | bfow, & !: now climatology of the west boundary barotropic stream function |
---|
145 | sshfow, & !: now climatology of the west boundary sea surface height |
---|
146 | ubtfow,vbtfow !: now climatology of the west boundary barotropic transport |
---|
147 | |
---|
148 | REAL(wp), DIMENSION(jpj,jpk) :: & !: |
---|
149 | ufow, vfow, & !: now climatology of the west velocities |
---|
150 | tfow, sfow, & !: now climatology of the west temperature and salinity |
---|
151 | ucliw !: baroclinic componant of the zonal velocity after the radiation |
---|
152 | ! ! in the obcdyn.F90 routine |
---|
153 | |
---|
154 | REAL(wp), DIMENSION(jpjwd:jpjwf,jpj) :: & !: |
---|
155 | sshfow_b !: west boundary ssh correction averaged over the barotropic loop |
---|
156 | !: (if Flather's algoritm applied at open boundary) |
---|
157 | |
---|
158 | !!------------------------------- |
---|
159 | !! Arrays for radiative West OBC |
---|
160 | !!------------------------------- |
---|
161 | REAL(wp), DIMENSION(jpj,jpk,3,3) :: & !: |
---|
162 | uwbnd, vwbnd !: baroclinic u & v components of the velocity over 3 rows |
---|
163 | ! ! and 3 time step (now, before, and before before) |
---|
164 | |
---|
165 | REAL(wp), DIMENSION(jpj,jpk,2,2) :: & !: |
---|
166 | twbnd, swbnd !: west boundary temperature and salinity over 2 rows and |
---|
167 | ! ! 2 time step (now and before) |
---|
168 | |
---|
169 | REAL(wp), DIMENSION(jpj,jpk) :: & !: |
---|
170 | u_cxwbnd, v_cxwbnd !: Zonal component of the phase speed ratio computed with |
---|
171 | ! ! radiation of zonal and meridional velocity (respectively) |
---|
172 | ! ! at the west open boundary (u_cxwbnd = cx rdt ) |
---|
173 | |
---|
174 | REAL(wp), DIMENSION(jpj,jpk) :: & !: |
---|
175 | uwmsk, vwmsk, twmsk !: 2D mask for the West OB |
---|
176 | |
---|
177 | ! Note that those arrays are optimized for mpp case |
---|
178 | ! (hence the dimension jpj is the size of one processor subdomain) |
---|
179 | |
---|
180 | !!--------------------- |
---|
181 | !! North open boundary |
---|
182 | !!--------------------- |
---|
183 | INTEGER :: nin0 , nin1 !: do loop index in mpp case for jpind, jpinf |
---|
184 | INTEGER :: nin0p1, nin1m1 !: do loop index in mpp case for jpindp1, jpinfm1 |
---|
185 | INTEGER :: nin1m2, nin0m1 !: do loop index in mpp case for jpinfm1-1,jpind |
---|
186 | INTEGER :: njn0 , njn1 !: do loop index in mpp case for jpnob |
---|
187 | INTEGER :: njn0p1, njn1p1 !: do loop index in mpp case for jpnob+1 |
---|
188 | INTEGER :: njn0m1, njn1m1 !: do loop index in mpp case for jpnob-1 |
---|
189 | |
---|
190 | REAL(wp), DIMENSION(jpi) :: & !: |
---|
191 | bsfnob !: now barotropic stream fuction computed at the OBC. The corres- |
---|
192 | ! ! ponding bsfn will be computed by the forward time step in dynspg. |
---|
193 | |
---|
194 | REAL(wp), DIMENSION(jpi,3,3) :: & !: |
---|
195 | bnbnd !: north boundary barotropic streamfunction over |
---|
196 | ! ! 3 rows and 3 time step (now, before, and before before) |
---|
197 | |
---|
198 | REAL(wp), DIMENSION(jpind:jpinf) :: & !: |
---|
199 | bfon, & !: now climatology of the north boundary barotropic stream function |
---|
200 | sshfon, & !: now climatology of the north boundary sea surface height |
---|
201 | ubtfon,vbtfon !: now climatology of the north boundary barotropic transport |
---|
202 | |
---|
203 | REAL(wp), DIMENSION(jpi,jpk) :: & !: |
---|
204 | ufon, vfon, & !: now climatology of the north boundary velocities |
---|
205 | tfon, sfon, & !: now climatology of the north boundary temperature and salinity |
---|
206 | vclin !: baroclinic componant of the meridian velocity after the radiation |
---|
207 | ! ! in yhe obcdyn.F90 routine |
---|
208 | |
---|
209 | REAL(wp), DIMENSION(jpind:jpinf,jpj) :: & !: |
---|
210 | sshfon_b !: north boundary ssh correction averaged over the barotropic loop |
---|
211 | !: (if Flather's algoritm applied at open boundary) |
---|
212 | |
---|
213 | !!-------------------------------- |
---|
214 | !! Arrays for radiative North OBC |
---|
215 | !!-------------------------------- |
---|
216 | !! |
---|
217 | REAL(wp), DIMENSION(jpi,jpk,3,3) :: & !: |
---|
218 | unbnd, vnbnd !: baroclinic u & v components of the velocity over 3 |
---|
219 | ! ! rows and 3 time step (now, before, and before before) |
---|
220 | |
---|
221 | REAL(wp), DIMENSION(jpi,jpk,2,2) :: & !: |
---|
222 | tnbnd, snbnd !: north boundary temperature and salinity over |
---|
223 | ! ! 2 rows and 2 time step (now and before) |
---|
224 | |
---|
225 | REAL(wp), DIMENSION(jpi,jpk) :: & !: |
---|
226 | u_cynbnd, v_cynbnd !: Meridional component of the phase speed ratio compu- |
---|
227 | ! ! ted with radiation of zonal and meridional velocity |
---|
228 | ! ! (respectively) at the north OB (u_cynbnd = cx rdt ) |
---|
229 | |
---|
230 | REAL(wp), DIMENSION(jpi,jpk) :: & !: |
---|
231 | unmsk, vnmsk, tnmsk !: 2D mask for the North OB |
---|
232 | |
---|
233 | ! Note that those arrays are optimized for mpp case |
---|
234 | ! (hence the dimension jpj is the size of one processor subdomain) |
---|
235 | |
---|
236 | !!--------------------- |
---|
237 | !! South open boundary |
---|
238 | !!--------------------- |
---|
239 | INTEGER :: nis0 , nis1 !: do loop index in mpp case for jpisd, jpisf |
---|
240 | INTEGER :: nis0p1, nis1m1 !: do loop index in mpp case for jpisdp1, jpisfm1 |
---|
241 | INTEGER :: nis1m2, nis0m1 !: do loop index in mpp case for jpisfm1-1,jpisd |
---|
242 | INTEGER :: njs0 , njs1 !: do loop index in mpp case for jpsob |
---|
243 | INTEGER :: njs0p1, njs1p1 !: do loop index in mpp case for jpsob+1 |
---|
244 | |
---|
245 | REAL(wp), DIMENSION(jpi) :: & !: |
---|
246 | bsfsob !: now barotropic stream fuction computed at the OBC.The corres- |
---|
247 | ! ! ponding bsfn will be computed by the forward time step in dynspg. |
---|
248 | REAL(wp), DIMENSION(jpi,3,3) :: & !: |
---|
249 | bsbnd !: south boundary barotropic stream function over |
---|
250 | ! ! 3 rows and 3 time step (now, before, and before before) |
---|
251 | |
---|
252 | REAL(wp), DIMENSION(jpisd:jpisf) :: & !: |
---|
253 | bfos, & !: now climatology of the south boundary barotropic stream function |
---|
254 | sshfos, & !: now climatology of the south boundary sea surface height |
---|
255 | ubtfos,vbtfos !: now climatology of the south boundary barotropic transport |
---|
256 | |
---|
257 | REAL(wp), DIMENSION(jpi,jpk) :: & !: |
---|
258 | ufos, vfos, & !: now climatology of the south boundary velocities |
---|
259 | tfos, sfos, & !: now climatology of the south boundary temperature and salinity |
---|
260 | vclis !: baroclinic componant of the meridian velocity after the radiation |
---|
261 | ! ! in the obcdyn.F90 routine |
---|
262 | |
---|
263 | REAL(wp), DIMENSION(jpisd:jpisf,jpj) :: & !: |
---|
264 | sshfos_b !: south boundary ssh correction averaged over the barotropic loop |
---|
265 | !: (if Flather's algoritm applied at open boundary) |
---|
266 | |
---|
267 | !!-------------------------------- |
---|
268 | !! Arrays for radiative South OBC |
---|
269 | !!-------------------------------- |
---|
270 | !! computed by the forward time step in dynspg. |
---|
271 | REAL(wp), DIMENSION(jpi,jpk,3,3) :: & !: |
---|
272 | usbnd, vsbnd !: baroclinic u & v components of the velocity over 3 |
---|
273 | ! ! rows and 3 time step (now, before, and before before) |
---|
274 | |
---|
275 | REAL(wp), DIMENSION(jpi,jpk,2,2) :: & !: |
---|
276 | tsbnd, ssbnd !: south boundary temperature and salinity over |
---|
277 | ! ! 2 rows and 2 time step (now and before) |
---|
278 | |
---|
279 | REAL(wp), DIMENSION(jpi,jpk) :: & !: |
---|
280 | u_cysbnd, v_cysbnd !: Meridional component of the phase speed ratio compu- |
---|
281 | ! ! ted with radiation of zonal and meridional velocity |
---|
282 | ! ! (repsectively) at the south OB (u_cynbnd = cx rdt ) |
---|
283 | |
---|
284 | REAL(wp), DIMENSION(jpi,jpk) :: & !: |
---|
285 | usmsk, vsmsk, tsmsk !: 2D mask for the South OB |
---|
286 | |
---|
287 | CHARACTER ( len=20 ) :: cffile |
---|
288 | |
---|
289 | #else |
---|
290 | !!---------------------------------------------------------------------- |
---|
291 | !! Default option : Empty module |
---|
292 | !!---------------------------------------------------------------------- |
---|
293 | #endif |
---|
294 | |
---|
295 | !!====================================================================== |
---|
296 | END MODULE obc_oce |
---|