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.
obcdyn_bt.F90 in trunk/NEMO/OPA_SRC/OBC – NEMO

source: trunk/NEMO/OPA_SRC/OBC/obcdyn_bt.F90 @ 2153

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

Correct USE obcdta in obcdyn_bt

  • Property svn:eol-style set to native
  • Property svn:keywords set to Id
File size: 14.0 KB
Line 
1MODULE obcdyn_bt
2#if ( defined key_dynspg_ts || defined key_dynspg_exp ) && defined key_obc
3   !!=================================================================================
4   !!                       ***  MODULE  obcdyn_bt  ***
5   !! Ocean dynamics:   Radiation/prescription of sea surface heights
6   !!                   on each open boundary
7   !!=================================================================================
8
9   !!---------------------------------------------------------------------------------
10   !!   obc_dyn_bt        : call the subroutine for each open boundary
11   !!   obc_dyn_bt_east   : Flather's algorithm at the east open boundary
12   !!   obc_dyn_bt_west   : Flather's algorithm at the west open boundary
13   !!   obc_dyn_bt_north  : Flather's algorithm at the north open boundary
14   !!   obc_dyn_bt_south  : Flather's algorithm at the south open boundary
15   !!----------------------------------------------------------------------------------
16
17   !!----------------------------------------------------------------------------------
18   !! * Modules used
19   USE oce             ! ocean dynamics and tracers
20   USE dom_oce         ! ocean space and time domain
21   USE phycst          ! physical constants
22   USE obc_oce         ! ocean open boundary conditions
23   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
24   USE lib_mpp         ! distributed memory computing
25   USE obcdta          ! ocean open boundary conditions
26   USE in_out_manager  ! I/O manager
27   USE dynspg_oce      ! surface pressure gradient     (free surface with time-splitting)
28
29   IMPLICIT NONE
30   PRIVATE
31
32   !! * Accessibility
33   PUBLIC obc_dyn_bt  ! routine called in dynnxt (explicit free surface case)
34
35   !!---------------------------------------------------------------------------------
36   !!  OPA 9.0 , LOCEAN-IPSL (2005)
37   !! $Id$
38   !! This software is governed by the CeCILL licence see modipsl/doc/NEMO_CeCILL.txt
39   !!----------------------------------------------------------------------
40
41CONTAINS
42
43   SUBROUTINE obc_dyn_bt ( kt )
44      !!------------------------------------------------------------------------------
45      !!                      SUBROUTINE obc_dyn_bt
46      !!                     ***********************
47      !! ** Purpose :
48      !!      Apply Flather's algorithm at open boundaries for the explicit
49      !!      free surface case and free surface case with time-splitting
50      !!
51      !!      This routine is called in dynnxt.F routine and updates ua, va and sshn.
52      !!
53      !!      The logical variable lp_obc_east, and/or lp_obc_west, and/or lp_obc_north,
54      !!      and/or lp_obc_south allow the user to determine which boundary is an
55      !!      open one (must be done in the param_obc.h90 file).
56      !!
57      !! ** Reference :
58      !!         Flather, R. A., 1976, Mem. Soc. R. Sci. Liege, Ser. 6, 10, 141-164
59      !!
60      !! History :
61      !!   9.0  !  05-12  (V. Garnier) original
62      !!----------------------------------------------------------------------
63      !! * Arguments
64      INTEGER, INTENT( in ) ::   kt
65
66      !!----------------------------------------------------------------------
67
68      IF( lp_obc_east  )   CALL obc_dyn_bt_east 
69      IF( lp_obc_west  )   CALL obc_dyn_bt_west 
70      IF( lp_obc_north )   CALL obc_dyn_bt_north
71      IF( lp_obc_south )   CALL obc_dyn_bt_south
72
73      IF( lk_mpp ) THEN
74         IF( kt >= nit000+3 .AND. ln_rstart ) THEN
75            CALL lbc_lnk( sshb, 'T',  1. )
76            CALL lbc_lnk( ub  , 'U', -1. )
77            CALL lbc_lnk( vb  , 'V', -1. )
78         END IF
79         CALL lbc_lnk( sshn, 'T',  1. )
80         CALL lbc_lnk( ua  , 'U', -1. )
81         CALL lbc_lnk( va  , 'V', -1. )
82      ENDIF
83
84   END SUBROUTINE obc_dyn_bt
85
86# if defined key_dynspg_exp
87   SUBROUTINE obc_dyn_bt_east 
88      !!------------------------------------------------------------------------------
89      !!                  ***  SUBROUTINE obc_dyn_bt_east  ***
90      !!             
91      !! ** Purpose :
92      !!      Apply Flather's algorithm on east OBC velocities ua, va
93      !!      Fix sea surface height (sshn) on east open boundary
94      !!      The logical lfbceast must be .TRUE.
95      !!
96      !!  History :
97      !!   9.0  !  05-12  (V. Garnier) original
98      !!------------------------------------------------------------------------------
99      !! * Local declaration
100      INTEGER ::   ji, jj, jk ! dummy loop indices
101      REAL(wp) ::   z05cx, ztau, zin
102      !!------------------------------------------------------------------------------
103
104      DO ji = nie0, nie1
105         DO jk = 1, jpkm1
106            DO jj = 1, jpj
107               ua(ji,jj,jk) = ua(ji,jj,jk) + sqrt( grav*hur (ji,jj) )               &
108                  &                      * ( ( sshn(ji,jj) + sshn(ji+1,jj) ) * 0.5  &
109                  &                          - sshfoe(jj) ) * uemsk(jj,jk)
110            END DO
111         END DO
112      END DO
113      DO ji = nie0p1, nie1p1
114         DO jj = 1, jpj
115            sshn(ji,jj) = sshn(ji,jj) * (1.-temsk(jj,1)) + temsk(jj,1)*sshfoe(jj)
116         END DO
117      END DO
118
119   END SUBROUTINE obc_dyn_bt_east
120
121
122   SUBROUTINE obc_dyn_bt_west 
123      !!------------------------------------------------------------------------------
124      !!                  ***  SUBROUTINE obc_dyn_bt_west  ***
125      !!                 
126      !! ** Purpose :
127      !!      Apply Flather algorithm on west OBC velocities ua, va
128      !!      Fix sea surface height (sshn) on west open boundary
129      !!      The logical lfbcwest must be .TRUE.
130      !!
131      !!  History :
132      !!   9.0  !  05-12  (V. Garnier) original
133      !!------------------------------------------------------------------------------
134      !! * Local declaration
135      INTEGER ::   ji, jj, jk ! dummy loop indices
136      REAL(wp) ::   z05cx, ztau, zin
137      !!------------------------------------------------------------------------------
138
139      DO ji = niw0, niw1
140         DO jk = 1, jpkm1
141            DO jj = 1, jpj
142               ua(ji,jj,jk) = ua(ji,jj,jk) - sqrt( grav*hur (ji,jj) )               &
143                  &                      * ( ( sshn(ji,jj) + sshn(ji+1,jj) ) * 0.5  &
144                  &                          - sshfow(jj) ) * uwmsk(jj,jk)
145            END DO
146         END DO
147         DO jj = 1, jpj
148            sshn(ji,jj) = sshn(ji,jj) * (1.-twmsk(jj,1)) + twmsk(jj,1)*sshfow(jj)
149         END DO
150      END DO
151
152   END SUBROUTINE obc_dyn_bt_west
153
154   SUBROUTINE obc_dyn_bt_north 
155      !!------------------------------------------------------------------------------
156      !!                ***  SUBROUTINE obc_dyn_bt_north  ***
157      !!
158      !! ** Purpose :
159      !!      Apply Flather algorithm on north OBC velocities ua, va
160      !!      Fix sea surface height (sshn) on north open boundary
161      !!      The logical lfbcnorth must be .TRUE.
162      !!
163      !!  History :
164      !!   9.0  !  05-12  (V. Garnier) original
165      !!------------------------------------------------------------------------------
166      !! * Local declaration
167      INTEGER ::   ji, jj, jk ! dummy loop indices
168      REAL(wp) ::   z05cx, ztau, zin
169      !!------------------------------------------------------------------------------
170
171      DO jj = njn0, njn1
172         DO jk = 1, jpkm1
173            DO ji = 1, jpi
174               va(ji,jj,jk) = va(ji,jj,jk) + sqrt( grav*hvr (ji,jj) )               &
175                  &                      * ( ( sshn(ji,jj) + sshn(ji,jj+1) ) * 0.5  &
176                  &                          - sshfon(ji) ) * vnmsk(ji,jk)
177            END DO
178         END DO
179      END DO
180      DO jj = njn0p1, njn1p1
181         DO ji = 1, jpi
182            sshn(ji,jj)= sshn(ji,jj) * (1.-tnmsk(ji,1)) + sshfon(ji)*tnmsk(ji,1)
183         END DO
184      END DO
185
186   END SUBROUTINE obc_dyn_bt_north
187
188   SUBROUTINE obc_dyn_bt_south 
189      !!------------------------------------------------------------------------------
190      !!                ***  SUBROUTINE obc_dyn_bt_south  ***
191      !!                   
192      !! ** Purpose :
193      !!      Apply Flather algorithm on south OBC velocities ua, va
194      !!      Fix sea surface height (sshn) on south open boundary
195      !!      The logical lfbcsouth must be .TRUE.
196      !!
197      !!  History :
198      !!   9.0  !  05-12  (V. Garnier) original
199      !!------------------------------------------------------------------------------
200      !! * Local declaration
201      INTEGER ::   ji, jj, jk ! dummy loop indices
202      REAL(wp) ::   z05cx, ztau, zin
203
204      !!------------------------------------------------------------------------------
205
206      DO jj = njs0, njs1
207         DO jk = 1, jpkm1
208            DO ji = 1, jpi
209               va(ji,jj,jk) = va(ji,jj,jk) - sqrt( grav*hvr (ji,jj) )               &
210                  &                       * ( ( sshn(ji,jj) + sshn(ji,jj+1) ) * 0.5 &
211                  &                           - sshfos(ji) ) * vsmsk(ji,jk)
212            END DO
213         END DO
214         DO ji = 1, jpi
215            sshn(ji,jj)= sshn(ji,jj) * (1.-tsmsk(ji,1)) + tsmsk(ji,1) * sshfos(ji)
216         END DO
217      END DO
218
219   END SUBROUTINE obc_dyn_bt_south
220
221# elif defined key_dynspg_ts
222
223   SUBROUTINE obc_dyn_bt_east 
224      !!------------------------------------------------------------------------------
225      !!                  ***  SUBROUTINE obc_dyn_bt_east  ***
226      !!
227      !! ** Purpose :
228      !!      Apply Flather's algorithm on east OBC velocities ua, va
229      !!      Fix sea surface height (sshn) on east open boundary
230      !!      The logical lfbceast must be .TRUE.
231      !!
232      !!  History :
233      !!   9.0  !  05-12  (V. Garnier) original
234      !!------------------------------------------------------------------------------
235      !! * Local declaration
236      INTEGER ::   ji, jj, jk ! dummy loop indices
237      REAL(wp) ::   z05cx, ztau, zin
238      !!------------------------------------------------------------------------------
239
240      DO ji = nie0, nie1
241         DO jk = 1, jpkm1
242            DO jj = 1, jpj
243               ua(ji,jj,jk) = ( ua(ji,jj,jk) + sshfoe_b(ji,jj) ) * uemsk(jj,jk)
244            END DO
245         END DO
246      END DO
247      DO ji = nie0p1, nie1p1
248         DO jj = 1, jpj
249            sshn(ji,jj) = sshn(ji,jj) * (1.-temsk(jj,1)) + temsk(jj,1)*sshn_b(ji,jj)
250         END DO
251      END DO
252
253   END SUBROUTINE obc_dyn_bt_east
254
255   SUBROUTINE obc_dyn_bt_west 
256      !!------------------------------------------------------------------------------
257      !!                  ***  SUBROUTINE obc_dyn_bt_west  ***
258      !!
259      !! ** Purpose :
260      !! ** Purpose :
261      !!      Apply Flather algorithm on west OBC velocities ua, va
262      !!      Fix sea surface height (sshn) on west open boundary
263      !!      The logical lfbcwest must be .TRUE.
264      !!
265      !!  History :
266      !!   9.0  !  05-12  (V. Garnier) original
267      !!------------------------------------------------------------------------------
268      !! * Local declaration
269      INTEGER ::   ji, jj, jk ! dummy loop indices
270      REAL(wp) ::   z05cx, ztau, zin
271      !!------------------------------------------------------------------------------
272
273      DO ji = niw0, niw1
274         DO jk = 1, jpkm1
275            DO jj = 1, jpj
276               ua(ji,jj,jk) = ( ua(ji,jj,jk) + sshfow_b(ji,jj) ) * uwmsk(jj,jk)
277            END DO
278         END DO
279         DO jj = 1, jpj
280            sshn(ji,jj) = sshn(ji,jj) * (1.-twmsk(jj,1)) + twmsk(jj,1)*sshn_b(ji,jj)
281         END DO
282      END DO
283
284   END SUBROUTINE obc_dyn_bt_west
285
286   SUBROUTINE obc_dyn_bt_north 
287      !!------------------------------------------------------------------------------
288      !!                     SUBROUTINE obc_dyn_bt_north
289      !!                    *************************
290      !! ** Purpose :
291      !!      Apply Flather algorithm on north OBC velocities ua, va
292      !!      Fix sea surface height (sshn) on north open boundary
293      !!      The logical lfbcnorth must be .TRUE.
294      !!
295      !!  History :
296      !!   9.0  !  05-12  (V. Garnier) original
297      !!------------------------------------------------------------------------------
298      !! * Local declaration
299      INTEGER ::   ji, jj, jk ! dummy loop indices
300      REAL(wp) ::   z05cx, ztau, zin
301      !!------------------------------------------------------------------------------
302
303      DO jj = njn0, njn1
304         DO jk = 1, jpkm1
305            DO ji = 1, jpi
306               va(ji,jj,jk) = ( va(ji,jj,jk) + sshfon_b(ji,jj) ) * vnmsk(jj,jk)
307            END DO
308         END DO
309      END DO
310      DO jj = njn0p1, njn1p1
311         DO ji = 1, jpi
312            sshn(ji,jj)= sshn(ji,jj) * (1.-tnmsk(ji,1)) + sshn_b(ji,jj)*tnmsk(ji,1)
313         END DO
314      END DO
315
316   END SUBROUTINE obc_dyn_bt_north
317
318   SUBROUTINE obc_dyn_bt_south 
319      !!------------------------------------------------------------------------------
320      !!                     SUBROUTINE obc_dyn_bt_south
321      !!                    *************************
322      !! ** Purpose :
323      !!      Apply Flather algorithm on south OBC velocities ua, va
324      !!      Fix sea surface height (sshn) on south open boundary
325      !!      The logical lfbcsouth must be .TRUE.
326      !!
327      !!  History :
328      !!   9.0  !  05-12  (V. Garnier) original
329      !!------------------------------------------------------------------------------
330      !! * Local declaration
331      INTEGER ::   ji, jj, jk ! dummy loop indices
332      REAL(wp) ::   z05cx, ztau, zin
333
334      !!------------------------------------------------------------------------------
335
336      DO jj = njs0, njs1
337         DO jk = 1, jpkm1
338            DO ji = 1, jpi
339               va(ji,jj,jk) = ( va(ji,jj,jk) + sshfos_b(ji,jj) ) * vsmsk(jj,jk)
340            END DO
341         END DO
342         DO ji = 1, jpi
343            sshn(ji,jj)= sshn(ji,jj) * (1.-tsmsk(ji,1)) + tsmsk(ji,1) * sshn_b(ji,jj)
344         END DO
345      END DO
346
347   END SUBROUTINE obc_dyn_bt_south
348
349# endif
350#else
351   !!=================================================================================
352   !!                       ***  MODULE  obcdyn_bt  ***
353   !! Ocean dynamics:   Radiation of velocities on each open boundary
354   !!=================================================================================
355CONTAINS
356
357   SUBROUTINE obc_dyn_bt
358                              ! No open boundaries ==> empty routine
359   END SUBROUTINE obc_dyn_bt
360#endif
361
362END MODULE obcdyn_bt
Note: See TracBrowser for help on using the repository browser.