New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
obc_oce.F90 in trunk/NEMO/OPA_SRC/OBC – NEMO

source: trunk/NEMO/OPA_SRC/OBC/obc_oce.F90 @ 1158

Last change on this file since 1158 was 1158, checked in by rblod, 16 years ago

OBC corrections, see ticket #224

  • Property svn:eol-style set to native
  • Property svn:keywords set to Id
File size: 16.8 KB
Line 
1MODULE 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   !! Rigid lid case:
73   !!----------------
74   INTEGER ::   nbic !: number of isolated coastlines ( 0 <= nbic <= 3 )
75         
76   INTEGER, DIMENSION(jpnic,0:4,3) ::   &  !:
77      miic, mjic     !: position of isolated coastlines points
78
79   INTEGER, DIMENSION(0:4,3) ::   &  !:
80      mnic           !: number of points on isolated coastlines
81
82   REAL(wp), DIMENSION(jpi,jpj) ::   &  !:
83      gcbob          !: right hand side of the barotropic elliptic equation associated
84      !              !  with the OBC
85                                             
86   REAL(wp), DIMENSION(jpi,jpj,3) ::   &  !:
87      gcfobc         !: coef. associated with the contribution of isolated coastlines
88      !              !  to the right hand side of the barotropic elliptic equation
89
90   REAL(wp), DIMENSION(3) ::   &  !:
91      gcbic          !: time variation of the barotropic stream function along the
92      !              !  isolated coastlines
93
94   REAL(wp), DIMENSION(1) ::   &  !:
95      bsfic0         !: barotropic stream function on isolated coastline
96         
97   REAL(wp), DIMENSION(3) ::   &  !:
98      bsfic          !: barotropic stream function on isolated coastline
99         
100   !!--------------------
101   !! East open boundary:
102   !!--------------------
103   INTEGER ::   nie0  , nie1      !: do loop index in mpp case for jpieob
104   INTEGER ::   nie0p1, nie1p1    !: do loop index in mpp case for jpieob+1
105   INTEGER ::   nie0m1, nie1m1    !: do loop index in mpp case for jpieob-1
106   INTEGER ::   nje0  , nje1      !: do loop index in mpp case for jpjed, jpjef
107   INTEGER ::   nje0p1, nje1m1    !: do loop index in mpp case for jpjedp1,jpjefm1
108   INTEGER ::   nje1m2, nje0m1    !: do loop index in mpp case for jpjefm1-1,jpjed
109
110   REAL(wp), DIMENSION(jpj) ::    &  !:
111      bsfeob              !: now barotropic stream fuction computed at the OBC. The corres-
112      !                   ! ponding bsfn will be computed by the forward time step in dynspg.
113
114   REAL(wp), DIMENSION(jpj,3,3) ::   &  !:
115      bebnd               !: east boundary barotropic streamfunction over 3 rows
116      !                   ! and 3 time step (now, before, and before before)
117
118   REAL(wp), DIMENSION(jpjed:jpjef) ::   &  !:
119      bfoe,             & !: now climatology of the east boundary barotropic stream function
120      sshfoe,           & !: now climatology of the east boundary sea surface height
121      ubtfoe,vbtfoe       !: now climatology of the east boundary barotropic transport
122     
123   REAL(wp), DIMENSION(jpj,jpk) ::   &  !:
124      ufoe, vfoe,       & !: now climatology of the east boundary velocities
125      tfoe, sfoe,       & !: now climatology of the east boundary temperature and salinity
126      uclie               !: baroclinic componant of the zonal velocity after radiation
127      !                   ! in the obcdyn.F90 routine
128
129   REAL(wp), DIMENSION(jpjed:jpjef,jpj) ::   &  !:
130      sshfoe_b            !: east boundary ssh correction averaged over the barotropic loop
131                          !: (if Flather's algoritm applied at open boundary)
132
133   !!-------------------------------
134   !! Arrays for radiative East OBC:
135   !!-------------------------------
136   REAL(wp), DIMENSION(jpj,jpk,3,3) ::   &  !:
137      uebnd, vebnd                  !: baroclinic u & v component of the velocity over 3 rows
138                                    ! and 3 time step (now, before, and before before)
139
140   REAL(wp), DIMENSION(jpj,jpk,2,2) ::   &  !:
141      tebnd, sebnd                  !: East boundary temperature and salinity over 2 rows
142                                    ! and 2 time step (now and before)
143
144   REAL(wp), DIMENSION(jpj,jpk) ::   &  !:
145      u_cxebnd, v_cxebnd            !: Zonal component of the phase speed ratio computed with
146                                    ! radiation of u and v velocity (respectively) at the
147                                    ! east open boundary (u_cxebnd = cx rdt )
148
149   REAL(wp), DIMENSION(jpj,jpk) ::   &  !:
150      uemsk, vemsk, temsk           !: 2D mask for the East OB
151
152   ! Note that those arrays are optimized for mpp case
153   ! (hence the dimension jpj is the size of one processor subdomain)
154
155   !!--------------------
156   !! West open boundary
157   !!--------------------
158   INTEGER ::   niw0  , niw1       !: do loop index in mpp case for jpiwob
159   INTEGER ::   niw0p1, niw1p1     !: do loop index in mpp case for jpiwob+1
160   INTEGER ::   njw0  , njw1       !: do loop index in mpp case for jpjwd, jpjwf
161   INTEGER ::   njw0p1, njw1m1     !: do loop index in mpp case for jpjwdp1,jpjwfm1
162   INTEGER ::   njw1m2, njw0m1     !: do loop index in mpp case for jpjwfm2,jpjwd
163
164   REAL(wp), DIMENSION(jpj) ::   &  !:
165      bsfwob              !: now barotropic stream fuction computed at the OBC. The corres-
166      !                   !  ponding bsfn will be computed by the forward time step in dynspg.
167
168   REAL(wp), DIMENSION(jpj,3,3) ::   &  !:
169      bwbnd               !: West boundary barotropic streamfunction over
170      !                   !  3 rows and 3 time step (now, before, and before before)
171
172   REAL(wp), DIMENSION(jpjwd:jpjwf) ::   &  !:
173      bfow,             & !: now climatology of the west boundary barotropic stream function
174      sshfow,           & !: now climatology of the west boundary sea surface height
175      ubtfow,vbtfow       !: now climatology of the west boundary barotropic transport
176
177   REAL(wp), DIMENSION(jpj,jpk) ::   &  !:
178      ufow, vfow,       & !: now climatology of the west velocities
179      tfow, sfow,       & !: now climatology of the west temperature and salinity
180      ucliw               !: baroclinic componant of the zonal velocity after the radiation
181      !                   !  in the obcdyn.F90 routine
182
183   REAL(wp), DIMENSION(jpjwd:jpjwf,jpj) ::   &  !:
184      sshfow_b            !: west boundary ssh correction averaged over the barotropic loop
185                          !: (if Flather's algoritm applied at open boundary)
186
187   !!-------------------------------
188   !! Arrays for radiative West OBC
189   !!-------------------------------
190   REAL(wp), DIMENSION(jpj,jpk,3,3) ::   &  !:
191      uwbnd, vwbnd                  !: baroclinic u & v components of the velocity over 3 rows
192      !                             !  and 3 time step (now, before, and before before)
193
194   REAL(wp), DIMENSION(jpj,jpk,2,2) ::   &  !:
195      twbnd, swbnd                  !: west boundary temperature and salinity over 2 rows and
196      !                             !  2 time step (now and before)
197
198   REAL(wp), DIMENSION(jpj,jpk) ::    &  !:
199      u_cxwbnd, v_cxwbnd            !: Zonal component of the phase speed ratio computed with
200      !                             !  radiation of zonal and meridional velocity (respectively)
201      !                             !  at the west open boundary (u_cxwbnd = cx rdt )
202
203   REAL(wp), DIMENSION(jpj,jpk) ::    &  !:
204      uwmsk, vwmsk, twmsk           !: 2D mask for the West OB
205
206   ! Note that those arrays are optimized for mpp case
207   ! (hence the dimension jpj is the size of one processor subdomain)
208
209   !!---------------------
210   !! North open boundary
211   !!---------------------
212   INTEGER ::   nin0  , nin1       !: do loop index in mpp case for jpind, jpinf
213   INTEGER ::   nin0p1, nin1m1     !: do loop index in mpp case for jpindp1, jpinfm1
214   INTEGER ::   nin1m2, nin0m1     !: do loop index in mpp case for jpinfm1-1,jpind
215   INTEGER ::   njn0  , njn1       !: do loop index in mpp case for jpnob
216   INTEGER ::   njn0p1, njn1p1     !: do loop index in mpp case for jpnob+1
217   INTEGER ::   njn0m1, njn1m1     !: do loop index in mpp case for jpnob-1
218
219   REAL(wp), DIMENSION(jpi) ::   &  !:
220      bsfnob              !: now barotropic stream fuction computed at the OBC. The corres-
221      !                   !  ponding bsfn will be computed by the forward time step in dynspg.
222
223   REAL(wp), DIMENSION(jpi,3,3) ::   &  !:
224      bnbnd               !: north boundary barotropic streamfunction over
225      !                   !  3 rows and 3 time step (now, before, and before before)
226
227   REAL(wp), DIMENSION(jpind:jpinf) ::   &  !:
228      bfon,             & !: now climatology of the north boundary barotropic stream function
229      sshfon,           & !: now climatology of the north boundary sea surface height
230      ubtfon,vbtfon       !: now climatology of the north boundary barotropic transport
231
232   REAL(wp), DIMENSION(jpi,jpk) ::   &    !:
233      ufon, vfon,       & !: now climatology of the north boundary velocities
234      tfon, sfon,       & !: now climatology of the north boundary temperature and salinity
235      vclin               !: baroclinic componant of the meridian velocity after the radiation
236      !                   !  in yhe obcdyn.F90 routine
237
238   REAL(wp), DIMENSION(jpind:jpinf,jpj) ::   &  !:
239      sshfon_b            !: north boundary ssh correction averaged over the barotropic loop
240                          !: (if Flather's algoritm applied at open boundary)
241
242   !!--------------------------------
243   !! Arrays for radiative North OBC
244   !!--------------------------------
245   !!   
246   REAL(wp), DIMENSION(jpi,jpk,3,3) ::   &   !:
247      unbnd, vnbnd                  !: baroclinic u & v components of the velocity over 3
248      !                             !  rows and 3 time step (now, before, and before before)
249
250   REAL(wp), DIMENSION(jpi,jpk,2,2) ::   &   !:
251      tnbnd, snbnd                  !: north boundary temperature and salinity over
252      !                             !  2 rows and 2 time step (now and before)
253
254   REAL(wp), DIMENSION(jpi,jpk) ::   &     !:
255      u_cynbnd, v_cynbnd            !: Meridional component of the phase speed ratio compu-
256      !                             !  ted with radiation of zonal and meridional velocity
257      !                             !  (respectively) at the north OB (u_cynbnd = cx rdt )
258
259   REAL(wp), DIMENSION(jpi,jpk) ::   &  !:
260      unmsk, vnmsk, tnmsk           !: 2D mask for the North OB
261
262   ! Note that those arrays are optimized for mpp case
263   ! (hence the dimension jpj is the size of one processor subdomain)
264   
265   !!---------------------
266   !! South open boundary
267   !!---------------------
268   INTEGER ::   nis0  , nis1       !: do loop index in mpp case for jpisd, jpisf
269   INTEGER ::   nis0p1, nis1m1     !: do loop index in mpp case for jpisdp1, jpisfm1
270   INTEGER ::   nis1m2, nis0m1     !: do loop index in mpp case for jpisfm1-1,jpisd
271   INTEGER ::   njs0  , njs1       !: do loop index in mpp case for jpsob
272   INTEGER ::   njs0p1, njs1p1     !: do loop index in mpp case for jpsob+1
273
274   REAL(wp), DIMENSION(jpi) ::    &   !:
275      bsfsob              !: now barotropic stream fuction computed at the OBC.The corres-
276      !                   !  ponding bsfn will be computed by the forward time step in dynspg.
277   REAL(wp), DIMENSION(jpi,3,3) ::   &   !:
278      bsbnd               !: south boundary barotropic stream function over
279      !                   !  3 rows and 3 time step (now, before, and before before)
280
281   REAL(wp), DIMENSION(jpisd:jpisf) ::    &   !:
282      bfos,             & !: now climatology of the south boundary barotropic stream function
283      sshfos,           & !: now climatology of the south boundary sea surface height
284      ubtfos,vbtfos       !: now climatology of the south boundary barotropic transport
285
286   REAL(wp), DIMENSION(jpi,jpk) ::    &   !:
287      ufos, vfos,       & !: now climatology of the south boundary velocities
288      tfos, sfos,       & !: now climatology of the south boundary temperature and salinity
289      vclis               !: baroclinic componant of the meridian velocity after the radiation
290      !                   !  in the obcdyn.F90 routine
291
292   REAL(wp), DIMENSION(jpisd:jpisf,jpj) ::   &  !:
293      sshfos_b            !: south boundary ssh correction averaged over the barotropic loop
294                          !: (if Flather's algoritm applied at open boundary)
295
296   !!--------------------------------
297   !! Arrays for radiative South OBC
298   !!--------------------------------
299   !!                        computed by the forward time step in dynspg.
300   REAL(wp), DIMENSION(jpi,jpk,3,3) ::   &   !:
301      usbnd, vsbnd                  !: baroclinic u & v components of the velocity over 3
302      !                             !  rows and 3 time step (now, before, and before before)
303
304   REAL(wp), DIMENSION(jpi,jpk,2,2) ::   &  !:
305      tsbnd, ssbnd                  !: south boundary temperature and salinity over
306      !                             !  2 rows and 2 time step (now and before)
307
308   REAL(wp), DIMENSION(jpi,jpk) ::   &  !:
309      u_cysbnd, v_cysbnd            !: Meridional component of the phase speed ratio compu-
310      !                             !  ted with radiation of zonal and meridional velocity
311      !                             !  (repsectively) at the south OB (u_cynbnd = cx rdt )
312
313   REAL(wp), DIMENSION(jpi,jpk) ::   &  !:
314      usmsk, vsmsk, tsmsk           !: 2D mask for the South OB
315
316   CHARACTER ( len=20 ) :: cffile
317
318#else
319   !!----------------------------------------------------------------------
320   !!   Default option :                                       Empty module
321   !!----------------------------------------------------------------------
322#endif
323
324   !!======================================================================
325END MODULE obc_oce
Note: See TracBrowser for help on using the repository browser.