diaar5.F90 14 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323
  1. MODULE diaar5
  2. !!======================================================================
  3. !! *** MODULE diaar5 ***
  4. !! AR5 diagnostics
  5. !!======================================================================
  6. !! History : 3.2 ! 2009-11 (S. Masson) Original code
  7. !! 3.3 ! 2010-10 (C. Ethe, G. Madec) reorganisation of initialisation phase + merge TRC-TRA
  8. !!----------------------------------------------------------------------
  9. #if defined key_diaar5 || defined key_esopa
  10. !!----------------------------------------------------------------------
  11. !! 'key_diaar5' : activate ar5 diagnotics
  12. !!----------------------------------------------------------------------
  13. !! dia_ar5 : AR5 diagnostics
  14. !! dia_ar5_init : initialisation of AR5 diagnostics
  15. !!----------------------------------------------------------------------
  16. USE oce ! ocean dynamics and active tracers
  17. USE dom_oce ! ocean space and time domain
  18. USE eosbn2 ! equation of state (eos_bn2 routine)
  19. USE lib_mpp ! distribued memory computing library
  20. USE iom ! I/O manager library
  21. USE timing ! preformance summary
  22. USE wrk_nemo ! working arrays
  23. USE fldread ! type FLD_N
  24. USE phycst ! physical constant
  25. USE in_out_manager ! I/O manager
  26. USE zdfddm
  27. USE zdf_oce
  28. IMPLICIT NONE
  29. PRIVATE
  30. PUBLIC dia_ar5 ! routine called in step.F90 module
  31. PUBLIC dia_ar5_init ! routine called in opa.F90 module
  32. PUBLIC dia_ar5_alloc ! routine called in nemogcm.F90 module
  33. LOGICAL, PUBLIC, PARAMETER :: lk_diaar5 = .TRUE. ! coupled flag
  34. REAL(wp) :: vol0 ! ocean volume (interior domain)
  35. REAL(wp) :: area_tot ! total ocean surface (interior domain)
  36. REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,: ) :: area ! cell surface (interior domain)
  37. REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,: ) :: thick0 ! ocean thickness (interior domain)
  38. REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: sn0 ! initial salinity
  39. !! * Substitutions
  40. # include "domzgr_substitute.h90"
  41. # include "zdfddm_substitute.h90"
  42. !!----------------------------------------------------------------------
  43. !! NEMO/OPA 3.3 , NEMO Consortium (2010)
  44. !! $Id: diaar5.F90 4990 2014-12-15 16:42:49Z timgraham $
  45. !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
  46. !!----------------------------------------------------------------------
  47. CONTAINS
  48. FUNCTION dia_ar5_alloc()
  49. !!----------------------------------------------------------------------
  50. !! *** ROUTINE dia_ar5_alloc ***
  51. !!----------------------------------------------------------------------
  52. INTEGER :: dia_ar5_alloc
  53. !!----------------------------------------------------------------------
  54. !
  55. ALLOCATE( area(jpi,jpj), thick0(jpi,jpj) , sn0(jpi,jpj,jpk) , STAT=dia_ar5_alloc )
  56. !
  57. IF( lk_mpp ) CALL mpp_sum ( dia_ar5_alloc )
  58. IF( dia_ar5_alloc /= 0 ) CALL ctl_warn('dia_ar5_alloc: failed to allocate arrays')
  59. !
  60. END FUNCTION dia_ar5_alloc
  61. SUBROUTINE dia_ar5( kt )
  62. !!----------------------------------------------------------------------
  63. !! *** ROUTINE dia_ar5 ***
  64. !!
  65. !! ** Purpose : compute and output some AR5 diagnostics
  66. !!----------------------------------------------------------------------
  67. !
  68. INTEGER, INTENT( in ) :: kt ! ocean time-step index
  69. !
  70. INTEGER :: ji, jj, jk ! dummy loop arguments
  71. REAL(wp) :: zvolssh, zvol, zssh_steric, zztmp, zarho, ztemp, zsal, zmass
  72. REAL(wp) :: zaw, zbw, zrw
  73. !
  74. REAL(wp), POINTER, DIMENSION(:,:) :: zarea_ssh , zbotpres ! 2D workspace
  75. REAL(wp), POINTER, DIMENSION(:,:) :: zpe ! 2D workspace
  76. REAL(wp), POINTER, DIMENSION(:,:,:) :: zrhd , zrhop ! 3D workspace
  77. REAL(wp), POINTER, DIMENSION(:,:,:,:) :: ztsn ! 4D workspace
  78. !!--------------------------------------------------------------------
  79. IF( nn_timing == 1 ) CALL timing_start('dia_ar5')
  80. !Call to init moved to here so that we can call iom_use in the
  81. !initialisation
  82. IF( kt == nit000 ) CALL dia_ar5_init
  83. CALL wrk_alloc( jpi , jpj , zarea_ssh , zbotpres, zpe )
  84. CALL wrk_alloc( jpi , jpj , jpk , zrhd , zrhop )
  85. CALL wrk_alloc( jpi , jpj , jpk , jpts , ztsn )
  86. zarea_ssh(:,:) = area(:,:) * sshn(:,:)
  87. ! ! total volume of liquid seawater
  88. zvolssh = SUM( zarea_ssh(:,:) )
  89. IF( lk_mpp ) CALL mpp_sum( zvolssh )
  90. zvol = vol0 + zvolssh
  91. CALL iom_put( 'voltot', zvol )
  92. CALL iom_put( 'sshtot', zvolssh / area_tot )
  93. CALL iom_put( 'sshdyn', sshn(:,:) - (zvolssh / area_tot) )
  94. !
  95. IF( iom_use('sshthster')) THEN
  96. ztsn(:,:,:,jp_tem) = tsn(:,:,:,jp_tem) ! thermosteric ssh
  97. ztsn(:,:,:,jp_sal) = sn0(:,:,:)
  98. CALL eos( ztsn, zrhd, fsdept_n(:,:,:) ) ! now in situ density using initial salinity
  99. !
  100. zbotpres(:,:) = 0._wp ! no atmospheric surface pressure, levitating sea-ice
  101. DO jk = 1, jpkm1
  102. zbotpres(:,:) = zbotpres(:,:) + fse3t(:,:,jk) * zrhd(:,:,jk)
  103. END DO
  104. IF( .NOT.lk_vvl ) THEN
  105. IF ( ln_isfcav ) THEN
  106. DO ji=1,jpi
  107. DO jj=1,jpj
  108. zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,mikt(ji,jj)) + riceload(ji,jj)
  109. END DO
  110. END DO
  111. ELSE
  112. zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1)
  113. END IF
  114. END IF
  115. !
  116. zarho = SUM( area(:,:) * zbotpres(:,:) )
  117. IF( lk_mpp ) CALL mpp_sum( zarho )
  118. zssh_steric = - zarho / area_tot
  119. CALL iom_put( 'sshthster', zssh_steric )
  120. ENDIF
  121. ! ! steric sea surface height
  122. CALL eos( tsn, zrhd, zrhop, fsdept_n(:,:,:) ) ! now in situ and potential density
  123. zrhop(:,:,jpk) = 0._wp
  124. CALL iom_put( 'rhop', zrhop )
  125. !
  126. zbotpres(:,:) = 0._wp ! no atmospheric surface pressure, levitating sea-ice
  127. DO jk = 1, jpkm1
  128. zbotpres(:,:) = zbotpres(:,:) + fse3t(:,:,jk) * zrhd(:,:,jk)
  129. END DO
  130. IF( .NOT.lk_vvl ) THEN
  131. IF ( ln_isfcav ) THEN
  132. DO ji=1,jpi
  133. DO jj=1,jpj
  134. zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,mikt(ji,jj)) + riceload(ji,jj)
  135. END DO
  136. END DO
  137. ELSE
  138. zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1)
  139. END IF
  140. END IF
  141. !
  142. zarho = SUM( area(:,:) * zbotpres(:,:) )
  143. IF( lk_mpp ) CALL mpp_sum( zarho )
  144. zssh_steric = - zarho / area_tot
  145. CALL iom_put( 'sshsteric', zssh_steric )
  146. ! ! ocean bottom pressure
  147. zztmp = rau0 * grav * 1.e-4_wp ! recover pressure from pressure anomaly and cover to dbar = 1.e4 Pa
  148. zbotpres(:,:) = zztmp * ( zbotpres(:,:) + sshn(:,:) + thick0(:,:) )
  149. CALL iom_put( 'botpres', zbotpres )
  150. ! ! Mean density anomalie, temperature and salinity
  151. ztemp = 0._wp
  152. zsal = 0._wp
  153. DO jk = 1, jpkm1
  154. DO jj = 1, jpj
  155. DO ji = 1, jpi
  156. zztmp = area(ji,jj) * fse3t(ji,jj,jk)
  157. ztemp = ztemp + zztmp * tsn(ji,jj,jk,jp_tem)
  158. zsal = zsal + zztmp * tsn(ji,jj,jk,jp_sal)
  159. END DO
  160. END DO
  161. END DO
  162. IF( .NOT.lk_vvl ) THEN
  163. IF ( ln_isfcav ) THEN
  164. DO ji=1,jpi
  165. DO jj=1,jpj
  166. ztemp = ztemp + zarea_ssh(ji,jj) * tsn(ji,jj,mikt(ji,jj),jp_tem)
  167. zsal = zsal + zarea_ssh(ji,jj) * tsn(ji,jj,mikt(ji,jj),jp_sal)
  168. END DO
  169. END DO
  170. ELSE
  171. ztemp = ztemp + SUM( zarea_ssh(:,:) * tsn(:,:,1,jp_tem) )
  172. zsal = zsal + SUM( zarea_ssh(:,:) * tsn(:,:,1,jp_sal) )
  173. END IF
  174. ENDIF
  175. IF( lk_mpp ) THEN
  176. CALL mpp_sum( ztemp )
  177. CALL mpp_sum( zsal )
  178. END IF
  179. !
  180. zmass = rau0 * ( zarho + zvol ) ! total mass of liquid seawater
  181. ztemp = ztemp / zvol ! potential temperature in liquid seawater
  182. zsal = zsal / zvol ! Salinity of liquid seawater
  183. !
  184. CALL iom_put( 'masstot', zmass )
  185. CALL iom_put( 'temptot', ztemp )
  186. CALL iom_put( 'saltot' , zsal )
  187. IF( iom_use( 'tnpeo' )) THEN
  188. ! Work done against stratification by vertical mixing
  189. ! Exclude points where rn2 is negative as convection kicks in here and
  190. ! work is not being done against stratification
  191. zpe(:,:) = 0._wp
  192. IF( lk_zdfddm ) THEN
  193. DO jk = 2, jpk
  194. DO jj = 1, jpj
  195. DO ji = 1, jpi
  196. IF( rn2(ji,jj,jk) > 0._wp ) THEN
  197. zrw = ( fsdepw(ji,jj,jk ) - fsdept(ji,jj,jk) ) &
  198. & / ( fsdept(ji,jj,jk-1) - fsdept(ji,jj,jk) )
  199. !
  200. zaw = rab_n(ji,jj,jk,jp_tem) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_tem)* zrw
  201. zbw = rab_n(ji,jj,jk,jp_sal) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_sal)* zrw
  202. !
  203. zpe(ji, jj) = zpe(ji, jj) &
  204. & - grav * ( avt(ji,jj,jk) * zaw * (tsn(ji,jj,jk-1,jp_tem) - tsn(ji,jj,jk,jp_tem) ) &
  205. & - fsavs(ji,jj,jk) * zbw * (tsn(ji,jj,jk-1,jp_sal) - tsn(ji,jj,jk,jp_sal) ) )
  206. ENDIF
  207. END DO
  208. END DO
  209. END DO
  210. ELSE
  211. DO jk = 1, jpk
  212. DO ji = 1, jpi
  213. DO jj = 1, jpj
  214. zpe(ji,jj) = zpe(ji,jj) + avt(ji, jj, jk) * MIN(0._wp,rn2(ji, jj, jk)) * rau0 * fse3w(ji, jj, jk)
  215. END DO
  216. END DO
  217. END DO
  218. ENDIF
  219. CALL lbc_lnk(zpe, 'T', 1._wp)
  220. CALL iom_put( 'tnpeo', zpe )
  221. ENDIF
  222. !
  223. CALL wrk_dealloc( jpi , jpj , zarea_ssh , zbotpres, zpe )
  224. CALL wrk_dealloc( jpi , jpj , jpk , zrhd , zrhop )
  225. CALL wrk_dealloc( jpi , jpj , jpk , jpts , ztsn )
  226. !
  227. IF( nn_timing == 1 ) CALL timing_stop('dia_ar5')
  228. !
  229. END SUBROUTINE dia_ar5
  230. SUBROUTINE dia_ar5_init
  231. !!----------------------------------------------------------------------
  232. !! *** ROUTINE dia_ar5_init ***
  233. !!
  234. !! ** Purpose : initialization for AR5 diagnostic computation
  235. !!----------------------------------------------------------------------
  236. INTEGER :: inum
  237. INTEGER :: ik
  238. INTEGER :: ji, jj, jk ! dummy loop indices
  239. REAL(wp) :: zztmp
  240. REAL(wp), POINTER, DIMENSION(:,:,:,:) :: zsaldta ! Jan/Dec levitus salinity
  241. !
  242. !!----------------------------------------------------------------------
  243. !
  244. IF( nn_timing == 1 ) CALL timing_start('dia_ar5_init')
  245. !
  246. CALL wrk_alloc( jpi, jpj, jpk, 2, zsaldta )
  247. ! ! allocate dia_ar5 arrays
  248. IF( dia_ar5_alloc() /= 0 ) CALL ctl_stop( 'STOP', 'dia_ar5_init : unable to allocate arrays' )
  249. area(:,:) = e1t(:,:) * e2t(:,:) * tmask_i(:,:)
  250. area_tot = SUM( area(:,:) ) ; IF( lk_mpp ) CALL mpp_sum( area_tot )
  251. vol0 = 0._wp
  252. thick0(:,:) = 0._wp
  253. DO jk = 1, jpkm1
  254. vol0 = vol0 + SUM( area (:,:) * tmask(:,:,jk) * e3t_0(:,:,jk) )
  255. thick0(:,:) = thick0(:,:) + tmask_i(:,:) * tmask(:,:,jk) * e3t_0(:,:,jk)
  256. END DO
  257. IF( lk_mpp ) CALL mpp_sum( vol0 )
  258. IF( iom_use('sshthster')) THEN
  259. CALL iom_open ( 'sali_ref_clim_monthly', inum )
  260. CALL iom_get ( inum, jpdom_data, 'vosaline' , zsaldta(:,:,:,1), 1 )
  261. CALL iom_get ( inum, jpdom_data, 'vosaline' , zsaldta(:,:,:,2), 12 )
  262. CALL iom_close( inum )
  263. sn0(:,:,:) = 0.5_wp * ( zsaldta(:,:,:,1) + zsaldta(:,:,:,2) )
  264. sn0(:,:,:) = sn0(:,:,:) * tmask(:,:,:)
  265. IF( ln_zps ) THEN ! z-coord. partial steps
  266. DO jj = 1, jpj ! interpolation of salinity at the last ocean level (i.e. the partial step)
  267. DO ji = 1, jpi
  268. ik = mbkt(ji,jj)
  269. IF( ik > 1 ) THEN
  270. zztmp = ( gdept_1d(ik) - gdept_0(ji,jj,ik) ) / ( gdept_1d(ik) - gdept_1d(ik-1) )
  271. sn0(ji,jj,ik) = ( 1._wp - zztmp ) * sn0(ji,jj,ik) + zztmp * sn0(ji,jj,ik-1)
  272. ENDIF
  273. END DO
  274. END DO
  275. ENDIF
  276. ENDIF
  277. !
  278. CALL wrk_dealloc( jpi, jpj, jpk, 2, zsaldta )
  279. !
  280. IF( nn_timing == 1 ) CALL timing_stop('dia_ar5_init')
  281. !
  282. END SUBROUTINE dia_ar5_init
  283. #else
  284. !!----------------------------------------------------------------------
  285. !! Default option : NO diaar5
  286. !!----------------------------------------------------------------------
  287. LOGICAL, PUBLIC, PARAMETER :: lk_diaar5 = .FALSE. ! coupled flag
  288. CONTAINS
  289. SUBROUTINE dia_ar5_init ! Dummy routine
  290. END SUBROUTINE dia_ar5_init
  291. SUBROUTINE dia_ar5( kt ) ! Empty routine
  292. INTEGER :: kt
  293. WRITE(*,*) 'dia_ar5: You should not have seen this print! error?', kt
  294. END SUBROUTINE dia_ar5
  295. #endif
  296. !!======================================================================
  297. END MODULE diaar5