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