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 @ 418

Last change on this file since 418 was 367, checked in by opalod, 18 years ago

nemo_v1_update_035 : CT : OBCs adapted to the new surface pressure gradient algorithms

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