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

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

source: branches/UKMO/dev_r5518_GO6_starthour_obsoper/NEMOGCM/NEMO/OPA_SRC/ASM/asminc.F90 @ 12555

Last change on this file since 12555 was 12555, checked in by charris, 4 years ago

Changes from GO6 package branch (GMED ticket 450):

svn merge -r 11035:11101 svn+ssh://charris@forge.ipsl.jussieu.fr/ipsl/forge/projets/nemo/svn/branches/UKMO/dev_r5518_GO6_package

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