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.
bdydyn3d.F90 in branches/UKMO/CO6_KD490_amm7oper/NEMOGCM/NEMO/OPA_SRC/BDY – NEMO

source: branches/UKMO/CO6_KD490_amm7oper/NEMOGCM/NEMO/OPA_SRC/BDY/bdydyn3d.F90 @ 6991

Last change on this file since 6991 was 6991, checked in by kingr, 8 years ago

Adding Momme's bdy changes for zerograd and neumann conditions.

File size: 17.0 KB
Line 
1MODULE bdydyn3d
2   !!======================================================================
3   !!                       ***  MODULE  bdydyn3d  ***
4   !! Unstructured Open Boundary Cond. :   Flow relaxation scheme on baroclinic velocities
5   !!======================================================================
6   !! History :  3.4  !  2011     (D. Storkey) new module as part of BDY rewrite
7   !!            3.5  !  2012     (S. Mocavero, I. Epicoco) Optimization of BDY communications
8   !!----------------------------------------------------------------------
9#if defined key_bdy 
10   !!----------------------------------------------------------------------
11   !!   'key_bdy' :                    Unstructured Open Boundary Condition
12   !!----------------------------------------------------------------------
13   !!   bdy_dyn3d        : apply open boundary conditions to baroclinic velocities
14   !!   bdy_dyn3d_frs    : apply Flow Relaxation Scheme
15   !!----------------------------------------------------------------------
16   USE timing          ! Timing
17   USE oce             ! ocean dynamics and tracers
18   USE dom_oce         ! ocean space and time domain
19   USE bdy_oce         ! ocean open boundary conditions
20   USE bdylib          ! for orlanski library routines
21   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
22   USE in_out_manager  !
23   Use phycst
24
25   IMPLICIT NONE
26   PRIVATE
27
28   PUBLIC   bdy_dyn3d     ! routine called by bdy_dyn
29   PUBLIC   bdy_dyn3d_dmp ! routine called by step
30
31   !! * Substitutions
32#  include "domzgr_substitute.h90"
33   !!----------------------------------------------------------------------
34   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
35   !! $Id$
36   !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
37   !!----------------------------------------------------------------------
38CONTAINS
39
40   SUBROUTINE bdy_dyn3d( kt )
41      !!----------------------------------------------------------------------
42      !!                  ***  SUBROUTINE bdy_dyn3d  ***
43      !!
44      !! ** Purpose : - Apply open boundary conditions for baroclinic velocities
45      !!
46      !!----------------------------------------------------------------------
47      INTEGER, INTENT( in ) :: kt     ! Main time step counter
48      !!
49      INTEGER               :: ib_bdy ! loop index
50      !!
51
52      DO ib_bdy=1, nb_bdy
53
54         SELECT CASE( cn_dyn3d(ib_bdy) )
55         CASE('none')
56            CYCLE
57         CASE('frs')
58            CALL bdy_dyn3d_frs( idx_bdy(ib_bdy), dta_bdy(ib_bdy), kt, ib_bdy )
59         CASE('specified')
60            CALL bdy_dyn3d_spe( idx_bdy(ib_bdy), dta_bdy(ib_bdy), kt, ib_bdy )
61         CASE('zerograd')
62            CALL bdy_dyn3d_zgrad( idx_bdy(ib_bdy), dta_bdy(ib_bdy), kt, ib_bdy )
63         CASE('zero')
64            CALL bdy_dyn3d_zro( idx_bdy(ib_bdy), dta_bdy(ib_bdy), kt, ib_bdy )
65         CASE('neumann')
66            CALL bdy_dyn3d_nmn( idx_bdy(ib_bdy), ib_bdy )
67         CASE('orlanski')
68            CALL bdy_dyn3d_orlanski( idx_bdy(ib_bdy), dta_bdy(ib_bdy), ib_bdy, ll_npo=.false. )
69         CASE('orlanski_npo')
70            CALL bdy_dyn3d_orlanski( idx_bdy(ib_bdy), dta_bdy(ib_bdy), ib_bdy, ll_npo=.true. )
71         CASE DEFAULT
72            CALL ctl_stop( 'bdy_dyn3d : unrecognised option for open boundaries for baroclinic velocities' )
73         END SELECT
74      ENDDO
75
76   END SUBROUTINE bdy_dyn3d
77
78   SUBROUTINE bdy_dyn3d_spe( idx, dta, kt , ib_bdy )
79      !!----------------------------------------------------------------------
80      !!                  ***  SUBROUTINE bdy_dyn3d_spe  ***
81      !!
82      !! ** Purpose : - Apply a specified value for baroclinic velocities
83      !!                at open boundaries.
84      !!
85      !!----------------------------------------------------------------------
86      INTEGER                     ::   kt
87      TYPE(OBC_INDEX), INTENT(in) ::   idx  ! OBC indices
88      TYPE(OBC_DATA),  INTENT(in) ::   dta  ! OBC external data
89      INTEGER,         INTENT(in) ::   ib_bdy  ! BDY set index
90      !!
91      INTEGER  ::   jb, jk         ! dummy loop indices
92      INTEGER  ::   ii, ij, igrd   ! local integers
93      REAL(wp) ::   zwgt           ! boundary weight
94      !!----------------------------------------------------------------------
95      !
96      IF( nn_timing == 1 ) CALL timing_start('bdy_dyn3d_spe')
97      !
98      igrd = 2                      ! Relaxation of zonal velocity
99      DO jb = 1, idx%nblenrim(igrd)
100         DO jk = 1, jpkm1
101            ii   = idx%nbi(jb,igrd)
102            ij   = idx%nbj(jb,igrd)
103            ua(ii,ij,jk) = dta%u3d(jb,jk) * umask(ii,ij,jk)
104         END DO
105      END DO
106      !
107      igrd = 3                      ! Relaxation of meridional velocity
108      DO jb = 1, idx%nblenrim(igrd)
109         DO jk = 1, jpkm1
110            ii   = idx%nbi(jb,igrd)
111            ij   = idx%nbj(jb,igrd)
112            va(ii,ij,jk) = dta%v3d(jb,jk) * vmask(ii,ij,jk)
113         END DO
114      END DO
115      CALL lbc_bdy_lnk( ua, 'U', -1., ib_bdy )   ! Boundary points should be updated 
116      CALL lbc_bdy_lnk( va, 'V', -1., ib_bdy )   
117      !
118      IF( kt .eq. nit000 ) CLOSE( unit = 102 )
119
120      IF( nn_timing == 1 ) CALL timing_stop('bdy_dyn3d_spe')
121
122   END SUBROUTINE bdy_dyn3d_spe
123
124   SUBROUTINE bdy_dyn3d_zgrad( idx, dta, kt , ib_bdy )
125      !!----------------------------------------------------------------------
126      !!                  ***  SUBROUTINE bdy_dyn3d_zgrad  ***
127      !!
128      !! ** Purpose : - Enforce a zero gradient of normal velocity
129      !!
130      !!----------------------------------------------------------------------
131      INTEGER                     ::   kt
132      TYPE(OBC_INDEX), INTENT(in) ::   idx  ! OBC indices
133      TYPE(OBC_DATA),  INTENT(in) ::   dta  ! OBC external data
134      INTEGER,         INTENT(in) ::   ib_bdy  ! BDY set index
135      !!
136      INTEGER  ::   jb, jk         ! dummy loop indices
137      INTEGER  ::   ii, ij, igrd   ! local integers
138      REAL(wp) ::   zwgt           ! boundary weight
139      INTEGER  ::   fu, fv
140      !!----------------------------------------------------------------------
141      !
142      IF( nn_timing == 1 ) CALL timing_start('bdy_dyn3d_zgrad')
143      !
144      igrd = 2                      ! Copying tangential velocity into bdy points
145      DO jb = 1, idx%nblenrim(igrd)
146         DO jk = 1, jpkm1
147            ii   = idx%nbi(jb,igrd)
148            ij   = idx%nbj(jb,igrd)
149            fu   = ABS( ABS (NINT( idx%flagu(jb,igrd) ) ) - 1 )
150            ua(ii,ij,jk) = ua(ii,ij,jk) * REAL( 1 - fu) + ( ua(ii,ij+fu,jk) * umask(ii,ij+fu,jk) &
151                        &+ ua(ii,ij-fu,jk) * umask(ii,ij-fu,jk) ) * umask(ii,ij,jk) * REAL( fu )
152         END DO
153      END DO
154      !
155      igrd = 3                      ! Copying tangential velocity into bdy points
156      DO jb = 1, idx%nblenrim(igrd)
157         DO jk = 1, jpkm1
158            ii   = idx%nbi(jb,igrd)
159            ij   = idx%nbj(jb,igrd)
160            fv   = ABS( ABS (NINT( idx%flagv(jb,igrd) ) ) - 1 )
161            va(ii,ij,jk) = va(ii,ij,jk) * REAL( 1 - fv ) + ( va(ii+fv,ij,jk) * vmask(ii+fv,ij,jk) &
162                        &+ va(ii-fv,ij,jk) * vmask(ii-fv,ij,jk) ) * vmask(ii,ij,jk) * REAL( fv )
163         END DO
164      END DO
165      CALL lbc_bdy_lnk( ua, 'U', -1., ib_bdy )   ! Boundary points should be updated
166      CALL lbc_bdy_lnk( va, 'V', -1., ib_bdy )   
167      !
168      IF( kt .eq. nit000 ) CLOSE( unit = 102 )
169      IF( nn_timing == 1 ) CALL timing_stop('bdy_dyn3d_zgrad')
170   END SUBROUTINE bdy_dyn3d_zgrad
171
172   SUBROUTINE bdy_dyn3d_zro( idx, dta, kt, ib_bdy )
173      !!----------------------------------------------------------------------
174      !!                  ***  SUBROUTINE bdy_dyn3d_zro  ***
175      !!
176      !! ** Purpose : - baroclinic velocities = 0. at open boundaries.
177      !!
178      !!----------------------------------------------------------------------
179      INTEGER                     ::   kt
180      TYPE(OBC_INDEX), INTENT(in) ::   idx  ! OBC indices
181      TYPE(OBC_DATA),  INTENT(in) ::   dta  ! OBC external data
182      INTEGER,         INTENT(in) ::   ib_bdy  ! BDY set index
183      !!
184      INTEGER  ::   ib, ik         ! dummy loop indices
185      INTEGER  ::   ii, ij, igrd, zcoef   ! local integers
186      REAL(wp) ::   zwgt           ! boundary weight
187      !!----------------------------------------------------------------------
188      !
189      IF( nn_timing == 1 ) CALL timing_start('bdy_dyn3d_zro')
190      !
191      igrd = 2                       ! Everything is at T-points here
192      DO ib = 1, idx%nblenrim(igrd)
193         ii = idx%nbi(ib,igrd)
194         ij = idx%nbj(ib,igrd)
195         DO ik = 1, jpkm1
196            ua(ii,ij,ik) = 0._wp
197         END DO
198      END DO
199
200      igrd = 3                       ! Everything is at T-points here
201      DO ib = 1, idx%nblenrim(igrd)
202         ii = idx%nbi(ib,igrd)
203         ij = idx%nbj(ib,igrd)
204         DO ik = 1, jpkm1
205            va(ii,ij,ik) = 0._wp
206         END DO
207      END DO
208      !
209      CALL lbc_bdy_lnk( ua, 'U', -1., ib_bdy )   ;   CALL lbc_bdy_lnk( va, 'V', -1.,ib_bdy )   ! Boundary points should be updated
210      !
211      IF( kt .eq. nit000 ) CLOSE( unit = 102 )
212
213      IF( nn_timing == 1 ) CALL timing_stop('bdy_dyn3d_zro')
214
215   END SUBROUTINE bdy_dyn3d_zro
216
217   SUBROUTINE bdy_dyn3d_frs( idx, dta, kt, ib_bdy )
218      !!----------------------------------------------------------------------
219      !!                  ***  SUBROUTINE bdy_dyn3d_frs  ***
220      !!
221      !! ** Purpose : - Apply the Flow Relaxation Scheme for baroclinic velocities
222      !!                at open boundaries.
223      !!
224      !! References :- Engedahl H., 1995: Use of the flow relaxation scheme in
225      !!               a three-dimensional baroclinic ocean model with realistic
226      !!               topography. Tellus, 365-382.
227      !!----------------------------------------------------------------------
228      INTEGER                     ::   kt
229      TYPE(OBC_INDEX), INTENT(in) ::   idx  ! OBC indices
230      TYPE(OBC_DATA),  INTENT(in) ::   dta  ! OBC external data
231      INTEGER,         INTENT(in) ::   ib_bdy  ! BDY set index
232      !!
233      INTEGER  ::   jb, jk         ! dummy loop indices
234      INTEGER  ::   ii, ij, igrd   ! local integers
235      REAL(wp) ::   zwgt           ! boundary weight
236      !!----------------------------------------------------------------------
237      !
238      IF( nn_timing == 1 ) CALL timing_start('bdy_dyn3d_frs')
239      !
240      igrd = 2                      ! Relaxation of zonal velocity
241      DO jb = 1, idx%nblen(igrd)
242         DO jk = 1, jpkm1
243            ii   = idx%nbi(jb,igrd)
244            ij   = idx%nbj(jb,igrd)
245            zwgt = idx%nbw(jb,igrd)
246            ua(ii,ij,jk) = ( ua(ii,ij,jk) + zwgt * ( dta%u3d(jb,jk) - ua(ii,ij,jk) ) ) * umask(ii,ij,jk)
247         END DO
248      END DO
249      !
250      igrd = 3                      ! Relaxation of meridional velocity
251      DO jb = 1, idx%nblen(igrd)
252         DO jk = 1, jpkm1
253            ii   = idx%nbi(jb,igrd)
254            ij   = idx%nbj(jb,igrd)
255            zwgt = idx%nbw(jb,igrd)
256            va(ii,ij,jk) = ( va(ii,ij,jk) + zwgt * ( dta%v3d(jb,jk) - va(ii,ij,jk) ) ) * vmask(ii,ij,jk)
257         END DO
258      END DO
259      CALL lbc_bdy_lnk( ua, 'U', -1., ib_bdy )    ! Boundary points should be updated
260      CALL lbc_bdy_lnk( va, 'V', -1., ib_bdy )   
261      !
262      IF( kt .eq. nit000 ) CLOSE( unit = 102 )
263
264      IF( nn_timing == 1 ) CALL timing_stop('bdy_dyn3d_frs')
265
266   END SUBROUTINE bdy_dyn3d_frs
267
268   SUBROUTINE bdy_dyn3d_orlanski( idx, dta, ib_bdy, ll_npo )
269      !!----------------------------------------------------------------------
270      !!                 ***  SUBROUTINE bdy_dyn3d_orlanski  ***
271      !!             
272      !!              - Apply Orlanski radiation to baroclinic velocities.
273      !!              - Wrapper routine for bdy_orlanski_3d
274      !!
275      !!
276      !! References:  Marchesiello, McWilliams and Shchepetkin, Ocean Modelling vol. 3 (2001)   
277      !!----------------------------------------------------------------------
278      TYPE(OBC_INDEX),              INTENT(in) ::   idx  ! OBC indices
279      TYPE(OBC_DATA),               INTENT(in) ::   dta  ! OBC external data
280      INTEGER,                      INTENT(in) ::   ib_bdy  ! BDY set index
281      LOGICAL,                      INTENT(in) ::   ll_npo  ! switch for NPO version
282
283      INTEGER  ::   jb, igrd                               ! dummy loop indices
284      !!----------------------------------------------------------------------
285
286      IF( nn_timing == 1 ) CALL timing_start('bdy_dyn3d_orlanski')
287      !
288      !! Note that at this stage the ub and ua arrays contain the baroclinic velocities.
289      !
290      igrd = 2      ! Orlanski bc on u-velocity;
291      !           
292      CALL bdy_orlanski_3d( idx, igrd, ub, ua, dta%u3d, ll_npo )
293
294      igrd = 3      ! Orlanski bc on v-velocity
295     
296      CALL bdy_orlanski_3d( idx, igrd, vb, va, dta%v3d, ll_npo )
297      !
298      CALL lbc_bdy_lnk( ua, 'U', -1., ib_bdy )    ! Boundary points should be updated
299      CALL lbc_bdy_lnk( va, 'V', -1., ib_bdy )   
300      !
301      IF( nn_timing == 1 ) CALL timing_stop('bdy_dyn3d_orlanski')
302      !
303   END SUBROUTINE bdy_dyn3d_orlanski
304
305
306   SUBROUTINE bdy_dyn3d_dmp( kt )
307      !!----------------------------------------------------------------------
308      !!                  ***  SUBROUTINE bdy_dyn3d_dmp  ***
309      !!
310      !! ** Purpose : Apply damping for baroclinic velocities at open boundaries.
311      !!
312      !!----------------------------------------------------------------------
313      INTEGER                     ::   kt
314      !!
315      INTEGER  ::   jb, jk         ! dummy loop indices
316      INTEGER  ::   ii, ij, igrd   ! local integers
317      REAL(wp) ::   zwgt           ! boundary weight
318      INTEGER  ::  ib_bdy          ! loop index
319      !!----------------------------------------------------------------------
320      !
321      IF( nn_timing == 1 ) CALL timing_start('bdy_dyn3d_dmp')
322      !
323      !-------------------------------------------------------
324
325      DO ib_bdy=1, nb_bdy
326         IF ( ln_dyn3d_dmp(ib_bdy) .and. cn_dyn3d(ib_bdy) /= 'none' ) THEN
327            igrd = 2                      ! Relaxation of zonal velocity
328            DO jb = 1, idx_bdy(ib_bdy)%nblen(igrd)
329               ii   = idx_bdy(ib_bdy)%nbi(jb,igrd)
330               ij   = idx_bdy(ib_bdy)%nbj(jb,igrd)
331               zwgt = idx_bdy(ib_bdy)%nbd(jb,igrd)
332               DO jk = 1, jpkm1
333                  ua(ii,ij,jk) = ( ua(ii,ij,jk) + zwgt * ( dta_bdy(ib_bdy)%u3d(jb,jk) - &
334                                   ub(ii,ij,jk) + ub_b(ii,ij)) ) * umask(ii,ij,jk)
335               END DO
336            END DO
337            !
338            igrd = 3                      ! Relaxation of meridional velocity
339            DO jb = 1, idx_bdy(ib_bdy)%nblen(igrd)
340               ii   = idx_bdy(ib_bdy)%nbi(jb,igrd)
341               ij   = idx_bdy(ib_bdy)%nbj(jb,igrd)
342               zwgt = idx_bdy(ib_bdy)%nbd(jb,igrd)
343               DO jk = 1, jpkm1
344                  va(ii,ij,jk) = ( va(ii,ij,jk) + zwgt * ( dta_bdy(ib_bdy)%v3d(jb,jk) -  &
345                                   vb(ii,ij,jk) + vb_b(ii,ij)) ) * vmask(ii,ij,jk)
346               END DO
347            END DO
348         ENDIF
349      ENDDO
350      !
351      CALL lbc_lnk( ua, 'U', -1. )   ;   CALL lbc_lnk( va, 'V', -1. )   ! Boundary points should be updated
352      !
353      IF( nn_timing == 1 ) CALL timing_stop('bdy_dyn3d_dmp')
354
355   END SUBROUTINE bdy_dyn3d_dmp
356
357   SUBROUTINE bdy_dyn3d_nmn( idx, ib_bdy )
358      !!----------------------------------------------------------------------
359      !!                 ***  SUBROUTINE bdy_dyn3d_nmn  ***
360      !!             
361      !!              - Apply Neumann condition to baroclinic velocities.
362      !!              - Wrapper routine for bdy_nmn
363      !!
364      !!
365      !!----------------------------------------------------------------------
366      TYPE(OBC_INDEX),              INTENT(in) ::   idx  ! OBC indices
367      INTEGER,                      INTENT(in) ::   ib_bdy  ! BDY set index
368      INTEGER  ::   jb, igrd                               ! dummy loop indices
369      !!----------------------------------------------------------------------
370      IF( nn_timing == 1 ) CALL timing_start('bdy_dyn3d_nmn')
371      !
372      !! Note that at this stage the ub and ua arrays contain the baroclinic velocities.
373      !
374      igrd = 2      ! Neumann bc on u-velocity;
375      !           
376      CALL bdy_nmn( idx, igrd, ua )
377      igrd = 3      ! Neumann bc on v-velocity
378      !
379      CALL bdy_nmn( idx, igrd, va )
380      !
381      CALL lbc_bdy_lnk( ua, 'U', -1., ib_bdy )    ! Boundary points should be updated
382      CALL lbc_bdy_lnk( va, 'V', -1., ib_bdy )
383      !
384      IF( nn_timing == 1 ) CALL timing_stop('bdy_dyn3d_nmn')
385      !
386   END SUBROUTINE bdy_dyn3d_nmn
387   
388#else
389   !!----------------------------------------------------------------------
390   !!   Dummy module                   NO Unstruct Open Boundary Conditions
391   !!----------------------------------------------------------------------
392CONTAINS
393   SUBROUTINE bdy_dyn3d( kt )      ! Empty routine
394      WRITE(*,*) 'bdy_dyn3d: You should not have seen this print! error?', kt
395   END SUBROUTINE bdy_dyn3d
396
397   SUBROUTINE bdy_dyn3d_dmp( kt )      ! Empty routine
398      WRITE(*,*) 'bdy_dyn3d_dmp: You should not have seen this print! error?', kt
399   END SUBROUTINE bdy_dyn3d_dmp
400
401#endif
402
403   !!======================================================================
404END MODULE bdydyn3d
Note: See TracBrowser for help on using the repository browser.