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/2015/dev_CMCC_merge_2015/NEMOGCM/NEMO/OPA_SRC/ASM – NEMO

source: branches/2015/dev_CMCC_merge_2015/NEMOGCM/NEMO/OPA_SRC/ASM/asminc.F90 @ 6051

Last change on this file since 6051 was 6051, checked in by lovato, 8 years ago

Merge branches/2015/dev_r5056_CMCC4_simplification (see ticket #1456)

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