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/dev_r4028_CNRS_LIM3_MV2014/NEMOGCM/NEMO/LIM_SRC_3 – NEMO

source: branches/dev_r4028_CNRS_LIM3_MV2014/NEMOGCM/NEMO/LIM_SRC_3/limwri.F90 @ 4566

Last change on this file since 4566 was 4566, checked in by flavoni, 10 years ago

remove control print in limwri.F90,see ticket #1259

  • Property svn:keywords set to Id
File size: 19.5 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 dom_ice
20   USE ice
21   USE limvar
22   USE in_out_manager
23   USE lbclnk
24   USE lib_mpp         ! MPP library
25   USE wrk_nemo        ! work arrays
26   USE par_ice
27   USE iom
28   USE timing          ! Timing
29   USE lib_fortran     ! Fortran utilities
30
31   IMPLICIT NONE
32   PRIVATE
33
34   PUBLIC lim_wri        ! routine called by lim_step.F90
35   PUBLIC lim_wri_state  ! called by dia_wri_state
36
37   REAL(wp)  ::   epsi06 = 1.e-6_wp
38   REAL(wp)  ::   zzero  = 0._wp
39   REAL(wp)  ::   zone   = 1._wp     
40   !!----------------------------------------------------------------------
41   !! NEMO/LIM3 4.0 , UCL - NEMO Consortium (2011)
42   !! $Id$
43   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
44   !!----------------------------------------------------------------------
45CONTAINS
46
47#if defined key_dimgout
48# include "limwri_dimg.h90"
49#else
50
51   SUBROUTINE lim_wri( kindic )
52      !!-------------------------------------------------------------------
53      !!  This routine computes the average of some variables and write it
54      !!  on the ouput files.
55      !!  ATTENTION cette routine n'est valable que si le pas de temps est
56      !!  egale a une fraction entiere de 1 jours.
57      !!  Diff 1-D 3-D : suppress common also included in etat
58      !!                 suppress cmoymo 11-18
59      !!  modif : 03/06/98
60      !!-------------------------------------------------------------------
61      INTEGER, INTENT(in) ::   kindic   ! if kindic < 0 there has been an error somewhere
62      !
63      INTEGER ::  ji, jj, jk, jl  ! dummy loop indices
64      REAL(wp) ::  zinda, zindb, z1_365
65      REAL(wp), POINTER, DIMENSION(:,:,:) ::   zoi, zei
66      REAL(wp), POINTER, DIMENSION(:,:)   :: z2d, z2da, z2db, zind    ! 2D workspace
67      !!-------------------------------------------------------------------
68
69      IF( nn_timing == 1 )  CALL timing_start('limwri')
70
71      CALL wrk_alloc( jpi, jpj, jpl, zoi, zei )
72      CALL wrk_alloc( jpi, jpj     , z2d, z2da, z2db, zind )
73
74
75      !-----------------------------
76      ! Mean category values
77      !-----------------------------
78
79      CALL lim_var_icetm      ! mean sea ice temperature
80
81      CALL lim_var_bv         ! brine volume
82
83      DO jj = 1, jpj          ! presence indicator of ice
84         DO ji = 1, jpi
85            zind(ji,jj)  = MAX( zzero , SIGN( zone , at_i(ji,jj) - epsi06 ) )
86         END DO
87      END DO
88      !
89      CALL iom_put( "iceconc"      , at_i               )        ! ice concentration
90      !                                             
91      IF ( iom_use( "icethic_cea" ) ) THEN                       ! mean ice thickness
92         DO jj = 1, jpj 
93            DO ji = 1, jpi
94               z2d(ji,jj)  = vt_i(ji,jj) / MAX( at_i(ji,jj), epsi06 ) * zind(ji,jj)
95            END DO
96         END DO
97         CALL iom_put( "icethic_cea"  , z2d              )
98      ENDIF
99      CALL iom_put( "icevolu"      , vt_i                )       ! ice volume = mean ice thickness over the cell
100      IF ( iom_use( "snowthic_cea" ) ) THEN                      ! snow thickness = mean snow thickness over the cell
101         DO jj = 1, jpj                                           
102            DO ji = 1, jpi
103               z2d(ji,jj)  = vt_s(ji,jj) / MAX( at_i(ji,jj), epsi06 ) * zind(ji,jj)
104            END DO
105         END DO
106         CALL iom_put( "snowthic_cea" , z2d              )       
107      ENDIF
108      IF ( iom_use( "ioceflxb" ) ) THEN
109         CALL iom_put( "ioceflxb"     , fbif * at_i         )       ! oceanic flux at the ice base
110      ENDIF
111      CALL iom_put( "isst"         , sst_m               )       ! sea surface temperature
112      CALL iom_put( "isss"         , sss_m               )       ! sea surface salinity
113      IF ( iom_use( "qt_oce" ) ) THEN
114         CALL iom_put( "qt_oce"       , qns + qsr           )       ! total flux at ocean surface
115      ENDIF
116      !
117      IF ( iom_use( "uice_ipa" ) .OR. iom_use( "vice_ipa" ) .OR. iom_use( "icevel" )) THEN
118         DO jj = 2 , jpjm1
119            DO ji = 2 , jpim1
120               z2da(ji,jj)  = (  u_ice(ji,jj) * tmu(ji,jj) + u_ice(ji-1,jj) * tmu(ji-1,jj) ) * 0.5_wp
121               z2db(ji,jj)  = (  v_ice(ji,jj) * tmv(ji,jj) + v_ice(ji,jj-1) * tmv(ji,jj-1) ) * 0.5_wp
122           END DO
123         END DO
124         CALL lbc_lnk( z2da, 'T', -1. )
125         CALL lbc_lnk( z2db, 'T', -1. )
126         CALL iom_put( "uice_ipa"     , z2da                )       ! ice velocity u component
127         CALL iom_put( "vice_ipa"     , z2db                )       ! ice velocity v component
128         DO jj = 1, jpj                                 
129            DO ji = 1, jpi
130               z2d(ji,jj)  = SQRT( z2da(ji,jj) * z2da(ji,jj) + z2db(ji,jj) * z2db(ji,jj) ) 
131            END DO
132         END DO
133         CALL iom_put( "icevel"       , z2d                 )       ! ice velocity module
134      ENDIF
135      IF ( iom_use( "icebopr" ) ) THEN
136         CALL iom_put( "icebopr"      , diag_bot_gr * rday  )       ! daily bottom thermodynamic ice production
137      ENDIF
138      IF ( iom_use( "icebopr" ) ) THEN
139         CALL iom_put( "icedypr"      , diag_dyn_gr * rday  )       ! daily dynamic ice production (rid/raft)
140      ENDIF
141!!SF BE CAREFUL : qsr_oce qnd qns_oce are after penetration over ice
142      CALL iom_put( "qsr_oce"      , qsr                 )       ! solar flux at ocean surface
143      CALL iom_put( "qns_oce"      , qns                 )       ! non-solar flux at ocean surface
144!!SF end be careful
145      CALL iom_put( "hfbri"        , fhbri               )       ! heat flux due to brine release
146      CALL iom_put( "utau_ice"     , utau_ice            )       ! wind stress over ice along i-axis at I-point
147      CALL iom_put( "vtau_ice"     , vtau_ice            )       ! wind stress over ice along j-axis at I-point
148!!SF commented because this computation is not ok
149 !SF because qsr is not qsr_ocean but it contains already qsr_ice
150      !SF
151      !SF DO jj = 1 , jpj
152      !SF    DO ji = 1 , jpi
153      !SF       z2d(ji,jj)  = ( 1._wp - at_i(ji,jj) ) * qsr(ji,jj)
154      !SF   END DO
155      !SF END DO
156      !SF CALL iom_put( "qsr_io"       , z2d                 )        ! solar flux at ice/ocean surface
157      !SF DO jj = 1 , jpj
158      !SF    DO ji = 1 , jpi
159      !SF       z2d(ji,jj)  = ( 1._wp - at_i(ji,jj) ) * qns(ji,jj)
160      !SF   END DO
161      !SF END DO
162      !SF CALL iom_put( "qns_io"       , z2d                 )        ! non-solar flux at ice/ocean surface
163      CALL iom_put( "snowpre"      , sprecip             )        ! snow precipitation
164      CALL iom_put( "micesalt"     , smt_i               )        ! mean ice salinity
165      !
166      IF ( iom_use( "miceage" ) ) THEN
167         z2d(:,:) = 0.e0
168         DO jl = 1, jpl
169            DO jj = 1, jpj
170               DO ji = 1, jpi
171                  z2d(ji,jj) = z2d(ji,jj) + zind(ji,jj) * oa_i(ji,jj,jl)
172               END DO
173            END DO
174         END DO
175         z1_365 = 1._wp / 365._wp
176         CALL iom_put( "miceage"     , z2d * z1_365         )        ! mean ice age
177      ENDIF
178      IF ( iom_use( "icelapr" ) ) THEN
179         CALL iom_put( "icelapr"     , diag_lat_gr * rday   )        ! daily lateral thermodynamic ice production
180      ENDIF
181      IF ( iom_use( "icesipr" ) ) THEN
182         CALL iom_put( "icesipr"     , diag_sni_gr * rday   )        ! daily snowice ice production
183      ENDIF
184      IF ( iom_use( "micet" ) ) THEN
185         DO jj = 1, jpj
186            DO ji = 1, jpi
187               z2d(ji,jj) = ( tm_i(ji,jj) - rtt ) * zind(ji,jj)
188            END DO
189         END DO
190         CALL iom_put( "micet"       , z2d                  )        ! mean ice temperature
191      ENDIF
192      CALL iom_put( "icehc"       , et_i                 )        ! ice total heat content
193      CALL iom_put( "isnowhc"     , et_s                 )        ! snow total heat content
194      !
195      IF ( iom_use( "icest" ) ) THEN
196         z2d(:,:) = 0.e0
197         DO jl = 1, jpl
198            DO jj = 1, jpj
199               DO ji = 1, jpi
200                  z2d(ji,jj) = z2d(ji,jj) + zind(ji,jj) * ( t_su(ji,jj,jl) - rtt ) * a_i(ji,jj,jl) / MAX( at_i(ji,jj) , epsi06 )
201               END DO
202            END DO
203         END DO
204         CALL iom_put( "icest"       , z2d                 )        ! ice surface temperature
205      ENDIF
206      IF ( iom_use( "sfxthd" ) ) THEN
207         CALL iom_put( "sfxthd"      , sfx_thd * rday      )        ! equivalent FW salt flux
208      ENDIF
209      IF ( iom_use( "ibrinv" ) ) THEN
210         CALL iom_put( "ibrinv"      , bv_i * 100._wp      )        ! brine volume
211      ENDIF
212      IF ( iom_use( "icecolf" ) ) THEN
213         DO jj = 1, jpj
214            DO ji = 1, jpi
215               zindb  = MAX( zzero , SIGN( zone , at_i(ji,jj) ) )
216               z2d(ji,jj) = hicol(ji,jj) * zindb
217            END DO
218         END DO
219         CALL iom_put( "icecolf"     , z2d                 )        ! frazil ice collection thickness
220      ENDIF
221      IF ( iom_use( "icestr" ) ) THEN
222         CALL iom_put( "icestr"      , strength * 0.001    )        ! ice strength
223      ENDIF
224      IF ( iom_use( "isume" ) ) THEN
225         CALL iom_put( "isume"       , diag_sur_me * rday  )        ! surface melt
226      ENDIF
227      IF ( iom_use( "ibome" ) ) THEN
228         CALL iom_put( "ibome"       , diag_bot_me * rday  )        ! bottom melt
229      ENDIF
230      IF ( iom_use( "idive" ) ) THEN
231         CALL iom_put( "idive"       , divu_i * 1.0e8      )        ! divergence
232      ENDIF
233      IF ( iom_use( "ishear" ) ) THEN
234         CALL iom_put( "ishear"      , shear_i * 1.0e8     )        ! shear
235      ENDIF
236      IF ( iom_use( "icerepr" ) ) THEN
237         CALL iom_put( "icerepr"     , diag_res_pr * rday  )        ! daily prod./melting due to limupdate
238      ENDIF
239      CALL iom_put( "snowvol"     , vt_s                )        ! snow volume
240      IF ( iom_use( "sfxmec" ) ) THEN
241         CALL iom_put( "sfxmec"      , sfx_mec * rday      )        ! salt flux from ridging rafting
242      ENDIF
243      IF ( iom_use( "sfxres" ) ) THEN
244         CALL iom_put( "sfxres"      , sfx_res * rday      )        ! salt flux from limupdate (resultant)
245      ENDIF
246      IF ( iom_use( "icetrp" ) ) THEN
247         CALL iom_put( "icetrp"      , diag_trp_vi * rday  )        ! ice volume transport
248      ENDIF
249
250
251      !--------------------------------
252      ! Output values for each category
253      !--------------------------------
254
255      IF ( iom_use( "iceage_cat" ) ) THEN
256         DO jl = 1, jpl 
257            CALL lbc_lnk( a_i(:,:,jl)  , 'T' ,  1. )
258            CALL lbc_lnk( sm_i(:,:,jl) , 'T' ,  1. )
259            CALL lbc_lnk( oa_i(:,:,jl) , 'T' ,  1. )
260            CALL lbc_lnk( ht_i(:,:,jl) , 'T' ,  1. )
261            CALL lbc_lnk( ht_s(:,:,jl) , 'T' ,  1. )
262         END DO
263
264         DO jl = 1, jpl 
265            DO jj = 1, jpj
266               DO ji = 1, jpi
267                  zinda = MAX( zzero , SIGN( zone , a_i(ji,jj,jl) - epsi06 ) )
268                  zoi(ji,jj,jl) = oa_i(ji,jj,jl)  / MAX( a_i(ji,jj,jl) , epsi06 ) * zinda
269               END DO
270            END DO
271         END DO
272 
273         CALL iom_put( "iceage_cat"     , zoi         )        ! ice age for categories
274      ENDIF
275
276      ! Compute brine volume
277      IF ( iom_use( "iceage_cat" ) ) THEN
278         zei(:,:,:) = 0._wp
279         DO jl = 1, jpl 
280            DO jk = 1, nlay_i
281               DO jj = 1, jpj
282                  DO ji = 1, jpi
283                     zinda = MAX( zzero , SIGN( zone , a_i(ji,jj,jl) - epsi06 ) )
284                     zei(ji,jj,jl) = zei(ji,jj,jl) + 100.0* &
285                        ( - tmut * s_i(ji,jj,jk,jl) / MIN( ( t_i(ji,jj,jk,jl) - rtt ), - epsi06 ) ) * &
286                        zinda / nlay_i
287                  END DO
288               END DO
289            END DO
290         END DO
291
292         DO jl = 1, jpl 
293            CALL lbc_lnk( zei(:,:,jl) , 'T' ,  1. )
294         END DO
295         CALL iom_put( "brinevol_cat"     , zei         )        ! brine volume for categories
296      ENDIF
297
298         CALL iom_put( "iceconc_cat"      , a_i         )        ! area for categories
299         CALL iom_put( "icethic_cat"      , ht_i        )        ! thickness for categories
300         CALL iom_put( "snowthic_cat"     , ht_s        )        ! snow depth for categories
301         CALL iom_put( "salinity_cat"     , sm_i        )        ! salinity for categories
302
303      !     !  Create an output files (output.lim.abort.nc) if S < 0 or u > 20
304      !     m/s
305      !     IF( kindic < 0 )   CALL lim_wri_state( 'output.abort' )
306      !     not yet implemented
307
308      CALL wrk_dealloc( jpi, jpj, jpl, zoi, zei )
309      CALL wrk_dealloc( jpi, jpj     , z2d, zind, z2da, z2db )
310
311      IF( nn_timing == 1 )  CALL timing_stop('limwri')
312     
313   END SUBROUTINE lim_wri
314#endif
315
316 
317   SUBROUTINE lim_wri_state( kt, kid, kh_i )
318      !!---------------------------------------------------------------------
319      !!                 ***  ROUTINE lim_wri_state  ***
320      !!       
321      !! ** Purpose :   create a NetCDF file named cdfile_name which contains
322      !!      the instantaneous ice state and forcing fields for ice model
323      !!        Used to find errors in the initial state or save the last
324      !!      ocean state in case of abnormal end of a simulation
325      !!
326      !! History :
327      !!   4.1  !  2013-06  (C. Rousset)
328      !!----------------------------------------------------------------------
329      INTEGER, INTENT( in ) ::   kt               ! ocean time-step index)
330      INTEGER, INTENT( in ) ::   kid , kh_i       
331      !!----------------------------------------------------------------------
332   
333
334   
335      CALL histdef( kid, "iicethic", "Ice thickness"           , "m"      , jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt )
336      CALL histdef( kid, "iiceconc", "Ice concentration"       , "%"      , jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt )
337      CALL histdef( kid, "iicetemp", "Ice temperature"         , "C"      , jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt )
338      CALL histdef( kid, "iicevelu", "i-Ice speed (I-point)"   , "m/s"    , jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt )
339      CALL histdef( kid, "iicevelv", "j-Ice speed (I-point)"   , "m/s"    , jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt ) 
340      CALL histdef( kid, "iicestru", "i-Wind stress over ice (I-pt)", "Pa", jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt )
341      CALL histdef( kid, "iicestrv", "j-Wind stress over ice (I-pt)", "Pa", jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt ) 
342      CALL histdef( kid, "iicesflx", "Solar flux over ocean"     , "w/m2"   , jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt ) 
343      CALL histdef( kid, "iicenflx", "Non-solar flux over ocean" , "w/m2"   , jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt )
344      CALL histdef( kid, "isnowpre", "Snow precipitation"      , "kg/m2/s", jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt ) 
345      CALL histdef( kid, "iicesali", "Ice salinity"            , "PSU"    , jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt ) 
346      CALL histdef( kid, "iicevolu", "Ice volume"              , "m"      , jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt ) 
347      CALL histdef( kid, "iicedive", "Ice divergence"          , "10-8s-1", jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt ) 
348      CALL histdef( kid, "iicebopr", "Ice bottom production"   , "m/s"      , jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt )
349      CALL histdef( kid, "iicedypr", "Ice dynamic production"  , "m/s"      , jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt )
350      CALL histdef( kid, "iicelapr", "Ice open water prod"     , "m/s"      , jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt )
351      CALL histdef( kid, "iicesipr", "Snow ice production "    , "m/s"      , jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt )
352      CALL histdef( kid, "iicerepr", "Ice prod from limupdate" , "m/s"      , jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt )
353      CALL histdef( kid, "iicebome", "Ice bottom melt"         , "m/s"      , jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt )
354      CALL histdef( kid, "iicesume", "Ice surface melt"        , "m/s"      , jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt )
355      CALL histdef( kid, "iisfxthd", "Salt flux from thermo"   , ""      , jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt )
356      CALL histdef( kid, "iisfxmec", "Salt flux from dynmics"  , ""      , jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt )
357      CALL histdef( kid, "iisfxres", "Salt flux from limupdate", ""      , jpi, jpj, kh_i, 1, 1, 1, -99, 32, "inst(x)", rdt, rdt )
358
359      CALL histend( kid, snc4set )   ! end of the file definition
360
361      CALL histwrite( kid, "iicethic", kt, icethi        , jpi*jpj, (/1/) )   
362      CALL histwrite( kid, "iiceconc", kt, at_i          , jpi*jpj, (/1/) )
363      CALL histwrite( kid, "iicetemp", kt, tm_i - rtt    , jpi*jpj, (/1/) )
364      CALL histwrite( kid, "iicevelu", kt, u_ice          , jpi*jpj, (/1/) )
365      CALL histwrite( kid, "iicevelv", kt, v_ice          , jpi*jpj, (/1/) )
366      CALL histwrite( kid, "iicestru", kt, utau_ice       , jpi*jpj, (/1/) )
367      CALL histwrite( kid, "iicestrv", kt, vtau_ice       , jpi*jpj, (/1/) )
368      CALL histwrite( kid, "iicesflx", kt, qsr , jpi*jpj, (/1/) )
369      CALL histwrite( kid, "iicenflx", kt, qns , jpi*jpj, (/1/) )
370      CALL histwrite( kid, "isnowpre", kt, sprecip        , jpi*jpj, (/1/) )
371      CALL histwrite( kid, "iicesali", kt, smt_i          , jpi*jpj, (/1/) )
372      CALL histwrite( kid, "iicevolu", kt, vt_i           , jpi*jpj, (/1/) )
373      CALL histwrite( kid, "iicedive", kt, divu_i*1.0e8   , jpi*jpj, (/1/) )
374
375      CALL histwrite( kid, "iicebopr", kt, diag_bot_gr        , jpi*jpj, (/1/) )
376      CALL histwrite( kid, "iicedypr", kt, diag_dyn_gr        , jpi*jpj, (/1/) )
377      CALL histwrite( kid, "iicelapr", kt, diag_lat_gr        , jpi*jpj, (/1/) )
378      CALL histwrite( kid, "iicesipr", kt, diag_sni_gr        , jpi*jpj, (/1/) )
379      CALL histwrite( kid, "iicerepr", kt, diag_res_pr        , jpi*jpj, (/1/) )
380      CALL histwrite( kid, "iicebome", kt, diag_bot_me        , jpi*jpj, (/1/) )
381      CALL histwrite( kid, "iicesume", kt, diag_sur_me        , jpi*jpj, (/1/) )
382      CALL histwrite( kid, "iisfxthd", kt, sfx_thd        , jpi*jpj, (/1/) )
383      CALL histwrite( kid, "iisfxmec", kt, sfx_mec        , jpi*jpj, (/1/) )
384      CALL histwrite( kid, "iisfxres", kt, sfx_res        , jpi*jpj, (/1/) )
385
386      ! Close the file
387      ! -----------------
388      CALL histclo( kid )
389
390    END SUBROUTINE lim_wri_state
391
392#else
393   !!----------------------------------------------------------------------
394   !!   Default option :         Empty module          NO LIM sea-ice model
395   !!----------------------------------------------------------------------
396CONTAINS
397   SUBROUTINE lim_wri          ! Empty routine
398   END SUBROUTINE lim_wri
399#endif
400
401   !!======================================================================
402END MODULE limwri
Note: See TracBrowser for help on using the repository browser.