diaar5.F90 17 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396
  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 lib_fortran ! Fortran routines library
  21. USE iom ! I/O manager library
  22. USE timing ! preformance summary
  23. USE wrk_nemo ! working arrays
  24. USE fldread ! type FLD_N
  25. USE phycst ! physical constant
  26. USE in_out_manager ! I/O manager
  27. USE zdfddm
  28. USE zdf_oce
  29. IMPLICIT NONE
  30. PRIVATE
  31. PUBLIC dia_ar5 ! routine called in step.F90 module
  32. PUBLIC dia_ar5_init ! routine called in opa.F90 module
  33. PUBLIC dia_ar5_alloc ! routine called in nemogcm.F90 module
  34. LOGICAL, PUBLIC, PARAMETER :: lk_diaar5 = .TRUE. ! coupled flag
  35. REAL(wp) :: vol0 ! ocean volume (interior domain)
  36. REAL(wp) :: area_tot ! total ocean surface (interior domain)
  37. REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,: ) :: area ! cell surface (interior domain)
  38. REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,: ) :: thick0 ! ocean thickness (interior domain)
  39. REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: sn0 ! initial salinity
  40. !! * Substitutions
  41. # include "domzgr_substitute.h90"
  42. # include "zdfddm_substitute.h90"
  43. !!----------------------------------------------------------------------
  44. !! NEMO/OPA 3.3 , NEMO Consortium (2010)
  45. !! $Id$
  46. !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
  47. !!----------------------------------------------------------------------
  48. CONTAINS
  49. FUNCTION dia_ar5_alloc()
  50. !!----------------------------------------------------------------------
  51. !! *** ROUTINE dia_ar5_alloc ***
  52. !!----------------------------------------------------------------------
  53. INTEGER :: dia_ar5_alloc
  54. !!----------------------------------------------------------------------
  55. !
  56. ALLOCATE( area(jpi,jpj), thick0(jpi,jpj), sn0(jpi,jpj,jpk) , STAT=dia_ar5_alloc )
  57. !
  58. IF( lk_mpp ) CALL mpp_sum ( dia_ar5_alloc )
  59. IF( dia_ar5_alloc /= 0 ) CALL ctl_warn('dia_ar5_alloc: failed to allocate arrays')
  60. !
  61. END FUNCTION dia_ar5_alloc
  62. SUBROUTINE dia_ar5( kt )
  63. !!----------------------------------------------------------------------
  64. !! *** ROUTINE dia_ar5 ***
  65. !!
  66. !! ** Purpose : compute and output some AR5 diagnostics
  67. !!----------------------------------------------------------------------
  68. !
  69. INTEGER, INTENT( in ) :: kt ! ocean time-step index
  70. !
  71. INTEGER :: ji, jj, jk, iks, ikb ! dummy loop arguments
  72. REAL(wp) :: zvolssh, zvol, zssh_steric, zztmp, zarho, ztemp, zsal, zmass
  73. REAL(wp) :: zaw, zbw, zrw, zsst, zsss
  74. !
  75. REAL(wp), POINTER, DIMENSION(:,:) :: zarea_ssh, zbotpres, zpe, z2d ! 2D workspace
  76. REAL(wp), POINTER, DIMENSION(:,:,:) :: zrhd , zrhop, z3d, ztpot ! 3D workspace
  77. REAL(wp), POINTER, DIMENSION(:,:,:,:) :: ztsn ! 4D workspace
  78. !!--------------------------------------------------------------------
  79. IF( nn_timing == 1 ) CALL timing_start('dia_ar5')
  80. CALL wrk_alloc( jpi, jpj , zarea_ssh, zbotpres, zpe, z2d )
  81. CALL wrk_alloc( jpi, jpj, jpk , zrhd , zrhop, z3d, ztpot )
  82. CALL wrk_alloc( jpi, jpj, jpk, jpts, ztsn )
  83. !Call to init moved to here so that we can call iom_use in the
  84. !initialisation
  85. IF( kt == nit000 ) CALL dia_ar5_init
  86. zarea_ssh(:,:) = area(:,:) * sshn(:,:)
  87. ! ! total volume of liquid seawater
  88. zvolssh = glob_sum( zarea_ssh(:,:) )
  89. zvol = vol0 + zvolssh
  90. !
  91. z3d(:,:,jpk) = 0._wp ! ocean volume
  92. DO jk = 1, jpkm1
  93. z3d (:,:,jk) = area(:,:) * fse3t(:,:,jk) * tmask(:,:,jk)
  94. END DO
  95. CALL iom_put( "bathy" , bathy(:,:) )
  96. CALL iom_put( "areacello" , area(:,:) )
  97. !
  98. CALL iom_put( "volcello" , z3d(:,:,:) ) ! WARNING not consistent with CMIP DR where volcello is at ca. 2000
  99. CALL iom_put( "masscello" , rau0 * fse3t(:,:,:) * tmask(:,:,:) ) ! ocean mass
  100. !
  101. CALL iom_put( 'voltot', zvol )
  102. CALL iom_put( 'sshtot', zvolssh / area_tot )
  103. CALL iom_put( 'sshdyn', sshn(:,:) - (zvolssh / area_tot) )
  104. !
  105. IF( iom_use('sshthster')) THEN
  106. ztsn(:,:,:,jp_tem) = tsn(:,:,:,jp_tem) ! thermosteric ssh
  107. ztsn(:,:,:,jp_sal) = sn0(:,:,:)
  108. CALL eos( ztsn, zrhd, fsdept_n(:,:,:) ) ! now in situ density using initial salinity
  109. !
  110. zbotpres(:,:) = 0._wp ! no atmospheric surface pressure, levitating sea-ice
  111. DO jk = 1, jpkm1
  112. zbotpres(:,:) = zbotpres(:,:) + fse3t(:,:,jk) * zrhd(:,:,jk)
  113. END DO
  114. IF( .NOT.lk_vvl ) THEN
  115. IF ( ln_isfcav ) THEN
  116. DO jj = 1, jpj
  117. DO ji = 1, jpi
  118. iks = mikt(ji,jj)
  119. zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,iks) + riceload(ji,jj)
  120. END DO
  121. END DO
  122. ELSE
  123. zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1)
  124. END IF
  125. END IF
  126. !
  127. zarho = glob_sum( area(:,:) * zbotpres(:,:) )
  128. zssh_steric = - zarho / area_tot
  129. CALL iom_put( 'sshthster', zssh_steric )
  130. !
  131. ENDIF
  132. ! ! steric sea surface height
  133. CALL eos( tsn, zrhd, zrhop, fsdept_n(:,:,:) ) ! now in situ and potential density
  134. zbotpres(:,:) = 0._wp ! no atmospheric surface pressure, levitating sea-ice
  135. DO jk = 1, jpkm1
  136. zbotpres(:,:) = zbotpres(:,:) + fse3t(:,:,jk) * zrhd(:,:,jk)
  137. END DO
  138. IF( .NOT.lk_vvl ) THEN
  139. IF ( ln_isfcav ) THEN
  140. DO jj = 1, jpj
  141. DO ji = 1, jpi
  142. iks = mikt(ji,jj)
  143. zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,iks) + riceload(ji,jj)
  144. END DO
  145. END DO
  146. ELSE
  147. zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1)
  148. END IF
  149. END IF
  150. !
  151. zrhop(:,:,jpk) = 0._wp
  152. CALL iom_put( 'rhop', zrhop )
  153. !
  154. zarho = glob_sum( area(:,:) * zbotpres(:,:) )
  155. zssh_steric = - zarho / area_tot
  156. CALL iom_put( 'sshsteric', zssh_steric )
  157. ! ! ocean bottom pressure
  158. zztmp = rau0 * grav * 1.e-4_wp ! recover pressure from pressure anomaly and cover to dbar = 1.e4 Pa
  159. zbotpres(:,:) = zztmp * ( zbotpres(:,:) + sshn(:,:) + thick0(:,:) )
  160. CALL iom_put( 'botpres', zbotpres )
  161. ! ! Mean density anomalie, temperature and salinity
  162. !
  163. ztsn(:,:,:,:) = 0._wp ! ztsn(:,:,1,jp_tem/sal) is used here as 2D Workspace for temperature & salinity
  164. DO jk = 1, jpkm1
  165. ztsn(:,:,1,jp_tem) = ztsn(:,:,1,jp_tem) + area(:,:) * fse3t(:,:,jk) * tsn(:,:,jk,jp_tem)
  166. ztsn(:,:,1,jp_sal) = ztsn(:,:,1,jp_sal) + area(:,:) * fse3t(:,:,jk) * tsn(:,:,jk,jp_sal)
  167. ENDDO
  168. IF( .NOT.lk_vvl ) THEN
  169. IF( ln_isfcav ) THEN
  170. DO jj = 1, jpj
  171. DO ji = 1, jpi
  172. iks = mikt(ji,jj)
  173. ztsn(ji,jj,1,jp_tem) = ztsn(ji,jj,1,jp_tem) + zarea_ssh(ji,jj) * tsn(ji,jj,iks,jp_tem)
  174. ztsn(ji,jj,1,jp_sal) = ztsn(ji,jj,1,jp_sal) + zarea_ssh(ji,jj) * tsn(ji,jj,iks,jp_sal)
  175. END DO
  176. END DO
  177. ELSE
  178. ztsn(:,:,1,jp_tem) = ztsn(:,:,1,jp_tem) + zarea_ssh(:,:) * tsn(:,:,1,jp_tem)
  179. ztsn(:,:,1,jp_sal) = ztsn(:,:,1,jp_sal) + zarea_ssh(:,:) * tsn(:,:,1,jp_sal)
  180. END IF
  181. ENDIF
  182. !
  183. zsss = glob_sum( area(:,:) * tsn(:,:,1,jp_sal) )
  184. ztemp = glob_sum( ztsn(:,:,1,jp_tem) )
  185. zsal = glob_sum( ztsn(:,:,1,jp_sal) )
  186. zmass = rau0 * ( zarho + zvol )
  187. !
  188. CALL iom_put( 'masstot', zmass ) ! total mass of liquid seawater
  189. CALL iom_put( 'temptot', ztemp / zvol ) ! potential temperature in liquid seawater
  190. CALL iom_put( 'saltot' , zsal / zvol ) ! Salinity of liquid seawater
  191. CALL iom_put( 'ssstot' , zsss / area_tot ) ! Salinity of liquid seawater at surface
  192. IF( iom_use( "e3tb" ) ) THEN
  193. DO jj = 1, jpj
  194. DO ji = 1, jpi
  195. ikb = mbkt(ji,jj)
  196. z2d(ji,jj) = fse3t(ji,jj,ikb)
  197. END DO
  198. END DO
  199. CALL iom_put( "e3tb", z2d )
  200. ENDIF
  201. IF( nn_eos == -1 ) THEN ! ! potential temperature (TEOS-10 case)
  202. ztpot(:,:,:) = eos_pt_from_ct( tsn(:,:,:,jp_tem), tsn(:,:,:,jp_sal) )
  203. ztpot(:,:,jpk) = 0._wp
  204. !
  205. CALL iom_put( "toce_pot", ztpot(:,:,:) ) ! potential temperature (TEOS-10 case)
  206. CALL iom_put( "sst_pot" , ztpot(:,:,1) ) ! surface temperature
  207. !
  208. IF( iom_use('temptot_pot') ) THEN ! Output potential temperature in case we use TEOS-10
  209. z2d(:,:) = 0._wp
  210. DO jk = 1, jpkm1
  211. z2d(:,:) = z2d(:,:) + area(:,:) * fse3t(:,:,jk) * ztpot(:,:,jk)
  212. END DO
  213. zsst = glob_sum( area(:,:) * ztpot(:,:,1) )
  214. ztemp = glob_sum( z2d(:,:) )
  215. CALL iom_put( 'temptot_pot', ztemp / zvol )
  216. CALL iom_put( 'ssttot' , zsst / area_tot )
  217. ENDIF
  218. ! Vertical integral of temperature
  219. IF( iom_use("tosmint_pot") ) THEN
  220. z2d(:,:) = 0._wp
  221. DO jk = 1, jpkm1
  222. DO jj = 1, jpj
  223. DO ji = 1, jpi ! vector opt.
  224. z2d(ji,jj) = z2d(ji,jj) + rau0 * fse3t(ji,jj,jk) * ztpot(ji,jj,jk)
  225. END DO
  226. END DO
  227. END DO
  228. CALL iom_put( "tosmint_pot", z2d )
  229. ENDIF
  230. ELSE
  231. zsst = glob_sum( area(:,:) * tsn(:,:,1,jp_tem) ) ! Case EOS-80 : compute sst anyway
  232. CALL iom_put('ssttot', zsst / area_tot )
  233. ENDIF
  234. IF( iom_use( 'tnpeo' )) THEN
  235. ! Work done against stratification by vertical mixing
  236. ! Exclude points where rn2 is negative as convection kicks in here and
  237. ! work is not being done against stratification
  238. zpe(:,:) = 0._wp
  239. IF( lk_zdfddm ) THEN
  240. DO jk = 2, jpk
  241. DO jj = 1, jpj
  242. DO ji = 1, jpi
  243. IF( rn2(ji,jj,jk) > 0._wp ) THEN
  244. zrw = ( fsdepw(ji,jj,jk ) - fsdept(ji,jj,jk) ) &
  245. & / ( fsdept(ji,jj,jk-1) - fsdept(ji,jj,jk) )
  246. !
  247. zaw = rab_n(ji,jj,jk,jp_tem) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_tem)* zrw
  248. zbw = rab_n(ji,jj,jk,jp_sal) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_sal)* zrw
  249. !
  250. zpe(ji, jj) = zpe(ji, jj) &
  251. & - grav * ( avt(ji,jj,jk) * zaw * (tsn(ji,jj,jk-1,jp_tem) - tsn(ji,jj,jk,jp_tem) ) &
  252. & - fsavs(ji,jj,jk) * zbw * (tsn(ji,jj,jk-1,jp_sal) - tsn(ji,jj,jk,jp_sal) ) )
  253. ENDIF
  254. END DO
  255. END DO
  256. END DO
  257. ELSE
  258. DO jk = 1, jpk
  259. DO ji = 1, jpi
  260. DO jj = 1, jpj
  261. zpe(ji,jj) = zpe(ji,jj) + avt(ji, jj, jk) * MIN(0._wp,rn2(ji, jj, jk)) * rau0 * fse3w(ji, jj, jk)
  262. END DO
  263. END DO
  264. END DO
  265. ENDIF
  266. CALL lbc_lnk(zpe, 'T', 1._wp)
  267. CALL iom_put( 'tnpeo', zpe )
  268. ENDIF
  269. !
  270. CALL wrk_dealloc( jpi, jpj , zarea_ssh, zbotpres, zpe, z2d )
  271. CALL wrk_dealloc( jpi, jpj, jpk , zrhd , zrhop, z3d, ztpot )
  272. CALL wrk_dealloc( jpi, jpj, jpk, jpts, ztsn )
  273. !
  274. IF( nn_timing == 1 ) CALL timing_stop('dia_ar5')
  275. !
  276. END SUBROUTINE dia_ar5
  277. SUBROUTINE dia_ar5_init
  278. !!----------------------------------------------------------------------
  279. !! *** ROUTINE dia_ar5_init ***
  280. !!
  281. !! ** Purpose : initialization for AR5 diagnostic computation
  282. !!----------------------------------------------------------------------
  283. INTEGER :: inum
  284. INTEGER :: ik, idep
  285. INTEGER :: ji, jj, jk ! dummy loop indices
  286. REAL(wp) :: zztmp
  287. REAL(wp), POINTER, DIMENSION(:,:,:,:) :: zsaldta ! Jan/Dec levitus salinity
  288. REAL(wp), POINTER, DIMENSION(:,:) :: zvol0
  289. ! reading initial file
  290. LOGICAL :: ln_tsd_init !: T & S data flag
  291. LOGICAL :: ln_tsd_tradmp !: internal damping toward input data flag
  292. CHARACTER(len=100) :: cn_dir
  293. TYPE(FLD_N) :: sn_tem,sn_sal
  294. INTEGER :: ios=0
  295. NAMELIST/namtsd/ ln_tsd_init,ln_tsd_tradmp,cn_dir,sn_tem,sn_sal
  296. !
  297. REWIND( numnam_ref ) ! Namelist namtsd in reference namelist :
  298. READ ( numnam_ref, namtsd, IOSTAT = ios, ERR = 901)
  299. 901 IF( ios /= 0 ) CALL ctl_nam ( ios , ' namtsd in reference namelist for dia_ar5', lwp )
  300. REWIND( numnam_cfg ) ! Namelist namtsd in configuration namelist : Parameters of the run
  301. READ ( numnam_cfg, namtsd, IOSTAT = ios, ERR = 902 )
  302. 902 IF( ios /= 0 ) CALL ctl_nam ( ios , ' namtsd in configuration namelist for dia_ar5', lwp )
  303. IF(lwm) WRITE ( numond, namtsd )
  304. !
  305. !!----------------------------------------------------------------------
  306. !
  307. IF( nn_timing == 1 ) CALL timing_start('dia_ar5_init')
  308. !
  309. CALL wrk_alloc( jpi, jpj, jpk, 2, zsaldta )
  310. CALL wrk_alloc( jpi, jpj, zvol0 )
  311. ! ! allocate dia_ar5 arrays
  312. IF( dia_ar5_alloc() /= 0 ) CALL ctl_stop( 'STOP', 'dia_ar5_init : unable to allocate arrays' )
  313. area(:,:) = e1t(:,:) * e2t(:,:)
  314. area_tot = glob_sum( area(:,:) )
  315. zvol0 (:,:) = 0._wp
  316. thick0(:,:) = 0._wp
  317. DO jk = 1, jpkm1
  318. DO jj = 1, jpj ! interpolation of salinity at the last ocean level (i.e. the partial step)
  319. DO ji = 1, jpi
  320. idep = tmask(ji,jj,jk) * e3t_0(ji,jj,jk)
  321. zvol0 (ji,jj) = zvol0 (ji,jj) + idep * area(ji,jj)
  322. thick0(ji,jj) = thick0(ji,jj) + idep
  323. END DO
  324. END DO
  325. END DO
  326. vol0 = glob_sum( zvol0 )
  327. IF( iom_use('sshthster')) THEN
  328. CALL iom_open ( 'sali_ref_clim_monthly', inum )
  329. CALL iom_get ( inum, jpdom_data, 'vosaline' , zsaldta(:,:,:,1), 1 )
  330. CALL iom_get ( inum, jpdom_data, 'vosaline' , zsaldta(:,:,:,2), 12 )
  331. CALL iom_close( inum )
  332. sn0(:,:,:) = 0.5_wp * ( zsaldta(:,:,:,1) + zsaldta(:,:,:,2) )
  333. sn0(:,:,:) = sn0(:,:,:) * tmask(:,:,:)
  334. IF( ln_zps ) THEN ! z-coord. partial steps
  335. DO jj = 1, jpj ! interpolation of salinity at the last ocean level (i.e. the partial step)
  336. DO ji = 1, jpi
  337. ik = mbkt(ji,jj)
  338. IF( ik > 1 ) THEN
  339. zztmp = ( gdept_1d(ik) - gdept_0(ji,jj,ik) ) / ( gdept_1d(ik) - gdept_1d(ik-1) )
  340. sn0(ji,jj,ik) = ( 1._wp - zztmp ) * sn0(ji,jj,ik) + zztmp * sn0(ji,jj,ik-1)
  341. ENDIF
  342. END DO
  343. END DO
  344. ENDIF
  345. ENDIF
  346. !
  347. CALL wrk_dealloc( jpi, jpj, jpk, jpts, zsaldta )
  348. CALL wrk_dealloc( jpi, jpj, zvol0 )
  349. !
  350. IF( nn_timing == 1 ) CALL timing_stop('dia_ar5_init')
  351. !
  352. END SUBROUTINE dia_ar5_init
  353. #else
  354. !!----------------------------------------------------------------------
  355. !! Default option : NO diaar5
  356. !!----------------------------------------------------------------------
  357. LOGICAL, PUBLIC, PARAMETER :: lk_diaar5 = .FALSE. ! coupled flag
  358. CONTAINS
  359. SUBROUTINE dia_ar5_init ! Dummy routine
  360. END SUBROUTINE dia_ar5_init
  361. SUBROUTINE dia_ar5( kt ) ! Empty routine
  362. INTEGER :: kt
  363. WRITE(*,*) 'dia_ar5: You should not have seen this print! error?', kt
  364. END SUBROUTINE dia_ar5
  365. #endif
  366. !!======================================================================
  367. END MODULE diaar5