New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
asminc.F90 in branches/2017/dev_r7881_ENHANCE09_RK3/NEMOGCM/NEMO/OPA_SRC/ASM – NEMO

source: branches/2017/dev_r7881_ENHANCE09_RK3/NEMOGCM/NEMO/OPA_SRC/ASM/asminc.F90 @ 8637

Last change on this file since 8637 was 8586, checked in by gm, 7 years ago

#1911 (ENHANCE-09): PART I.3 - phasing with branch dev_r8183_ICEMODEL revision 8575

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