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_c1d.F90 in branches/CMIP5_IPSL/NEMO/C1D_SRC – NEMO

source: branches/CMIP5_IPSL/NEMO/C1D_SRC/diawri_c1d.F90 @ 5600

Last change on this file since 5600 was 1846, checked in by mafoipsl, 14 years ago

Increase length of characters string to allow long name of experiments.

  • Property svn:eol-style set to native
  • Property svn:keywords set to Id
File size: 21.1 KB
RevLine 
[900]1MODULE diawri_c1d
[253]2   !!======================================================================
[900]3   !!                     ***  MODULE  diawri_c1d  ***
[253]4   !! Ocean diagnostics :  write ocean output files
5   !!=====================================================================
[900]6   !! History :   2.0  !  2004-10  (C. Ethe)  1D Configuration
7   !!             3.0  !  2008-04  (G. Madec)  adaptation to SBC
8   !!----------------------------------------------------------------------
[899]9#if defined key_c1d
[253]10   !!----------------------------------------------------------------------
[900]11   !!   'key_c1d'                                          1D Configuration
[253]12   !!---------------------------------------------------------------------- 
[900]13   !!   dia_wri_c1d       : create the standart NetCDF output files
[253]14   !!----------------------------------------------------------------------
15   USE oce             ! ocean dynamics and tracers
16   USE dom_oce         ! ocean space and time domain
17   USE zdf_oce         ! ocean vertical physics
[888]18   USE sbc_oce         ! surface boundary condition: ocean
19   USE sbc_ice         ! surface boundary condition: ice
[900]20   USE sbcmod          ! surface Boundary Codition
21   USE sbcssr          ! surface boundary condition: restauring to SSS and or SST
[1533]22   USE zdftke_old      ! TKE vertical mixing (old scheme)
[253]23   USE zdftke          ! TKE vertical mixing
[255]24   USE zdfkpp          ! KPP vertical mixing
[253]25   USE sol_oce         ! solver variables
26   USE phycst          ! physical constants
27   USE zdfmxl          ! mixed layer
28   USE dianam          ! build name of file (routine)
29   USE diawri
30   USE zdfddm          ! vertical  physics: double diffusion
31   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
32   USE in_out_manager  ! I/O manager
[900]33   USE ioipsl
[253]34
35   IMPLICIT NONE
36   PRIVATE
37
[900]38   PUBLIC dia_wri_c1d                 ! routines called by step.F90
[253]39
[900]40   INTEGER ::   nid_T, nz_T, nh_T, ndim_T, ndim_hT, ndex(1)   ! grid_T file
41   INTEGER, DIMENSION(jpi*jpj)     ::   ndex_hT
42   INTEGER, DIMENSION(jpi*jpj*jpk) ::   ndex_T
43
[253]44   !! * Substitutions
45#  include "zdfddm_substitute.h90"
46   !!----------------------------------------------------------------------
[900]47   !! NEMO/C1D 3.0 , LOCEAN-IPSL  (2008)
[888]48   !! $Id$
[900]49   !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt)
[253]50   !!----------------------------------------------------------------------
51
52CONTAINS
53
[900]54   SUBROUTINE dia_wri_c1d( kt, kindic )
[253]55      !!---------------------------------------------------------------------
[900]56      !!                  ***  ROUTINE dia_wri_c1d  ***
[253]57      !!                   
58      !! ** Purpose :   Standard output of opa: dynamics and tracer fields
59      !!      NETCDF format is used by default
60      !!
61      !! ** Method  :   At the beginning of the first time step (nit000),
62      !!      define all the NETCDF files and fields
63      !!      At each time step call histdef to compute the mean if ncessary
64      !!      Each nwrite time step, output the instantaneous or mean fields
65      !!      IF kindic <0, output of fields before the model interruption.
66      !!      IF kindic =0, time step loop
67      !!      IF kindic >0, output of fields before the time step loop
68      !!----------------------------------------------------------------------
69      INTEGER, INTENT( in ) ::   kt      ! ocean time-step index
70      INTEGER, INTENT( in ) ::   kindic  !
[900]71      !!
72      LOGICAL ::   ll_print = .FALSE.                ! =T print and flush numout
[1846]73      CHARACTER (len=72) ::   clhstnam, clop, clmx   ! temporary names
[900]74      INTEGER ::   inum = 11                         ! temporary logical unit
75      INTEGER ::   ji, jj, ik                        ! dummy loop indices
76      INTEGER ::   iimi, iima, ipk, it, ijmi, ijma   ! temporary integers
[1334]77      INTEGER ::   itmod                             !
[900]78      REAL(wp) ::   zsto, zout, zmax, zjulian, zdt   ! temporary scalars
79      REAL(wp), DIMENSION(jpi,jpj) ::   zw2d         ! temporary workspace
[253]80      !!----------------------------------------------------------------------
81     
82      ! 0. Initialisation
83      ! -----------------
84     
85      ! local variable for debugging
86      ll_print = .FALSE.
87      ll_print = ll_print .AND. lwp
88
89      ! Define frequency of output and means
90      zdt = rdt
91      IF( nacc == 1 ) zdt = rdtmin
[1312]92      IF( ln_mskland )   THEN   ;   clop = "only(x)"   ! put 1.e+20 on land (very expensive!!)
93      ELSE                      ;   clop = "x"         ! no use of the mask value (require less cpu time)
94      ENDIF
[253]95#if defined key_diainstant
96      zsto = nwrite * zdt
[1312]97      clop = "inst("//TRIM(clop)//")"
[253]98#else
99      zsto=zdt
[1312]100      clop = "ave("//TRIM(clop)//")"
[253]101#endif
102      zout = nwrite * zdt
103      zmax = ( nitend - nit000 + 1 ) * zdt
104
105      ! Define indices of the horizontal output zoom and vertical limit storage
106      iimi = 1      ;      iima = jpi
107      ijmi = 1      ;      ijma = jpj
108      ipk = jpk
109
110      ! define time axis
[1334]111      it = kt
112      itmod = kt - nit000 + 1
[253]113
114
115      ! 1. Define NETCDF files and fields at beginning of first time step
116      ! -----------------------------------------------------------------
117
[900]118      IF(ll_print) WRITE(numout,*) 'dia_wri_c1d kt = ', kt, ' kindic ', kindic
[253]119
120      IF( kt == nit000 ) THEN
121
122         ! Define the NETCDF files (one per grid)
123         
124         ! Compute julian date from starting date of the run
[1310]125         CALL ymds2ju( nyear, nmonth, nday, rdt, zjulian )
126         zjulian = zjulian - adatrj   !   set calendar origin to the beginning of the experiment
[253]127         IF(lwp)WRITE(numout,*)
128         IF(lwp)WRITE(numout,*) 'Date 0 used :', nit000, ' YEAR ', nyear,   &
129            &                    ' MONTH ', nmonth, ' DAY ', nday, 'Julian day : ', zjulian
130         IF(lwp)WRITE(numout,*) ' indexes of zoom = ', iimi, iima, ijmi, ijma,   &
131                                 ' limit storage in depth = ', ipk
132
133         ! WRITE root name in date.file for use by postpro
134         CALL dia_nam( clhstnam, nwrite,' ' )
[1581]135         CALL ctl_opn( inum, 'date.file', 'REPLACE', 'FORMATTED', 'SEQUENTIAL', 1, numout, lwp )
[253]136         WRITE(inum,*) clhstnam
137         CLOSE(inum)
138         
139         ! Define the T grid FILE ( nid_T )
140         
141         CALL dia_nam( clhstnam, nwrite, 'grid_T' )
142         IF(lwp) WRITE(numout,*) " Name of NETCDF file ", clhstnam    ! filename
143         CALL histbeg( clhstnam, jpi, glamt, jpj, gphit,           &  ! Horizontal grid: glamt and gphit
144            &          iimi, iima-iimi+1, ijmi, ijma-ijmi+1,       &
[1334]145            &          nit000-1, zjulian, zdt, nh_T, nid_T, domain_id=nidom )
[253]146         CALL histvert( nid_T, "deptht", "Vertical T levels",      &  ! Vertical grid: gdept
[1334]147            &           "m", ipk, gdept_0, nz_T, "down" )
[253]148         !                                                            ! Index of ocean points
149         CALL wheneq( jpi*jpj*ipk, tmask, 1, 1., ndex_T , ndim_T  )      ! volume
150         CALL wheneq( jpi*jpj    , tmask, 1, 1., ndex_hT, ndim_hT )      ! surface
151
152
153         ! Declare all the output fields as NETCDF variables
154
155         !                                                                                      !!! nid_T : 3D
156         CALL histdef( nid_T, "votemper", "Temperature"                        , "C"      ,   &  ! tn
157            &          jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout )
158         CALL histdef( nid_T, "vosaline", "Salinity"                           , "PSU"    ,   &  ! sn
159            &          jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout )
160         !                                                                                      !!! nid_T : 2D
161         CALL histdef( nid_T, "sosstsst", "Sea Surface temperature"            , "C"      ,   &  ! sst
162            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
163         CALL histdef( nid_T, "sosaline", "Sea Surface Salinity"               , "PSU"    ,   &  ! sss
164            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
165         CALL histdef( nid_T, "sowaflup", "Net Upward Water Flux"              , "Kg/m2/s",   &  ! emp
166            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
[900]167!!       CALL histdef( nid_T, "sorunoff", "Runoffs"                            , "Kg/m2/s",   &  ! runoffs
168!!          &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
[253]169         CALL histdef( nid_T, "sowaflcd", "concentration/dilution water flux"  , "kg/m2/s",   &  ! emps
170            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
171         CALL histdef( nid_T, "sosalflx", "Surface Salt Flux"                  , "Kg/m2/s",   &  ! emps * sn
172            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
[888]173         CALL histdef( nid_T, "sohefldo", "Net Downward Heat Flux"             , "W/m2"   ,   &  ! qsr + qns
[253]174            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
175         CALL histdef( nid_T, "soshfldo", "Shortwave Radiation"                , "W/m2"   ,   &  ! qsr
176            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
177         CALL histdef( nid_T, "somxl010", "Mixed Layer Depth 0.01"             , "m"      ,   &  ! hmlp
178            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
[255]179#if defined key_zdfkpp
180         CALL histdef( nid_T, "sokppekd", "Ekman depth                     "   , "m"      ,   &  ! sokppekd
181            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
182         CALL histdef( nid_T, "sokppbld", "Boundary Layer Depth            "   , "m"      ,   &  ! sokppbld
183            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
184#endif
[253]185         CALL histdef( nid_T, "somxlavt", "AVT : bottom of the mixed layer    ", "m"      ,   &  ! avt_mxl
186            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
187         CALL histdef( nid_T, "somixhgt", "Turbocline Depth"                   , "m"      ,   &  ! hmld
188            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
[1037]189         CALL histdef( nid_T, "soicecov", "Ice Fraction"                          , "[0,1]"  ,   &  ! fr_i
[253]190            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
[900]191         IF( ln_ssr ) THEN
192            CALL histdef( nid_T, "sohefldp", "Surface Heat Flux: Damping"         , "W/m2"   ,   &  ! qrp
193               &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
194            CALL histdef( nid_T, "sowafldp", "Surface Water Flux: Damping"        , "Kg/m2/s",   &  ! erp
195               &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
196            CALL histdef( nid_T, "sosafldp", "Surface salt flux: damping"         , "Kg/m2/s",   &  ! erp * sn
197               &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
198         ENDIF
[253]199         clmx ="l_max(only(x))"    ! max index on a period
200         CALL histdef( nid_T, "sobowlin", "Bowl Index"                         , "W-point",   &  ! bowl INDEX
201            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clmx, zsto, zout )
202         !                                                                                      !!! nid_U : 3D
203         CALL histdef( nid_T, "vozocrtx", "Zonal Current"                      , "m/s"    ,   &  ! un
204            &          jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout )
205#if defined key_diaeiv
206         CALL histdef( nid_T, "vozoeivu", "Zonal EIV Current"                  , "m/s"    ,   &  ! u_eiv
207            &          jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout )
208#endif
209         !                                                                                      !!! nid_U : 2D
[888]210         CALL histdef( nid_T, "sozotaux", "Wind Stress along i-axis"           , "N/m2"   ,   &  ! utau
[253]211            &          jpi, jpj, nh_T, 1  , 1, 1  , - 99, 32, clop, zsto, zout )
212
213         !                                                                                      !!! nid_V : 3D
214         CALL histdef( nid_T, "vomecrty", "Meridional Current"                 , "m/s"    ,   &  ! vn
215            &          jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout )
216#if defined key_diaeiv
217         CALL histdef( nid_T, "vomeeivv", "Meridional EIV Current"             , "m/s"    ,   &  ! v_eiv
218            &          jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout )
219#endif
220         !                                                                                      !!! nid_V : 2D
[888]221         CALL histdef( nid_T, "sometauy", "Wind Stress along j-axis"           , "N/m2"   ,   &  ! vtau
[253]222            &          jpi, jpj, nh_T, 1  , 1, 1  , - 99, 32, clop, zsto, zout )
[1531]223#if defined key_zdftke_old
[253]224         CALL histdef( nid_T, "votlsdis", " Dissipation Turbulent Lenght Scale", "m"      ,   &  ! e_dis
225            &          jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout )
226         CALL histdef( nid_T, "votlsmix", " Mixing Turbulent Lenght Scale"     , "m"      ,   &  ! e_mix
227            &          jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout )
228         CALL histdef( nid_T, "votlspdl", " Prandl Number",                      "-"       ,   &  ! e_pdl
229            &          jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout )
230         CALL histdef( nid_T, "votlsric", " Local Richardson Number",            "-"       ,   &  ! e_ric
231            &          jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout )
232         CALL histdef( nid_T, "votkeend", "TKE: Turbulent kinetic energy"       , "m2/s"   ,   &  ! TKE
233            &          jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout )
[1221]234         CALL histdef( nid_T, "sotkehlc", "TKE: Langmuir Circulation depth"     , "m"      ,   &  ! hlc
235            &          jpi, jpj, nh_T, 1  , 1, 1  , -99 , 32, clop, zsto, zout )
[253]236#endif
[255]237#if defined key_zdfkpp
238         CALL histdef( nid_T, "vokpprig", " Gradient Richardson Number"        ,  "-"      ,   &  ! rig
239            &          jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout )
240         CALL histdef( nid_T, "vokpprib", " Bulk Richardson Number    "        ,  "-"      ,   &   ! rib
241            &          jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout )
242         CALL histdef( nid_T, "vokppbsf", " Buoyancy forcing          "        , "N/m2"    ,   &  ! sokppbsf
243            &          jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout )
244         CALL histdef( nid_T, "vokppmol", "Moning Obukhov length scale     "   , "m"       ,   &  ! sokppmol
245            &          jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout )
246#endif
247         !
[253]248         CALL histdef( nid_T, "voeosbn2", "Brunt-Vaisala Frequency"             , "m2/s2"  ,   &  ! rn2
249            &          jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout )
250         CALL histdef( nid_T, "votkeavt", "Vertical Eddy Diffusivity"          , "m2/s"   ,   &  ! avt
251            &          jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout )
252         CALL histdef( nid_T, "votkeavm", "Vertical Eddy Viscosity",             "m2/s"   ,   &  ! avmu
253            &          jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout )
254         !
255         IF( lk_zdfddm ) THEN
256            CALL histdef( nid_T,"voddmavs","Salt Vertical Eddy Diffusivity"    , "m2/s"   ,   &  ! avs
257               &          jpi, jpj, nh_T, ipk, 1, ipk, nz_T, 32, clop, zsto, zout )
258         ENDIF
259
260         CALL histend( nid_T )
261
262         IF(lwp) WRITE(numout,*)
263         IF(lwp) WRITE(numout,*) 'End of NetCDF Initialization'
264         IF(ll_print) CALL FLUSH(numout )
265
266      ENDIF
267
268      ! 2. Start writing data
269      ! ---------------------
270
271      ! ndex(1) est utilise ssi l'avant dernier argument est diffferent de
272      ! la taille du tableau en sortie. Dans ce cas , l'avant dernier argument
273      ! donne le nombre d'elements, et ndex la liste des indices a sortir
274
[1334]275      IF( lwp .AND. MOD( itmod, nwrite ) == 0 ) THEN
[253]276         WRITE(numout,*) 'dia_wri : write model outputs in NetCDF files at ', kt, 'time-step'
277         WRITE(numout,*) '~~~~~~ '
278      ENDIF
279
280      ! Write fields on T grid
281      CALL histwrite( nid_T, "votemper", it, tn            , ndim_T , ndex_T  )   ! temperature
282      CALL histwrite( nid_T, "vosaline", it, sn            , ndim_T , ndex_T  )   ! salinity
283      CALL histwrite( nid_T, "sosstsst", it, tn(:,:,1)     , ndim_hT, ndex_hT )   ! sea surface temperature
284      CALL histwrite( nid_T, "sosaline", it, sn(:,:,1)     , ndim_hT, ndex_hT )   ! sea surface salinity
285      CALL histwrite( nid_T, "sowaflup", it, emp           , ndim_hT, ndex_hT )   ! upward water flux
[900]286!!    CALL histwrite( nid_T, "sorunoff", it, runoff        , ndim_hT, ndex_hT )   ! runoff
[253]287      CALL histwrite( nid_T, "sowaflcd", it, emps          , ndim_hT, ndex_hT )   ! c/d water flux
288      zw2d(:,:) = emps(:,:) * sn(:,:,1) * tmask(:,:,1)
289      CALL histwrite( nid_T, "sosalflx", it, zw2d          , ndim_hT, ndex_hT )   ! c/d salt flux
[888]290      CALL histwrite( nid_T, "sohefldo", it, qsr + qns     , ndim_hT, ndex_hT )   ! total heat flux
[253]291      CALL histwrite( nid_T, "soshfldo", it, qsr           , ndim_hT, ndex_hT )   ! solar heat flux
292      CALL histwrite( nid_T, "somxl010", it, hmlp          , ndim_hT, ndex_hT )   ! mixed layer depth
[255]293#if defined key_zdfkpp
294      CALL histwrite( nid_T, "sokppekd", it, ekdp          , ndim_hT, ndex_hT )   ! Ekman depht
295      CALL histwrite( nid_T, "sokppbld", it, hkpp          , ndim_hT, ndex_hT )   ! boundary layer depth
296#endif 
[253]297      ! store the vertical eddy diffusivity coef. at the bottom of the mixed layer
298      DO jj = 1, jpj
299         DO ji = 1, jpi
300            ik = nmln(ji,jj)
301            zw2d(ji,jj) = avt(ji,jj,ik) * tmask(ji,jj,1)
302         END DO
303      END DO
304      CALL histwrite( nid_T, "somxlavt", it, zw2d          , ndim_hT, ndex_hT )   ! Kz at bottom of mixed layer
305      CALL histwrite( nid_T, "somixhgt", it, hmld          , ndim_hT, ndex_hT )   ! turbocline depth
[1037]306      CALL histwrite( nid_T, "soicecov", it, fr_i          , ndim_hT, ndex_hT )   ! ice fraction
[900]307      IF( ln_ssr ) THEN
308         CALL histwrite( nid_T, "sohefldp", it, qrp           , ndim_hT, ndex_hT )   ! heat flux damping
309         CALL histwrite( nid_T, "sowafldp", it, erp           , ndim_hT, ndex_hT )   ! freshwater flux damping
[253]310         zw2d(:,:) = erp(:,:) * sn(:,:,1) * tmask(:,:,1)
[900]311         CALL histwrite( nid_T, "sosafldp", it, zw2d          , ndim_hT, ndex_hT )   ! salt flux damping
312      ENDIF
313      zw2d(:,:) = FLOAT( nmln(:,:) ) * tmask(:,:,1)
[253]314      CALL histwrite( nid_T, "sobowlin", it, zw2d          , ndim_hT, ndex_hT )   ! ???
315      CALL histwrite( nid_T, "vozocrtx", it, un            , ndim_T , ndex_T )    ! i-current
[888]316      CALL histwrite( nid_T, "sozotaux", it, utau          , ndim_hT, ndex_hT )   ! i-wind stress
[253]317      CALL histwrite( nid_T, "vomecrty", it, vn            , ndim_T , ndex_T  )   ! j-current
[888]318      CALL histwrite( nid_T, "sometauy", it, vtau          , ndim_hT, ndex_hT )   ! j-wind stress
[1531]319#if defined key_zdftke_old
[253]320      CALL histwrite( nid_T, "votlsdis", it, e_dis         , ndim_T , ndex_T )    ! Diss. Turb. lenght scale
321      CALL histwrite( nid_T, "votlsmix", it, e_mix         , ndim_T , ndex_T )    ! Mixing Turb. lenght scale
322      CALL histwrite( nid_T, "votlspdl", it, e_pdl         , ndim_T , ndex_T )    ! Prandl number
323      CALL histwrite( nid_T, "votlsric", it, e_ric         , ndim_T , ndex_T )    ! local Richardson number
324      CALL histwrite( nid_T, "votkeend", it, en            , ndim_T , ndex_T )    ! TKE
[1221]325      CALL histwrite( nid_T, "sotkehlc", it, hlc           , ndim_hT, ndex_hT )   ! TKE Langmuir Circulation depth
[253]326#endif
[255]327#if defined key_zdfkpp
328      CALL histwrite( nid_T, "vokpprig", it, rig           , ndim_T , ndex_T )    ! gradient Richardson number
329      CALL histwrite( nid_T, "vokpprib", it, rib           , ndim_T , ndex_T )    ! bulk Richardson number
330      CALL histwrite( nid_T, "vokppbsf", it, buof          , ndim_T , ndex_T )    ! buoyancy forcing
331      CALL histwrite( nid_T, "vokppmol", it, mols          , ndim_T , ndex_T )    ! Moning-Obukov length scale
332#endif
[253]333      CALL histwrite( nid_T, "voeosbn2", it, rn2           , ndim_T , ndex_T )    ! Brunt-Vaisala Frequency
334      CALL histwrite( nid_T, "votkeavt", it, avt           , ndim_T , ndex_T )    ! T vert. eddy diff. coef.
335      CALL histwrite( nid_T, "votkeavm", it, avmu          , ndim_T , ndex_T )    ! T vert. eddy visc. coef.
336      IF( lk_zdfddm ) THEN
337         CALL histwrite( nid_T, "voddmavs", it, fsavs(:,:,:), ndim_T, ndex_T )    ! S vert. eddy diff. coef.
338      ENDIF
339
340      !  Create an output files (output.abort.nc) if S < 0 or u > 20 m/s
[1334]341      IF( kindic < 0 )   CALL dia_wri_state( 'output.abort', kt )
[253]342
[1318]343      ! 3. Close all files
344      ! ---------------------------------------
[900]345      IF( kt == nitend .OR. kindic < 0 )   CALL histclo( nid_T )
346      !
347   END SUBROUTINE dia_wri_c1d
[253]348
349#else
350   !!----------------------------------------------------------------------
[900]351   !!   Default key                                     NO 1D Configuration
[253]352   !!----------------------------------------------------------------------
353CONTAINS
[900]354   SUBROUTINE dia_wri_c1d ( kt, kindic )      ! dummy routine
355      WRITE(*,*) 'dia_wri_c1d: You should not have seen this print! error?', kt, kindic
356   END SUBROUTINE dia_wri_c1d
[253]357#endif
358
359   !!======================================================================
[900]360END MODULE diawri_c1d
Note: See TracBrowser for help on using the repository browser.