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

source: branches/UKMO/dev_r5518_obs_oper_strthr/NEMOGCM/NEMO/OPA_SRC/ASM/asminc.F90 @ 15258

Last change on this file since 15258 was 15256, checked in by jroberts, 3 years ago

Changes manually merged from https://forge.ipsl.jussieu.fr/nemo/browser/branches/UKMO/dev_r5518_GO6_package_STARTHOUR?

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