[3] | 1 | MODULE diawri |
---|
| 2 | !!====================================================================== |
---|
| 3 | !! *** MODULE diawri *** |
---|
| 4 | !! Ocean diagnostics : write ocean output files |
---|
| 5 | !!===================================================================== |
---|
[2528] | 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 |
---|
[5836] | 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 |
---|
[2528] | 21 | !!---------------------------------------------------------------------- |
---|
[3] | 22 | |
---|
| 23 | !!---------------------------------------------------------------------- |
---|
[2528] | 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 | !!---------------------------------------------------------------------- |
---|
[3] | 27 | USE oce ! ocean dynamics and tracers |
---|
| 28 | USE dom_oce ! ocean space and time domain |
---|
[4292] | 29 | USE dynadv, ONLY: ln_dynadv_vec |
---|
[3] | 30 | USE zdf_oce ! ocean vertical physics |
---|
[5836] | 31 | USE ldftra ! lateral physics: eddy diffusivity coef. |
---|
[6140] | 32 | USE ldfdyn ! lateral physics: eddy viscosity coef. |
---|
[888] | 33 | USE sbc_oce ! Surface boundary condition: ocean fields |
---|
| 34 | USE sbc_ice ! Surface boundary condition: ice fields |
---|
[3609] | 35 | USE icb_oce ! Icebergs |
---|
| 36 | USE icbdia ! Iceberg budgets |
---|
[888] | 37 | USE sbcssr ! restoring term toward SST/SSS climatology |
---|
[3] | 38 | USE phycst ! physical constants |
---|
| 39 | USE zdfmxl ! mixed layer |
---|
| 40 | USE dianam ! build name of file (routine) |
---|
| 41 | USE zdfddm ! vertical physics: double diffusion |
---|
| 42 | USE diahth ! thermocline diagnostics |
---|
[6140] | 43 | ! |
---|
[3] | 44 | USE lbclnk ! ocean lateral boundary conditions (or mpp link) |
---|
| 45 | USE in_out_manager ! I/O manager |
---|
[6140] | 46 | USE diatmb ! Top,middle,bottom output |
---|
| 47 | USE dia25h ! 25h Mean output |
---|
[1359] | 48 | USE iom |
---|
[460] | 49 | USE ioipsl |
---|
[5463] | 50 | |
---|
[1482] | 51 | #if defined key_lim2 |
---|
| 52 | USE limwri_2 |
---|
[4161] | 53 | #elif defined key_lim3 |
---|
| 54 | USE limwri |
---|
[1482] | 55 | #endif |
---|
[2715] | 56 | USE lib_mpp ! MPP library |
---|
[3294] | 57 | USE timing ! preformance summary |
---|
[6140] | 58 | USE diurnal_bulk ! diurnal warm layer |
---|
| 59 | USE cool_skin ! Cool skin |
---|
[3294] | 60 | USE wrk_nemo ! working array |
---|
[2528] | 61 | |
---|
[3] | 62 | IMPLICIT NONE |
---|
| 63 | PRIVATE |
---|
| 64 | |
---|
[2528] | 65 | PUBLIC dia_wri ! routines called by step.F90 |
---|
| 66 | PUBLIC dia_wri_state |
---|
[2715] | 67 | PUBLIC dia_wri_alloc ! Called by nemogcm module |
---|
[3] | 68 | |
---|
[2528] | 69 | INTEGER :: nid_T, nz_T, nh_T, ndim_T, ndim_hT ! grid_T file |
---|
[3609] | 70 | INTEGER :: nb_T , ndim_bT ! grid_T file |
---|
[2528] | 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) ! ??? |
---|
[2715] | 75 | INTEGER, SAVE, ALLOCATABLE, DIMENSION(:) :: ndex_hT, ndex_hU, ndex_hV |
---|
| 76 | INTEGER, SAVE, ALLOCATABLE, DIMENSION(:) :: ndex_T, ndex_U, ndex_V |
---|
[3609] | 77 | INTEGER, SAVE, ALLOCATABLE, DIMENSION(:) :: ndex_bT |
---|
[3] | 78 | |
---|
| 79 | !! * Substitutions |
---|
| 80 | # include "zdfddm_substitute.h90" |
---|
[1756] | 81 | # include "vectopt_loop_substitute.h90" |
---|
[3] | 82 | !!---------------------------------------------------------------------- |
---|
[2528] | 83 | !! NEMO/OPA 3.3 , NEMO Consortium (2010) |
---|
[5217] | 84 | !! $Id$ |
---|
[2528] | 85 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
[3] | 86 | !!---------------------------------------------------------------------- |
---|
| 87 | CONTAINS |
---|
| 88 | |
---|
[2715] | 89 | INTEGER FUNCTION dia_wri_alloc() |
---|
| 90 | !!---------------------------------------------------------------------- |
---|
| 91 | INTEGER, DIMENSION(2) :: ierr |
---|
| 92 | !!---------------------------------------------------------------------- |
---|
| 93 | ierr = 0 |
---|
| 94 | ALLOCATE( ndex_hT(jpi*jpj) , ndex_T(jpi*jpj*jpk) , & |
---|
| 95 | & ndex_hU(jpi*jpj) , ndex_U(jpi*jpj*jpk) , & |
---|
| 96 | & ndex_hV(jpi*jpj) , ndex_V(jpi*jpj*jpk) , STAT=ierr(1) ) |
---|
| 97 | ! |
---|
| 98 | dia_wri_alloc = MAXVAL(ierr) |
---|
| 99 | IF( lk_mpp ) CALL mpp_sum( dia_wri_alloc ) |
---|
| 100 | ! |
---|
| 101 | END FUNCTION dia_wri_alloc |
---|
| 102 | |
---|
[3] | 103 | !!---------------------------------------------------------------------- |
---|
| 104 | !! Default option NetCDF output file |
---|
| 105 | !!---------------------------------------------------------------------- |
---|
[6140] | 106 | #if defined key_iomput |
---|
[3] | 107 | !!---------------------------------------------------------------------- |
---|
[2528] | 108 | !! 'key_iomput' use IOM library |
---|
| 109 | !!---------------------------------------------------------------------- |
---|
[2715] | 110 | |
---|
[1561] | 111 | SUBROUTINE dia_wri( kt ) |
---|
[1482] | 112 | !!--------------------------------------------------------------------- |
---|
| 113 | !! *** ROUTINE dia_wri *** |
---|
| 114 | !! |
---|
| 115 | !! ** Purpose : Standard output of opa: dynamics and tracer fields |
---|
| 116 | !! NETCDF format is used by default |
---|
| 117 | !! |
---|
| 118 | !! ** Method : use iom_put |
---|
| 119 | !!---------------------------------------------------------------------- |
---|
[1756] | 120 | !! |
---|
[1482] | 121 | INTEGER, INTENT( in ) :: kt ! ocean time-step index |
---|
[1756] | 122 | !! |
---|
| 123 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
[5460] | 124 | INTEGER :: jkbot ! |
---|
[1756] | 125 | REAL(wp) :: zztmp, zztmpx, zztmpy ! |
---|
[3294] | 126 | !! |
---|
[4761] | 127 | REAL(wp), POINTER, DIMENSION(:,:) :: z2d ! 2D workspace |
---|
[3294] | 128 | REAL(wp), POINTER, DIMENSION(:,:,:) :: z3d ! 3D workspace |
---|
[1482] | 129 | !!---------------------------------------------------------------------- |
---|
| 130 | ! |
---|
[3294] | 131 | IF( nn_timing == 1 ) CALL timing_start('dia_wri') |
---|
| 132 | ! |
---|
[4990] | 133 | CALL wrk_alloc( jpi , jpj , z2d ) |
---|
[3294] | 134 | CALL wrk_alloc( jpi , jpj, jpk , z3d ) |
---|
[2715] | 135 | ! |
---|
[1482] | 136 | ! Output the initial state and forcings |
---|
| 137 | IF( ninist == 1 ) THEN |
---|
| 138 | CALL dia_wri_state( 'output.init', kt ) |
---|
| 139 | ninist = 0 |
---|
| 140 | ENDIF |
---|
[3] | 141 | |
---|
[6351] | 142 | ! Output of initial vertical scale factor |
---|
| 143 | CALL iom_put("e3t_0", e3t_0(:,:,:) ) |
---|
| 144 | CALL iom_put("e3u_0", e3t_0(:,:,:) ) |
---|
| 145 | CALL iom_put("e3v_0", e3t_0(:,:,:) ) |
---|
| 146 | ! |
---|
[6387] | 147 | CALL iom_put( "e3t" , e3t_n(:,:,:) ) |
---|
| 148 | CALL iom_put( "e3u" , e3u_n(:,:,:) ) |
---|
| 149 | CALL iom_put( "e3v" , e3v_n(:,:,:) ) |
---|
| 150 | CALL iom_put( "e3w" , e3w_n(:,:,:) ) |
---|
[6351] | 151 | IF( iom_use("e3tdef") ) & |
---|
[6387] | 152 | CALL iom_put( "e3tdef" , ( ( e3t_n(:,:,:) - e3t_0(:,:,:) ) / e3t_0(:,:,:) * 100 * tmask(:,:,:) ) ** 2 ) |
---|
[5461] | 153 | |
---|
| 154 | CALL iom_put( "ssh" , sshn ) ! sea surface height |
---|
[5107] | 155 | |
---|
| 156 | CALL iom_put( "toce", tsn(:,:,:,jp_tem) ) ! 3D temperature |
---|
| 157 | CALL iom_put( "sst", tsn(:,:,1,jp_tem) ) ! surface temperature |
---|
| 158 | IF ( iom_use("sbt") ) THEN |
---|
[4990] | 159 | DO jj = 1, jpj |
---|
| 160 | DO ji = 1, jpi |
---|
[5460] | 161 | jkbot = mbkt(ji,jj) |
---|
| 162 | z2d(ji,jj) = tsn(ji,jj,jkbot,jp_tem) |
---|
[4990] | 163 | END DO |
---|
[5107] | 164 | END DO |
---|
| 165 | CALL iom_put( "sbt", z2d ) ! bottom temperature |
---|
| 166 | ENDIF |
---|
| 167 | |
---|
| 168 | CALL iom_put( "soce", tsn(:,:,:,jp_sal) ) ! 3D salinity |
---|
| 169 | CALL iom_put( "sss", tsn(:,:,1,jp_sal) ) ! surface salinity |
---|
| 170 | IF ( iom_use("sbs") ) THEN |
---|
[4990] | 171 | DO jj = 1, jpj |
---|
| 172 | DO ji = 1, jpi |
---|
[5460] | 173 | jkbot = mbkt(ji,jj) |
---|
| 174 | z2d(ji,jj) = tsn(ji,jj,jkbot,jp_sal) |
---|
[4990] | 175 | END DO |
---|
[5107] | 176 | END DO |
---|
| 177 | CALL iom_put( "sbs", z2d ) ! bottom salinity |
---|
| 178 | ENDIF |
---|
[5463] | 179 | |
---|
| 180 | IF ( iom_use("taubot") ) THEN ! bottom stress |
---|
| 181 | z2d(:,:) = 0._wp |
---|
| 182 | DO jj = 2, jpjm1 |
---|
| 183 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 184 | zztmpx = ( bfrua(ji ,jj) * un(ji ,jj,mbku(ji ,jj)) & |
---|
| 185 | & + bfrua(ji-1,jj) * un(ji-1,jj,mbku(ji-1,jj)) ) |
---|
| 186 | zztmpy = ( bfrva(ji, jj) * vn(ji,jj ,mbkv(ji,jj )) & |
---|
| 187 | & + bfrva(ji,jj-1) * vn(ji,jj-1,mbkv(ji,jj-1)) ) |
---|
| 188 | z2d(ji,jj) = rau0 * SQRT( zztmpx * zztmpx + zztmpy * zztmpy ) * tmask(ji,jj,1) |
---|
| 189 | ! |
---|
| 190 | ENDDO |
---|
| 191 | ENDDO |
---|
| 192 | CALL lbc_lnk( z2d, 'T', 1. ) |
---|
| 193 | CALL iom_put( "taubot", z2d ) |
---|
| 194 | ENDIF |
---|
[5107] | 195 | |
---|
| 196 | CALL iom_put( "uoce", un(:,:,:) ) ! 3D i-current |
---|
| 197 | CALL iom_put( "ssu", un(:,:,1) ) ! surface i-current |
---|
| 198 | IF ( iom_use("sbu") ) THEN |
---|
[4990] | 199 | DO jj = 1, jpj |
---|
| 200 | DO ji = 1, jpi |
---|
[5460] | 201 | jkbot = mbku(ji,jj) |
---|
| 202 | z2d(ji,jj) = un(ji,jj,jkbot) |
---|
[4990] | 203 | END DO |
---|
[5107] | 204 | END DO |
---|
| 205 | CALL iom_put( "sbu", z2d ) ! bottom i-current |
---|
| 206 | ENDIF |
---|
| 207 | |
---|
| 208 | CALL iom_put( "voce", vn(:,:,:) ) ! 3D j-current |
---|
| 209 | CALL iom_put( "ssv", vn(:,:,1) ) ! surface j-current |
---|
| 210 | IF ( iom_use("sbv") ) THEN |
---|
[4990] | 211 | DO jj = 1, jpj |
---|
| 212 | DO ji = 1, jpi |
---|
[5460] | 213 | jkbot = mbkv(ji,jj) |
---|
| 214 | z2d(ji,jj) = vn(ji,jj,jkbot) |
---|
[4990] | 215 | END DO |
---|
[5107] | 216 | END DO |
---|
| 217 | CALL iom_put( "sbv", z2d ) ! bottom j-current |
---|
[4990] | 218 | ENDIF |
---|
[1482] | 219 | |
---|
[5461] | 220 | CALL iom_put( "woce", wn ) ! vertical velocity |
---|
| 221 | IF( iom_use('w_masstr') .OR. iom_use('w_masstr2') ) THEN ! vertical mass transport & its square value |
---|
| 222 | ! Caution: in the VVL case, it only correponds to the baroclinic mass transport. |
---|
[5836] | 223 | z2d(:,:) = rau0 * e1e2t(:,:) |
---|
[5461] | 224 | DO jk = 1, jpk |
---|
| 225 | z3d(:,:,jk) = wn(:,:,jk) * z2d(:,:) |
---|
| 226 | END DO |
---|
| 227 | CALL iom_put( "w_masstr" , z3d ) |
---|
| 228 | IF( iom_use('w_masstr2') ) CALL iom_put( "w_masstr2", z3d(:,:,:) * z3d(:,:,:) ) |
---|
| 229 | ENDIF |
---|
| 230 | |
---|
[5107] | 231 | CALL iom_put( "avt" , avt ) ! T vert. eddy diff. coef. |
---|
| 232 | CALL iom_put( "avm" , avmu ) ! T vert. eddy visc. coef. |
---|
| 233 | CALL iom_put( "avs" , fsavs(:,:,:) ) ! S vert. eddy diff. coef. (useful only with key_zdfddm) |
---|
| 234 | |
---|
[6351] | 235 | IF( iom_use('logavt') ) CALL iom_put( "logavt", LOG( MAX( 1.e-20_wp, avt (:,:,:) ) ) ) |
---|
| 236 | IF( iom_use('logavs') ) CALL iom_put( "logavs", LOG( MAX( 1.e-20_wp, fsavs(:,:,:) ) ) ) |
---|
| 237 | |
---|
[5107] | 238 | IF ( iom_use("sstgrad") .OR. iom_use("sstgrad2") ) THEN |
---|
[4990] | 239 | DO jj = 2, jpjm1 ! sst gradient |
---|
| 240 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[5836] | 241 | zztmp = tsn(ji,jj,1,jp_tem) |
---|
| 242 | 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) |
---|
| 243 | 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) |
---|
[4990] | 244 | z2d(ji,jj) = 0.25 * ( zztmpx * zztmpx + zztmpy * zztmpy ) & |
---|
| 245 | & * umask(ji,jj,1) * umask(ji-1,jj,1) * vmask(ji,jj,1) * umask(ji,jj-1,1) |
---|
| 246 | END DO |
---|
[1756] | 247 | END DO |
---|
[4990] | 248 | CALL lbc_lnk( z2d, 'T', 1. ) |
---|
| 249 | CALL iom_put( "sstgrad2", z2d ) ! square of module of sst gradient |
---|
| 250 | z2d(:,:) = SQRT( z2d(:,:) ) |
---|
| 251 | CALL iom_put( "sstgrad" , z2d ) ! module of sst gradient |
---|
| 252 | ENDIF |
---|
| 253 | |
---|
[4761] | 254 | ! clem: heat and salt content |
---|
[4990] | 255 | IF( iom_use("heatc") ) THEN |
---|
| 256 | z2d(:,:) = 0._wp |
---|
| 257 | DO jk = 1, jpkm1 |
---|
[5107] | 258 | DO jj = 1, jpj |
---|
| 259 | DO ji = 1, jpi |
---|
[6140] | 260 | z2d(ji,jj) = z2d(ji,jj) + e3t_n(ji,jj,jk) * tsn(ji,jj,jk,jp_tem) * tmask(ji,jj,jk) |
---|
[4990] | 261 | END DO |
---|
[4761] | 262 | END DO |
---|
| 263 | END DO |
---|
[4990] | 264 | CALL iom_put( "heatc", (rau0 * rcp) * z2d ) ! vertically integrated heat content (J/m2) |
---|
| 265 | ENDIF |
---|
| 266 | |
---|
| 267 | IF( iom_use("saltc") ) THEN |
---|
| 268 | z2d(:,:) = 0._wp |
---|
| 269 | DO jk = 1, jpkm1 |
---|
[5107] | 270 | DO jj = 1, jpj |
---|
| 271 | DO ji = 1, jpi |
---|
[6140] | 272 | z2d(ji,jj) = z2d(ji,jj) + e3t_n(ji,jj,jk) * tsn(ji,jj,jk,jp_sal) * tmask(ji,jj,jk) |
---|
[4990] | 273 | END DO |
---|
| 274 | END DO |
---|
| 275 | END DO |
---|
| 276 | CALL iom_put( "saltc", rau0 * z2d ) ! vertically integrated salt content (PSU*kg/m2) |
---|
| 277 | ENDIF |
---|
[4840] | 278 | ! |
---|
[4990] | 279 | IF ( iom_use("eken") ) THEN |
---|
| 280 | rke(:,:,jk) = 0._wp ! kinetic energy |
---|
| 281 | DO jk = 1, jpkm1 |
---|
| 282 | DO jj = 2, jpjm1 |
---|
| 283 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[6140] | 284 | zztmp = 1._wp / ( e1e2t(ji,jj) * e3t_n(ji,jj,jk) ) |
---|
| 285 | zztmpx = 0.5 * ( un(ji-1,jj,jk) * un(ji-1,jj,jk) * e2u(ji-1,jj) * e3u_n(ji-1,jj,jk) & |
---|
| 286 | & + un(ji ,jj,jk) * un(ji ,jj,jk) * e2u(ji ,jj) * e3u_n(ji ,jj,jk) ) & |
---|
[4990] | 287 | & * zztmp |
---|
| 288 | ! |
---|
[6140] | 289 | zztmpy = 0.5 * ( vn(ji,jj-1,jk) * vn(ji,jj-1,jk) * e1v(ji,jj-1) * e3v_n(ji,jj-1,jk) & |
---|
| 290 | & + vn(ji,jj ,jk) * vn(ji,jj ,jk) * e1v(ji,jj ) * e3v_n(ji,jj ,jk) ) & |
---|
[4990] | 291 | & * zztmp |
---|
| 292 | ! |
---|
| 293 | rke(ji,jj,jk) = 0.5_wp * ( zztmpx + zztmpy ) |
---|
| 294 | ! |
---|
| 295 | ENDDO |
---|
[4840] | 296 | ENDDO |
---|
| 297 | ENDDO |
---|
[4990] | 298 | CALL lbc_lnk( rke, 'T', 1. ) |
---|
| 299 | CALL iom_put( "eken", rke ) |
---|
| 300 | ENDIF |
---|
[6351] | 301 | ! |
---|
| 302 | CALL iom_put( "hdiv", hdivn ) ! Horizontal divergence |
---|
| 303 | ! |
---|
[7240] | 304 | IF( iom_use("u_masstr") .OR. iom_use("u_masstr_vint") .OR. iom_use("u_heattr") .OR. iom_use("u_salttr") ) THEN |
---|
[1756] | 305 | z3d(:,:,jpk) = 0.e0 |
---|
[7240] | 306 | z2d(:,:) = 0.e0 |
---|
[1756] | 307 | DO jk = 1, jpkm1 |
---|
[7240] | 308 | z3d(:,:,jk) = rau0 * un(:,:,jk) * e2u(:,:) * fse3u(:,:,jk) * umask(:,:,jk) |
---|
| 309 | z2d(:,:) = z2d(:,:) + z3d(:,:,jk) |
---|
[1756] | 310 | END DO |
---|
| 311 | CALL iom_put( "u_masstr", z3d ) ! mass transport in i-direction |
---|
[7240] | 312 | CALL iom_put( "u_masstr_vint", z2d ) ! mass transport in i-direction vertical sum |
---|
[4990] | 313 | ENDIF |
---|
| 314 | |
---|
| 315 | IF( iom_use("u_heattr") ) THEN |
---|
| 316 | z2d(:,:) = 0.e0 |
---|
| 317 | DO jk = 1, jpkm1 |
---|
| 318 | DO jj = 2, jpjm1 |
---|
| 319 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 320 | z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( tsn(ji,jj,jk,jp_tem) + tsn(ji+1,jj,jk,jp_tem) ) |
---|
| 321 | END DO |
---|
| 322 | END DO |
---|
| 323 | END DO |
---|
| 324 | CALL lbc_lnk( z2d, 'U', -1. ) |
---|
| 325 | CALL iom_put( "u_heattr", (0.5 * rcp) * z2d ) ! heat transport in i-direction |
---|
| 326 | ENDIF |
---|
[4761] | 327 | |
---|
[4990] | 328 | IF( iom_use("u_salttr") ) THEN |
---|
[1756] | 329 | z2d(:,:) = 0.e0 |
---|
| 330 | DO jk = 1, jpkm1 |
---|
| 331 | DO jj = 2, jpjm1 |
---|
| 332 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[4990] | 333 | z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( tsn(ji,jj,jk,jp_sal) + tsn(ji+1,jj,jk,jp_sal) ) |
---|
[1756] | 334 | END DO |
---|
| 335 | END DO |
---|
| 336 | END DO |
---|
| 337 | CALL lbc_lnk( z2d, 'U', -1. ) |
---|
[4990] | 338 | CALL iom_put( "u_salttr", 0.5 * z2d ) ! heat transport in i-direction |
---|
| 339 | ENDIF |
---|
[4761] | 340 | |
---|
[4990] | 341 | |
---|
| 342 | IF( iom_use("v_masstr") .OR. iom_use("v_heattr") .OR. iom_use("v_salttr") ) THEN |
---|
[4761] | 343 | z3d(:,:,jpk) = 0.e0 |
---|
[1756] | 344 | DO jk = 1, jpkm1 |
---|
[6140] | 345 | z3d(:,:,jk) = rau0 * vn(:,:,jk) * e1v(:,:) * e3v_n(:,:,jk) * vmask(:,:,jk) |
---|
[1756] | 346 | END DO |
---|
| 347 | CALL iom_put( "v_masstr", z3d ) ! mass transport in j-direction |
---|
[4990] | 348 | ENDIF |
---|
| 349 | |
---|
| 350 | IF( iom_use("v_heattr") ) THEN |
---|
| 351 | z2d(:,:) = 0.e0 |
---|
| 352 | DO jk = 1, jpkm1 |
---|
| 353 | DO jj = 2, jpjm1 |
---|
| 354 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 355 | z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( tsn(ji,jj,jk,jp_tem) + tsn(ji,jj+1,jk,jp_tem) ) |
---|
| 356 | END DO |
---|
| 357 | END DO |
---|
| 358 | END DO |
---|
| 359 | CALL lbc_lnk( z2d, 'V', -1. ) |
---|
| 360 | CALL iom_put( "v_heattr", (0.5 * rcp) * z2d ) ! heat transport in j-direction |
---|
| 361 | ENDIF |
---|
[4761] | 362 | |
---|
[4990] | 363 | IF( iom_use("v_salttr") ) THEN |
---|
[1756] | 364 | z2d(:,:) = 0.e0 |
---|
| 365 | DO jk = 1, jpkm1 |
---|
| 366 | DO jj = 2, jpjm1 |
---|
| 367 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
[4990] | 368 | z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( tsn(ji,jj,jk,jp_sal) + tsn(ji,jj+1,jk,jp_sal) ) |
---|
[1756] | 369 | END DO |
---|
| 370 | END DO |
---|
| 371 | END DO |
---|
| 372 | CALL lbc_lnk( z2d, 'V', -1. ) |
---|
[4990] | 373 | CALL iom_put( "v_salttr", 0.5 * z2d ) ! heat transport in j-direction |
---|
[1756] | 374 | ENDIF |
---|
[7240] | 375 | |
---|
| 376 | ! Vertical integral of temperature |
---|
| 377 | IF( iom_use("tosmint") ) THEN |
---|
| 378 | 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) + rau0 * e3t_n(ji,jj,jk) * tsn(ji,jj,jk,jp_tem) |
---|
| 383 | END DO |
---|
| 384 | END DO |
---|
| 385 | END DO |
---|
| 386 | CALL lbc_lnk( z2d, 'T', -1. ) |
---|
| 387 | CALL iom_put( "tosmint", z2d ) |
---|
| 388 | ENDIF |
---|
| 389 | |
---|
| 390 | ! Vertical integral of salinity |
---|
| 391 | IF( iom_use("somint") ) THEN |
---|
| 392 | z2d(:,:)=0._wp |
---|
| 393 | DO jk = 1, jpkm1 |
---|
| 394 | DO jj = 2, jpjm1 |
---|
| 395 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
| 396 | z2d(ji,jj) = z2d(ji,jj) + rau0 * e3t_n(ji,jj,jk) * tsn(ji,jj,jk,jp_sal) |
---|
| 397 | END DO |
---|
| 398 | END DO |
---|
| 399 | END DO |
---|
| 400 | CALL lbc_lnk( z2d, 'T', -1. ) |
---|
| 401 | CALL iom_put( "somint", z2d ) |
---|
| 402 | ENDIF |
---|
| 403 | |
---|
| 404 | CALL iom_put( "bn2", rn2 ) !Brunt-Vaisala buoyancy frequency (N^2) |
---|
[2528] | 405 | ! |
---|
[4990] | 406 | CALL wrk_dealloc( jpi , jpj , z2d ) |
---|
[3294] | 407 | CALL wrk_dealloc( jpi , jpj, jpk , z3d ) |
---|
[2715] | 408 | ! |
---|
[6140] | 409 | ! If we want tmb values |
---|
| 410 | |
---|
| 411 | IF (ln_diatmb) THEN |
---|
| 412 | CALL dia_tmb |
---|
| 413 | ENDIF |
---|
| 414 | IF (ln_dia25h) THEN |
---|
| 415 | CALL dia_25h( kt ) |
---|
| 416 | ENDIF |
---|
| 417 | |
---|
[3294] | 418 | IF( nn_timing == 1 ) CALL timing_stop('dia_wri') |
---|
| 419 | ! |
---|
[1482] | 420 | END SUBROUTINE dia_wri |
---|
| 421 | |
---|
| 422 | #else |
---|
[2528] | 423 | !!---------------------------------------------------------------------- |
---|
| 424 | !! Default option use IOIPSL library |
---|
| 425 | !!---------------------------------------------------------------------- |
---|
| 426 | |
---|
[1561] | 427 | SUBROUTINE dia_wri( kt ) |
---|
[3] | 428 | !!--------------------------------------------------------------------- |
---|
| 429 | !! *** ROUTINE dia_wri *** |
---|
| 430 | !! |
---|
| 431 | !! ** Purpose : Standard output of opa: dynamics and tracer fields |
---|
| 432 | !! NETCDF format is used by default |
---|
| 433 | !! |
---|
| 434 | !! ** Method : At the beginning of the first time step (nit000), |
---|
| 435 | !! define all the NETCDF files and fields |
---|
| 436 | !! At each time step call histdef to compute the mean if ncessary |
---|
| 437 | !! Each nwrite time step, output the instantaneous or mean fields |
---|
| 438 | !!---------------------------------------------------------------------- |
---|
[5836] | 439 | INTEGER, INTENT( in ) :: kt ! ocean time-step index |
---|
| 440 | ! |
---|
[2528] | 441 | LOGICAL :: ll_print = .FALSE. ! =T print and flush numout |
---|
| 442 | CHARACTER (len=40) :: clhstnam, clop, clmx ! local names |
---|
| 443 | INTEGER :: inum = 11 ! temporary logical unit |
---|
[3294] | 444 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
| 445 | INTEGER :: ierr ! error code return from allocation |
---|
[2528] | 446 | INTEGER :: iimi, iima, ipk, it, itmod, ijmi, ijma ! local integers |
---|
[3609] | 447 | INTEGER :: jn, ierror ! local integers |
---|
[6140] | 448 | REAL(wp) :: zsto, zout, zmax, zjulian ! local scalars |
---|
[5836] | 449 | ! |
---|
[3294] | 450 | REAL(wp), POINTER, DIMENSION(:,:) :: zw2d ! 2D workspace |
---|
| 451 | REAL(wp), POINTER, DIMENSION(:,:,:) :: zw3d ! 3D workspace |
---|
[3] | 452 | !!---------------------------------------------------------------------- |
---|
[3294] | 453 | ! |
---|
| 454 | IF( nn_timing == 1 ) CALL timing_start('dia_wri') |
---|
[1482] | 455 | ! |
---|
[6140] | 456 | CALL wrk_alloc( jpi,jpj , zw2d ) |
---|
| 457 | IF( .NOT.ln_linssh ) CALL wrk_alloc( jpi,jpj,jpk , zw3d ) |
---|
[2715] | 458 | ! |
---|
[1482] | 459 | ! Output the initial state and forcings |
---|
| 460 | IF( ninist == 1 ) THEN |
---|
| 461 | CALL dia_wri_state( 'output.init', kt ) |
---|
| 462 | ninist = 0 |
---|
| 463 | ENDIF |
---|
| 464 | ! |
---|
[3] | 465 | ! 0. Initialisation |
---|
| 466 | ! ----------------- |
---|
[632] | 467 | |
---|
[3] | 468 | ! local variable for debugging |
---|
| 469 | ll_print = .FALSE. |
---|
| 470 | ll_print = ll_print .AND. lwp |
---|
| 471 | |
---|
| 472 | ! Define frequency of output and means |
---|
[5566] | 473 | clop = "x" ! no use of the mask value (require less cpu time and otherwise the model crashes) |
---|
[3] | 474 | #if defined key_diainstant |
---|
[6140] | 475 | zsto = nwrite * rdt |
---|
[1312] | 476 | clop = "inst("//TRIM(clop)//")" |
---|
[3] | 477 | #else |
---|
[6140] | 478 | zsto=rdt |
---|
[1312] | 479 | clop = "ave("//TRIM(clop)//")" |
---|
[3] | 480 | #endif |
---|
[6140] | 481 | zout = nwrite * rdt |
---|
| 482 | zmax = ( nitend - nit000 + 1 ) * rdt |
---|
[3] | 483 | |
---|
| 484 | ! Define indices of the horizontal output zoom and vertical limit storage |
---|
| 485 | iimi = 1 ; iima = jpi |
---|
| 486 | ijmi = 1 ; ijma = jpj |
---|
| 487 | ipk = jpk |
---|
| 488 | |
---|
| 489 | ! define time axis |
---|
[1334] | 490 | it = kt |
---|
| 491 | itmod = kt - nit000 + 1 |
---|
[3] | 492 | |
---|
| 493 | |
---|
| 494 | ! 1. Define NETCDF files and fields at beginning of first time step |
---|
| 495 | ! ----------------------------------------------------------------- |
---|
| 496 | |
---|
| 497 | IF( kt == nit000 ) THEN |
---|
| 498 | |
---|
| 499 | ! Define the NETCDF files (one per grid) |
---|
[632] | 500 | |
---|
[3] | 501 | ! Compute julian date from starting date of the run |
---|
[1309] | 502 | CALL ymds2ju( nyear, nmonth, nday, rdt, zjulian ) |
---|
| 503 | zjulian = zjulian - adatrj ! set calendar origin to the beginning of the experiment |
---|
[3] | 504 | IF(lwp)WRITE(numout,*) |
---|
| 505 | IF(lwp)WRITE(numout,*) 'Date 0 used :', nit000, ' YEAR ', nyear, & |
---|
| 506 | & ' MONTH ', nmonth, ' DAY ', nday, 'Julian day : ', zjulian |
---|
| 507 | IF(lwp)WRITE(numout,*) ' indexes of zoom = ', iimi, iima, ijmi, ijma, & |
---|
| 508 | ' limit storage in depth = ', ipk |
---|
| 509 | |
---|
| 510 | ! WRITE root name in date.file for use by postpro |
---|
[1581] | 511 | IF(lwp) THEN |
---|
[895] | 512 | CALL dia_nam( clhstnam, nwrite,' ' ) |
---|
[1581] | 513 | CALL ctl_opn( inum, 'date.file', 'REPLACE', 'FORMATTED', 'SEQUENTIAL', -1, numout, lwp, narea ) |
---|
[895] | 514 | WRITE(inum,*) clhstnam |
---|
| 515 | CLOSE(inum) |
---|
| 516 | ENDIF |
---|
[632] | 517 | |
---|
[3] | 518 | ! Define the T grid FILE ( nid_T ) |
---|
[632] | 519 | |
---|
[3] | 520 | CALL dia_nam( clhstnam, nwrite, 'grid_T' ) |
---|
| 521 | IF(lwp) WRITE(numout,*) " Name of NETCDF file ", clhstnam ! filename |
---|
| 522 | CALL histbeg( clhstnam, jpi, glamt, jpj, gphit, & ! Horizontal grid: glamt and gphit |
---|
| 523 | & iimi, iima-iimi+1, ijmi, ijma-ijmi+1, & |
---|
[6140] | 524 | & nit000-1, zjulian, rdt, nh_T, nid_T, domain_id=nidom, snc4chunks=snc4set ) |
---|
[3] | 525 | CALL histvert( nid_T, "deptht", "Vertical T levels", & ! Vertical grid: gdept |
---|
[4292] | 526 | & "m", ipk, gdept_1d, nz_T, "down" ) |
---|
[3] | 527 | ! ! Index of ocean points |
---|
| 528 | CALL wheneq( jpi*jpj*ipk, tmask, 1, 1., ndex_T , ndim_T ) ! volume |
---|
| 529 | CALL wheneq( jpi*jpj , tmask, 1, 1., ndex_hT, ndim_hT ) ! surface |
---|
[3609] | 530 | ! |
---|
| 531 | IF( ln_icebergs ) THEN |
---|
| 532 | ! |
---|
| 533 | !! allocation cant go in dia_wri_alloc because ln_icebergs is only set after |
---|
| 534 | !! that routine is called from nemogcm, so do it here immediately before its needed |
---|
| 535 | ALLOCATE( ndex_bT(jpi*jpj*nclasses), STAT=ierror ) |
---|
| 536 | IF( lk_mpp ) CALL mpp_sum( ierror ) |
---|
| 537 | IF( ierror /= 0 ) THEN |
---|
| 538 | CALL ctl_stop('dia_wri: failed to allocate iceberg diagnostic array') |
---|
| 539 | RETURN |
---|
| 540 | ENDIF |
---|
| 541 | ! |
---|
| 542 | !! iceberg vertical coordinate is class number |
---|
| 543 | CALL histvert( nid_T, "class", "Iceberg class", & ! Vertical grid: class |
---|
| 544 | & "number", nclasses, class_num, nb_T ) |
---|
| 545 | ! |
---|
| 546 | !! each class just needs the surface index pattern |
---|
| 547 | ndim_bT = 3 |
---|
| 548 | DO jn = 1,nclasses |
---|
| 549 | ndex_bT((jn-1)*jpi*jpj+1:jn*jpi*jpj) = ndex_hT(1:jpi*jpj) |
---|
| 550 | ENDDO |
---|
| 551 | ! |
---|
| 552 | ENDIF |
---|
[3] | 553 | |
---|
| 554 | ! Define the U grid FILE ( nid_U ) |
---|
| 555 | |
---|
| 556 | CALL dia_nam( clhstnam, nwrite, 'grid_U' ) |
---|
| 557 | IF(lwp) WRITE(numout,*) " Name of NETCDF file ", clhstnam ! filename |
---|
| 558 | CALL histbeg( clhstnam, jpi, glamu, jpj, gphiu, & ! Horizontal grid: glamu and gphiu |
---|
| 559 | & iimi, iima-iimi+1, ijmi, ijma-ijmi+1, & |
---|
[6140] | 560 | & nit000-1, zjulian, rdt, nh_U, nid_U, domain_id=nidom, snc4chunks=snc4set ) |
---|
[3] | 561 | CALL histvert( nid_U, "depthu", "Vertical U levels", & ! Vertical grid: gdept |
---|
[4292] | 562 | & "m", ipk, gdept_1d, nz_U, "down" ) |
---|
[3] | 563 | ! ! Index of ocean points |
---|
| 564 | CALL wheneq( jpi*jpj*ipk, umask, 1, 1., ndex_U , ndim_U ) ! volume |
---|
| 565 | CALL wheneq( jpi*jpj , umask, 1, 1., ndex_hU, ndim_hU ) ! surface |
---|
| 566 | |
---|
| 567 | ! Define the V grid FILE ( nid_V ) |
---|
| 568 | |
---|
| 569 | CALL dia_nam( clhstnam, nwrite, 'grid_V' ) ! filename |
---|
| 570 | IF(lwp) WRITE(numout,*) " Name of NETCDF file ", clhstnam |
---|
| 571 | CALL histbeg( clhstnam, jpi, glamv, jpj, gphiv, & ! Horizontal grid: glamv and gphiv |
---|
| 572 | & iimi, iima-iimi+1, ijmi, ijma-ijmi+1, & |
---|
[6140] | 573 | & nit000-1, zjulian, rdt, nh_V, nid_V, domain_id=nidom, snc4chunks=snc4set ) |
---|
[3] | 574 | CALL histvert( nid_V, "depthv", "Vertical V levels", & ! Vertical grid : gdept |
---|
[4292] | 575 | & "m", ipk, gdept_1d, nz_V, "down" ) |
---|
[3] | 576 | ! ! Index of ocean points |
---|
| 577 | CALL wheneq( jpi*jpj*ipk, vmask, 1, 1., ndex_V , ndim_V ) ! volume |
---|
| 578 | CALL wheneq( jpi*jpj , vmask, 1, 1., ndex_hV, ndim_hV ) ! surface |
---|
| 579 | |
---|
| 580 | ! Define the W grid FILE ( nid_W ) |
---|
| 581 | |
---|
| 582 | CALL dia_nam( clhstnam, nwrite, 'grid_W' ) ! filename |
---|
| 583 | IF(lwp) WRITE(numout,*) " Name of NETCDF file ", clhstnam |
---|
| 584 | CALL histbeg( clhstnam, jpi, glamt, jpj, gphit, & ! Horizontal grid: glamt and gphit |
---|
| 585 | & iimi, iima-iimi+1, ijmi, ijma-ijmi+1, & |
---|
[6140] | 586 | & nit000-1, zjulian, rdt, nh_W, nid_W, domain_id=nidom, snc4chunks=snc4set ) |
---|
[3] | 587 | CALL histvert( nid_W, "depthw", "Vertical W levels", & ! Vertical grid: gdepw |
---|
[4292] | 588 | & "m", ipk, gdepw_1d, nz_W, "down" ) |
---|
[3] | 589 | |
---|
[632] | 590 | |
---|
[3] | 591 | ! Declare all the output fields as NETCDF variables |
---|
| 592 | |
---|
| 593 | ! !!! nid_T : 3D |
---|
| 594 | CALL histdef( nid_T, "votemper", "Temperature" , "C" , & ! tn |
---|
| 595 | & jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout ) |
---|
| 596 | CALL histdef( nid_T, "vosaline", "Salinity" , "PSU" , & ! sn |
---|
| 597 | & jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout ) |
---|
[6140] | 598 | IF( .NOT.ln_linssh ) THEN |
---|
[4292] | 599 | CALL histdef( nid_T, "vovvle3t", "Level thickness" , "m" ,& ! e3t_n |
---|
| 600 | & jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout ) |
---|
| 601 | CALL histdef( nid_T, "vovvldep", "T point depth" , "m" ,& ! e3t_n |
---|
| 602 | & jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout ) |
---|
| 603 | CALL histdef( nid_T, "vovvldef", "Squared level deformation" , "%^2" ,& ! e3t_n |
---|
| 604 | & jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout ) |
---|
| 605 | ENDIF |
---|
[3] | 606 | ! !!! nid_T : 2D |
---|
| 607 | CALL histdef( nid_T, "sosstsst", "Sea Surface temperature" , "C" , & ! sst |
---|
| 608 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 609 | CALL histdef( nid_T, "sosaline", "Sea Surface Salinity" , "PSU" , & ! sss |
---|
| 610 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
[359] | 611 | CALL histdef( nid_T, "sossheig", "Sea Surface Height" , "m" , & ! ssh |
---|
[3] | 612 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
[2528] | 613 | CALL histdef( nid_T, "sowaflup", "Net Upward Water Flux" , "Kg/m2/s", & ! (emp-rnf) |
---|
[3] | 614 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
[4570] | 615 | CALL histdef( nid_T, "sorunoff", "River runoffs" , "Kg/m2/s", & ! runoffs |
---|
| 616 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
[3625] | 617 | CALL histdef( nid_T, "sosfldow", "downward salt flux" , "PSU/m2/s", & ! sfx |
---|
[3] | 618 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
[6140] | 619 | IF( ln_linssh ) THEN |
---|
[4292] | 620 | CALL histdef( nid_T, "sosst_cd", "Concentration/Dilution term on temperature" & ! emp * tsn(:,:,1,jp_tem) |
---|
[3625] | 621 | & , "KgC/m2/s", & ! sosst_cd |
---|
[4292] | 622 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 623 | CALL histdef( nid_T, "sosss_cd", "Concentration/Dilution term on salinity" & ! emp * tsn(:,:,1,jp_sal) |
---|
[3625] | 624 | & , "KgPSU/m2/s",& ! sosss_cd |
---|
[4292] | 625 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 626 | ENDIF |
---|
[888] | 627 | CALL histdef( nid_T, "sohefldo", "Net Downward Heat Flux" , "W/m2" , & ! qns + qsr |
---|
[3] | 628 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 629 | CALL histdef( nid_T, "soshfldo", "Shortwave Radiation" , "W/m2" , & ! qsr |
---|
| 630 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
[1585] | 631 | CALL histdef( nid_T, "somixhgt", "Turbocline Depth" , "m" , & ! hmld |
---|
| 632 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
[3] | 633 | CALL histdef( nid_T, "somxl010", "Mixed Layer Depth 0.01" , "m" , & ! hmlp |
---|
| 634 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
[1037] | 635 | CALL histdef( nid_T, "soicecov", "Ice fraction" , "[0,1]" , & ! fr_i |
---|
[3] | 636 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
[1649] | 637 | CALL histdef( nid_T, "sowindsp", "wind speed at 10m" , "m/s" , & ! wndm |
---|
| 638 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
[3609] | 639 | ! |
---|
| 640 | IF( ln_icebergs ) THEN |
---|
| 641 | CALL histdef( nid_T, "calving" , "calving mass input" , "kg/s" , & |
---|
| 642 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 643 | CALL histdef( nid_T, "calving_heat" , "calving heat flux" , "XXXX" , & |
---|
| 644 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 645 | CALL histdef( nid_T, "berg_floating_melt" , "Melt rate of icebergs + bits" , "kg/m2/s", & |
---|
| 646 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 647 | CALL histdef( nid_T, "berg_stored_ice" , "Accumulated ice mass by class" , "kg" , & |
---|
| 648 | & jpi, jpj, nh_T, nclasses , 1, nclasses , nb_T , 32, clop, zsto, zout ) |
---|
| 649 | IF( ln_bergdia ) THEN |
---|
| 650 | CALL histdef( nid_T, "berg_melt" , "Melt rate of icebergs" , "kg/m2/s", & |
---|
| 651 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 652 | CALL histdef( nid_T, "berg_buoy_melt" , "Buoyancy component of iceberg melt rate" , "kg/m2/s", & |
---|
| 653 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 654 | CALL histdef( nid_T, "berg_eros_melt" , "Erosion component of iceberg melt rate" , "kg/m2/s", & |
---|
| 655 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 656 | CALL histdef( nid_T, "berg_conv_melt" , "Convective component of iceberg melt rate", "kg/m2/s", & |
---|
| 657 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 658 | CALL histdef( nid_T, "berg_virtual_area" , "Virtual coverage by icebergs" , "m2" , & |
---|
| 659 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 660 | CALL histdef( nid_T, "bits_src" , "Mass source of bergy bits" , "kg/m2/s", & |
---|
| 661 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 662 | CALL histdef( nid_T, "bits_melt" , "Melt rate of bergy bits" , "kg/m2/s", & |
---|
| 663 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 664 | CALL histdef( nid_T, "bits_mass" , "Bergy bit density field" , "kg/m2" , & |
---|
| 665 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 666 | CALL histdef( nid_T, "berg_mass" , "Iceberg density field" , "kg/m2" , & |
---|
| 667 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 668 | CALL histdef( nid_T, "berg_real_calving" , "Calving into iceberg class" , "kg/s" , & |
---|
| 669 | & jpi, jpj, nh_T, nclasses , 1, nclasses , nb_T , 32, clop, zsto, zout ) |
---|
| 670 | ENDIF |
---|
| 671 | ENDIF |
---|
| 672 | |
---|
[5407] | 673 | IF( .NOT. ln_cpl ) THEN |
---|
[4990] | 674 | CALL histdef( nid_T, "sohefldp", "Surface Heat Flux: Damping" , "W/m2" , & ! qrp |
---|
| 675 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 676 | CALL histdef( nid_T, "sowafldp", "Surface Water Flux: Damping" , "Kg/m2/s", & ! erp |
---|
| 677 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 678 | CALL histdef( nid_T, "sosafldp", "Surface salt flux: damping" , "Kg/m2/s", & ! erp * sn |
---|
| 679 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 680 | ENDIF |
---|
[3] | 681 | |
---|
[5407] | 682 | IF( ln_cpl .AND. nn_ice <= 1 ) THEN |
---|
[4990] | 683 | CALL histdef( nid_T, "sohefldp", "Surface Heat Flux: Damping" , "W/m2" , & ! qrp |
---|
| 684 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 685 | CALL histdef( nid_T, "sowafldp", "Surface Water Flux: Damping" , "Kg/m2/s", & ! erp |
---|
| 686 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 687 | CALL histdef( nid_T, "sosafldp", "Surface salt flux: Damping" , "Kg/m2/s", & ! erp * sn |
---|
| 688 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 689 | ENDIF |
---|
| 690 | |
---|
[3] | 691 | clmx ="l_max(only(x))" ! max index on a period |
---|
[5836] | 692 | ! CALL histdef( nid_T, "sobowlin", "Bowl Index" , "W-point", & ! bowl INDEX |
---|
| 693 | ! & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clmx, zsto, zout ) |
---|
[3] | 694 | #if defined key_diahth |
---|
| 695 | CALL histdef( nid_T, "sothedep", "Thermocline Depth" , "m" , & ! hth |
---|
| 696 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 697 | CALL histdef( nid_T, "so20chgt", "Depth of 20C isotherm" , "m" , & ! hd20 |
---|
| 698 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 699 | CALL histdef( nid_T, "so28chgt", "Depth of 28C isotherm" , "m" , & ! hd28 |
---|
| 700 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 701 | CALL histdef( nid_T, "sohtc300", "Heat content 300 m" , "W" , & ! htc3 |
---|
| 702 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 703 | #endif |
---|
| 704 | |
---|
[5407] | 705 | IF( ln_cpl .AND. nn_ice == 2 ) THEN |
---|
[4990] | 706 | CALL histdef( nid_T,"soicetem" , "Ice Surface Temperature" , "K" , & ! tn_ice |
---|
| 707 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 708 | CALL histdef( nid_T,"soicealb" , "Ice Albedo" , "[0,1]" , & ! alb_ice |
---|
| 709 | & jpi, jpj, nh_T, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 710 | ENDIF |
---|
[3] | 711 | |
---|
[2528] | 712 | CALL histend( nid_T, snc4chunks=snc4set ) |
---|
[3] | 713 | |
---|
| 714 | ! !!! nid_U : 3D |
---|
| 715 | CALL histdef( nid_U, "vozocrtx", "Zonal Current" , "m/s" , & ! un |
---|
| 716 | & jpi, jpj, nh_U, ipk, 1, ipk, nz_U, 32, clop, zsto, zout ) |
---|
| 717 | ! !!! nid_U : 2D |
---|
[888] | 718 | CALL histdef( nid_U, "sozotaux", "Wind Stress along i-axis" , "N/m2" , & ! utau |
---|
[3] | 719 | & jpi, jpj, nh_U, 1 , 1, 1 , - 99, 32, clop, zsto, zout ) |
---|
| 720 | |
---|
[2528] | 721 | CALL histend( nid_U, snc4chunks=snc4set ) |
---|
[3] | 722 | |
---|
| 723 | ! !!! nid_V : 3D |
---|
| 724 | CALL histdef( nid_V, "vomecrty", "Meridional Current" , "m/s" , & ! vn |
---|
| 725 | & jpi, jpj, nh_V, ipk, 1, ipk, nz_V, 32, clop, zsto, zout ) |
---|
| 726 | ! !!! nid_V : 2D |
---|
[888] | 727 | CALL histdef( nid_V, "sometauy", "Wind Stress along j-axis" , "N/m2" , & ! vtau |
---|
[3] | 728 | & jpi, jpj, nh_V, 1 , 1, 1 , - 99, 32, clop, zsto, zout ) |
---|
| 729 | |
---|
[2528] | 730 | CALL histend( nid_V, snc4chunks=snc4set ) |
---|
[3] | 731 | |
---|
| 732 | ! !!! nid_W : 3D |
---|
| 733 | CALL histdef( nid_W, "vovecrtz", "Vertical Velocity" , "m/s" , & ! wn |
---|
| 734 | & jpi, jpj, nh_W, ipk, 1, ipk, nz_W, 32, clop, zsto, zout ) |
---|
| 735 | CALL histdef( nid_W, "votkeavt", "Vertical Eddy Diffusivity" , "m2/s" , & ! avt |
---|
| 736 | & jpi, jpj, nh_W, ipk, 1, ipk, nz_W, 32, clop, zsto, zout ) |
---|
[255] | 737 | CALL histdef( nid_W, "votkeavm", "Vertical Eddy Viscosity" , "m2/s" , & ! avmu |
---|
| 738 | & jpi, jpj, nh_W, ipk, 1, ipk, nz_W, 32, clop, zsto, zout ) |
---|
| 739 | |
---|
[3] | 740 | IF( lk_zdfddm ) THEN |
---|
| 741 | CALL histdef( nid_W,"voddmavs","Salt Vertical Eddy Diffusivity" , "m2/s" , & ! avs |
---|
| 742 | & jpi, jpj, nh_W, ipk, 1, ipk, nz_W, 32, clop, zsto, zout ) |
---|
| 743 | ENDIF |
---|
| 744 | ! !!! nid_W : 2D |
---|
[2528] | 745 | CALL histend( nid_W, snc4chunks=snc4set ) |
---|
[3] | 746 | |
---|
| 747 | IF(lwp) WRITE(numout,*) |
---|
| 748 | IF(lwp) WRITE(numout,*) 'End of NetCDF Initialization' |
---|
| 749 | IF(ll_print) CALL FLUSH(numout ) |
---|
| 750 | |
---|
| 751 | ENDIF |
---|
| 752 | |
---|
| 753 | ! 2. Start writing data |
---|
| 754 | ! --------------------- |
---|
| 755 | |
---|
[4292] | 756 | ! ndex(1) est utilise ssi l'avant dernier argument est different de |
---|
[3] | 757 | ! la taille du tableau en sortie. Dans ce cas , l'avant dernier argument |
---|
| 758 | ! donne le nombre d'elements, et ndex la liste des indices a sortir |
---|
| 759 | |
---|
[1334] | 760 | IF( lwp .AND. MOD( itmod, nwrite ) == 0 ) THEN |
---|
[3] | 761 | WRITE(numout,*) 'dia_wri : write model outputs in NetCDF files at ', kt, 'time-step' |
---|
| 762 | WRITE(numout,*) '~~~~~~ ' |
---|
| 763 | ENDIF |
---|
| 764 | |
---|
[6140] | 765 | IF( .NOT.ln_linssh ) THEN |
---|
| 766 | CALL histwrite( nid_T, "votemper", it, tsn(:,:,:,jp_tem) * e3t_n(:,:,:) , ndim_T , ndex_T ) ! heat content |
---|
| 767 | CALL histwrite( nid_T, "vosaline", it, tsn(:,:,:,jp_sal) * e3t_n(:,:,:) , ndim_T , ndex_T ) ! salt content |
---|
| 768 | CALL histwrite( nid_T, "sosstsst", it, tsn(:,:,1,jp_tem) * e3t_n(:,:,1) , ndim_hT, ndex_hT ) ! sea surface heat content |
---|
| 769 | CALL histwrite( nid_T, "sosaline", it, tsn(:,:,1,jp_sal) * e3t_n(:,:,1) , ndim_hT, ndex_hT ) ! sea surface salinity content |
---|
[4292] | 770 | ELSE |
---|
| 771 | CALL histwrite( nid_T, "votemper", it, tsn(:,:,:,jp_tem) , ndim_T , ndex_T ) ! temperature |
---|
| 772 | CALL histwrite( nid_T, "vosaline", it, tsn(:,:,:,jp_sal) , ndim_T , ndex_T ) ! salinity |
---|
| 773 | CALL histwrite( nid_T, "sosstsst", it, tsn(:,:,1,jp_tem) , ndim_hT, ndex_hT ) ! sea surface temperature |
---|
| 774 | CALL histwrite( nid_T, "sosaline", it, tsn(:,:,1,jp_sal) , ndim_hT, ndex_hT ) ! sea surface salinity |
---|
| 775 | ENDIF |
---|
[6140] | 776 | IF( .NOT.ln_linssh ) THEN |
---|
| 777 | zw3d(:,:,:) = ( ( e3t_n(:,:,:) - e3t_0(:,:,:) ) / e3t_0(:,:,:) * 100 * tmask(:,:,:) ) ** 2 |
---|
| 778 | CALL histwrite( nid_T, "vovvle3t", it, e3t_n (:,:,:) , ndim_T , ndex_T ) ! level thickness |
---|
| 779 | CALL histwrite( nid_T, "vovvldep", it, gdept_n(:,:,:) , ndim_T , ndex_T ) ! t-point depth |
---|
[4292] | 780 | CALL histwrite( nid_T, "vovvldef", it, zw3d , ndim_T , ndex_T ) ! level thickness deformation |
---|
| 781 | ENDIF |
---|
[3] | 782 | CALL histwrite( nid_T, "sossheig", it, sshn , ndim_hT, ndex_hT ) ! sea surface height |
---|
[2528] | 783 | CALL histwrite( nid_T, "sowaflup", it, ( emp-rnf ) , ndim_hT, ndex_hT ) ! upward water flux |
---|
[4570] | 784 | CALL histwrite( nid_T, "sorunoff", it, rnf , ndim_hT, ndex_hT ) ! river runoffs |
---|
[3625] | 785 | CALL histwrite( nid_T, "sosfldow", it, sfx , ndim_hT, ndex_hT ) ! downward salt flux |
---|
| 786 | ! (includes virtual salt flux beneath ice |
---|
| 787 | ! in linear free surface case) |
---|
[6140] | 788 | IF( ln_linssh ) THEN |
---|
[4292] | 789 | zw2d(:,:) = emp (:,:) * tsn(:,:,1,jp_tem) |
---|
| 790 | CALL histwrite( nid_T, "sosst_cd", it, zw2d, ndim_hT, ndex_hT ) ! c/d term on sst |
---|
| 791 | zw2d(:,:) = emp (:,:) * tsn(:,:,1,jp_sal) |
---|
| 792 | CALL histwrite( nid_T, "sosss_cd", it, zw2d, ndim_hT, ndex_hT ) ! c/d term on sss |
---|
| 793 | ENDIF |
---|
[888] | 794 | CALL histwrite( nid_T, "sohefldo", it, qns + qsr , ndim_hT, ndex_hT ) ! total heat flux |
---|
[3] | 795 | CALL histwrite( nid_T, "soshfldo", it, qsr , ndim_hT, ndex_hT ) ! solar heat flux |
---|
[1585] | 796 | CALL histwrite( nid_T, "somixhgt", it, hmld , ndim_hT, ndex_hT ) ! turbocline depth |
---|
[3] | 797 | CALL histwrite( nid_T, "somxl010", it, hmlp , ndim_hT, ndex_hT ) ! mixed layer depth |
---|
[1037] | 798 | CALL histwrite( nid_T, "soicecov", it, fr_i , ndim_hT, ndex_hT ) ! ice fraction |
---|
[1649] | 799 | CALL histwrite( nid_T, "sowindsp", it, wndm , ndim_hT, ndex_hT ) ! wind speed |
---|
[3609] | 800 | ! |
---|
| 801 | IF( ln_icebergs ) THEN |
---|
| 802 | ! |
---|
| 803 | CALL histwrite( nid_T, "calving" , it, berg_grid%calving , ndim_hT, ndex_hT ) |
---|
| 804 | CALL histwrite( nid_T, "calving_heat" , it, berg_grid%calving_hflx , ndim_hT, ndex_hT ) |
---|
| 805 | CALL histwrite( nid_T, "berg_floating_melt" , it, berg_grid%floating_melt, ndim_hT, ndex_hT ) |
---|
| 806 | ! |
---|
| 807 | CALL histwrite( nid_T, "berg_stored_ice" , it, berg_grid%stored_ice , ndim_bT, ndex_bT ) |
---|
| 808 | ! |
---|
| 809 | IF( ln_bergdia ) THEN |
---|
| 810 | CALL histwrite( nid_T, "berg_melt" , it, berg_melt , ndim_hT, ndex_hT ) |
---|
| 811 | CALL histwrite( nid_T, "berg_buoy_melt" , it, buoy_melt , ndim_hT, ndex_hT ) |
---|
| 812 | CALL histwrite( nid_T, "berg_eros_melt" , it, eros_melt , ndim_hT, ndex_hT ) |
---|
| 813 | CALL histwrite( nid_T, "berg_conv_melt" , it, conv_melt , ndim_hT, ndex_hT ) |
---|
| 814 | CALL histwrite( nid_T, "berg_virtual_area" , it, virtual_area , ndim_hT, ndex_hT ) |
---|
| 815 | CALL histwrite( nid_T, "bits_src" , it, bits_src , ndim_hT, ndex_hT ) |
---|
| 816 | CALL histwrite( nid_T, "bits_melt" , it, bits_melt , ndim_hT, ndex_hT ) |
---|
| 817 | CALL histwrite( nid_T, "bits_mass" , it, bits_mass , ndim_hT, ndex_hT ) |
---|
| 818 | CALL histwrite( nid_T, "berg_mass" , it, berg_mass , ndim_hT, ndex_hT ) |
---|
| 819 | ! |
---|
| 820 | CALL histwrite( nid_T, "berg_real_calving" , it, real_calving , ndim_bT, ndex_bT ) |
---|
| 821 | ENDIF |
---|
| 822 | ENDIF |
---|
| 823 | |
---|
[5407] | 824 | IF( .NOT. ln_cpl ) THEN |
---|
[4990] | 825 | CALL histwrite( nid_T, "sohefldp", it, qrp , ndim_hT, ndex_hT ) ! heat flux damping |
---|
| 826 | CALL histwrite( nid_T, "sowafldp", it, erp , ndim_hT, ndex_hT ) ! freshwater flux damping |
---|
[3294] | 827 | IF( ln_ssr ) zw2d(:,:) = erp(:,:) * tsn(:,:,1,jp_sal) * tmask(:,:,1) |
---|
[4990] | 828 | CALL histwrite( nid_T, "sosafldp", it, zw2d , ndim_hT, ndex_hT ) ! salt flux damping |
---|
| 829 | ENDIF |
---|
[5407] | 830 | IF( ln_cpl .AND. nn_ice <= 1 ) THEN |
---|
[4990] | 831 | CALL histwrite( nid_T, "sohefldp", it, qrp , ndim_hT, ndex_hT ) ! heat flux damping |
---|
| 832 | CALL histwrite( nid_T, "sowafldp", it, erp , ndim_hT, ndex_hT ) ! freshwater flux damping |
---|
| 833 | IF( ln_ssr ) zw2d(:,:) = erp(:,:) * tsn(:,:,1,jp_sal) * tmask(:,:,1) |
---|
| 834 | CALL histwrite( nid_T, "sosafldp", it, zw2d , ndim_hT, ndex_hT ) ! salt flux damping |
---|
| 835 | ENDIF |
---|
| 836 | ! zw2d(:,:) = FLOAT( nmln(:,:) ) * tmask(:,:,1) |
---|
| 837 | ! CALL histwrite( nid_T, "sobowlin", it, zw2d , ndim_hT, ndex_hT ) ! ??? |
---|
[3] | 838 | |
---|
| 839 | #if defined key_diahth |
---|
| 840 | CALL histwrite( nid_T, "sothedep", it, hth , ndim_hT, ndex_hT ) ! depth of the thermocline |
---|
| 841 | CALL histwrite( nid_T, "so20chgt", it, hd20 , ndim_hT, ndex_hT ) ! depth of the 20 isotherm |
---|
| 842 | CALL histwrite( nid_T, "so28chgt", it, hd28 , ndim_hT, ndex_hT ) ! depth of the 28 isotherm |
---|
| 843 | CALL histwrite( nid_T, "sohtc300", it, htc3 , ndim_hT, ndex_hT ) ! first 300m heaat content |
---|
| 844 | #endif |
---|
[888] | 845 | |
---|
[5407] | 846 | IF( ln_cpl .AND. nn_ice == 2 ) THEN |
---|
[4990] | 847 | CALL histwrite( nid_T, "soicetem", it, tn_ice(:,:,1) , ndim_hT, ndex_hT ) ! surf. ice temperature |
---|
| 848 | CALL histwrite( nid_T, "soicealb", it, alb_ice(:,:,1), ndim_hT, ndex_hT ) ! ice albedo |
---|
| 849 | ENDIF |
---|
| 850 | |
---|
[3] | 851 | CALL histwrite( nid_U, "vozocrtx", it, un , ndim_U , ndex_U ) ! i-current |
---|
[888] | 852 | CALL histwrite( nid_U, "sozotaux", it, utau , ndim_hU, ndex_hU ) ! i-wind stress |
---|
[3] | 853 | |
---|
| 854 | CALL histwrite( nid_V, "vomecrty", it, vn , ndim_V , ndex_V ) ! j-current |
---|
[888] | 855 | CALL histwrite( nid_V, "sometauy", it, vtau , ndim_hV, ndex_hV ) ! j-wind stress |
---|
[3] | 856 | |
---|
| 857 | CALL histwrite( nid_W, "vovecrtz", it, wn , ndim_T, ndex_T ) ! vert. current |
---|
| 858 | CALL histwrite( nid_W, "votkeavt", it, avt , ndim_T, ndex_T ) ! T vert. eddy diff. coef. |
---|
[255] | 859 | CALL histwrite( nid_W, "votkeavm", it, avmu , ndim_T, ndex_T ) ! T vert. eddy visc. coef. |
---|
[3] | 860 | IF( lk_zdfddm ) THEN |
---|
| 861 | CALL histwrite( nid_W, "voddmavs", it, fsavs(:,:,:), ndim_T, ndex_T ) ! S vert. eddy diff. coef. |
---|
| 862 | ENDIF |
---|
| 863 | |
---|
[1318] | 864 | ! 3. Close all files |
---|
| 865 | ! --------------------------------------- |
---|
[1561] | 866 | IF( kt == nitend ) THEN |
---|
[3] | 867 | CALL histclo( nid_T ) |
---|
| 868 | CALL histclo( nid_U ) |
---|
| 869 | CALL histclo( nid_V ) |
---|
| 870 | CALL histclo( nid_W ) |
---|
| 871 | ENDIF |
---|
[2528] | 872 | ! |
---|
[6140] | 873 | CALL wrk_dealloc( jpi , jpj , zw2d ) |
---|
| 874 | IF( .NOT.ln_linssh ) CALL wrk_dealloc( jpi , jpj , jpk , zw3d ) |
---|
[2715] | 875 | ! |
---|
[3294] | 876 | IF( nn_timing == 1 ) CALL timing_stop('dia_wri') |
---|
| 877 | ! |
---|
[3] | 878 | END SUBROUTINE dia_wri |
---|
[1567] | 879 | #endif |
---|
| 880 | |
---|
[1334] | 881 | SUBROUTINE dia_wri_state( cdfile_name, kt ) |
---|
[3] | 882 | !!--------------------------------------------------------------------- |
---|
| 883 | !! *** ROUTINE dia_wri_state *** |
---|
| 884 | !! |
---|
| 885 | !! ** Purpose : create a NetCDF file named cdfile_name which contains |
---|
| 886 | !! the instantaneous ocean state and forcing fields. |
---|
| 887 | !! Used to find errors in the initial state or save the last |
---|
| 888 | !! ocean state in case of abnormal end of a simulation |
---|
| 889 | !! |
---|
| 890 | !! ** Method : NetCDF files using ioipsl |
---|
| 891 | !! File 'output.init.nc' is created if ninist = 1 (namelist) |
---|
| 892 | !! File 'output.abort.nc' is created in case of abnormal job end |
---|
| 893 | !!---------------------------------------------------------------------- |
---|
[1334] | 894 | CHARACTER (len=* ), INTENT( in ) :: cdfile_name ! name of the file created |
---|
| 895 | INTEGER , INTENT( in ) :: kt ! ocean time-step index |
---|
[2528] | 896 | !! |
---|
[648] | 897 | CHARACTER (len=32) :: clname |
---|
[3] | 898 | CHARACTER (len=40) :: clop |
---|
[2528] | 899 | INTEGER :: id_i , nz_i, nh_i |
---|
| 900 | INTEGER, DIMENSION(1) :: idex ! local workspace |
---|
[6140] | 901 | REAL(wp) :: zsto, zout, zmax, zjulian |
---|
[3] | 902 | !!---------------------------------------------------------------------- |
---|
[3294] | 903 | ! |
---|
[3] | 904 | ! 0. Initialisation |
---|
| 905 | ! ----------------- |
---|
[632] | 906 | |
---|
[648] | 907 | ! Define name, frequency of output and means |
---|
| 908 | clname = cdfile_name |
---|
[1792] | 909 | IF( .NOT. Agrif_Root() ) clname = TRIM(Agrif_CFixed())//'_'//TRIM(clname) |
---|
[3] | 910 | zsto = rdt |
---|
| 911 | clop = "inst(x)" ! no use of the mask value (require less cpu time) |
---|
| 912 | zout = rdt |
---|
[6140] | 913 | zmax = ( nitend - nit000 + 1 ) * rdt |
---|
[3] | 914 | |
---|
[648] | 915 | IF(lwp) WRITE(numout,*) |
---|
| 916 | IF(lwp) WRITE(numout,*) 'dia_wri_state : single instantaneous ocean state' |
---|
| 917 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~ and forcing fields file created ' |
---|
| 918 | IF(lwp) WRITE(numout,*) ' and named :', clname, '.nc' |
---|
| 919 | |
---|
| 920 | |
---|
[3] | 921 | ! 1. Define NETCDF files and fields at beginning of first time step |
---|
| 922 | ! ----------------------------------------------------------------- |
---|
| 923 | |
---|
| 924 | ! Compute julian date from starting date of the run |
---|
[1310] | 925 | CALL ymds2ju( nyear, nmonth, nday, rdt, zjulian ) ! time axis |
---|
| 926 | zjulian = zjulian - adatrj ! set calendar origin to the beginning of the experiment |
---|
[648] | 927 | CALL histbeg( clname, jpi, glamt, jpj, gphit, & |
---|
[6140] | 928 | 1, jpi, 1, jpj, nit000-1, zjulian, rdt, nh_i, id_i, domain_id=nidom, snc4chunks=snc4set ) ! Horizontal grid : glamt and gphit |
---|
[3] | 929 | CALL histvert( id_i, "deptht", "Vertical T levels", & ! Vertical grid : gdept |
---|
[4292] | 930 | "m", jpk, gdept_1d, nz_i, "down") |
---|
[3] | 931 | |
---|
| 932 | ! Declare all the output fields as NetCDF variables |
---|
| 933 | |
---|
| 934 | CALL histdef( id_i, "vosaline", "Salinity" , "PSU" , & ! salinity |
---|
| 935 | & jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout ) |
---|
| 936 | CALL histdef( id_i, "votemper", "Temperature" , "C" , & ! temperature |
---|
| 937 | & jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout ) |
---|
[359] | 938 | CALL histdef( id_i, "sossheig", "Sea Surface Height" , "m" , & ! ssh |
---|
| 939 | & jpi, jpj, nh_i, 1 , 1, 1 , nz_i, 32, clop, zsto, zout ) |
---|
[3] | 940 | CALL histdef( id_i, "vozocrtx", "Zonal Current" , "m/s" , & ! zonal current |
---|
| 941 | & jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout ) |
---|
| 942 | CALL histdef( id_i, "vomecrty", "Meridional Current" , "m/s" , & ! meridonal current |
---|
| 943 | & jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout ) |
---|
| 944 | CALL histdef( id_i, "vovecrtz", "Vertical Velocity" , "m/s" , & ! vertical current |
---|
| 945 | & jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout ) |
---|
[6140] | 946 | ! |
---|
| 947 | CALL histdef( id_i, "ahtu" , "u-eddy diffusivity" , "m2/s" , & ! zonal current |
---|
| 948 | & jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout ) |
---|
| 949 | CALL histdef( id_i, "ahtv" , "v-eddy diffusivity" , "m2/s" , & ! meridonal current |
---|
| 950 | & jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout ) |
---|
| 951 | CALL histdef( id_i, "ahmt" , "t-eddy viscosity" , "m2/s" , & ! zonal current |
---|
| 952 | & jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout ) |
---|
| 953 | CALL histdef( id_i, "ahmf" , "f-eddy viscosity" , "m2/s" , & ! meridonal current |
---|
| 954 | & jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout ) |
---|
| 955 | ! |
---|
[3] | 956 | CALL histdef( id_i, "sowaflup", "Net Upward Water Flux" , "Kg/m2/S", & ! net freshwater |
---|
| 957 | & jpi, jpj, nh_i, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 958 | CALL histdef( id_i, "sohefldo", "Net Downward Heat Flux", "W/m2" , & ! net heat flux |
---|
| 959 | & jpi, jpj, nh_i, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 960 | CALL histdef( id_i, "soshfldo", "Shortwave Radiation" , "W/m2" , & ! solar flux |
---|
| 961 | & jpi, jpj, nh_i, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
[1037] | 962 | CALL histdef( id_i, "soicecov", "Ice fraction" , "[0,1]" , & ! fr_i |
---|
[3] | 963 | & jpi, jpj, nh_i, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 964 | CALL histdef( id_i, "sozotaux", "Zonal Wind Stress" , "N/m2" , & ! i-wind stress |
---|
| 965 | & jpi, jpj, nh_i, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
| 966 | CALL histdef( id_i, "sometauy", "Meridional Wind Stress", "N/m2" , & ! j-wind stress |
---|
| 967 | & jpi, jpj, nh_i, 1 , 1, 1 , -99 , 32, clop, zsto, zout ) |
---|
[6140] | 968 | IF( .NOT.ln_linssh ) THEN |
---|
| 969 | CALL histdef( id_i, "vovvldep", "T point depth" , "m" , & ! t-point depth |
---|
[4292] | 970 | & jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout ) |
---|
[6140] | 971 | CALL histdef( id_i, "vovvle3t", "T point thickness" , "m" , & ! t-point depth |
---|
| 972 | & jpi, jpj, nh_i, jpk, 1, jpk, nz_i, 32, clop, zsto, zout ) |
---|
[5836] | 973 | ENDIF |
---|
[3] | 974 | |
---|
[1482] | 975 | #if defined key_lim2 |
---|
| 976 | CALL lim_wri_state_2( kt, id_i, nh_i ) |
---|
[4161] | 977 | #elif defined key_lim3 |
---|
| 978 | CALL lim_wri_state( kt, id_i, nh_i ) |
---|
[1482] | 979 | #else |
---|
[2528] | 980 | CALL histend( id_i, snc4chunks=snc4set ) |
---|
[1482] | 981 | #endif |
---|
[3] | 982 | |
---|
| 983 | ! 2. Start writing data |
---|
| 984 | ! --------------------- |
---|
| 985 | ! idex(1) est utilise ssi l'avant dernier argument est diffferent de |
---|
| 986 | ! la taille du tableau en sortie. Dans ce cas , l'avant dernier argument |
---|
| 987 | ! donne le nombre d'elements, et idex la liste des indices a sortir |
---|
| 988 | idex(1) = 1 ! init to avoid compil warning |
---|
[632] | 989 | |
---|
[3] | 990 | ! Write all fields on T grid |
---|
[3294] | 991 | CALL histwrite( id_i, "votemper", kt, tsn(:,:,:,jp_tem), jpi*jpj*jpk, idex ) ! now temperature |
---|
| 992 | CALL histwrite( id_i, "vosaline", kt, tsn(:,:,:,jp_sal), jpi*jpj*jpk, idex ) ! now salinity |
---|
| 993 | CALL histwrite( id_i, "sossheig", kt, sshn , jpi*jpj , idex ) ! sea surface height |
---|
| 994 | CALL histwrite( id_i, "vozocrtx", kt, un , jpi*jpj*jpk, idex ) ! now i-velocity |
---|
| 995 | CALL histwrite( id_i, "vomecrty", kt, vn , jpi*jpj*jpk, idex ) ! now j-velocity |
---|
| 996 | CALL histwrite( id_i, "vovecrtz", kt, wn , jpi*jpj*jpk, idex ) ! now k-velocity |
---|
[6140] | 997 | ! |
---|
| 998 | CALL histwrite( id_i, "ahtu" , kt, ahtu , jpi*jpj*jpk, idex ) ! aht at u-point |
---|
| 999 | CALL histwrite( id_i, "ahtv" , kt, ahtv , jpi*jpj*jpk, idex ) ! - at v-point |
---|
| 1000 | CALL histwrite( id_i, "ahmt" , kt, ahmt , jpi*jpj*jpk, idex ) ! ahm at t-point |
---|
| 1001 | CALL histwrite( id_i, "ahmf" , kt, ahmf , jpi*jpj*jpk, idex ) ! - at f-point |
---|
| 1002 | ! |
---|
[5836] | 1003 | CALL histwrite( id_i, "sowaflup", kt, emp-rnf , jpi*jpj , idex ) ! freshwater budget |
---|
[3294] | 1004 | CALL histwrite( id_i, "sohefldo", kt, qsr + qns , jpi*jpj , idex ) ! total heat flux |
---|
| 1005 | CALL histwrite( id_i, "soshfldo", kt, qsr , jpi*jpj , idex ) ! solar heat flux |
---|
| 1006 | CALL histwrite( id_i, "soicecov", kt, fr_i , jpi*jpj , idex ) ! ice fraction |
---|
| 1007 | CALL histwrite( id_i, "sozotaux", kt, utau , jpi*jpj , idex ) ! i-wind stress |
---|
| 1008 | CALL histwrite( id_i, "sometauy", kt, vtau , jpi*jpj , idex ) ! j-wind stress |
---|
[3] | 1009 | |
---|
[6140] | 1010 | IF( .NOT.ln_linssh ) THEN |
---|
| 1011 | CALL histwrite( id_i, "vovvldep", kt, gdept_n(:,:,:), jpi*jpj*jpk, idex )! T-cell depth |
---|
| 1012 | CALL histwrite( id_i, "vovvle3t", kt, e3t_n (:,:,:) , jpi*jpj*jpk, idex )! T-cell thickness |
---|
| 1013 | END IF |
---|
[3] | 1014 | ! 3. Close the file |
---|
| 1015 | ! ----------------- |
---|
| 1016 | CALL histclo( id_i ) |
---|
[6140] | 1017 | #if ! defined key_iomput |
---|
[1561] | 1018 | IF( ninist /= 1 ) THEN |
---|
| 1019 | CALL histclo( nid_T ) |
---|
| 1020 | CALL histclo( nid_U ) |
---|
| 1021 | CALL histclo( nid_V ) |
---|
| 1022 | CALL histclo( nid_W ) |
---|
| 1023 | ENDIF |
---|
| 1024 | #endif |
---|
[3294] | 1025 | ! |
---|
[3] | 1026 | END SUBROUTINE dia_wri_state |
---|
[6140] | 1027 | |
---|
[3] | 1028 | !!====================================================================== |
---|
| 1029 | END MODULE diawri |
---|