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 NEMO/branches/UKMO/NEMO_4.0_mirror_text_diagnostics/src/OCE/ASM – NEMO

source: NEMO/branches/UKMO/NEMO_4.0_mirror_text_diagnostics/src/OCE/ASM/asminc.F90 @ 10986

Last change on this file since 10986 was 10986, checked in by andmirek, 5 years ago

GMED 462 add flush

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