1 | MODULE obcdyn2d |
---|
2 | !!====================================================================== |
---|
3 | !! *** MODULE obcdyn *** |
---|
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_obc |
---|
14 | !!---------------------------------------------------------------------- |
---|
15 | !! 'key_obc' : Unstructured Open Boundary Condition |
---|
16 | !!---------------------------------------------------------------------- |
---|
17 | !! obc_dyn2d : Apply open boundary conditions to barotropic variables. |
---|
18 | !! obc_dyn2d_fla : Apply Flather condition |
---|
19 | !!---------------------------------------------------------------------- |
---|
20 | USE oce ! ocean dynamics and tracers |
---|
21 | USE dom_oce ! ocean space and time domain |
---|
22 | USE obc_oce ! ocean open boundary conditions |
---|
23 | USE dynspg_oce ! for barotropic variables |
---|
24 | USE phycst ! physical constants |
---|
25 | USE lbclnk ! ocean lateral boundary conditions (or mpp link) |
---|
26 | USE obctides ! for tidal harmonic forcing at boundary |
---|
27 | USE in_out_manager ! |
---|
28 | |
---|
29 | IMPLICIT NONE |
---|
30 | PRIVATE |
---|
31 | |
---|
32 | PUBLIC obc_dyn2d ! routine called in dynspg_ts (time splitting case ONLY) |
---|
33 | |
---|
34 | !!---------------------------------------------------------------------- |
---|
35 | !! NEMO/OPA 3.3 , NEMO Consortium (2010) |
---|
36 | !! $Id: obcdyn.F90 2528 2010-12-27 17:33:53Z rblod $ |
---|
37 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
38 | !!---------------------------------------------------------------------- |
---|
39 | CONTAINS |
---|
40 | |
---|
41 | SUBROUTINE obc_dyn2d( kt, pssh ) |
---|
42 | !!---------------------------------------------------------------------- |
---|
43 | !! *** SUBROUTINE obc_dyn2d *** |
---|
44 | !! |
---|
45 | !! ** Purpose : - Apply open boundary conditions for barotropic variables |
---|
46 | !! |
---|
47 | !!---------------------------------------------------------------------- |
---|
48 | INTEGER, INTENT(in) :: kt ! Main time step counter |
---|
49 | REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: pssh |
---|
50 | !! |
---|
51 | INTEGER :: ib_obc ! Loop counter |
---|
52 | |
---|
53 | DO ib_obc=1, nb_obc |
---|
54 | |
---|
55 | SELECT CASE( nn_dyn2d(ib_obc) ) |
---|
56 | CASE(jp_none) |
---|
57 | CYCLE |
---|
58 | CASE(jp_flather) |
---|
59 | CALL obc_dyn2d_fla( idx_obc(ib_obc), dta_obc(ib_obc), pssh ) |
---|
60 | CASE DEFAULT |
---|
61 | CALL ctl_stop( 'obc_dyn3d : unrecognised option for open boundaries for barotropic variables' ) |
---|
62 | END SELECT |
---|
63 | ENDDO |
---|
64 | |
---|
65 | END SUBROUTINE obc_dyn2d |
---|
66 | |
---|
67 | SUBROUTINE obc_dyn2d_fla( idx, dta, pssh ) |
---|
68 | !!---------------------------------------------------------------------- |
---|
69 | !! *** SUBROUTINE obc_dyn_fla *** |
---|
70 | !! |
---|
71 | !! - Apply Flather boundary conditions on normal barotropic velocities |
---|
72 | !! (ln_dyn_fla=.true. or ln_tides=.true.) |
---|
73 | !! |
---|
74 | !! ** WARNINGS about FLATHER implementation: |
---|
75 | !!1. According to Palma and Matano, 1998 "after ssh" is used. |
---|
76 | !! In ROMS and POM implementations, it is "now ssh". In the current |
---|
77 | !! implementation (tested only in the EEL-R5 conf.), both cases were unstable. |
---|
78 | !! So I use "before ssh" in the following. |
---|
79 | !! |
---|
80 | !!2. We assume that the normal ssh gradient at the obc is zero. As a matter of |
---|
81 | !! fact, the model ssh just inside the dynamical boundary is used (the outside |
---|
82 | !! ssh in the code is not updated). |
---|
83 | !! |
---|
84 | !! References: Flather, R. A., 1976: A tidal model of the northwest European |
---|
85 | !! continental shelf. Mem. Soc. R. Sci. Liege, Ser. 6,10, 141-164. |
---|
86 | !!---------------------------------------------------------------------- |
---|
87 | TYPE(OBC_INDEX), INTENT(in) :: idx ! OBC indices |
---|
88 | TYPE(OBC_DATA), INTENT(in) :: dta ! OBC external data |
---|
89 | REAL(wp), DIMENSION(jpi,jpj), INTENT(in) :: pssh |
---|
90 | |
---|
91 | INTEGER :: jb, igrd ! dummy loop indices |
---|
92 | INTEGER :: ii, ij, iim1, iip1, ijm1, ijp1 ! 2D addresses |
---|
93 | REAL(wp) :: zcorr ! Flather correction |
---|
94 | REAL(wp) :: zforc ! temporary scalar |
---|
95 | !!---------------------------------------------------------------------- |
---|
96 | |
---|
97 | ! ---------------------------------! |
---|
98 | ! Flather boundary conditions :! |
---|
99 | ! ---------------------------------! |
---|
100 | |
---|
101 | !!!!!!!!!!!! SOME WORK TO DO ON THE TIDES HERE !!!!!!!!!!!!! |
---|
102 | |
---|
103 | ! Fill temporary array with ssh data (here spgu): |
---|
104 | igrd = 1 |
---|
105 | spgu(:,:) = 0.0 |
---|
106 | DO jb = 1, idx%nblenrim(igrd) |
---|
107 | ii = idx%nbi(jb,igrd) |
---|
108 | ij = idx%nbj(jb,igrd) |
---|
109 | spgu(ii, ij) = dta%ssh(jb) |
---|
110 | !!$ IF( ln_tides ) spgu(ii, ij) = spgu(ii, ij) + sshtide(jb) |
---|
111 | END DO |
---|
112 | ! |
---|
113 | igrd = 2 ! Flather bc on u-velocity; |
---|
114 | ! ! remember that flagu=-1 if normal velocity direction is outward |
---|
115 | ! ! I think we should rather use after ssh ? |
---|
116 | DO jb = 1, idx%nblenrim(igrd) |
---|
117 | ii = idx%nbi(jb,igrd) |
---|
118 | ij = idx%nbj(jb,igrd) |
---|
119 | iim1 = ii + MAX( 0, INT( idx%flagu(jb) ) ) ! T pts i-indice inside the boundary |
---|
120 | iip1 = ii - MIN( 0, INT( idx%flagu(jb) ) ) ! T pts i-indice outside the boundary |
---|
121 | ! |
---|
122 | zcorr = - idx%flagu(jb) * SQRT( grav * hur_e(ii, ij) ) * ( pssh(iim1, ij) - spgu(iip1,ij) ) |
---|
123 | !!$ zforc = dta%u2d(jb) + utide(jb) |
---|
124 | zforc = dta%u2d(jb) |
---|
125 | ua_e(ii,ij) = zforc + zcorr * umask(ii,ij,1) |
---|
126 | END DO |
---|
127 | ! |
---|
128 | igrd = 3 ! Flather bc on v-velocity |
---|
129 | ! ! remember that flagv=-1 if normal velocity direction is outward |
---|
130 | DO jb = 1, idx%nblenrim(igrd) |
---|
131 | ii = idx%nbi(jb,igrd) |
---|
132 | ij = idx%nbj(jb,igrd) |
---|
133 | ijm1 = ij + MAX( 0, INT( idx%flagv(jb) ) ) ! T pts j-indice inside the boundary |
---|
134 | ijp1 = ij - MIN( 0, INT( idx%flagv(jb) ) ) ! T pts j-indice outside the boundary |
---|
135 | ! |
---|
136 | zcorr = - idx%flagv(jb) * SQRT( grav * hvr_e(ii, ij) ) * ( pssh(ii, ijm1) - spgu(ii,ijp1) ) |
---|
137 | !!$ zforc = dta%v2d(jb) + vtide(jb) |
---|
138 | zforc = dta%v2d(jb) |
---|
139 | va_e(ii,ij) = zforc + zcorr * vmask(ii,ij,1) |
---|
140 | END DO |
---|
141 | CALL lbc_lnk( ua_e, 'U', -1. ) ! Boundary points should be updated |
---|
142 | CALL lbc_lnk( va_e, 'V', -1. ) ! |
---|
143 | ! |
---|
144 | ! |
---|
145 | END SUBROUTINE obc_dyn2d_fla |
---|
146 | #else |
---|
147 | !!---------------------------------------------------------------------- |
---|
148 | !! Dummy module NO Unstruct Open Boundary Conditions |
---|
149 | !!---------------------------------------------------------------------- |
---|
150 | CONTAINS |
---|
151 | SUBROUTINE obc_dyn2d( kt ) ! Empty routine |
---|
152 | WRITE(*,*) 'obc_dyn_frs: You should not have seen this print! error?', kt |
---|
153 | END SUBROUTINE obc_dyn2d |
---|
154 | #endif |
---|
155 | |
---|
156 | !!====================================================================== |
---|
157 | END MODULE obcdyn2d |
---|