cla.F90 42 KB

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