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.
limwri.F90 in branches/2017/dev_r7881_HPC09_ZDF/NEMOGCM/NEMO/LIM_SRC_3 – NEMO

source: branches/2017/dev_r7881_HPC09_ZDF/NEMOGCM/NEMO/LIM_SRC_3/limwri.F90 @ 7899

Last change on this file since 7899 was 7753, checked in by mocavero, 7 years ago

Reverting trunk to remove OpenMP

  • Property svn:keywords set to Id
File size: 21.2 KB
Line 
1MODULE limwri
2   !!======================================================================
3   !!                     ***  MODULE  limwri  ***
4   !!         Ice diagnostics :  write ice output files
5   !!======================================================================
6#if defined key_lim3
7   !!----------------------------------------------------------------------
8   !!   'key_lim3'                                      LIM3 sea-ice model
9   !!----------------------------------------------------------------------
10   !!   lim_wri      : write of the diagnostics variables in ouput file
11   !!   lim_wri_state : write for initial state or/and abandon
12   !!----------------------------------------------------------------------
13   USE ioipsl
14   USE dianam          ! build name of file (routine)
15   USE phycst
16   USE dom_oce
17   USE sbc_oce         ! Surface boundary condition: ocean fields
18   USE sbc_ice         ! Surface boundary condition: ice fields
19   USE ice
20   USE limvar
21   USE in_out_manager
22   USE lbclnk
23   USE lib_mpp         ! MPP library
24   USE wrk_nemo        ! work arrays
25   USE iom
26   USE timing          ! Timing
27   USE lib_fortran     ! Fortran utilities
28
29   IMPLICIT NONE
30   PRIVATE
31
32   PUBLIC lim_wri        ! routine called by lim_step.F90
33   PUBLIC lim_wri_state  ! called by dia_wri_state
34
35   !!----------------------------------------------------------------------
36   !! NEMO/LIM3 4.0 , UCL - NEMO Consortium (2011)
37   !! $Id$
38   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
39   !!----------------------------------------------------------------------
40CONTAINS
41
42
43   SUBROUTINE lim_wri( kindic )
44      !!-------------------------------------------------------------------
45      !!  This routine computes the average of some variables and write it
46      !!  on the ouput files.
47      !!  ATTENTION cette routine n'est valable que si le pas de temps est
48      !!  egale a une fraction entiere de 1 jours.
49      !!  Diff 1-D 3-D : suppress common also included in etat
50      !!                 suppress cmoymo 11-18
51      !!  modif : 03/06/98
52      !!-------------------------------------------------------------------
53      INTEGER, INTENT(in) ::   kindic   ! if kindic < 0 there has been an error somewhere
54      !
55      INTEGER  ::  ji, jj, jk, jl  ! dummy loop indices
56      REAL(wp) ::  z1_365
57      REAL(wp) ::  z2da, z2db, ztmp
58      REAL(wp), POINTER, DIMENSION(:,:,:) ::  zswi2
59      REAL(wp), POINTER, DIMENSION(:,:)   ::  z2d, zswi    ! 2D workspace
60      !!-------------------------------------------------------------------
61
62      IF( nn_timing == 1 )  CALL timing_start('limwri')
63
64      CALL wrk_alloc( jpi,jpj,jpl, zswi2 )
65      CALL wrk_alloc( jpi,jpj    , z2d, zswi )
66
67      !-----------------------------
68      ! Mean category values
69      !-----------------------------
70      z1_365 = 1._wp / 365._wp
71
72      ! brine volume
73      CALL lim_var_bv 
74
75      ! tresholds for outputs
76      DO jj = 1, jpj
77         DO ji = 1, jpi
78            zswi(ji,jj)  = MAX( 0._wp , SIGN( 1._wp , at_i(ji,jj) - epsi06 ) )
79         END DO
80      END DO
81      DO jl = 1, jpl
82         DO jj = 1, jpj
83            DO ji = 1, jpi
84               zswi2(ji,jj,jl)  = MAX( 0._wp , SIGN( 1._wp , a_i(ji,jj,jl) - epsi06 ) )
85            END DO
86         END DO
87      END DO
88      !
89      ! fluxes
90      ! pfrld is the lead fraction at the previous time step (actually between TRP and THD)
91      IF( iom_use('qsr_oce') )   CALL iom_put( "qsr_oce" , qsr_oce(:,:) * pfrld(:,:) )                                   !     solar flux at ocean surface
92      IF( iom_use('qns_oce') )   CALL iom_put( "qns_oce" , qns_oce(:,:) * pfrld(:,:) + qemp_oce(:,:) )                   ! non-solar flux at ocean surface
93      IF( iom_use('qsr_ice') )   CALL iom_put( "qsr_ice" , SUM( qsr_ice(:,:,:) * a_i_b(:,:,:), dim=3 ) )                 !     solar flux at ice surface
94      IF( iom_use('qns_ice') )   CALL iom_put( "qns_ice" , SUM( qns_ice(:,:,:) * a_i_b(:,:,:), dim=3 ) + qemp_ice(:,:) ) ! non-solar flux at ice surface
95      IF( iom_use('qtr_ice') )   CALL iom_put( "qtr_ice" , SUM( ftr_ice(:,:,:) * a_i_b(:,:,:), dim=3 ) )                 !     solar flux transmitted thru ice
96      IF( iom_use('qt_oce' ) )   CALL iom_put( "qt_oce"  , ( qsr_oce(:,:) + qns_oce(:,:) ) * pfrld(:,:) + qemp_oce(:,:) ) 
97      IF( iom_use('qt_ice' ) )   CALL iom_put( "qt_ice"  , SUM( ( qns_ice(:,:,:) + qsr_ice(:,:,:) )   &
98         &                                                      * a_i_b(:,:,:), dim=3 ) + qemp_ice(:,:) )
99      IF( iom_use('qemp_oce') )  CALL iom_put( "qemp_oce" , qemp_oce(:,:) ) 
100      IF( iom_use('qemp_ice') )  CALL iom_put( "qemp_ice" , qemp_ice(:,:) ) 
101      IF( iom_use('emp_oce' ) )  CALL iom_put( "emp_oce"  , emp_oce(:,:) )   ! emp over ocean (taking into account the snow blown away from the ice)
102      IF( iom_use('emp_ice' ) )  CALL iom_put( "emp_ice"  , emp_ice(:,:) )   ! emp over ice   (taking into account the snow blown away from the ice)
103
104      ! velocity
105      IF ( iom_use( "uice_ipa" ) .OR. iom_use( "vice_ipa" ) .OR. iom_use( "icevel" ) ) THEN
106         DO jj = 2 , jpjm1
107            DO ji = 2 , jpim1
108               z2da  = ( u_ice(ji,jj) * umask(ji,jj,1) + u_ice(ji-1,jj) * umask(ji-1,jj,1) ) * 0.5_wp
109               z2db  = ( v_ice(ji,jj) * vmask(ji,jj,1) + v_ice(ji,jj-1) * vmask(ji,jj-1,1) ) * 0.5_wp
110               z2d(ji,jj) = SQRT( z2da * z2da + z2db * z2db )
111           END DO
112         END DO
113         CALL lbc_lnk( z2d, 'T', 1. )
114         CALL iom_put( "uice_ipa"     , u_ice      )       ! ice velocity u component
115         CALL iom_put( "vice_ipa"     , v_ice      )       ! ice velocity v component
116         CALL iom_put( "icevel"       , z2d        )       ! ice velocity module
117      ENDIF
118
119      IF ( iom_use( "tau_icebfr" ) )    CALL iom_put( "tau_icebfr"  , tau_icebfr             )  ! ice friction with ocean bottom (landfast ice) 
120      !
121      IF ( iom_use( "miceage" ) )       CALL iom_put( "miceage"     , om_i * zswi * z1_365   )  ! mean ice age
122      IF ( iom_use( "icethic_cea" ) )   CALL iom_put( "icethic_cea" , htm_i * zswi           )  ! ice thickness mean
123      IF ( iom_use( "snowthic_cea" ) )  CALL iom_put( "snowthic_cea", htm_s * zswi           )  ! snow thickness mean
124      IF ( iom_use( "micet" ) )         CALL iom_put( "micet"       , ( tm_i  - rt0 ) * zswi )  ! ice mean    temperature
125      IF ( iom_use( "icest" ) )         CALL iom_put( "icest"       , ( tm_su - rt0 ) * zswi )  ! ice surface temperature
126      IF ( iom_use( "icecolf" ) )       CALL iom_put( "icecolf"     , hicol                  )  ! frazil ice collection thickness
127      !
128      CALL iom_put( "isst"        , sst_m               )        ! sea surface temperature
129      CALL iom_put( "isss"        , sss_m               )        ! sea surface salinity
130      CALL iom_put( "iceconc"     , at_i  * zswi        )        ! ice concentration
131      CALL iom_put( "icevolu"     , vt_i  * zswi        )        ! ice volume = mean ice thickness over the cell
132      CALL iom_put( "icehc"       , et_i  * zswi        )        ! ice total heat content
133      CALL iom_put( "isnowhc"     , et_s  * zswi        )        ! snow total heat content
134      CALL iom_put( "ibrinv"      , bvm_i * zswi * 100. )        ! brine volume
135      CALL iom_put( "utau_ice"    , utau_ice            )        ! wind stress over ice along i-axis at I-point
136      CALL iom_put( "vtau_ice"    , vtau_ice            )        ! wind stress over ice along j-axis at I-point
137      CALL iom_put( "snowpre"     , sprecip * 86400.    )        ! snow precipitation
138      CALL iom_put( "micesalt"    , smt_i   * zswi      )        ! mean ice salinity
139
140      CALL iom_put( "icestr"      , strength * zswi )    ! ice strength
141      CALL iom_put( "idive"       , divu_i * 1.0e8      )    ! divergence
142      CALL iom_put( "ishear"      , shear_i * 1.0e8     )    ! shear
143      CALL iom_put( "snowvol"     , vt_s   * zswi       )        ! snow volume
144     
145      CALL iom_put( "icetrp"      , diag_trp_vi * rday  )        ! ice volume transport
146      CALL iom_put( "snwtrp"      , diag_trp_vs * rday  )        ! snw volume transport
147      CALL iom_put( "saltrp"      , diag_trp_smv * rday * rhoic ) ! salt content transport
148      CALL iom_put( "deitrp"      , diag_trp_ei         )        ! advected ice enthalpy (W/m2)
149      CALL iom_put( "destrp"      , diag_trp_es         )        ! advected snw enthalpy (W/m2)
150
151      CALL iom_put( "sfxbog"      , sfx_bog * rday      )        ! salt flux from bottom growth
152      CALL iom_put( "sfxbom"      , sfx_bom * rday      )        ! salt flux from bottom melting
153      CALL iom_put( "sfxsum"      , sfx_sum * rday      )        ! salt flux from surface melting
154      CALL iom_put( "sfxlam"      , sfx_lam * rday      )        ! salt flux from lateral melting
155      CALL iom_put( "sfxsni"      , sfx_sni * rday      )        ! salt flux from snow ice formation
156      CALL iom_put( "sfxopw"      , sfx_opw * rday      )        ! salt flux from open water formation
157      CALL iom_put( "sfxdyn"      , sfx_dyn * rday      )        ! salt flux from ridging rafting
158      CALL iom_put( "sfxres"      , sfx_res * rday      )        ! salt flux from limupdate (resultant)
159      CALL iom_put( "sfxbri"      , sfx_bri * rday      )        ! salt flux from brines
160      CALL iom_put( "sfxsub"      , sfx_sub * rday      )        ! salt flux from sublimation
161      CALL iom_put( "sfx"         , sfx     * rday      )        ! total salt flux
162
163      ztmp = rday / rhoic
164      CALL iom_put( "vfxres"     , wfx_res * ztmp       )        ! daily prod./melting due to limupdate
165      CALL iom_put( "vfxopw"     , wfx_opw * ztmp       )        ! daily lateral thermodynamic ice production
166      CALL iom_put( "vfxsni"     , wfx_sni * ztmp       )        ! daily snowice ice production
167      CALL iom_put( "vfxbog"     , wfx_bog * ztmp       )        ! daily bottom thermodynamic ice production
168      CALL iom_put( "vfxdyn"     , wfx_dyn * ztmp       )        ! daily dynamic ice production (rid/raft)
169      CALL iom_put( "vfxsum"     , wfx_sum * ztmp       )        ! surface melt
170      CALL iom_put( "vfxbom"     , wfx_bom * ztmp       )        ! bottom melt
171      CALL iom_put( "vfxlam"     , wfx_lam * ztmp       )        ! lateral melt
172      CALL iom_put( "vfxice"     , wfx_ice * ztmp       )        ! total ice growth/melt
173
174      IF ( iom_use( "vfxthin" ) ) THEN   ! ice production for open water + thin ice (<20cm) => comparable to observations 
175         WHERE( htm_i(:,:) < 0.2 .AND. htm_i(:,:) > 0. ) ; z2d = wfx_bog
176         ELSEWHERE                                       ; z2d = 0._wp
177         END WHERE
178         CALL iom_put( "vfxthin", ( wfx_opw + z2d ) * ztmp )
179      ENDIF
180
181      ztmp = rday / rhosn
182      CALL iom_put( "vfxspr"     , wfx_spr * ztmp       )        ! precip (snow)
183      CALL iom_put( "vfxsnw"     , wfx_snw * ztmp       )        ! total snw growth/melt
184      CALL iom_put( "vfxsub"     , wfx_sub * ztmp       )        ! sublimation (snow/ice)
185      CALL iom_put( "vfxsub_err" , wfx_err_sub * ztmp   )        ! "excess" of sublimation sent to ocean
186     
187      CALL iom_put( "afxtot"     , afx_tot * rday       )        ! concentration tendency (total)
188      CALL iom_put( "afxdyn"     , afx_dyn * rday       )        ! concentration tendency (dynamics)
189      CALL iom_put( "afxthd"     , afx_thd * rday       )        ! concentration tendency (thermo)
190
191      CALL iom_put ('hfxthd'     , hfx_thd(:,:)         )   
192      CALL iom_put ('hfxdyn'     , hfx_dyn(:,:)         )   
193      CALL iom_put ('hfxres'     , hfx_res(:,:)         )   
194      CALL iom_put ('hfxout'     , hfx_out(:,:)         )   
195      CALL iom_put ('hfxin'      , hfx_in(:,:)          )   
196      CALL iom_put ('hfxsnw'     , hfx_snw(:,:)         )   
197      CALL iom_put ('hfxsub'     , hfx_sub(:,:)         )   
198      CALL iom_put ('hfxerr'     , hfx_err(:,:)         )   
199      CALL iom_put ('hfxerr_rem' , hfx_err_rem(:,:)     )   
200     
201      CALL iom_put ('hfxsum'     , hfx_sum(:,:)         )   
202      CALL iom_put ('hfxbom'     , hfx_bom(:,:)         )   
203      CALL iom_put ('hfxbog'     , hfx_bog(:,:)         )   
204      CALL iom_put ('hfxdif'     , hfx_dif(:,:)         )   
205      CALL iom_put ('hfxopw'     , hfx_opw(:,:)         )   
206      CALL iom_put ('hfxtur'     , fhtur(:,:) * at_i_b(:,:) ) ! turbulent heat flux at ice base
207      CALL iom_put ('hfxdhc'     , diag_heat(:,:)       )   ! Heat content variation in snow and ice
208      CALL iom_put ('hfxspr'     , hfx_spr(:,:)         )   ! Heat content of snow precip
209
210     
211      !--------------------------------
212      ! Output values for each category
213      !--------------------------------
214      IF ( iom_use( "iceconc_cat"  ) )  CALL iom_put( "iceconc_cat"      , a_i   * zswi2   )        ! area for categories
215      IF ( iom_use( "icethic_cat"  ) )  CALL iom_put( "icethic_cat"      , ht_i  * zswi2   )        ! thickness for categories
216      IF ( iom_use( "snowthic_cat" ) )  CALL iom_put( "snowthic_cat"     , ht_s  * zswi2   )        ! snow depth for categories
217      IF ( iom_use( "salinity_cat" ) )  CALL iom_put( "salinity_cat"     , sm_i  * zswi2   )        ! salinity for categories
218      ! ice temperature
219      IF ( iom_use( "icetemp_cat"  ) )  CALL iom_put( "icetemp_cat", ( SUM( t_i(:,:,:,:), dim=3 ) * r1_nlay_i - rt0 ) * zswi2 )
220      ! snow temperature
221      IF ( iom_use( "snwtemp_cat"  ) )  CALL iom_put( "snwtemp_cat", ( SUM( t_s(:,:,:,:), dim=3 ) * r1_nlay_s - rt0 ) * zswi2 )
222      ! ice age
223      IF ( iom_use( "iceage_cat"   ) )  CALL iom_put( "iceage_cat" , o_i * zswi2 * z1_365 )
224      ! brine volume
225      IF ( iom_use( "brinevol_cat" ) )  CALL iom_put( "brinevol_cat", bv_i * 100. * zswi2 )
226
227      !     !  Create an output files (output.lim.abort.nc) if S < 0 or u > 20 m/s
228      !     IF( kindic < 0 )   CALL lim_wri_state( 'output.abort' )
229      !     not yet implemented
230     
231      CALL wrk_dealloc( jpi, jpj, jpl, zswi2 )
232      CALL wrk_dealloc( jpi, jpj     , z2d, zswi )
233
234      IF( nn_timing == 1 )  CALL timing_stop('limwri')
235     
236   END SUBROUTINE lim_wri
237
238 
239   SUBROUTINE lim_wri_state( kt, kid, kh_i )
240      !!---------------------------------------------------------------------
241      !!                 ***  ROUTINE lim_wri_state  ***
242      !!       
243      !! ** Purpose :   create a NetCDF file named cdfile_name which contains
244      !!      the instantaneous ice state and forcing fields for ice model
245      !!        Used to find errors in the initial state or save the last
246      !!      ocean state in case of abnormal end of a simulation
247      !!
248      !! History :
249      !!   4.0  !  2013-06  (C. Rousset)
250      !!----------------------------------------------------------------------
251      INTEGER, INTENT( in )   ::   kt               ! ocean time-step index)
252      INTEGER, INTENT( in )   ::   kid , kh_i
253      INTEGER                 ::   nz_i, jl
254      REAL(wp), DIMENSION(jpl) :: jcat
255      !!----------------------------------------------------------------------
256      DO jl = 1, jpl
257         jcat(jl) = REAL(jl)
258      ENDDO
259     
260      CALL histvert( kid, "ncatice", "Ice Categories","", jpl, jcat, nz_i, "up")
261
262      CALL histdef( kid, "sithic", "Ice thickness"           , "m"      ,   &
263      &      jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt )
264      CALL histdef( kid, "siconc", "Ice concentration"       , "%"      ,   &
265      &      jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt )
266      CALL histdef( kid, "sitemp", "Ice temperature"         , "C"      ,   &
267      &      jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt )
268      CALL histdef( kid, "sivelu", "i-Ice speed "            , "m/s"    ,   &
269      &      jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt )
270      CALL histdef( kid, "sivelv", "j-Ice speed "            , "m/s"    ,   &
271      &      jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt ) 
272      CALL histdef( kid, "sistru", "i-Wind stress over ice " , "Pa"     ,   &
273      &      jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt )
274      CALL histdef( kid, "sistrv", "j-Wind stress over ice " , "Pa"     ,   &
275      &      jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt ) 
276      CALL histdef( kid, "sisflx", "Solar flux over ocean"     , "w/m2" ,   &
277      &      jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt ) 
278      CALL histdef( kid, "sinflx", "Non-solar flux over ocean" , "w/m2" ,   &
279      &      jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt )
280      CALL histdef( kid, "isnowpre", "Snow precipitation"      , "kg/m2/s",   &
281      &      jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt ) 
282      CALL histdef( kid, "sisali", "Ice salinity"            , "PSU"    ,   &
283      &      jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt ) 
284      CALL histdef( kid, "sivolu", "Ice volume"              , "m"      ,   &
285      &      jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt ) 
286      CALL histdef( kid, "sidive", "Ice divergence"          , "10-8s-1",   &
287      &      jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt ) 
288
289      CALL histdef( kid, "vfxbog", "Ice bottom production"   , "m/s"    ,   &
290      &      jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt )
291      CALL histdef( kid, "vfxdyn", "Ice dynamic production"  , "m/s"    ,   &
292      &      jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt )
293      CALL histdef( kid, "vfxopw", "Ice open water prod"     , "m/s"    ,   &
294      &       jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt )
295      CALL histdef( kid, "vfxsni", "Snow ice production "    , "m/s"    ,   &
296      &      jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt )
297      CALL histdef( kid, "vfxres", "Ice prod from limupdate" , "m/s"    ,   &
298      &      jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt )
299      CALL histdef( kid, "vfxbom", "Ice bottom melt"         , "m/s"    ,   &
300      &      jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt )
301      CALL histdef( kid, "vfxsum", "Ice surface melt"        , "m/s"    ,   &
302      &      jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt )
303
304      CALL histdef( kid, "sithicat", "Ice thickness"         , "m"      ,   &
305      &      jpi, jpj, kh_i, jpl, 1, jpl, nz_i, 32, "inst(x)", rdt, rdt )
306      CALL histdef( kid, "siconcat", "Ice concentration"     , "%"      ,   &
307      &      jpi, jpj, kh_i, jpl, 1, jpl, nz_i, 32, "inst(x)", rdt, rdt )
308      CALL histdef( kid, "sisalcat", "Ice salinity"           , ""      ,   &
309      &      jpi, jpj, kh_i, jpl, 1, jpl, nz_i, 32, "inst(x)", rdt, rdt )
310      CALL histdef( kid, "sitemcat", "Ice temperature"       , "C"      ,   &
311      &      jpi, jpj, kh_i, jpl, 1, jpl, nz_i, 32, "inst(x)", rdt, rdt )
312      CALL histdef( kid, "snthicat", "Snw thickness"         , "m"      ,   &
313      &      jpi, jpj, kh_i, jpl, 1, jpl, nz_i, 32, "inst(x)", rdt, rdt )
314      CALL histdef( kid, "sntemcat", "Snw temperature"       , "C"      ,   &
315      &      jpi, jpj, kh_i, jpl, 1, jpl, nz_i, 32, "inst(x)", rdt, rdt )
316
317      CALL histend( kid, snc4set )   ! end of the file definition
318
319      CALL histwrite( kid, "sithic", kt, htm_i         , jpi*jpj, (/1/) )   
320      CALL histwrite( kid, "siconc", kt, at_i          , jpi*jpj, (/1/) )
321      CALL histwrite( kid, "sitemp", kt, tm_i - rt0    , jpi*jpj, (/1/) )
322      CALL histwrite( kid, "sivelu", kt, u_ice          , jpi*jpj, (/1/) )
323      CALL histwrite( kid, "sivelv", kt, v_ice          , jpi*jpj, (/1/) )
324      CALL histwrite( kid, "sistru", kt, utau_ice       , jpi*jpj, (/1/) )
325      CALL histwrite( kid, "sistrv", kt, vtau_ice       , jpi*jpj, (/1/) )
326      CALL histwrite( kid, "sisflx", kt, qsr , jpi*jpj, (/1/) )
327      CALL histwrite( kid, "sinflx", kt, qns , jpi*jpj, (/1/) )
328      CALL histwrite( kid, "isnowpre", kt, sprecip        , jpi*jpj, (/1/) )
329      CALL histwrite( kid, "sisali", kt, smt_i          , jpi*jpj, (/1/) )
330      CALL histwrite( kid, "sivolu", kt, vt_i           , jpi*jpj, (/1/) )
331      CALL histwrite( kid, "sidive", kt, divu_i*1.0e8   , jpi*jpj, (/1/) )
332
333      CALL histwrite( kid, "vfxbog", kt, wfx_bog        , jpi*jpj, (/1/) )
334      CALL histwrite( kid, "vfxdyn", kt, wfx_dyn        , jpi*jpj, (/1/) )
335      CALL histwrite( kid, "vfxopw", kt, wfx_opw        , jpi*jpj, (/1/) )
336      CALL histwrite( kid, "vfxsni", kt, wfx_sni        , jpi*jpj, (/1/) )
337      CALL histwrite( kid, "vfxres", kt, wfx_res        , jpi*jpj, (/1/) )
338      CALL histwrite( kid, "vfxbom", kt, wfx_bom        , jpi*jpj, (/1/) )
339      CALL histwrite( kid, "vfxsum", kt, wfx_sum        , jpi*jpj, (/1/) )
340
341      CALL histwrite( kid, "sithicat", kt, ht_i        , jpi*jpj*jpl, (/1/) )   
342      CALL histwrite( kid, "siconcat", kt, a_i         , jpi*jpj*jpl, (/1/) )   
343      CALL histwrite( kid, "sisalcat", kt, sm_i        , jpi*jpj*jpl, (/1/) )   
344      CALL histwrite( kid, "sitemcat", kt, tm_i - rt0  , jpi*jpj*jpl, (/1/) )   
345      CALL histwrite( kid, "snthicat", kt, ht_s        , jpi*jpj*jpl, (/1/) )   
346      CALL histwrite( kid, "sntemcat", kt, tm_su - rt0 , jpi*jpj*jpl, (/1/) )   
347
348      ! Close the file
349      ! -----------------
350      !CALL histclo( kid )
351
352    END SUBROUTINE lim_wri_state
353
354#else
355   !!----------------------------------------------------------------------
356   !!   Default option :         Empty module          NO LIM sea-ice model
357   !!----------------------------------------------------------------------
358CONTAINS
359   SUBROUTINE lim_wri          ! Empty routine
360   END SUBROUTINE lim_wri
361#endif
362
363   !!======================================================================
364END MODULE limwri
Note: See TracBrowser for help on using the repository browser.