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.
bdydyn2d.F90 in branches/2011/dev_NEMO_MERGE_2011/NEMOGCM/NEMO/OPA_SRC/BDY – NEMO

source: branches/2011/dev_NEMO_MERGE_2011/NEMOGCM/NEMO/OPA_SRC/BDY/bdydyn2d.F90 @ 3186

Last change on this file since 3186 was 3182, checked in by davestorkey, 13 years ago

Change dynamic allocation and add timing to BDY module.

File size: 9.0 KB
Line 
1MODULE bdydyn2d
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   !!            3.3  !  2010-09  (E.O'Dea) modifications for Shelf configurations
11   !!            3.3  !  2010-09  (D.Storkey) add ice boundary conditions
12   !!----------------------------------------------------------------------
13#if defined key_bdy 
14   !!----------------------------------------------------------------------
15   !!   'key_bdy' :                    Unstructured Open Boundary Condition
16   !!----------------------------------------------------------------------
17   !!   bdy_dyn2d      : Apply open boundary conditions to barotropic variables.
18   !!   bdy_dyn2d_fla    : Apply Flather condition
19   !!----------------------------------------------------------------------
20   USE timing          ! Timing
21   USE oce             ! ocean dynamics and tracers
22   USE dom_oce         ! ocean space and time domain
23   USE bdy_oce         ! ocean open boundary conditions
24   USE dynspg_oce      ! for barotropic variables
25   USE phycst          ! physical constants
26   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
27   USE in_out_manager  !
28
29   IMPLICIT NONE
30   PRIVATE
31
32   PUBLIC   bdy_dyn2d     ! routine called in dynspg_ts and dyn_nxt and dynspg_flt
33
34   !!----------------------------------------------------------------------
35   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
36   !! $Id: bdydyn.F90 2528 2010-12-27 17:33:53Z rblod $
37   !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
38   !!----------------------------------------------------------------------
39CONTAINS
40
41   SUBROUTINE bdy_dyn2d( kt )
42      !!----------------------------------------------------------------------
43      !!                  ***  SUBROUTINE bdy_dyn2d  ***
44      !!
45      !! ** Purpose : - Apply open boundary conditions for barotropic variables
46      !!
47      !!----------------------------------------------------------------------
48      INTEGER,                      INTENT(in) ::   kt   ! Main time step counter
49      !!
50      INTEGER                                  ::   ib_bdy ! Loop counter
51
52      DO ib_bdy=1, nb_bdy
53
54         SELECT CASE( nn_dyn2d(ib_bdy) )
55         CASE(jp_none)
56            CYCLE
57         CASE(jp_frs)
58            CALL bdy_dyn2d_frs( idx_bdy(ib_bdy), dta_bdy(ib_bdy) )
59         CASE(jp_flather)
60            CALL bdy_dyn2d_fla( idx_bdy(ib_bdy), dta_bdy(ib_bdy) )
61         CASE DEFAULT
62            CALL ctl_stop( 'bdy_dyn3d : unrecognised option for open boundaries for barotropic variables' )
63         END SELECT
64      ENDDO
65
66   END SUBROUTINE bdy_dyn2d
67
68   SUBROUTINE bdy_dyn2d_frs( idx, dta )
69      !!----------------------------------------------------------------------
70      !!                  ***  SUBROUTINE bdy_dyn2d_frs  ***
71      !!
72      !! ** Purpose : - Apply the Flow Relaxation Scheme for barotropic velocities
73      !!                at open boundaries.
74      !!
75      !! References :- Engedahl H., 1995: Use of the flow relaxation scheme in
76      !!               a three-dimensional baroclinic ocean model with realistic
77      !!               topography. Tellus, 365-382.
78      !!----------------------------------------------------------------------
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_dyn2d_frs')
88      !
89      igrd = 2                      ! Relaxation of zonal velocity
90      DO jb = 1, idx%nblen(igrd)
91         ii   = idx%nbi(jb,igrd)
92         ij   = idx%nbj(jb,igrd)
93         zwgt = idx%nbw(jb,igrd)
94         pu2d(ii,ij) = ( pu2d(ii,ij) + zwgt * ( dta%u2d(jb) - pu2d(ii,ij) ) ) * umask(ii,ij,1)
95      END DO
96      !
97      igrd = 3                      ! Relaxation of meridional velocity
98      DO jb = 1, idx%nblen(igrd)
99         ii   = idx%nbi(jb,igrd)
100         ij   = idx%nbj(jb,igrd)
101         zwgt = idx%nbw(jb,igrd)
102         pv2d(ii,ij) = ( pv2d(ii,ij) + zwgt * ( dta%v2d(jb) - pv2d(ii,ij) ) ) * vmask(ii,ij,1)
103      END DO
104      CALL lbc_lnk( pu2d, 'U', -1. ) 
105      CALL lbc_lnk( pv2d, 'V', -1. )   ! Boundary points should be updated
106      !
107      IF( nn_timing == 1 ) CALL timing_stop('bdy_dyn2d_frs')
108      !
109
110   END SUBROUTINE bdy_dyn2d_frs
111
112
113   SUBROUTINE bdy_dyn2d_fla( idx, dta )
114      !!----------------------------------------------------------------------
115      !!                 ***  SUBROUTINE bdy_dyn2d_fla  ***
116      !!             
117      !!              - Apply Flather boundary conditions on normal barotropic velocities
118      !!
119      !! ** WARNINGS about FLATHER implementation:
120      !!1. According to Palma and Matano, 1998 "after ssh" is used.
121      !!   In ROMS and POM implementations, it is "now ssh". In the current
122      !!   implementation (tested only in the EEL-R5 conf.), both cases were unstable.
123      !!   So I use "before ssh" in the following.
124      !!
125      !!2. We assume that the normal ssh gradient at the bdy is zero. As a matter of
126      !!   fact, the model ssh just inside the dynamical boundary is used (the outside 
127      !!   ssh in the code is not updated).
128      !!
129      !! References:  Flather, R. A., 1976: A tidal model of the northwest European
130      !!              continental shelf. Mem. Soc. R. Sci. Liege, Ser. 6,10, 141-164.     
131      !!----------------------------------------------------------------------
132      TYPE(OBC_INDEX),              INTENT(in) ::   idx  ! OBC indices
133      TYPE(OBC_DATA),               INTENT(in) ::   dta  ! OBC external data
134
135      INTEGER  ::   jb, igrd                         ! dummy loop indices
136      INTEGER  ::   ii, ij, iim1, iip1, ijm1, ijp1   ! 2D addresses
137      REAL(wp) ::   zcorr                            ! Flather correction
138      REAL(wp) ::   zforc                            ! temporary scalar
139      !!----------------------------------------------------------------------
140
141      IF( nn_timing == 1 ) CALL timing_start('bdy_dyn2d_fla')
142
143      ! ---------------------------------!
144      ! Flather boundary conditions     :!
145      ! ---------------------------------!
146     
147!!! REPLACE spgu with nemo_wrk work space
148
149      ! Fill temporary array with ssh data (here spgu):
150      igrd = 1
151      spgu(:,:) = 0.0
152      DO jb = 1, idx%nblenrim(igrd)
153         ii = idx%nbi(jb,igrd)
154         ij = idx%nbj(jb,igrd)
155         spgu(ii, ij) = dta%ssh(jb)
156      END DO
157      !
158      igrd = 2      ! Flather bc on u-velocity;
159      !             ! remember that flagu=-1 if normal velocity direction is outward
160      !             ! I think we should rather use after ssh ?
161      DO jb = 1, idx%nblenrim(igrd)
162         ii  = idx%nbi(jb,igrd)
163         ij  = idx%nbj(jb,igrd) 
164         iim1 = ii + MAX( 0, INT( idx%flagu(jb) ) )   ! T pts i-indice inside the boundary
165         iip1 = ii - MIN( 0, INT( idx%flagu(jb) ) )   ! T pts i-indice outside the boundary
166         !
167         zcorr = - idx%flagu(jb) * SQRT( grav * phur(ii, ij) ) * ( pssh(iim1, ij) - spgu(iip1,ij) )
168         zforc = dta%u2d(jb)
169         pu2d(ii,ij) = zforc + zcorr * umask(ii,ij,1) 
170      END DO
171      !
172      igrd = 3      ! Flather bc on v-velocity
173      !             ! remember that flagv=-1 if normal velocity direction is outward
174      DO jb = 1, idx%nblenrim(igrd)
175         ii  = idx%nbi(jb,igrd)
176         ij  = idx%nbj(jb,igrd) 
177         ijm1 = ij + MAX( 0, INT( idx%flagv(jb) ) )   ! T pts j-indice inside the boundary
178         ijp1 = ij - MIN( 0, INT( idx%flagv(jb) ) )   ! T pts j-indice outside the boundary
179         !
180         zcorr = - idx%flagv(jb) * SQRT( grav * phvr(ii, ij) ) * ( pssh(ii, ijm1) - spgu(ii,ijp1) )
181         zforc = dta%v2d(jb)
182         pv2d(ii,ij) = zforc + zcorr * vmask(ii,ij,1)
183      END DO
184      CALL lbc_lnk( pu2d, 'U', -1. )   ! Boundary points should be updated
185      CALL lbc_lnk( pv2d, 'V', -1. )   !
186      !
187      IF( nn_timing == 1 ) CALL timing_stop('bdy_dyn2d_fla')
188      !
189   END SUBROUTINE bdy_dyn2d_fla
190#else
191   !!----------------------------------------------------------------------
192   !!   Dummy module                   NO Unstruct Open Boundary Conditions
193   !!----------------------------------------------------------------------
194CONTAINS
195   SUBROUTINE bdy_dyn2d( kt )      ! Empty routine
196      WRITE(*,*) 'bdy_dyn_frs: You should not have seen this print! error?', kt
197   END SUBROUTINE bdy_dyn2d
198#endif
199
200   !!======================================================================
201END MODULE bdydyn2d
Note: See TracBrowser for help on using the repository browser.