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

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
asminc.F90 in NEMO/branches/UKMO/NEMO_4.0.4_CO9_OBS_ASM/src/OCE/ASM – NEMO

source: NEMO/branches/UKMO/NEMO_4.0.4_CO9_OBS_ASM/src/OCE/ASM/asminc.F90 @ 15183

Last change on this file since 15183 was 15183, checked in by kingr, 3 years ago

Adding changes required to write out time-averaged assimilation background. Required moving call to asm_bkg_wri from asminc.F90 to nemogcm.F90 to avoid circular USE statements.

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