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_CNRS0_NOC1_LDF/NEMOGCM/NEMO/OPA_SRC/LBC – NEMO

source: branches/2014/dev_CNRS0_NOC1_LDF/NEMOGCM/NEMO/OPA_SRC/LBC/cla.F90 @ 4616

Last change on this file since 4616 was 4616, checked in by gm, 10 years ago

#1260 : see the associated wiki page for explanation

  • Property svn:keywords set to Id
File size: 41.5 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      ! Cross land advection hard coded only for ORCA_R2 with 31 levels linear filtred free surface
185      !
186      IF( cp_cfg /= 'orca'   )   CALL ctl_stop( 'cla_init: Cross Land Advection hard coded for ORCA_R2_L31' )
187      IF( jp_cfg /= 2        )   CALL ctl_stop( 'cla_init: Cross Land Advection hard coded for ORCA_R2_L31' )
188      IF( .NOT.lk_dynspg_flt )   CALL ctl_stop( 'cla_init: Cross Land Advection works only with lk_dynspg_flt=T ' )
189      IF( lk_vvl             )   CALL ctl_stop( 'cla_init: Cross Land Advection does not work with lk_vvl=T option' )
190      IF( jpk /= 31          )   CALL ctl_stop( 'cla_init: Cross Land Advection hard coded for ORCA_R2_L31' )         
191      !
192      !                           ! Allocate arrays for this module
193      ALLOCATE( hdiv_139_101(jpk) , hdiv_139_101_kt(jpk) ,     &    ! Gibraltar
194         &      hdiv_139_102(jpk) ,                            &
195         &      hdiv_141_102(jpk) , hdiv_141_102_kt(jpk) ,     &
196         &      hdiv_161_88 (jpk) , hdiv_161_88_kt (jpk) ,     &    ! Bab-el-Mandeb
197         &      hdiv_161_87 (jpk) ,                            &                     
198         &      hdiv_160_89 (jpk) , hdiv_160_89_kt (jpk) ,     &     ! Hormuz
199         &      hdiv_172_94 (jpk) ,                            &
200         &      t_171_94_hor(jpk) , s_171_94_hor   (jpk) , STAT=ierr )
201      IF( lk_mpp    )   CALL mpp_sum( ierr )
202      IF( ierr /= 0 )   CALL ctl_stop( 'STOP', 'cla_init: unable to allocate arrays' )
203      !
204      !                                        _|_______|_______|_
205      !                                     89  |       |///////| 
206      !                                        _|_______|_______|_
207      ! ------------------------ !          88  |///////|       |
208      !   Bab el Mandeb strait   !             _|_______|_______|_
209      ! ------------------------ !          87  |///////|       |
210      !                                        _|_______|_______|_
211      !                                         |  160  |  161  | 
212      !
213      ! The 6 Bab el Mandeb grid-points must be inside one of the interior of the
214      ! local domain for the cla advection to work properly (i.e. (2:jpim1,2:jpjm1)
215      nbab = 0
216      IF(  ( 1 <= mj0( 88) .AND. mj1( 89) <= jpj ) .AND.    &  !* (161,89), (161,88) and (161,88) on the local pocessor
217         & ( 1 <= mi0(160) .AND. mi1(161) <= jpi )       )    nbab = 1 
218      !
219      ! test if there is no local domain that includes all required grid-points
220      ztemp = REAL( nbab )
221      IF( lk_mpp )   CALL mpp_sum( ztemp )      ! sum with other processors value
222      IF( ztemp == 0 ) THEN                     ! Only 2 points in each direction, this should never be a problem
223         CALL ctl_stop( ' cross land advection at Bab-el_Mandeb does not work with your processor cutting: change it' )
224      ENDIF
225      !                                        ___________________________
226      ! ------------------------ !         102  |       |///////|       |
227      !     Gibraltar strait     !             _|_______|_______|_______|_
228      ! ------------------------ !         101  |       |///////|       |
229      !                                        _|_______|_______|_______|_
230      !                                         |  139  |  140  |  141  |
231      !
232      ! The 6 Gibraltar grid-points must be inside one of the interior of the
233      ! local domain for the cla advection to work properly (i.e. (2:jpim1,2:jpjm1)
234      ngib = 0
235      IF(  ( 2 <= mj0(101) .AND. mj1(102) <= jpjm1 ) .AND.    &  !* (139:141,101:102) on the local pocessor
236         & ( 2 <= mi0(139) .AND. mi1(141) <= jpim1 )       )    ngib = 1 
237      !
238      ! test if there is no local domain that includes all required grid-points
239      ztemp = REAL( ngib )
240      IF( lk_mpp )   CALL mpp_sum( ztemp )      ! sum with other processors value
241      IF( ztemp == 0 ) THEN                     ! 3 points in i-direction, this may be a problem with some cutting
242           CALL ctl_stop( ' cross land advection at Gibraltar does not work with your processor cutting: change it' )
243      ENDIF
244      !                                        _______________
245      ! ------------------------ !          94  |/////|     |
246      !       Hormuz strait      !             _|_____|_____|_
247      ! ------------------------ !                171   172     
248      !           
249      ! The 2 Hormuz grid-points must be inside one of the interior of the
250      ! local domain for the cla advection to work properly (i.e. (2:jpim1,2:jpjm1)
251      nhor = 0
252      IF(    2 <= mj0( 94) .AND. mj1( 94) <= jpjm1  .AND.  & 
253         &   2 <= mi0(171) .AND. mi1(172) <= jpim1         )   nhor = 1 
254      !
255      ! test if there is no local domain that includes all required grid-points
256      ztemp = REAL( nhor )
257      IF( lk_mpp )   CALL mpp_sum( ztemp )      ! sum with other processors value
258      IF( ztemp == 0 ) THEN                     ! 3 points in i-direction, this may be a problem with some cutting
259           CALL ctl_stop( ' cross land advection at Hormuz does not work with your processor cutting: change it' )
260      ENDIF
261      !
262   END SUBROUTINE cla_init
263
264
265   SUBROUTINE cla_bab_el_mandeb( cd_td )
266      !!----------------------------------------------------------------------
267      !!                ***  ROUTINE cla_bab_el_mandeb  ***
268      !!       
269      !! ** Purpose :   update the now horizontal divergence, the tracer tendancy
270      !!              and the after velocity in vicinity of Bab el Mandeb ( Red Sea - Indian ocean).
271      !!
272      !! ** Method  :   compute the exchanges at each side of the strait :
273      !!
274      !!       surf. zio_flow
275      !! (+ balance of emp) /\  |\\\\\\\\\\\|
276      !!                    ||  |\\\\\\\\\\\| 
277      !!    deep zio_flow   ||  |\\\\\\\\\\\| 
278      !!            |  ||   ||  |\\\\\\\\\\\| 
279      !!        89  |  ||   ||  |\\\\\\\\\\\| 
280      !!            |__\/_v_||__|____________
281      !!            !\\\\\\\\\\\|          surf. zio_flow
282      !!            |\\\\\\\\\\\|<===    (+ balance of emp)
283      !!            |\\\\\\\\\\\u
284      !!        88  |\\\\\\\\\\\|<---      deep  zrecirc (upper+deep at 2 different levels)
285      !!            |___________|__________   
286      !!            !\\\\\\\\\\\|         
287      !!            |\\\\\\\\\\\| ---\     deep  zrecirc (upper+deep)
288      !!        87  !\\\\\\\\\\\u ===/   + deep  zio_flow   (all at the same level)
289      !!            !\\\\\\\\\\\| 
290      !!            !___________|__________
291      !!                160         161
292      !!
293      !!----------------------------------------------------------------------
294      CHARACTER(len=1), INTENT(in) ::   cd_td   ! ='div' update the divergence
295      !                                         ! ='tra' update the tracers
296      !                                         ! ='spg' update after velocity
297      INTEGER  ::   ji, jj, jk   ! dummy loop indices
298      REAL(wp) ::   zemp_red     ! temporary scalar
299      REAL(wp) ::   zio_flow, zrecirc_upp, zrecirc_mid, zrecirc_bot
300      !!---------------------------------------------------------------------
301      !
302      SELECT CASE( cd_td ) 
303      !                     ! ---------------- !
304      CASE( 'ini' )         !  initialisation  !
305         !                  ! ---------------- !
306         !                                   
307         zio_flow    = 0.4e6                       ! imposed in/out flow
308         zrecirc_upp = 0.2e6                       ! imposed upper recirculation water
309         zrecirc_bot = 0.5e6                       ! imposed bottom  recirculation water
310
311         hdiv_161_88(:) = 0.e0                     ! (161,88) Gulf of Aden side, north point
312         hdiv_161_87(:) = 0.e0                     ! (161,87) Gulf of Aden side, south point
313         hdiv_160_89(:) = 0.e0                     ! (160,89) Red sea side
314
315         DO jj = mj0(88), mj1(88)              !** profile of hdiv at (161,88)   (Gulf of Aden side, north point)
316            DO ji = mi0(161), mi1(161)         !------------------------------
317               DO jk = 1, 8                        ! surface in/out flow   (Ind -> Red)   (div >0)
318                  hdiv_161_88(jk) = + zio_flow / ( 8. * e1e2t(ji,jj) * fse3t(ji,jj,jk) )
319               END DO
320               !                                   ! recirculation water   (Ind -> Red)   (div >0)
321               hdiv_161_88(20) =                 + zrecirc_upp   / ( e1e2t(ji,jj) * fse3t(ji,jj,20) )
322               hdiv_161_88(21) = + ( zrecirc_bot - zrecirc_upp ) / ( e1e2t(ji,jj) * fse3t(ji,jj,21) )
323            END DO
324         END DO
325         !
326         DO jj = mj0(87), mj1(87)              !** profile of hdiv at (161,88)   (Gulf of Aden side, south point)
327            DO ji = mi0(161), mi1(161)         !------------------------------
328               !                                   ! deep out flow + recirculation   (Red -> Ind)   (div <0)
329               hdiv_161_87(21) = - ( zio_flow + zrecirc_bot ) / ( e1e2t(ji,jj) * fse3t(ji,jj,21) )
330            END DO
331         END DO
332         !
333         DO jj = mj0(89), mj1(89)              !** profile of hdiv at (161,88)   (Red sea side)
334            DO ji = mi0(160), mi1(160)         !------------------------------
335               DO jk = 1, 8                        ! surface inflow    (Ind -> Red)   (div <0)
336                  hdiv_160_89(jk) = - zio_flow / ( 8. * e1e2t(ji,jj) * fse3t(ji,jj,jk) )
337               END DO
338               !                                   ! deep    outflow   (Red -> Ind)   (div >0)
339               hdiv_160_89(16)    = + zio_flow / (      e1e2t(ji,jj) * fse3t(ji,jj,16) )
340            END DO
341         END DO
342         !                  ! ---------------- !
343      CASE( 'div' )         !   update hdivn   ! (call by divcur module)
344         !                  ! ---------=====-- !
345         !                                     !** emp on the Red Sea   (div >0)
346         zemp_red = 0.e0                       !---------------------
347         DO jj = mj0(87), mj1(96)                  ! sum over the Red sea
348            DO ji = mi0(148), mi1(160) 
349               zemp_red = zemp_red + emp(ji,jj) * e1e2t(ji,jj) * tmask_i(ji,jj)
350            END DO
351         END DO
352         IF( lk_mpp )   CALL mpp_sum( zemp_red )   ! sum with other processors value
353         zemp_red = zemp_red * 1.e-3               ! convert in m3
354         !
355         !                                     !** Correct hdivn (including emp adjustment)
356         !                                     !-------------------------------------------
357         DO jj = mj0(88), mj1(88)                  !* profile of hdiv at (161,88)   (Gulf of Aden side, north point)
358            DO ji = mi0(161), mi1(161) 
359               hdiv_161_88_kt(:) = hdiv_161_88(:)
360               DO jk = 1, 8                              ! increase the inflow from the Indian   (div >0)
361                  hdiv_161_88_kt(jk) = hdiv_161_88(jk) + zemp_red / (8. * e2u(ji,jj) * fse3u(ji,jj,jk) )
362               END DO
363               hdivn(ji,jj,:) = hdivn(ji,jj,:) + hdiv_161_88_kt(:)
364            END DO
365         END DO
366         DO jj = mj0(87), mj1(87)                  !* profile of divergence at (161,87)   (Gulf of Aden side, south point)
367            DO ji = mi0(161), mi1(161) 
368               hdivn(ji,jj,:) = hdivn(ji,jj,:) + hdiv_161_87(:)
369            END DO
370         END DO
371         DO jj = mj0(89), mj1(89)                  !* profile of divergence at (160,89)   (Red sea side)
372            DO ji = mi0(160), mi1(160) 
373               hdiv_160_89_kt(:) = hdiv_160_89(:)
374               DO jk = 1, 18                              ! increase the inflow from the Indian   (div <0)
375                  hdiv_160_89_kt(jk) = hdiv_160_89(jk) - zemp_red / (10. * e1v(ji,jj) * fse3v(ji,jj,jk) )
376               END DO
377               hdivn(ji, jj,:) = hdivn(ji, jj,:) + hdiv_160_89_kt(:)
378            END DO
379         END DO
380         !                  ! ---------------- !
381      CASE( 'tra' )         !  update (ta,sa)  ! (call by traadv module)
382         !                  ! --------=======- !
383         !
384         DO jj = mj0(88), mj1(88)              !** (161,88)   (Gulf of Aden side, north point)
385            DO ji = mi0(161), mi1(161) 
386               DO jk = 1, jpkm1                         ! surf inflow + reciculation (from Gulf of Aden)
387                  tsa(ji,jj,jk,jp_tem) = tsa(ji,jj,jk,jp_tem) - hdiv_161_88_kt(jk) * tsn(ji,jj,jk,jp_tem)
388                  tsa(ji,jj,jk,jp_sal) = tsa(ji,jj,jk,jp_sal) - hdiv_161_88_kt(jk) * tsn(ji,jj,jk,jp_sal)
389               END DO
390            END DO
391         END DO
392         DO jj = mj0(87), mj1(87)              !** (161,87)   (Gulf of Aden side, south point)
393            DO ji = mi0(161), mi1(161) 
394               jk =  21                                 ! deep outflow + recirulation (combined flux)
395               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
396                  &                        + hdiv_161_88(21) * tsn(ji  ,jj+1,21,jp_tem)   &  ! deep  recirculation from Gulf of Aden
397                  &                        + hdiv_160_89(16) * tsn(ji-1,jj+2,16,jp_tem)      ! deep inflow from Red sea
398               tsa(ji,jj,jk,jp_sal) = tsa(ji,jj,jk,jp_sal) + hdiv_161_88(20) * tsn(ji  ,jj+1,20,jp_sal)   &
399                  &                        + hdiv_161_88(21) * tsn(ji  ,jj+1,21,jp_sal)   &
400                  &                        + hdiv_160_89(16) * tsn(ji-1,jj+2,16,jp_sal)   
401            END DO
402         END DO
403         DO jj = mj0(89), mj1(89)              !** (161,88)   (Red sea side)
404            DO ji = mi0(160), mi1(160)
405               DO jk = 1, 14                            ! surface inflow (from Gulf of Aden)
406                  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)
407                  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)
408               END DO
409               !                                  ! deep    outflow (from Red sea)
410               tsa(ji,jj,16,jp_tem) = tsa(ji,jj,16,jp_tem) - hdiv_160_89(16) * tsn(ji,jj,16,jp_tem)
411               tsa(ji,jj,16,jp_sal) = tsa(ji,jj,16,jp_sal) - hdiv_160_89(16) * tsn(ji,jj,16,jp_sal)
412            END DO
413         END DO
414         !
415         !                  ! ---------------- !
416      CASE( 'spg' )         !  update (ua,va)  ! (call by dynspg module)
417         !                  ! --------=======- !
418         ! at this stage, (ua,va) are the after velocity, not the tendancy
419         ! compute the velocity from the divergence at T-point
420         !
421         DO jj = mj0(88), mj1(88)              !** (160,88)   (Gulf of Aden side, north point)
422            DO ji = mi0(160), mi1(160)                   ! 160, not 161 as it is a U-point)
423               ua(ji,jj,:) = - hdiv_161_88_kt(:) / ( e1e2t(ji+1,jj) * fse3t(ji+1,jj,:) ) * e2u(ji,jj) * fse3u(ji,jj,:)
424            END DO
425         END DO
426         DO jj = mj0(87), mj1(87)              !** (160,87)   (Gulf of Aden side, south point)
427            DO ji = mi0(160), mi1(160)                   ! 160, not 161 as it is a U-point)
428               ua(ji,jj,:) = - hdiv_161_87(:) / ( e1e2t(ji+1,jj) * fse3t(ji+1,jj,:) ) * 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(:) / ( e1e2t(ji,jj+1) * fse3t(ji,jj+1,:) ) * e1v(ji,jj) * fse3v(ji,jj,:)
434            END DO
435         END DO
436      END SELECT
437      !
438   END SUBROUTINE cla_bab_el_mandeb
439   
440
441   SUBROUTINE cla_gibraltar( cd_td )
442      !! -------------------------------------------------------------------
443      !!                 ***  ROUTINE cla_gibraltar  ***
444      !!       
445      !! ** Purpose :   update the now horizontal divergence, the tracer
446      !!              tendancyand the after velocity in vicinity of Gibraltar
447      !!              strait ( Persian Gulf - Indian ocean ).
448      !!
449      !! ** Method :
450      !!                     _______________________           
451      !!      deep  zio_flow /====|///////|====> surf. zio_flow
452      !!    + deep  zrecirc  \----|///////|     (+balance of emp)
453      !! 102                      u///////u
454      !!      mid.  recicul    <--|///////|<==== deep  zio_flow
455      !!                     _____|_______|_____ 
456      !!      surf. zio_flow ====>|///////|       
457      !!    (+balance of emp)     |///////|
458      !! 101                      u///////|             
459      !!      mid.  recicul    -->|///////|               Caution: zrecirc split into
460      !!      deep  zrecirc  ---->|///////|                  upper & bottom recirculation
461      !!                   _______|_______|_______
462      !!                     139     140     141 
463      !!
464      !!---------------------------------------------------------------------
465      CHARACTER(len=1), INTENT(in) ::   cd_td   ! ='div' update the divergence
466      !                                         ! ='tra' update the tracers
467      !                                         ! ='spg' update after velocity
468      INTEGER  ::   ji, jj, jk   ! dummy loop indices
469      REAL(wp) ::   zemp_med     ! temporary scalar
470      REAL(wp) ::   zio_flow, zrecirc_upp, zrecirc_mid, zrecirc_bot
471      !!---------------------------------------------------------------------
472      !
473      SELECT CASE( cd_td ) 
474      !                     ! ---------------- !
475      CASE( 'ini' )         !  initialisation  !
476         !                  ! ---------------- !
477         !                                     !** initialization of the velocity
478         hdiv_139_101(:) = 0.e0                     !  139,101 (Atlantic side, south point)
479         hdiv_139_102(:) = 0.e0                     !  139,102 (Atlantic side, north point)
480         hdiv_141_102(:) = 0.e0                     !  141,102 (Med sea  side)
481           
482         !                                     !** imposed transport
483         zio_flow    = 0.8e6                        ! inflow surface  water
484         zrecirc_mid = 0.7e6                        ! middle recirculation water
485         zrecirc_upp = 2.5e6                        ! upper  recirculation water
486         zrecirc_bot = 3.5e6                        ! bottom recirculation water
487         !
488         DO jj = mj0(101), mj1(101)            !** profile of hdiv at 139,101 (Atlantic side, south point)
489            DO ji = mi0(139), mi1(139)         !-----------------------------
490               DO jk = 1, 14                        ! surface in/out flow (Atl -> Med)   (div >0)
491                  hdiv_139_101(jk) = + zio_flow / ( 14. * e1e2t(ji,jj) * fse3t(ji,jj,jk) )
492               END DO
493               DO jk = 15, 20                       ! middle  reciculation (Atl 101 -> Atl 102)   (div >0)   
494                  hdiv_139_101(jk) = + zrecirc_mid / ( 6. * e1e2t(ji,jj) * fse3t(ji,jj,jk) )
495               END DO
496               !                                    ! upper reciculation (Atl 101 -> Atl 101)   (div >0)
497               hdiv_139_101(21) =               + zrecirc_upp / ( e1e2t(ji,jj) * fse3t(ji,jj,jk) )
498               !
499               !                                    ! upper & bottom reciculation (Atl 101 -> Atl 101 & 102)   (div >0)
500               hdiv_139_101(22) = ( zrecirc_bot - zrecirc_upp ) / ( e1e2t(ji,jj) * fse3t(ji,jj,jk) )
501            END DO
502         END DO
503         DO jj = mj0(102), mj1(102)            !** profile of hdiv at 139,102 (Atlantic side, north point)
504            DO ji = mi0(139), mi1(139)         !-----------------------------
505               DO jk = 15, 20                       ! middle reciculation (Atl 101 -> Atl 102)   (div <0)               
506                  hdiv_139_102(jk) = - zrecirc_mid / ( 6. * e1e2t(ji,jj) * fse3t(ji,jj,jk) )
507               END DO
508               !                                    ! outflow of Mediterranean sea + deep recirculation   (div <0)
509               hdiv_139_102(22) = - ( zio_flow + zrecirc_bot ) / ( e1e2t(ji,jj) * fse3t(ji,jj,jk) )
510            END DO
511         END DO
512         DO jj = mj0(102), mj1(102)            !** velocity profile at 141,102  (Med sea side)
513            DO ji = mi0(141), mi1(141)         !------------------------------
514               DO  jk = 1, 14                       ! surface inflow in the Med     (div <0)
515                  hdiv_141_102(jk) = - zio_flow / ( 14. * e1e2t(ji,jj) * fse3t(ji,jj,jk) )
516               END DO
517               !                                    ! deep    outflow toward the Atlantic    (div >0)
518               hdiv_141_102(21)    = + zio_flow / ( e1e2t(ji,jj) * fse3t(ji,jj,jk) )
519            END DO
520         END DO
521         !                  ! ---------------- !
522      CASE( 'div' )         !   update hdivn   ! (call by divcur module)
523         !                  ! ---------=====-- !
524         !                                     !** emp on the Mediterranean Sea  (div >0)
525         zemp_med = 0.e0                       !-------------------------------
526         DO jj = mj0(96), mj1(110)                  ! sum over the Med sea
527            DO ji = mi0(141),mi1(181)
528               zemp_med = zemp_med + emp(ji,jj) * e1e2t(ji,jj) * tmask_i(ji,jj) 
529            END DO
530         END DO
531         DO jj = mj0(96), mj1(96)                   ! minus 2 points in Red Sea
532            DO ji = mi0(148),mi1(148)
533               zemp_med = zemp_med - emp(ji,jj) * e1e2t(ji,jj) * tmask_i(ji,jj)
534            END DO
535            DO ji = mi0(149),mi1(149)
536               zemp_med = zemp_med - emp(ji,jj) * e1e2t(ji,jj) * tmask_i(ji,jj)
537            END DO
538         END DO
539         IF( lk_mpp )   CALL mpp_sum( zemp_med )    ! sum with other processors value
540         zemp_med = zemp_med * 1.e-3                ! convert in m3
541         !
542         !                                     !** Correct hdivn (including emp adjustment)
543         !                                     !-------------------------------------------
544         DO jj = mj0(101), mj1(101)                 !* 139,101 (Atlantic side, south point)
545            DO ji = mi0(139), mi1(139) 
546               hdiv_139_101_kt(:) = hdiv_139_101(:)     
547               DO jk = 1, 14                              ! increase the inflow from the Atlantic   (div >0)
548                  hdiv_139_101_kt(jk) = hdiv_139_101(jk) + zemp_med / ( 14. * e1e2t(ji,jj) * fse3t(ji,jj,jk) )
549               END DO
550               hdivn(ji, jj,:) = hdivn(ji, jj,:) + hdiv_139_101_kt(:)
551            END DO
552         END DO
553         DO jj = mj0(102), mj1(102)                 !* 139,102 (Atlantic side, north point)
554            DO ji = mi0(139), mi1(139) 
555               hdivn(ji,jj,:) = hdivn(ji,jj,:) + hdiv_139_102(:)
556            END DO
557         END DO
558         DO jj = mj0(102), mj1(102)                 !* 141,102 (Med side)
559            DO ji = mi0(141), mi1(141) 
560               hdiv_141_102(:) = hdiv_141_102(:)
561               DO jk = 1, 14                              ! increase the inflow from the Atlantic   (div <0)
562                  hdiv_141_102_kt(jk) = hdiv_141_102(jk) - zemp_med / ( 14. * e1e2t(ji,jj) * fse3t(ji,jj,jk) )
563               END DO
564               hdivn(ji, jj,:) = hdivn(ji, jj,:) + hdiv_141_102_kt(:)
565            END DO
566         END DO
567         !                  ! ---------------- !
568      CASE( 'tra' )         !  update (ta,sa)  ! (call by traadv module)
569         !                  ! --------=======- !
570         !
571         DO jj = mj0(101), mj1(101)            !** 139,101 (Atlantic side, south point)   (div >0)
572            DO ji = mi0(139), mi1(139) 
573               DO jk = 1, jpkm1                         ! surf inflow + mid. & bottom reciculation (from Atlantic)   
574                  tsa(ji,jj,jk,jp_tem) = tsa(ji,jj,jk,jp_tem) - hdiv_139_101_kt(jk) * tsn(ji,jj,jk,jp_tem)
575                  tsa(ji,jj,jk,jp_sal) = tsa(ji,jj,jk,jp_sal) - hdiv_139_101_kt(jk) * tsn(ji,jj,jk,jp_sal)
576               END DO
577            END DO
578         END DO
579         !
580         DO jj = mj0(102), mj1(102)            !** 139,102 (Atlantic side, north point)   (div <0)
581            DO ji = mi0(139), mi1(139) 
582               DO jk = 15, 20                            ! middle  reciculation (Atl 101 -> Atl 102)   (div <0)
583                  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
584                  tsa(ji,jj,jk,jp_sal) = tsa(ji,jj,jk,jp_sal) - hdiv_139_102(jk) * tsn(ji,jj-1,jk,jp_sal)
585               END DO
586               !                                         ! upper & bottom Atl. reciculation (Atl 101 -> Atl 102) - (div <0)
587               !                                         ! deep Med flow                    (Med 102 -> Atl 102) - (div <0)
588               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 
589                  &                        + hdiv_139_101(21) * tsn(ji,jj-1,21,jp_tem)   &  ! upper  Atlantic recirculation 
590                  &                        + hdiv_139_101(22) * tsn(ji,jj-1,22,jp_tem)      ! bottom Atlantic recirculation 
591               tsa(ji,jj,22,jp_sal) = tsa(ji,jj,22,jp_sal) + hdiv_141_102(21) * tsn(ji+2,jj,21,jp_sal)   &
592                  &                        + hdiv_139_101(21) * tsn(ji,jj-1,21,jp_sal)   &
593                  &                        + hdiv_139_101(22) * tsn(ji,jj-1,22,jp_sal) 
594            END DO
595         END DO
596         DO jj = mj0(102), mj1(102)                 !* 141,102 (Med side)   (div <0)
597            DO ji = mi0(141), mi1(141) 
598               DO jk = 1, 14                             ! surface flow from Atlantic to Med sea
599                  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)
600                  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)
601               END DO
602               !                                         ! deeper flow from Med sea to Atlantic
603               tsa(ji,jj,21,jp_tem) = tsa(ji,jj,21,jp_tem) - hdiv_141_102(21) * tsn(ji,jj,21,jp_tem)
604               tsa(ji,jj,21,jp_sal) = tsa(ji,jj,21,jp_sal) - hdiv_141_102(21) * tsn(ji,jj,21,jp_sal)
605            END DO
606         END DO
607         !                  ! ---------------- !
608      CASE( 'spg' )         !  update (ua,va)  ! (call by dynspg module)
609         !                  ! --------=======- !
610         ! at this stage, (ua,va) are the after velocity, not the tendancy
611         ! compute the velocity from the divergence at T-point
612         !
613         DO jj = mj0(101), mj1(101)            !** 139,101 (Atlantic side, south point)
614            DO ji = mi0(139), mi1(139)                    ! div >0 => ua >0, same sign
615               ua(ji,jj,:) = hdiv_139_101_kt(:) / ( e1e2t(ji,jj) * fse3t(ji,jj,:) ) * e2u(ji,jj) * fse3u(ji,jj,:)
616            END DO
617         END DO
618         DO jj = mj0(102), mj1(102)            !** 139,102 (Atlantic side, north point)
619            DO ji = mi0(139), mi1(139)                    ! div <0 => ua <0, same sign
620               ua(ji,jj,:) = hdiv_139_102(:) / ( e1e2t(ji,jj) * fse3t(ji,jj,:) ) * e2u(ji,jj) * fse3u(ji,jj,:)   
621            END DO
622         END DO
623         DO jj = mj0(102), mj1(102)            !** 140,102 (Med side) (140 not 141 as it is a U-point)
624            DO ji = mi0(140), mi1(140)                    ! div >0 => ua <0, opposite sign
625               ua(ji,jj,:) = - hdiv_141_102(:) / ( e1e2t(ji+1,jj) * fse3t(ji+1,jj,:) ) * e2u(ji,jj) * fse3u(ji,jj,:)
626            END DO
627         END DO
628         !
629      END SELECT
630      !
631   END SUBROUTINE cla_gibraltar
632
633
634   SUBROUTINE cla_hormuz( cd_td )
635      !! -------------------------------------------------------------------
636      !!                   ***  ROUTINE div_hormuz  ***
637      !!             
638      !! ** Purpose :   update the now horizontal divergence, the tracer
639      !!              tendancyand the after velocity in vicinity of Hormuz
640      !!              strait ( Persian Gulf - Indian ocean ).
641      !!
642      !! ** Method  :   Hormuz strait
643      !!            ______________   
644      !!            |/////|<==      surface inflow
645      !!        94  |/////|     
646      !!            |/////|==>      deep    outflow
647      !!            |_____|_______
648      !!              171    172     
649      !!---------------------------------------------------------------------
650      CHARACTER(len=1), INTENT(in) ::   cd_td   ! ='ini' initialisation
651      !!                                        ! ='div' update the divergence
652      !!                                        ! ='tra' update the tracers
653      !!                                        ! ='spg' update after velocity
654      !!
655      INTEGER  ::   ji, jj, jk   ! dummy loop indices
656      REAL(wp) ::   zio_flow     ! temporary scalar
657      !!---------------------------------------------------------------------
658      !
659      SELECT CASE( cd_td ) 
660      !                     ! ---------------- !
661      CASE( 'ini' )         !  initialisation  !
662         !                  ! ---------------- !
663         !                                     !** profile of horizontal divergence due to cross-land advection
664         zio_flow  = 1.e6                          ! imposed in/out flow
665         !
666         hdiv_172_94(:) = 0.e0         
667         !
668         DO jj = mj0(94), mj1(94)                  ! in/out flow at (i,j) = (172,94)
669            DO ji = mi0(172), mi1(172) 
670               DO jk = 1, 8                            ! surface inflow  (Indian ocean to Persian Gulf) (div<0)
671                  hdiv_172_94(jk) = - ( zio_flow / 8.e0 * e1e2t(ji,jj) * fse3t(ji,jj,jk) )
672               END DO
673               DO jk = 16, 18                          ! deep    outflow (Persian Gulf to Indian ocean) (div>0)
674                  hdiv_172_94(jk) = + ( zio_flow / 3.e0 * e1e2t(ji,jj) * fse3t(ji,jj,jk) )
675               END DO
676            END DO
677         END DO
678         !                                     !** T & S profile in the Hormuz strait (use in deep outflow)
679         !      Temperature       and         Salinity
680         t_171_94_hor(:)  = 0.e0   ;   s_171_94_hor(:)  = 0.e0
681         t_171_94_hor(16) = 18.4   ;   s_171_94_hor(16) = 36.27
682         t_171_94_hor(17) = 17.8   ;   s_171_94_hor(17) = 36.4
683         t_171_94_hor(18) = 16.    ;   s_171_94_hor(18) = 36.27
684         !
685         !                  ! ---------------- !
686      CASE( 'div' )         !   update hdivn   ! (call by divcur module)
687         !                  ! ---------=====-- !
688         !                                   
689         DO jj = mj0(94), mj1(94)              !** 172,94 (Indian ocean side)
690            DO ji = mi0(172), mi1(172) 
691               hdivn(ji,jj,:) = hdivn(ji,jj,:) + hdiv_172_94(:)
692            END DO
693         END DO
694         !                  ! ---------------- !
695      CASE( 'tra' )         !  update (ta,sa)  ! (call by traadv module)
696         !                  ! --------=======- !
697         !                         
698         DO jj = mj0(94), mj1(94)              !** 172,94 (Indian ocean side)
699            DO ji = mi0(172), mi1(172) 
700               DO jk = 1, 8                          ! surface inflow   (Indian ocean to Persian Gulf) (div<0)
701                  tsa(ji,jj,jk,jp_tem) = tsa(ji,jj,jk,jp_tem) - hdiv_172_94(jk) * tsn(ji,jj,jk,jp_tem) 
702                  tsa(ji,jj,jk,jp_sal) = tsa(ji,jj,jk,jp_sal) - hdiv_172_94(jk) * tsn(ji,jj,jk,jp_sal) 
703               END DO
704               DO jk = 16, 18                        ! deep outflow     (Persian Gulf to Indian ocean) (div>0)
705                  tsa(ji,jj,jk,jp_tem) = tsa(ji,jj,jk,jp_tem) - hdiv_172_94(jk) * t_171_94_hor(jk)
706                  tsa(ji,jj,jk,jp_sal) = tsa(ji,jj,jk,jp_sal) - hdiv_172_94(jk) * s_171_94_hor(jk)
707               END DO
708            END DO
709         END DO
710         !                  ! ---------------- !
711      CASE( 'spg' )         !  update (ua,va)  ! (call by dynspg module)
712         !                  ! --------=======- !
713         ! No barotropic flow through Hormuz strait
714         ! at this stage, (ua,va) are the after velocity, not the tendancy
715         ! compute the velocity from the divergence at T-point
716         DO jj = mj0(94), mj1(94)              !** 171,94 (Indian ocean side) (171 not 172 as it is the western U-point)
717            DO ji = mi0(171), mi1(171)                ! div >0 => ua >0, opposite sign
718               ua(ji,jj,:) = - hdiv_172_94(:) / ( e1e2t(ji+1,jj) * fse3t(ji+1,jj,:) ) * e2u(ji,jj) * fse3u(ji,jj,:)
719            END DO
720         END DO
721         !
722      END SELECT
723      !
724   END SUBROUTINE cla_hormuz
725   
726   !!======================================================================
727END MODULE cla
Note: See TracBrowser for help on using the repository browser.