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 @ 911

Last change on this file since 911 was 911, checked in by ctlod, 16 years ago

Implementation of the BDY package, see ticket: #126

  • Property svn:executable set to *
File size: 8.4 KB
Line 
1MODULE bdydyn
2   !!=================================================================================
3   !!                       ***  MODULE  bdydyn  ***
4   !! Ocean dynamics:   Flow relaxation scheme of velocities on unstruc. open boundary
5   !!=================================================================================
6
7   !!---------------------------------------------------------------------------------
8   !!   bdy_dyn_frs    : relaxation of velocities on unstructured open boundary
9   !!   bdy_dyn_fla    : Flather condition for barotropic solution
10   !!---------------------------------------------------------------------------------
11#if defined key_bdy || defined key_bdy_tides
12   !!---------------------------------------------------------------------------------
13   !! * Modules used
14   USE oce             ! ocean dynamics and tracers
15   USE dom_oce         ! ocean space and time domain
16   USE bdy_oce         ! ocean open boundary conditions
17   USE dynspg_oce      ! for barotropic variables
18   USE phycst          ! physical constants
19   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
20   USE bdytides        ! for tidal harmonic forcing at boundary
21   USE in_out_manager
22
23   IMPLICIT NONE
24   PRIVATE
25
26   !! * Accessibility
27   PUBLIC bdy_dyn_frs  ! routine called in dynspg_flt (free surface case ONLY)
28#if defined key_dynspg_exp || defined key_dynspg_ts
29   PUBLIC bdy_dyn_fla  ! routine called in dynspg_flt (free surface case ONLY)
30#endif
31
32   !!---------------------------------------------------------------------------------
33
34CONTAINS
35
36   SUBROUTINE bdy_dyn_frs ( kt )
37      !!------------------------------------------------------------------------------
38      !!                      SUBROUTINE bdy_dyn_frs
39      !!                     ***********************
40      !! ** Purpose : - Apply the Flow Relaxation Scheme for dynamic in the 
41      !!                case of unstructured open boundaries.
42      !!
43      !! References :- Engedahl H., 1995: Use of the flow relaxation scheme in
44      !!               a three-dimensional baroclinic ocean model with realistic
45      !!               topography. Tellus, 365-382.
46      !!
47      !! History :
48      !!   9.0  !  05-02 (J. Chanut, A. Sellar) Original
49      !!        !  07-07 (D. Storkey) Move Flather implementation to separate routine.
50      !!------------------------------------------------------------------------------
51      !! * Arguments
52      INTEGER, INTENT( in ) ::   kt    ! Main time step counter
53
54      !! * Local declarations
55      INTEGER ::   jb, jk, jgrd        ! dummy loop indices
56      INTEGER ::   ii, ij              ! 2D addresses
57      REAL(wp) :: zwgt                 ! boundary weight
58
59      !!------------------------------------------------------------------------------
60      !!  OPA 9.0, LODYC-IPSL (2003)
61      !!------------------------------------------------------------------------------
62
63      ! ---------------------------!
64      ! FRS on the total velocity :!
65      ! ---------------------------! 
66 
67      jgrd=2 !: Relaxation of zonal velocity
68      DO jb = 1, nblen(jgrd)
69        DO jk = 1, jpkm1
70          ii = nbi(jb,jgrd)
71          ij = nbj(jb,jgrd)
72          zwgt = nbw(jb,jgrd)
73
74          ua(ii,ij,jk) = ( ua(ii,ij,jk)*(1.-zwgt) + ubdy(jb,jk)*zwgt ) &
75                                                        * umask(ii,ij,jk)
76        END DO
77      END DO
78
79      jgrd=3 !: Relaxation of meridional velocity
80      DO jb = 1, nblen(jgrd)
81        DO jk = 1, jpkm1
82          ii = nbi(jb,jgrd)
83          ij = nbj(jb,jgrd)
84          zwgt = nbw(jb,jgrd)
85
86          va(ii,ij,jk) = ( va(ii,ij,jk)*(1.-zwgt) + vbdy(jb,jk)*zwgt ) &
87                                                        * vmask(ii,ij,jk)
88        END DO
89      END DO
90
91      CALL lbc_lnk( ua, 'U', 1. ) ! Boundary points should be updated
92      CALL lbc_lnk( va, 'V', 1. ) !
93
94   END SUBROUTINE bdy_dyn_frs
95
96#if defined key_dynspg_exp || defined key_dynspg_ts
97!! Option to use Flather with dynspg_flt not coded yet...
98   SUBROUTINE bdy_dyn_fla
99      !!------------------------------------------------------------------------------
100      !!                      SUBROUTINE bdy_dyn_fla
101      !!                     ***********************
102      !!              - Apply Flather boundary conditions on normal barotropic velocities
103      !!                (ln_bdy_fla=.true.)
104      !!
105      !! ** WARNINGS about FLATHER implementation:
106      !!1. According to Palma and Matano, 1998 "after ssh" is used.
107      !!   In ROMS and POM implementations, it is "now ssh". In the current
108      !!   implementation (tested only in the EEL-R5 conf.), both cases were unstable.
109      !!   So I use "before ssh" in the following.
110      !!
111      !!2. We assume that the normal ssh gradient at the bdy is zero. As a matter of
112      !!   fact, the model ssh just inside the dynamical boundary is used (the outside 
113      !!   ssh in the code is not updated).
114      !!
115      !!             - Flather, R. A., 1976: A tidal model of the northwest European
116      !!               continental shelf. Mem. Soc. R. Sci. Liege, Ser. 6,10, 141-164.     
117      !! History :
118      !!   9.0  !  05-02 (J. Chanut, A. Sellar) Original
119      !!        !  07-07 (D. Storkey) Flather algorithm in separate routine.
120      !!------------------------------------------------------------------------------
121      !! * Local declarations
122      INTEGER ::   jb, jk, jgrd, ji, jj, jim1, jip1, jjm1, jjp1   ! dummy loop indices
123      INTEGER ::   ii, ij                    ! 2D addresses
124      REAL(wp) ::  corr ! Flather correction
125      REAL(wp) :: zwgt                       ! boundary weight
126      REAL(wp) :: elapsed
127
128      !!------------------------------------------------------------------------------
129      !!  OPA 9.0, LODYC-IPSL (2003)
130      !!------------------------------------------------------------------------------
131
132      ! ---------------------------------!
133      ! Flather boundary conditions     :!
134      ! ---------------------------------!
135     
136      ! Fill temporary array with ssh data (here spgu):
137      jgrd = 1
138      DO jb = 1, nblenrim(jgrd)
139        ji = nbi(jb,jgrd)
140        jj = nbj(jb,jgrd)
141        spgu(ji, jj) = sshbdy(jb)
142#if defined key_bdy_tides
143        spgu(ji, jj) = spgu(ji, jj) + sshtide(jb)
144#endif
145      END DO
146
147      jgrd = 2  !: Flather bc on u-velocity;
148                ! remember that flagu=-1 if normal velocity direction is outward
149                ! I think we should rather use after ssh ?
150
151      DO jb = 1, nblenrim(jgrd)
152        ji  = nbi(jb,jgrd)
153        jj  = nbj(jb,jgrd) 
154        jim1 = ji+MAX(0, INT(flagu(jb))) ! T pts i-indice inside the boundary
155        jip1 = ji-MIN(0, INT(flagu(jb))) ! T pts i-indice outside the boundary
156     
157        corr = - flagu(jb) * sqrt (grav / (hu_e(ji, jj) + 1.e-20) )  &
158                                      * ( sshn_e(jim1, jj) - spgu(jip1,jj) )
159        ua_e(ji, jj) = ( ubtbdy(jb) + utide(jb) ) * hu_e(ji,jj)
160        if ( ln_bdy_fla ) then
161          ua_e(ji,jj) = ua_e(ji,jj) + corr * umask(ji,jj,1) * hu_e(ji,jj)
162        endif
163
164      END DO
165       
166      jgrd = 3  !: Flather bc on v-velocity
167                ! remember that flagv=-1 if normal velocity direction is outward
168                 
169      DO jb = 1, nblenrim(jgrd)
170        ji  = nbi(jb,jgrd)
171        jj  = nbj(jb,jgrd) 
172        jjm1 = jj+MAX(0, INT(flagv(jb))) ! T pts j-indice inside the boundary
173        jjp1 = jj-MIN(0, INT(flagv(jb))) ! T pts j-indice outside the boundary
174     
175        corr = - flagv(jb) * sqrt (grav / (hv_e(ji, jj) + 1.e-20) )  &
176                                      * ( sshn_e(ji, jjm1) - spgu(ji,jjp1) )
177        va_e(ji, jj) = ( vbtbdy(jb) + vtide(jb) ) * hv_e(ji,jj)
178        if ( ln_bdy_fla ) then
179          va_e(ji,jj) = va_e(ji,jj) + corr * vmask(ji,jj,1) * hv_e(ji,jj)
180        endif
181
182      END DO
183
184      CALL lbc_lnk( ua_e, 'U', 1. ) ! Boundary points should be updated
185      CALL lbc_lnk( va_e, 'V', 1. ) !
186
187   END SUBROUTINE bdy_dyn_fla
188#endif
189
190#else
191   !!=================================================================================
192   !!                       ***  MODULE  bdydyn  ***
193   !! Ocean dynamics:   Flow relaxation scheme of velocities on unstruc. open boundary
194   !!=================================================================================
195CONTAINS
196
197   SUBROUTINE bdy_dyn_frs
198                              ! No Unstructured open boundaries ==> empty routine
199   END SUBROUTINE bdy_dyn_frs
200
201   SUBROUTINE bdy_dyn_fla
202                              ! No Unstructured open boundaries ==> empty routine
203   END SUBROUTINE bdy_dyn_fla
204
205#endif
206
207END MODULE bdydyn
Note: See TracBrowser for help on using the repository browser.