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 branches/TAM_V3_0/NEMO/OPA_SRC/OBC – NEMO

source: branches/TAM_V3_0/NEMO/OPA_SRC/OBC/obc_oce.F90 @ 1884

Last change on this file since 1884 was 1884, checked in by rblod, 14 years ago

Light adaptation of NEMO direct model routine to handle TAM

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