[9762] | 1 | MODULE diawri |
---|
| 2 | !!====================================================================== |
---|
| 3 | !! *** MODULE diawri *** |
---|
| 4 | !! Ocean diagnostics : write ocean output files |
---|
| 5 | !!===================================================================== |
---|
| 6 | !! History : OPA ! 1991-03 (M.-A. Foujols) Original code |
---|
| 7 | !! 4.0 ! 1991-11 (G. Madec) |
---|
| 8 | !! ! 1992-06 (M. Imbard) correction restart file |
---|
| 9 | !! ! 1992-07 (M. Imbard) split into diawri and rstwri |
---|
| 10 | !! ! 1993-03 (M. Imbard) suppress writibm |
---|
| 11 | !! ! 1998-01 (C. Levy) NETCDF format using ioipsl INTERFACE |
---|
| 12 | !! ! 1999-02 (E. Guilyardi) name of netCDF files + variables |
---|
| 13 | !! 8.2 ! 2000-06 (M. Imbard) Original code (diabort.F) |
---|
| 14 | !! NEMO 1.0 ! 2002-06 (A.Bozec, E. Durand) Original code (diainit.F) |
---|
| 15 | !! - ! 2002-09 (G. Madec) F90: Free form and module |
---|
| 16 | !! - ! 2002-12 (G. Madec) merge of diabort and diainit, F90 |
---|
| 17 | !! ! 2005-11 (V. Garnier) Surface pressure gradient organization |
---|
| 18 | !! 3.2 ! 2008-11 (B. Lemaire) creation from old diawri |
---|
| 19 | !! 3.7 ! 2014-01 (G. Madec) remove eddy induced velocity from no-IOM output |
---|
| 20 | !! ! change name of output variables in dia_wri_state |
---|
| 21 | !!---------------------------------------------------------------------- |
---|
| 22 | |
---|
| 23 | !!---------------------------------------------------------------------- |
---|
| 24 | !! dia_wri : create the standart output files |
---|
| 25 | !! dia_wri_state : create an output NetCDF file for a single instantaeous ocean state and forcing fields |
---|
| 26 | !!---------------------------------------------------------------------- |
---|
| 27 | USE oce ! ocean dynamics and tracers |
---|
| 28 | USE dom_oce ! ocean space and time domain |
---|
| 29 | USE phycst ! physical constants |
---|
| 30 | USE dianam ! build name of file (routine) |
---|
| 31 | USE diahth ! thermocline diagnostics |
---|
| 32 | USE dynadv , ONLY: ln_dynadv_vec |
---|
| 33 | USE icb_oce ! Icebergs |
---|
| 34 | USE icbdia ! Iceberg budgets |
---|
| 35 | USE ldftra ! lateral physics: eddy diffusivity coef. |
---|
| 36 | USE ldfdyn ! lateral physics: eddy viscosity coef. |
---|
| 37 | USE sbc_oce ! Surface boundary condition: ocean fields |
---|
| 38 | USE sbc_ice ! Surface boundary condition: ice fields |
---|
| 39 | USE sbcssr ! restoring term toward SST/SSS climatology |
---|
| 40 | USE sbcwave ! wave parameters |
---|
| 41 | USE wet_dry ! wetting and drying |
---|
| 42 | USE zdf_oce ! ocean vertical physics |
---|
| 43 | USE zdfdrg ! ocean vertical physics: top/bottom friction |
---|
| 44 | USE zdfmxl ! mixed layer |
---|
| 45 | ! |
---|
| 46 | USE lbclnk ! ocean lateral boundary conditions (or mpp link) |
---|
| 47 | USE in_out_manager ! I/O manager |
---|
| 48 | USE diatmb ! Top,middle,bottom output |
---|
| 49 | USE dia25h ! 25h Mean output |
---|
| 50 | USE iom ! |
---|
| 51 | USE ioipsl ! |
---|
| 52 | |
---|
| 53 | #if defined key_si3 |
---|
[11758] | 54 | USE ice |
---|
[9762] | 55 | USE icewri |
---|
| 56 | #endif |
---|
| 57 | USE lib_mpp ! MPP library |
---|
| 58 | USE timing ! preformance summary |
---|
[11758] | 59 | USE diu_bulk ! diurnal warm layer |
---|
| 60 | USE diu_coolskin ! Cool skin |
---|
[9762] | 61 | |
---|
| 62 | IMPLICIT NONE |
---|
| 63 | PRIVATE |
---|
| 64 | |
---|
| 65 | PUBLIC dia_wri ! routines called by step.F90 |
---|
| 66 | PUBLIC dia_wri_state |
---|
| 67 | PUBLIC dia_wri_alloc ! Called by nemogcm module |
---|
| 68 | |
---|
| 69 | INTEGER :: nid_T, nz_T, nh_T, ndim_T, ndim_hT ! grid_T file |
---|
| 70 | INTEGER :: nb_T , ndim_bT ! grid_T file |
---|
| 71 | INTEGER :: nid_U, nz_U, nh_U, ndim_U, ndim_hU ! grid_U file |
---|
| 72 | INTEGER :: nid_V, nz_V, nh_V, ndim_V, ndim_hV ! grid_V file |
---|
| 73 | INTEGER :: nid_W, nz_W, nh_W ! grid_W file |
---|
| 74 | INTEGER :: ndex(1) ! ??? |
---|
| 75 | INTEGER, SAVE, ALLOCATABLE, DIMENSION(:) :: ndex_hT, ndex_hU, ndex_hV |
---|
| 76 | INTEGER, SAVE, ALLOCATABLE, DIMENSION(:) :: ndex_T, ndex_U, ndex_V |
---|
| 77 | INTEGER, SAVE, ALLOCATABLE, DIMENSION(:) :: ndex_bT |
---|
| 78 | |
---|
| 79 | !! * Substitutions |
---|
| 80 | # include "vectopt_loop_substitute.h90" |
---|
| 81 | !!---------------------------------------------------------------------- |
---|
| 82 | !! NEMO/OCE 4.0 , NEMO Consortium (2018) |
---|
| 83 | !! $Id: diawri.F90 9598 2018-05-15 22:47:16Z nicolasmartin $ |
---|
| 84 | !! Software governed by the CeCILL licence (./LICENSE) |
---|
| 85 | !!---------------------------------------------------------------------- |
---|
| 86 | CONTAINS |
---|
| 87 | |
---|
| 88 | #if defined key_iomput |
---|
| 89 | !!---------------------------------------------------------------------- |
---|
| 90 | !! 'key_iomput' use IOM library |
---|
| 91 | !!---------------------------------------------------------------------- |
---|
| 92 | INTEGER FUNCTION dia_wri_alloc() |
---|
| 93 | ! |
---|
| 94 | dia_wri_alloc = 0 |
---|
| 95 | ! |
---|
| 96 | END FUNCTION dia_wri_alloc |
---|
| 97 | |
---|
| 98 | |
---|
[11758] | 99 | SUBROUTINE dia_wri( kt, Kmm ) |
---|
[9762] | 100 | !!--------------------------------------------------------------------- |
---|
| 101 | !! *** ROUTINE dia_wri *** |
---|
| 102 | !! |
---|
| 103 | !! ** Purpose : Standard output of opa: dynamics and tracer fields |
---|
| 104 | !! NETCDF format is used by default |
---|
| 105 | !! |
---|
| 106 | !! ** Method : use iom_put |
---|
| 107 | !!---------------------------------------------------------------------- |
---|
| 108 | INTEGER, INTENT( in ) :: kt ! ocean time-step index |
---|
[11758] | 109 | INTEGER, INTENT( in ) :: Kmm ! ocean time level index |
---|
[9762] | 110 | !! |
---|
| 111 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
| 112 | INTEGER :: ikbot ! local integer |
---|
| 113 | REAL(wp):: zztmp , zztmpx ! local scalar |
---|
| 114 | REAL(wp):: zztmp2, zztmpy ! - - |
---|
| 115 | REAL(wp), DIMENSION(jpi,jpj) :: z2d ! 2D workspace |
---|
| 116 | REAL(wp), DIMENSION(jpi,jpj,jpk) :: z3d ! 3D workspace |
---|
| 117 | !!---------------------------------------------------------------------- |
---|
| 118 | ! |
---|
| 119 | IF( ln_timing ) CALL timing_start('dia_wri') |
---|
| 120 | ! |
---|
| 121 | ! Output the initial state and forcings |
---|
| 122 | IF( ninist == 1 ) THEN |
---|
[11758] | 123 | CALL dia_wri_state( Kmm, 'output.init' ) |
---|
[9762] | 124 | ninist = 0 |
---|
| 125 | ENDIF |
---|
| 126 | |
---|
| 127 | ! Output of initial vertical scale factor |
---|
| 128 | CALL iom_put("e3t_0", e3t_0(:,:,:) ) |
---|
[10166] | 129 | CALL iom_put("e3u_0", e3u_0(:,:,:) ) |
---|
| 130 | CALL iom_put("e3v_0", e3v_0(:,:,:) ) |
---|
[9762] | 131 | ! |
---|
[11758] | 132 | CALL iom_put( "e3t" , e3t(:,:,:,Kmm) ) |
---|
| 133 | CALL iom_put( "e3u" , e3u(:,:,:,Kmm) ) |
---|
| 134 | CALL iom_put( "e3v" , e3v(:,:,:,Kmm) ) |
---|
| 135 | CALL iom_put( "e3w" , e3w(:,:,:,Kmm) ) |
---|
[9762] | 136 | IF( iom_use("e3tdef") ) & |
---|
[11758] | 137 | CALL iom_put( "e3tdef" , ( ( e3t(:,:,:,Kmm) - e3t_0(:,:,:) ) / e3t_0(:,:,:) * 100 * tmask(:,:,:) ) ** 2 ) |
---|
[9762] | 138 | |
---|
| 139 | IF( ll_wd ) THEN |
---|
[11758] | 140 | CALL iom_put( "ssh" , (ssh(:,:,Kmm)+ssh_ref)*tmask(:,:,1) ) ! sea surface height (brought back to the reference used for wetting and drying) |
---|
[9762] | 141 | ELSE |
---|
[11758] | 142 | CALL iom_put( "ssh" , ssh(:,:,Kmm) ) ! sea surface height |
---|
[9762] | 143 | ENDIF |
---|
| 144 | |
---|
| 145 | IF( iom_use("wetdep") ) & ! wet depth |
---|
[11758] | 146 | CALL iom_put( "wetdep" , ht_0(:,:) + ssh(:,:,Kmm) ) |
---|
[9762] | 147 | |
---|
[11758] | 148 | CALL iom_put( "toce", ts(:,:,:,jp_tem,Kmm) ) ! 3D temperature |
---|
| 149 | CALL iom_put( "sst", ts(:,:,1,jp_tem,Kmm) ) ! surface temperature |
---|
[9762] | 150 | IF ( iom_use("sbt") ) THEN |
---|
| 151 | DO jj = 1, jpj |
---|
| 152 | DO ji = 1, jpi |
---|
| 153 | ikbot = mbkt(ji,jj) |
---|
[11758] | 154 | z2d(ji,jj) = ts(ji,jj,ikbot,jp_tem,Kmm) |
---|
[9762] | 155 | END DO |
---|
| 156 | END DO |
---|
| 157 | CALL iom_put( "sbt", z2d ) ! bottom temperature |
---|
| 158 | ENDIF |
---|
| 159 | |
---|
[11758] | 160 | CALL iom_put( "soce", ts(:,:,:,jp_sal,Kmm) ) ! 3D salinity |
---|
| 161 | CALL iom_put( "sss", ts(:,:,1,jp_sal,Kmm) ) ! surface salinity |
---|
[9762] | 162 | IF ( iom_use("sbs") ) THEN |
---|
| 163 | DO jj = 1, jpj |
---|
| 164 | DO ji = 1, jpi |
---|
| 165 | ikbot = mbkt(ji,jj) |
---|
[11758] | 166 | z2d(ji,jj) = ts(ji,jj,ikbot,jp_sal,Kmm) |
---|
[9762] | 167 | END DO |
---|
| 168 | END DO |
---|
| 169 | CALL iom_put( "sbs", z2d ) ! bottom salinity |
---|
| 170 | ENDIF |
---|
| 171 | |
---|
| 172 | IF ( iom_use("taubot") ) THEN ! bottom stress |
---|
| 173 | zztmp = rau0 * 0.25 |
---|
| 174 | z2d(:,:) = 0._wp |
---|
| 175 | DO jj = 2, jpjm1 |
---|
| 176 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[11758] | 177 | zztmp2 = ( ( rCdU_bot(ji+1,jj)+rCdU_bot(ji ,jj) ) * uu(ji ,jj,mbku(ji ,jj),Kmm) )**2 & |
---|
| 178 | & + ( ( rCdU_bot(ji ,jj)+rCdU_bot(ji-1,jj) ) * uu(ji-1,jj,mbku(ji-1,jj),Kmm) )**2 & |
---|
| 179 | & + ( ( rCdU_bot(ji,jj+1)+rCdU_bot(ji,jj ) ) * vv(ji,jj ,mbkv(ji,jj ),Kmm) )**2 & |
---|
| 180 | & + ( ( rCdU_bot(ji,jj )+rCdU_bot(ji,jj-1) ) * vv(ji,jj-1,mbkv(ji,jj-1),Kmm) )**2 |
---|
[9762] | 181 | z2d(ji,jj) = zztmp * SQRT( zztmp2 ) * tmask(ji,jj,1) |
---|
| 182 | ! |
---|
| 183 | END DO |
---|
| 184 | END DO |
---|
[10170] | 185 | CALL lbc_lnk( 'diawri', z2d, 'T', 1. ) |
---|
[9762] | 186 | CALL iom_put( "taubot", z2d ) |
---|
| 187 | ENDIF |
---|
| 188 | |
---|
[11758] | 189 | CALL iom_put( "uoce", uu(:,:,:,Kmm) ) ! 3D i-current |
---|
| 190 | CALL iom_put( "ssu", uu(:,:,1,Kmm) ) ! surface i-current |
---|
[9762] | 191 | IF ( iom_use("sbu") ) THEN |
---|
| 192 | DO jj = 1, jpj |
---|
| 193 | DO ji = 1, jpi |
---|
| 194 | ikbot = mbku(ji,jj) |
---|
[11758] | 195 | z2d(ji,jj) = uu(ji,jj,ikbot,Kmm) |
---|
[9762] | 196 | END DO |
---|
| 197 | END DO |
---|
| 198 | CALL iom_put( "sbu", z2d ) ! bottom i-current |
---|
| 199 | ENDIF |
---|
| 200 | |
---|
[11758] | 201 | CALL iom_put( "voce", vv(:,:,:,Kmm) ) ! 3D j-current |
---|
| 202 | CALL iom_put( "ssv", vv(:,:,1,Kmm) ) ! surface j-current |
---|
[9762] | 203 | IF ( iom_use("sbv") ) THEN |
---|
| 204 | DO jj = 1, jpj |
---|
| 205 | DO ji = 1, jpi |
---|
| 206 | ikbot = mbkv(ji,jj) |
---|
[11758] | 207 | z2d(ji,jj) = vv(ji,jj,ikbot,Kmm) |
---|
[9762] | 208 | END DO |
---|
| 209 | END DO |
---|
| 210 | CALL iom_put( "sbv", z2d ) ! bottom j-current |
---|
| 211 | ENDIF |
---|
| 212 | |
---|
[11822] | 213 | IF( ln_zad_Aimp ) ww = ww + wi ! Recombine explicit and implicit parts of vertical velocity for diagnostic output |
---|
| 214 | ! |
---|
[11758] | 215 | CALL iom_put( "woce", ww ) ! vertical velocity |
---|
[9762] | 216 | IF( iom_use('w_masstr') .OR. iom_use('w_masstr2') ) THEN ! vertical mass transport & its square value |
---|
| 217 | ! Caution: in the VVL case, it only correponds to the baroclinic mass transport. |
---|
| 218 | z2d(:,:) = rau0 * e1e2t(:,:) |
---|
| 219 | DO jk = 1, jpk |
---|
[11758] | 220 | z3d(:,:,jk) = ww(:,:,jk) * z2d(:,:) |
---|
[9762] | 221 | END DO |
---|
| 222 | CALL iom_put( "w_masstr" , z3d ) |
---|
| 223 | IF( iom_use('w_masstr2') ) CALL iom_put( "w_masstr2", z3d(:,:,:) * z3d(:,:,:) ) |
---|
| 224 | ENDIF |
---|
[11822] | 225 | ! |
---|
| 226 | IF( ln_zad_Aimp ) ww = ww - wi ! Remove implicit part of vertical velocity that was added for diagnostic output |
---|
[9762] | 227 | |
---|
| 228 | CALL iom_put( "avt" , avt ) ! T vert. eddy diff. coef. |
---|
| 229 | CALL iom_put( "avs" , avs ) ! S vert. eddy diff. coef. |
---|
| 230 | CALL iom_put( "avm" , avm ) ! T vert. eddy visc. coef. |
---|
| 231 | |
---|
| 232 | IF( iom_use('logavt') ) CALL iom_put( "logavt", LOG( MAX( 1.e-20_wp, avt(:,:,:) ) ) ) |
---|
| 233 | IF( iom_use('logavs') ) CALL iom_put( "logavs", LOG( MAX( 1.e-20_wp, avs(:,:,:) ) ) ) |
---|
| 234 | |
---|
| 235 | IF ( iom_use("sstgrad") .OR. iom_use("sstgrad2") ) THEN |
---|
| 236 | DO jj = 2, jpjm1 ! sst gradient |
---|
| 237 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[11758] | 238 | zztmp = ts(ji,jj,1,jp_tem,Kmm) |
---|
| 239 | 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) |
---|
| 240 | 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) |
---|
[9762] | 241 | z2d(ji,jj) = 0.25 * ( zztmpx * zztmpx + zztmpy * zztmpy ) & |
---|
| 242 | & * umask(ji,jj,1) * umask(ji-1,jj,1) * vmask(ji,jj,1) * umask(ji,jj-1,1) |
---|
| 243 | END DO |
---|
| 244 | END DO |
---|
[10170] | 245 | CALL lbc_lnk( 'diawri', z2d, 'T', 1. ) |
---|
[9762] | 246 | CALL iom_put( "sstgrad2", z2d ) ! square of module of sst gradient |
---|
| 247 | z2d(:,:) = SQRT( z2d(:,:) ) |
---|
| 248 | CALL iom_put( "sstgrad" , z2d ) ! module of sst gradient |
---|
| 249 | ENDIF |
---|
| 250 | |
---|
| 251 | ! heat and salt contents |
---|
| 252 | IF( iom_use("heatc") ) THEN |
---|
| 253 | z2d(:,:) = 0._wp |
---|
| 254 | DO jk = 1, jpkm1 |
---|
| 255 | DO jj = 1, jpj |
---|
| 256 | DO ji = 1, jpi |
---|
[11758] | 257 | z2d(ji,jj) = z2d(ji,jj) + e3t(ji,jj,jk,Kmm) * ts(ji,jj,jk,jp_tem,Kmm) * tmask(ji,jj,jk) |
---|
[9762] | 258 | END DO |
---|
| 259 | END DO |
---|
| 260 | END DO |
---|
| 261 | CALL iom_put( "heatc", rau0_rcp * z2d ) ! vertically integrated heat content (J/m2) |
---|
| 262 | ENDIF |
---|
| 263 | |
---|
| 264 | IF( iom_use("saltc") ) THEN |
---|
| 265 | z2d(:,:) = 0._wp |
---|
| 266 | DO jk = 1, jpkm1 |
---|
| 267 | DO jj = 1, jpj |
---|
| 268 | DO ji = 1, jpi |
---|
[11758] | 269 | z2d(ji,jj) = z2d(ji,jj) + e3t(ji,jj,jk,Kmm) * ts(ji,jj,jk,jp_sal,Kmm) * tmask(ji,jj,jk) |
---|
[9762] | 270 | END DO |
---|
| 271 | END DO |
---|
| 272 | END DO |
---|
| 273 | CALL iom_put( "saltc", rau0 * z2d ) ! vertically integrated salt content (PSU*kg/m2) |
---|
| 274 | ENDIF |
---|
| 275 | ! |
---|
| 276 | IF ( iom_use("eken") ) THEN |
---|
| 277 | z3d(:,:,jpk) = 0._wp |
---|
| 278 | DO jk = 1, jpkm1 |
---|
| 279 | DO jj = 2, jpjm1 |
---|
| 280 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[11758] | 281 | zztmp = 0.25_wp * r1_e1e2t(ji,jj) / e3t(ji,jj,jk,Kmm) |
---|
| 282 | z3d(ji,jj,jk) = zztmp * ( uu(ji-1,jj,jk,Kmm)**2 * e2u(ji-1,jj) * e3u(ji-1,jj,jk,Kmm) & |
---|
| 283 | & + uu(ji ,jj,jk,Kmm)**2 * e2u(ji ,jj) * e3u(ji ,jj,jk,Kmm) & |
---|
| 284 | & + vv(ji,jj-1,jk,Kmm)**2 * e1v(ji,jj-1) * e3v(ji,jj-1,jk,Kmm) & |
---|
| 285 | & + vv(ji,jj ,jk,Kmm)**2 * e1v(ji,jj ) * e3v(ji,jj ,jk,Kmm) ) |
---|
[9762] | 286 | END DO |
---|
| 287 | END DO |
---|
| 288 | END DO |
---|
[10170] | 289 | CALL lbc_lnk( 'diawri', z3d, 'T', 1. ) |
---|
[9762] | 290 | CALL iom_put( "eken", z3d ) ! kinetic energy |
---|
| 291 | ENDIF |
---|
| 292 | ! |
---|
[11758] | 293 | CALL iom_put( "hdiv", hdiv ) ! Horizontal divergence |
---|
[9762] | 294 | ! |
---|
| 295 | IF( iom_use("u_masstr") .OR. iom_use("u_masstr_vint") .OR. iom_use("u_heattr") .OR. iom_use("u_salttr") ) THEN |
---|
| 296 | z3d(:,:,jpk) = 0.e0 |
---|
| 297 | z2d(:,:) = 0.e0 |
---|
| 298 | DO jk = 1, jpkm1 |
---|
[11758] | 299 | z3d(:,:,jk) = rau0 * uu(:,:,jk,Kmm) * e2u(:,:) * e3u(:,:,jk,Kmm) * umask(:,:,jk) |
---|
[9762] | 300 | z2d(:,:) = z2d(:,:) + z3d(:,:,jk) |
---|
| 301 | END DO |
---|
| 302 | CALL iom_put( "u_masstr" , z3d ) ! mass transport in i-direction |
---|
| 303 | CALL iom_put( "u_masstr_vint", z2d ) ! mass transport in i-direction vertical sum |
---|
| 304 | ENDIF |
---|
| 305 | |
---|
| 306 | IF( iom_use("u_heattr") ) THEN |
---|
| 307 | z2d(:,:) = 0._wp |
---|
| 308 | DO jk = 1, jpkm1 |
---|
| 309 | DO jj = 2, jpjm1 |
---|
| 310 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[11758] | 311 | 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) ) |
---|
[9762] | 312 | END DO |
---|
| 313 | END DO |
---|
| 314 | END DO |
---|
[10170] | 315 | CALL lbc_lnk( 'diawri', z2d, 'U', -1. ) |
---|
[9762] | 316 | CALL iom_put( "u_heattr", 0.5*rcp * z2d ) ! heat transport in i-direction |
---|
| 317 | ENDIF |
---|
| 318 | |
---|
| 319 | IF( iom_use("u_salttr") ) THEN |
---|
| 320 | z2d(:,:) = 0.e0 |
---|
| 321 | DO jk = 1, jpkm1 |
---|
| 322 | DO jj = 2, jpjm1 |
---|
| 323 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[11758] | 324 | 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) ) |
---|
[9762] | 325 | END DO |
---|
| 326 | END DO |
---|
| 327 | END DO |
---|
[10170] | 328 | CALL lbc_lnk( 'diawri', z2d, 'U', -1. ) |
---|
[9762] | 329 | CALL iom_put( "u_salttr", 0.5 * z2d ) ! heat transport in i-direction |
---|
| 330 | ENDIF |
---|
| 331 | |
---|
| 332 | |
---|
| 333 | IF( iom_use("v_masstr") .OR. iom_use("v_heattr") .OR. iom_use("v_salttr") ) THEN |
---|
| 334 | z3d(:,:,jpk) = 0.e0 |
---|
| 335 | DO jk = 1, jpkm1 |
---|
[11758] | 336 | z3d(:,:,jk) = rau0 * vv(:,:,jk,Kmm) * e1v(:,:) * e3v(:,:,jk,Kmm) * vmask(:,:,jk) |
---|
[9762] | 337 | END DO |
---|
| 338 | CALL iom_put( "v_masstr", z3d ) ! mass transport in j-direction |
---|
| 339 | ENDIF |
---|
| 340 | |
---|
| 341 | IF( iom_use("v_heattr") ) THEN |
---|
| 342 | z2d(:,:) = 0.e0 |
---|
| 343 | DO jk = 1, jpkm1 |
---|
| 344 | DO jj = 2, jpjm1 |
---|
| 345 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[11758] | 346 | 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) ) |
---|
[9762] | 347 | END DO |
---|
| 348 | END DO |
---|
| 349 | END DO |
---|
[10170] | 350 | CALL lbc_lnk( 'diawri', z2d, 'V', -1. ) |
---|
[9762] | 351 | CALL iom_put( "v_heattr", 0.5*rcp * z2d ) ! heat transport in j-direction |
---|
| 352 | ENDIF |
---|
| 353 | |
---|
| 354 | IF( iom_use("v_salttr") ) THEN |
---|
| 355 | z2d(:,:) = 0._wp |
---|
| 356 | DO jk = 1, jpkm1 |
---|
| 357 | DO jj = 2, jpjm1 |
---|
| 358 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[11758] | 359 | 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) ) |
---|
[9762] | 360 | END DO |
---|
| 361 | END DO |
---|
| 362 | END DO |
---|
[10170] | 363 | CALL lbc_lnk( 'diawri', z2d, 'V', -1. ) |
---|
[9762] | 364 | CALL iom_put( "v_salttr", 0.5 * z2d ) ! heat transport in j-direction |
---|
| 365 | ENDIF |
---|
| 366 | |
---|
| 367 | IF( iom_use("tosmint") ) THEN |
---|
| 368 | z2d(:,:) = 0._wp |
---|
| 369 | DO jk = 1, jpkm1 |
---|
| 370 | DO jj = 2, jpjm1 |
---|
| 371 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[11758] | 372 | z2d(ji,jj) = z2d(ji,jj) + e3t(ji,jj,jk,Kmm) * ts(ji,jj,jk,jp_tem,Kmm) |
---|
[9762] | 373 | END DO |
---|
| 374 | END DO |
---|
| 375 | END DO |
---|
[10170] | 376 | CALL lbc_lnk( 'diawri', z2d, 'T', -1. ) |
---|
[9762] | 377 | CALL iom_put( "tosmint", rau0 * z2d ) ! Vertical integral of temperature |
---|
| 378 | ENDIF |
---|
| 379 | IF( iom_use("somint") ) THEN |
---|
| 380 | z2d(:,:)=0._wp |
---|
| 381 | DO jk = 1, jpkm1 |
---|
| 382 | DO jj = 2, jpjm1 |
---|
| 383 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[11758] | 384 | z2d(ji,jj) = z2d(ji,jj) + e3t(ji,jj,jk,Kmm) * ts(ji,jj,jk,jp_sal,Kmm) |
---|
[9762] | 385 | END DO |
---|
| 386 | END DO |
---|
| 387 | END DO |
---|
[10170] | 388 | CALL lbc_lnk( 'diawri', z2d, 'T', -1. ) |
---|
[9762] | 389 | CALL iom_put( "somint", rau0 * z2d ) ! Vertical integral of salinity |
---|
| 390 | ENDIF |
---|
| 391 | |
---|
| 392 | CALL iom_put( "bn2", rn2 ) ! Brunt-Vaisala buoyancy frequency (N^2) |
---|
| 393 | ! |
---|
| 394 | |
---|
[11758] | 395 | IF (ln_diatmb) CALL dia_tmb( Kmm ) ! tmb values |
---|
[9762] | 396 | |
---|
[11758] | 397 | IF (ln_dia25h) CALL dia_25h( kt, Kmm ) ! 25h averaging |
---|
[9762] | 398 | |
---|
| 399 | IF( ln_timing ) CALL timing_stop('dia_wri') |
---|
| 400 | ! |
---|
| 401 | END SUBROUTINE dia_wri |
---|
| 402 | |
---|
| 403 | #else |
---|
| 404 | !!---------------------------------------------------------------------- |
---|
| 405 | !! Default option use IOIPSL library |
---|
| 406 | !!---------------------------------------------------------------------- |
---|
| 407 | |
---|
| 408 | INTEGER FUNCTION dia_wri_alloc() |
---|
| 409 | ! |
---|
| 410 | dia_wri_alloc = 0 |
---|
| 411 | ! |
---|
| 412 | END FUNCTION dia_wri_alloc |
---|
| 413 | |
---|
[11758] | 414 | |
---|
| 415 | SUBROUTINE dia_wri( kt, Kmm ) |
---|
[10179] | 416 | |
---|
[11758] | 417 | INTEGER, INTENT( in ) :: kt ! ocean time-step index |
---|
| 418 | INTEGER, INTENT( in ) :: Kmm ! ocean time level index |
---|
[10179] | 419 | IF( ninist == 1 ) THEN !== Output the initial state and forcings ==! |
---|
[11758] | 420 | CALL dia_wri_state( Kmm, 'output.init' ) |
---|
[10179] | 421 | ninist = 0 |
---|
| 422 | ENDIF |
---|
[11758] | 423 | ! |
---|
| 424 | ! 0. Initialisation |
---|
| 425 | ! ----------------- |
---|
[10179] | 426 | |
---|
[9762] | 427 | END SUBROUTINE dia_wri |
---|
| 428 | |
---|
| 429 | #endif |
---|
| 430 | |
---|
[11758] | 431 | SUBROUTINE dia_wri_state( Kmm, cdfile_name ) |
---|
[9762] | 432 | !!--------------------------------------------------------------------- |
---|
| 433 | !! *** ROUTINE dia_wri_state *** |
---|
| 434 | !! |
---|
| 435 | !! ** Purpose : create a NetCDF file named cdfile_name which contains |
---|
| 436 | !! the instantaneous ocean state and forcing fields. |
---|
| 437 | !! Used to find errors in the initial state or save the last |
---|
| 438 | !! ocean state in case of abnormal end of a simulation |
---|
| 439 | !! |
---|
| 440 | !! ** Method : NetCDF files using ioipsl |
---|
| 441 | !! File 'output.init.nc' is created if ninist = 1 (namelist) |
---|
| 442 | !! File 'output.abort.nc' is created in case of abnormal job end |
---|
| 443 | !!---------------------------------------------------------------------- |
---|
[11758] | 444 | INTEGER , INTENT( in ) :: Kmm ! time level index |
---|
[9762] | 445 | CHARACTER (len=* ), INTENT( in ) :: cdfile_name ! name of the file created |
---|
[10358] | 446 | !! |
---|
| 447 | INTEGER :: inum |
---|
[9762] | 448 | !!---------------------------------------------------------------------- |
---|
| 449 | ! |
---|
| 450 | IF(lwp) WRITE(numout,*) |
---|
| 451 | IF(lwp) WRITE(numout,*) 'dia_wri_state : single instantaneous ocean state' |
---|
| 452 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~ and forcing fields file created ' |
---|
[10358] | 453 | IF(lwp) WRITE(numout,*) ' and named :', cdfile_name, '...nc' |
---|
[9762] | 454 | |
---|
| 455 | #if defined key_si3 |
---|
[10410] | 456 | CALL iom_open( TRIM(cdfile_name), inum, ldwrt = .TRUE., kdlev = jpl ) |
---|
[9762] | 457 | #else |
---|
[10410] | 458 | CALL iom_open( TRIM(cdfile_name), inum, ldwrt = .TRUE. ) |
---|
[9762] | 459 | #endif |
---|
| 460 | |
---|
[11758] | 461 | CALL iom_rstput( 0, 0, inum, 'votemper', ts(:,:,:,jp_tem,Kmm) ) ! now temperature |
---|
| 462 | CALL iom_rstput( 0, 0, inum, 'vosaline', ts(:,:,:,jp_sal,Kmm) ) ! now salinity |
---|
| 463 | CALL iom_rstput( 0, 0, inum, 'sossheig', ssh(:,:,Kmm) ) ! sea surface height |
---|
| 464 | CALL iom_rstput( 0, 0, inum, 'vozocrtx', uu(:,:,:,Kmm) ) ! now i-velocity |
---|
| 465 | CALL iom_rstput( 0, 0, inum, 'vomecrty', vv(:,:,:,Kmm) ) ! now j-velocity |
---|
[11822] | 466 | IF( ln_zad_Aimp ) THEN |
---|
| 467 | CALL iom_rstput( 0, 0, inum, 'vovecrtz', ww + wi ) ! now k-velocity |
---|
| 468 | ELSE |
---|
| 469 | CALL iom_rstput( 0, 0, inum, 'vovecrtz', ww ) ! now k-velocity |
---|
| 470 | ENDIF |
---|
[9762] | 471 | IF( ALLOCATED(ahtu) ) THEN |
---|
[10358] | 472 | CALL iom_rstput( 0, 0, inum, 'ahtu', ahtu ) ! aht at u-point |
---|
| 473 | CALL iom_rstput( 0, 0, inum, 'ahtv', ahtv ) ! aht at v-point |
---|
[9762] | 474 | ENDIF |
---|
| 475 | IF( ALLOCATED(ahmt) ) THEN |
---|
[10358] | 476 | CALL iom_rstput( 0, 0, inum, 'ahmt', ahmt ) ! ahmt at u-point |
---|
| 477 | CALL iom_rstput( 0, 0, inum, 'ahmf', ahmf ) ! ahmf at v-point |
---|
[9762] | 478 | ENDIF |
---|
[10358] | 479 | CALL iom_rstput( 0, 0, inum, 'sowaflup', emp - rnf ) ! freshwater budget |
---|
| 480 | CALL iom_rstput( 0, 0, inum, 'sohefldo', qsr + qns ) ! total heat flux |
---|
| 481 | CALL iom_rstput( 0, 0, inum, 'soshfldo', qsr ) ! solar heat flux |
---|
| 482 | CALL iom_rstput( 0, 0, inum, 'soicecov', fr_i ) ! ice fraction |
---|
| 483 | CALL iom_rstput( 0, 0, inum, 'sozotaux', utau ) ! i-wind stress |
---|
| 484 | CALL iom_rstput( 0, 0, inum, 'sometauy', vtau ) ! j-wind stress |
---|
[9762] | 485 | IF( .NOT.ln_linssh ) THEN |
---|
[11758] | 486 | CALL iom_rstput( 0, 0, inum, 'vovvldep', gdept(:,:,:,Kmm) ) ! T-cell depth |
---|
| 487 | CALL iom_rstput( 0, 0, inum, 'vovvle3t', e3t(:,:,:,Kmm) ) ! T-cell thickness |
---|
[10358] | 488 | END IF |
---|
[9762] | 489 | IF( ln_wave .AND. ln_sdw ) THEN |
---|
[10358] | 490 | CALL iom_rstput( 0, 0, inum, 'sdzocrtx', usd ) ! now StokesDrift i-velocity |
---|
| 491 | CALL iom_rstput( 0, 0, inum, 'sdmecrty', vsd ) ! now StokesDrift j-velocity |
---|
| 492 | CALL iom_rstput( 0, 0, inum, 'sdvecrtz', wsd ) ! now StokesDrift k-velocity |
---|
[9762] | 493 | ENDIF |
---|
[10358] | 494 | |
---|
| 495 | #if defined key_si3 |
---|
| 496 | IF( nn_ice == 2 ) THEN ! condition needed in case agrif + ice-model but no-ice in child grid |
---|
| 497 | CALL ice_wri_state( inum ) |
---|
[9762] | 498 | ENDIF |
---|
| 499 | #endif |
---|
[10358] | 500 | ! |
---|
| 501 | CALL iom_close( inum ) |
---|
[9762] | 502 | ! |
---|
| 503 | END SUBROUTINE dia_wri_state |
---|
| 504 | |
---|
| 505 | !!====================================================================== |
---|
| 506 | END MODULE diawri |
---|