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.
bdydyn.F90 in trunk/NEMO/OPA_SRC/BDY – NEMO

source: trunk/NEMO/OPA_SRC/BDY/bdydyn.F90 @ 1502

Last change on this file since 1502 was 1502, checked in by rblod, 15 years ago

Update dynnxt and dynspg_ts for variable volume, see ticket #474

  • Property svn:executable set to *
  • Property svn:keywords set to Id
File size: 8.7 KB
Line 
1MODULE bdydyn
2   !!======================================================================
3   !!                       ***  MODULE  bdydyn  ***
4   !! Unstructured Open Boundary Cond. :   Flow relaxation scheme on velocities
5   !!======================================================================
6   !! History :  1.0  !  2005-02  (J. Chanut, A. Sellar)  Original code
7   !!             -   !  2007-07  (D. Storkey) Move Flather implementation to separate routine.
8   !!            3.0  !  2008-04  (NEMO team)  add in the reference version
9   !!            3.2  !  2008-04  (R. Benshila) consider velocity instead of transport
10   !!----------------------------------------------------------------------
11#if defined key_bdy 
12   !!----------------------------------------------------------------------
13   !!   'key_bdy' :                    Unstructured Open Boundary Condition
14   !!----------------------------------------------------------------------
15   !!   bdy_dyn_frs    : relaxation of velocities on unstructured open boundary
16   !!   bdy_dyn_fla    : Flather condition for barotropic solution
17   !!----------------------------------------------------------------------
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 dynspg_oce      ! for barotropic variables
22   USE phycst          ! physical constants
23   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
24   USE bdytides        ! for tidal harmonic forcing at boundary
25   USE in_out_manager  !
26
27   IMPLICIT NONE
28   PRIVATE
29
30   PUBLIC   bdy_dyn_frs   ! routine called in dynspg_flt (free surface case ONLY)
31# if defined key_dynspg_exp || defined key_dynspg_ts
32   PUBLIC   bdy_dyn_fla   ! routine called in dynspg_flt (free surface case ONLY)
33# endif
34
35   !!----------------------------------------------------------------------
36   !! NEMO/OPA 3.0 , LOCEAN-IPSL (2008)
37   !! $Id$
38   !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt)
39   !!----------------------------------------------------------------------
40
41CONTAINS
42
43   SUBROUTINE bdy_dyn_frs( kt )
44      !!----------------------------------------------------------------------
45      !!                  ***  SUBROUTINE bdy_dyn_frs  ***
46      !!
47      !! ** Purpose : - Apply the Flow Relaxation Scheme for dynamic in the 
48      !!                case of unstructured open boundaries.
49      !!
50      !! References :- Engedahl H., 1995: Use of the flow relaxation scheme in
51      !!               a three-dimensional baroclinic ocean model with realistic
52      !!               topography. Tellus, 365-382.
53      !!----------------------------------------------------------------------
54      INTEGER, INTENT( in ) ::   kt   ! Main time step counter
55      !!
56      INTEGER  ::   ib, ik, igrd      ! dummy loop indices
57      INTEGER  ::   ii, ij            ! 2D addresses
58      REAL(wp) ::   zwgt              ! boundary weight
59      !!----------------------------------------------------------------------
60      !
61      IF(ln_bdy_dyn_frs) THEN ! If this is false, then this routine does nothing.
62
63         IF( kt == nit000 ) THEN
64            IF(lwp) WRITE(numout,*)
65            IF(lwp) WRITE(numout,*) 'bdy_dyn : Flow Relaxation Scheme on momentum'
66            IF(lwp) WRITE(numout,*) '~~~~~~~'
67         ENDIF
68         !
69         igrd = 2                      ! Relaxation of zonal velocity
70         DO ib = 1, nblen(igrd)
71            DO ik = 1, jpkm1
72               ii = nbi(ib,igrd)
73               ij = nbj(ib,igrd)
74               zwgt = nbw(ib,igrd)
75               ua(ii,ij,ik) = ( ua(ii,ij,ik) * ( 1.- zwgt ) + ubdy(ib,ik) * zwgt ) * umask(ii,ij,ik)
76            END DO
77         END DO
78         !
79         igrd = 3                      ! Relaxation of meridional velocity
80         DO ib = 1, nblen(igrd)
81            DO ik = 1, jpkm1
82               ii = nbi(ib,igrd)
83               ij = nbj(ib,igrd)
84               zwgt = nbw(ib,igrd)
85               va(ii,ij,ik) = ( va(ii,ij,ik) * ( 1.- zwgt ) + vbdy(ib,ik) * zwgt ) * vmask(ii,ij,ik)
86            END DO
87         END DO 
88         !
89         CALL lbc_lnk( ua, 'U', -1. )   ! Boundary points should be updated
90         CALL lbc_lnk( va, 'V', -1. )   !
91         !
92      ENDIF ! ln_bdy_dyn_frs
93
94   END SUBROUTINE bdy_dyn_frs
95
96
97#if defined key_dynspg_exp || defined key_dynspg_ts
98!! Option to use Flather with dynspg_flt not coded yet...
99   SUBROUTINE bdy_dyn_fla( pssh )
100      !!----------------------------------------------------------------------
101      !!                 ***  SUBROUTINE bdy_dyn_fla  ***
102      !!             
103      !!              - Apply Flather boundary conditions on normal barotropic velocities
104      !!                (ln_bdy_dyn_fla=.true. or ln_bdy_tides=.true.)
105      !!
106      !! ** WARNINGS about FLATHER implementation:
107      !!1. According to Palma and Matano, 1998 "after ssh" is used.
108      !!   In ROMS and POM implementations, it is "now ssh". In the current
109      !!   implementation (tested only in the EEL-R5 conf.), both cases were unstable.
110      !!   So I use "before ssh" in the following.
111      !!
112      !!2. We assume that the normal ssh gradient at the bdy is zero. As a matter of
113      !!   fact, the model ssh just inside the dynamical boundary is used (the outside 
114      !!   ssh in the code is not updated).
115      !!
116      !! References:  Flather, R. A., 1976: A tidal model of the northwest European
117      !!              continental shelf. Mem. Soc. R. Sci. Liege, Ser. 6,10, 141-164.     
118      !!----------------------------------------------------------------------
119      REAL(wp), DIMENSION(jpi,jpj), INTENT(in) ::   pssh
120
121      INTEGER  ::   ib, igrd                         ! dummy loop indices
122      INTEGER  ::   ii, ij, iim1, iip1, ijm1, ijp1   ! 2D addresses
123      REAL(wp) ::   zcorr                            ! Flather correction
124      REAL(wp) ::   zforc                            ! temporary scalar
125      !!----------------------------------------------------------------------
126
127      ! ---------------------------------!
128      ! Flather boundary conditions     :!
129      ! ---------------------------------!
130     
131      IF(ln_bdy_dyn_fla .OR. ln_bdy_tides) THEN ! If these are both false, then this routine does nothing.
132
133         ! Fill temporary array with ssh data (here spgu):
134         igrd = 1
135         spgu(:,:) = 0.0
136         DO ib = 1, nblenrim(igrd)
137            ii = nbi(ib,igrd)
138            ij = nbj(ib,igrd)
139            IF( ln_bdy_dyn_fla ) spgu(ii, ij) = sshbdy(ib)
140            IF( ln_bdy_tides )   spgu(ii, ij) = spgu(ii, ij) + sshtide(ib)
141         END DO
142         !
143         igrd = 2      ! Flather bc on u-velocity;
144         !             ! remember that flagu=-1 if normal velocity direction is outward
145         !             ! I think we should rather use after ssh ?
146         DO ib = 1, nblenrim(igrd)
147            ii  = nbi(ib,igrd)
148            ij  = nbj(ib,igrd) 
149            iim1 = ii + MAX( 0, INT( flagu(ib) ) )   ! T pts i-indice inside the boundary
150            iip1 = ii - MIN( 0, INT( flagu(ib) ) )   ! T pts i-indice outside the boundary
151            !
152            zcorr = - flagu(ib) * SQRT( grav * hur_e(ii, ij) ) * ( pssh(iim1, ij) - spgu(iip1,ij) )
153            zforc = ubtbdy(ib) + utide(ib)
154            ua_e(ii,ij) = zforc + zcorr * umask(ii,ij,1) 
155         END DO
156         !
157         igrd = 3      ! Flather bc on v-velocity
158         !             ! remember that flagv=-1 if normal velocity direction is outward
159         DO ib = 1, nblenrim(igrd)
160            ii  = nbi(ib,igrd)
161            ij  = nbj(ib,igrd) 
162            ijm1 = ij + MAX( 0, INT( flagv(ib) ) )   ! T pts j-indice inside the boundary
163            ijp1 = ij - MIN( 0, INT( flagv(ib) ) )   ! T pts j-indice outside the boundary
164            !
165            zcorr = - flagv(ib) * SQRT( grav * hvr_e(ii, ij) ) * ( pssh(ii, ijm1) - spgu(ii,ijp1) )
166            zforc = vbtbdy(ib) + vtide(ib)
167            va_e(ii,ij) = zforc + zcorr * vmask(ii,ij,1)
168         END DO
169         CALL lbc_lnk( ua_e, 'U', 1. )   ! Boundary points should be updated
170         CALL lbc_lnk( va_e, 'V', 1. )   !
171         !
172      ENDIF ! ln_bdy_dyn_fla .or. ln_bdy_tides
173      !
174   END SUBROUTINE bdy_dyn_fla
175#endif
176
177#else
178   !!----------------------------------------------------------------------
179   !!   Dummy module                   NO Unstruct Open Boundary Conditions
180   !!----------------------------------------------------------------------
181CONTAINS
182   SUBROUTINE bdy_dyn_frs( kt )      ! Empty routine
183      WRITE(*,*) 'bdy_dyn_frs: You should not have seen this print! error?', kt
184   END SUBROUTINE bdy_dyn_frs
185   SUBROUTINE bdy_dyn_fla            ! Empty routine
186      WRITE(*,*) 'bdy_dyn_fla: You should not have seen this print! error?'
187   END SUBROUTINE bdy_dyn_fla
188#endif
189
190   !!======================================================================
191END MODULE bdydyn
Note: See TracBrowser for help on using the repository browser.