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_r11265_ASINTER-01_Guillaume_ABL1D/src/OCE/DIA – NEMO

source: NEMO/branches/2019/dev_r11265_ASINTER-01_Guillaume_ABL1D/src/OCE/DIA/diawri.F90 @ 11413

Last change on this file since 11413 was 11413, checked in by gsamson, 5 years ago

dev_r11265_ABL : see #2131

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