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

Last change on this file since 3 was 3, checked in by opalod, 20 years ago

Initial revision

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