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

source: branches/2012/dev_r3327_MERCATOR1_BDY/NEMOGCM/NEMO/OPA_SRC/BDY/bdydyn3d.F90 @ 3367

Last change on this file since 3367 was 3367, checked in by greffray, 12 years ago

Merge tidal packages
Merge OBC-BDY packages
Add atmospheric pressure at the open boundaries
Modification of the namelist for AMM12 configuration

File size: 10.8 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   !!----------------------------------------------------------------------
8#if defined key_bdy 
9   !!----------------------------------------------------------------------
10   !!   'key_bdy' :                    Unstructured Open Boundary Condition
11   !!----------------------------------------------------------------------
12   !!   bdy_dyn3d        : apply open boundary conditions to baroclinic velocities
13   !!   bdy_dyn3d_frs    : apply Flow Relaxation Scheme
14   !!----------------------------------------------------------------------
15   USE timing          ! Timing
16   USE oce             ! ocean dynamics and tracers
17   USE dom_oce         ! ocean space and time domain
18   USE bdy_oce         ! ocean open boundary conditions
19   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
20   USE in_out_manager  !
21   Use phycst
22
23   IMPLICIT NONE
24   PRIVATE
25
26   PUBLIC   bdy_dyn3d     ! routine called by bdy_dyn
27   PUBLIC   bdy_dyn3d_dmp ! routine called by step
28
29   !!----------------------------------------------------------------------
30   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
31   !! $Id: bdydyn.F90 2528 2010-12-27 17:33:53Z rblod $
32   !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
33   !!----------------------------------------------------------------------
34CONTAINS
35
36   SUBROUTINE bdy_dyn3d( kt )
37      !!----------------------------------------------------------------------
38      !!                  ***  SUBROUTINE bdy_dyn3d  ***
39      !!
40      !! ** Purpose : - Apply open boundary conditions for baroclinic velocities
41      !!
42      !!----------------------------------------------------------------------
43      INTEGER, INTENT( in ) :: kt     ! Main time step counter
44      !!
45      INTEGER               :: ib_bdy ! loop index
46      !!
47
48      DO ib_bdy=1, nb_bdy
49
50!!$         IF ( using Orlanski radiation conditions ) THEN
51!!$            CALL bdy_rad( kt,  bdyidx(ib_bdy) )
52!!$         ENDIF
53
54         SELECT CASE( nn_dyn3d(ib_bdy) )
55         CASE(jp_none)
56            CYCLE
57         CASE(jp_frs)
58            CALL bdy_dyn3d_frs( idx_bdy(ib_bdy), dta_bdy(ib_bdy), kt )
59         CASE(2)
60            CALL bdy_dyn3d_spe( idx_bdy(ib_bdy), dta_bdy(ib_bdy), kt )
61         CASE(3)
62            CALL bdy_dyn3d_zro( idx_bdy(ib_bdy), dta_bdy(ib_bdy), kt )
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 )
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      !!
82      INTEGER  ::   jb, jk         ! dummy loop indices
83      INTEGER  ::   ii, ij, igrd   ! local integers
84      REAL(wp) ::   zwgt           ! boundary weight
85      !!----------------------------------------------------------------------
86      !
87      IF( nn_timing == 1 ) CALL timing_start('bdy_dyn3d_spe')
88      !
89      igrd = 2                      ! Relaxation of zonal velocity
90      DO jb = 1, idx%nblenrim(igrd)
91         DO jk = 1, jpkm1
92            ii   = idx%nbi(jb,igrd)
93            ij   = idx%nbj(jb,igrd)
94            ua(ii,ij,jk) = dta%u3d(jb,jk) * umask(ii,ij,jk)
95         END DO
96      END DO
97      !
98      igrd = 3                      ! Relaxation of meridional 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            va(ii,ij,jk) = dta%v3d(jb,jk) * vmask(ii,ij,jk)
104         END DO
105      END DO
106      CALL lbc_lnk( ua, 'U', -1. )   ;   CALL lbc_lnk( va, 'V', -1. )   ! Boundary points should be updated
107      !
108      IF( kt .eq. nit000 ) CLOSE( unit = 102 )
109
110      IF( nn_timing == 1 ) CALL timing_stop('bdy_dyn3d_spe')
111
112   END SUBROUTINE bdy_dyn3d_spe
113
114   SUBROUTINE bdy_dyn3d_zro( idx, dta, kt )
115      !!----------------------------------------------------------------------
116      !!                  ***  SUBROUTINE bdy_dyn3d_zro  ***
117      !!
118      !! ** Purpose : - baroclinic velocities = 0. at open boundaries.
119      !!
120      !!----------------------------------------------------------------------
121      INTEGER                     ::   kt
122      TYPE(OBC_INDEX), INTENT(in) ::   idx  ! OBC indices
123      TYPE(OBC_DATA),  INTENT(in) ::   dta  ! OBC external data
124      !!
125      INTEGER  ::   ib, ik         ! dummy loop indices
126      INTEGER  ::   ii, ij, igrd, zcoef   ! 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_lnk( ua, 'U', -1. )   ;   CALL lbc_lnk( va, 'V', -1. )   ! Boundary points should be updated
151      !
152      IF( kt .eq. 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   SUBROUTINE bdy_dyn3d_frs( idx, dta, kt )
159      !!----------------------------------------------------------------------
160      !!                  ***  SUBROUTINE bdy_dyn3d_frs  ***
161      !!
162      !! ** Purpose : - Apply the Flow Relaxation Scheme for baroclinic velocities
163      !!                at open boundaries.
164      !!
165      !! References :- Engedahl H., 1995: Use of the flow relaxation scheme in
166      !!               a three-dimensional baroclinic ocean model with realistic
167      !!               topography. Tellus, 365-382.
168      !!----------------------------------------------------------------------
169      INTEGER                     ::   kt
170      TYPE(OBC_INDEX), INTENT(in) ::   idx  ! OBC indices
171      TYPE(OBC_DATA),  INTENT(in) ::   dta  ! OBC external data
172      !!
173      INTEGER  ::   jb, jk         ! dummy loop indices
174      INTEGER  ::   ii, ij, igrd   ! local integers
175      REAL(wp) ::   zwgt           ! boundary weight
176      !!----------------------------------------------------------------------
177      !
178      IF( nn_timing == 1 ) CALL timing_start('bdy_dyn3d_frs')
179      !
180      igrd = 2                      ! Relaxation of zonal velocity
181      DO jb = 1, idx%nblen(igrd)
182         DO jk = 1, jpkm1
183            ii   = idx%nbi(jb,igrd)
184            ij   = idx%nbj(jb,igrd)
185            zwgt = idx%nbw(jb,igrd)
186            ua(ii,ij,jk) = ( ua(ii,ij,jk) + zwgt * ( dta%u3d(jb,jk) - ua(ii,ij,jk) ) ) * umask(ii,ij,jk)
187         END DO
188      END DO
189      !
190      igrd = 3                      ! Relaxation of meridional velocity
191      DO jb = 1, idx%nblen(igrd)
192         DO jk = 1, jpkm1
193            ii   = idx%nbi(jb,igrd)
194            ij   = idx%nbj(jb,igrd)
195            zwgt = idx%nbw(jb,igrd)
196            va(ii,ij,jk) = ( va(ii,ij,jk) + zwgt * ( dta%v3d(jb,jk) - va(ii,ij,jk) ) ) * vmask(ii,ij,jk)
197         END DO
198      END DO
199      CALL lbc_lnk( ua, 'U', -1. )   ;   CALL lbc_lnk( va, 'V', -1. )   ! Boundary points should be updated
200      !
201      IF( kt .eq. nit000 ) CLOSE( unit = 102 )
202
203      IF( nn_timing == 1 ) CALL timing_stop('bdy_dyn3d_frs')
204
205   END SUBROUTINE bdy_dyn3d_frs
206
207   SUBROUTINE bdy_dyn3d_dmp( kt )
208      !!----------------------------------------------------------------------
209      !!                  ***  SUBROUTINE bdy_dyn3d_dmp  ***
210      !!
211      !! ** Purpose : Apply damping for baroclinic velocities at open boundaries.
212      !!
213      !!----------------------------------------------------------------------
214      INTEGER                     ::   kt
215      !!
216      INTEGER  ::   jb, jk         ! dummy loop indices
217      INTEGER  ::   ii, ij, igrd   ! local integers
218      REAL(wp) ::   zwgt           ! boundary weight
219      INTEGER  ::  ib_bdy          ! loop index
220      !!----------------------------------------------------------------------
221      !
222      IF( nn_timing == 1 ) CALL timing_start('bdy_dyn3d_dmp')
223      !
224      DO ib_bdy=1, nb_bdy
225         IF ( ln_dyn3d_dmp(ib_bdy).and.nn_dyn3d(ib_bdy).gt.0 ) THEN
226            igrd = 2                      ! Relaxation of zonal velocity
227            DO jb = 1, idx_bdy(ib_bdy)%nblen(igrd)
228               ii   = idx_bdy(ib_bdy)%nbi(jb,igrd)
229               ij   = idx_bdy(ib_bdy)%nbj(jb,igrd)
230               zwgt = idx_bdy(ib_bdy)%nbw(jb,igrd) / ( rn_time_dmp(ib_bdy) * rday )
231               DO jk = 1, jpkm1
232                  ua(ii,ij,jk) = ( ua(ii,ij,jk) + zwgt * ( dta_bdy(ib_bdy)%u3d(jb,jk) - ub(ii,ij,jk) ) ) * umask(ii,ij,jk)
233               END DO
234            END DO
235            !
236            igrd = 3                      ! Relaxation of meridional velocity
237            DO jb = 1, idx_bdy(ib_bdy)%nblen(igrd)
238               ii   = idx_bdy(ib_bdy)%nbi(jb,igrd)
239               ij   = idx_bdy(ib_bdy)%nbj(jb,igrd)
240               zwgt = idx_bdy(ib_bdy)%nbw(jb,igrd) / ( rn_time_dmp(ib_bdy) * rday )
241               DO jk = 1, jpkm1
242                  va(ii,ij,jk) = ( va(ii,ij,jk) + zwgt * ( dta_bdy(ib_bdy)%v3d(jb,jk) - vb(ii,ij,jk) ) ) * vmask(ii,ij,jk)
243               END DO
244            END DO
245         ENDIF
246      ENDDO
247      !
248      CALL lbc_lnk( ua, 'U', -1. )   ;   CALL lbc_lnk( va, 'V', -1. )   ! Boundary points should be updated
249      !
250      IF( nn_timing == 1 ) CALL timing_stop('bdy_dyn3d_dmp')
251
252   END SUBROUTINE bdy_dyn3d_dmp
253
254#else
255   !!----------------------------------------------------------------------
256   !!   Dummy module                   NO Unstruct Open Boundary Conditions
257   !!----------------------------------------------------------------------
258CONTAINS
259   SUBROUTINE bdy_dyn3d( kt )      ! Empty routine
260      WRITE(*,*) 'bdy_dyn_frs: You should not have seen this print! error?', kt
261   END SUBROUTINE bdy_dyn3d
262#endif
263
264   !!======================================================================
265END MODULE bdydyn3d
Note: See TracBrowser for help on using the repository browser.