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/2020/dev_r12512_HPC-04_mcastril_Mixed_Precision_implementation/src/OCE/DIA – NEMO

source: NEMO/branches/2020/dev_r12512_HPC-04_mcastril_Mixed_Precision_implementation/src/OCE/DIA/diawri.F90 @ 13135

Last change on this file since 13135 was 13135, checked in by orioltp, 4 years ago

dev_r12512_HPC-04_mcastril_Mixed_Precision_implementation: merge with trunk@13134, see #2364

  • Property svn:keywords set to Id
File size: 55.5 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      CALL iom_put( "rhop", rhop(:,:,:) )          ! 3D potential density (sigma0)
175
176      IF ( iom_use("taubot") ) THEN                ! bottom stress
177         zztmp = rho0 * 0.25
178         z2d(:,:) = 0._wp
179         DO_2D_00_00
180            zztmp2 = (  ( rCdU_bot(ji+1,jj)+rCdU_bot(ji  ,jj) ) * uu(ji  ,jj,mbku(ji  ,jj),Kmm)  )**2   &
181               &   + (  ( rCdU_bot(ji  ,jj)+rCdU_bot(ji-1,jj) ) * uu(ji-1,jj,mbku(ji-1,jj),Kmm)  )**2   &
182               &   + (  ( rCdU_bot(ji,jj+1)+rCdU_bot(ji,jj  ) ) * vv(ji,jj  ,mbkv(ji,jj  ),Kmm)  )**2   &
183               &   + (  ( rCdU_bot(ji,jj  )+rCdU_bot(ji,jj-1) ) * vv(ji,jj-1,mbkv(ji,jj-1),Kmm)  )**2
184            z2d(ji,jj) = zztmp * SQRT( zztmp2 ) * tmask(ji,jj,1) 
185            !
186         END_2D
187         CALL lbc_lnk( 'diawri', z2d, 'T', 1.0_wp )
188         CALL iom_put( "taubot", z2d )           
189      ENDIF
190         
191      CALL iom_put( "uoce", uu(:,:,:,Kmm) )            ! 3D i-current
192      CALL iom_put(  "ssu", uu(:,:,1,Kmm) )            ! surface i-current
193      IF ( iom_use("sbu") ) THEN
194         DO_2D_11_11
195            ikbot = mbku(ji,jj)
196            z2d(ji,jj) = uu(ji,jj,ikbot,Kmm)
197         END_2D
198         CALL iom_put( "sbu", z2d )                ! bottom i-current
199      ENDIF
200     
201      CALL iom_put( "voce", vv(:,:,:,Kmm) )            ! 3D j-current
202      CALL iom_put(  "ssv", vv(:,:,1,Kmm) )            ! surface j-current
203      IF ( iom_use("sbv") ) THEN
204         DO_2D_11_11
205            ikbot = mbkv(ji,jj)
206            z2d(ji,jj) = vv(ji,jj,ikbot,Kmm)
207         END_2D
208         CALL iom_put( "sbv", z2d )                ! bottom j-current
209      ENDIF
210
211      IF( ln_zad_Aimp ) ww = ww + wi               ! Recombine explicit and implicit parts of vertical velocity for diagnostic output
212      !
213      CALL iom_put( "woce", ww )                   ! vertical velocity
214      IF( iom_use('w_masstr') .OR. iom_use('w_masstr2') ) THEN   ! vertical mass transport & its square value
215         ! Caution: in the VVL case, it only correponds to the baroclinic mass transport.
216         z2d(:,:) = rho0 * e1e2t(:,:)
217         DO jk = 1, jpk
218            z3d(:,:,jk) = ww(:,:,jk) * z2d(:,:)
219         END DO
220         CALL iom_put( "w_masstr" , z3d ) 
221         IF( iom_use('w_masstr2') )   CALL iom_put( "w_masstr2", z3d(:,:,:) * z3d(:,:,:) )
222      ENDIF
223      !
224      IF( ln_zad_Aimp ) ww = ww - wi               ! Remove implicit part of vertical velocity that was added for diagnostic output
225
226      CALL iom_put( "avt" , avt )                  ! T vert. eddy diff. coef.
227      CALL iom_put( "avs" , avs )                  ! S vert. eddy diff. coef.
228      CALL iom_put( "avm" , avm )                  ! T vert. eddy visc. coef.
229
230      IF( iom_use('logavt') )   CALL iom_put( "logavt", LOG( MAX( 1.e-20_wp, avt(:,:,:) ) ) )
231      IF( iom_use('logavs') )   CALL iom_put( "logavs", LOG( MAX( 1.e-20_wp, avs(:,:,:) ) ) )
232
233      IF ( iom_use("sstgrad") .OR. iom_use("sstgrad2") ) THEN
234         DO_2D_00_00
235            zztmp  = ts(ji,jj,1,jp_tem,Kmm)
236            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)
237            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)
238            z2d(ji,jj) = 0.25 * ( zztmpx * zztmpx + zztmpy * zztmpy )   &
239               &              * umask(ji,jj,1) * umask(ji-1,jj,1) * vmask(ji,jj,1) * umask(ji,jj-1,1)
240         END_2D
241         CALL lbc_lnk( 'diawri', z2d, 'T', 1.0_wp )
242         CALL iom_put( "sstgrad2",  z2d )          ! square of module of sst gradient
243         z2d(:,:) = SQRT( z2d(:,:) )
244         CALL iom_put( "sstgrad" ,  z2d )          ! module of sst gradient
245      ENDIF
246         
247      ! heat and salt contents
248      IF( iom_use("heatc") ) THEN
249         z2d(:,:)  = 0._wp 
250         DO_3D_11_11( 1, jpkm1 )
251            z2d(ji,jj) = z2d(ji,jj) + e3t(ji,jj,jk,Kmm) * ts(ji,jj,jk,jp_tem,Kmm) * tmask(ji,jj,jk)
252         END_3D
253         CALL iom_put( "heatc", rho0_rcp * z2d )   ! vertically integrated heat content (J/m2)
254      ENDIF
255
256      IF( iom_use("saltc") ) THEN
257         z2d(:,:)  = 0._wp 
258         DO_3D_11_11( 1, jpkm1 )
259            z2d(ji,jj) = z2d(ji,jj) + e3t(ji,jj,jk,Kmm) * ts(ji,jj,jk,jp_sal,Kmm) * tmask(ji,jj,jk)
260         END_3D
261         CALL iom_put( "saltc", rho0 * z2d )          ! vertically integrated salt content (PSU*kg/m2)
262      ENDIF
263      !
264      IF ( iom_use("eken") ) THEN
265         z3d(:,:,jpk) = 0._wp 
266         DO_3D_00_00( 1, jpkm1 )
267            zztmp  = 0.25_wp * r1_e1e2t(ji,jj) / e3t(ji,jj,jk,Kmm)
268            z3d(ji,jj,jk) = zztmp * (  uu(ji-1,jj,jk,Kmm)**2 * e2u(ji-1,jj) * e3u(ji-1,jj,jk,Kmm)   &
269               &                     + uu(ji  ,jj,jk,Kmm)**2 * e2u(ji  ,jj) * e3u(ji  ,jj,jk,Kmm)   &
270               &                     + vv(ji,jj-1,jk,Kmm)**2 * e1v(ji,jj-1) * e3v(ji,jj-1,jk,Kmm)   &
271               &                     + vv(ji,jj  ,jk,Kmm)**2 * e1v(ji,jj  ) * e3v(ji,jj  ,jk,Kmm)   )
272         END_3D
273         CALL lbc_lnk( 'diawri', z3d, 'T', 1.0_wp )
274         CALL iom_put( "eken", z3d )                 ! kinetic energy
275      ENDIF
276      !
277      CALL iom_put( "hdiv", hdiv )                  ! Horizontal divergence
278      !
279      IF( iom_use("u_masstr") .OR. iom_use("u_masstr_vint") .OR. iom_use("u_heattr") .OR. iom_use("u_salttr") ) THEN
280         z3d(:,:,jpk) = 0.e0
281         z2d(:,:) = 0.e0
282         DO jk = 1, jpkm1
283            z3d(:,:,jk) = rho0 * uu(:,:,jk,Kmm) * e2u(:,:) * e3u(:,:,jk,Kmm) * umask(:,:,jk)
284            z2d(:,:) = z2d(:,:) + z3d(:,:,jk)
285         END DO
286         CALL iom_put( "u_masstr"     , z3d )         ! mass transport in i-direction
287         CALL iom_put( "u_masstr_vint", z2d )         ! mass transport in i-direction vertical sum
288      ENDIF
289     
290      IF( iom_use("u_heattr") ) THEN
291         z2d(:,:) = 0._wp 
292         DO_3D_00_00( 1, jpkm1 )
293            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) )
294         END_3D
295         CALL lbc_lnk( 'diawri', z2d, 'U', -1.0_wp )
296         CALL iom_put( "u_heattr", 0.5*rcp * z2d )    ! heat transport in i-direction
297      ENDIF
298
299      IF( iom_use("u_salttr") ) THEN
300         z2d(:,:) = 0.e0 
301         DO_3D_00_00( 1, jpkm1 )
302            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) )
303         END_3D
304         CALL lbc_lnk( 'diawri', z2d, 'U', -1.0_wp )
305         CALL iom_put( "u_salttr", 0.5 * z2d )        ! heat transport in i-direction
306      ENDIF
307
308     
309      IF( iom_use("v_masstr") .OR. iom_use("v_heattr") .OR. iom_use("v_salttr") ) THEN
310         z3d(:,:,jpk) = 0.e0
311         DO jk = 1, jpkm1
312            z3d(:,:,jk) = rho0 * vv(:,:,jk,Kmm) * e1v(:,:) * e3v(:,:,jk,Kmm) * vmask(:,:,jk)
313         END DO
314         CALL iom_put( "v_masstr", z3d )              ! mass transport in j-direction
315      ENDIF
316     
317      IF( iom_use("v_heattr") ) THEN
318         z2d(:,:) = 0.e0 
319         DO_3D_00_00( 1, jpkm1 )
320            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) )
321         END_3D
322         CALL lbc_lnk( 'diawri', z2d, 'V', -1.0_wp )
323         CALL iom_put( "v_heattr", 0.5*rcp * z2d )    !  heat transport in j-direction
324      ENDIF
325
326      IF( iom_use("v_salttr") ) THEN
327         z2d(:,:) = 0._wp 
328         DO_3D_00_00( 1, jpkm1 )
329            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) )
330         END_3D
331         CALL lbc_lnk( 'diawri', z2d, 'V', -1.0_wp )
332         CALL iom_put( "v_salttr", 0.5 * z2d )        !  heat transport in j-direction
333      ENDIF
334
335      IF( iom_use("tosmint") ) THEN
336         z2d(:,:) = 0._wp
337         DO_3D_00_00( 1, jpkm1 )
338            z2d(ji,jj) = z2d(ji,jj) + e3t(ji,jj,jk,Kmm) *  ts(ji,jj,jk,jp_tem,Kmm)
339         END_3D
340         CALL lbc_lnk( 'diawri', z2d, 'T', -1.0_wp )
341         CALL iom_put( "tosmint", rho0 * z2d )        ! Vertical integral of temperature
342      ENDIF
343      IF( iom_use("somint") ) THEN
344         z2d(:,:)=0._wp
345         DO_3D_00_00( 1, jpkm1 )
346            z2d(ji,jj) = z2d(ji,jj) + e3t(ji,jj,jk,Kmm) * ts(ji,jj,jk,jp_sal,Kmm)
347         END_3D
348         CALL lbc_lnk( 'diawri', z2d, 'T', -1.0_wp )
349         CALL iom_put( "somint", rho0 * z2d )         ! Vertical integral of salinity
350      ENDIF
351
352      CALL iom_put( "bn2", rn2 )                      ! Brunt-Vaisala buoyancy frequency (N^2)
353      !
354     
355      IF (ln_dia25h)   CALL dia_25h( kt, Kmm )        ! 25h averaging
356
357      IF( ln_timing )   CALL timing_stop('dia_wri')
358      !
359   END SUBROUTINE dia_wri
360
361#else
362   !!----------------------------------------------------------------------
363   !!   Default option                                  use IOIPSL  library
364   !!----------------------------------------------------------------------
365
366   INTEGER FUNCTION dia_wri_alloc()
367      !!----------------------------------------------------------------------
368      INTEGER, DIMENSION(2) :: ierr
369      !!----------------------------------------------------------------------
370      IF( nn_write == -1 ) THEN
371         dia_wri_alloc = 0
372      ELSE   
373         ierr = 0
374         ALLOCATE( ndex_hT(jpi*jpj) , ndex_T(jpi*jpj*jpk) ,     &
375            &      ndex_hU(jpi*jpj) , ndex_U(jpi*jpj*jpk) ,     &
376            &      ndex_hV(jpi*jpj) , ndex_V(jpi*jpj*jpk) , STAT=ierr(1) )
377         !
378         dia_wri_alloc = MAXVAL(ierr)
379         CALL mpp_sum( 'diawri', dia_wri_alloc )
380         !
381      ENDIF
382      !
383   END FUNCTION dia_wri_alloc
384 
385   INTEGER FUNCTION dia_wri_alloc_abl()
386      !!----------------------------------------------------------------------
387     ALLOCATE(   ndex_hA(jpi*jpj), ndex_A (jpi*jpj*jpkam1), STAT=dia_wri_alloc_abl)
388      CALL mpp_sum( 'diawri', dia_wri_alloc_abl )
389      !
390   END FUNCTION dia_wri_alloc_abl
391
392   
393   SUBROUTINE dia_wri( kt, Kmm )
394      !!---------------------------------------------------------------------
395      !!                  ***  ROUTINE dia_wri  ***
396      !!                   
397      !! ** Purpose :   Standard output of opa: dynamics and tracer fields
398      !!      NETCDF format is used by default
399      !!
400      !! ** Method  :   At the beginning of the first time step (nit000),
401      !!      define all the NETCDF files and fields
402      !!      At each time step call histdef to compute the mean if ncessary
403      !!      Each nn_write time step, output the instantaneous or mean fields
404      !!----------------------------------------------------------------------
405      INTEGER, INTENT( in ) ::   kt   ! ocean time-step index
406      INTEGER, INTENT( in ) ::   Kmm  ! ocean time level index
407      !
408      LOGICAL ::   ll_print = .FALSE.                        ! =T print and flush numout
409      CHARACTER (len=40) ::   clhstnam, clop, clmx           ! local names
410      INTEGER  ::   inum = 11                                ! temporary logical unit
411      INTEGER  ::   ji, jj, jk                               ! dummy loop indices
412      INTEGER  ::   ierr                                     ! error code return from allocation
413      INTEGER  ::   iimi, iima, ipk, it, itmod, ijmi, ijma   ! local integers
414      INTEGER  ::   ipka                                     ! ABL
415      INTEGER  ::   jn, ierror                               ! local integers
416      REAL(wp) ::   zsto, zout, zmax, zjulian                ! local scalars
417      !
418      REAL(wp), DIMENSION(jpi,jpj)   :: zw2d       ! 2D workspace
419      REAL(wp), DIMENSION(jpi,jpj,jpk) :: zw3d       ! 3D workspace
420      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE :: zw3d_abl   ! ABL 3D workspace
421      !!----------------------------------------------------------------------
422      !
423      IF( ninist == 1 ) THEN     !==  Output the initial state and forcings  ==!
424         CALL dia_wri_state( Kmm, 'output.init' )
425         ninist = 0
426      ENDIF
427      !
428      IF( nn_write == -1 )   RETURN   ! we will never do any output
429      !
430      IF( ln_timing )   CALL timing_start('dia_wri')
431      !
432      ! 0. Initialisation
433      ! -----------------
434
435      ll_print = .FALSE.                  ! local variable for debugging
436      ll_print = ll_print .AND. lwp
437
438      ! Define frequency of output and means
439      clop = "x"         ! no use of the mask value (require less cpu time and otherwise the model crashes)
440#if defined key_diainstant
441      zsto = nn_write * rn_Dt
442      clop = "inst("//TRIM(clop)//")"
443#else
444      zsto=rn_Dt
445      clop = "ave("//TRIM(clop)//")"
446#endif
447      zout = nn_write * rn_Dt
448      zmax = ( nitend - nit000 + 1 ) * rn_Dt
449
450      ! Define indices of the horizontal output zoom and vertical limit storage
451      iimi = 1      ;      iima = jpi
452      ijmi = 1      ;      ijma = jpj
453      ipk = jpk
454      IF(ln_abl) ipka = jpkam1
455
456      ! define time axis
457      it = kt
458      itmod = kt - nit000 + 1
459
460
461      ! 1. Define NETCDF files and fields at beginning of first time step
462      ! -----------------------------------------------------------------
463
464      IF( kt == nit000 ) THEN
465
466         ! Define the NETCDF files (one per grid)
467
468         ! Compute julian date from starting date of the run
469         CALL ymds2ju( nyear, nmonth, nday, rn_Dt, zjulian )
470         zjulian = zjulian - adatrj   !   set calendar origin to the beginning of the experiment
471         IF(lwp)WRITE(numout,*)
472         IF(lwp)WRITE(numout,*) 'Date 0 used :', nit000, ' YEAR ', nyear,   &
473            &                    ' MONTH ', nmonth, ' DAY ', nday, 'Julian day : ', zjulian
474         IF(lwp)WRITE(numout,*) ' indexes of zoom = ', iimi, iima, ijmi, ijma,   &
475                                 ' limit storage in depth = ', ipk
476
477         ! WRITE root name in date.file for use by postpro
478         IF(lwp) THEN
479            CALL dia_nam( clhstnam, nn_write,' ' )
480            CALL ctl_opn( inum, 'date.file', 'REPLACE', 'FORMATTED', 'SEQUENTIAL', -1, numout, lwp, narea )
481            WRITE(inum,*) clhstnam
482            CLOSE(inum)
483         ENDIF
484
485         ! Define the T grid FILE ( nid_T )
486
487         CALL dia_nam( clhstnam, nn_write, 'grid_T' )
488         IF(lwp) WRITE(numout,*) " Name of NETCDF file ", clhstnam    ! filename
489         CALL histbeg( clhstnam, jpi, glamt, jpj, gphit,           &  ! Horizontal grid: glamt and gphit
490            &          iimi, iima-iimi+1, ijmi, ijma-ijmi+1,       &
491            &          nit000-1, zjulian, rn_Dt, nh_T, nid_T, domain_id=nidom, snc4chunks=snc4set )
492         CALL histvert( nid_T, "deptht", "Vertical T levels",      &  ! Vertical grid: gdept
493            &           "m", ipk, gdept_1d, nz_T, "down" )
494         !                                                            ! Index of ocean points
495         CALL wheneq( jpi*jpj*ipk, tmask, 1, 1., ndex_T , ndim_T  )      ! volume
496         CALL wheneq( jpi*jpj    , tmask, 1, 1., ndex_hT, ndim_hT )      ! surface
497         !
498         IF( ln_icebergs ) THEN
499            !
500            !! allocation cant go in dia_wri_alloc because ln_icebergs is only set after
501            !! that routine is called from nemogcm, so do it here immediately before its needed
502            ALLOCATE( ndex_bT(jpi*jpj*nclasses), STAT=ierror )
503            CALL mpp_sum( 'diawri', ierror )
504            IF( ierror /= 0 ) THEN
505               CALL ctl_stop('dia_wri: failed to allocate iceberg diagnostic array')
506               RETURN
507            ENDIF
508            !
509            !! iceberg vertical coordinate is class number
510            CALL histvert( nid_T, "class", "Iceberg class",      &  ! Vertical grid: class
511               &           "number", nclasses, class_num, nb_T )
512            !
513            !! each class just needs the surface index pattern
514            ndim_bT = 3
515            DO jn = 1,nclasses
516               ndex_bT((jn-1)*jpi*jpj+1:jn*jpi*jpj) = ndex_hT(1:jpi*jpj)
517            ENDDO
518            !
519         ENDIF
520
521         ! Define the U grid FILE ( nid_U )
522
523         CALL dia_nam( clhstnam, nn_write, 'grid_U' )
524         IF(lwp) WRITE(numout,*) " Name of NETCDF file ", clhstnam    ! filename
525         CALL histbeg( clhstnam, jpi, glamu, jpj, gphiu,           &  ! Horizontal grid: glamu and gphiu
526            &          iimi, iima-iimi+1, ijmi, ijma-ijmi+1,       &
527            &          nit000-1, zjulian, rn_Dt, nh_U, nid_U, domain_id=nidom, snc4chunks=snc4set )
528         CALL histvert( nid_U, "depthu", "Vertical U levels",      &  ! Vertical grid: gdept
529            &           "m", ipk, gdept_1d, nz_U, "down" )
530         !                                                            ! Index of ocean points
531         CALL wheneq( jpi*jpj*ipk, umask, 1, 1., ndex_U , ndim_U  )      ! volume
532         CALL wheneq( jpi*jpj    , umask, 1, 1., ndex_hU, ndim_hU )      ! surface
533
534         ! Define the V grid FILE ( nid_V )
535
536         CALL dia_nam( clhstnam, nn_write, 'grid_V' )                   ! filename
537         IF(lwp) WRITE(numout,*) " Name of NETCDF file ", clhstnam
538         CALL histbeg( clhstnam, jpi, glamv, jpj, gphiv,           &  ! Horizontal grid: glamv and gphiv
539            &          iimi, iima-iimi+1, ijmi, ijma-ijmi+1,       &
540            &          nit000-1, zjulian, rn_Dt, nh_V, nid_V, domain_id=nidom, snc4chunks=snc4set )
541         CALL histvert( nid_V, "depthv", "Vertical V levels",      &  ! Vertical grid : gdept
542            &          "m", ipk, gdept_1d, nz_V, "down" )
543         !                                                            ! Index of ocean points
544         CALL wheneq( jpi*jpj*ipk, vmask, 1, 1., ndex_V , ndim_V  )      ! volume
545         CALL wheneq( jpi*jpj    , vmask, 1, 1., ndex_hV, ndim_hV )      ! surface
546
547         ! Define the W grid FILE ( nid_W )
548
549         CALL dia_nam( clhstnam, nn_write, 'grid_W' )                   ! filename
550         IF(lwp) WRITE(numout,*) " Name of NETCDF file ", clhstnam
551         CALL histbeg( clhstnam, jpi, glamt, jpj, gphit,           &  ! Horizontal grid: glamt and gphit
552            &          iimi, iima-iimi+1, ijmi, ijma-ijmi+1,       &
553            &          nit000-1, zjulian, rn_Dt, nh_W, nid_W, domain_id=nidom, snc4chunks=snc4set )
554         CALL histvert( nid_W, "depthw", "Vertical W levels",      &  ! Vertical grid: gdepw
555            &          "m", ipk, gdepw_1d, nz_W, "down" )
556
557         IF( ln_abl ) THEN 
558         ! Define the ABL grid FILE ( nid_A )
559            CALL dia_nam( clhstnam, nn_write, 'grid_ABL' )
560            IF(lwp) WRITE(numout,*) " Name of NETCDF file ", clhstnam    ! filename
561            CALL histbeg( clhstnam, jpi, glamt, jpj, gphit,           &  ! Horizontal grid: glamt and gphit
562               &          iimi, iima-iimi+1, ijmi, ijma-ijmi+1,       &
563               &          nit000-1, zjulian, rn_Dt, nh_A, nid_A, domain_id=nidom, snc4chunks=snc4set )
564            CALL histvert( nid_A, "ght_abl", "Vertical T levels",      &  ! Vertical grid: gdept
565               &           "m", ipka, ght_abl(2:jpka), nz_A, "up" )
566            !                                                            ! Index of ocean points
567         ALLOCATE( zw3d_abl(jpi,jpj,ipka) ) 
568         zw3d_abl(:,:,:) = 1._wp 
569         CALL wheneq( jpi*jpj*ipka, zw3d_abl, 1, 1., ndex_A , ndim_A  )      ! volume
570            CALL wheneq( jpi*jpj     , zw3d_abl, 1, 1., ndex_hA, ndim_hA )      ! surface
571         DEALLOCATE(zw3d_abl)
572         ENDIF
573
574         ! Declare all the output fields as NETCDF variables
575
576         !                                                                                      !!! nid_T : 3D
577         CALL histdef( nid_T, "votemper", "Temperature"                        , "C"      ,   &  ! tn
578            &          jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout )
579         CALL histdef( nid_T, "vosaline", "Salinity"                           , "PSU"    ,   &  ! sn
580            &          jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout )
581         IF(  .NOT.ln_linssh  ) THEN
582            CALL histdef( nid_T, "vovvle3t", "Level thickness"                    , "m"      ,&  ! e3t(:,:,:,Kmm)
583            &             jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout )
584            CALL histdef( nid_T, "vovvldep", "T point depth"                      , "m"      ,&  ! e3t(:,:,:,Kmm)
585            &             jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout )
586            CALL histdef( nid_T, "vovvldef", "Squared level deformation"          , "%^2"    ,&  ! e3t(:,:,:,Kmm)
587            &             jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout )
588         ENDIF
589         !                                                                                      !!! nid_T : 2D
590         CALL histdef( nid_T, "sosstsst", "Sea Surface temperature"            , "C"      ,   &  ! sst
591            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
592         CALL histdef( nid_T, "sosaline", "Sea Surface Salinity"               , "PSU"    ,   &  ! sss
593            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
594         CALL histdef( nid_T, "sossheig", "Sea Surface Height"                 , "m"      ,   &  ! ssh
595            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
596         CALL histdef( nid_T, "sowaflup", "Net Upward Water Flux"              , "Kg/m2/s",   &  ! (emp-rnf)
597            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
598         CALL histdef( nid_T, "sorunoff", "River runoffs"                      , "Kg/m2/s",   &  ! runoffs
599            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
600         CALL histdef( nid_T, "sosfldow", "downward salt flux"                 , "PSU/m2/s",  &  ! sfx
601            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
602         IF(  ln_linssh  ) THEN
603            CALL histdef( nid_T, "sosst_cd", "Concentration/Dilution term on temperature"     &  ! emp * ts(:,:,1,jp_tem,Kmm)
604            &                                                                  , "KgC/m2/s",  &  ! sosst_cd
605            &             jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
606            CALL histdef( nid_T, "sosss_cd", "Concentration/Dilution term on salinity"        &  ! emp * ts(:,:,1,jp_sal,Kmm)
607            &                                                                  , "KgPSU/m2/s",&  ! sosss_cd
608            &             jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
609         ENDIF
610         CALL histdef( nid_T, "sohefldo", "Net Downward Heat Flux"             , "W/m2"   ,   &  ! qns + qsr
611            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
612         CALL histdef( nid_T, "soshfldo", "Shortwave Radiation"                , "W/m2"   ,   &  ! qsr
613            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
614         CALL histdef( nid_T, "somixhgt", "Turbocline Depth"                   , "m"      ,   &  ! hmld
615            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
616         CALL histdef( nid_T, "somxl010", "Mixed Layer Depth 0.01"             , "m"      ,   &  ! hmlp
617            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
618         CALL histdef( nid_T, "soicecov", "Ice fraction"                       , "[0,1]"  ,   &  ! fr_i
619            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
620         CALL histdef( nid_T, "sowindsp", "wind speed at 10m"                  , "m/s"    ,   &  ! wndm
621            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
622         !
623         IF( ln_abl ) THEN
624            CALL histdef( nid_A, "t_abl", "Potential Temperature"     , "K"        ,       &  ! t_abl
625               &          jpi, jpj, nh_A, ipka, 1, ipka, nz_A, 32, clop, zsto, zout )
626            CALL histdef( nid_A, "q_abl", "Humidity"                  , "kg/kg"    ,       &  ! q_abl
627               &          jpi, jpj, nh_A, ipka, 1, ipka, nz_A, 32, clop, zsto, zout ) 
628            CALL histdef( nid_A, "u_abl", "Atmospheric U-wind   "     , "m/s"        ,     &  ! u_abl
629               &          jpi, jpj, nh_A, ipka, 1, ipka, nz_A, 32, clop, zsto, zout )
630            CALL histdef( nid_A, "v_abl", "Atmospheric V-wind   "     , "m/s"    ,         &  ! v_abl
631               &          jpi, jpj, nh_A, ipka, 1, ipka, nz_A, 32, clop, zsto, zout ) 
632            CALL histdef( nid_A, "tke_abl", "Atmospheric TKE   "     , "m2/s2"    ,        &  ! tke_abl
633               &          jpi, jpj, nh_A, ipka, 1, ipka, nz_A, 32, clop, zsto, zout ) 
634            CALL histdef( nid_A, "avm_abl", "Atmospheric turbulent viscosity", "m2/s"   ,  &  ! avm_abl
635               &          jpi, jpj, nh_A, ipka, 1, ipka, nz_A, 32, clop, zsto, zout ) 
636            CALL histdef( nid_A, "avt_abl", "Atmospheric turbulent diffusivity", "m2/s2",  &  ! avt_abl
637               &          jpi, jpj, nh_A, ipka, 1, ipka, nz_A, 32, clop, zsto, zout ) 
638            CALL histdef( nid_A, "pblh", "Atmospheric boundary layer height "  , "m",      &  ! pblh
639               &          jpi, jpj, nh_A,  1  , 1, 1   , -99 , 32, clop, zsto, zout )                 
640#if defined key_si3
641            CALL histdef( nid_A, "oce_frac", "Fraction of open ocean"  , " ",      &  ! ato_i
642               &          jpi, jpj, nh_A,  1  , 1, 1   , -99 , 32, clop, zsto, zout )
643#endif
644            CALL histend( nid_A, snc4chunks=snc4set )
645         ENDIF
646         !
647         IF( ln_icebergs ) THEN
648            CALL histdef( nid_T, "calving"             , "calving mass input"                       , "kg/s"   , &
649               &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
650            CALL histdef( nid_T, "calving_heat"        , "calving heat flux"                        , "XXXX"   , &
651               &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
652            CALL histdef( nid_T, "berg_floating_melt"  , "Melt rate of icebergs + bits"             , "kg/m2/s", &
653               &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
654            CALL histdef( nid_T, "berg_stored_ice"     , "Accumulated ice mass by class"            , "kg"     , &
655               &          jpi, jpj, nh_T, nclasses  , 1, nclasses  , nb_T , 32, clop, zsto, zout )
656            IF( ln_bergdia ) THEN
657               CALL histdef( nid_T, "berg_melt"           , "Melt rate of icebergs"                    , "kg/m2/s", &
658                  &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
659               CALL histdef( nid_T, "berg_buoy_melt"      , "Buoyancy component of iceberg melt rate"  , "kg/m2/s", &
660                  &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
661               CALL histdef( nid_T, "berg_eros_melt"      , "Erosion component of iceberg melt rate"   , "kg/m2/s", &
662                  &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
663               CALL histdef( nid_T, "berg_conv_melt"      , "Convective component of iceberg melt rate", "kg/m2/s", &
664                  &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
665               CALL histdef( nid_T, "berg_virtual_area"   , "Virtual coverage by icebergs"             , "m2"     , &
666                  &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
667               CALL histdef( nid_T, "bits_src"           , "Mass source of bergy bits"                , "kg/m2/s", &
668                  &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
669               CALL histdef( nid_T, "bits_melt"          , "Melt rate of bergy bits"                  , "kg/m2/s", &
670                  &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
671               CALL histdef( nid_T, "bits_mass"          , "Bergy bit density field"                  , "kg/m2"  , &
672                  &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
673               CALL histdef( nid_T, "berg_mass"           , "Iceberg density field"                    , "kg/m2"  , &
674                  &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
675               CALL histdef( nid_T, "berg_real_calving"   , "Calving into iceberg class"               , "kg/s"   , &
676                  &          jpi, jpj, nh_T, nclasses  , 1, nclasses  , nb_T , 32, clop, zsto, zout )
677            ENDIF
678         ENDIF
679
680         IF( ln_ssr ) THEN
681            CALL histdef( nid_T, "sohefldp", "Surface Heat Flux: Damping"         , "W/m2"   ,   &  ! qrp
682               &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
683            CALL histdef( nid_T, "sowafldp", "Surface Water Flux: Damping"        , "Kg/m2/s",   &  ! erp
684               &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
685            CALL histdef( nid_T, "sosafldp", "Surface salt flux: damping"         , "Kg/m2/s",   &  ! erp * sn
686               &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
687         ENDIF
688       
689         clmx ="l_max(only(x))"    ! max index on a period
690!         CALL histdef( nid_T, "sobowlin", "Bowl Index"                         , "W-point",   &  ! bowl INDEX
691!            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clmx, zsto, zout )
692#if defined key_diahth
693         CALL histdef( nid_T, "sothedep", "Thermocline Depth"                  , "m"      ,   & ! hth
694            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
695         CALL histdef( nid_T, "so20chgt", "Depth of 20C isotherm"              , "m"      ,   & ! hd20
696            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
697         CALL histdef( nid_T, "so28chgt", "Depth of 28C isotherm"              , "m"      ,   & ! hd28
698            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
699         CALL histdef( nid_T, "sohtc300", "Heat content 300 m"                 , "J/m2"   ,   & ! htc3
700            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
701#endif
702
703         CALL histend( nid_T, snc4chunks=snc4set )
704
705         !                                                                                      !!! nid_U : 3D
706         CALL histdef( nid_U, "vozocrtx", "Zonal Current"                      , "m/s"    ,   &  ! uu(:,:,:,Kmm)
707            &          jpi, jpj, nh_U, ipk, 1, ipk, nz_U, 32, clop, zsto, zout )
708         IF( ln_wave .AND. ln_sdw) THEN
709            CALL histdef( nid_U, "sdzocrtx", "Stokes Drift Zonal Current"         , "m/s"    ,   &  ! usd
710               &          jpi, jpj, nh_U, ipk, 1, ipk, nz_U, 32, clop, zsto, zout )
711         ENDIF
712         !                                                                                      !!! nid_U : 2D
713         CALL histdef( nid_U, "sozotaux", "Wind Stress along i-axis"           , "N/m2"   ,   &  ! utau
714            &          jpi, jpj, nh_U, 1  , 1, 1  , - 99, 32, clop, zsto, zout )
715
716         CALL histend( nid_U, snc4chunks=snc4set )
717
718         !                                                                                      !!! nid_V : 3D
719         CALL histdef( nid_V, "vomecrty", "Meridional Current"                 , "m/s"    ,   &  ! vv(:,:,:,Kmm)
720            &          jpi, jpj, nh_V, ipk, 1, ipk, nz_V, 32, clop, zsto, zout )
721         IF( ln_wave .AND. ln_sdw) THEN
722            CALL histdef( nid_V, "sdmecrty", "Stokes Drift Meridional Current"    , "m/s"    ,   &  ! vsd
723               &          jpi, jpj, nh_V, ipk, 1, ipk, nz_V, 32, clop, zsto, zout )
724         ENDIF
725         !                                                                                      !!! nid_V : 2D
726         CALL histdef( nid_V, "sometauy", "Wind Stress along j-axis"           , "N/m2"   ,   &  ! vtau
727            &          jpi, jpj, nh_V, 1  , 1, 1  , - 99, 32, clop, zsto, zout )
728
729         CALL histend( nid_V, snc4chunks=snc4set )
730
731         !                                                                                      !!! nid_W : 3D
732         CALL histdef( nid_W, "vovecrtz", "Vertical Velocity"                  , "m/s"    ,   &  ! ww
733            &          jpi, jpj, nh_W, ipk, 1, ipk, nz_W, 32, clop, zsto, zout )
734         CALL histdef( nid_W, "votkeavt", "Vertical Eddy Diffusivity"          , "m2/s"   ,   &  ! avt
735            &          jpi, jpj, nh_W, ipk, 1, ipk, nz_W, 32, clop, zsto, zout )
736         CALL histdef( nid_W, "votkeavm", "Vertical Eddy Viscosity"             , "m2/s"  ,   &  ! avm
737            &          jpi, jpj, nh_W, ipk, 1, ipk, nz_W, 32, clop, zsto, zout )
738
739         IF( ln_zdfddm ) THEN
740            CALL histdef( nid_W,"voddmavs","Salt Vertical Eddy Diffusivity"    , "m2/s"   ,   &  ! avs
741               &          jpi, jpj, nh_W, ipk, 1, ipk, nz_W, 32, clop, zsto, zout )
742         ENDIF
743         
744         IF( ln_wave .AND. ln_sdw) THEN
745            CALL histdef( nid_W, "sdvecrtz", "Stokes Drift Vertical Current"   , "m/s"    ,   &  ! wsd
746               &          jpi, jpj, nh_W, ipk, 1, ipk, nz_W, 32, clop, zsto, zout )
747         ENDIF
748         !                                                                                      !!! nid_W : 2D
749         CALL histend( nid_W, snc4chunks=snc4set )
750
751         IF(lwp) WRITE(numout,*)
752         IF(lwp) WRITE(numout,*) 'End of NetCDF Initialization'
753         IF(ll_print) CALL FLUSH(numout )
754
755      ENDIF
756
757      ! 2. Start writing data
758      ! ---------------------
759
760      ! ndex(1) est utilise ssi l'avant dernier argument est different de
761      ! la taille du tableau en sortie. Dans ce cas , l'avant dernier argument
762      ! donne le nombre d'elements, et ndex la liste des indices a sortir
763
764      IF( lwp .AND. MOD( itmod, nn_write ) == 0 ) THEN
765         WRITE(numout,*) 'dia_wri : write model outputs in NetCDF files at ', kt, 'time-step'
766         WRITE(numout,*) '~~~~~~ '
767      ENDIF
768
769      IF( .NOT.ln_linssh ) THEN
770         CALL histwrite( nid_T, "votemper", it, ts(:,:,:,jp_tem,Kmm) * e3t(:,:,:,Kmm) , ndim_T , ndex_T  )   ! heat content
771         CALL histwrite( nid_T, "vosaline", it, ts(:,:,:,jp_sal,Kmm) * e3t(:,:,:,Kmm) , ndim_T , ndex_T  )   ! salt content
772         CALL histwrite( nid_T, "sosstsst", it, ts(:,:,1,jp_tem,Kmm) * e3t(:,:,1,Kmm) , ndim_hT, ndex_hT )   ! sea surface heat content
773         CALL histwrite( nid_T, "sosaline", it, ts(:,:,1,jp_sal,Kmm) * e3t(:,:,1,Kmm) , ndim_hT, ndex_hT )   ! sea surface salinity content
774      ELSE
775         CALL histwrite( nid_T, "votemper", it, ts(:,:,:,jp_tem,Kmm) , ndim_T , ndex_T  )   ! temperature
776         CALL histwrite( nid_T, "vosaline", it, ts(:,:,:,jp_sal,Kmm) , ndim_T , ndex_T  )   ! salinity
777         CALL histwrite( nid_T, "sosstsst", it, ts(:,:,1,jp_tem,Kmm) , ndim_hT, ndex_hT )   ! sea surface temperature
778         CALL histwrite( nid_T, "sosaline", it, ts(:,:,1,jp_sal,Kmm) , ndim_hT, ndex_hT )   ! sea surface salinity
779      ENDIF
780      IF( .NOT.ln_linssh ) THEN
781         zw3d(:,:,:) = ( ( e3t(:,:,:,Kmm) - e3t_0(:,:,:) ) / e3t_0(:,:,:) * 100 * tmask(:,:,:) ) ** 2
782         CALL histwrite( nid_T, "vovvle3t", it, e3t (:,:,:,Kmm) , ndim_T , ndex_T  )   ! level thickness
783         CALL histwrite( nid_T, "vovvldep", it, gdept(:,:,:,Kmm) , ndim_T , ndex_T  )   ! t-point depth
784         CALL histwrite( nid_T, "vovvldef", it, zw3d             , ndim_T , ndex_T  )   ! level thickness deformation
785      ENDIF
786      CALL histwrite( nid_T, "sossheig", it, ssh(:,:,Kmm)          , ndim_hT, ndex_hT )   ! sea surface height
787      CALL histwrite( nid_T, "sowaflup", it, ( emp-rnf )   , ndim_hT, ndex_hT )   ! upward water flux
788      CALL histwrite( nid_T, "sorunoff", it, rnf           , ndim_hT, ndex_hT )   ! river runoffs
789      CALL histwrite( nid_T, "sosfldow", it, sfx           , ndim_hT, ndex_hT )   ! downward salt flux
790                                                                                  ! (includes virtual salt flux beneath ice
791                                                                                  ! in linear free surface case)
792      IF( ln_linssh ) THEN
793         zw2d(:,:) = emp (:,:) * ts(:,:,1,jp_tem,Kmm)
794         CALL histwrite( nid_T, "sosst_cd", it, zw2d, ndim_hT, ndex_hT )          ! c/d term on sst
795         zw2d(:,:) = emp (:,:) * ts(:,:,1,jp_sal,Kmm)
796         CALL histwrite( nid_T, "sosss_cd", it, zw2d, ndim_hT, ndex_hT )          ! c/d term on sss
797      ENDIF
798      CALL histwrite( nid_T, "sohefldo", it, qns + qsr     , ndim_hT, ndex_hT )   ! total heat flux
799      CALL histwrite( nid_T, "soshfldo", it, qsr           , ndim_hT, ndex_hT )   ! solar heat flux
800      CALL histwrite( nid_T, "somixhgt", it, hmld          , ndim_hT, ndex_hT )   ! turbocline depth
801      CALL histwrite( nid_T, "somxl010", it, hmlp          , ndim_hT, ndex_hT )   ! mixed layer depth
802      CALL histwrite( nid_T, "soicecov", it, fr_i          , ndim_hT, ndex_hT )   ! ice fraction   
803      CALL histwrite( nid_T, "sowindsp", it, wndm          , ndim_hT, ndex_hT )   ! wind speed   
804      !
805      IF( ln_abl ) THEN
806         ALLOCATE( zw3d_abl(jpi,jpj,jpka) )
807         IF( ln_mskland )   THEN
808            DO jk=1,jpka
809               zw3d_abl(:,:,jk) = tmask(:,:,1)
810            END DO       
811         ELSE
812            zw3d_abl(:,:,:) = 1._wp     
813         ENDIF       
814         CALL histwrite( nid_A,  "pblh"   , it, pblh(:,:)                  *zw3d_abl(:,:,1     ), ndim_hA, ndex_hA )   ! pblh
815         CALL histwrite( nid_A,  "u_abl"  , it, u_abl   (:,:,2:jpka,nt_n  )*zw3d_abl(:,:,2:jpka), ndim_A , ndex_A  )   ! u_abl
816         CALL histwrite( nid_A,  "v_abl"  , it, v_abl   (:,:,2:jpka,nt_n  )*zw3d_abl(:,:,2:jpka), ndim_A , ndex_A  )   ! v_abl
817         CALL histwrite( nid_A,  "t_abl"  , it, tq_abl  (:,:,2:jpka,nt_n,1)*zw3d_abl(:,:,2:jpka), ndim_A , ndex_A  )   ! t_abl
818         CALL histwrite( nid_A,  "q_abl"  , it, tq_abl  (:,:,2:jpka,nt_n,2)*zw3d_abl(:,:,2:jpka), ndim_A , ndex_A  )   ! q_abl       
819         CALL histwrite( nid_A,  "tke_abl", it, tke_abl (:,:,2:jpka,nt_n  )*zw3d_abl(:,:,2:jpka), ndim_A , ndex_A  )   ! tke_abl
820         CALL histwrite( nid_A,  "avm_abl", it, avm_abl (:,:,2:jpka       )*zw3d_abl(:,:,2:jpka), ndim_A , ndex_A  )   ! avm_abl
821         CALL histwrite( nid_A,  "avt_abl", it, avt_abl (:,:,2:jpka       )*zw3d_abl(:,:,2:jpka), ndim_A , ndex_A  )   ! avt_abl
822#if defined key_si3
823         CALL histwrite( nid_A,  "oce_frac"   , it, ato_i(:,:)                                  , ndim_hA, ndex_hA )   ! ato_i
824#endif
825         DEALLOCATE(zw3d_abl)
826      ENDIF
827      !
828      IF( ln_icebergs ) THEN
829         !
830         CALL histwrite( nid_T, "calving"             , it, berg_grid%calving      , ndim_hT, ndex_hT ) 
831         CALL histwrite( nid_T, "calving_heat"        , it, berg_grid%calving_hflx , ndim_hT, ndex_hT )         
832         CALL histwrite( nid_T, "berg_floating_melt"  , it, berg_grid%floating_melt, ndim_hT, ndex_hT ) 
833         !
834         CALL histwrite( nid_T, "berg_stored_ice"     , it, berg_grid%stored_ice   , ndim_bT, ndex_bT )
835         !
836         IF( ln_bergdia ) THEN
837            CALL histwrite( nid_T, "berg_melt"           , it, berg_melt        , ndim_hT, ndex_hT   ) 
838            CALL histwrite( nid_T, "berg_buoy_melt"      , it, buoy_melt        , ndim_hT, ndex_hT   ) 
839            CALL histwrite( nid_T, "berg_eros_melt"      , it, eros_melt        , ndim_hT, ndex_hT   ) 
840            CALL histwrite( nid_T, "berg_conv_melt"      , it, conv_melt        , ndim_hT, ndex_hT   ) 
841            CALL histwrite( nid_T, "berg_virtual_area"   , it, virtual_area     , ndim_hT, ndex_hT   ) 
842            CALL histwrite( nid_T, "bits_src"            , it, bits_src         , ndim_hT, ndex_hT   ) 
843            CALL histwrite( nid_T, "bits_melt"           , it, bits_melt        , ndim_hT, ndex_hT   ) 
844            CALL histwrite( nid_T, "bits_mass"           , it, bits_mass        , ndim_hT, ndex_hT   ) 
845            CALL histwrite( nid_T, "berg_mass"           , it, berg_mass        , ndim_hT, ndex_hT   ) 
846            !
847            CALL histwrite( nid_T, "berg_real_calving"   , it, real_calving     , ndim_bT, ndex_bT   )
848         ENDIF
849      ENDIF
850
851      IF( ln_ssr ) THEN
852         CALL histwrite( nid_T, "sohefldp", it, qrp           , ndim_hT, ndex_hT )   ! heat flux damping
853         CALL histwrite( nid_T, "sowafldp", it, erp           , ndim_hT, ndex_hT )   ! freshwater flux damping
854         zw2d(:,:) = erp(:,:) * ts(:,:,1,jp_sal,Kmm) * tmask(:,:,1)
855         CALL histwrite( nid_T, "sosafldp", it, zw2d          , ndim_hT, ndex_hT )   ! salt flux damping
856      ENDIF
857!      zw2d(:,:) = FLOAT( nmln(:,:) ) * tmask(:,:,1)
858!      CALL histwrite( nid_T, "sobowlin", it, zw2d          , ndim_hT, ndex_hT )   ! ???
859
860#if defined key_diahth
861      CALL histwrite( nid_T, "sothedep", it, hth           , ndim_hT, ndex_hT )   ! depth of the thermocline
862      CALL histwrite( nid_T, "so20chgt", it, hd20          , ndim_hT, ndex_hT )   ! depth of the 20 isotherm
863      CALL histwrite( nid_T, "so28chgt", it, hd28          , ndim_hT, ndex_hT )   ! depth of the 28 isotherm
864      CALL histwrite( nid_T, "sohtc300", it, htc3          , ndim_hT, ndex_hT )   ! first 300m heaat content
865#endif
866
867      CALL histwrite( nid_U, "vozocrtx", it, uu(:,:,:,Kmm)            , ndim_U , ndex_U )    ! i-current
868      CALL histwrite( nid_U, "sozotaux", it, utau          , ndim_hU, ndex_hU )   ! i-wind stress
869
870      CALL histwrite( nid_V, "vomecrty", it, vv(:,:,:,Kmm)            , ndim_V , ndex_V  )   ! j-current
871      CALL histwrite( nid_V, "sometauy", it, vtau          , ndim_hV, ndex_hV )   ! j-wind stress
872
873      IF( ln_zad_Aimp ) THEN
874         CALL histwrite( nid_W, "vovecrtz", it, ww + wi     , ndim_T, ndex_T )    ! vert. current
875      ELSE
876         CALL histwrite( nid_W, "vovecrtz", it, ww          , ndim_T, ndex_T )    ! vert. current
877      ENDIF
878      CALL histwrite( nid_W, "votkeavt", it, avt            , ndim_T, ndex_T )    ! T vert. eddy diff. coef.
879      CALL histwrite( nid_W, "votkeavm", it, avm            , ndim_T, ndex_T )    ! T vert. eddy visc. coef.
880      IF( ln_zdfddm ) THEN
881         CALL histwrite( nid_W, "voddmavs", it, avs         , ndim_T, ndex_T )    ! S vert. eddy diff. coef.
882      ENDIF
883
884      IF( ln_wave .AND. ln_sdw ) THEN
885         CALL histwrite( nid_U, "sdzocrtx", it, usd         , ndim_U , ndex_U )    ! i-StokesDrift-current
886         CALL histwrite( nid_V, "sdmecrty", it, vsd         , ndim_V , ndex_V )    ! j-StokesDrift-current
887         CALL histwrite( nid_W, "sdvecrtz", it, wsd         , ndim_T , ndex_T )    ! StokesDrift vert. current
888      ENDIF
889
890      ! 3. Close all files
891      ! ---------------------------------------
892      IF( kt == nitend ) THEN
893         CALL histclo( nid_T )
894         CALL histclo( nid_U )
895         CALL histclo( nid_V )
896         CALL histclo( nid_W )
897         IF(ln_abl) CALL histclo( nid_A )
898      ENDIF
899      !
900      IF( ln_timing )   CALL timing_stop('dia_wri')
901      !
902   END SUBROUTINE dia_wri
903#endif
904
905   SUBROUTINE dia_wri_state( Kmm, cdfile_name )
906      !!---------------------------------------------------------------------
907      !!                 ***  ROUTINE dia_wri_state  ***
908      !!       
909      !! ** Purpose :   create a NetCDF file named cdfile_name which contains
910      !!      the instantaneous ocean state and forcing fields.
911      !!        Used to find errors in the initial state or save the last
912      !!      ocean state in case of abnormal end of a simulation
913      !!
914      !! ** Method  :   NetCDF files using ioipsl
915      !!      File 'output.init.nc'  is created if ninist = 1 (namelist)
916      !!      File 'output.abort.nc' is created in case of abnormal job end
917      !!----------------------------------------------------------------------
918      INTEGER           , INTENT( in ) ::   Kmm              ! time level index
919      CHARACTER (len=* ), INTENT( in ) ::   cdfile_name      ! name of the file created
920      !!
921      INTEGER :: inum, jk
922      !!----------------------------------------------------------------------
923      !
924      IF(lwp) WRITE(numout,*)
925      IF(lwp) WRITE(numout,*) 'dia_wri_state : single instantaneous ocean state'
926      IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~   and forcing fields file created '
927      IF(lwp) WRITE(numout,*) '                and named :', cdfile_name, '...nc'
928      !
929      CALL iom_open( TRIM(cdfile_name), inum, ldwrt = .TRUE. )
930      !
931      CALL iom_rstput( 0, 0, inum, 'votemper', ts(:,:,:,jp_tem,Kmm) )    ! now temperature
932      CALL iom_rstput( 0, 0, inum, 'vosaline', ts(:,:,:,jp_sal,Kmm) )    ! now salinity
933      CALL iom_rstput( 0, 0, inum, 'sossheig', ssh(:,:,Kmm)              )    ! sea surface height
934      CALL iom_rstput( 0, 0, inum, 'vozocrtx', uu(:,:,:,Kmm)                )    ! now i-velocity
935      CALL iom_rstput( 0, 0, inum, 'vomecrty', vv(:,:,:,Kmm)                )    ! now j-velocity
936      IF( ln_zad_Aimp ) THEN
937         CALL iom_rstput( 0, 0, inum, 'vovecrtz', ww + wi        )    ! now k-velocity
938      ELSE
939         CALL iom_rstput( 0, 0, inum, 'vovecrtz', ww             )    ! now k-velocity
940      ENDIF
941      CALL iom_rstput( 0, 0, inum, 'risfdep', risfdep            )    ! now k-velocity
942      CALL iom_rstput( 0, 0, inum, 'ht'     , ht                 )    ! now water column height
943      !
944      IF ( ln_isf ) THEN
945         IF (ln_isfcav_mlt) THEN
946            CALL iom_rstput( 0, 0, inum, 'fwfisf_cav', fwfisf_cav          )    ! now k-velocity
947            CALL iom_rstput( 0, 0, inum, 'rhisf_cav_tbl', rhisf_tbl_cav    )    ! now k-velocity
948            CALL iom_rstput( 0, 0, inum, 'rfrac_cav_tbl', rfrac_tbl_cav    )    ! now k-velocity
949            CALL iom_rstput( 0, 0, inum, 'misfkb_cav', REAL(misfkb_cav,wp) )    ! now k-velocity
950            CALL iom_rstput( 0, 0, inum, 'misfkt_cav', REAL(misfkt_cav,wp) )    ! now k-velocity
951            CALL iom_rstput( 0, 0, inum, 'mskisf_cav', REAL(mskisf_cav,wp), ktype = jp_i1 )
952         END IF
953         IF (ln_isfpar_mlt) THEN
954            CALL iom_rstput( 0, 0, inum, 'isfmsk_par', REAL(mskisf_par,wp) )    ! now k-velocity
955            CALL iom_rstput( 0, 0, inum, 'fwfisf_par', fwfisf_par          )    ! now k-velocity
956            CALL iom_rstput( 0, 0, inum, 'rhisf_par_tbl', rhisf_tbl_par    )    ! now k-velocity
957            CALL iom_rstput( 0, 0, inum, 'rfrac_par_tbl', rfrac_tbl_par    )    ! now k-velocity
958            CALL iom_rstput( 0, 0, inum, 'misfkb_par', REAL(misfkb_par,wp) )    ! now k-velocity
959            CALL iom_rstput( 0, 0, inum, 'misfkt_par', REAL(misfkt_par,wp) )    ! now k-velocity
960            CALL iom_rstput( 0, 0, inum, 'mskisf_par', REAL(mskisf_par,wp), ktype = jp_i1 )
961         END IF
962      END IF
963      !
964      IF( ALLOCATED(ahtu) ) THEN
965         CALL iom_rstput( 0, 0, inum,  'ahtu', ahtu              )    ! aht at u-point
966         CALL iom_rstput( 0, 0, inum,  'ahtv', ahtv              )    ! aht at v-point
967      ENDIF
968      IF( ALLOCATED(ahmt) ) THEN
969         CALL iom_rstput( 0, 0, inum,  'ahmt', ahmt              )    ! ahmt at u-point
970         CALL iom_rstput( 0, 0, inum,  'ahmf', ahmf              )    ! ahmf at v-point
971      ENDIF
972      CALL iom_rstput( 0, 0, inum, 'sowaflup', emp - rnf         )    ! freshwater budget
973      CALL iom_rstput( 0, 0, inum, 'sohefldo', qsr + qns         )    ! total heat flux
974      CALL iom_rstput( 0, 0, inum, 'soshfldo', qsr               )    ! solar heat flux
975      CALL iom_rstput( 0, 0, inum, 'soicecov', fr_i              )    ! ice fraction
976      CALL iom_rstput( 0, 0, inum, 'sozotaux', utau              )    ! i-wind stress
977      CALL iom_rstput( 0, 0, inum, 'sometauy', vtau              )    ! j-wind stress
978      IF(  .NOT.ln_linssh  ) THEN             
979         CALL iom_rstput( 0, 0, inum, 'vovvldep', gdept(:,:,:,Kmm)        )    !  T-cell depth
980         CALL iom_rstput( 0, 0, inum, 'vovvle3t', e3t(:,:,:,Kmm)          )    !  T-cell thickness 
981      END IF
982      IF( ln_wave .AND. ln_sdw ) THEN
983         CALL iom_rstput( 0, 0, inum, 'sdzocrtx', usd            )    ! now StokesDrift i-velocity
984         CALL iom_rstput( 0, 0, inum, 'sdmecrty', vsd            )    ! now StokesDrift j-velocity
985         CALL iom_rstput( 0, 0, inum, 'sdvecrtz', wsd            )    ! now StokesDrift k-velocity
986      ENDIF
987      IF ( ln_abl ) THEN
988         CALL iom_rstput ( 0, 0, inum, "uz1_abl",   u_abl(:,:,2,nt_a  ) )   ! now first level i-wind
989         CALL iom_rstput ( 0, 0, inum, "vz1_abl",   v_abl(:,:,2,nt_a  ) )   ! now first level j-wind
990         CALL iom_rstput ( 0, 0, inum, "tz1_abl",  tq_abl(:,:,2,nt_a,1) )   ! now first level temperature
991         CALL iom_rstput ( 0, 0, inum, "qz1_abl",  tq_abl(:,:,2,nt_a,2) )   ! now first level humidity
992      ENDIF
993      !
994      CALL iom_close( inum )
995      !
996#if defined key_si3
997      IF( nn_ice == 2 ) THEN   ! condition needed in case agrif + ice-model but no-ice in child grid
998         CALL iom_open( TRIM(cdfile_name)//'_ice', inum, ldwrt = .TRUE., kdlev = jpl, cdcomp = 'ICE' )
999         CALL ice_wri_state( inum )
1000         CALL iom_close( inum )
1001      ENDIF
1002      !
1003#endif
1004   END SUBROUTINE dia_wri_state
1005
1006   !!======================================================================
1007END MODULE diawri
Note: See TracBrowser for help on using the repository browser.