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.
asminc.F90 in branches/2016/dev_v3_6_STABLE_r6506_AGRIF_LIM3/NEMOGCM/NEMO/OPA_SRC/ASM – NEMO

source: branches/2016/dev_v3_6_STABLE_r6506_AGRIF_LIM3/NEMOGCM/NEMO/OPA_SRC/ASM/asminc.F90 @ 7158

Last change on this file since 7158 was 7158, checked in by clem, 7 years ago

debug branch

  • Property svn:keywords set to Id
File size: 48.8 KB
Line 
1MODULE asminc
2   !!======================================================================
3   !!                       ***  MODULE asminc  ***
4   !! Assimilation increment : Apply an increment generated by data
5   !!                          assimilation
6   !!======================================================================
7   !! History :       ! 2007-03  (M. Martin)  Met Office version
8   !!                 ! 2007-04  (A. Weaver)  calc_date original code
9   !!                 ! 2007-04  (A. Weaver)  Merge with OPAVAR/NEMOVAR
10   !!   NEMO     3.3  ! 2010-05  (D. Lea)  Update to work with NEMO v3.2
11   !!             -   ! 2010-05  (D. Lea)  add calc_month_len routine based on day_init
12   !!            3.4  ! 2012-10  (A. Weaver and K. Mogensen) Fix for direct initialization
13   !!----------------------------------------------------------------------
14
15   !!----------------------------------------------------------------------
16   !!   'key_asminc'   : Switch on the assimilation increment interface
17   !!----------------------------------------------------------------------
18   !!   asm_inc_init   : Initialize the increment arrays and IAU weights
19   !!   calc_date      : Compute the calendar date YYYYMMDD on a given step
20   !!   tra_asm_inc    : Apply the tracer (T and S) increments
21   !!   dyn_asm_inc    : Apply the dynamic (u and v) increments
22   !!   ssh_asm_inc    : Apply the SSH increment
23   !!   seaice_asm_inc : Apply the seaice increment
24   !!----------------------------------------------------------------------
25   USE wrk_nemo         ! Memory Allocation
26   USE par_oce          ! Ocean space and time domain variables
27   USE dom_oce          ! Ocean space and time domain
28   USE domvvl           ! domain: variable volume level
29   USE oce              ! Dynamics and active tracers defined in memory
30   USE ldfdyn_oce       ! ocean dynamics: lateral physics
31   USE eosbn2           ! Equation of state - in situ and potential density
32   USE zpshde           ! Partial step : Horizontal Derivative
33   USE iom              ! Library to read input files
34   USE asmpar           ! Parameters for the assmilation interface
35   USE c1d              ! 1D initialization
36   USE in_out_manager   ! I/O manager
37   USE lib_mpp          ! MPP library
38#if defined key_lim2
39   USE ice_2            ! LIM2
40#endif
41   USE sbc_oce          ! Surface boundary condition variables.
42
43   IMPLICIT NONE
44   PRIVATE
45   
46   PUBLIC   asm_inc_init   !: Initialize the increment arrays and IAU weights
47   PUBLIC   calc_date      !: Compute the calendar date YYYYMMDD on a given step
48   PUBLIC   tra_asm_inc    !: Apply the tracer (T and S) increments
49   PUBLIC   dyn_asm_inc    !: Apply the dynamic (u and v) increments
50   PUBLIC   ssh_asm_inc    !: Apply the SSH increment
51   PUBLIC   seaice_asm_inc !: Apply the seaice increment
52
53#if defined key_asminc
54    LOGICAL, PUBLIC, PARAMETER :: lk_asminc = .TRUE.   !: Logical switch for assimilation increment interface
55#else
56    LOGICAL, PUBLIC, PARAMETER :: lk_asminc = .FALSE.  !: No assimilation increments
57#endif
58   LOGICAL, PUBLIC :: ln_bkgwri = .FALSE.      !: No output of the background state fields
59   LOGICAL, PUBLIC :: ln_asmiau = .FALSE.      !: No applying forcing with an assimilation increment
60   LOGICAL, PUBLIC :: ln_asmdin = .FALSE.      !: No direct initialization
61   LOGICAL, PUBLIC :: ln_trainc = .FALSE.      !: No tracer (T and S) assimilation increments
62   LOGICAL, PUBLIC :: ln_dyninc = .FALSE.      !: No dynamics (u and v) assimilation increments
63   LOGICAL, PUBLIC :: ln_sshinc = .FALSE.      !: No sea surface height assimilation increment
64   LOGICAL, PUBLIC :: ln_seaiceinc             !: No sea ice concentration increment
65   LOGICAL, PUBLIC :: ln_salfix = .FALSE.      !: Apply minimum salinity check
66   LOGICAL, PUBLIC :: ln_temnofreeze = .FALSE. !: Don't allow the temperature to drop below freezing
67   INTEGER, PUBLIC :: nn_divdmp                !: Apply divergence damping filter nn_divdmp times
68
69   REAL(wp), PUBLIC, DIMENSION(:,:,:), ALLOCATABLE ::   t_bkg   , s_bkg      !: Background temperature and salinity
70   REAL(wp), PUBLIC, DIMENSION(:,:,:), ALLOCATABLE ::   u_bkg   , v_bkg      !: Background u- & v- velocity components
71   REAL(wp), PUBLIC, DIMENSION(:,:,:), ALLOCATABLE ::   t_bkginc, s_bkginc   !: Increment to the background T & S
72   REAL(wp), PUBLIC, DIMENSION(:,:,:), ALLOCATABLE ::   u_bkginc, v_bkginc   !: Increment to the u- & v-components
73   REAL(wp), PUBLIC, DIMENSION(:)    , ALLOCATABLE ::   wgtiau               !: IAU weights for each time step
74#if defined key_asminc
75   REAL(wp), PUBLIC, DIMENSION(:,:), ALLOCATABLE ::   ssh_iau           !: IAU-weighted sea surface height increment
76#endif
77   !                                !!! time steps relative to the cycle interval [0,nitend-nit000-1]
78   INTEGER , PUBLIC ::   nitbkg      !: Time step of the background state used in the Jb term
79   INTEGER , PUBLIC ::   nitdin      !: Time step of the background state for direct initialization
80   INTEGER , PUBLIC ::   nitiaustr   !: Time step of the start of the IAU interval
81   INTEGER , PUBLIC ::   nitiaufin   !: Time step of the end of the IAU interval
82   !
83   INTEGER , PUBLIC ::   niaufn      !: Type of IAU weighing function: = 0   Constant weighting
84   !                                 !: = 1   Linear hat-like, centred in middle of IAU interval
85   REAL(wp), PUBLIC ::   salfixmin   !: Ensure that the salinity is larger than this  value if (ln_salfix)
86
87   REAL(wp), DIMENSION(:,:), ALLOCATABLE ::   ssh_bkg, ssh_bkginc   ! Background sea surface height and its increment
88   REAL(wp), DIMENSION(:,:), ALLOCATABLE ::   seaice_bkginc         ! Increment to the background sea ice conc
89
90   !! * Substitutions
91#  include "domzgr_substitute.h90"
92#  include "ldfdyn_substitute.h90"
93#  include "vectopt_loop_substitute.h90"
94   !!----------------------------------------------------------------------
95   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
96   !! $Id$
97   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
98   !!----------------------------------------------------------------------
99CONTAINS
100
101   SUBROUTINE asm_inc_init
102      !!----------------------------------------------------------------------
103      !!                    ***  ROUTINE asm_inc_init  ***
104      !!         
105      !! ** Purpose : Initialize the assimilation increment and IAU weights.
106      !!
107      !! ** Method  : Initialize the assimilation increment and IAU weights.
108      !!
109      !! ** Action  :
110      !!----------------------------------------------------------------------
111      INTEGER :: ji, jj, jk, jt  ! dummy loop indices
112      INTEGER :: imid, inum      ! local integers
113      INTEGER :: ios             ! Local integer output status for namelist read
114      INTEGER :: iiauper         ! Number of time steps in the IAU period
115      INTEGER :: icycper         ! Number of time steps in the cycle
116      INTEGER :: iitend_date     ! Date YYYYMMDD of final time step
117      INTEGER :: iitbkg_date     ! Date YYYYMMDD of background time step for Jb term
118      INTEGER :: iitdin_date     ! Date YYYYMMDD of background time step for DI
119      INTEGER :: iitiaustr_date  ! Date YYYYMMDD of IAU interval start time step
120      INTEGER :: iitiaufin_date  ! Date YYYYMMDD of IAU interval final time step
121      !
122      REAL(wp) :: znorm        ! Normalization factor for IAU weights
123      REAL(wp) :: ztotwgt      ! Value of time-integrated IAU weights (should be equal to one)
124      REAL(wp) :: z_inc_dateb  ! Start date of interval on which increment is valid
125      REAL(wp) :: z_inc_datef  ! End date of interval on which increment is valid
126      REAL(wp) :: zdate_bkg    ! Date in background state file for DI
127      REAL(wp) :: zdate_inc    ! Time axis in increments file
128      !
129      REAL(wp), POINTER, DIMENSION(:,:) ::   hdiv   ! 2D workspace
130      !!
131      NAMELIST/nam_asminc/ ln_bkgwri,                                      &
132         &                 ln_trainc, ln_dyninc, ln_sshinc,                &
133         &                 ln_asmdin, ln_asmiau,                           &
134         &                 nitbkg, nitdin, nitiaustr, nitiaufin, niaufn,   &
135         &                 ln_salfix, salfixmin, nn_divdmp
136      !!----------------------------------------------------------------------
137
138      !-----------------------------------------------------------------------
139      ! Read Namelist nam_asminc : assimilation increment interface
140      !-----------------------------------------------------------------------
141      ln_seaiceinc = .FALSE.
142      ln_temnofreeze = .FALSE.
143
144      REWIND( numnam_ref )              ! Namelist nam_asminc in reference namelist : Assimilation increment
145      READ  ( numnam_ref, nam_asminc, IOSTAT = ios, ERR = 901)
146901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'nam_asminc in reference namelist', lwp )
147
148      REWIND( numnam_cfg )              ! Namelist nam_asminc in configuration namelist : Assimilation increment
149      READ  ( numnam_cfg, nam_asminc, IOSTAT = ios, ERR = 902 )
150902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'nam_asminc in configuration namelist', lwp )
151      IF(lwm) WRITE ( numond, nam_asminc )
152
153      ! Control print
154      IF(lwp) THEN
155         WRITE(numout,*)
156         WRITE(numout,*) 'asm_inc_init : Assimilation increment initialization :'
157         WRITE(numout,*) '~~~~~~~~~~~~'
158         WRITE(numout,*) '   Namelist namasm : set assimilation increment parameters'
159         WRITE(numout,*) '      Logical switch for writing out background state          ln_bkgwri = ', ln_bkgwri
160         WRITE(numout,*) '      Logical switch for applying tracer increments            ln_trainc = ', ln_trainc
161         WRITE(numout,*) '      Logical switch for applying velocity increments          ln_dyninc = ', ln_dyninc
162         WRITE(numout,*) '      Logical switch for applying SSH increments               ln_sshinc = ', ln_sshinc
163         WRITE(numout,*) '      Logical switch for Direct Initialization (DI)            ln_asmdin = ', ln_asmdin
164         WRITE(numout,*) '      Logical switch for applying sea ice increments        ln_seaiceinc = ', ln_seaiceinc
165         WRITE(numout,*) '      Logical switch for Incremental Analysis Updating (IAU)   ln_asmiau = ', ln_asmiau
166         WRITE(numout,*) '      Timestep of background in [0,nitend-nit000-1]            nitbkg    = ', nitbkg
167         WRITE(numout,*) '      Timestep of background for DI in [0,nitend-nit000-1]     nitdin    = ', nitdin
168         WRITE(numout,*) '      Timestep of start of IAU interval in [0,nitend-nit000-1] nitiaustr = ', nitiaustr
169         WRITE(numout,*) '      Timestep of end of IAU interval in [0,nitend-nit000-1]   nitiaufin = ', nitiaufin
170         WRITE(numout,*) '      Type of IAU weighting function                           niaufn    = ', niaufn
171         WRITE(numout,*) '      Logical switch for ensuring that the sa > salfixmin      ln_salfix = ', ln_salfix
172         WRITE(numout,*) '      Minimum salinity after applying the increments           salfixmin = ', salfixmin
173      ENDIF
174
175      nitbkg_r    = nitbkg    + nit000 - 1  ! Background time referenced to nit000
176      nitdin_r    = nitdin    + nit000 - 1  ! Background time for DI referenced to nit000
177      nitiaustr_r = nitiaustr + nit000 - 1  ! Start of IAU interval referenced to nit000
178      nitiaufin_r = nitiaufin + nit000 - 1  ! End of IAU interval referenced to nit000
179
180      iiauper = nitiaufin_r - nitiaustr_r + 1  ! IAU interval length
181      icycper = nitend      - nit000      + 1  ! Cycle interval length
182
183      CALL calc_date( nit000, nitend     , ndate0, iitend_date    )     ! Date of final time step
184      CALL calc_date( nit000, nitbkg_r   , ndate0, iitbkg_date    )     ! Background time for Jb referenced to ndate0
185      CALL calc_date( nit000, nitdin_r   , ndate0, iitdin_date    )     ! Background time for DI referenced to ndate0
186      CALL calc_date( nit000, nitiaustr_r, ndate0, iitiaustr_date )     ! IAU start time referenced to ndate0
187      CALL calc_date( nit000, nitiaufin_r, ndate0, iitiaufin_date )     ! IAU end time referenced to ndate0
188      !
189      IF(lwp) THEN
190         WRITE(numout,*)
191         WRITE(numout,*) '   Time steps referenced to current cycle:'
192         WRITE(numout,*) '       iitrst      = ', nit000 - 1
193         WRITE(numout,*) '       nit000      = ', nit000
194         WRITE(numout,*) '       nitend      = ', nitend
195         WRITE(numout,*) '       nitbkg_r    = ', nitbkg_r
196         WRITE(numout,*) '       nitdin_r    = ', nitdin_r
197         WRITE(numout,*) '       nitiaustr_r = ', nitiaustr_r
198         WRITE(numout,*) '       nitiaufin_r = ', nitiaufin_r
199         WRITE(numout,*)
200         WRITE(numout,*) '   Dates referenced to current cycle:'
201         WRITE(numout,*) '       ndastp         = ', ndastp
202         WRITE(numout,*) '       ndate0         = ', ndate0
203         WRITE(numout,*) '       iitend_date    = ', iitend_date
204         WRITE(numout,*) '       iitbkg_date    = ', iitbkg_date
205         WRITE(numout,*) '       iitdin_date    = ', iitdin_date
206         WRITE(numout,*) '       iitiaustr_date = ', iitiaustr_date
207         WRITE(numout,*) '       iitiaufin_date = ', iitiaufin_date
208      ENDIF
209
210      IF ( ( ln_asmdin ).AND.( ln_asmiau ) )   &
211         & CALL ctl_stop( ' ln_asmdin and ln_asmiau :', &
212         &                ' Choose Direct Initialization OR Incremental Analysis Updating')
213
214      IF (      ( ( .NOT. ln_asmdin ).AND.( .NOT. ln_asmiau ) ) &
215           .AND.( ( ln_trainc ).OR.( ln_dyninc ).OR.( ln_sshinc ) .OR. ( ln_seaiceinc) )) &
216         & CALL ctl_stop( ' One or more of ln_trainc, ln_dyninc, ln_sshinc and ln_seaiceinc is set to .true.', &
217         &                ' but ln_asmdin and ln_asmiau are both set to .false. :', &
218         &                ' Inconsistent options')
219
220      IF ( ( niaufn /= 0 ).AND.( niaufn /= 1 ) ) &
221         & CALL ctl_stop( ' niaufn /= 0 or niaufn /=1 :',  &
222         &                ' Type IAU weighting function is invalid')
223
224      IF ( ( .NOT. ln_trainc ).AND.( .NOT. ln_dyninc ).AND.( .NOT. ln_sshinc ).AND.( .NOT. ln_seaiceinc ) &
225         &                     )  &
226         & CALL ctl_warn( ' ln_trainc, ln_dyninc, ln_sshinc and ln_seaiceinc are set to .false. :', &
227         &                ' The assimilation increments are not applied')
228
229      IF ( ( ln_asmiau ).AND.( nitiaustr == nitiaufin ) ) &
230         & CALL ctl_stop( ' nitiaustr = nitiaufin :',  &
231         &                ' IAU interval is of zero length')
232
233      IF ( ( ln_asmiau ).AND.( ( nitiaustr_r < nit000 ).OR.( nitiaufin_r > nitend ) ) ) &
234         & CALL ctl_stop( ' nitiaustr or nitiaufin :',  &
235         &                ' IAU starting or final time step is outside the cycle interval', &
236         &                 ' Valid range nit000 to nitend')
237
238      IF ( ( nitbkg_r < nit000 - 1 ).OR.( nitbkg_r > nitend ) ) &
239         & CALL ctl_stop( ' nitbkg :',  &
240         &                ' Background time step is outside the cycle interval')
241
242      IF ( ( nitdin_r < nit000 - 1 ).OR.( nitdin_r > nitend ) ) &
243         & CALL ctl_stop( ' nitdin :',  &
244         &                ' Background time step for Direct Initialization is outside', &
245         &                ' the cycle interval')
246
247      IF ( nstop > 0 ) RETURN       ! if there are any errors then go no further
248
249      !--------------------------------------------------------------------
250      ! Initialize the Incremental Analysis Updating weighting function
251      !--------------------------------------------------------------------
252
253      IF ( ln_asmiau ) THEN
254
255         ALLOCATE( wgtiau( icycper ) )
256
257         wgtiau(:) = 0.0
258
259         IF ( niaufn == 0 ) THEN
260
261            !---------------------------------------------------------
262            ! Constant IAU forcing
263            !---------------------------------------------------------
264
265            DO jt = 1, iiauper
266               wgtiau(jt+nitiaustr-1) = 1.0 / REAL( iiauper )
267            END DO
268
269         ELSEIF ( niaufn == 1 ) THEN
270
271            !---------------------------------------------------------
272            ! Linear hat-like, centred in middle of IAU interval
273            !---------------------------------------------------------
274
275            ! Compute the normalization factor
276            znorm = 0.0
277            IF ( MOD( iiauper, 2 ) == 0 ) THEN  ! Even number of time steps in IAU interval
278               imid = iiauper / 2 
279               DO jt = 1, imid
280                  znorm = znorm + REAL( jt )
281               END DO
282               znorm = 2.0 * znorm
283            ELSE                               ! Odd number of time steps in IAU interval
284               imid = ( iiauper + 1 ) / 2       
285               DO jt = 1, imid - 1
286                  znorm = znorm + REAL( jt )
287               END DO
288               znorm = 2.0 * znorm + REAL( imid )
289            ENDIF
290            znorm = 1.0 / znorm
291
292            DO jt = 1, imid - 1
293               wgtiau(jt+nitiaustr-1) = REAL( jt ) * znorm
294            END DO
295            DO jt = imid, iiauper
296               wgtiau(jt+nitiaustr-1) = REAL( iiauper - jt + 1 ) * znorm
297            END DO
298
299         ENDIF
300
301         ! Test that the integral of the weights over the weighting interval equals 1
302          IF(lwp) THEN
303             WRITE(numout,*)
304             WRITE(numout,*) 'asm_inc_init : IAU weights'
305             WRITE(numout,*) '~~~~~~~~~~~~'
306             WRITE(numout,*) '             time step         IAU  weight'
307             WRITE(numout,*) '             =========     ====================='
308             ztotwgt = 0.0
309             DO jt = 1, icycper
310                ztotwgt = ztotwgt + wgtiau(jt)
311                WRITE(numout,*) '         ', jt, '       ', wgtiau(jt) 
312             END DO   
313             WRITE(numout,*) '         ==================================='
314             WRITE(numout,*) '         Time-integrated weight = ', ztotwgt
315             WRITE(numout,*) '         ==================================='
316          ENDIF
317         
318      ENDIF
319
320      !--------------------------------------------------------------------
321      ! Allocate and initialize the increment arrays
322      !--------------------------------------------------------------------
323
324      ALLOCATE( t_bkginc(jpi,jpj,jpk) )
325      ALLOCATE( s_bkginc(jpi,jpj,jpk) )
326      ALLOCATE( u_bkginc(jpi,jpj,jpk) )
327      ALLOCATE( v_bkginc(jpi,jpj,jpk) )
328      ALLOCATE( ssh_bkginc(jpi,jpj)   )
329      ALLOCATE( seaice_bkginc(jpi,jpj))
330#if defined key_asminc
331      ALLOCATE( ssh_iau(jpi,jpj)      )
332#endif
333      t_bkginc(:,:,:) = 0.0
334      s_bkginc(:,:,:) = 0.0
335      u_bkginc(:,:,:) = 0.0
336      v_bkginc(:,:,:) = 0.0
337      ssh_bkginc(:,:) = 0.0
338      seaice_bkginc(:,:) = 0.0
339#if defined key_asminc
340      ssh_iau(:,:)    = 0.0
341#endif
342      IF ( ( ln_trainc ).OR.( ln_dyninc ).OR.( ln_sshinc ).OR.( ln_seaiceinc ) ) THEN
343
344         !--------------------------------------------------------------------
345         ! Read the increments from file
346         !--------------------------------------------------------------------
347
348         CALL iom_open( c_asminc, inum )
349
350         CALL iom_get( inum, 'time', zdate_inc ) 
351
352         CALL iom_get( inum, 'z_inc_dateb', z_inc_dateb )
353         CALL iom_get( inum, 'z_inc_datef', z_inc_datef )
354         z_inc_dateb = zdate_inc
355         z_inc_datef = zdate_inc
356
357         IF(lwp) THEN
358            WRITE(numout,*) 
359            WRITE(numout,*) 'asm_inc_init : Assimilation increments valid ', &
360               &            ' between dates ', NINT( z_inc_dateb ),' and ',  &
361               &            NINT( z_inc_datef )
362            WRITE(numout,*) '~~~~~~~~~~~~'
363         ENDIF
364
365         IF (     ( NINT( z_inc_dateb ) < ndastp      ) &
366            & .OR.( NINT( z_inc_datef ) > iitend_date ) ) &
367            & CALL ctl_warn( ' Validity time of assimilation increments is ', &
368            &                ' outside the assimilation interval' )
369
370         IF ( ( ln_asmdin ).AND.( NINT( zdate_inc ) /= iitdin_date ) ) &
371            & CALL ctl_warn( ' Validity time of assimilation increments does ', &
372            &                ' not agree with Direct Initialization time' )
373
374         IF ( ln_trainc ) THEN   
375            CALL iom_get( inum, jpdom_autoglo, 'bckint', t_bkginc, 1 )
376            CALL iom_get( inum, jpdom_autoglo, 'bckins', s_bkginc, 1 )
377            ! Apply the masks
378            t_bkginc(:,:,:) = t_bkginc(:,:,:) * tmask(:,:,:)
379            s_bkginc(:,:,:) = s_bkginc(:,:,:) * tmask(:,:,:)
380            ! Set missing increments to 0.0 rather than 1e+20
381            ! to allow for differences in masks
382            WHERE( ABS( t_bkginc(:,:,:) ) > 1.0e+10 ) t_bkginc(:,:,:) = 0.0
383            WHERE( ABS( s_bkginc(:,:,:) ) > 1.0e+10 ) s_bkginc(:,:,:) = 0.0
384         ENDIF
385
386         IF ( ln_dyninc ) THEN   
387            CALL iom_get( inum, jpdom_autoglo, 'bckinu', u_bkginc, 1 )             
388            CALL iom_get( inum, jpdom_autoglo, 'bckinv', v_bkginc, 1 )             
389            ! Apply the masks
390            u_bkginc(:,:,:) = u_bkginc(:,:,:) * umask(:,:,:)
391            v_bkginc(:,:,:) = v_bkginc(:,:,:) * vmask(:,:,:)
392            ! Set missing increments to 0.0 rather than 1e+20
393            ! to allow for differences in masks
394            WHERE( ABS( u_bkginc(:,:,:) ) > 1.0e+10 ) u_bkginc(:,:,:) = 0.0
395            WHERE( ABS( v_bkginc(:,:,:) ) > 1.0e+10 ) v_bkginc(:,:,:) = 0.0
396         ENDIF
397       
398         IF ( ln_sshinc ) THEN
399            CALL iom_get( inum, jpdom_autoglo, 'bckineta', ssh_bkginc, 1 )
400            ! Apply the masks
401            ssh_bkginc(:,:) = ssh_bkginc(:,:) * tmask(:,:,1)
402            ! Set missing increments to 0.0 rather than 1e+20
403            ! to allow for differences in masks
404            WHERE( ABS( ssh_bkginc(:,:) ) > 1.0e+10 ) ssh_bkginc(:,:) = 0.0
405         ENDIF
406
407         IF ( ln_seaiceinc ) THEN
408            CALL iom_get( inum, jpdom_autoglo, 'bckinseaice', seaice_bkginc, 1 )
409            ! Apply the masks
410            seaice_bkginc(:,:) = seaice_bkginc(:,:) * tmask(:,:,1)
411            ! Set missing increments to 0.0 rather than 1e+20
412            ! to allow for differences in masks
413            WHERE( ABS( seaice_bkginc(:,:) ) > 1.0e+10 ) seaice_bkginc(:,:) = 0.0
414         ENDIF
415
416         CALL iom_close( inum )
417 
418      ENDIF
419
420      !-----------------------------------------------------------------------
421      ! Apply divergence damping filter
422      !-----------------------------------------------------------------------
423
424      IF ( ln_dyninc .AND. nn_divdmp > 0 ) THEN
425
426         CALL wrk_alloc(jpi,jpj,hdiv) 
427
428         DO  jt = 1, nn_divdmp
429
430            DO jk = 1, jpkm1
431
432               hdiv(:,:) = 0._wp
433
434               DO jj = 2, jpjm1
435                  DO ji = fs_2, fs_jpim1   ! vector opt.
436                     hdiv(ji,jj) =   &
437                        (  e2u(ji  ,jj  ) * fse3u(ji  ,jj  ,jk) * u_bkginc(ji  ,jj  ,jk)     &
438                         - e2u(ji-1,jj  ) * fse3u(ji-1,jj  ,jk) * u_bkginc(ji-1,jj  ,jk)     &
439                         + e1v(ji  ,jj  ) * fse3v(ji  ,jj  ,jk) * v_bkginc(ji  ,jj  ,jk)     &
440                         - e1v(ji  ,jj-1) * fse3v(ji  ,jj-1,jk) * v_bkginc(ji  ,jj-1,jk)  )  &
441                         / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) )
442                  END DO
443               END DO
444
445               CALL lbc_lnk( hdiv, 'T', 1. )   ! lateral boundary cond. (no sign change)
446
447               DO jj = 2, jpjm1
448                  DO ji = fs_2, fs_jpim1   ! vector opt.
449                     u_bkginc(ji,jj,jk) = u_bkginc(ji,jj,jk) + 0.2_wp * ( e1t(ji+1,jj)*e2t(ji+1,jj) * hdiv(ji+1,jj)   &
450                                                                        - e1t(ji  ,jj)*e2t(ji  ,jj) * hdiv(ji  ,jj) ) &
451                                                                      / e1u(ji,jj) * umask(ji,jj,jk) 
452                     v_bkginc(ji,jj,jk) = v_bkginc(ji,jj,jk) + 0.2_wp * ( e1t(ji,jj+1)*e2t(ji,jj+1) * hdiv(ji,jj+1)   &
453                                                                        - e1t(ji,jj  )*e2t(ji,jj  ) * hdiv(ji,jj  ) ) &
454                                                                      / e2v(ji,jj) * vmask(ji,jj,jk) 
455                  END DO
456               END DO
457
458            END DO
459
460         END DO
461
462         CALL wrk_dealloc(jpi,jpj,hdiv) 
463
464      ENDIF
465
466
467
468      !-----------------------------------------------------------------------
469      ! Allocate and initialize the background state arrays
470      !-----------------------------------------------------------------------
471
472      IF ( ln_asmdin ) THEN
473
474         ALLOCATE( t_bkg(jpi,jpj,jpk) )
475         ALLOCATE( s_bkg(jpi,jpj,jpk) )
476         ALLOCATE( u_bkg(jpi,jpj,jpk) )
477         ALLOCATE( v_bkg(jpi,jpj,jpk) )
478         ALLOCATE( ssh_bkg(jpi,jpj)   )
479
480         t_bkg(:,:,:) = 0.0
481         s_bkg(:,:,:) = 0.0
482         u_bkg(:,:,:) = 0.0
483         v_bkg(:,:,:) = 0.0
484         ssh_bkg(:,:) = 0.0
485
486         !--------------------------------------------------------------------
487         ! Read from file the background state at analysis time
488         !--------------------------------------------------------------------
489
490         CALL iom_open( c_asmdin, inum )
491
492         CALL iom_get( inum, 'rdastp', zdate_bkg ) 
493       
494         IF(lwp) THEN
495            WRITE(numout,*) 
496            WRITE(numout,*) 'asm_inc_init : Assimilation background state valid at : ', &
497               &  NINT( zdate_bkg )
498            WRITE(numout,*) '~~~~~~~~~~~~'
499         ENDIF
500
501         IF ( NINT( zdate_bkg ) /= iitdin_date ) &
502            & CALL ctl_warn( ' Validity time of assimilation background state does', &
503            &                ' not agree with Direct Initialization time' )
504
505         IF ( ln_trainc ) THEN   
506            CALL iom_get( inum, jpdom_autoglo, 'tn', t_bkg )
507            CALL iom_get( inum, jpdom_autoglo, 'sn', s_bkg )
508            t_bkg(:,:,:) = t_bkg(:,:,:) * tmask(:,:,:)
509            s_bkg(:,:,:) = s_bkg(:,:,:) * tmask(:,:,:)
510         ENDIF
511
512         IF ( ln_dyninc ) THEN   
513            CALL iom_get( inum, jpdom_autoglo, 'un', u_bkg )
514            CALL iom_get( inum, jpdom_autoglo, 'vn', v_bkg )
515            u_bkg(:,:,:) = u_bkg(:,:,:) * umask(:,:,:)
516            v_bkg(:,:,:) = v_bkg(:,:,:) * vmask(:,:,:)
517         ENDIF
518       
519         IF ( ln_sshinc ) THEN
520            CALL iom_get( inum, jpdom_autoglo, 'sshn', ssh_bkg )
521            ssh_bkg(:,:) = ssh_bkg(:,:) * tmask(:,:,1)
522         ENDIF
523
524         CALL iom_close( inum )
525
526      ENDIF
527      !
528   END SUBROUTINE asm_inc_init
529
530
531   SUBROUTINE calc_date( kit000, kt, kdate0, kdate )
532      !!----------------------------------------------------------------------
533      !!                    ***  ROUTINE calc_date  ***
534      !!         
535      !! ** Purpose : Compute the calendar date YYYYMMDD at a given time step.
536      !!
537      !! ** Method  : Compute the calendar date YYYYMMDD at a given time step.
538      !!
539      !! ** Action  :
540      !!----------------------------------------------------------------------
541      INTEGER, INTENT(IN) :: kit000  ! Initial time step
542      INTEGER, INTENT(IN) :: kt      ! Current time step referenced to kit000
543      INTEGER, INTENT(IN) :: kdate0  ! Initial date
544      INTEGER, INTENT(OUT) :: kdate  ! Current date reference to kdate0
545      !
546      INTEGER :: iyea0    ! Initial year
547      INTEGER :: imon0    ! Initial month
548      INTEGER :: iday0    ! Initial day
549      INTEGER :: iyea     ! Current year
550      INTEGER :: imon     ! Current month
551      INTEGER :: iday     ! Current day
552      INTEGER :: idaystp  ! Number of days between initial and current date
553      INTEGER :: idaycnt  ! Day counter
554
555      INTEGER, DIMENSION(12) ::   imonth_len    !: length in days of the months of the current year
556
557      !-----------------------------------------------------------------------
558      ! Compute the calendar date YYYYMMDD
559      !-----------------------------------------------------------------------
560
561      ! Initial date
562      iyea0 =   kdate0 / 10000
563      imon0 = ( kdate0 - ( iyea0 * 10000 ) ) / 100
564      iday0 =   kdate0 - ( iyea0 * 10000 ) - ( imon0 * 100 ) 
565
566      ! Check that kt >= kit000 - 1
567      IF ( kt < kit000 - 1 ) CALL ctl_stop( ' kt must be >= kit000 - 1')
568
569      ! If kt = kit000 - 1 then set the date to the restart date
570      IF ( kt == kit000 - 1 ) THEN
571
572         kdate = ndastp
573         RETURN
574
575      ENDIF
576
577      ! Compute the number of days from the initial date
578      idaystp = INT( REAL( kt - kit000 ) * rdt / 86400. )
579   
580      iday    = iday0
581      imon    = imon0
582      iyea    = iyea0
583      idaycnt = 0
584
585      CALL calc_month_len( iyea, imonth_len )
586
587      DO WHILE ( idaycnt < idaystp )
588         iday = iday + 1
589         IF ( iday > imonth_len(imon) )  THEN
590            iday = 1
591            imon = imon + 1
592         ENDIF
593         IF ( imon > 12 ) THEN
594            imon = 1
595            iyea = iyea + 1
596            CALL calc_month_len( iyea, imonth_len )  ! update month lengths
597         ENDIF                 
598         idaycnt = idaycnt + 1
599      END DO
600      !
601      kdate = iyea * 10000 + imon * 100 + iday
602      !
603   END SUBROUTINE
604
605
606   SUBROUTINE calc_month_len( iyear, imonth_len )
607      !!----------------------------------------------------------------------
608      !!                    ***  ROUTINE calc_month_len  ***
609      !!         
610      !! ** Purpose : Compute the number of days in a months given a year.
611      !!
612      !! ** Method  :
613      !!----------------------------------------------------------------------
614      INTEGER, DIMENSION(12) ::   imonth_len    !: length in days of the months of the current year
615      INTEGER :: iyear         !: year
616      !!----------------------------------------------------------------------
617      !
618      ! length of the month of the current year (from nleapy, read in namelist)
619      IF ( nleapy < 2 ) THEN
620         imonth_len(:) = (/ 31, 28, 31, 30, 31, 30, 31, 31, 30, 31, 30, 31 /)
621         IF ( nleapy == 1 ) THEN   ! we are using calendar with leap years
622            IF ( MOD(iyear, 4) == 0 .AND. ( MOD(iyear, 400) == 0 .OR. MOD(iyear, 100) /= 0 ) ) THEN
623               imonth_len(2) = 29
624            ENDIF
625         ENDIF
626      ELSE
627         imonth_len(:) = nleapy   ! all months with nleapy days per year
628      ENDIF
629      !
630   END SUBROUTINE
631
632
633   SUBROUTINE tra_asm_inc( kt )
634      !!----------------------------------------------------------------------
635      !!                    ***  ROUTINE tra_asm_inc  ***
636      !!         
637      !! ** Purpose : Apply the tracer (T and S) assimilation increments
638      !!
639      !! ** Method  : Direct initialization or Incremental Analysis Updating
640      !!
641      !! ** Action  :
642      !!----------------------------------------------------------------------
643      INTEGER, INTENT(IN) :: kt               ! Current time step
644      !
645      INTEGER :: ji,jj,jk
646      INTEGER :: it
647      REAL(wp) :: zincwgt  ! IAU weight for current time step
648      REAL (wp), DIMENSION(jpi,jpj,jpk) :: fzptnz ! 3d freezing point values
649      !!----------------------------------------------------------------------
650
651      ! freezing point calculation taken from oc_fz_pt (but calculated for all depths)
652      ! used to prevent the applied increments taking the temperature below the local freezing point
653
654      DO jk = 1, jpkm1
655        CALL eos_fzp( tsn(:,:,jk,jp_sal), fzptnz(:,:,jk), fsdept(:,:,jk) )
656      END DO
657
658      IF ( ln_asmiau ) THEN
659
660         !--------------------------------------------------------------------
661         ! Incremental Analysis Updating
662         !--------------------------------------------------------------------
663
664         IF ( ( kt >= nitiaustr_r ).AND.( kt <= nitiaufin_r ) ) THEN
665
666            it = kt - nit000 + 1
667            zincwgt = wgtiau(it) / rdt   ! IAU weight for the current time step
668
669            IF(lwp) THEN
670               WRITE(numout,*) 
671               WRITE(numout,*) 'tra_asm_inc : Tracer IAU at time step = ', kt,' with IAU weight = ', wgtiau(it)
672               WRITE(numout,*) '~~~~~~~~~~~~'
673            ENDIF
674
675            ! Update the tracer tendencies
676            DO jk = 1, jpkm1
677               IF (ln_temnofreeze) THEN
678                  ! Do not apply negative increments if the temperature will fall below freezing
679                  WHERE(t_bkginc(:,:,jk) > 0.0_wp .OR. &
680                     &   tsn(:,:,jk,jp_tem) + tsa(:,:,jk,jp_tem) + t_bkginc(:,:,jk) * wgtiau(it) > fzptnz(:,:,jk) ) 
681                     tsa(:,:,jk,jp_tem) = tsa(:,:,jk,jp_tem) + t_bkginc(:,:,jk) * zincwgt 
682                  END WHERE
683               ELSE
684                  tsa(:,:,jk,jp_tem) = tsa(:,:,jk,jp_tem) + t_bkginc(:,:,jk) * zincwgt 
685               ENDIF
686               IF (ln_salfix) THEN
687                  ! Do not apply negative increments if the salinity will fall below a specified
688                  ! minimum value salfixmin
689                  WHERE(s_bkginc(:,:,jk) > 0.0_wp .OR. &
690                     &   tsn(:,:,jk,jp_sal) + tsa(:,:,jk,jp_sal) + s_bkginc(:,:,jk) * wgtiau(it) > salfixmin ) 
691                     tsa(:,:,jk,jp_sal) = tsa(:,:,jk,jp_sal) + s_bkginc(:,:,jk) * zincwgt
692                  END WHERE
693               ELSE
694                  tsa(:,:,jk,jp_sal) = tsa(:,:,jk,jp_sal) + s_bkginc(:,:,jk) * zincwgt
695               ENDIF
696            END DO
697
698         ENDIF
699
700         IF ( kt == nitiaufin_r + 1  ) THEN   ! For bias crcn to work
701            DEALLOCATE( t_bkginc )
702            DEALLOCATE( s_bkginc )
703         ENDIF
704
705
706      ELSEIF ( ln_asmdin ) THEN
707
708         !--------------------------------------------------------------------
709         ! Direct Initialization
710         !--------------------------------------------------------------------
711           
712         IF ( kt == nitdin_r ) THEN
713
714            neuler = 0  ! Force Euler forward step
715
716            ! Initialize the now fields with the background + increment
717            IF (ln_temnofreeze) THEN
718               ! Do not apply negative increments if the temperature will fall below freezing
719               WHERE( t_bkginc(:,:,:) > 0.0_wp .OR. tsn(:,:,:,jp_tem) + t_bkginc(:,:,:) > fzptnz(:,:,:) ) 
720                  tsn(:,:,:,jp_tem) = t_bkg(:,:,:) + t_bkginc(:,:,:)   
721               END WHERE
722            ELSE
723               tsn(:,:,:,jp_tem) = t_bkg(:,:,:) + t_bkginc(:,:,:)   
724            ENDIF
725            IF (ln_salfix) THEN
726               ! Do not apply negative increments if the salinity will fall below a specified
727               ! minimum value salfixmin
728               WHERE( s_bkginc(:,:,:) > 0.0_wp .OR. tsn(:,:,:,jp_sal) + s_bkginc(:,:,:) > salfixmin ) 
729                  tsn(:,:,:,jp_sal) = s_bkg(:,:,:) + s_bkginc(:,:,:)   
730               END WHERE
731            ELSE
732               tsn(:,:,:,jp_sal) = s_bkg(:,:,:) + s_bkginc(:,:,:)   
733            ENDIF
734
735            tsb(:,:,:,:) = tsn(:,:,:,:)                 ! Update before fields
736
737            CALL eos( tsb, rhd, rhop, gdept_0(:,:,:) )  ! Before potential and in situ densities
738!!gm  fabien
739!            CALL eos( tsb, rhd, rhop )                ! Before potential and in situ densities
740!!gm
741
742
743            IF( ln_zps .AND. .NOT. lk_c1d .AND. .NOT. ln_isfcav)      &
744               &  CALL zps_hde    ( kt, jpts, tsb, gtsu, gtsv,        &  ! Partial steps: before horizontal gradient
745               &                              rhd, gru , grv          )  ! of t, s, rd at the last ocean level
746            IF( ln_zps .AND. .NOT. lk_c1d .AND.       ln_isfcav)      &
747               &  CALL zps_hde_isf( nit000, jpts, tsb, gtsu, gtsv,    &    ! Partial steps for top cell (ISF)
748               &                                  rhd, gru , grv , aru , arv , gzu , gzv , ge3ru , ge3rv ,   &
749               &                           gtui, gtvi, grui, grvi, arui, arvi, gzui, gzvi, ge3rui, ge3rvi    ) ! of t, s, rd at the last ocean level
750
751#if defined key_zdfkpp
752            CALL eos( tsn, rhd, fsdept_n(:,:,:) )                      ! Compute rhd
753!!gm fabien            CALL eos( tsn, rhd )                      ! Compute rhd
754#endif
755
756            DEALLOCATE( t_bkginc )
757            DEALLOCATE( s_bkginc )
758            DEALLOCATE( t_bkg    )
759            DEALLOCATE( s_bkg    )
760         ENDIF
761         
762      ENDIF
763      ! Perhaps the following call should be in step
764      IF   ( ln_seaiceinc  )   CALL seaice_asm_inc ( kt )   ! apply sea ice concentration increment
765      !
766   END SUBROUTINE tra_asm_inc
767
768
769   SUBROUTINE dyn_asm_inc( kt )
770      !!----------------------------------------------------------------------
771      !!                    ***  ROUTINE dyn_asm_inc  ***
772      !!         
773      !! ** Purpose : Apply the dynamics (u and v) assimilation increments.
774      !!
775      !! ** Method  : Direct initialization or Incremental Analysis Updating.
776      !!
777      !! ** Action  :
778      !!----------------------------------------------------------------------
779      INTEGER, INTENT(IN) :: kt   ! Current time step
780      !
781      INTEGER :: jk
782      INTEGER :: it
783      REAL(wp) :: zincwgt  ! IAU weight for current time step
784      !!----------------------------------------------------------------------
785
786      IF ( ln_asmiau ) THEN
787
788         !--------------------------------------------------------------------
789         ! Incremental Analysis Updating
790         !--------------------------------------------------------------------
791
792         IF ( ( kt >= nitiaustr_r ).AND.( kt <= nitiaufin_r ) ) THEN
793
794            it = kt - nit000 + 1
795            zincwgt = wgtiau(it) / rdt   ! IAU weight for the current time step
796
797            IF(lwp) THEN
798               WRITE(numout,*) 
799               WRITE(numout,*) 'dyn_asm_inc : Dynamics IAU at time step = ', &
800                  &  kt,' with IAU weight = ', wgtiau(it)
801               WRITE(numout,*) '~~~~~~~~~~~~'
802            ENDIF
803
804            ! Update the dynamic tendencies
805            DO jk = 1, jpkm1
806               ua(:,:,jk) = ua(:,:,jk) + u_bkginc(:,:,jk) * zincwgt
807               va(:,:,jk) = va(:,:,jk) + v_bkginc(:,:,jk) * zincwgt
808            END DO
809           
810            IF ( kt == nitiaufin_r ) THEN
811               DEALLOCATE( u_bkginc )
812               DEALLOCATE( v_bkginc )
813            ENDIF
814
815         ENDIF
816
817      ELSEIF ( ln_asmdin ) THEN 
818
819         !--------------------------------------------------------------------
820         ! Direct Initialization
821         !--------------------------------------------------------------------
822         
823         IF ( kt == nitdin_r ) THEN
824
825            neuler = 0                    ! Force Euler forward step
826
827            ! Initialize the now fields with the background + increment
828            un(:,:,:) = u_bkg(:,:,:) + u_bkginc(:,:,:)
829            vn(:,:,:) = v_bkg(:,:,:) + v_bkginc(:,:,:) 
830
831            ub(:,:,:) = un(:,:,:)         ! Update before fields
832            vb(:,:,:) = vn(:,:,:)
833 
834            DEALLOCATE( u_bkg    )
835            DEALLOCATE( v_bkg    )
836            DEALLOCATE( u_bkginc )
837            DEALLOCATE( v_bkginc )
838         ENDIF
839         !
840      ENDIF
841      !
842   END SUBROUTINE dyn_asm_inc
843
844
845   SUBROUTINE ssh_asm_inc( kt )
846      !!----------------------------------------------------------------------
847      !!                    ***  ROUTINE ssh_asm_inc  ***
848      !!         
849      !! ** Purpose : Apply the sea surface height assimilation increment.
850      !!
851      !! ** Method  : Direct initialization or Incremental Analysis Updating.
852      !!
853      !! ** Action  :
854      !!----------------------------------------------------------------------
855      INTEGER, INTENT(IN) :: kt   ! Current time step
856      !
857      INTEGER :: it
858      INTEGER :: jk
859      REAL(wp) :: zincwgt  ! IAU weight for current time step
860      !!----------------------------------------------------------------------
861
862      IF ( ln_asmiau ) THEN
863
864         !--------------------------------------------------------------------
865         ! Incremental Analysis Updating
866         !--------------------------------------------------------------------
867
868         IF ( ( kt >= nitiaustr_r ).AND.( kt <= nitiaufin_r ) ) THEN
869
870            it = kt - nit000 + 1
871            zincwgt = wgtiau(it) / rdt   ! IAU weight for the current time step
872
873            IF(lwp) THEN
874               WRITE(numout,*) 
875               WRITE(numout,*) 'ssh_asm_inc : SSH IAU at time step = ', &
876                  &  kt,' with IAU weight = ', wgtiau(it)
877               WRITE(numout,*) '~~~~~~~~~~~~'
878            ENDIF
879
880            ! Save the tendency associated with the IAU weighted SSH increment
881            ! (applied in dynspg.*)
882#if defined key_asminc
883            ssh_iau(:,:) = ssh_bkginc(:,:) * zincwgt
884#endif
885            IF ( kt == nitiaufin_r ) THEN
886               DEALLOCATE( ssh_bkginc )
887            ENDIF
888
889         ENDIF
890
891      ELSEIF ( ln_asmdin ) THEN
892
893         !--------------------------------------------------------------------
894         ! Direct Initialization
895         !--------------------------------------------------------------------
896
897         IF ( kt == nitdin_r ) THEN
898
899            neuler = 0                    ! Force Euler forward step
900
901            ! Initialize the now fields the background + increment
902            sshn(:,:) = ssh_bkg(:,:) + ssh_bkginc(:,:) 
903
904            ! Update before fields
905            sshb(:,:) = sshn(:,:)         
906
907            IF( lk_vvl ) THEN
908               DO jk = 1, jpk
909                  fse3t_b(:,:,jk) = fse3t_n(:,:,jk)
910               END DO
911            ENDIF
912
913            DEALLOCATE( ssh_bkg    )
914            DEALLOCATE( ssh_bkginc )
915
916         ENDIF
917         !
918      ENDIF
919      !
920   END SUBROUTINE ssh_asm_inc
921
922
923   SUBROUTINE seaice_asm_inc( kt, kindic )
924      !!----------------------------------------------------------------------
925      !!                    ***  ROUTINE seaice_asm_inc  ***
926      !!         
927      !! ** Purpose : Apply the sea ice assimilation increment.
928      !!
929      !! ** Method  : Direct initialization or Incremental Analysis Updating.
930      !!
931      !! ** Action  :
932      !!
933      !!----------------------------------------------------------------------
934      IMPLICIT NONE
935      !
936      INTEGER, INTENT(in)           ::   kt   ! Current time step
937      INTEGER, INTENT(in), OPTIONAL ::   kindic   ! flag for disabling the deallocation
938      !
939      INTEGER  ::   it
940      REAL(wp) ::   zincwgt   ! IAU weight for current time step
941#if defined key_lim2
942      REAL(wp), DIMENSION(jpi,jpj) ::   zofrld, zohicif, zseaicendg, zhicifinc  ! LIM
943      REAL(wp) ::   zhicifmin = 0.5_wp      ! ice minimum depth in metres
944#endif
945      !!----------------------------------------------------------------------
946
947      IF ( ln_asmiau ) THEN
948
949         !--------------------------------------------------------------------
950         ! Incremental Analysis Updating
951         !--------------------------------------------------------------------
952
953         IF ( ( kt >= nitiaustr_r ).AND.( kt <= nitiaufin_r ) ) THEN
954
955            it = kt - nit000 + 1
956            zincwgt = wgtiau(it)      ! IAU weight for the current time step
957            ! note this is not a tendency so should not be divided by rdt (as with the tracer and other increments)
958
959            IF(lwp) THEN
960               WRITE(numout,*) 
961               WRITE(numout,*) 'seaice_asm_inc : sea ice conc IAU at time step = ', &
962                  &  kt,' with IAU weight = ', wgtiau(it)
963               WRITE(numout,*) '~~~~~~~~~~~~'
964            ENDIF
965
966            ! Sea-ice : LIM-3 case (to add)
967
968#if defined key_lim2
969            ! Sea-ice : LIM-2 case
970            zofrld (:,:) = frld(:,:)
971            zohicif(:,:) = hicif(:,:)
972            !
973            frld  = MIN( MAX( frld (:,:) - seaice_bkginc(:,:) * zincwgt, 0.0_wp), 1.0_wp)
974            pfrld = MIN( MAX( pfrld(:,:) - seaice_bkginc(:,:) * zincwgt, 0.0_wp), 1.0_wp)
975            fr_i(:,:) = 1.0_wp - frld(:,:)        ! adjust ice fraction
976            !
977            zseaicendg(:,:) = zofrld(:,:) - frld(:,:)   ! find out actual sea ice nudge applied
978            !
979            ! Nudge sea ice depth to bring it up to a required minimum depth
980            WHERE( zseaicendg(:,:) > 0.0_wp .AND. hicif(:,:) < zhicifmin ) 
981               zhicifinc(:,:) = (zhicifmin - hicif(:,:)) * zincwgt   
982            ELSEWHERE
983               zhicifinc(:,:) = 0.0_wp
984            END WHERE
985            !
986            ! nudge ice depth
987            hicif (:,:) = hicif (:,:) + zhicifinc(:,:)
988            phicif(:,:) = phicif(:,:) + zhicifinc(:,:)       
989            !
990            ! seaice salinity balancing (to add)
991#endif
992
993#if defined key_cice && defined key_asminc
994            ! Sea-ice : CICE case. Pass ice increment tendency into CICE
995            ndaice_da(:,:) = seaice_bkginc(:,:) * zincwgt / rdt
996#endif
997
998            IF ( kt == nitiaufin_r ) THEN
999               DEALLOCATE( seaice_bkginc )
1000            ENDIF
1001
1002         ELSE
1003
1004#if defined key_cice && defined key_asminc
1005            ! Sea-ice : CICE case. Zero ice increment tendency into CICE
1006            ndaice_da(:,:) = 0.0_wp
1007#endif
1008
1009         ENDIF
1010
1011      ELSEIF ( ln_asmdin ) THEN
1012
1013         !--------------------------------------------------------------------
1014         ! Direct Initialization
1015         !--------------------------------------------------------------------
1016
1017         IF ( kt == nitdin_r ) THEN
1018
1019            neuler = 0                    ! Force Euler forward step
1020
1021            ! Sea-ice : LIM-3 case (to add)
1022
1023#if defined key_lim2
1024            ! Sea-ice : LIM-2 case.
1025            zofrld(:,:)=frld(:,:)
1026            zohicif(:,:)=hicif(:,:)
1027            !
1028            ! Initialize the now fields the background + increment
1029            frld (:,:) = MIN( MAX( frld(:,:) - seaice_bkginc(:,:), 0.0_wp), 1.0_wp)
1030            pfrld(:,:) = frld(:,:) 
1031            fr_i (:,:) = 1.0_wp - frld(:,:)                ! adjust ice fraction
1032            zseaicendg(:,:) = zofrld(:,:) - frld(:,:)      ! find out actual sea ice nudge applied
1033            !
1034            ! Nudge sea ice depth to bring it up to a required minimum depth
1035            WHERE( zseaicendg(:,:) > 0.0_wp .AND. hicif(:,:) < zhicifmin ) 
1036               zhicifinc(:,:) = (zhicifmin - hicif(:,:)) * zincwgt   
1037            ELSEWHERE
1038               zhicifinc(:,:) = 0.0_wp
1039            END WHERE
1040            !
1041            ! nudge ice depth
1042            hicif (:,:) = hicif (:,:) + zhicifinc(:,:)
1043            phicif(:,:) = phicif(:,:)       
1044            !
1045            ! seaice salinity balancing (to add)
1046#endif
1047 
1048#if defined key_cice && defined key_asminc
1049            ! Sea-ice : CICE case. Pass ice increment tendency into CICE
1050           ndaice_da(:,:) = seaice_bkginc(:,:) / rdt
1051#endif
1052           IF ( .NOT. PRESENT(kindic) ) THEN
1053              DEALLOCATE( seaice_bkginc )
1054           END IF
1055
1056         ELSE
1057
1058#if defined key_cice && defined key_asminc
1059            ! Sea-ice : CICE case. Zero ice increment tendency into CICE
1060            ndaice_da(:,:) = 0.0_wp
1061#endif
1062         
1063         ENDIF
1064
1065!#if defined defined key_lim2 || defined key_cice
1066!
1067!            IF (ln_seaicebal ) THEN       
1068!             !! balancing salinity increments
1069!             !! simple case from limflx.F90 (doesn't include a mass flux)
1070!             !! assumption is that as ice concentration is reduced or increased
1071!             !! the snow and ice depths remain constant
1072!             !! note that snow is being created where ice concentration is being increased
1073!             !! - could be more sophisticated and
1074!             !! not do this (but would need to alter h_snow)
1075!
1076!             usave(:,:,:)=sb(:,:,:)   ! use array as a temporary store
1077!
1078!             DO jj = 1, jpj
1079!               DO ji = 1, jpi
1080!           ! calculate change in ice and snow mass per unit area
1081!           ! positive values imply adding salt to the ocean (results from ice formation)
1082!           ! fwf : ice formation and melting
1083!
1084!                 zfons = ( -nfresh_da(ji,jj)*soce + nfsalt_da(ji,jj) )*rdt
1085!
1086!           ! change salinity down to mixed layer depth
1087!                 mld=hmld_kara(ji,jj)
1088!
1089!           ! prevent small mld
1090!           ! less than 10m can cause salinity instability
1091!                 IF (mld < 10) mld=10
1092!
1093!           ! set to bottom of a level
1094!                 DO jk = jpk-1, 2, -1
1095!                   IF ((mld > gdepw(ji,jj,jk)) .and. (mld < gdepw(ji,jj,jk+1))) THEN
1096!                     mld=gdepw(ji,jj,jk+1)
1097!                     jkmax=jk
1098!                   ENDIF
1099!                 ENDDO
1100!
1101!            ! avoid applying salinity balancing in shallow water or on land
1102!            !
1103!
1104!            ! dsal_ocn (psu kg m^-2) / (kg m^-3 * m)
1105!
1106!                 dsal_ocn=0.0_wp
1107!                 sal_thresh=5.0_wp        ! minimum salinity threshold for salinity balancing
1108!
1109!                 if (tmask(ji,jj,1) > 0 .AND. tmask(ji,jj,jkmax) > 0 ) &
1110!                              dsal_ocn = zfons / (rhop(ji,jj,1) * mld)
1111!
1112!           ! put increments in for levels in the mixed layer
1113!           ! but prevent salinity below a threshold value
1114!
1115!                   DO jk = 1, jkmax             
1116!
1117!                     IF (dsal_ocn > 0.0_wp .or. sb(ji,jj,jk)+dsal_ocn > sal_thresh) THEN
1118!                           sb(ji,jj,jk) = sb(ji,jj,jk) + dsal_ocn
1119!                           sn(ji,jj,jk) = sn(ji,jj,jk) + dsal_ocn
1120!                     ENDIF
1121!
1122!                   ENDDO
1123!
1124!      !            !  salt exchanges at the ice/ocean interface
1125!      !            zpmess         = zfons / rdt_ice    ! rdt_ice is ice timestep
1126!      !
1127!      !! Adjust fsalt. A +ve fsalt means adding salt to ocean
1128!      !!           fsalt(ji,jj) =  fsalt(ji,jj) + zpmess     ! adjust fsalt 
1129!      !!               
1130!      !!           emps(ji,jj) = emps(ji,jj) + zpmess        ! or adjust emps (see icestp1d)
1131!      !!                                                     ! E-P (kg m-2 s-2)
1132!      !            emp(ji,jj) = emp(ji,jj) + zpmess          ! E-P (kg m-2 s-2)
1133!               ENDDO !ji
1134!             ENDDO !jj!
1135!
1136!            ENDIF !ln_seaicebal
1137!
1138!#endif
1139
1140      ENDIF
1141
1142   END SUBROUTINE seaice_asm_inc
1143   
1144   !!======================================================================
1145END MODULE asminc
Note: See TracBrowser for help on using the repository browser.