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 branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/OBC – NEMO

source: branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/OBC/obcdyn_bt.F90 @ 4400

Last change on this file since 4400 was 3211, checked in by spickles2, 12 years ago

Stephen Pickles, 11 Dec 2011

Commit to bring the rest of the DCSE NEMO development branch
in line with the latest development version. This includes
array index re-ordering of all OPA_SRC/.

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