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

Last change on this file since 1596 was 1528, checked in by rblod, 15 years ago

Suppress rigid-lid option, see ticket #486

  • Property svn:eol-style set to native
  • Property svn:keywords set to Id
File size: 15.7 KB
RevLine 
[3]1MODULE obc_oce
2   !!==============================================================================
3   !!                       ***  MODULE obc_oce   ***
4   !! Open Boundary Cond. :   define related variables
5   !!==============================================================================
[465]6   !!  OPA 9.0 , LOCEAN-IPSL (2005)
[1152]7   !! $Id$
[465]8   !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt)
9   !!----------------------------------------------------------------------
[3]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
[353]16   !!        06/04   (F. Durand) ORCA2_ZIND config
17   !!                 
[3]18   !!----------------------------------------------------------------------
19   !! * Modules used
20   USE par_oce         ! ocean parameters
21   USE obc_par         ! open boundary condition parameters
22
[1158]23#if defined key_obc
24
[3]25   IMPLICIT NONE
26   PUBLIC
27
28   !!----------------------------------------------------------------------
29   !! open boundary variables
30   !!----------------------------------------------------------------------
31   !!
32   !!General variables for open boundaries:
33   !!--------------------------------------
[32]34   INTEGER ::              & !: * namelist ??? *
[367]35      nbobc    = 2  ,      & !: number of open boundaries ( 1=< nbobc =< 4 )
36      nobc_dta = 0           !:  = 0 use the initial state as obc data
[32]37      !                      !   = 1 read obc data in obcxxx.dta files
[3]38
[353]39   LOGICAL ::  ln_obc_clim = .true.  !:  obc data files are climatological
[367]40   LOGICAL ::  ln_obc_fla  = .false. !:  Flather open boundary condition not used
41   LOGICAL ::  ln_vol_cst  = .true.  !:  Conservation of the whole volume
[353]42
[32]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.
[3]55
[32]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).
[3]61
[465]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
[32]67      !                      !  it used to calculate the cumulate flux E-P in the obcvol.F90 routine
[465]68      obcumask, obcvmask     !: u-, v- Force filtering mask for the open
69      !                      !  boundary condition on grad D
70
[32]71   !!--------------------
[3]72   !! East open boundary:
73   !!--------------------
[32]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
[3]80
[32]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.
[3]84
[32]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)
[3]88
[32]89   REAL(wp), DIMENSION(jpjed:jpjef) ::   &  !:
[367]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
[3]93     
[32]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
[3]99
[367]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
[32]104   !!-------------------------------
[3]105   !! Arrays for radiative East OBC:
106   !!-------------------------------
[32]107   REAL(wp), DIMENSION(jpj,jpk,3,3) ::   &  !:
108      uebnd, vebnd                  !: baroclinic u & v component of the velocity over 3 rows
[3]109                                    ! and 3 time step (now, before, and before before)
110
[32]111   REAL(wp), DIMENSION(jpj,jpk,2,2) ::   &  !:
112      tebnd, sebnd                  !: East boundary temperature and salinity over 2 rows
[3]113                                    ! and 2 time step (now and before)
114
[32]115   REAL(wp), DIMENSION(jpj,jpk) ::   &  !:
116      u_cxebnd, v_cxebnd            !: Zonal component of the phase speed ratio computed with
[3]117                                    ! radiation of u and v velocity (respectively) at the
118                                    ! east open boundary (u_cxebnd = cx rdt )
119
[32]120   REAL(wp), DIMENSION(jpj,jpk) ::   &  !:
121      uemsk, vemsk, temsk           !: 2D mask for the East OB
[3]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   !!--------------------
[32]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
[3]134
[32]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.
[3]138
[32]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)
[3]142
[32]143   REAL(wp), DIMENSION(jpjwd:jpjwf) ::   &  !:
[367]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
[3]147
[32]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
[3]153
[367]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
[3]158   !!-------------------------------
[32]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)
[3]164
[32]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)
[3]168
[32]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 )
[3]173
[32]174   REAL(wp), DIMENSION(jpj,jpk) ::    &  !:
175      uwmsk, vwmsk, twmsk           !: 2D mask for the West OB
[3]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   !!---------------------
[32]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
[3]189
[32]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.
[3]193
[32]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)
[3]197
[32]198   REAL(wp), DIMENSION(jpind:jpinf) ::   &  !:
[367]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
[3]202
[32]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
[3]208
[367]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
[3]213   !!--------------------------------
[32]214   !! Arrays for radiative North OBC
215   !!--------------------------------
[3]216   !!   
[32]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)
[3]220
[32]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)
[3]224
[32]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 )
[3]229
[32]230   REAL(wp), DIMENSION(jpi,jpk) ::   &  !:
231      unmsk, vnmsk, tnmsk           !: 2D mask for the North OB
[3]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   !!---------------------
[32]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
[3]244
[32]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)
[3]251
[32]252   REAL(wp), DIMENSION(jpisd:jpisf) ::    &   !:
[367]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
[3]256
[32]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
[3]262
[367]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
[3]267   !!--------------------------------
[32]268   !! Arrays for radiative South OBC
269   !!--------------------------------
[3]270   !!                        computed by the forward time step in dynspg.
[32]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)
[3]274
[32]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)
[3]278
[32]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 )
[3]283
[32]284   REAL(wp), DIMENSION(jpi,jpk) ::   &  !:
285      usmsk, vsmsk, tsmsk           !: 2D mask for the South OB
[3]286
[1151]287   CHARACTER ( len=20 ) :: cffile
[3]288
289#else
290   !!----------------------------------------------------------------------
291   !!   Default option :                                       Empty module
292   !!----------------------------------------------------------------------
293#endif
[32]294
[3]295   !!======================================================================
296END MODULE obc_oce
Note: See TracBrowser for help on using the repository browser.