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

source: branches/UKMO/r6232_tracer_advection/NEMOGCM/NEMO/OPA_SRC/BDY/bdydyn3d.F90 @ 9295

Last change on this file since 9295 was 9295, checked in by jcastill, 6 years ago

Remove svn keywords

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