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/2015/dev_merge_2015/NEMOGCM/NEMO/OPA_SRC/BDY – NEMO

source: branches/2015/dev_merge_2015/NEMOGCM/NEMO/OPA_SRC/BDY/bdydyn3d.F90 @ 6060

Last change on this file since 6060 was 6060, checked in by timgraham, 8 years ago

Merged dev_r5836_noc2_VVL_BY_DEFAULT into branch

  • Property svn:keywords set to Id
File size: 13.4 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   !!----------------------------------------------------------------------
32   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
33   !! $Id$
34   !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
35   !!----------------------------------------------------------------------
36CONTAINS
37
38   SUBROUTINE bdy_dyn3d( kt )
39      !!----------------------------------------------------------------------
40      !!                  ***  SUBROUTINE bdy_dyn3d  ***
41      !!
42      !! ** Purpose : - Apply open boundary conditions for baroclinic velocities
43      !!
44      !!----------------------------------------------------------------------
45      INTEGER, INTENT(in) ::   kt   ! Main time step counter
46      !
47      INTEGER ::   ib_bdy   ! loop index
48      !!----------------------------------------------------------------------
49      !
50      DO ib_bdy=1, nb_bdy
51         !
52         SELECT CASE( cn_dyn3d(ib_bdy) )
53         CASE('none')        ;   CYCLE
54         CASE('frs' )        ;   CALL bdy_dyn3d_frs( idx_bdy(ib_bdy), dta_bdy(ib_bdy), kt, ib_bdy )
55         CASE('specified')   ;   CALL bdy_dyn3d_spe( idx_bdy(ib_bdy), dta_bdy(ib_bdy), kt, ib_bdy )
56         CASE('zero')        ;   CALL bdy_dyn3d_zro( idx_bdy(ib_bdy), dta_bdy(ib_bdy), kt, ib_bdy )
57         CASE('orlanski' )   ;   CALL bdy_dyn3d_orlanski( idx_bdy(ib_bdy), dta_bdy(ib_bdy), ib_bdy, ll_npo=.false. )
58         CASE('orlanski_npo');   CALL bdy_dyn3d_orlanski( idx_bdy(ib_bdy), dta_bdy(ib_bdy), ib_bdy, ll_npo=.true. )
59         CASE DEFAULT        ;   CALL ctl_stop( 'bdy_dyn3d : unrecognised option for open boundaries for baroclinic velocities' )
60         END SELECT
61      END DO
62      !
63   END SUBROUTINE bdy_dyn3d
64
65
66   SUBROUTINE bdy_dyn3d_spe( idx, dta, kt , ib_bdy )
67      !!----------------------------------------------------------------------
68      !!                  ***  SUBROUTINE bdy_dyn3d_spe  ***
69      !!
70      !! ** Purpose : - Apply a specified value for baroclinic velocities
71      !!                at open boundaries.
72      !!
73      !!----------------------------------------------------------------------
74      INTEGER        , INTENT(in) ::   kt      ! time step index
75      TYPE(OBC_INDEX), INTENT(in) ::   idx     ! OBC indices
76      TYPE(OBC_DATA) , INTENT(in) ::   dta     ! OBC external data
77      INTEGER        , INTENT(in) ::   ib_bdy  ! BDY set index
78      !
79      INTEGER  ::   jb, jk         ! dummy loop indices
80      INTEGER  ::   ii, ij, igrd   ! local integers
81      REAL(wp) ::   zwgt           ! boundary weight
82      !!----------------------------------------------------------------------
83      !
84      IF( nn_timing == 1 ) CALL timing_start('bdy_dyn3d_spe')
85      !
86      igrd = 2                      ! Relaxation of zonal velocity
87      DO jb = 1, idx%nblenrim(igrd)
88         DO jk = 1, jpkm1
89            ii   = idx%nbi(jb,igrd)
90            ij   = idx%nbj(jb,igrd)
91            ua(ii,ij,jk) = dta%u3d(jb,jk) * umask(ii,ij,jk)
92         END DO
93      END DO
94      !
95      igrd = 3                      ! Relaxation of meridional velocity
96      DO jb = 1, idx%nblenrim(igrd)
97         DO jk = 1, jpkm1
98            ii   = idx%nbi(jb,igrd)
99            ij   = idx%nbj(jb,igrd)
100            va(ii,ij,jk) = dta%v3d(jb,jk) * vmask(ii,ij,jk)
101         END DO
102      END DO
103      CALL lbc_bdy_lnk( ua, 'U', -1., ib_bdy )   ! Boundary points should be updated 
104      CALL lbc_bdy_lnk( va, 'V', -1., ib_bdy )   
105      !
106      IF( kt == nit000 )   CLOSE( unit = 102 )
107      !
108      IF( nn_timing == 1 ) CALL timing_stop('bdy_dyn3d_spe')
109      !
110   END SUBROUTINE bdy_dyn3d_spe
111
112
113   SUBROUTINE bdy_dyn3d_zro( idx, dta, kt, ib_bdy )
114      !!----------------------------------------------------------------------
115      !!                  ***  SUBROUTINE bdy_dyn3d_zro  ***
116      !!
117      !! ** Purpose : - baroclinic velocities = 0. at open boundaries.
118      !!
119      !!----------------------------------------------------------------------
120      INTEGER        , INTENT(in) ::   kt      ! time step index
121      TYPE(OBC_INDEX), INTENT(in) ::   idx     ! OBC indices
122      TYPE(OBC_DATA) , INTENT(in) ::   dta     ! OBC external data
123      INTEGER,         INTENT(in) ::   ib_bdy  ! BDY set index
124      !
125      INTEGER  ::   ib, ik         ! dummy loop indices
126      INTEGER  ::   ii, ij, igrd   ! local integers
127      REAL(wp) ::   zwgt           ! boundary weight
128      !!----------------------------------------------------------------------
129      !
130      IF( nn_timing == 1 ) CALL timing_start('bdy_dyn3d_zro')
131      !
132      igrd = 2                       ! Everything is at T-points here
133      DO ib = 1, idx%nblenrim(igrd)
134         ii = idx%nbi(ib,igrd)
135         ij = idx%nbj(ib,igrd)
136         DO ik = 1, jpkm1
137            ua(ii,ij,ik) = 0._wp
138         END DO
139      END DO
140
141      igrd = 3                       ! Everything is at T-points here
142      DO ib = 1, idx%nblenrim(igrd)
143         ii = idx%nbi(ib,igrd)
144         ij = idx%nbj(ib,igrd)
145         DO ik = 1, jpkm1
146            va(ii,ij,ik) = 0._wp
147         END DO
148      END DO
149      !
150      CALL lbc_bdy_lnk( ua, 'U', -1., ib_bdy )   ;   CALL lbc_bdy_lnk( va, 'V', -1.,ib_bdy )   ! Boundary points should be updated
151      !
152      IF( kt == nit000 )   CLOSE( unit = 102 )
153      !
154      IF( nn_timing == 1 )   CALL timing_stop('bdy_dyn3d_zro')
155      !
156   END SUBROUTINE bdy_dyn3d_zro
157
158
159   SUBROUTINE bdy_dyn3d_frs( idx, dta, kt, ib_bdy )
160      !!----------------------------------------------------------------------
161      !!                  ***  SUBROUTINE bdy_dyn3d_frs  ***
162      !!
163      !! ** Purpose : - Apply the Flow Relaxation Scheme for baroclinic velocities
164      !!                at open boundaries.
165      !!
166      !! References :- Engedahl H., 1995: Use of the flow relaxation scheme in
167      !!               a three-dimensional baroclinic ocean model with realistic
168      !!               topography. Tellus, 365-382.
169      !!----------------------------------------------------------------------
170      INTEGER        , INTENT(in) ::   kt      ! time step index
171      TYPE(OBC_INDEX), INTENT(in) ::   idx     ! OBC indices
172      TYPE(OBC_DATA) , INTENT(in) ::   dta     ! OBC external data
173      INTEGER,         INTENT(in) ::   ib_bdy  ! BDY set index
174      !
175      INTEGER  ::   jb, jk         ! dummy loop indices
176      INTEGER  ::   ii, ij, igrd   ! local integers
177      REAL(wp) ::   zwgt           ! boundary weight
178      !!----------------------------------------------------------------------
179      !
180      IF( nn_timing == 1 ) CALL timing_start('bdy_dyn3d_frs')
181      !
182      igrd = 2                      ! Relaxation of zonal velocity
183      DO jb = 1, idx%nblen(igrd)
184         DO jk = 1, jpkm1
185            ii   = idx%nbi(jb,igrd)
186            ij   = idx%nbj(jb,igrd)
187            zwgt = idx%nbw(jb,igrd)
188            ua(ii,ij,jk) = ( ua(ii,ij,jk) + zwgt * ( dta%u3d(jb,jk) - ua(ii,ij,jk) ) ) * umask(ii,ij,jk)
189         END DO
190      END DO
191      !
192      igrd = 3                      ! Relaxation of meridional velocity
193      DO jb = 1, idx%nblen(igrd)
194         DO jk = 1, jpkm1
195            ii   = idx%nbi(jb,igrd)
196            ij   = idx%nbj(jb,igrd)
197            zwgt = idx%nbw(jb,igrd)
198            va(ii,ij,jk) = ( va(ii,ij,jk) + zwgt * ( dta%v3d(jb,jk) - va(ii,ij,jk) ) ) * vmask(ii,ij,jk)
199         END DO
200      END DO
201      CALL lbc_bdy_lnk( ua, 'U', -1., ib_bdy )    ! Boundary points should be updated
202      CALL lbc_bdy_lnk( va, 'V', -1., ib_bdy )   
203      !
204      IF( kt == nit000 )   CLOSE( unit = 102 )
205      !
206      IF( nn_timing == 1 )   CALL timing_stop('bdy_dyn3d_frs')
207      !
208   END SUBROUTINE bdy_dyn3d_frs
209
210
211   SUBROUTINE bdy_dyn3d_orlanski( idx, dta, ib_bdy, ll_npo )
212      !!----------------------------------------------------------------------
213      !!                 ***  SUBROUTINE bdy_dyn3d_orlanski  ***
214      !!             
215      !!              - Apply Orlanski radiation to baroclinic velocities.
216      !!              - Wrapper routine for bdy_orlanski_3d
217      !!
218      !!
219      !! References:  Marchesiello, McWilliams and Shchepetkin, Ocean Modelling vol. 3 (2001)   
220      !!----------------------------------------------------------------------
221      TYPE(OBC_INDEX),              INTENT(in) ::   idx  ! OBC indices
222      TYPE(OBC_DATA),               INTENT(in) ::   dta  ! OBC external data
223      INTEGER,                      INTENT(in) ::   ib_bdy  ! BDY set index
224      LOGICAL,                      INTENT(in) ::   ll_npo  ! switch for NPO version
225
226      INTEGER  ::   jb, igrd                               ! dummy loop indices
227      !!----------------------------------------------------------------------
228
229      IF( nn_timing == 1 ) CALL timing_start('bdy_dyn3d_orlanski')
230      !
231      !! Note that at this stage the ub and ua arrays contain the baroclinic velocities.
232      !
233      igrd = 2      ! Orlanski bc on u-velocity;
234      !           
235      CALL bdy_orlanski_3d( idx, igrd, ub, ua, dta%u3d, ll_npo )
236
237      igrd = 3      ! Orlanski bc on v-velocity
238     
239      CALL bdy_orlanski_3d( idx, igrd, vb, va, dta%v3d, ll_npo )
240      !
241      CALL lbc_bdy_lnk( ua, 'U', -1., ib_bdy )    ! Boundary points should be updated
242      CALL lbc_bdy_lnk( va, 'V', -1., ib_bdy )   
243      !
244      IF( nn_timing == 1 ) CALL timing_stop('bdy_dyn3d_orlanski')
245      !
246   END SUBROUTINE bdy_dyn3d_orlanski
247
248
249   SUBROUTINE bdy_dyn3d_dmp( kt )
250      !!----------------------------------------------------------------------
251      !!                  ***  SUBROUTINE bdy_dyn3d_dmp  ***
252      !!
253      !! ** Purpose : Apply damping for baroclinic velocities at open boundaries.
254      !!
255      !!----------------------------------------------------------------------
256      INTEGER, INTENT(in) ::   kt   ! time step index
257      !
258      INTEGER  ::   jb, jk         ! dummy loop indices
259      INTEGER  ::   ib_bdy         ! loop index
260      INTEGER  ::   ii, ij, igrd   ! local integers
261      REAL(wp) ::   zwgt           ! boundary weight
262      !!----------------------------------------------------------------------
263      !
264      IF( nn_timing == 1 )   CALL timing_start('bdy_dyn3d_dmp')
265      !
266      DO ib_bdy=1, nb_bdy
267         IF ( ln_dyn3d_dmp(ib_bdy) .and. cn_dyn3d(ib_bdy) /= 'none' ) THEN
268            igrd = 2                      ! Relaxation of zonal velocity
269            DO jb = 1, idx_bdy(ib_bdy)%nblen(igrd)
270               ii   = idx_bdy(ib_bdy)%nbi(jb,igrd)
271               ij   = idx_bdy(ib_bdy)%nbj(jb,igrd)
272               zwgt = idx_bdy(ib_bdy)%nbd(jb,igrd)
273               DO jk = 1, jpkm1
274                  ua(ii,ij,jk) = ( ua(ii,ij,jk) + zwgt * ( dta_bdy(ib_bdy)%u3d(jb,jk) - &
275                                   ub(ii,ij,jk) + ub_b(ii,ij)) ) * umask(ii,ij,jk)
276               END DO
277            END DO
278            !
279            igrd = 3                      ! Relaxation of meridional velocity
280            DO jb = 1, idx_bdy(ib_bdy)%nblen(igrd)
281               ii   = idx_bdy(ib_bdy)%nbi(jb,igrd)
282               ij   = idx_bdy(ib_bdy)%nbj(jb,igrd)
283               zwgt = idx_bdy(ib_bdy)%nbd(jb,igrd)
284               DO jk = 1, jpkm1
285                  va(ii,ij,jk) = ( va(ii,ij,jk) + zwgt * ( dta_bdy(ib_bdy)%v3d(jb,jk) -  &
286                                   vb(ii,ij,jk) + vb_b(ii,ij)) ) * vmask(ii,ij,jk)
287               END DO
288            END DO
289         ENDIF
290      END DO
291      !
292      CALL lbc_lnk( ua, 'U', -1. )   ;   CALL lbc_lnk( va, 'V', -1. )   ! Boundary points should be updated
293      !
294      IF( nn_timing == 1 )   CALL timing_stop('bdy_dyn3d_dmp')
295      !
296   END SUBROUTINE bdy_dyn3d_dmp
297
298#else
299   !!----------------------------------------------------------------------
300   !!   Dummy module                   NO Unstruct Open Boundary Conditions
301   !!----------------------------------------------------------------------
302CONTAINS
303   SUBROUTINE bdy_dyn3d( kt )      ! Empty routine
304      WRITE(*,*) 'bdy_dyn3d: You should not have seen this print! error?', kt
305   END SUBROUTINE bdy_dyn3d
306   SUBROUTINE bdy_dyn3d_dmp( kt )      ! Empty routine
307      WRITE(*,*) 'bdy_dyn3d_dmp: You should not have seen this print! error?', kt
308   END SUBROUTINE bdy_dyn3d_dmp
309#endif
310
311   !!======================================================================
312END MODULE bdydyn3d
Note: See TracBrowser for help on using the repository browser.