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_r12563_ASINTER-06_ABL_improvement/src/SAS – NEMO

source: NEMO/branches/2020/dev_r12563_ASINTER-06_ABL_improvement/src/SAS/diawri.F90 @ 12808

Last change on this file since 12808 was 12489, checked in by davestorkey, 4 years ago

Preparation for new timestepping scheme #2390.
Main changes:

  1. Initial euler timestep now handled in stp and not in TRA/DYN routines.
  2. Renaming of all timestep parameters. In summary, the namelist parameter is now rn_Dt and the current timestep is rDt (and rDt_ice, rDt_trc etc).
  3. Renaming of a few miscellaneous parameters, eg. atfp -> rn_atfp (namelist parameter used everywhere) and rau0 -> rho0.

This version gives bit-comparable results to the previous version of the trunk.

  • Property svn:keywords set to Id
File size: 24.2 KB
Line 
1MODULE diawri
2   !!======================================================================
3   !!                     ***  MODULE  diawri  ***
4   !! Ocean diagnostics :  write ocean output files
5   !!=====================================================================
6   !! History :  OPA  ! 1991-03  (M.-A. Foujols)  Original code
7   !!            4.0  ! 1991-11  (G. Madec)
8   !!                 ! 1992-06  (M. Imbard)  correction restart file
9   !!                 ! 1992-07  (M. Imbard)  split into diawri and rstwri
10   !!                 ! 1993-03  (M. Imbard)  suppress writibm
11   !!                 ! 1998-01  (C. Levy)  NETCDF format using ioipsl INTERFACE
12   !!                 ! 1999-02  (E. Guilyardi)  name of netCDF files + variables
13   !!            8.2  ! 2000-06  (M. Imbard)  Original code (diabort.F)
14   !!   NEMO     1.0  ! 2002-06  (A.Bozec, E. Durand)  Original code (diainit.F)
15   !!             -   ! 2002-09  (G. Madec)  F90: Free form and module
16   !!             -   ! 2002-12  (G. Madec)  merge of diabort and diainit, F90
17   !!                 ! 2005-11  (V. Garnier) Surface pressure gradient organization
18   !!            3.2  ! 2008-11  (B. Lemaire) creation from old diawri
19   !!----------------------------------------------------------------------
20
21   !!----------------------------------------------------------------------
22   !!   dia_wri       : create the standart output files
23   !!   dia_wri_state : create an output NetCDF file for a single instantaeous ocean state and forcing fields
24   !!----------------------------------------------------------------------
25   USE oce             ! ocean dynamics and tracers
26   USE abl            ! abl variables in case ln_abl = .true.
27   USE dom_oce         ! ocean space and time domain
28   USE zdf_oce         ! ocean vertical physics
29   USE sbc_oce         ! Surface boundary condition: ocean fields
30   USE sbc_ice         ! Surface boundary condition: ice fields
31   USE sbcssr          ! restoring term toward SST/SSS climatology
32   USE phycst          ! physical constants
33   USE zdfmxl          ! mixed layer
34   USE dianam          ! build name of file (routine)
35   USE zdfddm          ! vertical  physics: double diffusion
36   USE diahth          ! thermocline diagnostics
37   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
38   USE in_out_manager  ! I/O manager
39   USE iom
40   USE ioipsl
41#if defined key_si3
42   USE ice
43   USE icewri
44#endif
45   USE lib_mpp         ! MPP library
46   USE timing          ! preformance summary
47
48   IMPLICIT NONE
49   PRIVATE
50
51   PUBLIC   dia_wri                 ! routines called by step.F90
52   PUBLIC   dia_wri_state
53   PUBLIC   dia_wri_alloc           ! Called by nemogcm module
54#if ! defined key_iomput   
55   PUBLIC   dia_wri_alloc_abl       ! Called by sbcabl  module (if ln_abl = .true.)
56#endif
57   INTEGER ::   nid_T, nz_T, nh_T, ndim_T, ndim_hT   ! grid_T file
58   INTEGER ::   nid_U, nz_U, nh_U, ndim_U, ndim_hU   ! grid_U file
59   INTEGER ::   nid_V, nz_V, nh_V, ndim_V, ndim_hV   ! grid_V file
60   INTEGER ::   ndim_A, ndim_hA                      ! ABL file   
61   INTEGER ::   nid_A, nz_A, nh_A                    ! grid_ABL file   
62   INTEGER ::   ndex(1)                              ! ???
63   INTEGER, SAVE, ALLOCATABLE, DIMENSION(:) :: ndex_hT, ndex_hU, ndex_hV
64   INTEGER, SAVE, ALLOCATABLE, DIMENSION(:) :: ndex_hA, ndex_A ! ABL
65
66   !!----------------------------------------------------------------------
67   !! NEMO/SAS 4.0 , NEMO Consortium (2018)
68   !! $Id$
69   !! Software governed by the CeCILL license (see ./LICENSE)
70   !!----------------------------------------------------------------------
71CONTAINS
72
73# if defined key_iomput
74   !!----------------------------------------------------------------------
75   !!   'key_iomput'                                        use IOM library
76   !!----------------------------------------------------------------------
77   INTEGER FUNCTION dia_wri_alloc()
78      !
79      dia_wri_alloc = 0
80      !
81   END FUNCTION dia_wri_alloc
82
83   
84   SUBROUTINE dia_wri( kt, Kmm )
85      !!---------------------------------------------------------------------
86      !!                  ***  ROUTINE dia_wri  ***
87      !!                   
88      !! ** Purpose :   Standard output of opa: dynamics and tracer fields
89      !!      NETCDF format is used by default
90      !!      Standalone surface scheme
91      !!
92      !! ** Method  :  use iom_put
93      !!----------------------------------------------------------------------
94      !!
95      INTEGER, INTENT( in ) ::   kt      ! ocean time-step index
96      INTEGER, INTENT( in ) ::   Kmm     ! ocean time levelindex
97      !!----------------------------------------------------------------------
98      !
99      ! Output the initial state and forcings
100      IF( ninist == 1 ) THEN
101         CALL dia_wri_state( 'output.init', Kmm )
102         ninist = 0
103      ENDIF
104      !
105   END SUBROUTINE dia_wri
106
107#else
108   !!----------------------------------------------------------------------
109   !!   Default option                                  use IOIPSL  library
110   !!----------------------------------------------------------------------
111   INTEGER FUNCTION dia_wri_alloc()
112      !!----------------------------------------------------------------------
113      INTEGER :: ierr
114      !!----------------------------------------------------------------------
115      !
116      ALLOCATE( ndex_hT(jpi*jpj), ndex_hU(jpi*jpj), ndex_hV(jpi*jpj), STAT=dia_wri_alloc )
117      CALL mpp_sum( 'diawri', dia_wri_alloc )
118      !
119   END FUNCTION dia_wri_alloc
120   
121   INTEGER FUNCTION dia_wri_alloc_abl()
122      !!----------------------------------------------------------------------
123     ALLOCATE(   ndex_hA(jpi*jpj), ndex_A (jpi*jpj*jpkam1), STAT=dia_wri_alloc_abl)
124      CALL mpp_sum( 'diawri', dia_wri_alloc_abl )
125      !
126   END FUNCTION dia_wri_alloc_abl
127 
128   SUBROUTINE dia_wri( kt )
129      !!---------------------------------------------------------------------
130      !!                  ***  ROUTINE dia_wri  ***
131      !!                   
132      !! ** Purpose :   Standard output of opa: dynamics and tracer fields
133      !!      NETCDF format is used by default
134      !!
135      !! ** Method  :   At the beginning of the first time step (nit000),
136      !!      define all the NETCDF files and fields
137      !!      At each time step call histdef to compute the mean if ncessary
138      !!      Each nn_write time step, output the instantaneous or mean fields
139      !!----------------------------------------------------------------------
140      !!
141      INTEGER, INTENT( in ) ::   kt      ! ocean time-step index
142      !!
143      LOGICAL ::   ll_print = .FALSE.                        ! =T print and flush numout
144      CHARACTER (len=40) ::   clhstnam, clop, clmx           ! local names
145      INTEGER  ::   inum = 11                                ! temporary logical unit
146      INTEGER  ::   ji, jj, jk                               ! dummy loop indices
147      INTEGER  ::   ierr                                     ! error code return from allocation
148      INTEGER  ::   iimi, iima, ipk, it, itmod, ijmi, ijma   ! local integers
149      INTEGER  ::   ipka                                     ! ABL
150      REAL(wp) ::   zsto, zout, zmax, zjulian                ! local scalars
151      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE :: zw3d_abl   ! ABL 3D workspace
152      !!----------------------------------------------------------------------
153      !
154      ! Output the initial state and forcings
155      IF( ninist == 1 ) THEN                       
156         CALL dia_wri_state( 'output.init' )
157         ninist = 0
158      ENDIF
159      !
160      IF( nn_write == -1 )   RETURN   ! we will never do any output
161      !
162      IF( ln_timing )   CALL timing_start('dia_wri')
163      !
164      ! 0. Initialisation
165      ! -----------------
166
167      ! local variable for debugging
168      ll_print = .FALSE.
169      ll_print = ll_print .AND. lwp
170
171      ! Define frequency of output and means
172      IF( ln_mskland )   THEN   ;   clop = "only(x)"   ! put 1.e+20 on land (very expensive!!)
173      ELSE                      ;   clop = "x"         ! no use of the mask value (require less cpu time)
174      ENDIF
175#if defined key_diainstant
176      zsto = nn_write * rn_Dt
177      clop = "inst("//TRIM(clop)//")"
178#else
179      zsto=rn_Dt
180      clop = "ave("//TRIM(clop)//")"
181#endif
182      zout = nn_write * rn_Dt
183      zmax = ( nitend - nit000 + 1 ) * rn_Dt
184
185      ! Define indices of the horizontal output zoom and vertical limit storage
186      iimi = 1      ;      iima = jpi
187      ijmi = 1      ;      ijma = jpj
188      ipk = jpk
189     IF(ln_abl) ipka = jpkam1
190
191      ! define time axis
192      it = kt
193      itmod = kt - nit000 + 1
194
195
196      ! 1. Define NETCDF files and fields at beginning of first time step
197      ! -----------------------------------------------------------------
198
199      IF( kt == nit000 ) THEN
200
201         ! Define the NETCDF files (one per grid)
202
203         ! Compute julian date from starting date of the run
204         CALL ymds2ju( nyear, nmonth, nday, rn_Dt, zjulian )
205         zjulian = zjulian - adatrj   !   set calendar origin to the beginning of the experiment
206         IF(lwp)WRITE(numout,*)
207         IF(lwp)WRITE(numout,*) 'Date 0 used :', nit000, ' YEAR ', nyear,   &
208            &                    ' MONTH ', nmonth, ' DAY ', nday, 'Julian day : ', zjulian
209         IF(lwp)WRITE(numout,*) ' indexes of zoom = ', iimi, iima, ijmi, ijma,   &
210                                 ' limit storage in depth = ', ipk
211
212         ! WRITE root name in date.file for use by postpro
213         IF(lwp) THEN
214            CALL dia_nam( clhstnam, nn_write,' ' )
215            CALL ctl_opn( inum, 'date.file', 'REPLACE', 'FORMATTED', 'SEQUENTIAL', -1, numout, lwp, narea )
216            WRITE(inum,*) clhstnam
217            CLOSE(inum)
218         ENDIF
219
220         ! Define the T grid FILE ( nid_T )
221
222         CALL dia_nam( clhstnam, nn_write, 'grid_T' )
223         IF(lwp) WRITE(numout,*) " Name of NETCDF file ", clhstnam    ! filename
224         CALL histbeg( clhstnam, jpi, glamt, jpj, gphit,           &  ! Horizontal grid: glamt and gphit
225            &          iimi, iima-iimi+1, ijmi, ijma-ijmi+1,       &
226            &          nit000-1, zjulian, rn_Dt, nh_T, nid_T, domain_id=nidom, snc4chunks=snc4set )
227         CALL histvert( nid_T, "deptht", "Vertical T levels",      &  ! Vertical grid: gdept
228            &           "m", ipk, gdept_1d, nz_T, "down" )
229         !                                                            ! Index of ocean points
230         CALL wheneq( jpi*jpj    , tmask, 1, 1., ndex_hT, ndim_hT )      ! surface
231
232         ! Define the U grid FILE ( nid_U )
233
234         CALL dia_nam( clhstnam, nn_write, 'grid_U' )
235         IF(lwp) WRITE(numout,*) " Name of NETCDF file ", clhstnam    ! filename
236         CALL histbeg( clhstnam, jpi, glamu, jpj, gphiu,           &  ! Horizontal grid: glamu and gphiu
237            &          iimi, iima-iimi+1, ijmi, ijma-ijmi+1,       &
238            &          nit000-1, zjulian, rn_Dt, nh_U, nid_U, domain_id=nidom, snc4chunks=snc4set )
239         CALL histvert( nid_U, "depthu", "Vertical U levels",      &  ! Vertical grid: gdept
240            &           "m", ipk, gdept_1d, nz_U, "down" )
241         !                                                            ! Index of ocean points
242         CALL wheneq( jpi*jpj    , umask, 1, 1., ndex_hU, ndim_hU )      ! surface
243
244         ! Define the V grid FILE ( nid_V )
245
246         CALL dia_nam( clhstnam, nn_write, 'grid_V' )                   ! filename
247         IF(lwp) WRITE(numout,*) " Name of NETCDF file ", clhstnam
248         CALL histbeg( clhstnam, jpi, glamv, jpj, gphiv,           &  ! Horizontal grid: glamv and gphiv
249            &          iimi, iima-iimi+1, ijmi, ijma-ijmi+1,       &
250            &          nit000-1, zjulian, rn_Dt, nh_V, nid_V, domain_id=nidom, snc4chunks=snc4set )
251         CALL histvert( nid_V, "depthv", "Vertical V levels",      &  ! Vertical grid : gdept
252            &          "m", ipk, gdept_1d, nz_V, "down" )
253         !                                                            ! Index of ocean points
254         CALL wheneq( jpi*jpj    , vmask, 1, 1., ndex_hV, ndim_hV )      ! surface
255
256         ! No W grid FILE
257         IF( ln_abl ) THEN 
258         ! Define the ABL grid FILE ( nid_A )
259            CALL dia_nam( clhstnam, nwrite, 'grid_ABL' )
260            IF(lwp) WRITE(numout,*) " Name of NETCDF file ", clhstnam    ! filename
261            CALL histbeg( clhstnam, jpi, glamt, jpj, gphit,           &  ! Horizontal grid: glamt and gphit
262               &          iimi, iima-iimi+1, ijmi, ijma-ijmi+1,       &
263               &          nit000-1, zjulian, rn_Dt, nh_A, nid_A, domain_id=nidom, snc4chunks=snc4set )
264            CALL histvert( nid_A, "ght_abl", "Vertical T levels",      &  ! Vertical grid: gdept
265               &           "m", ipka, ght_abl(2:jpka), nz_A, "up" )
266            !                                                            ! Index of ocean points
267         ALLOCATE( zw3d_abl(jpi,jpj,ipka) ) 
268         zw3d_abl(:,:,:) = 1._wp 
269         CALL wheneq( jpi*jpj*ipka, zw3d_abl, 1, 1., ndex_A , ndim_A  )      ! volume
270            CALL wheneq( jpi*jpj     , zw3d_abl, 1, 1., ndex_hA, ndim_hA )      ! surface
271         DEALLOCATE(zw3d_abl)
272         ENDIF
273
274         ! Declare all the output fields as NETCDF variables
275
276         !                                                                                      !!! nid_T : 3D
277         CALL histdef( nid_T, "sst_m", "Sea Surface temperature"            , "C"      ,   &  ! sst
278            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
279         CALL histdef( nid_T, "sss_m", "Sea Surface Salinity"               , "PSU"    ,   &  ! sss
280            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
281         CALL histdef( nid_T, "sowaflup", "Net Upward Water Flux"              , "Kg/m2/s",   &  ! (emp-rnf)
282            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
283         CALL histdef( nid_T, "sosfldow", "downward salt flux"                 , "PSU/m2/s",  &  ! (sfx)
284             &         jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
285         CALL histdef( nid_T, "sohefldo", "Net Downward Heat Flux"             , "W/m2"   ,   &  ! qns + qsr
286            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
287         CALL histdef( nid_T, "soshfldo", "Shortwave Radiation"                , "W/m2"   ,   &  ! qsr
288            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
289         CALL histdef( nid_T, "soicecov", "Ice fraction"                       , "[0,1]"  ,   &  ! fr_i
290            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
291         CALL histdef( nid_T, "sowindsp", "wind speed at 10m"                  , "m/s"    ,   &  ! wndm
292            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
293!
294         IF( ln_abl ) THEN
295         !                                                                                      !!! nid_A : 3D
296         CALL histdef( nid_A, "t_abl", "Potential Temperature"     , "K"        ,       &  ! t_abl
297               &          jpi, jpj, nh_A, ipka, 1, ipka, nz_A, 32, clop, zsto, zout )
298            CALL histdef( nid_A, "q_abl", "Humidity"                  , "kg/kg"    ,       &  ! q_abl
299               &          jpi, jpj, nh_A, ipka, 1, ipka, nz_A, 32, clop, zsto, zout ) 
300            CALL histdef( nid_A, "u_abl", "Atmospheric U-wind   "     , "m/s"        ,     &  ! u_abl
301               &          jpi, jpj, nh_A, ipka, 1, ipka, nz_A, 32, clop, zsto, zout )
302            CALL histdef( nid_A, "v_abl", "Atmospheric V-wind   "     , "m/s"    ,         &  ! v_abl
303               &          jpi, jpj, nh_A, ipka, 1, ipka, nz_A, 32, clop, zsto, zout ) 
304            CALL histdef( nid_A, "tke_abl", "Atmospheric TKE   "     , "m2/s2"    ,        &  ! tke_abl
305               &          jpi, jpj, nh_A, ipka, 1, ipka, nz_A, 32, clop, zsto, zout ) 
306            CALL histdef( nid_A, "avm_abl", "Atmospheric turbulent viscosity", "m2/s"   ,  &  ! avm_abl
307               &          jpi, jpj, nh_A, ipka, 1, ipka, nz_A, 32, clop, zsto, zout ) 
308            CALL histdef( nid_A, "avt_abl", "Atmospheric turbulent diffusivity", "m2/s2",  &  ! avt_abl
309               &          jpi, jpj, nh_A, ipka, 1, ipka, nz_A, 32, clop, zsto, zout ) 
310            CALL histdef( nid_A, "pblh", "Atmospheric boundary layer height "  , "m",      &  ! pblh
311               &          jpi, jpj, nh_A,  1  , 1, 1   , -99 , 32, clop, zsto, zout )                 
312#if defined key_si3
313            CALL histdef( nid_A, "oce_frac", "Fraction of open ocean"  , " ",      &  ! ato_i
314               &          jpi, jpj, nh_A,  1  , 1, 1   , -99 , 32, clop, zsto, zout )
315#endif
316          CALL histend( nid_A, snc4chunks=snc4set )
317       !
318       ENDIF
319!
320
321         CALL histend( nid_T, snc4chunks=snc4set )
322
323         !                                                                                      !!! nid_U : 3D
324         CALL histdef( nid_U, "ssu_m", "Velocity component in x-direction", "m/s"   ,         &  ! ssu
325            &          jpi, jpj, nh_U, 1  , 1, 1  , - 99, 32, clop, zsto, zout )
326         CALL histdef( nid_U, "sozotaux", "Wind Stress along i-axis"           , "N/m2"   ,   &  ! utau
327            &          jpi, jpj, nh_U, 1  , 1, 1  , - 99, 32, clop, zsto, zout )
328
329         CALL histend( nid_U, snc4chunks=snc4set )
330
331         !                                                                                      !!! nid_V : 3D
332         CALL histdef( nid_V, "ssv_m", "Velocity component in y-direction", "m/s",            &  ! ssv_m
333            &          jpi, jpj, nh_V, 1  , 1, 1  , - 99, 32, clop, zsto, zout )
334         CALL histdef( nid_V, "sometauy", "Wind Stress along j-axis"           , "N/m2"   ,   &  ! vtau
335            &          jpi, jpj, nh_V, 1  , 1, 1  , - 99, 32, clop, zsto, zout )
336
337         CALL histend( nid_V, snc4chunks=snc4set )
338
339         IF(lwp) WRITE(numout,*)
340         IF(lwp) WRITE(numout,*) 'End of NetCDF Initialization'
341         IF(ll_print) CALL FLUSH(numout )
342
343      ENDIF
344
345      ! 2. Start writing data
346      ! ---------------------
347
348      ! ndex(1) est utilise ssi l'avant dernier argument est diffferent de
349      ! la taille du tableau en sortie. Dans ce cas , l'avant dernier argument
350      ! donne le nombre d'elements, et ndex la liste des indices a sortir
351
352      IF( lwp .AND. MOD( itmod, nn_write ) == 0 ) THEN
353         WRITE(numout,*) 'dia_wri : write model outputs in NetCDF files at ', kt, 'time-step'
354         WRITE(numout,*) '~~~~~~ '
355      ENDIF
356
357      ! Write fields on T grid
358      CALL histwrite( nid_T, "sst_m", it, sst_m, ndim_hT, ndex_hT )   ! sea surface temperature
359      CALL histwrite( nid_T, "sss_m", it, sss_m, ndim_hT, ndex_hT )   ! sea surface salinity
360      CALL histwrite( nid_T, "sowaflup", it, (emp - rnf )  , ndim_hT, ndex_hT )   ! upward water flux
361      CALL histwrite( nid_T, "sosfldow", it, sfx           , ndim_hT, ndex_hT )   ! downward salt flux
362                                                                                  ! (includes virtual salt flux beneath ice
363                                                                                  ! in linear free surface case)
364
365      CALL histwrite( nid_T, "sohefldo", it, qns + qsr     , ndim_hT, ndex_hT )   ! total heat flux
366      CALL histwrite( nid_T, "soshfldo", it, qsr           , ndim_hT, ndex_hT )   ! solar heat flux
367      CALL histwrite( nid_T, "soicecov", it, fr_i          , ndim_hT, ndex_hT )   ! ice fraction   
368      CALL histwrite( nid_T, "sowindsp", it, wndm          , ndim_hT, ndex_hT )   ! wind speed   
369!
370      IF( ln_abl ) THEN
371        ALLOCATE( zw3d_abl(jpi,jpj,jpka) )
372        IF( ln_mskland )   THEN
373          DO jk=1,jpka
374             zw3d_abl(:,:,jk) = tmask(:,:,1)
375            END DO
376       ELSE
377            zw3d_abl(:,:,:) = 1._wp     
378         ENDIF       
379       CALL histwrite( nid_A,  "pblh"   , it, pblh(:,:)                  *zw3d_abl(:,:,1     ), ndim_hA, ndex_hA )   ! pblh
380        CALL histwrite( nid_A,  "u_abl"  , it, u_abl   (:,:,2:jpka,nt_n  )*zw3d_abl(:,:,2:jpka), ndim_A , ndex_A  )   ! u_abl
381        CALL histwrite( nid_A,  "v_abl"  , it, v_abl   (:,:,2:jpka,nt_n  )*zw3d_abl(:,:,2:jpka), ndim_A , ndex_A  )   ! v_abl
382        CALL histwrite( nid_A,  "t_abl"  , it, tq_abl  (:,:,2:jpka,nt_n,1)*zw3d_abl(:,:,2:jpka), ndim_A , ndex_A  )   ! t_abl
383        CALL histwrite( nid_A,  "q_abl"  , it, tq_abl  (:,:,2:jpka,nt_n,2)*zw3d_abl(:,:,2:jpka), ndim_A , ndex_A  )   ! q_abl     
384        CALL histwrite( nid_A,  "tke_abl", it, tke_abl (:,:,2:jpka,nt_n  )*zw3d_abl(:,:,2:jpka), ndim_A , ndex_A  )   ! tke_abl
385        CALL histwrite( nid_A,  "avm_abl", it, avm_abl (:,:,2:jpka       )*zw3d_abl(:,:,2:jpka), ndim_A , ndex_A  )   ! avm_abl
386        CALL histwrite( nid_A,  "avt_abl", it, avt_abl (:,:,2:jpka       )*zw3d_abl(:,:,2:jpka), ndim_A , ndex_A  )   ! avt_abl 
387#if defined key_si3
388         CALL histwrite( nid_A,  "oce_frac"   , it, ato_i(:,:)                                  , ndim_hA, ndex_hA )   ! ato_i
389#endif
390       DEALLOCATE(zw3d_abl)
391     ENDIF
392!
393
394         ! Write fields on U grid
395      CALL histwrite( nid_U, "ssu_m"   , it, ssu_m         , ndim_hU, ndex_hU )   ! i-current speed
396      CALL histwrite( nid_U, "sozotaux", it, utau          , ndim_hU, ndex_hU )   ! i-wind stress
397
398         ! Write fields on V grid
399      CALL histwrite( nid_V, "ssv_m"   , it, ssv_m         , ndim_hV, ndex_hV )   ! j-current speed
400      CALL histwrite( nid_V, "sometauy", it, vtau          , ndim_hV, ndex_hV )   ! j-wind stress
401
402      ! 3. Close all files
403      ! ---------------------------------------
404      IF( kt == nitend ) THEN
405         CALL histclo( nid_T )
406         CALL histclo( nid_U )
407         CALL histclo( nid_V )
408         IF(ln_abl) CALL histclo( nid_A )
409      ENDIF
410      !
411      IF( ln_timing )   CALL timing_stop('dia_wri')
412      !
413   END SUBROUTINE dia_wri
414#endif
415
416   SUBROUTINE dia_wri_state( cdfile_name, Kmm )
417      !!---------------------------------------------------------------------
418      !!                 ***  ROUTINE dia_wri_state  ***
419      !!       
420      !! ** Purpose :   create a NetCDF file named cdfile_name which contains
421      !!      the instantaneous ocean state and forcing fields.
422      !!        Used to find errors in the initial state or save the last
423      !!      ocean state in case of abnormal end of a simulation
424      !!
425      !! ** Method  :   NetCDF files using ioipsl
426      !!      File 'output.init.nc'  is created if ninist = 1 (namelist)
427      !!      File 'output.abort.nc' is created in case of abnormal job end
428      !!----------------------------------------------------------------------
429      CHARACTER (len=* ), INTENT( in ) ::   cdfile_name      ! name of the file created
430      INTEGER           , INTENT( in ) ::   Kmm              ! ocean time levelindex
431      !!
432      INTEGER :: inum
433      !!----------------------------------------------------------------------
434      !
435      IF(lwp) WRITE(numout,*)
436      IF(lwp) WRITE(numout,*) 'dia_wri_state : single instantaneous ocean state'
437      IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~   and forcing fields file created '
438      IF(lwp) WRITE(numout,*) '                and named :', cdfile_name, '...nc'
439
440#if defined key_si3
441     CALL iom_open( TRIM(cdfile_name), inum, ldwrt = .TRUE., kdlev = jpl )
442#else
443     CALL iom_open( TRIM(cdfile_name), inum, ldwrt = .TRUE. )
444#endif
445
446      CALL iom_rstput( 0, 0, inum, 'votemper', ts (:,:,:,jp_tem,Kmm) )    ! now temperature
447      CALL iom_rstput( 0, 0, inum, 'vosaline', ts (:,:,:,jp_sal,Kmm) )    ! now salinity
448      CALL iom_rstput( 0, 0, inum, 'sossheig', ssh(:,:,         Kmm) )    ! sea surface height
449      CALL iom_rstput( 0, 0, inum, 'vozocrtx', uu (:,:,:,       Kmm) )    ! now i-velocity
450      CALL iom_rstput( 0, 0, inum, 'vomecrty', vv (:,:,:,       Kmm) )    ! now j-velocity
451      CALL iom_rstput( 0, 0, inum, 'vovecrtz', ww                    )    ! now k-velocity
452      CALL iom_rstput( 0, 0, inum, 'sowaflup', emp - rnf             )    ! freshwater budget
453      CALL iom_rstput( 0, 0, inum, 'sohefldo', qsr + qns             )    ! total heat flux
454      CALL iom_rstput( 0, 0, inum, 'soshfldo', qsr                   )    ! solar heat flux
455      CALL iom_rstput( 0, 0, inum, 'soicecov', fr_i                  )    ! ice fraction
456      CALL iom_rstput( 0, 0, inum, 'sozotaux', utau                  )    ! i-wind stress
457      CALL iom_rstput( 0, 0, inum, 'sometauy', vtau                  )    ! j-wind stress
458 
459#if defined key_si3
460      IF( nn_ice == 2 ) THEN   ! condition needed in case agrif + ice-model but no-ice in child grid
461         CALL ice_wri_state( inum )
462      ENDIF
463#endif
464      !
465      CALL iom_close( inum )
466      !
467   END SUBROUTINE dia_wri_state
468
469   !!======================================================================
470END MODULE diawri
Note: See TracBrowser for help on using the repository browser.