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