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/2020/dev_r13296_HPC-07_mocavero_mpi3/src/SWE – NEMO

source: NEMO/branches/2020/dev_r13296_HPC-07_mocavero_mpi3/src/SWE/asminc.F90 @ 13630

Last change on this file since 13630 was 13630, checked in by mocavero, 4 years ago

Add neighborhood collectives calls in the NEMO src - ticket #2496

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