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

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

source: NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/ASM/asminc.F90 @ 11960

Last change on this file since 11960 was 11960, checked in by acc, 4 years ago

Branch 2019/dev_r11943_MERGE_2019. Merge in changes from 2019/dev_r11613_ENHANCE-04_namelists_as_internalfiles. (svn merge -r 11614:11954). Resolved tree conflicts and one actual conflict. Sette tested(these changes alter the ext/AGRIF reference; remember to update). See ticket #2341

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