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/OCE/ASM – NEMO

source: NEMO/branches/2020/dev_r13296_HPC-07_mocavero_mpi3/src/OCE/ASM/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

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