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/2012/dev_CMCC_2012/NEMOGCM/NEMO/OPA_SRC/OBC – NEMO

source: branches/2012/dev_CMCC_2012/NEMOGCM/NEMO/OPA_SRC/OBC/obcdyn_bt.F90 @ 3593

Last change on this file since 3593 was 3593, checked in by vichi, 11 years ago

Add in branch 2012/dev_CMCC_2012 changes from dev_r3365_CMCC1_BDYOBCopt & dev_r3379_CMCC6_topbfm, see ticket 1002

File size: 13.1 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   !!            3.5  ! 2012     (S. Mocavero, I. Epicoco) Updates for the
8   !!                             optimization of OBC communications
9   !!----------------------------------------------------------------------
10#if ( defined key_dynspg_ts || defined key_dynspg_exp ) && defined key_obc
11   !!----------------------------------------------------------------------
12   !!   'key_dynspg_ts'     OR                   time spliting free surface
13   !!   'key_dynspg_exp'    AND                       explicit free surface
14   !!   'key_obc'                                   Open Boundary Condition
15   !!----------------------------------------------------------------------
16   !!   obc_dyn_bt        : call the subroutine for each open boundary
17   !!   obc_dyn_bt_east   : Flather's algorithm at the east open boundary
18   !!   obc_dyn_bt_west   : Flather's algorithm at the west open boundary
19   !!   obc_dyn_bt_north  : Flather's algorithm at the north open boundary
20   !!   obc_dyn_bt_south  : Flather's algorithm at the south open boundary
21   !!----------------------------------------------------------------------
22   USE oce             ! ocean dynamics and tracers
23   USE dom_oce         ! ocean space and time domain
24   USE phycst          ! physical constants
25   USE obc_oce         ! ocean open boundary conditions
26   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
27   USE lib_mpp         ! distributed memory computing
28   USE obcdta          ! ocean open boundary conditions
29   USE in_out_manager  ! I/O manager
30   USE dynspg_oce      ! surface pressure gradient     (free surface with time-splitting)
31
32   IMPLICIT NONE
33   PRIVATE
34
35   PUBLIC   obc_dyn_bt  ! routine called in dynnxt (explicit free surface case)
36
37   !!----------------------------------------------------------------------
38   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
39   !! $Id: obcdyn_bt.F90 2715 2011-03-30 15:58:35Z rblod $
40   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
41   !!----------------------------------------------------------------------
42CONTAINS
43
44   SUBROUTINE obc_dyn_bt( kt )
45      !!----------------------------------------------------------------------
46      !!                 ***  SUBROUTINE obc_dyn_bt  ***
47      !!
48      !! ** Purpose :   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 :   Flather, R. A., 1976, Mem. Soc. R. Sci. Liege, Ser. 6, 10, 141-164
58      !!----------------------------------------------------------------------
59      INTEGER, INTENT(in) ::   kt
60      !!----------------------------------------------------------------------
61
62      IF( lp_obc_east  )   CALL obc_dyn_bt_east 
63      IF( lp_obc_west  )   CALL obc_dyn_bt_west 
64      IF( lp_obc_north )   CALL obc_dyn_bt_north
65      IF( lp_obc_south )   CALL obc_dyn_bt_south
66
67      IF( lk_mpp ) THEN
68         IF( kt >= nit000+3 .AND. ln_rstart ) THEN
69            CALL lbc_obc_lnk( sshb, 'T',  1. )
70            CALL lbc_obc_lnk( ub  , 'U', -1. )
71            CALL lbc_obc_lnk( vb  , 'V', -1. )
72         END IF
73         CALL lbc_obc_lnk( sshn, 'T',  1. )
74         CALL lbc_obc_lnk( ua  , 'U', -1. )
75         CALL lbc_obc_lnk( va  , 'V', -1. )
76      ENDIF
77
78   END SUBROUTINE obc_dyn_bt
79
80# if defined key_dynspg_exp
81
82   SUBROUTINE obc_dyn_bt_east 
83      !!----------------------------------------------------------------------
84      !!                  ***  SUBROUTINE obc_dyn_bt_east  ***
85      !!             
86      !! ** Purpose :
87      !!      Apply Flather's algorithm on east OBC velocities ua, va
88      !!      Fix sea surface height (sshn) on east open boundary
89      !!      The logical lfbceast must be .TRUE.
90      !!----------------------------------------------------------------------
91      INTEGER, INTENT(in) ::   kt
92      !!----------------------------------------------------------------------
93      INTEGER ::   ji, jj, jk   ! dummy loop indices
94      !!----------------------------------------------------------------------
95
96      DO ji = nie0, nie1
97         DO jk = 1, jpkm1
98            DO jj = 1, jpj
99               ua(ji,jj,jk) = ua(ji,jj,jk) + sqrt( grav*hur (ji,jj) )               &
100                  &                      * ( ( sshn(ji,jj) + sshn(ji+1,jj) ) * 0.5  &
101                  &                          - sshfoe(jj) ) * uemsk(jj,jk)
102            END DO
103         END DO
104      END DO
105      DO ji = nie0p1, nie1p1
106         DO jj = 1, jpj
107            sshn(ji,jj) = sshn(ji,jj) * (1.-temsk(jj,1)) + temsk(jj,1)*sshfoe(jj)
108         END DO
109      END DO
110
111   END SUBROUTINE obc_dyn_bt_east
112
113
114   SUBROUTINE obc_dyn_bt_west 
115      !!----------------------------------------------------------------------
116      !!                  ***  SUBROUTINE obc_dyn_bt_west  ***
117      !!                 
118      !! ** Purpose :
119      !!      Apply Flather algorithm on west OBC velocities ua, va
120      !!      Fix sea surface height (sshn) on west open boundary
121      !!      The logical lfbcwest must be .TRUE.
122      !!----------------------------------------------------------------------
123      INTEGER ::   ji, jj, jk   ! dummy loop indices
124      !!----------------------------------------------------------------------
125      !
126      DO ji = niw0, niw1
127         DO jk = 1, jpkm1
128            DO jj = 1, jpj
129               ua(ji,jj,jk) = ua(ji,jj,jk) - sqrt( grav*hur (ji,jj) )               &
130                  &                      * ( ( sshn(ji,jj) + sshn(ji+1,jj) ) * 0.5  &
131                  &                          - sshfow(jj) ) * uwmsk(jj,jk)
132            END DO
133         END DO
134         DO jj = 1, jpj
135            sshn(ji,jj) = sshn(ji,jj) * (1.-twmsk(jj,1)) + twmsk(jj,1)*sshfow(jj)
136         END DO
137      END DO
138      !
139   END SUBROUTINE obc_dyn_bt_west
140
141
142   SUBROUTINE obc_dyn_bt_north 
143      !!------------------------------------------------------------------------------
144      !!                ***  SUBROUTINE obc_dyn_bt_north  ***
145      !!
146      !! ** Purpose :
147      !!      Apply Flather algorithm on north OBC velocities ua, va
148      !!      Fix sea surface height (sshn) on north open boundary
149      !!      The logical lfbcnorth must be .TRUE.
150      !!----------------------------------------------------------------------
151      INTEGER ::   ji, jj, jk   ! dummy loop indices
152      !!----------------------------------------------------------------------
153      !
154      DO jj = njn0, njn1
155         DO jk = 1, jpkm1
156            DO ji = 1, jpi
157               va(ji,jj,jk) = va(ji,jj,jk) + sqrt( grav*hvr (ji,jj) )               &
158                  &                      * ( ( sshn(ji,jj) + sshn(ji,jj+1) ) * 0.5  &
159                  &                          - sshfon(ji) ) * vnmsk(ji,jk)
160            END DO
161         END DO
162      END DO
163      DO jj = njn0p1, njn1p1
164         DO ji = 1, jpi
165            sshn(ji,jj)= sshn(ji,jj) * (1.-tnmsk(ji,1)) + sshfon(ji)*tnmsk(ji,1)
166         END DO
167      END DO
168      !
169   END SUBROUTINE obc_dyn_bt_north
170
171
172   SUBROUTINE obc_dyn_bt_south 
173      !!----------------------------------------------------------------------
174      !!                ***  SUBROUTINE obc_dyn_bt_south  ***
175      !!                   
176      !! ** Purpose :
177      !!      Apply Flather algorithm on south OBC velocities ua, va
178      !!      Fix sea surface height (sshn) on south open boundary
179      !!      The logical lfbcsouth must be .TRUE.
180      !!----------------------------------------------------------------------
181      INTEGER ::   ji, jj, jk   ! dummy loop indices
182      !!----------------------------------------------------------------------
183      !
184      DO jj = njs0, njs1
185         DO jk = 1, jpkm1
186            DO ji = 1, jpi
187               va(ji,jj,jk) = va(ji,jj,jk) - sqrt( grav*hvr (ji,jj) )               &
188                  &                       * ( ( sshn(ji,jj) + sshn(ji,jj+1) ) * 0.5 &
189                  &                           - sshfos(ji) ) * vsmsk(ji,jk)
190            END DO
191         END DO
192         DO ji = 1, jpi
193            sshn(ji,jj)= sshn(ji,jj) * (1.-tsmsk(ji,1)) + tsmsk(ji,1) * sshfos(ji)
194         END DO
195      END DO
196      !
197   END SUBROUTINE obc_dyn_bt_south
198
199# elif defined key_dynspg_ts
200
201   SUBROUTINE obc_dyn_bt_east 
202      !!------------------------------------------------------------------------------
203      !!                  ***  SUBROUTINE obc_dyn_bt_east  ***
204      !!
205      !! ** Purpose :
206      !!      Apply Flather's algorithm on east OBC velocities ua, va
207      !!      Fix sea surface height (sshn) on east open boundary
208      !!      The logical lfbceast must be .TRUE.
209      !!----------------------------------------------------------------------
210      INTEGER ::   ji, jj, jk   ! dummy loop indices
211      !!----------------------------------------------------------------------
212      !
213      DO ji = nie0, nie1
214         DO jk = 1, jpkm1
215            DO jj = 1, jpj
216               ua(ji,jj,jk) = ( ua(ji,jj,jk) + sshfoe_b(ji,jj) ) * uemsk(jj,jk)
217            END DO
218         END DO
219      END DO
220      DO ji = nie0p1, nie1p1
221         DO jj = 1, jpj
222            sshn(ji,jj) = sshn(ji,jj) * (1.-temsk(jj,1)) + temsk(jj,1)*sshn_b(ji,jj)
223         END DO
224      END DO
225      !
226   END SUBROUTINE obc_dyn_bt_east
227
228
229   SUBROUTINE obc_dyn_bt_west 
230      !!---------------------------------------------------------------------
231      !!                  ***  SUBROUTINE obc_dyn_bt_west  ***
232      !!
233      !! ** Purpose :   Apply Flather algorithm on west OBC velocities ua, va
234      !!      Fix sea surface height (sshn) on west open boundary
235      !!      The logical lfbcwest must be .TRUE.
236      !!----------------------------------------------------------------------
237      INTEGER ::   ji, jj, jk   ! dummy loop indices
238      !!----------------------------------------------------------------------
239      !
240      DO ji = niw0, niw1
241         DO jk = 1, jpkm1
242            DO jj = 1, jpj
243               ua(ji,jj,jk) = ( ua(ji,jj,jk) + sshfow_b(ji,jj) ) * uwmsk(jj,jk)
244            END DO
245         END DO
246         DO jj = 1, jpj
247            sshn(ji,jj) = sshn(ji,jj) * (1.-twmsk(jj,1)) + twmsk(jj,1)*sshn_b(ji,jj)
248         END DO
249      END DO
250      !
251   END SUBROUTINE obc_dyn_bt_west
252
253
254   SUBROUTINE obc_dyn_bt_north 
255      !!------------------------------------------------------------------------------
256      !!                ***  SUBROUTINE obc_dyn_bt_north  ***
257      !!               
258      !! ** Purpose :
259      !!      Apply Flather algorithm on north OBC velocities ua, va
260      !!      Fix sea surface height (sshn) on north open boundary
261      !!      The logical lfbcnorth must be .TRUE.
262      !!----------------------------------------------------------------------
263      INTEGER ::   ji, jj, jk   ! dummy loop indices
264      !!----------------------------------------------------------------------
265      !
266      DO jj = njn0, njn1
267         DO jk = 1, jpkm1
268            DO ji = 1, jpi
269               va(ji,jj,jk) = ( va(ji,jj,jk) + sshfon_b(ji,jj) ) * vnmsk(jj,jk)
270            END DO
271         END DO
272      END DO
273      DO jj = njn0p1, njn1p1
274         DO ji = 1, jpi
275            sshn(ji,jj)= sshn(ji,jj) * (1.-tnmsk(ji,1)) + sshn_b(ji,jj)*tnmsk(ji,1)
276         END DO
277      END DO
278      !
279   END SUBROUTINE obc_dyn_bt_north
280
281
282   SUBROUTINE obc_dyn_bt_south 
283      !!------------------------------------------------------------------------------
284      !!                ***  SUBROUTINE obc_dyn_bt_south  ***
285      !!                 
286      !! ** Purpose :
287      !!      Apply Flather algorithm on south OBC velocities ua, va
288      !!      Fix sea surface height (sshn) on south open boundary
289      !!      The logical lfbcsouth must be .TRUE.
290      !!----------------------------------------------------------------------
291      INTEGER ::   ji, jj, jk   ! dummy loop indices
292      !!----------------------------------------------------------------------
293      !
294      DO jj = njs0, njs1
295         DO jk = 1, jpkm1
296            DO ji = 1, jpi
297               va(ji,jj,jk) = ( va(ji,jj,jk) + sshfos_b(ji,jj) ) * vsmsk(jj,jk)
298            END DO
299         END DO
300         DO ji = 1, jpi
301            sshn(ji,jj)= sshn(ji,jj) * (1.-tsmsk(ji,1)) + tsmsk(ji,1) * sshn_b(ji,jj)
302         END DO
303      END DO
304      !
305   END SUBROUTINE obc_dyn_bt_south
306
307# endif
308
309#else
310   !!----------------------------------------------------------------------
311   !!   Default option       No Open Boundaries or not explicit fre surface
312   !!----------------------------------------------------------------------
313CONTAINS
314   SUBROUTINE obc_dyn_bt      ! Dummy routine
315   END SUBROUTINE obc_dyn_bt
316#endif
317
318   !!======================================================================
319END MODULE obcdyn_bt
Note: See TracBrowser for help on using the repository browser.