- Timestamp:
- 2018-07-23T11:33:03+02:00 (6 years ago)
- Location:
- branches/UKMO/dev_r5518_obs_oper_update_icethick/NEMOGCM/NEMO/OPA_SRC/DIA
- Files:
-
- 7 edited
- 2 copied
Legend:
- Unmodified
- Added
- Removed
-
branches/UKMO/dev_r5518_obs_oper_update_icethick/NEMOGCM/NEMO/OPA_SRC/DIA/diaar5.F90
r7960 r9987 24 24 USE phycst ! physical constant 25 25 USE in_out_manager ! I/O manager 26 USE zdfddm 27 USE zdf_oce 26 28 27 29 IMPLICIT NONE … … 42 44 !! * Substitutions 43 45 # include "domzgr_substitute.h90" 46 # include "zdfddm_substitute.h90" 44 47 !!---------------------------------------------------------------------- 45 48 !! NEMO/OPA 3.3 , NEMO Consortium (2010) … … 75 78 INTEGER :: ji, jj, jk ! dummy loop arguments 76 79 REAL(wp) :: zvolssh, zvol, zssh_steric, zztmp, zarho, ztemp, zsal, zmass 80 REAL(wp) :: zaw, zbw, zrw 77 81 ! 78 82 REAL(wp), POINTER, DIMENSION(:,:) :: zarea_ssh , zbotpres ! 2D workspace 83 REAL(wp), POINTER, DIMENSION(:,:) :: pe ! 2D workspace 79 84 REAL(wp), POINTER, DIMENSION(:,:,:) :: zrhd , zrhop ! 3D workspace 80 85 REAL(wp), POINTER, DIMENSION(:,:,:,:) :: ztsn ! 4D workspace 81 86 !!-------------------------------------------------------------------- 82 87 IF( nn_timing == 1 ) CALL timing_start('dia_ar5') 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 83 92 84 CALL wrk_alloc( jpi , jpj , zarea_ssh , zbotpres )93 CALL wrk_alloc( jpi , jpj , zarea_ssh , zbotpres, pe ) 85 94 CALL wrk_alloc( jpi , jpj , jpk , zrhd , zrhop ) 86 95 CALL wrk_alloc( jpi , jpj , jpk , jpts , ztsn ) … … 95 104 CALL iom_put( 'voltot', zvol ) 96 105 CALL iom_put( 'sshtot', zvolssh / area_tot ) 106 CALL iom_put( 'sshdyn', sshn(:,:) - (zvolssh / area_tot) ) 97 107 98 108 ! 99 ztsn(:,:,:,jp_tem) = tsn(:,:,:,jp_tem) ! thermosteric ssh 100 ztsn(:,:,:,jp_sal) = sn0(:,:,:) 101 CALL eos( ztsn, zrhd, fsdept_n(:,:,:) ) ! now in situ density using initial salinity 102 ! 103 zbotpres(:,:) = 0._wp ! no atmospheric surface pressure, levitating sea-ice 104 DO jk = 1, jpkm1 105 zbotpres(:,:) = zbotpres(:,:) + fse3t(:,:,jk) * zrhd(:,:,jk) 106 END DO 107 IF( .NOT.lk_vvl ) THEN 108 IF ( ln_isfcav ) THEN 109 DO ji=1,jpi 110 DO jj=1,jpj 111 zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,mikt(ji,jj)) + riceload(ji,jj) 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 112 124 END DO 113 E ND DO114 ELSE115 zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1)125 ELSE 126 zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1) 127 END IF 116 128 END IF 117 END IF118 129 ! 119 zarho = SUM( area(:,:) * zbotpres(:,:) ) 120 IF( lk_mpp ) CALL mpp_sum( zarho ) 121 zssh_steric = - zarho / area_tot 122 CALL iom_put( 'sshthster', zssh_steric ) 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 123 135 124 136 ! ! steric sea surface height … … 190 202 CALL iom_put( 'temptot', ztemp ) 191 203 CALL iom_put( 'saltot' , zsal ) 192 ! 193 CALL wrk_dealloc( jpi , jpj , zarea_ssh , zbotpres ) 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 238 ! 239 CALL wrk_dealloc( jpi , jpj , zarea_ssh , zbotpres, pe ) 194 240 CALL wrk_dealloc( jpi , jpj , jpk , zrhd , zrhop ) 195 241 CALL wrk_dealloc( jpi , jpj , jpk , jpts , ztsn ) … … 211 257 REAL(wp) :: zztmp 212 258 REAL(wp), POINTER, DIMENSION(:,:,:,:) :: zsaldta ! Jan/Dec levitus salinity 213 ! reading initial file214 LOGICAL :: ln_tsd_init !: T & S data flag215 LOGICAL :: ln_tsd_tradmp !: internal damping toward input data flag216 CHARACTER(len=100) :: cn_dir217 TYPE(FLD_N) :: sn_tem,sn_sal218 INTEGER :: ios=0219 220 NAMELIST/namtsd/ ln_tsd_init,ln_tsd_tradmp,cn_dir,sn_tem,sn_sal221 !222 223 REWIND( numnam_ref ) ! Namelist namtsd in reference namelist :224 READ ( numnam_ref, namtsd, IOSTAT = ios, ERR = 901)225 901 IF( ios /= 0 ) CALL ctl_nam ( ios , ' namtsd in reference namelist for dia_ar5', lwp )226 REWIND( numnam_cfg ) ! Namelist namtsd in configuration namelist : Parameters of the run227 READ ( numnam_cfg, namtsd, IOSTAT = ios, ERR = 902 )228 902 IF( ios /= 0 ) CALL ctl_nam ( ios , ' namtsd in configuration namelist for dia_ar5', lwp )229 IF(lwm) WRITE ( numond, namtsd )230 259 ! 231 260 !!---------------------------------------------------------------------- … … 233 262 IF( nn_timing == 1 ) CALL timing_start('dia_ar5_init') 234 263 ! 235 CALL wrk_alloc( jpi , jpj , jpk, jpts, zsaldta )264 CALL wrk_alloc( jpi, jpj, jpk, 2, zsaldta ) 236 265 ! ! allocate dia_ar5 arrays 237 266 IF( dia_ar5_alloc() /= 0 ) CALL ctl_stop( 'STOP', 'dia_ar5_init : unable to allocate arrays' ) … … 249 278 IF( lk_mpp ) CALL mpp_sum( vol0 ) 250 279 251 CALL iom_open ( TRIM( cn_dir )//TRIM(sn_sal%clname), inum ) 252 CALL iom_get ( inum, jpdom_data, TRIM(sn_sal%clvar), zsaldta(:,:,:,1), 1 ) 253 CALL iom_get ( inum, jpdom_data, TRIM(sn_sal%clvar), zsaldta(:,:,:,2), 12 ) 254 CALL iom_close( inum ) 255 sn0(:,:,:) = 0.5_wp * ( zsaldta(:,:,:,1) + zsaldta(:,:,:,2) ) 256 sn0(:,:,:) = sn0(:,:,:) * tmask(:,:,:) 257 IF( ln_zps ) THEN ! z-coord. partial steps 258 DO jj = 1, jpj ! interpolation of salinity at the last ocean level (i.e. the partial step) 259 DO ji = 1, jpi 260 ik = mbkt(ji,jj) 261 IF( ik > 1 ) THEN 262 zztmp = ( gdept_1d(ik) - gdept_0(ji,jj,ik) ) / ( gdept_1d(ik) - gdept_1d(ik-1) ) 263 sn0(ji,jj,ik) = ( 1._wp - zztmp ) * sn0(ji,jj,ik) + zztmp * sn0(ji,jj,ik-1) 264 ENDIF 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 ) 285 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 265 297 END DO 266 END DO298 ENDIF 267 299 ENDIF 268 300 ! 269 CALL wrk_dealloc( jpi , jpj , jpk, jpts, zsaldta )301 CALL wrk_dealloc( jpi, jpj, jpk, 2, zsaldta ) 270 302 ! 271 303 IF( nn_timing == 1 ) CALL timing_stop('dia_ar5_init') -
branches/UKMO/dev_r5518_obs_oper_update_icethick/NEMOGCM/NEMO/OPA_SRC/DIA/diadimg.F90
r7960 r9987 124 124 125 125 CASE DEFAULT 126 IF(lwp) WRITE(numout,*) ' E R R O R : bad cd_type in dia_wri_dimg ' 127 STOP 'dia_wri_dimg' 126 127 WRITE(numout,*) 'dia_wri_dimg : E R R O R : bad cd_type in dia_wri_dimg' 128 CALL ctl_stop( 'STOP', 'dia_wri_dimg :bad cd_type in dia_wri_dimg ' ) 128 129 129 130 END SELECT -
branches/UKMO/dev_r5518_obs_oper_update_icethick/NEMOGCM/NEMO/OPA_SRC/DIA/diaharm.F90
r7960 r9987 196 196 DO ji = 1,jpi 197 197 ! Elevation 198 ana_temp(ji,jj,nhc,1) = ana_temp(ji,jj,nhc,1) + ztemp*sshn(ji,jj) *tmask_i(ji,jj) 199 #if defined key_dynspg_ts 200 ana_temp(ji,jj,nhc,2) = ana_temp(ji,jj,nhc,2) + ztemp*un_b(ji,jj)*hur(ji,jj)*umask_i(ji,jj) 201 ana_temp(ji,jj,nhc,3) = ana_temp(ji,jj,nhc,3) + ztemp*vn_b(ji,jj)*hvr(ji,jj)*vmask_i(ji,jj) 202 #endif 198 ana_temp(ji,jj,nhc,1) = ana_temp(ji,jj,nhc,1) + ztemp*sshn(ji,jj)*tmask_i(ji,jj) 199 ana_temp(ji,jj,nhc,2) = ana_temp(ji,jj,nhc,2) + ztemp*un_b(ji,jj)*umask_i(ji,jj) 200 ana_temp(ji,jj,nhc,3) = ana_temp(ji,jj,nhc,3) + ztemp*vn_b(ji,jj)*vmask_i(ji,jj) 203 201 END DO 204 202 END DO -
branches/UKMO/dev_r5518_obs_oper_update_icethick/NEMOGCM/NEMO/OPA_SRC/DIA/diahsb.F90
r7960 r9987 93 93 ! 1 - Trends due to forcing ! 94 94 ! ------------------------- ! 95 z_frc_trd_v = r1_rau0 * glob_sum( - ( emp(:,:) - rnf(:,:) + rdivisf *fwfisf(:,:) ) * surf(:,:) ) ! volume fluxes95 z_frc_trd_v = r1_rau0 * glob_sum( - ( emp(:,:) - rnf(:,:) + fwfisf(:,:) ) * surf(:,:) ) ! volume fluxes 96 96 z_frc_trd_t = glob_sum( sbc_tsc(:,:,jp_tem) * surf(:,:) ) ! heat fluxes 97 97 z_frc_trd_s = glob_sum( sbc_tsc(:,:,jp_sal) * surf(:,:) ) ! salt fluxes … … 101 101 ! Add ice shelf heat & salt input 102 102 IF( nn_isf .GE. 1 ) THEN 103 z_frc_trd_t = z_frc_trd_t & 104 & + glob_sum( ( risf_tsc(:,:,jp_tem) - rdivisf * fwfisf(:,:) * (-1.9) * r1_rau0 ) * surf(:,:) ) 105 z_frc_trd_s = z_frc_trd_s + (1.0_wp - rdivisf) * glob_sum( risf_tsc(:,:,jp_sal) * surf(:,:) ) 103 z_frc_trd_t = z_frc_trd_t + glob_sum( risf_tsc(:,:,jp_tem) * surf(:,:) ) 104 z_frc_trd_s = z_frc_trd_s + glob_sum( risf_tsc(:,:,jp_sal) * surf(:,:) ) 106 105 ENDIF 107 106 … … 200 199 ! ENDIF 201 200 !!gm end 202 203 201 204 202 IF( lk_vvl ) THEN -
branches/UKMO/dev_r5518_obs_oper_update_icethick/NEMOGCM/NEMO/OPA_SRC/DIA/diaptr.F90
r7960 r9987 9 9 !! 3.3 ! 2010-10 (G. Madec) dynamical allocation 10 10 !! 3.6 ! 2014-12 (C. Ethe) use of IOM 11 !! 3.6 ! 2016-06 (T. Graham) Addition of diagnostics for CMIP6 11 12 !!---------------------------------------------------------------------- 12 13 … … 21 22 USE dom_oce ! ocean space and time domain 22 23 USE phycst ! physical constants 24 USE ldftra_oce 23 25 ! 24 26 USE iom ! IOM library … … 38 40 PUBLIC dia_ptr_init ! call in step module 39 41 PUBLIC dia_ptr ! call in step module 42 PUBLIC dia_ptr_ohst_components ! called from tra_ldf/tra_adv routines 40 43 41 44 ! !!** namelist namptr ** 42 REAL(wp), ALLOCATABLE, SAVE, PUBLIC, DIMENSION(:) :: htr_adv, htr_ldf !: Heat TRansports (adv, diff, overturn.) 43 REAL(wp), ALLOCATABLE, SAVE, PUBLIC, DIMENSION(:) :: str_adv, str_ldf !: Salt TRansports (adv, diff, overturn.) 44 45 REAL(wp), ALLOCATABLE, SAVE, PUBLIC, DIMENSION(:,:) :: htr_adv, htr_ldf, htr_eiv, htr_vt !: Heat TRansports (adv, diff, Bolus.) 46 REAL(wp), ALLOCATABLE, SAVE, PUBLIC, DIMENSION(:,:) :: str_adv, str_ldf, str_eiv, str_vs !: Salt TRansports (adv, diff, Bolus.) 47 REAL(wp), ALLOCATABLE, SAVE, PUBLIC, DIMENSION(:,:) :: htr_ove, str_ove !: heat Salt TRansports ( overturn.) 48 REAL(wp), ALLOCATABLE, SAVE, PUBLIC, DIMENSION(:,:) :: htr_btr, str_btr !: heat Salt TRansports ( barotropic ) 45 49 46 50 LOGICAL, PUBLIC :: ln_diaptr ! Poleward transport flag (T) or not (F) 47 51 LOGICAL, PUBLIC :: ln_subbas ! Atlantic/Pacific/Indian basins calculation 48 INTEGER 52 INTEGER, PUBLIC :: nptr ! = 1 (l_subbas=F) or = 5 (glo, atl, pac, ind, ipc) (l_subbas=T) 49 53 50 54 REAL(wp) :: rc_sv = 1.e-6_wp ! conversion from m3/s to Sverdrup … … 77 81 ! 78 82 INTEGER :: ji, jj, jk, jn ! dummy loop indices 79 REAL(wp) :: z v, zsfc ! local scalar83 REAL(wp) :: zsfc,zvfc ! local scalar 80 84 REAL(wp), DIMENSION(jpi,jpj) :: z2d ! 2D workspace 81 85 REAL(wp), DIMENSION(jpi,jpj,jpk) :: z3d ! 3D workspace 82 86 REAL(wp), DIMENSION(jpi,jpj,jpk) :: zmask ! 3D workspace 83 87 REAL(wp), DIMENSION(jpi,jpj,jpk,jpts) :: zts ! 3D workspace 84 CHARACTER( len = 10 ) :: cl1 88 REAL(wp), DIMENSION(jpj) :: vsum ! 1D workspace 89 REAL(wp), DIMENSION(jpj,jpts) :: tssum ! 1D workspace 90 91 ! 92 !overturning calculation 93 REAL(wp), DIMENSION(jpj,jpk,nptr) :: sjk , r1_sjk ! i-mean i-k-surface and its inverse 94 REAL(wp), DIMENSION(jpj,jpk,nptr) :: v_msf, sn_jk , tn_jk ! i-mean T and S, j-Stream-Function 95 REAL(wp), DIMENSION(jpi,jpj,jpk) :: zvn ! 3D workspace 96 97 98 CHARACTER( len = 12 ) :: cl1 85 99 !!---------------------------------------------------------------------- 86 100 ! … … 88 102 89 103 ! 104 z3d(:,:,:) = 0._wp 90 105 IF( PRESENT( pvtr ) ) THEN 91 106 IF( iom_use("zomsfglo") ) THEN ! effective MSF 92 107 z3d(1,:,:) = ptr_sjk( pvtr(:,:,:) ) ! zonal cumulative effective transport 93 DO jk = 2, jpkm194 z3d(1,:,jk) = z3d(1,:,jk -1) +z3d(1,:,jk) ! effective j-Stream-Function (MSF)108 DO jk = jpkm1,1,-1 !Integrate from bottom up to get 109 z3d(1,:,jk) = z3d(1,:,jk+1) - z3d(1,:,jk) ! effective j-Stream-Function (MSF) 95 110 END DO 96 111 DO ji = 1, jpi … … 100 115 CALL iom_put( cl1, z3d * rc_sv ) 101 116 DO jn = 2, nptr ! by sub-basins 102 z3d(1,:,:) = ptr_sjk( pvtr(:,:,:), btmsk(:,:,jn) *btm30(:,:))103 DO jk = 2, jpkm1104 z3d(1,:,jk) = z3d(1,:,jk -1) +z3d(1,:,jk) ! effective j-Stream-Function (MSF)117 z3d(1,:,:) = ptr_sjk( pvtr(:,:,:), btmsk(:,:,jn) ) 118 DO jk = jpkm1,1,-1 119 z3d(1,:,jk) = z3d(1,:,jk+1) - z3d(1,:,jk) ! effective j-Stream-Function (MSF) 105 120 END DO 106 121 DO ji = 1, jpi … … 111 126 END DO 112 127 ENDIF 128 IF( iom_use("sopstove") .OR. iom_use("sophtove") .OR. iom_use("sopstbtr") .OR. iom_use("sophtbtr") ) THEN 129 ! define fields multiplied by scalar 130 zmask(:,:,:) = 0._wp 131 zts(:,:,:,:) = 0._wp 132 zvn(:,:,:) = 0._wp 133 DO jk = 1, jpkm1 134 DO jj = 1, jpjm1 135 DO ji = 1, jpi 136 zvfc = e1v(ji,jj) * fse3v(ji,jj,jk) 137 zmask(ji,jj,jk) = vmask(ji,jj,jk) * zvfc 138 zts(ji,jj,jk,jp_tem) = (tsn(ji,jj,jk,jp_tem)+tsn(ji,jj+1,jk,jp_tem)) * 0.5 * zvfc !Tracers averaged onto V grid 139 zts(ji,jj,jk,jp_sal) = (tsn(ji,jj,jk,jp_sal)+tsn(ji,jj+1,jk,jp_sal)) * 0.5 * zvfc 140 zvn(ji,jj,jk) = vn(ji,jj,jk) * zvfc 141 ENDDO 142 ENDDO 143 ENDDO 144 ENDIF 145 IF( iom_use("sopstove") .OR. iom_use("sophtove") ) THEN 146 sjk(:,:,1) = ptr_sjk( zmask(:,:,:), btmsk(:,:,1) ) 147 r1_sjk(:,:,1) = 0._wp 148 WHERE( sjk(:,:,1) /= 0._wp ) r1_sjk(:,:,1) = 1._wp / sjk(:,:,1) 149 150 ! i-mean T and S, j-Stream-Function, global 151 tn_jk(:,:,1) = ptr_sjk( zts(:,:,:,jp_tem) ) * r1_sjk(:,:,1) 152 sn_jk(:,:,1) = ptr_sjk( zts(:,:,:,jp_sal) ) * r1_sjk(:,:,1) 153 v_msf(:,:,1) = ptr_sjk( zvn(:,:,:) ) 154 155 htr_ove(:,1) = SUM( v_msf(:,:,1)*tn_jk(:,:,1) ,2 ) 156 str_ove(:,1) = SUM( v_msf(:,:,1)*sn_jk(:,:,1) ,2 ) 157 158 z2d(1,:) = htr_ove(:,1) * rc_pwatt ! (conversion in PW) 159 DO ji = 1, jpi 160 z2d(ji,:) = z2d(1,:) 161 ENDDO 162 cl1 = 'sophtove' 163 CALL iom_put( TRIM(cl1), z2d ) 164 z2d(1,:) = str_ove(:,1) * rc_ggram ! (conversion in Gg) 165 DO ji = 1, jpi 166 z2d(ji,:) = z2d(1,:) 167 ENDDO 168 cl1 = 'sopstove' 169 CALL iom_put( TRIM(cl1), z2d ) 170 IF( ln_subbas ) THEN 171 DO jn = 2, nptr 172 sjk(:,:,jn) = ptr_sjk( zmask(:,:,:), btmsk(:,:,jn) ) 173 r1_sjk(:,:,jn) = 0._wp 174 WHERE( sjk(:,:,jn) /= 0._wp ) r1_sjk(:,:,jn) = 1._wp / sjk(:,:,jn) 175 176 ! i-mean T and S, j-Stream-Function, basin 177 tn_jk(:,:,jn) = ptr_sjk( zts(:,:,:,jp_tem), btmsk(:,:,jn) ) * r1_sjk(:,:,jn) 178 sn_jk(:,:,jn) = ptr_sjk( zts(:,:,:,jp_sal), btmsk(:,:,jn) ) * r1_sjk(:,:,jn) 179 v_msf(:,:,jn) = ptr_sjk( zvn(:,:,:), btmsk(:,:,jn) ) 180 htr_ove(:,jn) = SUM( v_msf(:,:,jn)*tn_jk(:,:,jn) ,2 ) 181 str_ove(:,jn) = SUM( v_msf(:,:,jn)*sn_jk(:,:,jn) ,2 ) 182 183 z2d(1,:) = htr_ove(:,jn) * rc_pwatt ! (conversion in PW) 184 DO ji = 1, jpi 185 z2d(ji,:) = z2d(1,:) 186 ENDDO 187 cl1 = TRIM('sophtove_'//clsubb(jn)) 188 CALL iom_put( cl1, z2d ) 189 z2d(1,:) = str_ove(:,jn) * rc_ggram ! (conversion in Gg) 190 DO ji = 1, jpi 191 z2d(ji,:) = z2d(1,:) 192 ENDDO 193 cl1 = TRIM('sopstove_'//clsubb(jn)) 194 CALL iom_put( cl1, z2d ) 195 END DO 196 ENDIF 197 ENDIF 198 IF( iom_use("sopstbtr") .OR. iom_use("sophtbtr") ) THEN 199 ! Calculate barotropic heat and salt transport here 200 sjk(:,1,1) = ptr_sj( zmask(:,:,:), btmsk(:,:,1) ) 201 r1_sjk(:,1,1) = 0._wp 202 WHERE( sjk(:,1,1) /= 0._wp ) r1_sjk(:,1,1) = 1._wp / sjk(:,1,1) 203 204 vsum = ptr_sj( zvn(:,:,:), btmsk(:,:,1)) 205 tssum(:,jp_tem) = ptr_sj( zts(:,:,:,jp_tem), btmsk(:,:,1) ) 206 tssum(:,jp_sal) = ptr_sj( zts(:,:,:,jp_sal), btmsk(:,:,1) ) 207 htr_btr(:,1) = vsum * tssum(:,jp_tem) * r1_sjk(:,1,1) 208 str_btr(:,1) = vsum * tssum(:,jp_sal) * r1_sjk(:,1,1) 209 z2d(1,:) = htr_btr(:,1) * rc_pwatt ! (conversion in PW) 210 DO ji = 2, jpi 211 z2d(ji,:) = z2d(1,:) 212 ENDDO 213 cl1 = 'sophtbtr' 214 CALL iom_put( TRIM(cl1), z2d ) 215 z2d(1,:) = str_btr(:,1) * rc_ggram ! (conversion in Gg) 216 DO ji = 2, jpi 217 z2d(ji,:) = z2d(1,:) 218 ENDDO 219 cl1 = 'sopstbtr' 220 CALL iom_put( TRIM(cl1), z2d ) 221 IF( ln_subbas ) THEN 222 DO jn = 2, nptr 223 sjk(:,1,jn) = ptr_sj( zmask(:,:,:), btmsk(:,:,jn) ) 224 r1_sjk(:,1,jn) = 0._wp 225 WHERE( sjk(:,1,jn) /= 0._wp ) r1_sjk(:,1,jn) = 1._wp / sjk(:,1,jn) 226 vsum = ptr_sj( zvn(:,:,:), btmsk(:,:,jn)) 227 tssum(:,jp_tem) = ptr_sj( zts(:,:,:,jp_tem), btmsk(:,:,jn) ) 228 tssum(:,jp_sal) = ptr_sj( zts(:,:,:,jp_sal), btmsk(:,:,jn) ) 229 htr_btr(:,jn) = vsum * tssum(:,jp_tem) * r1_sjk(:,1,jn) 230 str_btr(:,jn) = vsum * tssum(:,jp_sal) * r1_sjk(:,1,jn) 231 z2d(1,:) = htr_btr(:,jn) * rc_pwatt ! (conversion in PW) 232 DO ji = 1, jpi 233 z2d(ji,:) = z2d(1,:) 234 ENDDO 235 cl1 = TRIM('sophtbtr_'//clsubb(jn)) 236 CALL iom_put( cl1, z2d ) 237 z2d(1,:) = str_btr(:,jn) * rc_ggram ! (conversion in Gg) 238 DO ji = 1, jpi 239 z2d(ji,:) = z2d(1,:) 240 ENDDO 241 cl1 = TRIM('sopstbtr_'//clsubb(jn)) 242 CALL iom_put( cl1, z2d ) 243 ENDDO 244 ENDIF !ln_subbas 245 ENDIF !iom_use("sopstbtr....) 113 246 ! 114 247 ELSE … … 150 283 ! ! Advective and diffusive heat and salt transport 151 284 IF( iom_use("sophtadv") .OR. iom_use("sopstadv") ) THEN 152 z2d(1,:) = htr_adv(: ) * rc_pwatt ! (conversion in PW)285 z2d(1,:) = htr_adv(:,1) * rc_pwatt ! (conversion in PW) 153 286 DO ji = 1, jpi 154 287 z2d(ji,:) = z2d(1,:) … … 156 289 cl1 = 'sophtadv' 157 290 CALL iom_put( TRIM(cl1), z2d ) 158 z2d(1,:) = str_adv(: ) * rc_ggram ! (conversion in Gg)291 z2d(1,:) = str_adv(:,1) * rc_ggram ! (conversion in Gg) 159 292 DO ji = 1, jpi 160 293 z2d(ji,:) = z2d(1,:) … … 162 295 cl1 = 'sopstadv' 163 296 CALL iom_put( TRIM(cl1), z2d ) 297 IF( ln_subbas ) THEN 298 DO jn=2,nptr 299 z2d(1,:) = htr_adv(:,jn) * rc_pwatt ! (conversion in PW) 300 DO ji = 1, jpi 301 z2d(ji,:) = z2d(1,:) 302 ENDDO 303 cl1 = TRIM('sophtadv_'//clsubb(jn)) 304 CALL iom_put( cl1, z2d ) 305 z2d(1,:) = str_adv(:,jn) * rc_ggram ! (conversion in Gg) 306 DO ji = 1, jpi 307 z2d(ji,:) = z2d(1,:) 308 ENDDO 309 cl1 = TRIM('sopstadv_'//clsubb(jn)) 310 CALL iom_put( cl1, z2d ) 311 ENDDO 312 ENDIF 164 313 ENDIF 165 314 ! 166 315 IF( iom_use("sophtldf") .OR. iom_use("sopstldf") ) THEN 167 z2d(1,:) = htr_ldf(: ) * rc_pwatt ! (conversion in PW)316 z2d(1,:) = htr_ldf(:,1) * rc_pwatt ! (conversion in PW) 168 317 DO ji = 1, jpi 169 318 z2d(ji,:) = z2d(1,:) … … 171 320 cl1 = 'sophtldf' 172 321 CALL iom_put( TRIM(cl1), z2d ) 173 z2d(1,:) = str_ldf(: ) * rc_ggram ! (conversion in Gg)322 z2d(1,:) = str_ldf(:,1) * rc_ggram ! (conversion in Gg) 174 323 DO ji = 1, jpi 175 324 z2d(ji,:) = z2d(1,:) … … 177 326 cl1 = 'sopstldf' 178 327 CALL iom_put( TRIM(cl1), z2d ) 179 ENDIF 328 IF( ln_subbas ) THEN 329 DO jn=2,nptr 330 z2d(1,:) = htr_ldf(:,jn) * rc_pwatt ! (conversion in PW) 331 DO ji = 1, jpi 332 z2d(ji,:) = z2d(1,:) 333 ENDDO 334 cl1 = TRIM('sophtldf_'//clsubb(jn)) 335 CALL iom_put( cl1, z2d ) 336 z2d(1,:) = str_ldf(:,jn) * rc_ggram ! (conversion in Gg) 337 DO ji = 1, jpi 338 z2d(ji,:) = z2d(1,:) 339 ENDDO 340 cl1 = TRIM('sopstldf_'//clsubb(jn)) 341 CALL iom_put( cl1, z2d ) 342 ENDDO 343 ENDIF 344 ENDIF 345 346 IF( iom_use("sopht_vt") .OR. iom_use("sopst_vs") ) THEN 347 z2d(1,:) = htr_vt(:,1) * rc_pwatt ! (conversion in PW) 348 DO ji = 1, jpi 349 z2d(ji,:) = z2d(1,:) 350 ENDDO 351 cl1 = 'sopht_vt' 352 CALL iom_put( TRIM(cl1), z2d ) 353 z2d(1,:) = str_vs(:,1) * rc_ggram ! (conversion in Gg) 354 DO ji = 1, jpi 355 z2d(ji,:) = z2d(1,:) 356 ENDDO 357 cl1 = 'sopst_vs' 358 CALL iom_put( TRIM(cl1), z2d ) 359 IF( ln_subbas ) THEN 360 DO jn=2,nptr 361 z2d(1,:) = htr_vt(:,jn) * rc_pwatt ! (conversion in PW) 362 DO ji = 1, jpi 363 z2d(ji,:) = z2d(1,:) 364 ENDDO 365 cl1 = TRIM('sopht_vt_'//clsubb(jn)) 366 CALL iom_put( cl1, z2d ) 367 z2d(1,:) = str_vs(:,jn) * rc_ggram ! (conversion in Gg) 368 DO ji = 1, jpi 369 z2d(ji,:) = z2d(1,:) 370 ENDDO 371 cl1 = TRIM('sopst_vs_'//clsubb(jn)) 372 CALL iom_put( cl1, z2d ) 373 ENDDO 374 ENDIF 375 ENDIF 376 377 #ifdef key_diaeiv 378 IF(lk_traldf_eiv) THEN 379 IF( iom_use("sophteiv") .OR. iom_use("sopsteiv") ) THEN 380 z2d(1,:) = htr_eiv(:,1) * rc_pwatt ! (conversion in PW) 381 DO ji = 1, jpi 382 z2d(ji,:) = z2d(1,:) 383 ENDDO 384 cl1 = 'sophteiv' 385 CALL iom_put( TRIM(cl1), z2d ) 386 z2d(1,:) = str_eiv(:,1) * rc_ggram ! (conversion in Gg) 387 DO ji = 1, jpi 388 z2d(ji,:) = z2d(1,:) 389 ENDDO 390 cl1 = 'sopsteiv' 391 CALL iom_put( TRIM(cl1), z2d ) 392 IF( ln_subbas ) THEN 393 DO jn=2,nptr 394 z2d(1,:) = htr_eiv(:,jn) * rc_pwatt ! (conversion in PW) 395 DO ji = 1, jpi 396 z2d(ji,:) = z2d(1,:) 397 ENDDO 398 cl1 = TRIM('sophteiv_'//clsubb(jn)) 399 CALL iom_put( cl1, z2d ) 400 z2d(1,:) = str_eiv(:,jn) * rc_ggram ! (conversion in Gg) 401 DO ji = 1, jpi 402 z2d(ji,:) = z2d(1,:) 403 ENDDO 404 cl1 = TRIM('sopsteiv_'//clsubb(jn)) 405 CALL iom_put( cl1, z2d ) 406 ENDDO 407 ENDIF 408 ENDIF 409 IF( iom_use("zomsfeivglo") ) THEN 410 z3d(1,:,:) = ptr_sjk( v_eiv(:,:,:) ) ! zonal cumulative effective transport 411 DO jk = jpkm1,1,-1 412 z3d(1,:,jk) = z3d(1,:,jk+1) - z3d(1,:,jk) ! effective j-Stream-Function (MSF) 413 END DO 414 DO ji = 1, jpi 415 z3d(ji,:,:) = z3d(1,:,:) 416 ENDDO 417 cl1 = TRIM('zomsfeiv'//clsubb(1) ) 418 CALL iom_put( cl1, z3d * rc_sv ) 419 IF( ln_subbas ) THEN 420 DO jn = 2, nptr ! by sub-basins 421 z3d(1,:,:) = ptr_sjk( v_eiv(:,:,:), btmsk(:,:,jn) ) 422 DO jk = jpkm1,1,-1 423 z3d(1,:,jk) = z3d(1,:,jk+1) - z3d(1,:,jk) ! effective j-Stream-Function (MSF) 424 END DO 425 DO ji = 1, jpi 426 z3d(ji,:,:) = z3d(1,:,:) 427 ENDDO 428 cl1 = TRIM('zomsfeiv'//clsubb(jn) ) 429 CALL iom_put( cl1, z3d * rc_sv ) 430 END DO 431 ENDIF 432 ENDIF 433 ENDIF 434 #endif 180 435 ! 181 436 ENDIF … … 256 511 ! Initialise arrays to zero because diatpr is called before they are first calculated 257 512 ! Note that this means diagnostics will not be exactly correct when model run is restarted. 258 htr_adv(:) = 0._wp ; str_adv(:) = 0._wp 259 htr_ldf(:) = 0._wp ; str_ldf(:) = 0._wp 513 htr_adv(:,:) = 0._wp ; str_adv(:,:) = 0._wp 514 htr_ldf(:,:) = 0._wp ; str_ldf(:,:) = 0._wp 515 htr_eiv(:,:) = 0._wp ; str_eiv(:,:) = 0._wp 516 htr_vt(:,:) = 0._wp ; str_vs(:,:) = 0._wp 517 htr_ove(:,:) = 0._wp ; str_ove(:,:) = 0._wp 518 htr_btr(:,:) = 0._wp ; str_btr(:,:) = 0._wp 260 519 ! 261 520 ENDIF … … 263 522 END SUBROUTINE dia_ptr_init 264 523 524 SUBROUTINE dia_ptr_ohst_components( ktra, cptr, pva ) 525 !!---------------------------------------------------------------------- 526 !! *** ROUTINE dia_ptr_ohst_components *** 527 !!---------------------------------------------------------------------- 528 !! Wrapper for heat and salt transport calculations to calculate them for each basin 529 !! Called from all advection and/or diffusion routines 530 !!---------------------------------------------------------------------- 531 INTEGER , INTENT(in ) :: ktra ! tracer index 532 CHARACTER(len=3) , INTENT(in) :: cptr ! transport type 'adv'/'ldf'/'eiv' 533 REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in) :: pva ! 3D input array of advection/diffusion 534 INTEGER :: jn ! 535 536 IF( cptr == 'adv' ) THEN 537 IF( ktra == jp_tem ) htr_adv(:,1) = ptr_sj( pva(:,:,:) ) 538 IF( ktra == jp_sal ) str_adv(:,1) = ptr_sj( pva(:,:,:) ) 539 ENDIF 540 IF( cptr == 'ldf' ) THEN 541 IF( ktra == jp_tem ) htr_ldf(:,1) = ptr_sj( pva(:,:,:) ) 542 IF( ktra == jp_sal ) str_ldf(:,1) = ptr_sj( pva(:,:,:) ) 543 ENDIF 544 IF( cptr == 'eiv' ) THEN 545 IF( ktra == jp_tem ) htr_eiv(:,1) = ptr_sj( pva(:,:,:) ) 546 IF( ktra == jp_sal ) str_eiv(:,1) = ptr_sj( pva(:,:,:) ) 547 ENDIF 548 IF( cptr == 'vts' ) THEN 549 IF( ktra == jp_tem ) htr_vt(:,1) = ptr_sj( pva(:,:,:) ) 550 IF( ktra == jp_sal ) str_vs(:,1) = ptr_sj( pva(:,:,:) ) 551 ENDIF 552 ! 553 IF( ln_subbas ) THEN 554 ! 555 IF( cptr == 'adv' ) THEN 556 IF( ktra == jp_tem ) THEN 557 DO jn = 2, nptr 558 htr_adv(:,jn) = ptr_sj( pva(:,:,:), btmsk(:,:,jn) ) 559 END DO 560 ENDIF 561 IF( ktra == jp_sal ) THEN 562 DO jn = 2, nptr 563 str_adv(:,jn) = ptr_sj( pva(:,:,:), btmsk(:,:,jn) ) 564 END DO 565 ENDIF 566 ENDIF 567 IF( cptr == 'ldf' ) THEN 568 IF( ktra == jp_tem ) THEN 569 DO jn = 2, nptr 570 htr_ldf(:,jn) = ptr_sj( pva(:,:,:), btmsk(:,:,jn) ) 571 END DO 572 ENDIF 573 IF( ktra == jp_sal ) THEN 574 DO jn = 2, nptr 575 str_ldf(:,jn) = ptr_sj( pva(:,:,:), btmsk(:,:,jn) ) 576 END DO 577 ENDIF 578 ENDIF 579 IF( cptr == 'eiv' ) THEN 580 IF( ktra == jp_tem ) THEN 581 DO jn = 2, nptr 582 htr_eiv(:,jn) = ptr_sj( pva(:,:,:), btmsk(:,:,jn) ) 583 END DO 584 ENDIF 585 IF( ktra == jp_sal ) THEN 586 DO jn = 2, nptr 587 str_eiv(:,jn) = ptr_sj( pva(:,:,:), btmsk(:,:,jn) ) 588 END DO 589 ENDIF 590 ENDIF 591 IF( cptr == 'vts' ) THEN 592 IF( ktra == jp_tem ) THEN 593 DO jn = 2, nptr 594 htr_vt(:,jn) = ptr_sj( pva(:,:,:), btmsk(:,:,jn) ) 595 END DO 596 ENDIF 597 IF( ktra == jp_sal ) THEN 598 DO jn = 2, nptr 599 str_vs(:,jn) = ptr_sj( pva(:,:,:), btmsk(:,:,jn) ) 600 END DO 601 ENDIF 602 ENDIF 603 ! 604 ENDIF 605 END SUBROUTINE dia_ptr_ohst_components 606 265 607 266 608 FUNCTION dia_ptr_alloc() … … 273 615 ierr(:) = 0 274 616 ! 275 ALLOCATE( btmsk(jpi,jpj,nptr) , & 276 & htr_adv(jpj) , str_adv(jpj) , & 277 & htr_ldf(jpj) , str_ldf(jpj) , STAT=ierr(1) ) 617 ALLOCATE( btmsk(jpi,jpj,nptr) , & 618 & htr_adv(jpj,nptr) , str_adv(jpj,nptr) , & 619 & htr_eiv(jpj,nptr) , str_eiv(jpj,nptr) , & 620 & htr_vt(jpj,nptr) , str_vs(jpj,nptr) , & 621 & htr_ove(jpj,nptr) , str_ove(jpj,nptr) , & 622 & htr_btr(jpj,nptr) , str_btr(jpj,nptr) , & 623 & htr_ldf(jpj,nptr) , str_ldf(jpj,nptr) , STAT=ierr(1) ) 278 624 ! 279 625 ALLOCATE( p_fval1d(jpj), p_fval2d(jpj,jpk), Stat=ierr(2)) … … 402 748 #endif 403 749 !!-------------------------------------------------------------------- 404 750 ! 405 751 p_fval => p_fval2d 406 752 … … 434 780 #endif 435 781 ! 782 436 783 END FUNCTION ptr_sjk 437 784 -
branches/UKMO/dev_r5518_obs_oper_update_icethick/NEMOGCM/NEMO/OPA_SRC/DIA/diawri.F90
r7960 r9987 39 39 USE zdfmxl ! mixed layer 40 40 USE dianam ! build name of file (routine) 41 USE zdftke ! vertical physics: one-equation scheme 41 42 USE zdfddm ! vertical physics: double diffusion 42 43 USE diahth ! thermocline diagnostics … … 46 47 USE iom 47 48 USE ioipsl 48 USE dynspg_oce, ONLY: un_adv, vn_adv ! barotropic velocities 49 49 USE dynspg_oce, ONLY: un_adv, vn_adv ! barotropic velocities 50 USE insitu_tem, ONLY: insitu_t, theta2t 50 51 #if defined key_lim2 51 52 USE limwri_2 … … 145 146 ENDIF 146 147 147 IF( .NOT.lk_vvl ) THEN 148 CALL iom_put( "e3t" , fse3t_n(:,:,:) ) 149 CALL iom_put( "e3u" , fse3u_n(:,:,:) ) 150 CALL iom_put( "e3v" , fse3v_n(:,:,:) ) 151 CALL iom_put( "e3w" , fse3w_n(:,:,:) ) 152 ENDIF 148 ! Output of initial vertical scale factor 149 CALL iom_put("e3t_0", e3t_0(:,:,:) ) 150 CALL iom_put("e3u_0", e3t_0(:,:,:) ) 151 CALL iom_put("e3v_0", e3t_0(:,:,:) ) 152 ! 153 CALL iom_put( "e3t" , fse3t_n(:,:,:) ) 154 CALL iom_put( "e3u" , fse3u_n(:,:,:) ) 155 CALL iom_put( "e3v" , fse3v_n(:,:,:) ) 156 CALL iom_put( "e3w" , fse3w_n(:,:,:) ) 157 IF( iom_use("e3tdef") ) & 158 CALL iom_put( "e3tdef" , ( ( fse3t_n(:,:,:) - e3t_0(:,:,:) ) / e3t_0(:,:,:) * 100 * tmask(:,:,:) ) ** 2 ) 159 CALL iom_put("tpt_dep", fsdept_n(:,:,:) ) 160 CALL iom_put("wpt_dep", fsdepw_n(:,:,:) ) 161 153 162 154 163 CALL iom_put( "ssh" , sshn ) ! sea surface height 155 if( iom_use('ssh2') ) CALL iom_put( "ssh2", sshn(:,:) * sshn(:,:) ) ! square of sea surface height156 164 157 165 CALL iom_put( "toce", tsn(:,:,:,jp_tem) ) ! 3D temperature 166 CALL theta2t ! in-situ temperature conversion 167 CALL iom_put( "tinsitu", insitu_t(:,:,:)) ! in-situ temperature 158 168 CALL iom_put( "sst", tsn(:,:,1,jp_tem) ) ! surface temperature 159 169 IF ( iom_use("sbt") ) THEN … … 194 204 CALL iom_put( "taubot", z2d ) 195 205 ENDIF 196 206 197 207 CALL iom_put( "uoce", un(:,:,:) ) ! 3D i-current 198 208 CALL iom_put( "ssu", un(:,:,1) ) ! surface i-current … … 242 252 CALL iom_put( "avt" , avt ) ! T vert. eddy diff. coef. 243 253 CALL iom_put( "avm" , avmu ) ! T vert. eddy visc. coef. 254 IF( lk_zdftke ) THEN 255 CALL iom_put( "tke" , en ) ! TKE budget: Turbulent Kinetic Energy 256 CALL iom_put( "tke_niw" , e_niw ) ! TKE budget: Near-inertial waves 257 ENDIF 244 258 CALL iom_put( "avs" , fsavs(:,:,:) ) ! S vert. eddy diff. coef. (useful only with key_zdfddm) 259 ! Log of eddy diff coef 260 IF( iom_use('logavt') ) CALL iom_put( "logavt", LOG( MAX( 1.e-20_wp, avt (:,:,:) ) ) ) 261 IF( iom_use('logavs') ) CALL iom_put( "logavs", LOG( MAX( 1.e-20_wp, fsavs(:,:,:) ) ) ) 245 262 246 263 IF ( iom_use("sstgrad") .OR. iom_use("sstgrad2") ) THEN … … 307 324 CALL iom_put( "eken", rke ) 308 325 ENDIF 309 310 IF( iom_use("u_masstr") .OR. iom_use("u_heattr") .OR. iom_use("u_salttr") ) THEN 326 ! 327 CALL iom_put( "hdiv", hdivn ) ! Horizontal divergence 328 ! 329 IF( iom_use("u_masstr") .OR. iom_use("u_masstr_vint") .OR. iom_use("u_heattr") .OR. iom_use("u_salttr") ) THEN 311 330 z3d(:,:,jpk) = 0.e0 331 z2d(:,:) = 0.e0 312 332 DO jk = 1, jpkm1 313 333 z3d(:,:,jk) = rau0 * un(:,:,jk) * e2u(:,:) * fse3u(:,:,jk) * umask(:,:,jk) 334 z2d(:,:) = z2d(:,:) + z3d(:,:,jk) 314 335 END DO 315 336 CALL iom_put( "u_masstr", z3d ) ! mass transport in i-direction 337 CALL iom_put( "u_masstr_vint", z2d ) ! mass transport in i-direction vertical sum 316 338 ENDIF 317 339 … … 376 398 CALL iom_put( "v_salttr", 0.5 * z2d ) ! heat transport in j-direction 377 399 ENDIF 400 401 ! Vertical integral of temperature 402 IF( iom_use("tosmint") ) THEN 403 z2d(:,:)=0._wp 404 DO jk = 1, jpkm1 405 DO jj = 2, jpjm1 406 DO ji = fs_2, fs_jpim1 ! vector opt. 407 z2d(ji,jj) = z2d(ji,jj) + rau0 * fse3t(ji,jj,jk) * tsn(ji,jj,jk,jp_tem) 408 END DO 409 END DO 410 END DO 411 CALL lbc_lnk( z2d, 'T', -1. ) 412 CALL iom_put( "tosmint", z2d ) 413 ENDIF 414 415 ! Vertical integral of salinity 416 IF( iom_use("somint") ) THEN 417 z2d(:,:)=0._wp 418 DO jk = 1, jpkm1 419 DO jj = 2, jpjm1 420 DO ji = fs_2, fs_jpim1 ! vector opt. 421 z2d(ji,jj) = z2d(ji,jj) + rau0 * fse3t(ji,jj,jk) * tsn(ji,jj,jk,jp_sal) 422 END DO 423 END DO 424 END DO 425 CALL lbc_lnk( z2d, 'T', -1. ) 426 CALL iom_put( "somint", z2d ) 427 ENDIF 428 429 CALL iom_put( "bn2", rn2 ) !Brunt-Vaisala buoyancy frequency (N^2) 378 430 ! 379 431 CALL wrk_dealloc( jpi , jpj , z2d ) … … 438 490 zdt = rdt 439 491 IF( nacc == 1 ) zdt = rdtmin 440 IF( ln_mskland ) THEN ; clop = "only(x)" ! put 1.e+20 on land (very expensive!!) 441 ELSE ; clop = "x" ! no use of the mask value (require less cpu time) 442 ENDIF 492 clop = "x" ! no use of the mask value (require less cpu time, and otherwise the model crashes) 443 493 #if defined key_diainstant 444 494 zsto = nwrite * zdt … … 1020 1070 CALL histdef( id_i, "vovvldep", "T point depth" , "m" , & ! t-point depth 1021 1071 & jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout ) 1072 CALL histdef( id_i, "vovvle3t", "T point thickness" , "m" , & ! t-point depth 1073 & jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout ) 1022 1074 END IF 1023 1075 … … 1050 1102 CALL histwrite( id_i, "sozotaux", kt, utau , jpi*jpj , idex ) ! i-wind stress 1051 1103 CALL histwrite( id_i, "sometauy", kt, vtau , jpi*jpj , idex ) ! j-wind stress 1104 IF( lk_vvl ) THEN 1105 CALL histwrite( id_i, "vovvldep", kt, fsdept_n(:,:,:), jpi*jpj*jpk, idex )! T-cell depth 1106 CALL histwrite( id_i, "vovvle3t", kt, fse3t_n (:,:,:), jpi*jpj*jpk, idex )! T-cell thickness 1107 END IF 1052 1108 1053 1109 ! 3. Close the file … … 1062 1118 ENDIF 1063 1119 #endif 1120 1121 IF (cdfile_name == "output.abort") THEN 1122 CALL ctl_stop('MPPSTOP', 'NEMO abort from dia_wri_state') 1123 END IF 1064 1124 1065 1125 ! IF( nn_timing == 1 ) CALL timing_stop('dia_wri_state') ! not sure this works for routines not called in first timestep -
branches/UKMO/dev_r5518_obs_oper_update_icethick/NEMOGCM/NEMO/OPA_SRC/DIA/diawri_dimg.h90
r7960 r9987 112 112 IF( inbsel > jpk ) THEN 113 113 IF(lwp) WRITE(numout,*) ' STOP inbsel =',inbsel,' is larger than jpk=',jpk 114 STOP114 CALL ctl_stop('STOP', 'NEMO aborted from dia_wri') 115 115 ENDIF 116 116
Note: See TracChangeset
for help on using the changeset viewer.