Changeset 12377 for NEMO/trunk/src/OCE/DIA/diawri.F90
- Timestamp:
- 2020-02-12T15:39:06+01:00 (4 years ago)
- Location:
- NEMO/trunk
- Files:
-
- 2 edited
Legend:
- Unmodified
- Added
- Removed
-
NEMO/trunk
- Property svn:externals
-
old new 3 3 ^/utils/build/mk@HEAD mk 4 4 ^/utils/tools@HEAD tools 5 ^/vendors/AGRIF/dev @HEAD ext/AGRIF5 ^/vendors/AGRIF/dev_r11615_ENHANCE-04_namelists_as_internalfiles_agrif@HEAD ext/AGRIF 6 6 ^/vendors/FCM@HEAD ext/FCM 7 7 ^/vendors/IOIPSL@HEAD ext/IOIPSL
-
- Property svn:externals
-
NEMO/trunk/src/OCE/DIA/diawri.F90
r12206 r12377 26 26 !!---------------------------------------------------------------------- 27 27 USE oce ! ocean dynamics and tracers 28 USE isf_oce 29 USE isfcpl 30 USE abl ! abl variables in case ln_abl = .true. 28 31 USE dom_oce ! ocean space and time domain 29 32 USE phycst ! physical constants … … 56 59 USE lib_mpp ! MPP library 57 60 USE timing ! preformance summary 58 USE diu rnal_bulk! diurnal warm layer59 USE cool_skin! Cool skin61 USE diu_bulk ! diurnal warm layer 62 USE diu_coolskin ! Cool skin 60 63 61 64 IMPLICIT NONE … … 65 68 PUBLIC dia_wri_state 66 69 PUBLIC dia_wri_alloc ! Called by nemogcm module 67 70 #if ! defined key_iomput 71 PUBLIC dia_wri_alloc_abl ! Called by sbcabl module (if ln_abl = .true.) 72 #endif 68 73 INTEGER :: nid_T, nz_T, nh_T, ndim_T, ndim_hT ! grid_T file 69 74 INTEGER :: nb_T , ndim_bT ! grid_T file … … 71 76 INTEGER :: nid_V, nz_V, nh_V, ndim_V, ndim_hV ! grid_V file 72 77 INTEGER :: nid_W, nz_W, nh_W ! grid_W file 78 INTEGER :: nid_A, nz_A, nh_A, ndim_A, ndim_hA ! grid_ABL file 73 79 INTEGER :: ndex(1) ! ??? 74 80 INTEGER, SAVE, ALLOCATABLE, DIMENSION(:) :: ndex_hT, ndex_hU, ndex_hV 81 INTEGER, SAVE, ALLOCATABLE, DIMENSION(:) :: ndex_hA, ndex_A ! ABL 75 82 INTEGER, SAVE, ALLOCATABLE, DIMENSION(:) :: ndex_T, ndex_U, ndex_V 76 83 INTEGER, SAVE, ALLOCATABLE, DIMENSION(:) :: ndex_bT 77 84 78 85 !! * Substitutions 79 # include " vectopt_loop_substitute.h90"86 # include "do_loop_substitute.h90" 80 87 !!---------------------------------------------------------------------- 81 88 !! NEMO/OCE 4.0 , NEMO Consortium (2018) … … 96 103 97 104 98 SUBROUTINE dia_wri( kt )105 SUBROUTINE dia_wri( kt, Kmm ) 99 106 !!--------------------------------------------------------------------- 100 107 !! *** ROUTINE dia_wri *** … … 106 113 !!---------------------------------------------------------------------- 107 114 INTEGER, INTENT( in ) :: kt ! ocean time-step index 115 INTEGER, INTENT( in ) :: Kmm ! ocean time level index 108 116 !! 109 117 INTEGER :: ji, jj, jk ! dummy loop indices … … 119 127 ! Output the initial state and forcings 120 128 IF( ninist == 1 ) THEN 121 CALL dia_wri_state( 'output.init' )129 CALL dia_wri_state( Kmm, 'output.init' ) 122 130 ninist = 0 123 131 ENDIF … … 128 136 CALL iom_put("e3v_0", e3v_0(:,:,:) ) 129 137 ! 130 CALL iom_put( "e3t" , e3t _n(:,:,:) )131 CALL iom_put( "e3u" , e3u _n(:,:,:) )132 CALL iom_put( "e3v" , e3v _n(:,:,:) )133 CALL iom_put( "e3w" , e3w _n(:,:,:) )138 CALL iom_put( "e3t" , e3t(:,:,:,Kmm) ) 139 CALL iom_put( "e3u" , e3u(:,:,:,Kmm) ) 140 CALL iom_put( "e3v" , e3v(:,:,:,Kmm) ) 141 CALL iom_put( "e3w" , e3w(:,:,:,Kmm) ) 134 142 IF( iom_use("e3tdef") ) & 135 CALL iom_put( "e3tdef" , ( ( e3t _n(:,:,:) - e3t_0(:,:,:) ) / e3t_0(:,:,:) * 100 * tmask(:,:,:) ) ** 2 )143 CALL iom_put( "e3tdef" , ( ( e3t(:,:,:,Kmm) - e3t_0(:,:,:) ) / e3t_0(:,:,:) * 100 * tmask(:,:,:) ) ** 2 ) 136 144 137 145 IF( ll_wd ) THEN 138 CALL iom_put( "ssh" , (ssh n+ssh_ref)*tmask(:,:,1) ) ! sea surface height (brought back to the reference used for wetting and drying)146 CALL iom_put( "ssh" , (ssh(:,:,Kmm)+ssh_ref)*tmask(:,:,1) ) ! sea surface height (brought back to the reference used for wetting and drying) 139 147 ELSE 140 CALL iom_put( "ssh" , ssh n) ! sea surface height148 CALL iom_put( "ssh" , ssh(:,:,Kmm) ) ! sea surface height 141 149 ENDIF 142 150 143 151 IF( iom_use("wetdep") ) & ! wet depth 144 CALL iom_put( "wetdep" , ht_0(:,:) + ssh n(:,:) )152 CALL iom_put( "wetdep" , ht_0(:,:) + ssh(:,:,Kmm) ) 145 153 146 CALL iom_put( "toce", ts n(:,:,:,jp_tem) ) ! 3D temperature147 CALL iom_put( "sst", ts n(:,:,1,jp_tem) ) ! surface temperature154 CALL iom_put( "toce", ts(:,:,:,jp_tem,Kmm) ) ! 3D temperature 155 CALL iom_put( "sst", ts(:,:,1,jp_tem,Kmm) ) ! surface temperature 148 156 IF ( iom_use("sbt") ) THEN 149 DO jj = 1, jpj 150 DO ji = 1, jpi 151 ikbot = mbkt(ji,jj) 152 z2d(ji,jj) = tsn(ji,jj,ikbot,jp_tem) 153 END DO 154 END DO 157 DO_2D_11_11 158 ikbot = mbkt(ji,jj) 159 z2d(ji,jj) = ts(ji,jj,ikbot,jp_tem,Kmm) 160 END_2D 155 161 CALL iom_put( "sbt", z2d ) ! bottom temperature 156 162 ENDIF 157 163 158 CALL iom_put( "soce", ts n(:,:,:,jp_sal) ) ! 3D salinity159 CALL iom_put( "sss", ts n(:,:,1,jp_sal) ) ! surface salinity164 CALL iom_put( "soce", ts(:,:,:,jp_sal,Kmm) ) ! 3D salinity 165 CALL iom_put( "sss", ts(:,:,1,jp_sal,Kmm) ) ! surface salinity 160 166 IF ( iom_use("sbs") ) THEN 161 DO jj = 1, jpj 162 DO ji = 1, jpi 163 ikbot = mbkt(ji,jj) 164 z2d(ji,jj) = tsn(ji,jj,ikbot,jp_sal) 165 END DO 166 END DO 167 DO_2D_11_11 168 ikbot = mbkt(ji,jj) 169 z2d(ji,jj) = ts(ji,jj,ikbot,jp_sal,Kmm) 170 END_2D 167 171 CALL iom_put( "sbs", z2d ) ! bottom salinity 168 172 ENDIF … … 171 175 zztmp = rau0 * 0.25 172 176 z2d(:,:) = 0._wp 173 DO jj = 2, jpjm1 174 DO ji = fs_2, fs_jpim1 ! vector opt. 175 zztmp2 = ( ( rCdU_bot(ji+1,jj)+rCdU_bot(ji ,jj) ) * un(ji ,jj,mbku(ji ,jj)) )**2 & 176 & + ( ( rCdU_bot(ji ,jj)+rCdU_bot(ji-1,jj) ) * un(ji-1,jj,mbku(ji-1,jj)) )**2 & 177 & + ( ( rCdU_bot(ji,jj+1)+rCdU_bot(ji,jj ) ) * vn(ji,jj ,mbkv(ji,jj )) )**2 & 178 & + ( ( rCdU_bot(ji,jj )+rCdU_bot(ji,jj-1) ) * vn(ji,jj-1,mbkv(ji,jj-1)) )**2 179 z2d(ji,jj) = zztmp * SQRT( zztmp2 ) * tmask(ji,jj,1) 180 ! 181 END DO 182 END DO 177 DO_2D_00_00 178 zztmp2 = ( ( rCdU_bot(ji+1,jj)+rCdU_bot(ji ,jj) ) * uu(ji ,jj,mbku(ji ,jj),Kmm) )**2 & 179 & + ( ( rCdU_bot(ji ,jj)+rCdU_bot(ji-1,jj) ) * uu(ji-1,jj,mbku(ji-1,jj),Kmm) )**2 & 180 & + ( ( rCdU_bot(ji,jj+1)+rCdU_bot(ji,jj ) ) * vv(ji,jj ,mbkv(ji,jj ),Kmm) )**2 & 181 & + ( ( rCdU_bot(ji,jj )+rCdU_bot(ji,jj-1) ) * vv(ji,jj-1,mbkv(ji,jj-1),Kmm) )**2 182 z2d(ji,jj) = zztmp * SQRT( zztmp2 ) * tmask(ji,jj,1) 183 ! 184 END_2D 183 185 CALL lbc_lnk( 'diawri', z2d, 'T', 1. ) 184 186 CALL iom_put( "taubot", z2d ) 185 187 ENDIF 186 188 187 CALL iom_put( "uoce", u n(:,:,:) ) ! 3D i-current188 CALL iom_put( "ssu", u n(:,:,1) ) ! surface i-current189 CALL iom_put( "uoce", uu(:,:,:,Kmm) ) ! 3D i-current 190 CALL iom_put( "ssu", uu(:,:,1,Kmm) ) ! surface i-current 189 191 IF ( iom_use("sbu") ) THEN 190 DO jj = 1, jpj 191 DO ji = 1, jpi 192 ikbot = mbku(ji,jj) 193 z2d(ji,jj) = un(ji,jj,ikbot) 194 END DO 195 END DO 192 DO_2D_11_11 193 ikbot = mbku(ji,jj) 194 z2d(ji,jj) = uu(ji,jj,ikbot,Kmm) 195 END_2D 196 196 CALL iom_put( "sbu", z2d ) ! bottom i-current 197 197 ENDIF 198 198 199 CALL iom_put( "voce", v n(:,:,:) ) ! 3D j-current200 CALL iom_put( "ssv", v n(:,:,1) ) ! surface j-current199 CALL iom_put( "voce", vv(:,:,:,Kmm) ) ! 3D j-current 200 CALL iom_put( "ssv", vv(:,:,1,Kmm) ) ! surface j-current 201 201 IF ( iom_use("sbv") ) THEN 202 DO jj = 1, jpj 203 DO ji = 1, jpi 204 ikbot = mbkv(ji,jj) 205 z2d(ji,jj) = vn(ji,jj,ikbot) 206 END DO 207 END DO 202 DO_2D_11_11 203 ikbot = mbkv(ji,jj) 204 z2d(ji,jj) = vv(ji,jj,ikbot,Kmm) 205 END_2D 208 206 CALL iom_put( "sbv", z2d ) ! bottom j-current 209 207 ENDIF 210 208 211 IF( ln_zad_Aimp ) w n = wn+ wi ! Recombine explicit and implicit parts of vertical velocity for diagnostic output212 ! 213 CALL iom_put( "woce", w n) ! vertical velocity209 IF( ln_zad_Aimp ) ww = ww + wi ! Recombine explicit and implicit parts of vertical velocity for diagnostic output 210 ! 211 CALL iom_put( "woce", ww ) ! vertical velocity 214 212 IF( iom_use('w_masstr') .OR. iom_use('w_masstr2') ) THEN ! vertical mass transport & its square value 215 213 ! Caution: in the VVL case, it only correponds to the baroclinic mass transport. 216 214 z2d(:,:) = rau0 * e1e2t(:,:) 217 215 DO jk = 1, jpk 218 z3d(:,:,jk) = w n(:,:,jk) * z2d(:,:)216 z3d(:,:,jk) = ww(:,:,jk) * z2d(:,:) 219 217 END DO 220 218 CALL iom_put( "w_masstr" , z3d ) … … 222 220 ENDIF 223 221 ! 224 IF( ln_zad_Aimp ) w n = wn- wi ! Remove implicit part of vertical velocity that was added for diagnostic output222 IF( ln_zad_Aimp ) ww = ww - wi ! Remove implicit part of vertical velocity that was added for diagnostic output 225 223 226 224 CALL iom_put( "avt" , avt ) ! T vert. eddy diff. coef. … … 232 230 233 231 IF ( iom_use("sstgrad") .OR. iom_use("sstgrad2") ) THEN 234 DO jj = 2, jpjm1 ! sst gradient 235 DO ji = fs_2, fs_jpim1 ! vector opt. 236 zztmp = tsn(ji,jj,1,jp_tem) 237 zztmpx = ( tsn(ji+1,jj,1,jp_tem) - zztmp ) * r1_e1u(ji,jj) + ( zztmp - tsn(ji-1,jj ,1,jp_tem) ) * r1_e1u(ji-1,jj) 238 zztmpy = ( tsn(ji,jj+1,1,jp_tem) - zztmp ) * r1_e2v(ji,jj) + ( zztmp - tsn(ji ,jj-1,1,jp_tem) ) * r1_e2v(ji,jj-1) 239 z2d(ji,jj) = 0.25 * ( zztmpx * zztmpx + zztmpy * zztmpy ) & 240 & * umask(ji,jj,1) * umask(ji-1,jj,1) * vmask(ji,jj,1) * umask(ji,jj-1,1) 241 END DO 242 END DO 232 DO_2D_00_00 233 zztmp = ts(ji,jj,1,jp_tem,Kmm) 234 zztmpx = ( ts(ji+1,jj,1,jp_tem,Kmm) - zztmp ) * r1_e1u(ji,jj) + ( zztmp - ts(ji-1,jj ,1,jp_tem,Kmm) ) * r1_e1u(ji-1,jj) 235 zztmpy = ( ts(ji,jj+1,1,jp_tem,Kmm) - zztmp ) * r1_e2v(ji,jj) + ( zztmp - ts(ji ,jj-1,1,jp_tem,Kmm) ) * r1_e2v(ji,jj-1) 236 z2d(ji,jj) = 0.25 * ( zztmpx * zztmpx + zztmpy * zztmpy ) & 237 & * umask(ji,jj,1) * umask(ji-1,jj,1) * vmask(ji,jj,1) * umask(ji,jj-1,1) 238 END_2D 243 239 CALL lbc_lnk( 'diawri', z2d, 'T', 1. ) 244 240 CALL iom_put( "sstgrad2", z2d ) ! square of module of sst gradient … … 250 246 IF( iom_use("heatc") ) THEN 251 247 z2d(:,:) = 0._wp 252 DO jk = 1, jpkm1 253 DO jj = 1, jpj 254 DO ji = 1, jpi 255 z2d(ji,jj) = z2d(ji,jj) + e3t_n(ji,jj,jk) * tsn(ji,jj,jk,jp_tem) * tmask(ji,jj,jk) 256 END DO 257 END DO 258 END DO 248 DO_3D_11_11( 1, jpkm1 ) 249 z2d(ji,jj) = z2d(ji,jj) + e3t(ji,jj,jk,Kmm) * ts(ji,jj,jk,jp_tem,Kmm) * tmask(ji,jj,jk) 250 END_3D 259 251 CALL iom_put( "heatc", rau0_rcp * z2d ) ! vertically integrated heat content (J/m2) 260 252 ENDIF … … 262 254 IF( iom_use("saltc") ) THEN 263 255 z2d(:,:) = 0._wp 264 DO jk = 1, jpkm1 265 DO jj = 1, jpj 266 DO ji = 1, jpi 267 z2d(ji,jj) = z2d(ji,jj) + e3t_n(ji,jj,jk) * tsn(ji,jj,jk,jp_sal) * tmask(ji,jj,jk) 268 END DO 269 END DO 270 END DO 256 DO_3D_11_11( 1, jpkm1 ) 257 z2d(ji,jj) = z2d(ji,jj) + e3t(ji,jj,jk,Kmm) * ts(ji,jj,jk,jp_sal,Kmm) * tmask(ji,jj,jk) 258 END_3D 271 259 CALL iom_put( "saltc", rau0 * z2d ) ! vertically integrated salt content (PSU*kg/m2) 272 260 ENDIF … … 274 262 IF ( iom_use("eken") ) THEN 275 263 z3d(:,:,jpk) = 0._wp 276 DO jk = 1, jpkm1 277 DO jj = 2, jpjm1 278 DO ji = fs_2, fs_jpim1 ! vector opt. 279 zztmp = 0.25_wp * r1_e1e2t(ji,jj) / e3t_n(ji,jj,jk) 280 z3d(ji,jj,jk) = zztmp * ( un(ji-1,jj,jk)**2 * e2u(ji-1,jj) * e3u_n(ji-1,jj,jk) & 281 & + un(ji ,jj,jk)**2 * e2u(ji ,jj) * e3u_n(ji ,jj,jk) & 282 & + vn(ji,jj-1,jk)**2 * e1v(ji,jj-1) * e3v_n(ji,jj-1,jk) & 283 & + vn(ji,jj ,jk)**2 * e1v(ji,jj ) * e3v_n(ji,jj ,jk) ) 284 END DO 285 END DO 286 END DO 264 DO_3D_00_00( 1, jpkm1 ) 265 zztmp = 0.25_wp * r1_e1e2t(ji,jj) / e3t(ji,jj,jk,Kmm) 266 z3d(ji,jj,jk) = zztmp * ( uu(ji-1,jj,jk,Kmm)**2 * e2u(ji-1,jj) * e3u(ji-1,jj,jk,Kmm) & 267 & + uu(ji ,jj,jk,Kmm)**2 * e2u(ji ,jj) * e3u(ji ,jj,jk,Kmm) & 268 & + vv(ji,jj-1,jk,Kmm)**2 * e1v(ji,jj-1) * e3v(ji,jj-1,jk,Kmm) & 269 & + vv(ji,jj ,jk,Kmm)**2 * e1v(ji,jj ) * e3v(ji,jj ,jk,Kmm) ) 270 END_3D 287 271 CALL lbc_lnk( 'diawri', z3d, 'T', 1. ) 288 272 CALL iom_put( "eken", z3d ) ! kinetic energy 289 273 ENDIF 290 274 ! 291 CALL iom_put( "hdiv", hdiv n) ! Horizontal divergence275 CALL iom_put( "hdiv", hdiv ) ! Horizontal divergence 292 276 ! 293 277 IF( iom_use("u_masstr") .OR. iom_use("u_masstr_vint") .OR. iom_use("u_heattr") .OR. iom_use("u_salttr") ) THEN … … 295 279 z2d(:,:) = 0.e0 296 280 DO jk = 1, jpkm1 297 z3d(:,:,jk) = rau0 * u n(:,:,jk) * e2u(:,:) * e3u_n(:,:,jk) * umask(:,:,jk)281 z3d(:,:,jk) = rau0 * uu(:,:,jk,Kmm) * e2u(:,:) * e3u(:,:,jk,Kmm) * umask(:,:,jk) 298 282 z2d(:,:) = z2d(:,:) + z3d(:,:,jk) 299 283 END DO … … 304 288 IF( iom_use("u_heattr") ) THEN 305 289 z2d(:,:) = 0._wp 306 DO jk = 1, jpkm1 307 DO jj = 2, jpjm1 308 DO ji = fs_2, fs_jpim1 ! vector opt. 309 z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( tsn(ji,jj,jk,jp_tem) + tsn(ji+1,jj,jk,jp_tem) ) 310 END DO 311 END DO 312 END DO 290 DO_3D_00_00( 1, jpkm1 ) 291 z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( ts(ji,jj,jk,jp_tem,Kmm) + ts(ji+1,jj,jk,jp_tem,Kmm) ) 292 END_3D 313 293 CALL lbc_lnk( 'diawri', z2d, 'U', -1. ) 314 294 CALL iom_put( "u_heattr", 0.5*rcp * z2d ) ! heat transport in i-direction … … 317 297 IF( iom_use("u_salttr") ) THEN 318 298 z2d(:,:) = 0.e0 319 DO jk = 1, jpkm1 320 DO jj = 2, jpjm1 321 DO ji = fs_2, fs_jpim1 ! vector opt. 322 z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( tsn(ji,jj,jk,jp_sal) + tsn(ji+1,jj,jk,jp_sal) ) 323 END DO 324 END DO 325 END DO 299 DO_3D_00_00( 1, jpkm1 ) 300 z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( ts(ji,jj,jk,jp_sal,Kmm) + ts(ji+1,jj,jk,jp_sal,Kmm) ) 301 END_3D 326 302 CALL lbc_lnk( 'diawri', z2d, 'U', -1. ) 327 303 CALL iom_put( "u_salttr", 0.5 * z2d ) ! heat transport in i-direction … … 332 308 z3d(:,:,jpk) = 0.e0 333 309 DO jk = 1, jpkm1 334 z3d(:,:,jk) = rau0 * v n(:,:,jk) * e1v(:,:) * e3v_n(:,:,jk) * vmask(:,:,jk)310 z3d(:,:,jk) = rau0 * vv(:,:,jk,Kmm) * e1v(:,:) * e3v(:,:,jk,Kmm) * vmask(:,:,jk) 335 311 END DO 336 312 CALL iom_put( "v_masstr", z3d ) ! mass transport in j-direction … … 339 315 IF( iom_use("v_heattr") ) THEN 340 316 z2d(:,:) = 0.e0 341 DO jk = 1, jpkm1 342 DO jj = 2, jpjm1 343 DO ji = fs_2, fs_jpim1 ! vector opt. 344 z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( tsn(ji,jj,jk,jp_tem) + tsn(ji,jj+1,jk,jp_tem) ) 345 END DO 346 END DO 347 END DO 317 DO_3D_00_00( 1, jpkm1 ) 318 z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( ts(ji,jj,jk,jp_tem,Kmm) + ts(ji,jj+1,jk,jp_tem,Kmm) ) 319 END_3D 348 320 CALL lbc_lnk( 'diawri', z2d, 'V', -1. ) 349 321 CALL iom_put( "v_heattr", 0.5*rcp * z2d ) ! heat transport in j-direction … … 352 324 IF( iom_use("v_salttr") ) THEN 353 325 z2d(:,:) = 0._wp 354 DO jk = 1, jpkm1 355 DO jj = 2, jpjm1 356 DO ji = fs_2, fs_jpim1 ! vector opt. 357 z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( tsn(ji,jj,jk,jp_sal) + tsn(ji,jj+1,jk,jp_sal) ) 358 END DO 359 END DO 360 END DO 326 DO_3D_00_00( 1, jpkm1 ) 327 z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( ts(ji,jj,jk,jp_sal,Kmm) + ts(ji,jj+1,jk,jp_sal,Kmm) ) 328 END_3D 361 329 CALL lbc_lnk( 'diawri', z2d, 'V', -1. ) 362 330 CALL iom_put( "v_salttr", 0.5 * z2d ) ! heat transport in j-direction … … 365 333 IF( iom_use("tosmint") ) THEN 366 334 z2d(:,:) = 0._wp 367 DO jk = 1, jpkm1 368 DO jj = 2, jpjm1 369 DO ji = fs_2, fs_jpim1 ! vector opt. 370 z2d(ji,jj) = z2d(ji,jj) + e3t_n(ji,jj,jk) * tsn(ji,jj,jk,jp_tem) 371 END DO 372 END DO 373 END DO 335 DO_3D_00_00( 1, jpkm1 ) 336 z2d(ji,jj) = z2d(ji,jj) + e3t(ji,jj,jk,Kmm) * ts(ji,jj,jk,jp_tem,Kmm) 337 END_3D 374 338 CALL lbc_lnk( 'diawri', z2d, 'T', -1. ) 375 339 CALL iom_put( "tosmint", rau0 * z2d ) ! Vertical integral of temperature … … 377 341 IF( iom_use("somint") ) THEN 378 342 z2d(:,:)=0._wp 379 DO jk = 1, jpkm1 380 DO jj = 2, jpjm1 381 DO ji = fs_2, fs_jpim1 ! vector opt. 382 z2d(ji,jj) = z2d(ji,jj) + e3t_n(ji,jj,jk) * tsn(ji,jj,jk,jp_sal) 383 END DO 384 END DO 385 END DO 343 DO_3D_00_00( 1, jpkm1 ) 344 z2d(ji,jj) = z2d(ji,jj) + e3t(ji,jj,jk,Kmm) * ts(ji,jj,jk,jp_sal,Kmm) 345 END_3D 386 346 CALL lbc_lnk( 'diawri', z2d, 'T', -1. ) 387 347 CALL iom_put( "somint", rau0 * z2d ) ! Vertical integral of salinity … … 390 350 CALL iom_put( "bn2", rn2 ) ! Brunt-Vaisala buoyancy frequency (N^2) 391 351 ! 392 393 IF (ln_dia25h) CALL dia_25h( kt )! 25h averaging352 353 IF (ln_dia25h) CALL dia_25h( kt, Kmm ) ! 25h averaging 394 354 395 355 IF( ln_timing ) CALL timing_stop('dia_wri') … … 411 371 & ndex_hV(jpi*jpj) , ndex_V(jpi*jpj*jpk) , STAT=ierr(1) ) 412 372 ! 413 373 dia_wri_alloc = MAXVAL(ierr) 414 374 CALL mpp_sum( 'diawri', dia_wri_alloc ) 415 375 ! 416 376 END FUNCTION dia_wri_alloc 377 378 INTEGER FUNCTION dia_wri_alloc_abl() 379 !!---------------------------------------------------------------------- 380 ALLOCATE( ndex_hA(jpi*jpj), ndex_A (jpi*jpj*jpkam1), STAT=dia_wri_alloc_abl) 381 CALL mpp_sum( 'diawri', dia_wri_alloc_abl ) 382 ! 383 END FUNCTION dia_wri_alloc_abl 417 384 418 385 419 SUBROUTINE dia_wri( kt )386 SUBROUTINE dia_wri( kt, Kmm ) 420 387 !!--------------------------------------------------------------------- 421 388 !! *** ROUTINE dia_wri *** … … 430 397 !!---------------------------------------------------------------------- 431 398 INTEGER, INTENT( in ) :: kt ! ocean time-step index 399 INTEGER, INTENT( in ) :: Kmm ! ocean time level index 432 400 ! 433 401 LOGICAL :: ll_print = .FALSE. ! =T print and flush numout … … 437 405 INTEGER :: ierr ! error code return from allocation 438 406 INTEGER :: iimi, iima, ipk, it, itmod, ijmi, ijma ! local integers 407 INTEGER :: ipka ! ABL 439 408 INTEGER :: jn, ierror ! local integers 440 409 REAL(wp) :: zsto, zout, zmax, zjulian ! local scalars … … 442 411 REAL(wp), DIMENSION(jpi,jpj) :: zw2d ! 2D workspace 443 412 REAL(wp), DIMENSION(jpi,jpj,jpk) :: zw3d ! 3D workspace 413 REAL(wp), DIMENSION(:,:,:), ALLOCATABLE :: zw3d_abl ! ABL 3D workspace 444 414 !!---------------------------------------------------------------------- 445 415 ! 446 416 IF( ninist == 1 ) THEN !== Output the initial state and forcings ==! 447 CALL dia_wri_state( 'output.init' )417 CALL dia_wri_state( Kmm, 'output.init' ) 448 418 ninist = 0 449 419 ENDIF … … 475 445 ijmi = 1 ; ijma = jpj 476 446 ipk = jpk 447 IF(ln_abl) ipka = jpkam1 477 448 478 449 ! define time axis … … 577 548 & "m", ipk, gdepw_1d, nz_W, "down" ) 578 549 550 IF( ln_abl ) THEN 551 ! Define the ABL grid FILE ( nid_A ) 552 CALL dia_nam( clhstnam, nn_write, 'grid_ABL' ) 553 IF(lwp) WRITE(numout,*) " Name of NETCDF file ", clhstnam ! filename 554 CALL histbeg( clhstnam, jpi, glamt, jpj, gphit, & ! Horizontal grid: glamt and gphit 555 & iimi, iima-iimi+1, ijmi, ijma-ijmi+1, & 556 & nit000-1, zjulian, rdt, nh_A, nid_A, domain_id=nidom, snc4chunks=snc4set ) 557 CALL histvert( nid_A, "ght_abl", "Vertical T levels", & ! Vertical grid: gdept 558 & "m", ipka, ght_abl(2:jpka), nz_A, "up" ) 559 ! ! Index of ocean points 560 ALLOCATE( zw3d_abl(jpi,jpj,ipka) ) 561 zw3d_abl(:,:,:) = 1._wp 562 CALL wheneq( jpi*jpj*ipka, zw3d_abl, 1, 1., ndex_A , ndim_A ) ! volume 563 CALL wheneq( jpi*jpj , zw3d_abl, 1, 1., ndex_hA, ndim_hA ) ! surface 564 DEALLOCATE(zw3d_abl) 565 ENDIF 579 566 580 567 ! Declare all the output fields as NETCDF variables … … 586 573 & jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout ) 587 574 IF( .NOT.ln_linssh ) THEN 588 CALL histdef( nid_T, "vovvle3t", "Level thickness" , "m" ,& ! e3t _n575 CALL histdef( nid_T, "vovvle3t", "Level thickness" , "m" ,& ! e3t(:,:,:,Kmm) 589 576 & jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout ) 590 CALL histdef( nid_T, "vovvldep", "T point depth" , "m" ,& ! e3t _n577 CALL histdef( nid_T, "vovvldep", "T point depth" , "m" ,& ! e3t(:,:,:,Kmm) 591 578 & jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout ) 592 CALL histdef( nid_T, "vovvldef", "Squared level deformation" , "%^2" ,& ! e3t _n579 CALL histdef( nid_T, "vovvldef", "Squared level deformation" , "%^2" ,& ! e3t(:,:,:,Kmm) 593 580 & jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout ) 594 581 ENDIF … … 607 594 & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) 608 595 IF( ln_linssh ) THEN 609 CALL histdef( nid_T, "sosst_cd", "Concentration/Dilution term on temperature" & ! emp * ts n(:,:,1,jp_tem)596 CALL histdef( nid_T, "sosst_cd", "Concentration/Dilution term on temperature" & ! emp * ts(:,:,1,jp_tem,Kmm) 610 597 & , "KgC/m2/s", & ! sosst_cd 611 598 & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) 612 CALL histdef( nid_T, "sosss_cd", "Concentration/Dilution term on salinity" & ! emp * ts n(:,:,1,jp_sal)599 CALL histdef( nid_T, "sosss_cd", "Concentration/Dilution term on salinity" & ! emp * ts(:,:,1,jp_sal,Kmm) 613 600 & , "KgPSU/m2/s",& ! sosss_cd 614 601 & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) … … 626 613 CALL histdef( nid_T, "sowindsp", "wind speed at 10m" , "m/s" , & ! wndm 627 614 & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) 628 ! 615 ! 616 IF( ln_abl ) THEN 617 CALL histdef( nid_A, "t_abl", "Potential Temperature" , "K" , & ! t_abl 618 & jpi, jpj, nh_A, ipka, 1, ipka, nz_A, 32, clop, zsto, zout ) 619 CALL histdef( nid_A, "q_abl", "Humidity" , "kg/kg" , & ! q_abl 620 & jpi, jpj, nh_A, ipka, 1, ipka, nz_A, 32, clop, zsto, zout ) 621 CALL histdef( nid_A, "u_abl", "Atmospheric U-wind " , "m/s" , & ! u_abl 622 & jpi, jpj, nh_A, ipka, 1, ipka, nz_A, 32, clop, zsto, zout ) 623 CALL histdef( nid_A, "v_abl", "Atmospheric V-wind " , "m/s" , & ! v_abl 624 & jpi, jpj, nh_A, ipka, 1, ipka, nz_A, 32, clop, zsto, zout ) 625 CALL histdef( nid_A, "tke_abl", "Atmospheric TKE " , "m2/s2" , & ! tke_abl 626 & jpi, jpj, nh_A, ipka, 1, ipka, nz_A, 32, clop, zsto, zout ) 627 CALL histdef( nid_A, "avm_abl", "Atmospheric turbulent viscosity", "m2/s" , & ! avm_abl 628 & jpi, jpj, nh_A, ipka, 1, ipka, nz_A, 32, clop, zsto, zout ) 629 CALL histdef( nid_A, "avt_abl", "Atmospheric turbulent diffusivity", "m2/s2", & ! avt_abl 630 & jpi, jpj, nh_A, ipka, 1, ipka, nz_A, 32, clop, zsto, zout ) 631 CALL histdef( nid_A, "pblh", "Atmospheric boundary layer height " , "m", & ! pblh 632 & jpi, jpj, nh_A, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) 633 #if defined key_si3 634 CALL histdef( nid_A, "oce_frac", "Fraction of open ocean" , " ", & ! ato_i 635 & jpi, jpj, nh_A, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) 636 #endif 637 CALL histend( nid_A, snc4chunks=snc4set ) 638 ENDIF 639 ! 629 640 IF( ln_icebergs ) THEN 630 641 CALL histdef( nid_T, "calving" , "calving mass input" , "kg/s" , & … … 686 697 687 698 ! !!! nid_U : 3D 688 CALL histdef( nid_U, "vozocrtx", "Zonal Current" , "m/s" , & ! u n699 CALL histdef( nid_U, "vozocrtx", "Zonal Current" , "m/s" , & ! uu(:,:,:,Kmm) 689 700 & jpi, jpj, nh_U, ipk, 1, ipk, nz_U, 32, clop, zsto, zout ) 690 701 IF( ln_wave .AND. ln_sdw) THEN … … 699 710 700 711 ! !!! nid_V : 3D 701 CALL histdef( nid_V, "vomecrty", "Meridional Current" , "m/s" , & ! v n712 CALL histdef( nid_V, "vomecrty", "Meridional Current" , "m/s" , & ! vv(:,:,:,Kmm) 702 713 & jpi, jpj, nh_V, ipk, 1, ipk, nz_V, 32, clop, zsto, zout ) 703 714 IF( ln_wave .AND. ln_sdw) THEN … … 712 723 713 724 ! !!! nid_W : 3D 714 CALL histdef( nid_W, "vovecrtz", "Vertical Velocity" , "m/s" , & ! w n725 CALL histdef( nid_W, "vovecrtz", "Vertical Velocity" , "m/s" , & ! ww 715 726 & jpi, jpj, nh_W, ipk, 1, ipk, nz_W, 32, clop, zsto, zout ) 716 727 CALL histdef( nid_W, "votkeavt", "Vertical Eddy Diffusivity" , "m2/s" , & ! avt … … 750 761 751 762 IF( .NOT.ln_linssh ) THEN 752 CALL histwrite( nid_T, "votemper", it, ts n(:,:,:,jp_tem) * e3t_n(:,:,:) , ndim_T , ndex_T ) ! heat content753 CALL histwrite( nid_T, "vosaline", it, ts n(:,:,:,jp_sal) * e3t_n(:,:,:) , ndim_T , ndex_T ) ! salt content754 CALL histwrite( nid_T, "sosstsst", it, ts n(:,:,1,jp_tem) * e3t_n(:,:,1) , ndim_hT, ndex_hT ) ! sea surface heat content755 CALL histwrite( nid_T, "sosaline", it, ts n(:,:,1,jp_sal) * e3t_n(:,:,1) , ndim_hT, ndex_hT ) ! sea surface salinity content763 CALL histwrite( nid_T, "votemper", it, ts(:,:,:,jp_tem,Kmm) * e3t(:,:,:,Kmm) , ndim_T , ndex_T ) ! heat content 764 CALL histwrite( nid_T, "vosaline", it, ts(:,:,:,jp_sal,Kmm) * e3t(:,:,:,Kmm) , ndim_T , ndex_T ) ! salt content 765 CALL histwrite( nid_T, "sosstsst", it, ts(:,:,1,jp_tem,Kmm) * e3t(:,:,1,Kmm) , ndim_hT, ndex_hT ) ! sea surface heat content 766 CALL histwrite( nid_T, "sosaline", it, ts(:,:,1,jp_sal,Kmm) * e3t(:,:,1,Kmm) , ndim_hT, ndex_hT ) ! sea surface salinity content 756 767 ELSE 757 CALL histwrite( nid_T, "votemper", it, ts n(:,:,:,jp_tem) , ndim_T , ndex_T ) ! temperature758 CALL histwrite( nid_T, "vosaline", it, ts n(:,:,:,jp_sal) , ndim_T , ndex_T ) ! salinity759 CALL histwrite( nid_T, "sosstsst", it, ts n(:,:,1,jp_tem) , ndim_hT, ndex_hT ) ! sea surface temperature760 CALL histwrite( nid_T, "sosaline", it, ts n(:,:,1,jp_sal) , ndim_hT, ndex_hT ) ! sea surface salinity768 CALL histwrite( nid_T, "votemper", it, ts(:,:,:,jp_tem,Kmm) , ndim_T , ndex_T ) ! temperature 769 CALL histwrite( nid_T, "vosaline", it, ts(:,:,:,jp_sal,Kmm) , ndim_T , ndex_T ) ! salinity 770 CALL histwrite( nid_T, "sosstsst", it, ts(:,:,1,jp_tem,Kmm) , ndim_hT, ndex_hT ) ! sea surface temperature 771 CALL histwrite( nid_T, "sosaline", it, ts(:,:,1,jp_sal,Kmm) , ndim_hT, ndex_hT ) ! sea surface salinity 761 772 ENDIF 762 773 IF( .NOT.ln_linssh ) THEN 763 zw3d(:,:,:) = ( ( e3t _n(:,:,:) - e3t_0(:,:,:) ) / e3t_0(:,:,:) * 100 * tmask(:,:,:) ) ** 2764 CALL histwrite( nid_T, "vovvle3t", it, e3t _n (:,:,:) , ndim_T , ndex_T ) ! level thickness765 CALL histwrite( nid_T, "vovvldep", it, gdept _n(:,:,:) , ndim_T , ndex_T ) ! t-point depth774 zw3d(:,:,:) = ( ( e3t(:,:,:,Kmm) - e3t_0(:,:,:) ) / e3t_0(:,:,:) * 100 * tmask(:,:,:) ) ** 2 775 CALL histwrite( nid_T, "vovvle3t", it, e3t (:,:,:,Kmm) , ndim_T , ndex_T ) ! level thickness 776 CALL histwrite( nid_T, "vovvldep", it, gdept(:,:,:,Kmm) , ndim_T , ndex_T ) ! t-point depth 766 777 CALL histwrite( nid_T, "vovvldef", it, zw3d , ndim_T , ndex_T ) ! level thickness deformation 767 778 ENDIF 768 CALL histwrite( nid_T, "sossheig", it, ssh n, ndim_hT, ndex_hT ) ! sea surface height779 CALL histwrite( nid_T, "sossheig", it, ssh(:,:,Kmm) , ndim_hT, ndex_hT ) ! sea surface height 769 780 CALL histwrite( nid_T, "sowaflup", it, ( emp-rnf ) , ndim_hT, ndex_hT ) ! upward water flux 770 781 CALL histwrite( nid_T, "sorunoff", it, rnf , ndim_hT, ndex_hT ) ! river runoffs … … 773 784 ! in linear free surface case) 774 785 IF( ln_linssh ) THEN 775 zw2d(:,:) = emp (:,:) * ts n(:,:,1,jp_tem)786 zw2d(:,:) = emp (:,:) * ts(:,:,1,jp_tem,Kmm) 776 787 CALL histwrite( nid_T, "sosst_cd", it, zw2d, ndim_hT, ndex_hT ) ! c/d term on sst 777 zw2d(:,:) = emp (:,:) * ts n(:,:,1,jp_sal)788 zw2d(:,:) = emp (:,:) * ts(:,:,1,jp_sal,Kmm) 778 789 CALL histwrite( nid_T, "sosss_cd", it, zw2d, ndim_hT, ndex_hT ) ! c/d term on sss 779 790 ENDIF … … 784 795 CALL histwrite( nid_T, "soicecov", it, fr_i , ndim_hT, ndex_hT ) ! ice fraction 785 796 CALL histwrite( nid_T, "sowindsp", it, wndm , ndim_hT, ndex_hT ) ! wind speed 786 ! 797 ! 798 IF( ln_abl ) THEN 799 ALLOCATE( zw3d_abl(jpi,jpj,jpka) ) 800 IF( ln_mskland ) THEN 801 DO jk=1,jpka 802 zw3d_abl(:,:,jk) = tmask(:,:,1) 803 END DO 804 ELSE 805 zw3d_abl(:,:,:) = 1._wp 806 ENDIF 807 CALL histwrite( nid_A, "pblh" , it, pblh(:,:) *zw3d_abl(:,:,1 ), ndim_hA, ndex_hA ) ! pblh 808 CALL histwrite( nid_A, "u_abl" , it, u_abl (:,:,2:jpka,nt_n )*zw3d_abl(:,:,2:jpka), ndim_A , ndex_A ) ! u_abl 809 CALL histwrite( nid_A, "v_abl" , it, v_abl (:,:,2:jpka,nt_n )*zw3d_abl(:,:,2:jpka), ndim_A , ndex_A ) ! v_abl 810 CALL histwrite( nid_A, "t_abl" , it, tq_abl (:,:,2:jpka,nt_n,1)*zw3d_abl(:,:,2:jpka), ndim_A , ndex_A ) ! t_abl 811 CALL histwrite( nid_A, "q_abl" , it, tq_abl (:,:,2:jpka,nt_n,2)*zw3d_abl(:,:,2:jpka), ndim_A , ndex_A ) ! q_abl 812 CALL histwrite( nid_A, "tke_abl", it, tke_abl (:,:,2:jpka,nt_n )*zw3d_abl(:,:,2:jpka), ndim_A , ndex_A ) ! tke_abl 813 CALL histwrite( nid_A, "avm_abl", it, avm_abl (:,:,2:jpka )*zw3d_abl(:,:,2:jpka), ndim_A , ndex_A ) ! avm_abl 814 CALL histwrite( nid_A, "avt_abl", it, avt_abl (:,:,2:jpka )*zw3d_abl(:,:,2:jpka), ndim_A , ndex_A ) ! avt_abl 815 #if defined key_si3 816 CALL histwrite( nid_A, "oce_frac" , it, ato_i(:,:) , ndim_hA, ndex_hA ) ! ato_i 817 #endif 818 DEALLOCATE(zw3d_abl) 819 ENDIF 820 ! 787 821 IF( ln_icebergs ) THEN 788 822 ! … … 811 845 CALL histwrite( nid_T, "sohefldp", it, qrp , ndim_hT, ndex_hT ) ! heat flux damping 812 846 CALL histwrite( nid_T, "sowafldp", it, erp , ndim_hT, ndex_hT ) ! freshwater flux damping 813 zw2d(:,:) = erp(:,:) * ts n(:,:,1,jp_sal) * tmask(:,:,1)847 zw2d(:,:) = erp(:,:) * ts(:,:,1,jp_sal,Kmm) * tmask(:,:,1) 814 848 CALL histwrite( nid_T, "sosafldp", it, zw2d , ndim_hT, ndex_hT ) ! salt flux damping 815 849 ENDIF … … 824 858 #endif 825 859 826 CALL histwrite( nid_U, "vozocrtx", it, u n, ndim_U , ndex_U ) ! i-current860 CALL histwrite( nid_U, "vozocrtx", it, uu(:,:,:,Kmm) , ndim_U , ndex_U ) ! i-current 827 861 CALL histwrite( nid_U, "sozotaux", it, utau , ndim_hU, ndex_hU ) ! i-wind stress 828 862 829 CALL histwrite( nid_V, "vomecrty", it, v n, ndim_V , ndex_V ) ! j-current863 CALL histwrite( nid_V, "vomecrty", it, vv(:,:,:,Kmm) , ndim_V , ndex_V ) ! j-current 830 864 CALL histwrite( nid_V, "sometauy", it, vtau , ndim_hV, ndex_hV ) ! j-wind stress 831 865 832 866 IF( ln_zad_Aimp ) THEN 833 CALL histwrite( nid_W, "vovecrtz", it, w n+ wi , ndim_T, ndex_T ) ! vert. current867 CALL histwrite( nid_W, "vovecrtz", it, ww + wi , ndim_T, ndex_T ) ! vert. current 834 868 ELSE 835 CALL histwrite( nid_W, "vovecrtz", it, w n, ndim_T, ndex_T ) ! vert. current869 CALL histwrite( nid_W, "vovecrtz", it, ww , ndim_T, ndex_T ) ! vert. current 836 870 ENDIF 837 871 CALL histwrite( nid_W, "votkeavt", it, avt , ndim_T, ndex_T ) ! T vert. eddy diff. coef. … … 854 888 CALL histclo( nid_V ) 855 889 CALL histclo( nid_W ) 890 IF(ln_abl) CALL histclo( nid_A ) 856 891 ENDIF 857 892 ! … … 861 896 #endif 862 897 863 SUBROUTINE dia_wri_state( cdfile_name )898 SUBROUTINE dia_wri_state( Kmm, cdfile_name ) 864 899 !!--------------------------------------------------------------------- 865 900 !! *** ROUTINE dia_wri_state *** … … 874 909 !! File 'output.abort.nc' is created in case of abnormal job end 875 910 !!---------------------------------------------------------------------- 911 INTEGER , INTENT( in ) :: Kmm ! time level index 876 912 CHARACTER (len=* ), INTENT( in ) :: cdfile_name ! name of the file created 877 913 !! 878 INTEGER :: inum 914 INTEGER :: inum, jk 879 915 !!---------------------------------------------------------------------- 880 916 ! … … 890 926 #endif 891 927 892 CALL iom_rstput( 0, 0, inum, 'votemper', ts n(:,:,:,jp_tem) ) ! now temperature893 CALL iom_rstput( 0, 0, inum, 'vosaline', ts n(:,:,:,jp_sal) ) ! now salinity894 CALL iom_rstput( 0, 0, inum, 'sossheig', ssh n) ! sea surface height895 CALL iom_rstput( 0, 0, inum, 'vozocrtx', u n) ! now i-velocity896 CALL iom_rstput( 0, 0, inum, 'vomecrty', v n) ! now j-velocity928 CALL iom_rstput( 0, 0, inum, 'votemper', ts(:,:,:,jp_tem,Kmm) ) ! now temperature 929 CALL iom_rstput( 0, 0, inum, 'vosaline', ts(:,:,:,jp_sal,Kmm) ) ! now salinity 930 CALL iom_rstput( 0, 0, inum, 'sossheig', ssh(:,:,Kmm) ) ! sea surface height 931 CALL iom_rstput( 0, 0, inum, 'vozocrtx', uu(:,:,:,Kmm) ) ! now i-velocity 932 CALL iom_rstput( 0, 0, inum, 'vomecrty', vv(:,:,:,Kmm) ) ! now j-velocity 897 933 IF( ln_zad_Aimp ) THEN 898 CALL iom_rstput( 0, 0, inum, 'vovecrtz', w n+ wi ) ! now k-velocity934 CALL iom_rstput( 0, 0, inum, 'vovecrtz', ww + wi ) ! now k-velocity 899 935 ELSE 900 CALL iom_rstput( 0, 0, inum, 'vovecrtz', wn ) ! now k-velocity 901 ENDIF 936 CALL iom_rstput( 0, 0, inum, 'vovecrtz', ww ) ! now k-velocity 937 ENDIF 938 CALL iom_rstput( 0, 0, inum, 'risfdep', risfdep ) ! now k-velocity 939 CALL iom_rstput( 0, 0, inum, 'ht' , ht ) ! now water column height 940 941 IF ( ln_isf ) THEN 942 IF (ln_isfcav_mlt) THEN 943 CALL iom_rstput( 0, 0, inum, 'fwfisf_cav', fwfisf_cav ) ! now k-velocity 944 CALL iom_rstput( 0, 0, inum, 'rhisf_cav_tbl', rhisf_tbl_cav ) ! now k-velocity 945 CALL iom_rstput( 0, 0, inum, 'rfrac_cav_tbl', rfrac_tbl_cav ) ! now k-velocity 946 CALL iom_rstput( 0, 0, inum, 'misfkb_cav', REAL(misfkb_cav,8) ) ! now k-velocity 947 CALL iom_rstput( 0, 0, inum, 'misfkt_cav', REAL(misfkt_cav,8) ) ! now k-velocity 948 CALL iom_rstput( 0, 0, inum, 'mskisf_cav', REAL(mskisf_cav,8), ktype = jp_i1 ) 949 END IF 950 IF (ln_isfpar_mlt) THEN 951 CALL iom_rstput( 0, 0, inum, 'isfmsk_par', REAL(mskisf_par,8) ) ! now k-velocity 952 CALL iom_rstput( 0, 0, inum, 'fwfisf_par', fwfisf_par ) ! now k-velocity 953 CALL iom_rstput( 0, 0, inum, 'rhisf_par_tbl', rhisf_tbl_par ) ! now k-velocity 954 CALL iom_rstput( 0, 0, inum, 'rfrac_par_tbl', rfrac_tbl_par ) ! now k-velocity 955 CALL iom_rstput( 0, 0, inum, 'misfkb_par', REAL(misfkb_par,8) ) ! now k-velocity 956 CALL iom_rstput( 0, 0, inum, 'misfkt_par', REAL(misfkt_par,8) ) ! now k-velocity 957 CALL iom_rstput( 0, 0, inum, 'mskisf_par', REAL(mskisf_par,8), ktype = jp_i1 ) 958 END IF 959 END IF 960 902 961 IF( ALLOCATED(ahtu) ) THEN 903 962 CALL iom_rstput( 0, 0, inum, 'ahtu', ahtu ) ! aht at u-point … … 915 974 CALL iom_rstput( 0, 0, inum, 'sometauy', vtau ) ! j-wind stress 916 975 IF( .NOT.ln_linssh ) THEN 917 CALL iom_rstput( 0, 0, inum, 'vovvldep', gdept _n) ! T-cell depth918 CALL iom_rstput( 0, 0, inum, 'vovvle3t', e3t _n) ! T-cell thickness976 CALL iom_rstput( 0, 0, inum, 'vovvldep', gdept(:,:,:,Kmm) ) ! T-cell depth 977 CALL iom_rstput( 0, 0, inum, 'vovvle3t', e3t(:,:,:,Kmm) ) ! T-cell thickness 919 978 END IF 920 979 IF( ln_wave .AND. ln_sdw ) THEN … … 923 982 CALL iom_rstput( 0, 0, inum, 'sdvecrtz', wsd ) ! now StokesDrift k-velocity 924 983 ENDIF 984 IF ( ln_abl ) THEN 985 CALL iom_rstput ( 0, 0, inum, "uz1_abl", u_abl(:,:,2,nt_a ) ) ! now first level i-wind 986 CALL iom_rstput ( 0, 0, inum, "vz1_abl", v_abl(:,:,2,nt_a ) ) ! now first level j-wind 987 CALL iom_rstput ( 0, 0, inum, "tz1_abl", tq_abl(:,:,2,nt_a,1) ) ! now first level temperature 988 CALL iom_rstput ( 0, 0, inum, "qz1_abl", tq_abl(:,:,2,nt_a,2) ) ! now first level humidity 989 ENDIF 925 990 926 991 #if defined key_si3
Note: See TracChangeset
for help on using the changeset viewer.