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

source: branches/2015/dev_merge_2015/NEMOGCM/NEMO/OPA_SRC/ASM/asminc.F90 @ 6060

Last change on this file since 6060 was 6060, checked in by timgraham, 8 years ago

Merged dev_r5836_noc2_VVL_BY_DEFAULT into branch

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