[911] | 1 | MODULE bdydyn |
---|
[1125] | 2 | !!====================================================================== |
---|
[911] | 3 | !! *** MODULE bdydyn *** |
---|
[1125] | 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 | !!---------------------------------------------------------------------- |
---|
| 10 | #if defined key_bdy |
---|
| 11 | !!---------------------------------------------------------------------- |
---|
| 12 | !! 'key_bdy' : Unstructured Open Boundary Condition |
---|
| 13 | !!---------------------------------------------------------------------- |
---|
[911] | 14 | !! bdy_dyn_frs : relaxation of velocities on unstructured open boundary |
---|
| 15 | !! bdy_dyn_fla : Flather condition for barotropic solution |
---|
[1125] | 16 | !!---------------------------------------------------------------------- |
---|
[911] | 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 dynspg_oce ! for barotropic variables |
---|
| 21 | USE phycst ! physical constants |
---|
| 22 | USE lbclnk ! ocean lateral boundary conditions (or mpp link) |
---|
| 23 | USE bdytides ! for tidal harmonic forcing at boundary |
---|
[1125] | 24 | USE in_out_manager ! |
---|
[911] | 25 | |
---|
| 26 | IMPLICIT NONE |
---|
| 27 | PRIVATE |
---|
| 28 | |
---|
[1125] | 29 | PUBLIC bdy_dyn_frs ! routine called in dynspg_flt (free surface case ONLY) |
---|
| 30 | # if defined key_dynspg_exp || defined key_dynspg_ts |
---|
| 31 | PUBLIC bdy_dyn_fla ! routine called in dynspg_flt (free surface case ONLY) |
---|
| 32 | # endif |
---|
[911] | 33 | |
---|
[1125] | 34 | !!---------------------------------------------------------------------- |
---|
| 35 | !! NEMO/OPA 3.0 , LOCEAN-IPSL (2008) |
---|
| 36 | !! $Id: $ |
---|
| 37 | !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt) |
---|
| 38 | !!---------------------------------------------------------------------- |
---|
[911] | 39 | |
---|
| 40 | CONTAINS |
---|
| 41 | |
---|
[1125] | 42 | SUBROUTINE bdy_dyn_frs( kt ) |
---|
| 43 | !!---------------------------------------------------------------------- |
---|
| 44 | !! *** SUBROUTINE bdy_dyn_frs *** |
---|
| 45 | !! |
---|
[911] | 46 | !! ** Purpose : - Apply the Flow Relaxation Scheme for dynamic in the |
---|
| 47 | !! case of unstructured open boundaries. |
---|
| 48 | !! |
---|
| 49 | !! References :- Engedahl H., 1995: Use of the flow relaxation scheme in |
---|
| 50 | !! a three-dimensional baroclinic ocean model with realistic |
---|
| 51 | !! topography. Tellus, 365-382. |
---|
[1125] | 52 | !!---------------------------------------------------------------------- |
---|
| 53 | INTEGER, INTENT( in ) :: kt ! Main time step counter |
---|
[911] | 54 | !! |
---|
[1125] | 55 | INTEGER :: ib, ik, igrd ! dummy loop indices |
---|
| 56 | INTEGER :: ii, ij ! 2D addresses |
---|
| 57 | REAL(wp) :: zwgt ! boundary weight |
---|
| 58 | !!---------------------------------------------------------------------- |
---|
| 59 | ! |
---|
| 60 | IF(ln_bdy_dyn_frs) THEN ! If this is false, then this routine does nothing. |
---|
[911] | 61 | |
---|
[1125] | 62 | IF( kt == nit000 ) THEN |
---|
| 63 | IF(lwp) WRITE(numout,*) |
---|
| 64 | IF(lwp) WRITE(numout,*) 'bdy_dyn : Flow Relaxation Scheme on momentum' |
---|
| 65 | IF(lwp) WRITE(numout,*) '~~~~~~~' |
---|
| 66 | ENDIF |
---|
| 67 | ! |
---|
| 68 | igrd = 2 ! Relaxation of zonal velocity |
---|
| 69 | DO ib = 1, nblen(igrd) |
---|
| 70 | DO ik = 1, jpkm1 |
---|
| 71 | ii = nbi(ib,igrd) |
---|
| 72 | ij = nbj(ib,igrd) |
---|
| 73 | zwgt = nbw(ib,igrd) |
---|
| 74 | ua(ii,ij,ik) = ( ua(ii,ij,ik) * ( 1.- zwgt ) + ubdy(ib,ik) * zwgt ) * umask(ii,ij,ik) |
---|
| 75 | END DO |
---|
| 76 | END DO |
---|
| 77 | ! |
---|
| 78 | igrd = 3 ! Relaxation of meridional velocity |
---|
| 79 | DO ib = 1, nblen(igrd) |
---|
| 80 | DO ik = 1, jpkm1 |
---|
| 81 | ii = nbi(ib,igrd) |
---|
| 82 | ij = nbj(ib,igrd) |
---|
| 83 | zwgt = nbw(ib,igrd) |
---|
| 84 | va(ii,ij,ik) = ( va(ii,ij,ik) * ( 1.- zwgt ) + vbdy(ib,ik) * zwgt ) * vmask(ii,ij,ik) |
---|
| 85 | END DO |
---|
| 86 | END DO |
---|
| 87 | ! |
---|
| 88 | CALL lbc_lnk( ua, 'U', 1. ) ! Boundary points should be updated |
---|
| 89 | CALL lbc_lnk( va, 'V', 1. ) ! |
---|
| 90 | ! |
---|
| 91 | ENDIF ! ln_bdy_dyn_frs |
---|
[911] | 92 | |
---|
[1125] | 93 | END SUBROUTINE bdy_dyn_frs |
---|
[911] | 94 | |
---|
| 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 |
---|
[1125] | 99 | !!---------------------------------------------------------------------- |
---|
| 100 | !! *** SUBROUTINE bdy_dyn_fla *** |
---|
| 101 | !! |
---|
[911] | 102 | !! - Apply Flather boundary conditions on normal barotropic velocities |
---|
[1125] | 103 | !! (ln_bdy_dyn_fla=.true. or ln_bdy_tides=.true.) |
---|
[911] | 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 | !! |
---|
[1125] | 115 | !! References: 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 | !!---------------------------------------------------------------------- |
---|
| 118 | INTEGER :: ib, igrd ! dummy loop indices |
---|
| 119 | INTEGER :: ii, ij, iim1, iip1, ijm1, ijp1 ! 2D addresses |
---|
| 120 | REAL(wp) :: zcorr ! Flather correction |
---|
| 121 | !!---------------------------------------------------------------------- |
---|
[911] | 122 | |
---|
| 123 | ! ---------------------------------! |
---|
| 124 | ! Flather boundary conditions :! |
---|
| 125 | ! ---------------------------------! |
---|
| 126 | |
---|
[1125] | 127 | IF(ln_bdy_dyn_fla .OR. ln_bdy_tides) THEN ! If these are both false, then this routine does nothing. |
---|
| 128 | |
---|
[911] | 129 | ! Fill temporary array with ssh data (here spgu): |
---|
[1125] | 130 | igrd = 1 |
---|
| 131 | spgu(:,:) = 0.0 |
---|
| 132 | DO ib = 1, nblenrim(igrd) |
---|
| 133 | ii = nbi(ib,igrd) |
---|
| 134 | ij = nbj(ib,igrd) |
---|
| 135 | IF( ln_bdy_dyn_fla ) spgu(ii, ij) = sshbdy(ib) |
---|
| 136 | IF( ln_bdy_tides ) spgu(ii, ij) = spgu(ii, ij) + sshtide(ib) |
---|
[911] | 137 | END DO |
---|
[1125] | 138 | ! |
---|
| 139 | igrd = 2 ! Flather bc on u-velocity; |
---|
| 140 | ! ! remember that flagu=-1 if normal velocity direction is outward |
---|
| 141 | ! ! I think we should rather use after ssh ? |
---|
| 142 | DO ib = 1, nblenrim(igrd) |
---|
| 143 | ii = nbi(ib,igrd) |
---|
| 144 | ij = nbj(ib,igrd) |
---|
| 145 | iim1 = ii + MAX( 0, INT( flagu(ib) ) ) ! T pts i-indice inside the boundary |
---|
| 146 | iip1 = ii - MIN( 0, INT( flagu(ib) ) ) ! T pts i-indice outside the boundary |
---|
| 147 | ! |
---|
| 148 | zcorr = - flagu(ib) * SQRT( grav / (hu_e(ii, ij) + 1.e-20) ) * ( sshn_e(iim1, ij) - spgu(iip1,ij) ) |
---|
| 149 | ua_e(ii, ij) = ( ubtbdy(ib) + utide(ib) ) * hu_e(ii,ij) |
---|
| 150 | ua_e(ii,ij) = ua_e(ii,ij) + zcorr * umask(ii,ij,1) * hu_e(ii,ij) |
---|
[911] | 151 | END DO |
---|
[1125] | 152 | ! |
---|
| 153 | igrd = 3 ! Flather bc on v-velocity |
---|
| 154 | ! ! remember that flagv=-1 if normal velocity direction is outward |
---|
| 155 | DO ib = 1, nblenrim(igrd) |
---|
| 156 | ii = nbi(ib,igrd) |
---|
| 157 | ij = nbj(ib,igrd) |
---|
| 158 | ijm1 = ij + MAX( 0, INT( flagv(ib) ) ) ! T pts j-indice inside the boundary |
---|
| 159 | ijp1 = ij - MIN( 0, INT( flagv(ib) ) ) ! T pts j-indice outside the boundary |
---|
| 160 | ! |
---|
| 161 | zcorr = - flagv(ib) * SQRT( grav / (hv_e(ii, ij) + 1.e-20) ) * ( sshn_e(ii, ijm1) - spgu(ii,ijp1) ) |
---|
| 162 | va_e(ii, ij) = ( vbtbdy(ib) + vtide(ib) ) * hv_e(ii,ij) |
---|
| 163 | va_e(ii,ij) = va_e(ii,ij) + zcorr * vmask(ii,ij,1) * hv_e(ii,ij) |
---|
[911] | 164 | END DO |
---|
[1125] | 165 | ! |
---|
[911] | 166 | CALL lbc_lnk( ua_e, 'U', 1. ) ! Boundary points should be updated |
---|
| 167 | CALL lbc_lnk( va_e, 'V', 1. ) ! |
---|
[1125] | 168 | ! |
---|
| 169 | ENDIF ! ln_bdy_dyn_fla .or. ln_bdy_tides |
---|
[911] | 170 | |
---|
| 171 | END SUBROUTINE bdy_dyn_fla |
---|
| 172 | #endif |
---|
| 173 | |
---|
| 174 | #else |
---|
[1125] | 175 | !!---------------------------------------------------------------------- |
---|
| 176 | !! Dummy module NO Unstruct Open Boundary Conditions |
---|
| 177 | !!---------------------------------------------------------------------- |
---|
[911] | 178 | CONTAINS |
---|
[1125] | 179 | SUBROUTINE bdy_dyn_frs( kt ) ! Empty routine |
---|
| 180 | WRITE(*,*) 'bdy_dyn_frs: You should not have seen this print! error?', kt |
---|
[911] | 181 | END SUBROUTINE bdy_dyn_frs |
---|
[1125] | 182 | SUBROUTINE bdy_dyn_fla ! Empty routine |
---|
| 183 | WRITE(*,*) 'bdy_dyn_fla: You should not have seen this print! error?' |
---|
[911] | 184 | END SUBROUTINE bdy_dyn_fla |
---|
| 185 | #endif |
---|
| 186 | |
---|
[1125] | 187 | !!====================================================================== |
---|
[911] | 188 | END MODULE bdydyn |
---|