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.
cla.F90 in branches/2014/dev_r4650_UKMO14.4_OBS_GENERAL_VINTERP/NEMOGCM/NEMO/OPA_SRC/LBC – NEMO

source: branches/2014/dev_r4650_UKMO14.4_OBS_GENERAL_VINTERP/NEMOGCM/NEMO/OPA_SRC/LBC/cla.F90 @ 5950

Last change on this file since 5950 was 5950, checked in by timgraham, 8 years ago

Reinstated Id keyword before merging

  • Property svn:keywords set to Id
File size: 41.8 KB
Line 
1MODULE cla
2   !!======================================================================
3   !!                    ***  MODULE  cla  ***
4   !! Cross Land Advection : specific update of the horizontal divergence,
5   !!                        tracer trends and after velocity
6   !!
7   !!                 ---   Specific to ORCA_R2   ---
8   !!
9   !!======================================================================
10   !! History :  1.0  ! 2002-11 (A. Bozec)  Original code
11   !!            3.2  ! 2009-07 (G. Madec)  merge cla, cla_div, tra_cla, cla_dynspg
12   !!                 !                     and correct a mpp bug reported by A.R. Porter
13   !!----------------------------------------------------------------------
14   !!   cla_div           : update of horizontal divergence at cla straits
15   !!   tra_cla           : update of tracers at cla straits
16   !!   cla_dynspg        : update of after horizontal velocities at cla straits
17   !!   cla_init          : initialisation - control check
18   !!   cla_bab_el_mandeb : cross land advection for Bab-el-mandeb strait
19   !!   cla_gibraltar     : cross land advection for Gibraltar strait
20   !!   cla_hormuz        : cross land advection for Hormuz strait
21   !!----------------------------------------------------------------------
22   USE oce            ! ocean dynamics and tracers
23   USE dom_oce        ! ocean space and time domain
24   USE sbc_oce        ! surface boundary condition: ocean
25   USE dynspg_oce     ! ocean dynamics: surface pressure gradient variables
26   USE in_out_manager ! I/O manager
27   USE lib_mpp        ! distributed memory computing library
28   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
29   USE lib_mpp        ! MPP library
30
31   IMPLICIT NONE
32   PRIVATE
33   
34   PUBLIC   cla_init     ! routine called by opa.F90
35   PUBLIC   cla_div      ! routine called by divcur.F90
36   PUBLIC   cla_traadv   ! routine called by traadv.F90
37   PUBLIC   cla_dynspg   ! routine called by dynspg_flt.F90
38
39   INTEGER  ::   nbab, ngib, nhor   ! presence or not of required grid-point on local domain
40   !                                ! for Bab-el-Mandeb, Gibraltar, and Hormuz straits
41   
42   !                                           !   fixed part  !  time evolving    !!! profile of hdiv for some straits
43   REAL(wp), ALLOCATABLE, SAVE, DIMENSION (:) ::   hdiv_139_101, hdiv_139_101_kt    ! Gibraltar     (i,j)=(172,101)
44   REAL(wp), ALLOCATABLE, SAVE, DIMENSION (:) ::   hdiv_139_102                     ! Gibraltar     (i,j)=(139,102)
45   REAL(wp), ALLOCATABLE, SAVE, DIMENSION (:) ::   hdiv_141_102, hdiv_141_102_kt    ! Gibraltar     (i,j)=(141,102)
46   REAL(wp), ALLOCATABLE, SAVE, DIMENSION (:) ::   hdiv_161_88 , hdiv_161_88_kt     ! Bab-el-Mandeb (i,j)=(161,88)
47   REAL(wp), ALLOCATABLE, SAVE, DIMENSION (:) ::   hdiv_161_87                      ! Bab-el-Mandeb (i,j)=(161,87)
48   REAL(wp), ALLOCATABLE, SAVE, DIMENSION (:) ::   hdiv_160_89 , hdiv_160_89_kt     ! Bab-el-Mandeb (i,j)=(160,89)
49   REAL(wp), ALLOCATABLE, SAVE, DIMENSION (:) ::   hdiv_172_94                      ! Hormuz        (i,j)=(172, 94)
50
51   REAL(wp), ALLOCATABLE, SAVE, DIMENSION (:) ::   t_171_94_hor, s_171_94_hor       ! Temperature, salinity in Hormuz strait
52   
53   !! * Substitutions
54#  include "domzgr_substitute.h90"
55   !!----------------------------------------------------------------------
56   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
57   !! $Id$
58   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
59   !!----------------------------------------------------------------------
60CONTAINS
61
62   SUBROUTINE cla_div( kt )
63      !!----------------------------------------------------------------------
64      !!                 ***  ROUTINE div_cla  ***
65      !!
66      !! ** Purpose :   update the horizontal divergence of the velocity field
67      !!              at some straits ( Gibraltar, Bab el Mandeb and Hormuz ).
68      !!
69      !! ** Method  : - first time-step: initialisation of cla
70      !!              - all   time-step: using imposed transport at each strait,
71      !!              the now horizontal divergence is updated
72      !!
73      !! ** Action  :   phdivn   updted now horizontal divergence at cla straits
74      !!----------------------------------------------------------------------
75      INTEGER, INTENT( in ) ::   kt      ! ocean time-step index
76      !!----------------------------------------------------------------------
77      !     
78      IF( kt == nit000 ) THEN
79         !
80         CALL cla_init                                        ! control check
81         !
82         IF(lwp) WRITE(numout,*)
83         IF(lwp) WRITE(numout,*) 'div_cla : cross land advection on hdiv '
84         IF(lwp) WRITE(numout,*) '~~~~~~~~'
85         !
86         IF( nbab == 1 )   CALL cla_bab_el_mandeb('ini')    ! Bab el Mandeb ( Red Sea - Indian ocean )
87         IF( ngib == 1 )   CALL cla_gibraltar    ('ini')    ! Gibraltar strait (Med Sea - Atlantic ocean)
88         IF( nhor == 1 )   CALL cla_hormuz       ('ini')    ! Hormuz Strait ( Persian Gulf - Indian ocean )
89         !
90      ENDIF                           
91      !
92      IF( nbab == 1    )   CALL cla_bab_el_mandeb('div')    ! Bab el Mandeb ( Red Sea - Indian ocean )
93      IF( ngib == 1    )   CALL cla_gibraltar    ('div')    ! Gibraltar strait (Med Sea - Atlantic ocean)
94      IF( nhor == 1    )   CALL cla_hormuz       ('div')    ! Hormuz Strait ( Persian Gulf - Indian ocean )
95      !
96!!gm  lbc useless here, no?
97!!gm      CALL lbc_lnk( hdivn, 'T', 1. )                    ! Lateral boundary conditions on hdivn
98      !
99   END SUBROUTINE cla_div
100   
101   
102   SUBROUTINE cla_traadv( kt )
103      !!----------------------------------------------------------------------
104      !!                 ***  ROUTINE tra_cla  ***
105      !!                   
106      !! ** Purpose :   Update the now trend due to the advection of tracers
107      !!      and add it to the general trend of passive tracer equations
108      !!      at some straits ( Bab el Mandeb, Gibraltar, Hormuz ).
109      !!
110      !! ** Method  :   using both imposed transport at each strait and T & S
111      !!              budget, the now tracer trends is updated
112      !!
113      !! ** Action  :   (ta,sa)   updated now tracer trends at cla straits
114      !!----------------------------------------------------------------------
115      INTEGER, INTENT( in ) ::   kt         ! ocean time-step index
116      !!----------------------------------------------------------------------
117      !
118      IF( kt == nit000 ) THEN
119         IF(lwp) WRITE(numout,*)
120         IF(lwp) WRITE(numout,*) 'tra_cla : cross land advection on tracers '
121         IF(lwp) WRITE(numout,*) '~~~~~~~~'
122      ENDIF
123      !
124      IF( nbab == 1    )   CALL cla_bab_el_mandeb('tra')      ! Bab el Mandeb strait
125      IF( ngib == 1    )   CALL cla_gibraltar    ('tra')      ! Gibraltar strait
126      IF( nhor == 1    )   CALL cla_hormuz       ('tra')      ! Hormuz Strait ( Persian Gulf)
127      !
128   END SUBROUTINE cla_traadv
129
130   
131   SUBROUTINE cla_dynspg( kt )
132      !!----------------------------------------------------------------------
133      !!                 ***  ROUTINE cla_dynspg  ***
134      !!                   
135      !! ** Purpose :   Update the after velocity at some straits
136      !!              (Bab el Mandeb, Gibraltar, Hormuz).
137      !!
138      !! ** Method  :   required to compute the filtered surface pressure gradient
139      !!
140      !! ** Action  :   (ua,va)   after velocity at the cla straits
141      !!----------------------------------------------------------------------
142      INTEGER, INTENT( in ) ::   kt         ! ocean time-step index
143      !!----------------------------------------------------------------------
144      !
145      IF( kt == nit000 ) THEN
146         IF(lwp) WRITE(numout,*)
147         IF(lwp) WRITE(numout,*) 'cla_dynspg : cross land advection on (ua,va) '
148         IF(lwp) WRITE(numout,*) '~~~~~~~~~~'
149      ENDIF
150      !
151      IF( nbab == 1    )   CALL cla_bab_el_mandeb('spg')      ! Bab el Mandeb strait
152      IF( ngib == 1    )   CALL cla_gibraltar    ('spg')      ! Gibraltar strait
153      IF( nhor == 1    )   CALL cla_hormuz       ('spg')      ! Hormuz Strait ( Persian Gulf)
154      !
155!!gm lbc is needed here, not?
156!!gm      CALL lbc_lnk( hdivn, 'U', -1. )   ;   CALL lbc_lnk( hdivn, 'V', -1. )      ! Lateral boundary conditions
157      !
158   END SUBROUTINE cla_dynspg
159
160
161   SUBROUTINE cla_init
162      !! -------------------------------------------------------------------
163      !!                   ***  ROUTINE cla_init  ***
164      !!           
165      !! ** Purpose :   control check for mpp computation 
166      !!
167      !! ** Method  : - All the strait grid-points must be inside one of the
168      !!              local domain interior for the cla advection to work
169      !!              properly in mpp (i.e. inside (2:jpim1,2:jpjm1) ).
170      !!              Define the corresponding indicators (nbab, ngib, nhor)
171      !!              - The profiles of cross-land fluxes are currently hard
172      !!              coded for L31 levels. Stop if jpk/=31
173      !!
174      !! ** Action  :   nbab, ngib, nhor   strait inside the local domain or not
175      !!---------------------------------------------------------------------
176      REAL(wp) ::   ztemp   ! local scalar
177      INTEGER  ::   ierr    ! local integer
178      !!---------------------------------------------------------------------
179      !
180      IF(lwp) WRITE(numout,*)
181      IF(lwp) WRITE(numout,*) 'cla_init : cross land advection initialisation '
182      IF(lwp) WRITE(numout,*) '~~~~~~~~~'
183      !
184      !                           ! Allocate arrays for this module
185      ALLOCATE( hdiv_139_101(jpk) , hdiv_139_101_kt(jpk) ,     &    ! Gibraltar
186         &      hdiv_139_102(jpk) ,                            &
187         &      hdiv_141_102(jpk) , hdiv_141_102_kt(jpk) ,     &
188         &      hdiv_161_88 (jpk) , hdiv_161_88_kt (jpk) ,     &    ! Bab-el-Mandeb
189         &      hdiv_161_87 (jpk) ,                            &                     
190         &      hdiv_160_89 (jpk) , hdiv_160_89_kt (jpk) ,     &     ! Hormuz
191         &      hdiv_172_94 (jpk) ,                            &
192         &      t_171_94_hor(jpk) , s_171_94_hor   (jpk) , STAT=ierr )
193      IF( lk_mpp    )   CALL mpp_sum( ierr )
194      IF( ierr /= 0 )   CALL ctl_stop( 'STOP', 'cla_init: unable to allocate arrays' )
195      !
196      IF( .NOT.lk_dynspg_flt )   CALL ctl_stop( 'cla_init: Cross Land Advection works only with lk_dynspg_flt=T ' )
197      !
198      IF( lk_vvl             )   CALL ctl_stop( 'cla_init: Cross Land Advection does not work with lk_vvl=T option' )
199      !
200      IF( jpk /= 31          )   CALL ctl_stop( 'cla_init: Cross Land Advection hard coded for ORCA_R2_L31' )
201      !
202      !                                        _|_______|_______|_
203      !                                     89  |       |///////| 
204      !                                        _|_______|_______|_
205      ! ------------------------ !          88  |///////|       |
206      !   Bab el Mandeb strait   !             _|_______|_______|_
207      ! ------------------------ !          87  |///////|       |
208      !                                        _|_______|_______|_
209      !                                         |  160  |  161  | 
210      !
211      ! The 6 Bab el Mandeb grid-points must be inside one of the interior of the
212      ! local domain for the cla advection to work properly (i.e. (2:jpim1,2:jpjm1)
213      nbab = 0
214      IF(  ( 1 <= mj0( 88) .AND. mj1( 89) <= jpj ) .AND.    &  !* (161,89), (161,88) and (161,88) on the local pocessor
215         & ( 1 <= mi0(160) .AND. mi1(161) <= jpi )       )    nbab = 1 
216      !
217      ! test if there is no local domain that includes all required grid-points
218      ztemp = REAL( nbab )
219      IF( lk_mpp )   CALL mpp_sum( ztemp )      ! sum with other processors value
220      IF( ztemp == 0 ) THEN                     ! Only 2 points in each direction, this should never be a problem
221         CALL ctl_stop( ' cross land advection at Bab-el_Mandeb does not work with your processor cutting: change it' )
222      ENDIF
223      !                                        ___________________________
224      ! ------------------------ !         102  |       |///////|       |
225      !     Gibraltar strait     !             _|_______|_______|_______|_
226      ! ------------------------ !         101  |       |///////|       |
227      !                                        _|_______|_______|_______|_
228      !                                         |  139  |  140  |  141  |
229      !
230      ! The 6 Gibraltar grid-points must be inside one of the interior of the
231      ! local domain for the cla advection to work properly (i.e. (2:jpim1,2:jpjm1)
232      ngib = 0
233      IF(  ( 2 <= mj0(101) .AND. mj1(102) <= jpjm1 ) .AND.    &  !* (139:141,101:102) on the local pocessor
234         & ( 2 <= mi0(139) .AND. mi1(141) <= jpim1 )       )    ngib = 1 
235      !
236      ! test if there is no local domain that includes all required grid-points
237      ztemp = REAL( ngib )
238      IF( lk_mpp )   CALL mpp_sum( ztemp )      ! sum with other processors value
239      IF( ztemp == 0 ) THEN                     ! 3 points in i-direction, this may be a problem with some cutting
240           CALL ctl_stop( ' cross land advection at Gibraltar does not work with your processor cutting: change it' )
241      ENDIF
242      !                                        _______________
243      ! ------------------------ !          94  |/////|     |
244      !       Hormuz strait      !             _|_____|_____|_
245      ! ------------------------ !                171   172     
246      !           
247      ! The 2 Hormuz grid-points must be inside one of the interior of the
248      ! local domain for the cla advection to work properly (i.e. (2:jpim1,2:jpjm1)
249      nhor = 0
250      IF(    2 <= mj0( 94) .AND. mj1( 94) <= jpjm1  .AND.  & 
251         &   2 <= mi0(171) .AND. mi1(172) <= jpim1         )   nhor = 1 
252      !
253      ! test if there is no local domain that includes all required grid-points
254      ztemp = REAL( nhor )
255      IF( lk_mpp )   CALL mpp_sum( ztemp )      ! sum with other processors value
256      IF( ztemp == 0 ) THEN                     ! 3 points in i-direction, this may be a problem with some cutting
257           CALL ctl_stop( ' cross land advection at Hormuz does not work with your processor cutting: change it' )
258      ENDIF
259      !
260   END SUBROUTINE cla_init
261
262
263   SUBROUTINE cla_bab_el_mandeb( cd_td )
264      !!----------------------------------------------------------------------
265      !!                ***  ROUTINE cla_bab_el_mandeb  ***
266      !!       
267      !! ** Purpose :   update the now horizontal divergence, the tracer tendancy
268      !!              and the after velocity in vicinity of Bab el Mandeb ( Red Sea - Indian ocean).
269      !!
270      !! ** Method  :   compute the exchanges at each side of the strait :
271      !!
272      !!       surf. zio_flow
273      !! (+ balance of emp) /\  |\\\\\\\\\\\|
274      !!                    ||  |\\\\\\\\\\\| 
275      !!    deep zio_flow   ||  |\\\\\\\\\\\| 
276      !!            |  ||   ||  |\\\\\\\\\\\| 
277      !!        89  |  ||   ||  |\\\\\\\\\\\| 
278      !!            |__\/_v_||__|____________
279      !!            !\\\\\\\\\\\|          surf. zio_flow
280      !!            |\\\\\\\\\\\|<===    (+ balance of emp)
281      !!            |\\\\\\\\\\\u
282      !!        88  |\\\\\\\\\\\|<---      deep  zrecirc (upper+deep at 2 different levels)
283      !!            |___________|__________   
284      !!            !\\\\\\\\\\\|         
285      !!            |\\\\\\\\\\\| ---\     deep  zrecirc (upper+deep)
286      !!        87  !\\\\\\\\\\\u ===/   + deep  zio_flow   (all at the same level)
287      !!            !\\\\\\\\\\\| 
288      !!            !___________|__________
289      !!                160         161
290      !!
291      !!----------------------------------------------------------------------
292      CHARACTER(len=1), INTENT(in) ::   cd_td   ! ='div' update the divergence
293      !                                         ! ='tra' update the tracers
294      !                                         ! ='spg' update after velocity
295      INTEGER  ::   ji, jj, jk   ! dummy loop indices
296      REAL(wp) ::   zemp_red     ! temporary scalar
297      REAL(wp) ::   zio_flow, zrecirc_upp, zrecirc_mid, zrecirc_bot
298      !!---------------------------------------------------------------------
299      !
300      SELECT CASE( cd_td ) 
301      !                     ! ---------------- !
302      CASE( 'ini' )         !  initialisation  !
303         !                  ! ---------------- !
304         !                                   
305         zio_flow    = 0.4e6                       ! imposed in/out flow
306         zrecirc_upp = 0.2e6                       ! imposed upper recirculation water
307         zrecirc_bot = 0.5e6                       ! imposed bottom  recirculation water
308
309         hdiv_161_88(:) = 0.e0                     ! (161,88) Gulf of Aden side, north point
310         hdiv_161_87(:) = 0.e0                     ! (161,87) Gulf of Aden side, south point
311         hdiv_160_89(:) = 0.e0                     ! (160,89) Red sea side
312
313         DO jj = mj0(88), mj1(88)              !** profile of hdiv at (161,88)   (Gulf of Aden side, north point)
314            DO ji = mi0(161), mi1(161)         !------------------------------
315               DO jk = 1, 8                        ! surface in/out flow   (Ind -> Red)   (div >0)
316                  hdiv_161_88(jk) = + zio_flow / ( 8. * e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) )
317               END DO
318               !                                   ! recirculation water   (Ind -> Red)   (div >0)
319               hdiv_161_88(20) =                 + zrecirc_upp   / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,20) )
320               hdiv_161_88(21) = + ( zrecirc_bot - zrecirc_upp ) / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,21) )
321            END DO
322         END DO
323         !
324         DO jj = mj0(87), mj1(87)              !** profile of hdiv at (161,88)   (Gulf of Aden side, south point)
325            DO ji = mi0(161), mi1(161)         !------------------------------
326               !                                   ! deep out flow + recirculation   (Red -> Ind)   (div <0)
327               hdiv_161_87(21) = - ( zio_flow + zrecirc_bot ) / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,21) )
328            END DO
329         END DO
330         !
331         DO jj = mj0(89), mj1(89)              !** profile of hdiv at (161,88)   (Red sea side)
332            DO ji = mi0(160), mi1(160)         !------------------------------
333               DO jk = 1, 8                        ! surface inflow    (Ind -> Red)   (div <0)
334                  hdiv_160_89(jk) = - zio_flow / ( 8. * e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) )
335               END DO
336               !                                   ! deep    outflow   (Red -> Ind)   (div >0)
337               hdiv_160_89(16)    = + zio_flow / (      e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,16) )
338            END DO
339         END DO
340         !                  ! ---------------- !
341      CASE( 'div' )         !   update hdivn   ! (call by divcur module)
342         !                  ! ---------=====-- !
343         !                                     !** emp on the Red Sea   (div >0)
344         zemp_red = 0.e0                       !---------------------
345         DO jj = mj0(87), mj1(96)                  ! sum over the Red sea
346            DO ji = mi0(148), mi1(160) 
347               zemp_red = zemp_red + emp(ji,jj) * e1t(ji,jj) * e2t(ji,jj) * tmask_i(ji,jj)
348            END DO
349         END DO
350         IF( lk_mpp )   CALL mpp_sum( zemp_red )   ! sum with other processors value
351         zemp_red = zemp_red * 1.e-3               ! convert in m3
352         !
353         !                                     !** Correct hdivn (including emp adjustment)
354         !                                     !-------------------------------------------
355         DO jj = mj0(88), mj1(88)                  !* profile of hdiv at (161,88)   (Gulf of Aden side, north point)
356            DO ji = mi0(161), mi1(161) 
357               hdiv_161_88_kt(:) = hdiv_161_88(:)
358               DO jk = 1, 8                              ! increase the inflow from the Indian   (div >0)
359                  hdiv_161_88_kt(jk) = hdiv_161_88(jk) + zemp_red / (8. * e2u(ji,jj) * fse3u(ji,jj,jk) )
360               END DO
361               hdivn(ji,jj,:) = hdivn(ji,jj,:) + hdiv_161_88_kt(:)
362            END DO
363         END DO
364         DO jj = mj0(87), mj1(87)                  !* profile of divergence at (161,87)   (Gulf of Aden side, south point)
365            DO ji = mi0(161), mi1(161) 
366               hdivn(ji,jj,:) = hdivn(ji,jj,:) + hdiv_161_87(:)
367            END DO
368         END DO
369         DO jj = mj0(89), mj1(89)                  !* profile of divergence at (160,89)   (Red sea side)
370            DO ji = mi0(160), mi1(160) 
371               hdiv_160_89_kt(:) = hdiv_160_89(:)
372               DO jk = 1, 18                              ! increase the inflow from the Indian   (div <0)
373                  hdiv_160_89_kt(jk) = hdiv_160_89(jk) - zemp_red / (10. * e1v(ji,jj) * fse3v(ji,jj,jk) )
374               END DO
375               hdivn(ji, jj,:) = hdivn(ji, jj,:) + hdiv_160_89_kt(:)
376            END DO
377         END DO
378         !                  ! ---------------- !
379      CASE( 'tra' )         !  update (ta,sa)  ! (call by traadv module)
380         !                  ! --------=======- !
381         !
382         DO jj = mj0(88), mj1(88)              !** (161,88)   (Gulf of Aden side, north point)
383            DO ji = mi0(161), mi1(161) 
384               DO jk = 1, jpkm1                         ! surf inflow + reciculation (from Gulf of Aden)
385                  tsa(ji,jj,jk,jp_tem) = tsa(ji,jj,jk,jp_tem) - hdiv_161_88_kt(jk) * tsn(ji,jj,jk,jp_tem)
386                  tsa(ji,jj,jk,jp_sal) = tsa(ji,jj,jk,jp_sal) - hdiv_161_88_kt(jk) * tsn(ji,jj,jk,jp_sal)
387               END DO
388            END DO
389         END DO
390         DO jj = mj0(87), mj1(87)              !** (161,87)   (Gulf of Aden side, south point)
391            DO ji = mi0(161), mi1(161) 
392               jk =  21                                 ! deep outflow + recirulation (combined flux)
393               tsa(ji,jj,jk,jp_tem) = tsa(ji,jj,jk,jp_tem) + hdiv_161_88(20) * tsn(ji  ,jj+1,20,jp_tem)   &  ! upper recirculation from Gulf of Aden
394                  &                        + hdiv_161_88(21) * tsn(ji  ,jj+1,21,jp_tem)   &  ! deep  recirculation from Gulf of Aden
395                  &                        + hdiv_160_89(16) * tsn(ji-1,jj+2,16,jp_tem)      ! deep inflow from Red sea
396               tsa(ji,jj,jk,jp_sal) = tsa(ji,jj,jk,jp_sal) + hdiv_161_88(20) * tsn(ji  ,jj+1,20,jp_sal)   &
397                  &                        + hdiv_161_88(21) * tsn(ji  ,jj+1,21,jp_sal)   &
398                  &                        + hdiv_160_89(16) * tsn(ji-1,jj+2,16,jp_sal)   
399            END DO
400         END DO
401         DO jj = mj0(89), mj1(89)              !** (161,88)   (Red sea side)
402            DO ji = mi0(160), mi1(160)
403               DO jk = 1, 14                            ! surface inflow (from Gulf of Aden)
404                  tsa(ji,jj,jk,jp_tem) = tsa(ji,jj,jk,jp_tem) - hdiv_160_89_kt(jk) * tsn(ji+1,jj-1,jk,jp_tem)
405                  tsa(ji,jj,jk,jp_sal) = tsa(ji,jj,jk,jp_sal) - hdiv_160_89_kt(jk) * tsn(ji+1,jj-1,jk,jp_sal)
406               END DO
407               !                                  ! deep    outflow (from Red sea)
408               tsa(ji,jj,16,jp_tem) = tsa(ji,jj,16,jp_tem) - hdiv_160_89(16) * tsn(ji,jj,16,jp_tem)
409               tsa(ji,jj,16,jp_sal) = tsa(ji,jj,16,jp_sal) - hdiv_160_89(16) * tsn(ji,jj,16,jp_sal)
410            END DO
411         END DO
412         !
413         !                  ! ---------------- !
414      CASE( 'spg' )         !  update (ua,va)  ! (call by dynspg module)
415         !                  ! --------=======- !
416         ! at this stage, (ua,va) are the after velocity, not the tendancy
417         ! compute the velocity from the divergence at T-point
418         !
419         DO jj = mj0(88), mj1(88)              !** (160,88)   (Gulf of Aden side, north point)
420            DO ji = mi0(160), mi1(160)                   ! 160, not 161 as it is a U-point)
421               ua(ji,jj,:) = - hdiv_161_88_kt(:) / ( e1t(ji+1,jj) * e2t(ji+1,jj) * fse3t(ji+1,jj,:) )   &
422                  &                              * e2u(ji,jj) * fse3u(ji,jj,:)
423            END DO
424         END DO
425         DO jj = mj0(87), mj1(87)              !** (160,87)   (Gulf of Aden side, south point)
426            DO ji = mi0(160), mi1(160)                   ! 160, not 161 as it is a U-point)
427               ua(ji,jj,:) = - hdiv_161_87(:) / ( e1t(ji+1,jj) * e2t(ji+1,jj) * fse3t(ji+1,jj,:) )   &
428                  &                           * e2u(ji,jj) * fse3u(ji,jj,:)
429            END DO
430         END DO
431         DO jj = mj0(88), mj1(88)              !** profile of divergence at (160,89)   (Red sea side)
432            DO ji = mi0(160), mi1(160)                   ! 88, not 89 as it is a V-point)
433               va(ji,jj,:) = - hdiv_160_89_kt(:) / ( e1t(ji,jj+1) * e2t(ji,jj+1) * fse3t(ji,jj+1,:) )   &
434                  &                              * e1v(ji,jj) * fse3v(ji,jj,:)
435            END DO
436         END DO
437      END SELECT
438      !
439   END SUBROUTINE cla_bab_el_mandeb
440   
441
442   SUBROUTINE cla_gibraltar( cd_td )
443      !! -------------------------------------------------------------------
444      !!                 ***  ROUTINE cla_gibraltar  ***
445      !!       
446      !! ** Purpose :   update the now horizontal divergence, the tracer
447      !!              tendancyand the after velocity in vicinity of Gibraltar
448      !!              strait ( Persian Gulf - Indian ocean ).
449      !!
450      !! ** Method :
451      !!                     _______________________           
452      !!      deep  zio_flow /====|///////|====> surf. zio_flow
453      !!    + deep  zrecirc  \----|///////|     (+balance of emp)
454      !! 102                      u///////u
455      !!      mid.  recicul    <--|///////|<==== deep  zio_flow
456      !!                     _____|_______|_____ 
457      !!      surf. zio_flow ====>|///////|       
458      !!    (+balance of emp)     |///////|
459      !! 101                      u///////|             
460      !!      mid.  recicul    -->|///////|               Caution: zrecirc split into
461      !!      deep  zrecirc  ---->|///////|                  upper & bottom recirculation
462      !!                   _______|_______|_______
463      !!                     139     140     141 
464      !!
465      !!---------------------------------------------------------------------
466      CHARACTER(len=1), INTENT(in) ::   cd_td   ! ='div' update the divergence
467      !                                         ! ='tra' update the tracers
468      !                                         ! ='spg' update after velocity
469      INTEGER  ::   ji, jj, jk   ! dummy loop indices
470      REAL(wp) ::   zemp_med     ! temporary scalar
471      REAL(wp) ::   zio_flow, zrecirc_upp, zrecirc_mid, zrecirc_bot
472      !!---------------------------------------------------------------------
473      !
474      SELECT CASE( cd_td ) 
475      !                     ! ---------------- !
476      CASE( 'ini' )         !  initialisation  !
477         !                  ! ---------------- !
478         !                                     !** initialization of the velocity
479         hdiv_139_101(:) = 0.e0                     !  139,101 (Atlantic side, south point)
480         hdiv_139_102(:) = 0.e0                     !  139,102 (Atlantic side, north point)
481         hdiv_141_102(:) = 0.e0                     !  141,102 (Med sea  side)
482           
483         !                                     !** imposed transport
484         zio_flow    = 0.8e6                        ! inflow surface  water
485         zrecirc_mid = 0.7e6                        ! middle recirculation water
486         zrecirc_upp = 2.5e6                        ! upper  recirculation water
487         zrecirc_bot = 3.5e6                        ! bottom recirculation water
488         !
489         DO jj = mj0(101), mj1(101)            !** profile of hdiv at 139,101 (Atlantic side, south point)
490            DO ji = mi0(139), mi1(139)         !-----------------------------
491               DO jk = 1, 14                        ! surface in/out flow (Atl -> Med)   (div >0)
492                  hdiv_139_101(jk) = + zio_flow / ( 14. * e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) )
493               END DO
494               DO jk = 15, 20                       ! middle  reciculation (Atl 101 -> Atl 102)   (div >0)   
495                  hdiv_139_101(jk) = + zrecirc_mid / ( 6. * e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) )
496               END DO
497               !                                    ! upper reciculation (Atl 101 -> Atl 101)   (div >0)
498               hdiv_139_101(21) =               + zrecirc_upp / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) )
499               !
500               !                                    ! upper & bottom reciculation (Atl 101 -> Atl 101 & 102)   (div >0)
501               hdiv_139_101(22) = ( zrecirc_bot - zrecirc_upp ) / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) )
502            END DO
503         END DO
504         DO jj = mj0(102), mj1(102)            !** profile of hdiv at 139,102 (Atlantic side, north point)
505            DO ji = mi0(139), mi1(139)         !-----------------------------
506               DO jk = 15, 20                       ! middle reciculation (Atl 101 -> Atl 102)   (div <0)               
507                  hdiv_139_102(jk) = - zrecirc_mid / ( 6. * e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) )
508               END DO
509               !                                    ! outflow of Mediterranean sea + deep recirculation   (div <0)
510               hdiv_139_102(22) = - ( zio_flow + zrecirc_bot ) / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) )
511            END DO
512         END DO
513         DO jj = mj0(102), mj1(102)            !** velocity profile at 141,102  (Med sea side)
514            DO ji = mi0(141), mi1(141)         !------------------------------
515               DO  jk = 1, 14                       ! surface inflow in the Med     (div <0)
516                  hdiv_141_102(jk) = - zio_flow / ( 14. * e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) )
517               END DO
518               !                                    ! deep    outflow toward the Atlantic    (div >0)
519               hdiv_141_102(21)    = + zio_flow / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) )
520            END DO
521         END DO
522         !                  ! ---------------- !
523      CASE( 'div' )         !   update hdivn   ! (call by divcur module)
524         !                  ! ---------=====-- !
525         !                                     !** emp on the Mediterranean Sea  (div >0)
526         zemp_med = 0.e0                       !-------------------------------
527         DO jj = mj0(96), mj1(110)                  ! sum over the Med sea
528            DO ji = mi0(141),mi1(181)
529               zemp_med = zemp_med + emp(ji,jj) * e1t(ji,jj) * e2t(ji,jj) * tmask_i(ji,jj) 
530            END DO
531         END DO
532         DO jj = mj0(96), mj1(96)                   ! minus 2 points in Red Sea
533            DO ji = mi0(148),mi1(148)
534               zemp_med = zemp_med - emp(ji,jj) * e1t(ji,jj) * e2t(ji,jj) * tmask_i(ji,jj)
535            END DO
536            DO ji = mi0(149),mi1(149)
537               zemp_med = zemp_med - emp(ji,jj) * e1t(ji,jj) * e2t(ji,jj) * tmask_i(ji,jj)
538            END DO
539         END DO
540         IF( lk_mpp )   CALL mpp_sum( zemp_med )    ! sum with other processors value
541         zemp_med = zemp_med * 1.e-3                ! convert in m3
542         !
543         !                                     !** Correct hdivn (including emp adjustment)
544         !                                     !-------------------------------------------
545         DO jj = mj0(101), mj1(101)                 !* 139,101 (Atlantic side, south point)
546            DO ji = mi0(139), mi1(139) 
547               hdiv_139_101_kt(:) = hdiv_139_101(:)     
548               DO jk = 1, 14                              ! increase the inflow from the Atlantic   (div >0)
549                  hdiv_139_101_kt(jk) = hdiv_139_101(jk) + zemp_med / ( 14. * e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) )
550               END DO
551               hdivn(ji, jj,:) = hdivn(ji, jj,:) + hdiv_139_101_kt(:)
552            END DO
553         END DO
554         DO jj = mj0(102), mj1(102)                 !* 139,102 (Atlantic side, north point)
555            DO ji = mi0(139), mi1(139) 
556               hdivn(ji,jj,:) = hdivn(ji,jj,:) + hdiv_139_102(:)
557            END DO
558         END DO
559         DO jj = mj0(102), mj1(102)                 !* 141,102 (Med side)
560            DO ji = mi0(141), mi1(141) 
561               hdiv_141_102(:) = hdiv_141_102(:)
562               DO jk = 1, 14                              ! increase the inflow from the Atlantic   (div <0)
563                  hdiv_141_102_kt(jk) = hdiv_141_102(jk) - zemp_med / ( 14. * e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) )
564               END DO
565               hdivn(ji, jj,:) = hdivn(ji, jj,:) + hdiv_141_102_kt(:)
566            END DO
567         END DO
568         !                  ! ---------------- !
569      CASE( 'tra' )         !  update (ta,sa)  ! (call by traadv module)
570         !                  ! --------=======- !
571         !
572         DO jj = mj0(101), mj1(101)            !** 139,101 (Atlantic side, south point)   (div >0)
573            DO ji = mi0(139), mi1(139) 
574               DO jk = 1, jpkm1                         ! surf inflow + mid. & bottom reciculation (from Atlantic)   
575                  tsa(ji,jj,jk,jp_tem) = tsa(ji,jj,jk,jp_tem) - hdiv_139_101_kt(jk) * tsn(ji,jj,jk,jp_tem)
576                  tsa(ji,jj,jk,jp_sal) = tsa(ji,jj,jk,jp_sal) - hdiv_139_101_kt(jk) * tsn(ji,jj,jk,jp_sal)
577               END DO
578            END DO
579         END DO
580         !
581         DO jj = mj0(102), mj1(102)            !** 139,102 (Atlantic side, north point)   (div <0)
582            DO ji = mi0(139), mi1(139) 
583               DO jk = 15, 20                            ! middle  reciculation (Atl 101 -> Atl 102)   (div <0)
584                  tsa(ji,jj,jk,jp_tem) = tsa(ji,jj,jk,jp_tem) - hdiv_139_102(jk) * tsn(ji,jj-1,jk,jp_tem)  ! middle Atlantic recirculation
585                  tsa(ji,jj,jk,jp_sal) = tsa(ji,jj,jk,jp_sal) - hdiv_139_102(jk) * tsn(ji,jj-1,jk,jp_sal)
586               END DO
587               !                                         ! upper & bottom Atl. reciculation (Atl 101 -> Atl 102) - (div <0)
588               !                                         ! deep Med flow                    (Med 102 -> Atl 102) - (div <0)
589               tsa(ji,jj,22,jp_tem) = tsa(ji,jj,22,jp_tem) + hdiv_141_102(21) * tsn(ji+2,jj,21,jp_tem)   &  ! deep Med flow 
590                  &                        + hdiv_139_101(21) * tsn(ji,jj-1,21,jp_tem)   &  ! upper  Atlantic recirculation 
591                  &                        + hdiv_139_101(22) * tsn(ji,jj-1,22,jp_tem)      ! bottom Atlantic recirculation 
592               tsa(ji,jj,22,jp_sal) = tsa(ji,jj,22,jp_sal) + hdiv_141_102(21) * tsn(ji+2,jj,21,jp_sal)   &
593                  &                        + hdiv_139_101(21) * tsn(ji,jj-1,21,jp_sal)   &
594                  &                        + hdiv_139_101(22) * tsn(ji,jj-1,22,jp_sal) 
595            END DO
596         END DO
597         DO jj = mj0(102), mj1(102)                 !* 141,102 (Med side)   (div <0)
598            DO ji = mi0(141), mi1(141) 
599               DO jk = 1, 14                             ! surface flow from Atlantic to Med sea
600                  tsa(ji,jj,jk,jp_tem) = tsa(ji,jj,jk,jp_tem) - hdiv_141_102_kt(jk) * tsn(ji-2,jj-1,jk,jp_tem)
601                  tsa(ji,jj,jk,jp_sal) = tsa(ji,jj,jk,jp_sal) - hdiv_141_102_kt(jk) * tsn(ji-2,jj-1,jk,jp_sal)
602               END DO
603               !                                         ! deeper flow from Med sea to Atlantic
604               tsa(ji,jj,21,jp_tem) = tsa(ji,jj,21,jp_tem) - hdiv_141_102(21) * tsn(ji,jj,21,jp_tem)
605               tsa(ji,jj,21,jp_sal) = tsa(ji,jj,21,jp_sal) - hdiv_141_102(21) * tsn(ji,jj,21,jp_sal)
606            END DO
607         END DO
608         !                  ! ---------------- !
609      CASE( 'spg' )         !  update (ua,va)  ! (call by dynspg module)
610         !                  ! --------=======- !
611         ! at this stage, (ua,va) are the after velocity, not the tendancy
612         ! compute the velocity from the divergence at T-point
613         !
614         DO jj = mj0(101), mj1(101)            !** 139,101 (Atlantic side, south point)
615            DO ji = mi0(139), mi1(139)                    ! div >0 => ua >0, same sign
616               ua(ji,jj,:) = hdiv_139_101_kt(:) / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,:) )   &
617                  &                             * e2u(ji,jj) * fse3u(ji,jj,:)
618            END DO
619         END DO
620         DO jj = mj0(102), mj1(102)            !** 139,102 (Atlantic side, north point)
621            DO ji = mi0(139), mi1(139)                    ! div <0 => ua <0, same sign
622               ua(ji,jj,:) = hdiv_139_102(:) / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,:) )   &
623                  &                          * e2u(ji,jj) * fse3u(ji,jj,:)   
624            END DO
625         END DO
626         DO jj = mj0(102), mj1(102)            !** 140,102 (Med side) (140 not 141 as it is a U-point)
627            DO ji = mi0(140), mi1(140)                    ! div >0 => ua <0, opposite sign
628               ua(ji,jj,:) = - hdiv_141_102(:) / ( e1t(ji+1,jj) * e2t(ji+1,jj) * fse3t(ji+1,jj,:) )   &
629                  &                            * e2u(ji,jj) * fse3u(ji,jj,:)
630            END DO
631         END DO
632         !
633      END SELECT
634      !
635   END SUBROUTINE cla_gibraltar
636
637
638   SUBROUTINE cla_hormuz( cd_td )
639      !! -------------------------------------------------------------------
640      !!                   ***  ROUTINE div_hormuz  ***
641      !!             
642      !! ** Purpose :   update the now horizontal divergence, the tracer
643      !!              tendancyand the after velocity in vicinity of Hormuz
644      !!              strait ( Persian Gulf - Indian ocean ).
645      !!
646      !! ** Method  :   Hormuz strait
647      !!            ______________   
648      !!            |/////|<==      surface inflow
649      !!        94  |/////|     
650      !!            |/////|==>      deep    outflow
651      !!            |_____|_______
652      !!              171    172     
653      !!---------------------------------------------------------------------
654      CHARACTER(len=1), INTENT(in) ::   cd_td   ! ='ini' initialisation
655      !!                                        ! ='div' update the divergence
656      !!                                        ! ='tra' update the tracers
657      !!                                        ! ='spg' update after velocity
658      !!
659      INTEGER  ::   ji, jj, jk   ! dummy loop indices
660      REAL(wp) ::   zio_flow     ! temporary scalar
661      !!---------------------------------------------------------------------
662      !
663      SELECT CASE( cd_td ) 
664      !                     ! ---------------- !
665      CASE( 'ini' )         !  initialisation  !
666         !                  ! ---------------- !
667         !                                     !** profile of horizontal divergence due to cross-land advection
668         zio_flow  = 1.e6                          ! imposed in/out flow
669         !
670         hdiv_172_94(:) = 0.e0         
671         !
672         DO jj = mj0(94), mj1(94)                  ! in/out flow at (i,j) = (172,94)
673            DO ji = mi0(172), mi1(172) 
674               DO jk = 1, 8                            ! surface inflow  (Indian ocean to Persian Gulf) (div<0)
675                  hdiv_172_94(jk) = - ( zio_flow / 8.e0 * e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) )
676               END DO
677               DO jk = 16, 18                          ! deep    outflow (Persian Gulf to Indian ocean) (div>0)
678                  hdiv_172_94(jk) = + ( zio_flow / 3.e0 * e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) )
679               END DO
680            END DO
681         END DO
682         !                                     !** T & S profile in the Hormuz strait (use in deep outflow)
683         !      Temperature       and         Salinity
684         t_171_94_hor(:)  = 0.e0   ;   s_171_94_hor(:)  = 0.e0
685         t_171_94_hor(16) = 18.4   ;   s_171_94_hor(16) = 36.27
686         t_171_94_hor(17) = 17.8   ;   s_171_94_hor(17) = 36.4
687         t_171_94_hor(18) = 16.    ;   s_171_94_hor(18) = 36.27
688         !
689         !                  ! ---------------- !
690      CASE( 'div' )         !   update hdivn   ! (call by divcur module)
691         !                  ! ---------=====-- !
692         !                                   
693         DO jj = mj0(94), mj1(94)              !** 172,94 (Indian ocean side)
694            DO ji = mi0(172), mi1(172) 
695               hdivn(ji,jj,:) = hdivn(ji,jj,:) + hdiv_172_94(:)
696            END DO
697         END DO
698         !                  ! ---------------- !
699      CASE( 'tra' )         !  update (ta,sa)  ! (call by traadv module)
700         !                  ! --------=======- !
701         !                         
702         DO jj = mj0(94), mj1(94)              !** 172,94 (Indian ocean side)
703            DO ji = mi0(172), mi1(172) 
704               DO jk = 1, 8                          ! surface inflow   (Indian ocean to Persian Gulf) (div<0)
705                  tsa(ji,jj,jk,jp_tem) = tsa(ji,jj,jk,jp_tem) - hdiv_172_94(jk) * tsn(ji,jj,jk,jp_tem) 
706                  tsa(ji,jj,jk,jp_sal) = tsa(ji,jj,jk,jp_sal) - hdiv_172_94(jk) * tsn(ji,jj,jk,jp_sal) 
707               END DO
708               DO jk = 16, 18                        ! deep outflow     (Persian Gulf to Indian ocean) (div>0)
709                  tsa(ji,jj,jk,jp_tem) = tsa(ji,jj,jk,jp_tem) - hdiv_172_94(jk) * t_171_94_hor(jk)
710                  tsa(ji,jj,jk,jp_sal) = tsa(ji,jj,jk,jp_sal) - hdiv_172_94(jk) * s_171_94_hor(jk)
711               END DO
712            END DO
713         END DO
714         !                  ! ---------------- !
715      CASE( 'spg' )         !  update (ua,va)  ! (call by dynspg module)
716         !                  ! --------=======- !
717         ! No barotropic flow through Hormuz strait
718         ! at this stage, (ua,va) are the after velocity, not the tendancy
719         ! compute the velocity from the divergence at T-point
720         DO jj = mj0(94), mj1(94)              !** 171,94 (Indian ocean side) (171 not 172 as it is the western U-point)
721            DO ji = mi0(171), mi1(171)                ! div >0 => ua >0, opposite sign
722               ua(ji,jj,:) = - hdiv_172_94(:) / ( e1t(ji+1,jj) * e2t(ji+1,jj) * fse3t(ji+1,jj,:) )   &
723                  &                           * e2u(ji,jj) * fse3u(ji,jj,:)
724            END DO
725         END DO
726         !
727      END SELECT
728      !
729   END SUBROUTINE cla_hormuz
730   
731   !!======================================================================
732END MODULE cla
Note: See TracBrowser for help on using the repository browser.