source: trunk/NEMOGCM/NEMO/OPA_SRC/OBC/obc_oce.F90 @ 23

Last change on this file since 23 was 23, checked in by rblod, 12 years ago

Add possibility to choose where to apply the barotropic correction

File size: 19.4 KB
Line 
1MODULE obc_oce
2   !!==============================================================================
3   !!                       ***  MODULE obc_oce   ***
4   !! Open Boundary Cond. :   define related variables
5   !!==============================================================================
6   !! history :  OPA  ! 1991-01 (CLIPPER)  Original code
7   !!   NEMO     1.0  ! 2002-02 (C. Talandier)  modules, F90
8   !!----------------------------------------------------------------------
9#if defined key_obc
10   !!----------------------------------------------------------------------
11   !!   'key_obc' :                                 Open Boundary Condition
12   !!----------------------------------------------------------------------
13   USE par_oce         ! ocean parameters
14   USE obc_par         ! open boundary condition parameters
15
16   IMPLICIT NONE
17   PUBLIC
18   
19   PUBLIC   obc_oce_alloc   ! called by obcini.F90 module
20
21   !!----------------------------------------------------------------------
22   !! open boundary variables
23   !!----------------------------------------------------------------------
24   !
25   !                                            !!* Namelist namobc: open boundary condition *
26   INTEGER           ::   nn_obcdta   = 0        !:  = 0 use the initial state as obc data
27   !                                             !   = 1 read obc data in obcxxx.dta files
28   CHARACTER(len=20) ::   cn_obcdta   = 'annual' !: set to annual  if obc datafile hold 1 year of data
29   !                                             !  set to monthly if obc datafile hold 1 month of data
30   LOGICAL           ::   ln_obc_clim = .true.   !:  obc data files are climatological
31   LOGICAL           ::   ln_obc_fla  = .false.  !:  Flather open boundary condition not used
32   LOGICAL           ::   ln_vol_cst  = .true.   !:  Conservation of the whole volume
33   REAL(wp)          ::   rn_dpein    =  1.      !: damping time scale for inflow at East open boundary
34   REAL(wp)          ::   rn_dpwin    =  1.      !:    "                      "   at West open boundary
35   REAL(wp)          ::   rn_dpsin    =  1.      !:    "                      "   at South open boundary
36   REAL(wp)          ::   rn_dpnin    =  1.      !:    "                      "   at North open boundary
37   REAL(wp)          ::   rn_dpeob    = 15.      !: damping time scale for the climatology at East open boundary
38   REAL(wp)          ::   rn_dpwob    = 15.      !:    "                           "       at West open boundary
39   REAL(wp)          ::   rn_dpsob    = 15.      !:    "                           "       at South open boundary
40   REAL(wp)          ::   rn_dpnob    = 15.      !:    "                           "       at North open boundary
41   REAL(wp)          ::   rn_volemp   =  1.      !: = 0 the total volume will have the variability of the
42   !                                             !      surface Flux E-P else (volemp = 1) the volume will be constant
43   !                                             !  = 1 the volume will be constant during all the integration.
44
45   !                                  !!! OLD non-DOCTOR name of namelist variables
46   INTEGER  ::   nbobc                 !: number of open boundaries ( 1=< nbobc =< 4 )
47   INTEGER  ::   nobc_dta              !:  = 0 use the initial state as obc data
48   REAL(wp) ::   rdpein                !: damping time scale for inflow at East open boundary
49   REAL(wp) ::   rdpwin                !:    "                      "   at West open boundary
50   REAL(wp) ::   rdpsin                !:    "                      "   at South open boundary
51   REAL(wp) ::   rdpnin                !:    "                      "   at North open boundary
52   REAL(wp) ::   rdpeob                !: damping time scale for the climatology at East open boundary
53   REAL(wp) ::   rdpwob                !:    "                           "       at West open boundary
54   REAL(wp) ::   rdpsob                !:    "                           "       at South open boundary
55   REAL(wp) ::   rdpnob                !:    "                           "       at North open boundary
56   REAL(wp) ::   volemp                !: = 0 the total volume will have the variability of the
57   CHARACTER(len=20) :: cffile
58
59
60   !!General variables for open boundaries:
61   !!--------------------------------------
62   LOGICAL ::   lfbceast, lfbcwest      !: logical flag for a fixed East and West open boundaries       
63   LOGICAL ::   lfbcnorth, lfbcsouth    !: logical flag for a fixed North and South open boundaries       
64   !                                    !  These logical flags are set to 'true' if damping time
65   !                                    !  scale are set to 0 in the namelist, for both inflow and outflow).
66
67!--> IND4 Begin
68!     Logical for restarting with radiative OBCs, but without an OBC restart restart.obc.output file.
69!     During the first 30 time steps, used FIXED boundary conditions.
70!     We modify : obcini, obctra, obcdyn
71   LOGICAL :: ln_obc_rstart = .TRUE. !: radiative OBCs, but do not read restart.obc.output
72!<-- IND4 End
73
74   REAL(wp), PUBLIC ::   obcsurftot       !: Total lateral surface of open boundaries
75!     Add computation of E/W/N/S lateral surface of open boundaries
76   REAL(wp), PUBLIC ::    &  !:
77      obcsurfeast     ,   &  !: East  lateral surface of open boundaries
78      obcsurfwest     ,   &  !: West  lateral surface of open boundaries
79      obcsurfnorth    ,   &  !: North lateral surface of open boundaries
80      obcsurfsouth           !: South lateral surface of open boundaries
81   
82   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: &  !:
83      obctmsk,            &  !: mask array identical to tmask, execpt along OBC where it is set to 0
84      !                      !  it used to calculate the cumulate flux E-P in the obcvol.F90 routine
85      obcumask, obcvmask     !: u-, v- Force filtering mask for the open
86      !                      !  boundary condition on grad D
87
88   !!--------------------
89   !! East open boundary:
90   !!--------------------
91   INTEGER ::   nie0  , nie1      !: do loop index in mpp case for jpieob
92   INTEGER ::   nie0p1, nie1p1    !: do loop index in mpp case for jpieob+1
93   INTEGER ::   nie0m1, nie1m1    !: do loop index in mpp case for jpieob-1
94   INTEGER ::   nje0  , nje1      !: do loop index in mpp case for jpjed, jpjef
95   INTEGER ::   nje0p1, nje1m1    !: do loop index in mpp case for jpjedp1,jpjefm1
96   INTEGER ::   nje1m2, nje0m1    !: do loop index in mpp case for jpjefm1-1,jpjed
97
98   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:) ::   &  !:
99      sshfoe,           & !: now climatology of the east boundary sea surface height
100      ubtfoe,vbtfoe       !: now climatology of the east boundary barotropic transport
101     
102   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   &  !:
103      ufoe, vfoe,       & !: now climatology of the east boundary velocities
104      tfoe, sfoe,       & !: now climatology of the east boundary temperature and salinity
105      uclie               !: baroclinic componant of the zonal velocity after radiation
106      !                   ! in the obcdyn.F90 routine
107
108   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   sshfoe_b   !: east boundary ssh correction averaged over the barotropic loop
109      !                                            !  (if Flather's algoritm applied at open boundary)
110
111   !!-------------------------------
112   !! Arrays for radiative East OBC:
113   !!-------------------------------
114   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::   uebnd, vebnd      !: baroclinic u & v component of the velocity over 3 rows
115      !                                                    !  and 3 time step (now, before, and before before)
116   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::   tebnd, sebnd      !: East boundary temperature and salinity over 2 rows
117      !                                                    !  and 2 time step (now and before)
118   REAL(wp), ALLOCATABLE, SAVE,     DIMENSION(:,:) ::   u_cxebnd, v_cxebnd    !: Zonal component of the phase speed ratio computed with
119      !                                                    !  radiation of u and v velocity (respectively) at the
120      !                                                    !  east open boundary (u_cxebnd = cx rdt )
121   REAL(wp), ALLOCATABLE, SAVE,     DIMENSION(:,:) ::   uemsk, vemsk, temsk   !: 2D mask for the East OB
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   !!--------------------
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
134
135   REAL(wp), ALLOCATABLE, SAVE,   DIMENSION(:) ::   &  !:
136      sshfow,           & !: now climatology of the west boundary sea surface height
137      ubtfow,vbtfow       !: now climatology of the west boundary barotropic transport
138
139   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   &  !:
140      ufow, vfow,       & !: now climatology of the west velocities
141      tfow, sfow,       & !: now climatology of the west temperature and salinity
142      ucliw               !: baroclinic componant of the zonal velocity after the radiation
143      !                   !  in the obcdyn.F90 routine
144
145   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   sshfow_b    !: west boundary ssh correction averaged over the barotropic loop
146      !                                          !  (if Flather's algoritm applied at open boundary)
147
148   !!-------------------------------
149   !! Arrays for radiative West OBC
150   !!-------------------------------
151   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::   uwbnd, vwbnd     !: baroclinic u & v components of the velocity over 3 rows
152      !                                                   !  and 3 time step (now, before, and before before)
153   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::   twbnd, swbnd     !: west boundary temperature and salinity over 2 rows and
154      !                                                   !  2 time step (now and before)
155   REAL(wp), ALLOCATABLE, SAVE,     DIMENSION(:,:) ::   u_cxwbnd, v_cxwbnd   !: Zonal component of the phase speed ratio computed with
156      !                                                   !  radiation of zonal and meridional velocity (respectively)
157      !                                                   !  at the west open boundary (u_cxwbnd = cx rdt )
158   REAL(wp), ALLOCATABLE, SAVE,     DIMENSION(:,:) ::   uwmsk, vwmsk, twmsk  !: 2D mask for the West 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   !! North open boundary
165   !!---------------------
166   INTEGER ::   nin0  , nin1       !: do loop index in mpp case for jpind, jpinf
167   INTEGER ::   nin0p1, nin1m1     !: do loop index in mpp case for jpindp1, jpinfm1
168   INTEGER ::   nin1m2, nin0m1     !: do loop index in mpp case for jpinfm1-1,jpind
169   INTEGER ::   njn0  , njn1       !: do loop index in mpp case for jpnob
170   INTEGER ::   njn0p1, njn1p1     !: do loop index in mpp case for jpnob+1
171   INTEGER ::   njn0m1, njn1m1     !: do loop index in mpp case for jpnob-1
172
173   REAL(wp), ALLOCATABLE, SAVE,   DIMENSION(:) ::   &  !:
174      sshfon,           & !: now climatology of the north boundary sea surface height
175      ubtfon,vbtfon       !: now climatology of the north boundary barotropic transport
176
177   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   &    !:
178      ufon, vfon,       & !: now climatology of the north boundary velocities
179      tfon, sfon,       & !: now climatology of the north boundary temperature and salinity
180      vclin               !: baroclinic componant of the meridian velocity after the radiation
181      !                   !  in yhe obcdyn.F90 routine
182
183   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   sshfon_b      !: north boundary ssh correction averaged over the barotropic loop
184      !                                            !  (if Flather's algoritm applied at open boundary)
185
186   !!--------------------------------
187   !! Arrays for radiative North OBC
188   !!--------------------------------
189   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::   unbnd, vnbnd      !: baroclinic u & v components of the velocity over 3
190      !                                                    !  rows and 3 time step (now, before, and before before)
191   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::   tnbnd, snbnd      !: north boundary temperature and salinity over
192      !                                                    !  2 rows and 2 time step (now and before)
193   REAL(wp), ALLOCATABLE, SAVE,     DIMENSION(:,:) ::   u_cynbnd, v_cynbnd    !: Meridional component of the phase speed ratio compu-
194      !                                                    !  ted with radiation of zonal and meridional velocity
195      !                                                    !  (respectively) at the north OB (u_cynbnd = cx rdt )
196   REAL(wp), ALLOCATABLE, SAVE,     DIMENSION(:,:) ::   unmsk, vnmsk, tnmsk   !: 2D mask for the North OB
197
198   ! Note that those arrays are optimized for mpp case
199   ! (hence the dimension jpj is the size of one processor subdomain)
200   
201   !!---------------------
202   !! South open boundary
203   !!---------------------
204   INTEGER ::   nis0  , nis1       !: do loop index in mpp case for jpisd, jpisf
205   INTEGER ::   nis0p1, nis1m1     !: do loop index in mpp case for jpisdp1, jpisfm1
206   INTEGER ::   nis1m2, nis0m1     !: do loop index in mpp case for jpisfm1-1,jpisd
207   INTEGER ::   njs0  , njs1       !: do loop index in mpp case for jpsob
208   INTEGER ::   njs0p1, njs1p1     !: do loop index in mpp case for jpsob+1
209
210   REAL(wp), ALLOCATABLE, SAVE,   DIMENSION(:) ::    &   !:
211      sshfos,           & !: now climatology of the south boundary sea surface height
212      ubtfos,vbtfos       !: now climatology of the south boundary barotropic transport
213
214   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::    &   !:
215      ufos, vfos,       & !: now climatology of the south boundary velocities
216      tfos, sfos,       & !: now climatology of the south boundary temperature and salinity
217      vclis               !: baroclinic componant of the meridian velocity after the radiation
218      !                   !  in the obcdyn.F90 routine
219
220   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) ::   sshfos_b     !: south boundary ssh correction averaged over the barotropic loop
221      !                                           !  (if Flather's algoritm applied at open boundary)
222
223   !!--------------------------------
224   !! Arrays for radiative South OBC   (computed by the forward time step in dynspg)
225   !!--------------------------------
226   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::   usbnd, vsbnd     !: baroclinic u & v components of the velocity over 3
227      !                                                   !  rows and 3 time step (now, before, and before before)
228   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) ::   tsbnd, ssbnd     !: south boundary temperature and salinity over
229      !                                                   !  2 rows and 2 time step (now and before)
230   REAL(wp), ALLOCATABLE, SAVE,     DIMENSION(:,:) ::   u_cysbnd, v_cysbnd   !: Meridional component of the phase speed ratio
231      !                                                   !  computed with radiation of zonal and meridional velocity
232      !                                                   !  (repsectively) at the south OB (u_cynbnd = cx rdt )
233   REAL(wp), ALLOCATABLE, SAVE,     DIMENSION(:,:) ::   usmsk, vsmsk, tsmsk  !: 2D mask for the South OB
234
235   !!----------------------------------------------------------------------
236   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
237   !! $Id: obc_oce.F90 2715 2011-03-30 15:58:35Z rblod $
238   !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
239   !!----------------------------------------------------------------------
240CONTAINS
241
242   INTEGER FUNCTION obc_oce_alloc()
243      !!----------------------------------------------------------------------
244      !!               ***  FUNCTION obc_oce_alloc  ***
245      !!----------------------------------------------------------------------
246
247      ALLOCATE(                                                               &
248              !! East open boundary
249              obctmsk(jpi,jpj), obcumask(jpi,jpj), obcvmask(jpi,jpj),        &
250              sshfoe(jpjed:jpjef), ubtfoe(jpjed:jpjef), vbtfoe(jpjed:jpjef), &
251              ufoe(jpj,jpk), vfoe(jpj,jpk), tfoe(jpj,jpk), sfoe(jpj,jpk),    &
252              uclie(jpj,jpk), sshfoe_b(jpjed:jpjef,jpj),                     &
253              !! Arrays for radiative East OBC
254              uebnd(jpj,jpk,3,3), vebnd(jpj,jpk,3,3) ,                       &
255              tebnd(jpj,jpk,2,2), sebnd(jpj,jpk,2,2),                        &
256              u_cxebnd(jpj,jpk), v_cxebnd(jpj,jpk),                          &
257              uemsk(jpj,jpk), vemsk(jpj,jpk), temsk(jpj,jpk),                &
258              !! West open boundary
259              sshfow(jpjwd:jpjwf), ubtfow(jpjwd:jpjwf), vbtfow(jpjwd:jpjwf), &
260              ufow(jpj,jpk), vfow(jpj,jpk), tfow(jpj,jpk),                   &
261              sfow(jpj,jpk), ucliw(jpj,jpk), sshfow_b(jpjwd:jpjwf,jpj),      &
262              !! Arrays for radiative West OBC
263              uwbnd(jpj,jpk,3,3), vwbnd(jpj,jpk,3,3),                        &
264              twbnd(jpj,jpk,2,2), swbnd(jpj,jpk,2,2),                        &
265              u_cxwbnd(jpj,jpk), v_cxwbnd(jpj,jpk),                          &
266              uwmsk(jpj,jpk), vwmsk(jpj,jpk), twmsk(jpj,jpk),                &
267              !! North open boundary
268              sshfon(jpind:jpinf), ubtfon(jpind:jpinf), vbtfon(jpind:jpinf), &
269              ufon(jpi,jpk), vfon(jpi,jpk), tfon(jpi,jpk),                   &
270              sfon(jpi,jpk), vclin(jpi,jpk), sshfon_b(jpind:jpinf,jpj),      &
271              !! Arrays for radiative North OBC
272              unbnd(jpi,jpk,3,3), vnbnd(jpi,jpk,3,3),                        &
273              tnbnd(jpi,jpk,2,2), snbnd(jpi,jpk,2,2),                        &
274              u_cynbnd(jpi,jpk), v_cynbnd(jpi,jpk),                          &
275              unmsk(jpi,jpk), vnmsk(jpi,jpk), tnmsk (jpi,jpk),               &
276              !! South open boundary
277              sshfos(jpisd:jpisf), ubtfos(jpisd:jpisf), vbtfos(jpisd:jpisf), &
278              ufos(jpi,jpk), vfos(jpi,jpk), tfos(jpi,jpk),                   &
279              sfos(jpi,jpk), vclis(jpi,jpk),                                 &
280              sshfos_b(jpisd:jpisf,jpj),                                     &
281              !! Arrays for radiative South OBC
282              usbnd(jpi,jpk,3,3), vsbnd(jpi,jpk,3,3),                        &
283              tsbnd(jpi,jpk,2,2), ssbnd(jpi,jpk,2,2),                        &
284              u_cysbnd(jpi,jpk), v_cysbnd(jpi,jpk),                          &
285              usmsk(jpi,jpk), vsmsk(jpi,jpk), tsmsk(jpi,jpk),                &
286              !!
287              STAT=obc_oce_alloc )
288      !
289   END FUNCTION obc_oce_alloc
290   
291#else
292   !!----------------------------------------------------------------------
293   !!   Default option         Empty module                          No OBC
294   !!----------------------------------------------------------------------
295#endif
296
297   !!======================================================================
298END MODULE obc_oce
Note: See TracBrowser for help on using the repository browser.