New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
diawri.F90 in NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/DIA – NEMO

source: NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/DIA/diawri.F90 @ 12353

Last change on this file since 12353 was 12353, checked in by acc, 4 years ago

Branch 2019/dev_r11943_MERGE_2019. Additions to the do loop macro implementation: converted a few loops previously missed because they used jpi-1 instead of jpim1 etc.; changed internal macro names in do_loop_substitute.h90 to strings that are much more unlikely to appear in any future code elsewhere and removed the key_vectopt_loop option (and all related code) since the do loop macros have suppressed this option. These changes have been fully SETTE-tested and this branch should now be ready to go back to the trunk.

  • Property svn:keywords set to Id
File size: 55.2 KB
Line 
1MODULE 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 isf_oce
29   USE isfcpl
30   USE abl            ! abl variables in case ln_abl = .true.
31   USE dom_oce        ! ocean space and time domain
32   USE phycst         ! physical constants
33   USE dianam         ! build name of file (routine)
34   USE diahth         ! thermocline diagnostics
35   USE dynadv   , ONLY: ln_dynadv_vec
36   USE icb_oce        ! Icebergs
37   USE icbdia         ! Iceberg budgets
38   USE ldftra         ! lateral physics: eddy diffusivity coef.
39   USE ldfdyn         ! lateral physics: eddy viscosity   coef.
40   USE sbc_oce        ! Surface boundary condition: ocean fields
41   USE sbc_ice        ! Surface boundary condition: ice fields
42   USE sbcssr         ! restoring term toward SST/SSS climatology
43   USE sbcwave        ! wave parameters
44   USE wet_dry        ! wetting and drying
45   USE zdf_oce        ! ocean vertical physics
46   USE zdfdrg         ! ocean vertical physics: top/bottom friction
47   USE zdfmxl         ! mixed layer
48   !
49   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
50   USE in_out_manager ! I/O manager
51   USE dia25h         ! 25h Mean output
52   USE iom            !
53   USE ioipsl         !
54
55#if defined key_si3
56   USE ice 
57   USE icewri 
58#endif
59   USE lib_mpp         ! MPP library
60   USE timing          ! preformance summary
61   USE diu_bulk        ! diurnal warm layer
62   USE diu_coolskin    ! Cool skin
63
64   IMPLICIT NONE
65   PRIVATE
66
67   PUBLIC   dia_wri                 ! routines called by step.F90
68   PUBLIC   dia_wri_state
69   PUBLIC   dia_wri_alloc           ! Called by nemogcm module
70#if ! defined key_iomput   
71   PUBLIC   dia_wri_alloc_abl       ! Called by sbcabl  module (if ln_abl = .true.)
72#endif
73   INTEGER ::   nid_T, nz_T, nh_T, ndim_T, ndim_hT   ! grid_T file
74   INTEGER ::          nb_T              , ndim_bT   ! grid_T file
75   INTEGER ::   nid_U, nz_U, nh_U, ndim_U, ndim_hU   ! grid_U file
76   INTEGER ::   nid_V, nz_V, nh_V, ndim_V, ndim_hV   ! grid_V file
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   
79   INTEGER ::   ndex(1)                              ! ???
80   INTEGER, SAVE, ALLOCATABLE, DIMENSION(:) :: ndex_hT, ndex_hU, ndex_hV
81   INTEGER, SAVE, ALLOCATABLE, DIMENSION(:) :: ndex_hA, ndex_A ! ABL
82   INTEGER, SAVE, ALLOCATABLE, DIMENSION(:) :: ndex_T, ndex_U, ndex_V
83   INTEGER, SAVE, ALLOCATABLE, DIMENSION(:) :: ndex_bT
84
85   !! * Substitutions
86#  include "do_loop_substitute.h90"
87   !!----------------------------------------------------------------------
88   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
89   !! $Id$
90   !! Software governed by the CeCILL license (see ./LICENSE)
91   !!----------------------------------------------------------------------
92CONTAINS
93
94#if defined key_iomput
95   !!----------------------------------------------------------------------
96   !!   'key_iomput'                                        use IOM library
97   !!----------------------------------------------------------------------
98   INTEGER FUNCTION dia_wri_alloc()
99      !
100      dia_wri_alloc = 0
101      !
102   END FUNCTION dia_wri_alloc
103
104   
105   SUBROUTINE dia_wri( kt, Kmm )
106      !!---------------------------------------------------------------------
107      !!                  ***  ROUTINE dia_wri  ***
108      !!                   
109      !! ** Purpose :   Standard output of opa: dynamics and tracer fields
110      !!      NETCDF format is used by default
111      !!
112      !! ** Method  :  use iom_put
113      !!----------------------------------------------------------------------
114      INTEGER, INTENT( in ) ::   kt      ! ocean time-step index
115      INTEGER, INTENT( in ) ::   Kmm     ! ocean time level index
116      !!
117      INTEGER ::   ji, jj, jk       ! dummy loop indices
118      INTEGER ::   ikbot            ! local integer
119      REAL(wp)::   zztmp , zztmpx   ! local scalar
120      REAL(wp)::   zztmp2, zztmpy   !   -      -
121      REAL(wp), DIMENSION(jpi,jpj)     ::   z2d   ! 2D workspace
122      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   z3d   ! 3D workspace
123      !!----------------------------------------------------------------------
124      !
125      IF( ln_timing )   CALL timing_start('dia_wri')
126      !
127      ! Output the initial state and forcings
128      IF( ninist == 1 ) THEN                       
129         CALL dia_wri_state( Kmm, 'output.init' )
130         ninist = 0
131      ENDIF
132
133      ! Output of initial vertical scale factor
134      CALL iom_put("e3t_0", e3t_0(:,:,:) )
135      CALL iom_put("e3u_0", e3u_0(:,:,:) )
136      CALL iom_put("e3v_0", e3v_0(:,:,:) )
137      !
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) )
142      IF( iom_use("e3tdef") )   &
143         CALL iom_put( "e3tdef"  , ( ( e3t(:,:,:,Kmm) - e3t_0(:,:,:) ) / e3t_0(:,:,:) * 100 * tmask(:,:,:) ) ** 2 )
144
145      IF( ll_wd ) THEN
146         CALL iom_put( "ssh" , (ssh(:,:,Kmm)+ssh_ref)*tmask(:,:,1) )   ! sea surface height (brought back to the reference used for wetting and drying)
147      ELSE
148         CALL iom_put( "ssh" , ssh(:,:,Kmm) )              ! sea surface height
149      ENDIF
150
151      IF( iom_use("wetdep") )   &                  ! wet depth
152         CALL iom_put( "wetdep" , ht_0(:,:) + ssh(:,:,Kmm) )
153     
154      CALL iom_put( "toce", ts(:,:,:,jp_tem,Kmm) )    ! 3D temperature
155      CALL iom_put(  "sst", ts(:,:,1,jp_tem,Kmm) )    ! surface temperature
156      IF ( iom_use("sbt") ) THEN
157         DO_2D_11_11
158            ikbot = mbkt(ji,jj)
159            z2d(ji,jj) = ts(ji,jj,ikbot,jp_tem,Kmm)
160         END_2D
161         CALL iom_put( "sbt", z2d )                ! bottom temperature
162      ENDIF
163     
164      CALL iom_put( "soce", ts(:,:,:,jp_sal,Kmm) )    ! 3D salinity
165      CALL iom_put(  "sss", ts(:,:,1,jp_sal,Kmm) )    ! surface salinity
166      IF ( iom_use("sbs") ) THEN
167         DO_2D_11_11
168            ikbot = mbkt(ji,jj)
169            z2d(ji,jj) = ts(ji,jj,ikbot,jp_sal,Kmm)
170         END_2D
171         CALL iom_put( "sbs", z2d )                ! bottom salinity
172      ENDIF
173
174      IF ( iom_use("taubot") ) THEN                ! bottom stress
175         zztmp = rau0 * 0.25
176         z2d(:,:) = 0._wp
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
185         CALL lbc_lnk( 'diawri', z2d, 'T', 1. )
186         CALL iom_put( "taubot", z2d )           
187      ENDIF
188         
189      CALL iom_put( "uoce", uu(:,:,:,Kmm) )            ! 3D i-current
190      CALL iom_put(  "ssu", uu(:,:,1,Kmm) )            ! surface i-current
191      IF ( iom_use("sbu") ) THEN
192         DO_2D_11_11
193            ikbot = mbku(ji,jj)
194            z2d(ji,jj) = uu(ji,jj,ikbot,Kmm)
195         END_2D
196         CALL iom_put( "sbu", z2d )                ! bottom i-current
197      ENDIF
198     
199      CALL iom_put( "voce", vv(:,:,:,Kmm) )            ! 3D j-current
200      CALL iom_put(  "ssv", vv(:,:,1,Kmm) )            ! surface j-current
201      IF ( iom_use("sbv") ) THEN
202         DO_2D_11_11
203            ikbot = mbkv(ji,jj)
204            z2d(ji,jj) = vv(ji,jj,ikbot,Kmm)
205         END_2D
206         CALL iom_put( "sbv", z2d )                ! bottom j-current
207      ENDIF
208
209      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
212      IF( iom_use('w_masstr') .OR. iom_use('w_masstr2') ) THEN   ! vertical mass transport & its square value
213         ! Caution: in the VVL case, it only correponds to the baroclinic mass transport.
214         z2d(:,:) = rau0 * e1e2t(:,:)
215         DO jk = 1, jpk
216            z3d(:,:,jk) = ww(:,:,jk) * z2d(:,:)
217         END DO
218         CALL iom_put( "w_masstr" , z3d ) 
219         IF( iom_use('w_masstr2') )   CALL iom_put( "w_masstr2", z3d(:,:,:) * z3d(:,:,:) )
220      ENDIF
221      !
222      IF( ln_zad_Aimp ) ww = ww - wi               ! Remove implicit part of vertical velocity that was added for diagnostic output
223
224      CALL iom_put( "avt" , avt )                  ! T vert. eddy diff. coef.
225      CALL iom_put( "avs" , avs )                  ! S vert. eddy diff. coef.
226      CALL iom_put( "avm" , avm )                  ! T vert. eddy visc. coef.
227
228      IF( iom_use('logavt') )   CALL iom_put( "logavt", LOG( MAX( 1.e-20_wp, avt(:,:,:) ) ) )
229      IF( iom_use('logavs') )   CALL iom_put( "logavs", LOG( MAX( 1.e-20_wp, avs(:,:,:) ) ) )
230
231      IF ( iom_use("sstgrad") .OR. iom_use("sstgrad2") ) THEN
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
239         CALL lbc_lnk( 'diawri', z2d, 'T', 1. )
240         CALL iom_put( "sstgrad2",  z2d )          ! square of module of sst gradient
241         z2d(:,:) = SQRT( z2d(:,:) )
242         CALL iom_put( "sstgrad" ,  z2d )          ! module of sst gradient
243      ENDIF
244         
245      ! heat and salt contents
246      IF( iom_use("heatc") ) THEN
247         z2d(:,:)  = 0._wp 
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
251         CALL iom_put( "heatc", rau0_rcp * z2d )   ! vertically integrated heat content (J/m2)
252      ENDIF
253
254      IF( iom_use("saltc") ) THEN
255         z2d(:,:)  = 0._wp 
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
259         CALL iom_put( "saltc", rau0 * z2d )          ! vertically integrated salt content (PSU*kg/m2)
260      ENDIF
261      !
262      IF ( iom_use("eken") ) THEN
263         z3d(:,:,jpk) = 0._wp 
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
271         CALL lbc_lnk( 'diawri', z3d, 'T', 1. )
272         CALL iom_put( "eken", z3d )                 ! kinetic energy
273      ENDIF
274      !
275      CALL iom_put( "hdiv", hdiv )                  ! Horizontal divergence
276      !
277      IF( iom_use("u_masstr") .OR. iom_use("u_masstr_vint") .OR. iom_use("u_heattr") .OR. iom_use("u_salttr") ) THEN
278         z3d(:,:,jpk) = 0.e0
279         z2d(:,:) = 0.e0
280         DO jk = 1, jpkm1
281            z3d(:,:,jk) = rau0 * uu(:,:,jk,Kmm) * e2u(:,:) * e3u(:,:,jk,Kmm) * umask(:,:,jk)
282            z2d(:,:) = z2d(:,:) + z3d(:,:,jk)
283         END DO
284         CALL iom_put( "u_masstr"     , z3d )         ! mass transport in i-direction
285         CALL iom_put( "u_masstr_vint", z2d )         ! mass transport in i-direction vertical sum
286      ENDIF
287     
288      IF( iom_use("u_heattr") ) THEN
289         z2d(:,:) = 0._wp 
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
293         CALL lbc_lnk( 'diawri', z2d, 'U', -1. )
294         CALL iom_put( "u_heattr", 0.5*rcp * z2d )    ! heat transport in i-direction
295      ENDIF
296
297      IF( iom_use("u_salttr") ) THEN
298         z2d(:,:) = 0.e0 
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
302         CALL lbc_lnk( 'diawri', z2d, 'U', -1. )
303         CALL iom_put( "u_salttr", 0.5 * z2d )        ! heat transport in i-direction
304      ENDIF
305
306     
307      IF( iom_use("v_masstr") .OR. iom_use("v_heattr") .OR. iom_use("v_salttr") ) THEN
308         z3d(:,:,jpk) = 0.e0
309         DO jk = 1, jpkm1
310            z3d(:,:,jk) = rau0 * vv(:,:,jk,Kmm) * e1v(:,:) * e3v(:,:,jk,Kmm) * vmask(:,:,jk)
311         END DO
312         CALL iom_put( "v_masstr", z3d )              ! mass transport in j-direction
313      ENDIF
314     
315      IF( iom_use("v_heattr") ) THEN
316         z2d(:,:) = 0.e0 
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
320         CALL lbc_lnk( 'diawri', z2d, 'V', -1. )
321         CALL iom_put( "v_heattr", 0.5*rcp * z2d )    !  heat transport in j-direction
322      ENDIF
323
324      IF( iom_use("v_salttr") ) THEN
325         z2d(:,:) = 0._wp 
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
329         CALL lbc_lnk( 'diawri', z2d, 'V', -1. )
330         CALL iom_put( "v_salttr", 0.5 * z2d )        !  heat transport in j-direction
331      ENDIF
332
333      IF( iom_use("tosmint") ) THEN
334         z2d(:,:) = 0._wp
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
338         CALL lbc_lnk( 'diawri', z2d, 'T', -1. )
339         CALL iom_put( "tosmint", rau0 * z2d )        ! Vertical integral of temperature
340      ENDIF
341      IF( iom_use("somint") ) THEN
342         z2d(:,:)=0._wp
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
346         CALL lbc_lnk( 'diawri', z2d, 'T', -1. )
347         CALL iom_put( "somint", rau0 * z2d )         ! Vertical integral of salinity
348      ENDIF
349
350      CALL iom_put( "bn2", rn2 )                      ! Brunt-Vaisala buoyancy frequency (N^2)
351      !
352     
353      IF (ln_dia25h)   CALL dia_25h( kt, Kmm )        ! 25h averaging
354
355      IF( ln_timing )   CALL timing_stop('dia_wri')
356      !
357   END SUBROUTINE dia_wri
358
359#else
360   !!----------------------------------------------------------------------
361   !!   Default option                                  use IOIPSL  library
362   !!----------------------------------------------------------------------
363
364   INTEGER FUNCTION dia_wri_alloc()
365      !!----------------------------------------------------------------------
366      INTEGER, DIMENSION(2) :: ierr
367      !!----------------------------------------------------------------------
368      ierr = 0
369      ALLOCATE( ndex_hT(jpi*jpj) , ndex_T(jpi*jpj*jpk) ,     &
370         &      ndex_hU(jpi*jpj) , ndex_U(jpi*jpj*jpk) ,     &
371         &      ndex_hV(jpi*jpj) , ndex_V(jpi*jpj*jpk) , STAT=ierr(1) )
372         !
373     dia_wri_alloc = MAXVAL(ierr)
374      CALL mpp_sum( 'diawri', dia_wri_alloc )
375      !
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
384
385   
386   SUBROUTINE dia_wri( kt, Kmm )
387      !!---------------------------------------------------------------------
388      !!                  ***  ROUTINE dia_wri  ***
389      !!                   
390      !! ** Purpose :   Standard output of opa: dynamics and tracer fields
391      !!      NETCDF format is used by default
392      !!
393      !! ** Method  :   At the beginning of the first time step (nit000),
394      !!      define all the NETCDF files and fields
395      !!      At each time step call histdef to compute the mean if ncessary
396      !!      Each nn_write time step, output the instantaneous or mean fields
397      !!----------------------------------------------------------------------
398      INTEGER, INTENT( in ) ::   kt   ! ocean time-step index
399      INTEGER, INTENT( in ) ::   Kmm  ! ocean time level index
400      !
401      LOGICAL ::   ll_print = .FALSE.                        ! =T print and flush numout
402      CHARACTER (len=40) ::   clhstnam, clop, clmx           ! local names
403      INTEGER  ::   inum = 11                                ! temporary logical unit
404      INTEGER  ::   ji, jj, jk                               ! dummy loop indices
405      INTEGER  ::   ierr                                     ! error code return from allocation
406      INTEGER  ::   iimi, iima, ipk, it, itmod, ijmi, ijma   ! local integers
407      INTEGER  ::   ipka                                     ! ABL
408      INTEGER  ::   jn, ierror                               ! local integers
409      REAL(wp) ::   zsto, zout, zmax, zjulian                ! local scalars
410      !
411      REAL(wp), DIMENSION(jpi,jpj)   :: zw2d       ! 2D workspace
412      REAL(wp), DIMENSION(jpi,jpj,jpk) :: zw3d       ! 3D workspace
413      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE :: zw3d_abl   ! ABL 3D workspace
414      !!----------------------------------------------------------------------
415      !
416      IF( ninist == 1 ) THEN     !==  Output the initial state and forcings  ==!
417         CALL dia_wri_state( Kmm, 'output.init' )
418         ninist = 0
419      ENDIF
420      !
421      IF( nn_write == -1 )   RETURN   ! we will never do any output
422      !
423      IF( ln_timing )   CALL timing_start('dia_wri')
424      !
425      ! 0. Initialisation
426      ! -----------------
427
428      ll_print = .FALSE.                  ! local variable for debugging
429      ll_print = ll_print .AND. lwp
430
431      ! Define frequency of output and means
432      clop = "x"         ! no use of the mask value (require less cpu time and otherwise the model crashes)
433#if defined key_diainstant
434      zsto = nn_write * rdt
435      clop = "inst("//TRIM(clop)//")"
436#else
437      zsto=rdt
438      clop = "ave("//TRIM(clop)//")"
439#endif
440      zout = nn_write * rdt
441      zmax = ( nitend - nit000 + 1 ) * rdt
442
443      ! Define indices of the horizontal output zoom and vertical limit storage
444      iimi = 1      ;      iima = jpi
445      ijmi = 1      ;      ijma = jpj
446      ipk = jpk
447      IF(ln_abl) ipka = jpkam1
448
449      ! define time axis
450      it = kt
451      itmod = kt - nit000 + 1
452
453
454      ! 1. Define NETCDF files and fields at beginning of first time step
455      ! -----------------------------------------------------------------
456
457      IF( kt == nit000 ) THEN
458
459         ! Define the NETCDF files (one per grid)
460
461         ! Compute julian date from starting date of the run
462         CALL ymds2ju( nyear, nmonth, nday, rdt, zjulian )
463         zjulian = zjulian - adatrj   !   set calendar origin to the beginning of the experiment
464         IF(lwp)WRITE(numout,*)
465         IF(lwp)WRITE(numout,*) 'Date 0 used :', nit000, ' YEAR ', nyear,   &
466            &                    ' MONTH ', nmonth, ' DAY ', nday, 'Julian day : ', zjulian
467         IF(lwp)WRITE(numout,*) ' indexes of zoom = ', iimi, iima, ijmi, ijma,   &
468                                 ' limit storage in depth = ', ipk
469
470         ! WRITE root name in date.file for use by postpro
471         IF(lwp) THEN
472            CALL dia_nam( clhstnam, nn_write,' ' )
473            CALL ctl_opn( inum, 'date.file', 'REPLACE', 'FORMATTED', 'SEQUENTIAL', -1, numout, lwp, narea )
474            WRITE(inum,*) clhstnam
475            CLOSE(inum)
476         ENDIF
477
478         ! Define the T grid FILE ( nid_T )
479
480         CALL dia_nam( clhstnam, nn_write, 'grid_T' )
481         IF(lwp) WRITE(numout,*) " Name of NETCDF file ", clhstnam    ! filename
482         CALL histbeg( clhstnam, jpi, glamt, jpj, gphit,           &  ! Horizontal grid: glamt and gphit
483            &          iimi, iima-iimi+1, ijmi, ijma-ijmi+1,       &
484            &          nit000-1, zjulian, rdt, nh_T, nid_T, domain_id=nidom, snc4chunks=snc4set )
485         CALL histvert( nid_T, "deptht", "Vertical T levels",      &  ! Vertical grid: gdept
486            &           "m", ipk, gdept_1d, nz_T, "down" )
487         !                                                            ! Index of ocean points
488         CALL wheneq( jpi*jpj*ipk, tmask, 1, 1., ndex_T , ndim_T  )      ! volume
489         CALL wheneq( jpi*jpj    , tmask, 1, 1., ndex_hT, ndim_hT )      ! surface
490         !
491         IF( ln_icebergs ) THEN
492            !
493            !! allocation cant go in dia_wri_alloc because ln_icebergs is only set after
494            !! that routine is called from nemogcm, so do it here immediately before its needed
495            ALLOCATE( ndex_bT(jpi*jpj*nclasses), STAT=ierror )
496            CALL mpp_sum( 'diawri', ierror )
497            IF( ierror /= 0 ) THEN
498               CALL ctl_stop('dia_wri: failed to allocate iceberg diagnostic array')
499               RETURN
500            ENDIF
501            !
502            !! iceberg vertical coordinate is class number
503            CALL histvert( nid_T, "class", "Iceberg class",      &  ! Vertical grid: class
504               &           "number", nclasses, class_num, nb_T )
505            !
506            !! each class just needs the surface index pattern
507            ndim_bT = 3
508            DO jn = 1,nclasses
509               ndex_bT((jn-1)*jpi*jpj+1:jn*jpi*jpj) = ndex_hT(1:jpi*jpj)
510            ENDDO
511            !
512         ENDIF
513
514         ! Define the U grid FILE ( nid_U )
515
516         CALL dia_nam( clhstnam, nn_write, 'grid_U' )
517         IF(lwp) WRITE(numout,*) " Name of NETCDF file ", clhstnam    ! filename
518         CALL histbeg( clhstnam, jpi, glamu, jpj, gphiu,           &  ! Horizontal grid: glamu and gphiu
519            &          iimi, iima-iimi+1, ijmi, ijma-ijmi+1,       &
520            &          nit000-1, zjulian, rdt, nh_U, nid_U, domain_id=nidom, snc4chunks=snc4set )
521         CALL histvert( nid_U, "depthu", "Vertical U levels",      &  ! Vertical grid: gdept
522            &           "m", ipk, gdept_1d, nz_U, "down" )
523         !                                                            ! Index of ocean points
524         CALL wheneq( jpi*jpj*ipk, umask, 1, 1., ndex_U , ndim_U  )      ! volume
525         CALL wheneq( jpi*jpj    , umask, 1, 1., ndex_hU, ndim_hU )      ! surface
526
527         ! Define the V grid FILE ( nid_V )
528
529         CALL dia_nam( clhstnam, nn_write, 'grid_V' )                   ! filename
530         IF(lwp) WRITE(numout,*) " Name of NETCDF file ", clhstnam
531         CALL histbeg( clhstnam, jpi, glamv, jpj, gphiv,           &  ! Horizontal grid: glamv and gphiv
532            &          iimi, iima-iimi+1, ijmi, ijma-ijmi+1,       &
533            &          nit000-1, zjulian, rdt, nh_V, nid_V, domain_id=nidom, snc4chunks=snc4set )
534         CALL histvert( nid_V, "depthv", "Vertical V levels",      &  ! Vertical grid : gdept
535            &          "m", ipk, gdept_1d, nz_V, "down" )
536         !                                                            ! Index of ocean points
537         CALL wheneq( jpi*jpj*ipk, vmask, 1, 1., ndex_V , ndim_V  )      ! volume
538         CALL wheneq( jpi*jpj    , vmask, 1, 1., ndex_hV, ndim_hV )      ! surface
539
540         ! Define the W grid FILE ( nid_W )
541
542         CALL dia_nam( clhstnam, nn_write, 'grid_W' )                   ! filename
543         IF(lwp) WRITE(numout,*) " Name of NETCDF file ", clhstnam
544         CALL histbeg( clhstnam, jpi, glamt, jpj, gphit,           &  ! Horizontal grid: glamt and gphit
545            &          iimi, iima-iimi+1, ijmi, ijma-ijmi+1,       &
546            &          nit000-1, zjulian, rdt, nh_W, nid_W, domain_id=nidom, snc4chunks=snc4set )
547         CALL histvert( nid_W, "depthw", "Vertical W levels",      &  ! Vertical grid: gdepw
548            &          "m", ipk, gdepw_1d, nz_W, "down" )
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
566
567         ! Declare all the output fields as NETCDF variables
568
569         !                                                                                      !!! nid_T : 3D
570         CALL histdef( nid_T, "votemper", "Temperature"                        , "C"      ,   &  ! tn
571            &          jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout )
572         CALL histdef( nid_T, "vosaline", "Salinity"                           , "PSU"    ,   &  ! sn
573            &          jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout )
574         IF(  .NOT.ln_linssh  ) THEN
575            CALL histdef( nid_T, "vovvle3t", "Level thickness"                    , "m"      ,&  ! e3t(:,:,:,Kmm)
576            &             jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout )
577            CALL histdef( nid_T, "vovvldep", "T point depth"                      , "m"      ,&  ! e3t(:,:,:,Kmm)
578            &             jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout )
579            CALL histdef( nid_T, "vovvldef", "Squared level deformation"          , "%^2"    ,&  ! e3t(:,:,:,Kmm)
580            &             jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout )
581         ENDIF
582         !                                                                                      !!! nid_T : 2D
583         CALL histdef( nid_T, "sosstsst", "Sea Surface temperature"            , "C"      ,   &  ! sst
584            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
585         CALL histdef( nid_T, "sosaline", "Sea Surface Salinity"               , "PSU"    ,   &  ! sss
586            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
587         CALL histdef( nid_T, "sossheig", "Sea Surface Height"                 , "m"      ,   &  ! ssh
588            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
589         CALL histdef( nid_T, "sowaflup", "Net Upward Water Flux"              , "Kg/m2/s",   &  ! (emp-rnf)
590            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
591         CALL histdef( nid_T, "sorunoff", "River runoffs"                      , "Kg/m2/s",   &  ! runoffs
592            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
593         CALL histdef( nid_T, "sosfldow", "downward salt flux"                 , "PSU/m2/s",  &  ! sfx
594            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
595         IF(  ln_linssh  ) THEN
596            CALL histdef( nid_T, "sosst_cd", "Concentration/Dilution term on temperature"     &  ! emp * ts(:,:,1,jp_tem,Kmm)
597            &                                                                  , "KgC/m2/s",  &  ! sosst_cd
598            &             jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
599            CALL histdef( nid_T, "sosss_cd", "Concentration/Dilution term on salinity"        &  ! emp * ts(:,:,1,jp_sal,Kmm)
600            &                                                                  , "KgPSU/m2/s",&  ! sosss_cd
601            &             jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
602         ENDIF
603         CALL histdef( nid_T, "sohefldo", "Net Downward Heat Flux"             , "W/m2"   ,   &  ! qns + qsr
604            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
605         CALL histdef( nid_T, "soshfldo", "Shortwave Radiation"                , "W/m2"   ,   &  ! qsr
606            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
607         CALL histdef( nid_T, "somixhgt", "Turbocline Depth"                   , "m"      ,   &  ! hmld
608            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
609         CALL histdef( nid_T, "somxl010", "Mixed Layer Depth 0.01"             , "m"      ,   &  ! hmlp
610            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
611         CALL histdef( nid_T, "soicecov", "Ice fraction"                       , "[0,1]"  ,   &  ! fr_i
612            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
613         CALL histdef( nid_T, "sowindsp", "wind speed at 10m"                  , "m/s"    ,   &  ! wndm
614            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
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         !
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
673         IF( ln_ssr ) THEN
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
681       
682         clmx ="l_max(only(x))"    ! max index on a period
683!         CALL histdef( nid_T, "sobowlin", "Bowl Index"                         , "W-point",   &  ! bowl INDEX
684!            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clmx, zsto, zout )
685#if defined key_diahth
686         CALL histdef( nid_T, "sothedep", "Thermocline Depth"                  , "m"      ,   & ! hth
687            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
688         CALL histdef( nid_T, "so20chgt", "Depth of 20C isotherm"              , "m"      ,   & ! hd20
689            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
690         CALL histdef( nid_T, "so28chgt", "Depth of 28C isotherm"              , "m"      ,   & ! hd28
691            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
692         CALL histdef( nid_T, "sohtc300", "Heat content 300 m"                 , "J/m2"   ,   & ! htc3
693            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
694#endif
695
696         CALL histend( nid_T, snc4chunks=snc4set )
697
698         !                                                                                      !!! nid_U : 3D
699         CALL histdef( nid_U, "vozocrtx", "Zonal Current"                      , "m/s"    ,   &  ! uu(:,:,:,Kmm)
700            &          jpi, jpj, nh_U, ipk, 1, ipk, nz_U, 32, clop, zsto, zout )
701         IF( ln_wave .AND. ln_sdw) THEN
702            CALL histdef( nid_U, "sdzocrtx", "Stokes Drift Zonal Current"         , "m/s"    ,   &  ! usd
703               &          jpi, jpj, nh_U, ipk, 1, ipk, nz_U, 32, clop, zsto, zout )
704         ENDIF
705         !                                                                                      !!! nid_U : 2D
706         CALL histdef( nid_U, "sozotaux", "Wind Stress along i-axis"           , "N/m2"   ,   &  ! utau
707            &          jpi, jpj, nh_U, 1  , 1, 1  , - 99, 32, clop, zsto, zout )
708
709         CALL histend( nid_U, snc4chunks=snc4set )
710
711         !                                                                                      !!! nid_V : 3D
712         CALL histdef( nid_V, "vomecrty", "Meridional Current"                 , "m/s"    ,   &  ! vv(:,:,:,Kmm)
713            &          jpi, jpj, nh_V, ipk, 1, ipk, nz_V, 32, clop, zsto, zout )
714         IF( ln_wave .AND. ln_sdw) THEN
715            CALL histdef( nid_V, "sdmecrty", "Stokes Drift Meridional Current"    , "m/s"    ,   &  ! vsd
716               &          jpi, jpj, nh_V, ipk, 1, ipk, nz_V, 32, clop, zsto, zout )
717         ENDIF
718         !                                                                                      !!! nid_V : 2D
719         CALL histdef( nid_V, "sometauy", "Wind Stress along j-axis"           , "N/m2"   ,   &  ! vtau
720            &          jpi, jpj, nh_V, 1  , 1, 1  , - 99, 32, clop, zsto, zout )
721
722         CALL histend( nid_V, snc4chunks=snc4set )
723
724         !                                                                                      !!! nid_W : 3D
725         CALL histdef( nid_W, "vovecrtz", "Vertical Velocity"                  , "m/s"    ,   &  ! ww
726            &          jpi, jpj, nh_W, ipk, 1, ipk, nz_W, 32, clop, zsto, zout )
727         CALL histdef( nid_W, "votkeavt", "Vertical Eddy Diffusivity"          , "m2/s"   ,   &  ! avt
728            &          jpi, jpj, nh_W, ipk, 1, ipk, nz_W, 32, clop, zsto, zout )
729         CALL histdef( nid_W, "votkeavm", "Vertical Eddy Viscosity"             , "m2/s"  ,   &  ! avm
730            &          jpi, jpj, nh_W, ipk, 1, ipk, nz_W, 32, clop, zsto, zout )
731
732         IF( ln_zdfddm ) THEN
733            CALL histdef( nid_W,"voddmavs","Salt Vertical Eddy Diffusivity"    , "m2/s"   ,   &  ! avs
734               &          jpi, jpj, nh_W, ipk, 1, ipk, nz_W, 32, clop, zsto, zout )
735         ENDIF
736         
737         IF( ln_wave .AND. ln_sdw) THEN
738            CALL histdef( nid_W, "sdvecrtz", "Stokes Drift Vertical Current"   , "m/s"    ,   &  ! wsd
739               &          jpi, jpj, nh_W, ipk, 1, ipk, nz_W, 32, clop, zsto, zout )
740         ENDIF
741         !                                                                                      !!! nid_W : 2D
742         CALL histend( nid_W, snc4chunks=snc4set )
743
744         IF(lwp) WRITE(numout,*)
745         IF(lwp) WRITE(numout,*) 'End of NetCDF Initialization'
746         IF(ll_print) CALL FLUSH(numout )
747
748      ENDIF
749
750      ! 2. Start writing data
751      ! ---------------------
752
753      ! ndex(1) est utilise ssi l'avant dernier argument est different de
754      ! la taille du tableau en sortie. Dans ce cas , l'avant dernier argument
755      ! donne le nombre d'elements, et ndex la liste des indices a sortir
756
757      IF( lwp .AND. MOD( itmod, nn_write ) == 0 ) THEN
758         WRITE(numout,*) 'dia_wri : write model outputs in NetCDF files at ', kt, 'time-step'
759         WRITE(numout,*) '~~~~~~ '
760      ENDIF
761
762      IF( .NOT.ln_linssh ) THEN
763         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
767      ELSE
768         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
772      ENDIF
773      IF( .NOT.ln_linssh ) THEN
774         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
777         CALL histwrite( nid_T, "vovvldef", it, zw3d             , ndim_T , ndex_T  )   ! level thickness deformation
778      ENDIF
779      CALL histwrite( nid_T, "sossheig", it, ssh(:,:,Kmm)          , ndim_hT, ndex_hT )   ! sea surface height
780      CALL histwrite( nid_T, "sowaflup", it, ( emp-rnf )   , ndim_hT, ndex_hT )   ! upward water flux
781      CALL histwrite( nid_T, "sorunoff", it, rnf           , ndim_hT, ndex_hT )   ! river runoffs
782      CALL histwrite( nid_T, "sosfldow", it, sfx           , ndim_hT, ndex_hT )   ! downward salt flux
783                                                                                  ! (includes virtual salt flux beneath ice
784                                                                                  ! in linear free surface case)
785      IF( ln_linssh ) THEN
786         zw2d(:,:) = emp (:,:) * ts(:,:,1,jp_tem,Kmm)
787         CALL histwrite( nid_T, "sosst_cd", it, zw2d, ndim_hT, ndex_hT )          ! c/d term on sst
788         zw2d(:,:) = emp (:,:) * ts(:,:,1,jp_sal,Kmm)
789         CALL histwrite( nid_T, "sosss_cd", it, zw2d, ndim_hT, ndex_hT )          ! c/d term on sss
790      ENDIF
791      CALL histwrite( nid_T, "sohefldo", it, qns + qsr     , ndim_hT, ndex_hT )   ! total heat flux
792      CALL histwrite( nid_T, "soshfldo", it, qsr           , ndim_hT, ndex_hT )   ! solar heat flux
793      CALL histwrite( nid_T, "somixhgt", it, hmld          , ndim_hT, ndex_hT )   ! turbocline depth
794      CALL histwrite( nid_T, "somxl010", it, hmlp          , ndim_hT, ndex_hT )   ! mixed layer depth
795      CALL histwrite( nid_T, "soicecov", it, fr_i          , ndim_hT, ndex_hT )   ! ice fraction   
796      CALL histwrite( nid_T, "sowindsp", it, wndm          , ndim_hT, ndex_hT )   ! wind speed   
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      !
821      IF( ln_icebergs ) THEN
822         !
823         CALL histwrite( nid_T, "calving"             , it, berg_grid%calving      , ndim_hT, ndex_hT ) 
824         CALL histwrite( nid_T, "calving_heat"        , it, berg_grid%calving_hflx , ndim_hT, ndex_hT )         
825         CALL histwrite( nid_T, "berg_floating_melt"  , it, berg_grid%floating_melt, ndim_hT, ndex_hT ) 
826         !
827         CALL histwrite( nid_T, "berg_stored_ice"     , it, berg_grid%stored_ice   , ndim_bT, ndex_bT )
828         !
829         IF( ln_bergdia ) THEN
830            CALL histwrite( nid_T, "berg_melt"           , it, berg_melt        , ndim_hT, ndex_hT   ) 
831            CALL histwrite( nid_T, "berg_buoy_melt"      , it, buoy_melt        , ndim_hT, ndex_hT   ) 
832            CALL histwrite( nid_T, "berg_eros_melt"      , it, eros_melt        , ndim_hT, ndex_hT   ) 
833            CALL histwrite( nid_T, "berg_conv_melt"      , it, conv_melt        , ndim_hT, ndex_hT   ) 
834            CALL histwrite( nid_T, "berg_virtual_area"   , it, virtual_area     , ndim_hT, ndex_hT   ) 
835            CALL histwrite( nid_T, "bits_src"            , it, bits_src         , ndim_hT, ndex_hT   ) 
836            CALL histwrite( nid_T, "bits_melt"           , it, bits_melt        , ndim_hT, ndex_hT   ) 
837            CALL histwrite( nid_T, "bits_mass"           , it, bits_mass        , ndim_hT, ndex_hT   ) 
838            CALL histwrite( nid_T, "berg_mass"           , it, berg_mass        , ndim_hT, ndex_hT   ) 
839            !
840            CALL histwrite( nid_T, "berg_real_calving"   , it, real_calving     , ndim_bT, ndex_bT   )
841         ENDIF
842      ENDIF
843
844      IF( ln_ssr ) THEN
845         CALL histwrite( nid_T, "sohefldp", it, qrp           , ndim_hT, ndex_hT )   ! heat flux damping
846         CALL histwrite( nid_T, "sowafldp", it, erp           , ndim_hT, ndex_hT )   ! freshwater flux damping
847         zw2d(:,:) = erp(:,:) * ts(:,:,1,jp_sal,Kmm) * tmask(:,:,1)
848         CALL histwrite( nid_T, "sosafldp", it, zw2d          , ndim_hT, ndex_hT )   ! salt flux damping
849      ENDIF
850!      zw2d(:,:) = FLOAT( nmln(:,:) ) * tmask(:,:,1)
851!      CALL histwrite( nid_T, "sobowlin", it, zw2d          , ndim_hT, ndex_hT )   ! ???
852
853#if defined key_diahth
854      CALL histwrite( nid_T, "sothedep", it, hth           , ndim_hT, ndex_hT )   ! depth of the thermocline
855      CALL histwrite( nid_T, "so20chgt", it, hd20          , ndim_hT, ndex_hT )   ! depth of the 20 isotherm
856      CALL histwrite( nid_T, "so28chgt", it, hd28          , ndim_hT, ndex_hT )   ! depth of the 28 isotherm
857      CALL histwrite( nid_T, "sohtc300", it, htc3          , ndim_hT, ndex_hT )   ! first 300m heaat content
858#endif
859
860      CALL histwrite( nid_U, "vozocrtx", it, uu(:,:,:,Kmm)            , ndim_U , ndex_U )    ! i-current
861      CALL histwrite( nid_U, "sozotaux", it, utau          , ndim_hU, ndex_hU )   ! i-wind stress
862
863      CALL histwrite( nid_V, "vomecrty", it, vv(:,:,:,Kmm)            , ndim_V , ndex_V  )   ! j-current
864      CALL histwrite( nid_V, "sometauy", it, vtau          , ndim_hV, ndex_hV )   ! j-wind stress
865
866      IF( ln_zad_Aimp ) THEN
867         CALL histwrite( nid_W, "vovecrtz", it, ww + wi     , ndim_T, ndex_T )    ! vert. current
868      ELSE
869         CALL histwrite( nid_W, "vovecrtz", it, ww          , ndim_T, ndex_T )    ! vert. current
870      ENDIF
871      CALL histwrite( nid_W, "votkeavt", it, avt            , ndim_T, ndex_T )    ! T vert. eddy diff. coef.
872      CALL histwrite( nid_W, "votkeavm", it, avm            , ndim_T, ndex_T )    ! T vert. eddy visc. coef.
873      IF( ln_zdfddm ) THEN
874         CALL histwrite( nid_W, "voddmavs", it, avs         , ndim_T, ndex_T )    ! S vert. eddy diff. coef.
875      ENDIF
876
877      IF( ln_wave .AND. ln_sdw ) THEN
878         CALL histwrite( nid_U, "sdzocrtx", it, usd         , ndim_U , ndex_U )    ! i-StokesDrift-current
879         CALL histwrite( nid_V, "sdmecrty", it, vsd         , ndim_V , ndex_V )    ! j-StokesDrift-current
880         CALL histwrite( nid_W, "sdvecrtz", it, wsd         , ndim_T , ndex_T )    ! StokesDrift vert. current
881      ENDIF
882
883      ! 3. Close all files
884      ! ---------------------------------------
885      IF( kt == nitend ) THEN
886         CALL histclo( nid_T )
887         CALL histclo( nid_U )
888         CALL histclo( nid_V )
889         CALL histclo( nid_W )
890         IF(ln_abl) CALL histclo( nid_A )
891      ENDIF
892      !
893      IF( ln_timing )   CALL timing_stop('dia_wri')
894      !
895   END SUBROUTINE dia_wri
896#endif
897
898   SUBROUTINE dia_wri_state( Kmm, cdfile_name )
899      !!---------------------------------------------------------------------
900      !!                 ***  ROUTINE dia_wri_state  ***
901      !!       
902      !! ** Purpose :   create a NetCDF file named cdfile_name which contains
903      !!      the instantaneous ocean state and forcing fields.
904      !!        Used to find errors in the initial state or save the last
905      !!      ocean state in case of abnormal end of a simulation
906      !!
907      !! ** Method  :   NetCDF files using ioipsl
908      !!      File 'output.init.nc'  is created if ninist = 1 (namelist)
909      !!      File 'output.abort.nc' is created in case of abnormal job end
910      !!----------------------------------------------------------------------
911      INTEGER           , INTENT( in ) ::   Kmm              ! time level index
912      CHARACTER (len=* ), INTENT( in ) ::   cdfile_name      ! name of the file created
913      !!
914      INTEGER :: inum, jk
915      !!----------------------------------------------------------------------
916      !
917      IF(lwp) WRITE(numout,*)
918      IF(lwp) WRITE(numout,*) 'dia_wri_state : single instantaneous ocean state'
919      IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~   and forcing fields file created '
920      IF(lwp) WRITE(numout,*) '                and named :', cdfile_name, '...nc'
921
922#if defined key_si3
923     CALL iom_open( TRIM(cdfile_name), inum, ldwrt = .TRUE., kdlev = jpl )
924#else
925     CALL iom_open( TRIM(cdfile_name), inum, ldwrt = .TRUE. )
926#endif
927
928      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
933      IF( ln_zad_Aimp ) THEN
934         CALL iom_rstput( 0, 0, inum, 'vovecrtz', ww + wi        )    ! now k-velocity
935      ELSE
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
961      IF( ALLOCATED(ahtu) ) THEN
962         CALL iom_rstput( 0, 0, inum,  'ahtu', ahtu              )    ! aht at u-point
963         CALL iom_rstput( 0, 0, inum,  'ahtv', ahtv              )    ! aht at v-point
964      ENDIF
965      IF( ALLOCATED(ahmt) ) THEN
966         CALL iom_rstput( 0, 0, inum,  'ahmt', ahmt              )    ! ahmt at u-point
967         CALL iom_rstput( 0, 0, inum,  'ahmf', ahmf              )    ! ahmf at v-point
968      ENDIF
969      CALL iom_rstput( 0, 0, inum, 'sowaflup', emp - rnf         )    ! freshwater budget
970      CALL iom_rstput( 0, 0, inum, 'sohefldo', qsr + qns         )    ! total heat flux
971      CALL iom_rstput( 0, 0, inum, 'soshfldo', qsr               )    ! solar heat flux
972      CALL iom_rstput( 0, 0, inum, 'soicecov', fr_i              )    ! ice fraction
973      CALL iom_rstput( 0, 0, inum, 'sozotaux', utau              )    ! i-wind stress
974      CALL iom_rstput( 0, 0, inum, 'sometauy', vtau              )    ! j-wind stress
975      IF(  .NOT.ln_linssh  ) THEN             
976         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 
978      END IF
979      IF( ln_wave .AND. ln_sdw ) THEN
980         CALL iom_rstput( 0, 0, inum, 'sdzocrtx', usd            )    ! now StokesDrift i-velocity
981         CALL iom_rstput( 0, 0, inum, 'sdmecrty', vsd            )    ! now StokesDrift j-velocity
982         CALL iom_rstput( 0, 0, inum, 'sdvecrtz', wsd            )    ! now StokesDrift k-velocity
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
990 
991#if defined key_si3
992      IF( nn_ice == 2 ) THEN   ! condition needed in case agrif + ice-model but no-ice in child grid
993         CALL ice_wri_state( inum )
994      ENDIF
995#endif
996      !
997      CALL iom_close( inum )
998      !
999   END SUBROUTINE dia_wri_state
1000
1001   !!======================================================================
1002END MODULE diawri
Note: See TracBrowser for help on using the repository browser.