[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 |
---|
[1756] | 23 | |
---|
| 24 | IMPLICIT NONE |
---|
| 25 | PRIVATE |
---|
| 26 | |
---|
[2528] | 27 | PUBLIC dia_ar5 ! routine called in step.F90 module |
---|
| 28 | PUBLIC dia_ar5_init ! routine called in opa.F90 module |
---|
[2715] | 29 | PUBLIC dia_ar5_alloc ! routine called in nemogcm.F90 module |
---|
[1756] | 30 | |
---|
| 31 | LOGICAL, PUBLIC, PARAMETER :: lk_diaar5 = .TRUE. ! coupled flag |
---|
| 32 | |
---|
[2528] | 33 | REAL(wp) :: vol0 ! ocean volume (interior domain) |
---|
| 34 | REAL(wp) :: area_tot ! total ocean surface (interior domain) |
---|
[2715] | 35 | REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,: ) :: area ! cell surface (interior domain) |
---|
| 36 | REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,: ) :: thick0 ! ocean thickness (interior domain) |
---|
| 37 | REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: sn0 ! initial salinity |
---|
[1756] | 38 | |
---|
| 39 | !! * Substitutions |
---|
| 40 | # include "domzgr_substitute.h90" |
---|
| 41 | !!---------------------------------------------------------------------- |
---|
[2528] | 42 | !! NEMO/OPA 3.3 , NEMO Consortium (2010) |
---|
| 43 | !! $Id$ |
---|
| 44 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
[1756] | 45 | !!---------------------------------------------------------------------- |
---|
| 46 | CONTAINS |
---|
| 47 | |
---|
[2715] | 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 | |
---|
| 62 | |
---|
[1756] | 63 | SUBROUTINE dia_ar5( kt ) |
---|
| 64 | !!---------------------------------------------------------------------- |
---|
| 65 | !! *** ROUTINE dia_ar5 *** |
---|
| 66 | !! |
---|
[2528] | 67 | !! ** Purpose : compute and output some AR5 diagnostics |
---|
[1756] | 68 | !!---------------------------------------------------------------------- |
---|
[2715] | 69 | ! |
---|
[1756] | 70 | INTEGER, INTENT( in ) :: kt ! ocean time-step index |
---|
[2715] | 71 | ! |
---|
[1756] | 72 | INTEGER :: ji, jj, jk ! dummy loop arguments |
---|
| 73 | REAL(wp) :: zvolssh, zvol, zssh_steric, zztmp, zarho, ztemp, zsal, zmass |
---|
[3294] | 74 | ! |
---|
| 75 | REAL(wp), POINTER, DIMENSION(:,:) :: zarea_ssh , zbotpres ! 2D workspace |
---|
| 76 | REAL(wp), POINTER, DIMENSION(:,:,:) :: zrhd , zrhop ! 3D workspace |
---|
| 77 | REAL(wp), POINTER, DIMENSION(:,:,:,:) :: ztsn ! 4D workspace |
---|
[1756] | 78 | !!-------------------------------------------------------------------- |
---|
[3294] | 79 | IF( nn_timing == 1 ) CALL timing_start('dia_ar5') |
---|
| 80 | |
---|
| 81 | CALL wrk_alloc( jpi , jpj , zarea_ssh , zbotpres ) |
---|
| 82 | CALL wrk_alloc( jpi , jpj , jpk , zrhd , zrhop ) |
---|
| 83 | CALL wrk_alloc( jpi , jpj , jpk , jpts , ztsn ) |
---|
[1756] | 84 | |
---|
| 85 | CALL iom_put( 'cellthc', fse3t(:,:,:) ) |
---|
| 86 | |
---|
| 87 | zarea_ssh(:,:) = area(:,:) * sshn(:,:) |
---|
| 88 | |
---|
| 89 | ! ! total volume of liquid seawater |
---|
| 90 | zvolssh = SUM( zarea_ssh(:,:) ) |
---|
| 91 | IF( lk_mpp ) CALL mpp_sum( zvolssh ) |
---|
| 92 | zvol = vol0 + zvolssh |
---|
| 93 | |
---|
| 94 | CALL iom_put( 'voltot', zvol ) |
---|
| 95 | CALL iom_put( 'sshtot', zvolssh / area_tot ) |
---|
| 96 | |
---|
[3294] | 97 | ! |
---|
| 98 | ztsn(:,:,:,jp_tem) = tsn(:,:,:,jp_tem) ! thermosteric ssh |
---|
[2528] | 99 | ztsn(:,:,:,jp_sal) = sn0(:,:,:) |
---|
[4313] | 100 | CALL eos( ztsn, zrhd, fsdept_n(:,:,:) ) ! now in situ density using initial salinity |
---|
[1756] | 101 | ! |
---|
[2528] | 102 | zbotpres(:,:) = 0._wp ! no atmospheric surface pressure, levitating sea-ice |
---|
[1756] | 103 | DO jk = 1, jpkm1 |
---|
| 104 | zbotpres(:,:) = zbotpres(:,:) + fse3t(:,:,jk) * zrhd(:,:,jk) |
---|
| 105 | END DO |
---|
[4990] | 106 | IF( .NOT.lk_vvl ) THEN |
---|
[5111] | 107 | IF ( ln_isfcav ) THEN |
---|
| 108 | DO ji=1,jpi |
---|
| 109 | DO jj=1,jpj |
---|
| 110 | zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,mikt(ji,jj)) + riceload(ji,jj) |
---|
| 111 | END DO |
---|
[4990] | 112 | END DO |
---|
[5111] | 113 | ELSE |
---|
| 114 | zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1) |
---|
| 115 | END IF |
---|
[4990] | 116 | END IF |
---|
[1756] | 117 | ! |
---|
| 118 | zarho = SUM( area(:,:) * zbotpres(:,:) ) |
---|
| 119 | IF( lk_mpp ) CALL mpp_sum( zarho ) |
---|
| 120 | zssh_steric = - zarho / area_tot |
---|
| 121 | CALL iom_put( 'sshthster', zssh_steric ) |
---|
| 122 | |
---|
| 123 | ! ! steric sea surface height |
---|
[4313] | 124 | CALL eos( tsn, zrhd, zrhop, fsdept_n(:,:,:) ) ! now in situ and potential density |
---|
[2528] | 125 | zrhop(:,:,jpk) = 0._wp |
---|
[1756] | 126 | CALL iom_put( 'rhop', zrhop ) |
---|
| 127 | ! |
---|
[2528] | 128 | zbotpres(:,:) = 0._wp ! no atmospheric surface pressure, levitating sea-ice |
---|
[1756] | 129 | DO jk = 1, jpkm1 |
---|
| 130 | zbotpres(:,:) = zbotpres(:,:) + fse3t(:,:,jk) * zrhd(:,:,jk) |
---|
| 131 | END DO |
---|
[4990] | 132 | IF( .NOT.lk_vvl ) THEN |
---|
[5111] | 133 | IF ( ln_isfcav ) THEN |
---|
| 134 | DO ji=1,jpi |
---|
| 135 | DO jj=1,jpj |
---|
| 136 | zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,mikt(ji,jj)) + riceload(ji,jj) |
---|
| 137 | END DO |
---|
[4990] | 138 | END DO |
---|
[5111] | 139 | ELSE |
---|
| 140 | zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1) |
---|
| 141 | END IF |
---|
[4990] | 142 | END IF |
---|
[1756] | 143 | ! |
---|
| 144 | zarho = SUM( area(:,:) * zbotpres(:,:) ) |
---|
| 145 | IF( lk_mpp ) CALL mpp_sum( zarho ) |
---|
| 146 | zssh_steric = - zarho / area_tot |
---|
| 147 | CALL iom_put( 'sshsteric', zssh_steric ) |
---|
| 148 | |
---|
| 149 | ! ! ocean bottom pressure |
---|
[2528] | 150 | zztmp = rau0 * grav * 1.e-4_wp ! recover pressure from pressure anomaly and cover to dbar = 1.e4 Pa |
---|
[1756] | 151 | zbotpres(:,:) = zztmp * ( zbotpres(:,:) + sshn(:,:) + thick0(:,:) ) |
---|
| 152 | CALL iom_put( 'botpres', zbotpres ) |
---|
| 153 | |
---|
| 154 | ! ! Mean density anomalie, temperature and salinity |
---|
[2528] | 155 | ztemp = 0._wp |
---|
| 156 | zsal = 0._wp |
---|
[1756] | 157 | DO jk = 1, jpkm1 |
---|
| 158 | DO jj = 1, jpj |
---|
| 159 | DO ji = 1, jpi |
---|
| 160 | zztmp = area(ji,jj) * fse3t(ji,jj,jk) |
---|
[3294] | 161 | ztemp = ztemp + zztmp * tsn(ji,jj,jk,jp_tem) |
---|
| 162 | zsal = zsal + zztmp * tsn(ji,jj,jk,jp_sal) |
---|
[1756] | 163 | END DO |
---|
| 164 | END DO |
---|
| 165 | END DO |
---|
[2528] | 166 | IF( .NOT.lk_vvl ) THEN |
---|
[5111] | 167 | IF ( ln_isfcav ) THEN |
---|
| 168 | DO ji=1,jpi |
---|
| 169 | DO jj=1,jpj |
---|
| 170 | ztemp = ztemp + zarea_ssh(ji,jj) * tsn(ji,jj,mikt(ji,jj),jp_tem) |
---|
| 171 | zsal = zsal + zarea_ssh(ji,jj) * tsn(ji,jj,mikt(ji,jj),jp_sal) |
---|
| 172 | END DO |
---|
[4990] | 173 | END DO |
---|
[5111] | 174 | ELSE |
---|
| 175 | ztemp = ztemp + zarea_ssh(:,:) * tsn(:,:,1,jp_tem) |
---|
| 176 | zsal = zsal + zarea_ssh(:,:) * tsn(:,:,1,jp_sal) |
---|
| 177 | END IF |
---|
[1756] | 178 | ENDIF |
---|
| 179 | IF( lk_mpp ) THEN |
---|
| 180 | CALL mpp_sum( ztemp ) |
---|
| 181 | CALL mpp_sum( zsal ) |
---|
| 182 | END IF |
---|
| 183 | ! |
---|
| 184 | zmass = rau0 * ( zarho + zvol ) ! total mass of liquid seawater |
---|
| 185 | ztemp = ztemp / zvol ! potential temperature in liquid seawater |
---|
| 186 | zsal = zsal / zvol ! Salinity of liquid seawater |
---|
| 187 | ! |
---|
| 188 | CALL iom_put( 'masstot', zmass ) |
---|
| 189 | CALL iom_put( 'temptot', ztemp ) |
---|
| 190 | CALL iom_put( 'saltot' , zsal ) |
---|
| 191 | ! |
---|
[3294] | 192 | CALL wrk_dealloc( jpi , jpj , zarea_ssh , zbotpres ) |
---|
| 193 | CALL wrk_dealloc( jpi , jpj , jpk , zrhd , zrhop ) |
---|
| 194 | CALL wrk_dealloc( jpi , jpj , jpk , jpts , ztsn ) |
---|
[2715] | 195 | ! |
---|
[3294] | 196 | IF( nn_timing == 1 ) CALL timing_stop('dia_ar5') |
---|
| 197 | ! |
---|
[1756] | 198 | END SUBROUTINE dia_ar5 |
---|
| 199 | |
---|
| 200 | |
---|
| 201 | SUBROUTINE dia_ar5_init |
---|
| 202 | !!---------------------------------------------------------------------- |
---|
| 203 | !! *** ROUTINE dia_ar5_init *** |
---|
| 204 | !! |
---|
[2528] | 205 | !! ** Purpose : initialization for AR5 diagnostic computation |
---|
[1756] | 206 | !!---------------------------------------------------------------------- |
---|
| 207 | INTEGER :: inum |
---|
| 208 | INTEGER :: ik |
---|
| 209 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
| 210 | REAL(wp) :: zztmp |
---|
[2715] | 211 | REAL(wp), POINTER, DIMENSION(:,:,:,:) :: zsaldta ! Jan/Dec levitus salinity |
---|
[1756] | 212 | !!---------------------------------------------------------------------- |
---|
| 213 | ! |
---|
[3294] | 214 | IF( nn_timing == 1 ) CALL timing_start('dia_ar5_init') |
---|
| 215 | ! |
---|
| 216 | CALL wrk_alloc( jpi , jpj , jpk, jpts, zsaldta ) |
---|
[2715] | 217 | ! ! allocate dia_ar5 arrays |
---|
| 218 | IF( dia_ar5_alloc() /= 0 ) CALL ctl_stop( 'STOP', 'dia_ar5_init : unable to allocate arrays' ) |
---|
| 219 | |
---|
[1756] | 220 | area(:,:) = e1t(:,:) * e2t(:,:) * tmask_i(:,:) |
---|
| 221 | |
---|
| 222 | area_tot = SUM( area(:,:) ) ; IF( lk_mpp ) CALL mpp_sum( area_tot ) |
---|
| 223 | |
---|
[2528] | 224 | vol0 = 0._wp |
---|
| 225 | thick0(:,:) = 0._wp |
---|
[1756] | 226 | DO jk = 1, jpkm1 |
---|
[4292] | 227 | vol0 = vol0 + SUM( area (:,:) * tmask(:,:,jk) * e3t_0(:,:,jk) ) |
---|
| 228 | thick0(:,:) = thick0(:,:) + tmask_i(:,:) * tmask(:,:,jk) * e3t_0(:,:,jk) |
---|
[1756] | 229 | END DO |
---|
[1948] | 230 | IF( lk_mpp ) CALL mpp_sum( vol0 ) |
---|
[1756] | 231 | |
---|
| 232 | CALL iom_open ( 'data_1m_salinity_nomask', inum ) |
---|
| 233 | CALL iom_get ( inum, jpdom_data, 'vosaline', zsaldta(:,:,:,1), 1 ) |
---|
| 234 | CALL iom_get ( inum, jpdom_data, 'vosaline', zsaldta(:,:,:,2), 12 ) |
---|
| 235 | CALL iom_close( inum ) |
---|
[2528] | 236 | sn0(:,:,:) = 0.5_wp * ( zsaldta(:,:,:,1) + zsaldta(:,:,:,2) ) |
---|
| 237 | sn0(:,:,:) = sn0(:,:,:) * tmask(:,:,:) |
---|
[1756] | 238 | IF( ln_zps ) THEN ! z-coord. partial steps |
---|
| 239 | DO jj = 1, jpj ! interpolation of salinity at the last ocean level (i.e. the partial step) |
---|
| 240 | DO ji = 1, jpi |
---|
[2528] | 241 | ik = mbkt(ji,jj) |
---|
| 242 | IF( ik > 1 ) THEN |
---|
[4292] | 243 | zztmp = ( gdept_1d(ik) - gdept_0(ji,jj,ik) ) / ( gdept_1d(ik) - gdept_1d(ik-1) ) |
---|
[2528] | 244 | sn0(ji,jj,ik) = ( 1._wp - zztmp ) * sn0(ji,jj,ik) + zztmp * sn0(ji,jj,ik-1) |
---|
[1756] | 245 | ENDIF |
---|
| 246 | END DO |
---|
| 247 | END DO |
---|
| 248 | ENDIF |
---|
| 249 | ! |
---|
[3294] | 250 | CALL wrk_dealloc( jpi , jpj , jpk, jpts, zsaldta ) |
---|
[2715] | 251 | ! |
---|
[3294] | 252 | IF( nn_timing == 1 ) CALL timing_stop('dia_ar5_init') |
---|
| 253 | ! |
---|
[1756] | 254 | END SUBROUTINE dia_ar5_init |
---|
| 255 | |
---|
| 256 | #else |
---|
| 257 | !!---------------------------------------------------------------------- |
---|
| 258 | !! Default option : NO diaar5 |
---|
| 259 | !!---------------------------------------------------------------------- |
---|
| 260 | LOGICAL, PUBLIC, PARAMETER :: lk_diaar5 = .FALSE. ! coupled flag |
---|
| 261 | CONTAINS |
---|
[2528] | 262 | SUBROUTINE dia_ar5_init ! Dummy routine |
---|
| 263 | END SUBROUTINE dia_ar5_init |
---|
[1756] | 264 | SUBROUTINE dia_ar5( kt ) ! Empty routine |
---|
[2528] | 265 | INTEGER :: kt |
---|
[1756] | 266 | WRITE(*,*) 'dia_ar5: You should not have seen this print! error?', kt |
---|
| 267 | END SUBROUTINE dia_ar5 |
---|
| 268 | #endif |
---|
| 269 | |
---|
| 270 | !!====================================================================== |
---|
| 271 | END MODULE diaar5 |
---|