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

source: branches/2017/dev_merge_2017/NEMOGCM/NEMO/OPA_SRC/ASM/asminc.F90 @ 9168

Last change on this file since 9168 was 9168, checked in by gm, 6 years ago

dev_merge_2017: OPA_SRC & CONFIG: remove useless warning when reading namelist_cfg

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