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.
obcfla.F90 in branches/DEV_r1879_FCM/NEMOGCM/NEMO/OPA_SRC/OBC – NEMO

source: branches/DEV_r1879_FCM/NEMOGCM/NEMO/OPA_SRC/OBC/obcfla.F90 @ 2007

Last change on this file since 2007 was 2007, checked in by smasson, 14 years ago

update branches/DEV_r1879_FCM/NEMOGCM/NEMO with tags/nemo_v3_2_1/NEMO

  • Property svn:eol-style set to native
  • Property svn:keywords set to Id
File size: 9.4 KB
Line 
1MODULE obcfla
2#if defined key_obc && defined key_dynspg_ts
3   !!=================================================================================
4   !!                       ***  MODULE  obcfla  ***
5   !! Ocean dynamics:   Flather's algorithm at open boundaries for the time-splitting
6   !!=================================================================================
7
8   !!---------------------------------------------------------------------------------
9   !!   obc_fla_ts        : call the subroutine for each open boundary
10   !!   obc_fla_ts_east   : Flather on the east  open boundary velocities & ssh
11   !!   obc_fla_ts_west   : Flather on the west  open boundary velocities & ssh
12   !!   obc_fla_ts_north  : Flather on the north open boundary velocities & ssh
13   !!   obc_fla_ts_south  : Flather on the south open boundary velocities & ssh
14   !!----------------------------------------------------------------------------------
15
16   !!----------------------------------------------------------------------------------
17   !! * Modules used
18   USE oce             ! ocean dynamics and tracers
19   USE dom_oce         ! ocean space and time domain
20   USE dynspg_oce      ! surface pressure gradient variables
21   USE phycst          ! physical constants
22   USE obc_oce         ! ocean open boundary conditions
23   USE obcdta          ! ocean open boundary conditions: climatology
24
25   IMPLICIT NONE
26   PRIVATE
27
28   !! * Accessibility
29   PUBLIC obc_fla_ts  ! routine called in dynspg_ts (free surface time splitting case)
30
31   !!---------------------------------------------------------------------------------
32   !!  OPA 9.0 , LOCEAN-IPSL (2005)
33   !! $Id$
34   !! This software is governed by the CeCILL licence see modipsl/doc/NEMO_CeCILL.txt
35   !!---------------------------------------------------------------------------------
36
37CONTAINS
38
39   SUBROUTINE obc_fla_ts
40      !!------------------------------------------------------------------------------
41      !!                      SUBROUTINE obc_fla_ts
42      !!                     **********************
43      !! ** Purpose :
44      !!      Apply Flather's algorithm at open boundaries for the time-splitting
45      !!      free surface case (barotropic variables)
46      !!
47      !!      This routine is called in dynspg_ts.F90 routine
48      !!
49      !!      The logical variable lp_obc_east, and/or lp_obc_west, and/or lp_obc_north,
50      !!      and/or lp_obc_south allow the user to determine which boundary is an
51      !!      open one (must be done in the obc_par.F90 file).
52      !!
53      !! ** Reference :
54      !!         Flather, R. A., 1976, Mem. Soc. R. Sci. Liege, Ser. 6, 10, 141-164
55      !!
56      !! History :
57      !!   9.0  !  05-12  (V. Garnier) original
58      !!------------------------------------------------------------------------------
59
60      IF( lp_obc_east  )   CALL obc_fla_ts_east 
61      IF( lp_obc_west  )   CALL obc_fla_ts_west 
62      IF( lp_obc_north )   CALL obc_fla_ts_north
63      IF( lp_obc_south )   CALL obc_fla_ts_south
64
65   END SUBROUTINE obc_fla_ts
66
67
68   SUBROUTINE obc_fla_ts_east
69      !!------------------------------------------------------------------------------
70      !!                  ***  SUBROUTINE obc_fla_ts_east  ***
71      !!
72      !! ** Purpose :
73      !!      Apply Flather's algorithm on east OBC velocities ua, va
74      !!      Fix sea surface height (sshn_e) on east open boundary
75      !!
76      !!  History :
77      !!   9.0  !  05-12  (V. Garnier) original
78      !!------------------------------------------------------------------------------
79      !! * Local declaration
80      INTEGER ::   ji, jj, jk ! dummy loop indices
81      !!------------------------------------------------------------------------------
82
83      DO ji = nie0, nie1
84         DO jj = 1, jpj
85            ua_e(ji,jj) = (  ubtfoe(jj) * hur(ji,jj) + sqrt( grav*hur(ji,jj) )   &
86               &            * ( ( sshn_e(ji,jj) + sshn_e(ji+1,jj) ) * 0.5        &
87               &            - sshfoe(jj) )  ) * uemsk(jj,jk)
88         END DO
89      END DO
90      DO ji = nie0p1, nie1p1
91         DO jj = 1, jpj
92            sshfoe_b(ji,jj) = sshfoe_b(ji,jj) + sqrt( grav*hur(ji,jj) )     &
93               &             * ( ( sshn_e(ji,jj) + sshn_e(ji+1,jj) ) * 0.5  &
94               &                 - sshfoe(jj) ) * uemsk(jj,1)
95            ssha_e(ji,jj) = ssha_e(ji,jj) * ( 1. - temsk(jj,1) ) &
96               &            + temsk(jj,1) * sshfoe(jj)
97            va_e(ji,jj) = vbtfoe(jj) * hvr(ji,jj) * uemsk(jj,jk)
98         END DO
99      END DO
100
101   END SUBROUTINE obc_fla_ts_east
102
103
104   SUBROUTINE obc_fla_ts_west
105      !!------------------------------------------------------------------------------
106      !!                  ***  SUBROUTINE obc_fla_ts_west  ***
107      !!
108      !! ** Purpose :
109      !!      Apply Flather's algorithm on west OBC velocities ua, va
110      !!      Fix sea surface height (sshn_e) on west open boundary
111      !!
112      !!  History :
113      !!   9.0  !  05-12  (V. Garnier) original
114      !!------------------------------------------------------------------------------
115      !! * Local declaration
116      INTEGER ::   ji, jj, jk ! dummy loop indices
117      !!------------------------------------------------------------------------------
118
119      DO ji = niw0, niw1
120         DO jj = 1, jpj
121            ua_e(ji,jj) = ( ubtfow(jj) * hur(ji,jj) - sqrt( grav * hur(ji,jj) )   &
122               &            * ( ( sshn_e(ji,jj) + sshn_e(ji+1,jj) ) * 0.5         &
123               &                - sshfow(jj) ) ) * uwmsk(jj,jk)
124            va_e(ji,jj) = vbtfow(jj) * hvr(ji,jj) * uwmsk(jj,jk)
125         END DO
126         DO jj = 1, jpj
127            sshfow_b(ji,jj) = sshfow_b(ji,jj) - sqrt( grav * hur(ji,jj) )     &
128                              * ( ( sshn_e(ji,jj) + sshn_e(ji+1,jj) ) * 0.5   &
129                                 - sshfow(jj) ) * uwmsk(jj,1)
130            ssha_e(ji,jj) = ssha_e(ji,jj) * ( 1. - twmsk(jj,1) ) &
131               &            + twmsk(jj,1)*sshfow(jj)
132         END DO
133      END DO
134
135   END SUBROUTINE obc_fla_ts_west
136
137   SUBROUTINE obc_fla_ts_north
138      !!------------------------------------------------------------------------------
139      !!                     SUBROUTINE obc_fla_ts_north
140      !!                    *************************
141      !! ** Purpose :
142      !!      Apply Flather's algorithm on north OBC velocities ua, va
143      !!      Fix sea surface height (sshn_e) on north open boundary
144      !!
145      !!  History :
146      !!   9.0  !  05-12  (V. Garnier) original
147      !!------------------------------------------------------------------------------
148      !! * Local declaration
149      INTEGER ::   ji, jj, jk ! dummy loop indices
150      !!------------------------------------------------------------------------------
151
152      DO jj = njn0, njn1
153         DO ji = 1, jpi
154            va_e(ji,jj) = ( vbtfon(ji) * hvr(ji,jj) + sqrt( grav * hvr(ji,jj) )   &
155               &            * ( ( sshn_e(ji,jj) + sshn_e(ji,jj+1) ) * 0.5         &
156               &                - sshfon(ji) ) ) * vnmsk(ji,jk)
157         END DO
158      END DO
159      DO jj = njn0p1, njn1p1
160         DO ji = 1, jpi
161            sshfon_b(ji,jj) = sshfon_b(ji,jj) + sqrt( grav * hvr(ji,jj) )  &
162               &              * ( ( sshn_e(ji,jj) + sshn_e(ji,jj+1) ) * 0.5    &
163               &                  - sshfon(ji) ) * vnmsk(ji,1)
164            ssha_e(ji,jj) = ssha_e(ji,jj) * ( 1. - tnmsk(ji,1) ) &
165               &            + sshfon(ji) * tnmsk(ji,1)
166            ua_e(ji,jj) = ubtfon(ji) * hur(ji,jj) * vnmsk(ji,jk)
167         END DO
168      END DO
169
170   END SUBROUTINE obc_fla_ts_north
171
172   SUBROUTINE obc_fla_ts_south
173      !!------------------------------------------------------------------------------
174      !!                     SUBROUTINE obc_fla_ts_south
175      !!                    *************************
176      !! ** Purpose :
177      !!      Apply Flather's algorithm on south OBC velocities ua, va
178      !!      Fix sea surface height (sshn_e) on south open boundary
179      !!
180      !!  History :
181      !!   9.0  !  05-12  (V. Garnier) original
182      !!------------------------------------------------------------------------------
183      !! * Local declaration
184      INTEGER ::   ji, jj, jk ! dummy loop indices
185
186      !!------------------------------------------------------------------------------
187
188      DO jj = njs0, njs1
189         DO ji = 1, jpi
190            va_e(ji,jj) = ( vbtfos(ji) * hvr(ji,jj) - sqrt( grav * hvr(ji,jj) )   &
191               &            * ( ( sshn_e(ji,jj) + sshn_e(ji,jj+1) ) * 0.5         &
192               &                - sshfos(ji) ) ) * vsmsk(ji,jk)
193            ua_e(ji,jj) = ubtfos(ji) * hur(ji,jj) * vsmsk(ji,jk)
194         END DO
195         DO ji = 1, jpi
196            sshfos_b(ji,jj) = sshfos_b(ji,jj) - sqrt( grav * hvr(ji,jj) )      &
197               &              * ( ( sshn_e(ji,jj) + sshn_e(ji,jj+1) ) * 0.5    &
198               &                  - sshfos(ji) ) * vsmsk(ji,1)
199            ssha_e(ji,jj) = ssha_e(ji,jj) * (1. - tsmsk(ji,1) ) &
200               &            + tsmsk(ji,1) * sshfos(ji)
201         END DO
202      END DO
203
204   END SUBROUTINE obc_fla_ts_south
205#else
206   !!=================================================================================
207   !!                       ***  MODULE  obcfla  ***
208   !! Ocean dynamics:   Flather's algorithm at open boundaries for the time-splitting
209   !!=================================================================================
210CONTAINS
211
212   SUBROUTINE obc_fla_ts
213      WRITE(*,*) 'obc_fla_ts: You should not have seen this print! error?'
214   END SUBROUTINE obc_fla_ts
215#endif
216
217END MODULE obcfla
Note: See TracBrowser for help on using the repository browser.