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/2019/dev_r12072_MERGE_OPTION2_2019/src/OCE/ASM – NEMO

source: NEMO/branches/2019/dev_r12072_MERGE_OPTION2_2019/src/OCE/ASM/asminc.F90 @ 12202

Last change on this file since 12202 was 12202, checked in by cetlod, 4 years ago

dev_merge_option2 : merge in dev_r11613_ENHANCE-04_namelists_as_internalfiles

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