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

source: branches/2014/dev_r4650_UKMO7_STARTHOUR/NEMOGCM/NEMO/OPA_SRC/ASM/asminc.F90 @ 5872

Last change on this file since 5872 was 5872, checked in by djlea, 8 years ago

Fixes to enable restartability for runs of less than one day.

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