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/2013/dev_r3891_METO1_MERCATOR6_OBC_BDY_merge/NEMOGCM/NEMO/OPA_SRC/BDY – NEMO

source: branches/2013/dev_r3891_METO1_MERCATOR6_OBC_BDY_merge/NEMOGCM/NEMO/OPA_SRC/BDY/bdydyn3d.F90 @ 3900

Last change on this file since 3900 was 3900, checked in by davestorkey, 11 years ago

First sketch of Orlanski implementation (untested).

File size: 12.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 wrk_nemo        ! Memory Allocation
18   USE oce             ! ocean dynamics and tracers
19   USE dom_oce         ! ocean space and time domain
20   USE bdy_oce         ! ocean open boundary conditions
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: bdydyn.F90 2528 2010-12-27 17:33:53Z rblod $
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('zero')
62            CALL bdy_dyn3d_zro( idx_bdy(ib_bdy), dta_bdy(ib_bdy), kt, ib_bdy )
63         CASE DEFAULT
64            CALL ctl_stop( 'bdy_dyn3d : unrecognised option for open boundaries for baroclinic velocities' )
65         END SELECT
66      ENDDO
67
68   END SUBROUTINE bdy_dyn3d
69
70   SUBROUTINE bdy_dyn3d_spe( idx, dta, kt , ib_bdy )
71      !!----------------------------------------------------------------------
72      !!                  ***  SUBROUTINE bdy_dyn3d_spe  ***
73      !!
74      !! ** Purpose : - Apply a specified value for baroclinic velocities
75      !!                at open boundaries.
76      !!
77      !!----------------------------------------------------------------------
78      INTEGER                     ::   kt
79      TYPE(OBC_INDEX), INTENT(in) ::   idx  ! OBC indices
80      TYPE(OBC_DATA),  INTENT(in) ::   dta  ! OBC external data
81      INTEGER,         INTENT(in) ::   ib_bdy  ! BDY set index
82      !!
83      INTEGER  ::   jb, jk         ! dummy loop indices
84      INTEGER  ::   ii, ij, igrd   ! local integers
85      REAL(wp) ::   zwgt           ! boundary weight
86      !!----------------------------------------------------------------------
87      !
88      IF( nn_timing == 1 ) CALL timing_start('bdy_dyn3d_spe')
89      !
90      igrd = 2                      ! Relaxation of zonal velocity
91      DO jb = 1, idx%nblenrim(igrd)
92         DO jk = 1, jpkm1
93            ii   = idx%nbi(jb,igrd)
94            ij   = idx%nbj(jb,igrd)
95            ua(ii,ij,jk) = dta%u3d(jb,jk) * umask(ii,ij,jk)
96         END DO
97      END DO
98      !
99      igrd = 3                      ! Relaxation of meridional velocity
100      DO jb = 1, idx%nblenrim(igrd)
101         DO jk = 1, jpkm1
102            ii   = idx%nbi(jb,igrd)
103            ij   = idx%nbj(jb,igrd)
104            va(ii,ij,jk) = dta%v3d(jb,jk) * vmask(ii,ij,jk)
105         END DO
106      END DO
107      CALL lbc_bdy_lnk( ua, 'U', -1., ib_bdy )   ;   CALL lbc_bdy_lnk( va, 'V', -1.,ib_bdy )   ! Boundary points should be updated
108      !
109      IF( kt .eq. nit000 ) CLOSE( unit = 102 )
110
111      IF( nn_timing == 1 ) CALL timing_stop('bdy_dyn3d_spe')
112
113   END SUBROUTINE bdy_dyn3d_spe
114
115   SUBROUTINE bdy_dyn3d_zro( idx, dta, kt, ib_bdy )
116      !!----------------------------------------------------------------------
117      !!                  ***  SUBROUTINE bdy_dyn3d_zro  ***
118      !!
119      !! ** Purpose : - baroclinic velocities = 0. at open boundaries.
120      !!
121      !!----------------------------------------------------------------------
122      INTEGER                     ::   kt
123      TYPE(OBC_INDEX), INTENT(in) ::   idx  ! OBC indices
124      TYPE(OBC_DATA),  INTENT(in) ::   dta  ! OBC external data
125      INTEGER,         INTENT(in) ::   ib_bdy  ! BDY set index
126      !!
127      INTEGER  ::   ib, ik         ! dummy loop indices
128      INTEGER  ::   ii, ij, igrd, zcoef   ! local integers
129      REAL(wp) ::   zwgt           ! boundary weight
130      !!----------------------------------------------------------------------
131      !
132      IF( nn_timing == 1 ) CALL timing_start('bdy_dyn3d_zro')
133      !
134      igrd = 2                       ! Everything is at T-points here
135      DO ib = 1, idx%nblenrim(igrd)
136         ii = idx%nbi(ib,igrd)
137         ij = idx%nbj(ib,igrd)
138         DO ik = 1, jpkm1
139            ua(ii,ij,ik) = 0._wp
140         END DO
141      END DO
142
143      igrd = 3                       ! Everything is at T-points here
144      DO ib = 1, idx%nblenrim(igrd)
145         ii = idx%nbi(ib,igrd)
146         ij = idx%nbj(ib,igrd)
147         DO ik = 1, jpkm1
148            va(ii,ij,ik) = 0._wp
149         END DO
150      END DO
151      !
152      CALL lbc_bdy_lnk( ua, 'U', -1., ib_bdy )   ;   CALL lbc_bdy_lnk( va, 'V', -1.,ib_bdy )   ! Boundary points should be updated
153      !
154      IF( kt .eq. nit000 ) CLOSE( unit = 102 )
155
156      IF( nn_timing == 1 ) CALL timing_stop('bdy_dyn3d_zro')
157
158   END SUBROUTINE bdy_dyn3d_zro
159
160   SUBROUTINE bdy_dyn3d_frs( idx, dta, kt, ib_bdy )
161      !!----------------------------------------------------------------------
162      !!                  ***  SUBROUTINE bdy_dyn3d_frs  ***
163      !!
164      !! ** Purpose : - Apply the Flow Relaxation Scheme for baroclinic velocities
165      !!                at open boundaries.
166      !!
167      !! References :- Engedahl H., 1995: Use of the flow relaxation scheme in
168      !!               a three-dimensional baroclinic ocean model with realistic
169      !!               topography. Tellus, 365-382.
170      !!----------------------------------------------------------------------
171      INTEGER                     ::   kt
172      TYPE(OBC_INDEX), INTENT(in) ::   idx  ! OBC indices
173      TYPE(OBC_DATA),  INTENT(in) ::   dta  ! OBC external data
174      INTEGER,         INTENT(in) ::   ib_bdy  ! BDY set index
175      !!
176      INTEGER  ::   jb, jk         ! dummy loop indices
177      INTEGER  ::   ii, ij, igrd   ! local integers
178      REAL(wp) ::   zwgt           ! boundary weight
179      !!----------------------------------------------------------------------
180      !
181      IF( nn_timing == 1 ) CALL timing_start('bdy_dyn3d_frs')
182      !
183      igrd = 2                      ! Relaxation of zonal velocity
184      DO jb = 1, idx%nblen(igrd)
185         DO jk = 1, jpkm1
186            ii   = idx%nbi(jb,igrd)
187            ij   = idx%nbj(jb,igrd)
188            zwgt = idx%nbw(jb,igrd)
189            ua(ii,ij,jk) = ( ua(ii,ij,jk) + zwgt * ( dta%u3d(jb,jk) - ua(ii,ij,jk) ) ) * umask(ii,ij,jk)
190         END DO
191      END DO
192      !
193      igrd = 3                      ! Relaxation of meridional velocity
194      DO jb = 1, idx%nblen(igrd)
195         DO jk = 1, jpkm1
196            ii   = idx%nbi(jb,igrd)
197            ij   = idx%nbj(jb,igrd)
198            zwgt = idx%nbw(jb,igrd)
199            va(ii,ij,jk) = ( va(ii,ij,jk) + zwgt * ( dta%v3d(jb,jk) - va(ii,ij,jk) ) ) * vmask(ii,ij,jk)
200         END DO
201      END DO
202      CALL lbc_bdy_lnk( ua, 'U', -1., ib_bdy )   ;   CALL lbc_bdy_lnk( va, 'V', -1.,ib_bdy )   ! Boundary points should be updated
203      !
204      IF( kt .eq. 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   SUBROUTINE bdy_dyn3d_dmp( kt )
211      !!----------------------------------------------------------------------
212      !!                  ***  SUBROUTINE bdy_dyn3d_dmp  ***
213      !!
214      !! ** Purpose : Apply damping for baroclinic velocities at open boundaries.
215      !!
216      !!----------------------------------------------------------------------
217      INTEGER                     ::   kt
218      !!
219      INTEGER  ::   jb, jk         ! dummy loop indices
220      INTEGER  ::   ii, ij, igrd   ! local integers
221      REAL(wp) ::   zwgt           ! boundary weight
222      INTEGER  ::  ib_bdy          ! loop index
223      !!----------------------------------------------------------------------
224      !
225      IF( nn_timing == 1 ) CALL timing_start('bdy_dyn3d_dmp')
226      !
227      !-------------------------------------------------------
228      ! Remove barotropic part from before velocity
229      !-------------------------------------------------------
230      CALL wrk_alloc(jpi,jpj,pu2d,pv2d) 
231
232      pu2d(:,:) = 0.e0
233      pv2d(:,:) = 0.e0
234
235      DO jk = 1, jpkm1
236#if defined key_vvl
237         pu2d(:,:) = pu2d(:,:) + fse3u_b(:,:,jk)* ub(:,:,jk)   *umask(:,:,jk) 
238         pv2d(:,:) = pv2d(:,:) + fse3v_b(:,:,jk)* vb(:,:,jk)   *vmask(:,:,jk)
239#else
240         pu2d(:,:) = pu2d(:,:) + fse3u_0(:,:,jk) * ub(:,:,jk)  * umask(:,:,jk)
241         pv2d(:,:) = pv2d(:,:) + fse3v_0(:,:,jk) * vb(:,:,jk)  * vmask(:,:,jk)
242#endif
243      END DO
244
245      IF( lk_vvl ) THEN
246         pu2d(:,:) = pu2d(:,:) * umask(:,:,1) / ( hu_0(:,:) + sshu_b(:,:) + 1._wp - umask(:,:,1) )
247         pv2d(:,:) = pv2d(:,:) * vmask(:,:,1) / ( hv_0(:,:) + sshv_b(:,:) + 1._wp - vmask(:,:,1) )
248      ELSE
249         pu2d(:,:) = pv2d(:,:) * hur(:,:)
250         pv2d(:,:) = pu2d(:,:) * hvr(:,:)
251      ENDIF
252
253      DO ib_bdy=1, nb_bdy
254         IF ( ln_dyn3d_dmp(ib_bdy).and.nn_dyn3d(ib_bdy).gt.0 ) THEN
255            igrd = 2                      ! Relaxation of zonal velocity
256            DO jb = 1, idx_bdy(ib_bdy)%nblen(igrd)
257               ii   = idx_bdy(ib_bdy)%nbi(jb,igrd)
258               ij   = idx_bdy(ib_bdy)%nbj(jb,igrd)
259               zwgt = idx_bdy(ib_bdy)%nbd(jb,igrd)
260               DO jk = 1, jpkm1
261                  ua(ii,ij,jk) = ( ua(ii,ij,jk) + zwgt * ( dta_bdy(ib_bdy)%u3d(jb,jk) - &
262                                   ub(ii,ij,jk) + pu2d(ii,ij)) ) * umask(ii,ij,jk)
263               END DO
264            END DO
265            !
266            igrd = 3                      ! Relaxation of meridional velocity
267            DO jb = 1, idx_bdy(ib_bdy)%nblen(igrd)
268               ii   = idx_bdy(ib_bdy)%nbi(jb,igrd)
269               ij   = idx_bdy(ib_bdy)%nbj(jb,igrd)
270               zwgt = idx_bdy(ib_bdy)%nbd(jb,igrd)
271               DO jk = 1, jpkm1
272                  va(ii,ij,jk) = ( va(ii,ij,jk) + zwgt * ( dta_bdy(ib_bdy)%v3d(jb,jk) -  &
273                                   vb(ii,ij,jk) + pv2d(ii,ij)) ) * vmask(ii,ij,jk)
274               END DO
275            END DO
276         ENDIF
277      ENDDO
278      !
279      CALL wrk_dealloc(jpi,jpj,pu2d,pv2d) 
280      !
281      CALL lbc_lnk( ua, 'U', -1. )   ;   CALL lbc_lnk( va, 'V', -1. )   ! Boundary points should be updated
282      !
283      IF( nn_timing == 1 ) CALL timing_stop('bdy_dyn3d_dmp')
284
285   END SUBROUTINE bdy_dyn3d_dmp
286
287#else
288   !!----------------------------------------------------------------------
289   !!   Dummy module                   NO Unstruct Open Boundary Conditions
290   !!----------------------------------------------------------------------
291CONTAINS
292   SUBROUTINE bdy_dyn3d( kt )      ! Empty routine
293      WRITE(*,*) 'bdy_dyn3d: You should not have seen this print! error?', kt
294   END SUBROUTINE bdy_dyn3d
295
296   SUBROUTINE bdy_dyn3d_dmp( kt )      ! Empty routine
297      WRITE(*,*) 'bdy_dyn3d_dmp: You should not have seen this print! error?', kt
298   END SUBROUTINE bdy_dyn3d_dmp
299
300#endif
301
302   !!======================================================================
303END MODULE bdydyn3d
Note: See TracBrowser for help on using the repository browser.