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

Last change on this file since 1601 was 1601, checked in by ctlod, 15 years ago

Doctor naming of OPA namelist variables , see ticket: #526

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