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

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

source: NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/ASM/asminc.F90 @ 10954

Last change on this file since 10954 was 10954, checked in by acc, 5 years ago

2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps : Convert TRA modules and all knock on effects of these conversions. SETTE tested

  • Property svn:keywords set to Id
File size: 47.9 KB
RevLine 
[2128]1MODULE asminc
2   !!======================================================================
3   !!                       ***  MODULE asminc  ***
4   !! Assimilation increment : Apply an increment generated by data
5   !!                          assimilation
6   !!======================================================================
[2392]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
[3764]12   !!            3.4  ! 2012-10  (A. Weaver and K. Mogensen) Fix for direct initialization
[6140]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
[2392]15   !!----------------------------------------------------------------------
[2128]16
17   !!----------------------------------------------------------------------
[3785]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
[9023]22   !!   ssh_asm_div    : Apply divergence associated with SSH increment
[3785]23   !!   seaice_asm_inc : Apply the seaice increment
[2128]24   !!----------------------------------------------------------------------
[9019]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
[9254]33   USE asmbkg          !
[9019]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
[9570]37#if defined key_si3
[9019]38   USE ice      , ONLY : hm_i, at_i, at_i_b
[3764]39#endif
[9019]40   !
41   USE in_out_manager  ! I/O manager
42   USE iom             ! Library to read input files
43   USE lib_mpp         ! MPP library
[2128]44
45   IMPLICIT NONE
46   PRIVATE
[2392]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
[9023]52   PUBLIC   ssh_asm_div    !: Apply the SSH divergence
[3764]53   PUBLIC   seaice_asm_inc !: Apply the seaice increment
[2128]54
55#if defined key_asminc
[2392]56    LOGICAL, PUBLIC, PARAMETER :: lk_asminc = .TRUE.   !: Logical switch for assimilation increment interface
[2128]57#else
[2392]58    LOGICAL, PUBLIC, PARAMETER :: lk_asminc = .FALSE.  !: No assimilation increments
[2128]59#endif
[5836]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
[3764]68   LOGICAL, PUBLIC :: ln_temnofreeze = .FALSE. !: Don't allow the temperature to drop below freezing
[5836]69   INTEGER, PUBLIC :: nn_divdmp     !: Apply divergence damping filter nn_divdmp times
[2128]70
[2392]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
[2128]76#if defined key_asminc
[9213]77   REAL(wp), PUBLIC, DIMENSION(:,:)  , ALLOCATABLE ::   ssh_iau              !: IAU-weighted sea surface height increment
[2128]78#endif
[2392]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)
[2128]88
[2392]89   REAL(wp), DIMENSION(:,:), ALLOCATABLE ::   ssh_bkg, ssh_bkginc   ! Background sea surface height and its increment
[3764]90   REAL(wp), DIMENSION(:,:), ALLOCATABLE ::   seaice_bkginc         ! Increment to the background sea ice conc
[8030]91#if defined key_cice && defined key_asminc
92   REAL(wp), DIMENSION(:,:), ALLOCATABLE ::   ndaice_da             ! ice increment tendency into CICE
93#endif
[2128]94
[3294]95   !! * Substitutions
96#  include "vectopt_loop_substitute.h90"
[2287]97   !!----------------------------------------------------------------------
[9598]98   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
[2287]99   !! $Id$
[10068]100   !! Software governed by the CeCILL license (see ./LICENSE)
[2287]101   !!----------------------------------------------------------------------
[2128]102CONTAINS
103
[10954]104   SUBROUTINE asm_inc_init( Kmm )
[2128]105      !!----------------------------------------------------------------------
106      !!                    ***  ROUTINE asm_inc_init  ***
107      !!         
108      !! ** Purpose : Initialize the assimilation increment and IAU weights.
109      !!
110      !! ** Method  : Initialize the assimilation increment and IAU weights.
111      !!
112      !! ** Action  :
113      !!----------------------------------------------------------------------
[10954]114      INTEGER, INTENT(in) ::   Kmm  ! time level index
[4990]115      INTEGER :: ji, jj, jk, jt  ! dummy loop indices
116      INTEGER :: imid, inum      ! local integers
117      INTEGER :: ios             ! Local integer output status for namelist read
[2128]118      INTEGER :: iiauper         ! Number of time steps in the IAU period
119      INTEGER :: icycper         ! Number of time steps in the cycle
[6140]120      REAL(KIND=dp) :: ditend_date     ! Date YYYYMMDD.HHMMSS of final time step
121      REAL(KIND=dp) :: ditbkg_date     ! Date YYYYMMDD.HHMMSS of background time step for Jb term
122      REAL(KIND=dp) :: ditdin_date     ! Date YYYYMMDD.HHMMSS of background time step for DI
123      REAL(KIND=dp) :: ditiaustr_date  ! Date YYYYMMDD.HHMMSS of IAU interval start time step
124      REAL(KIND=dp) :: ditiaufin_date  ! Date YYYYMMDD.HHMMSS of IAU interval final time step
125
[2128]126      REAL(wp) :: znorm        ! Normalization factor for IAU weights
[4990]127      REAL(wp) :: ztotwgt      ! Value of time-integrated IAU weights (should be equal to one)
[2128]128      REAL(wp) :: z_inc_dateb  ! Start date of interval on which increment is valid
129      REAL(wp) :: z_inc_datef  ! End date of interval on which increment is valid
130      REAL(wp) :: zdate_bkg    ! Date in background state file for DI
131      REAL(wp) :: zdate_inc    ! Time axis in increments file
[4990]132      !
[9019]133      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zhdiv   ! 2D workspace
[2392]134      !!
[3764]135      NAMELIST/nam_asminc/ ln_bkgwri,                                      &
[2392]136         &                 ln_trainc, ln_dyninc, ln_sshinc,                &
137         &                 ln_asmdin, ln_asmiau,                           &
138         &                 nitbkg, nitdin, nitiaustr, nitiaufin, niaufn,   &
[4990]139         &                 ln_salfix, salfixmin, nn_divdmp
[2392]140      !!----------------------------------------------------------------------
[2128]141
142      !-----------------------------------------------------------------------
143      ! Read Namelist nam_asminc : assimilation increment interface
144      !-----------------------------------------------------------------------
[5836]145      ln_seaiceinc   = .FALSE.
[3764]146      ln_temnofreeze = .FALSE.
[2128]147
[4147]148      REWIND( numnam_ref )              ! Namelist nam_asminc in reference namelist : Assimilation increment
149      READ  ( numnam_ref, nam_asminc, IOSTAT = ios, ERR = 901)
[9168]150901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'nam_asminc in reference namelist', lwp )
[4147]151      REWIND( numnam_cfg )              ! Namelist nam_asminc in configuration namelist : Assimilation increment
152      READ  ( numnam_cfg, nam_asminc, IOSTAT = ios, ERR = 902 )
[9168]153902   IF( ios >  0 )   CALL ctl_nam ( ios , 'nam_asminc in configuration namelist', lwp )
[4624]154      IF(lwm) WRITE ( numond, nam_asminc )
[4147]155
[2128]156      ! Control print
157      IF(lwp) THEN
158         WRITE(numout,*)
159         WRITE(numout,*) 'asm_inc_init : Assimilation increment initialization :'
160         WRITE(numout,*) '~~~~~~~~~~~~'
[2392]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
[3764]167         WRITE(numout,*) '      Logical switch for applying sea ice increments        ln_seaiceinc = ', ln_seaiceinc
[2392]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
[2128]176      ENDIF
177
[9019]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
[2128]182
[9019]183      iiauper     = nitiaufin_r - nitiaustr_r + 1     ! IAU interval length
184      icycper     = nitend      - nit000      + 1     ! Cycle interval length
[2128]185
[9019]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
[6140]191
[2128]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
[6140]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
[2128]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 ) ) &
[3764]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.', &
[2128]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
[3764]229      IF ( ( .NOT. ln_trainc ).AND.( .NOT. ln_dyninc ).AND.( .NOT. ln_sshinc ).AND.( .NOT. ln_seaiceinc ) &
[2128]230         &                     )  &
[3764]231         & CALL ctl_warn( ' ln_trainc, ln_dyninc, ln_sshinc and ln_seaiceinc are set to .false. :', &
[2128]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
[9213]258      IF( ln_asmiau ) THEN
259         !
[2128]260         ALLOCATE( wgtiau( icycper ) )
[9213]261         !
[9019]262         wgtiau(:) = 0._wp
[9213]263         !
264         !                                !---------------------------------------------------------
265         IF( niaufn == 0 ) THEN           ! Constant IAU forcing
266            !                             !---------------------------------------------------------
[2128]267            DO jt = 1, iiauper
268               wgtiau(jt+nitiaustr-1) = 1.0 / REAL( iiauper )
269            END DO
[9213]270            !                             !---------------------------------------------------------
271         ELSEIF ( niaufn == 1 ) THEN      ! Linear hat-like, centred in middle of IAU interval
272            !                             !---------------------------------------------------------
[2128]273            ! Compute the normalization factor
[9213]274            znorm = 0._wp
275            IF( MOD( iiauper, 2 ) == 0 ) THEN   ! Even number of time steps in IAU interval
[2128]276               imid = iiauper / 2 
277               DO jt = 1, imid
278                  znorm = znorm + REAL( jt )
279               END DO
280               znorm = 2.0 * znorm
[9213]281            ELSE                                ! Odd number of time steps in IAU interval
[2128]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
[9213]289            !
[2128]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
[9213]296            !
[2128]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
[9213]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
[2128]328#if defined key_asminc
[9213]329      ALLOCATE( ssh_iau      (jpi,jpj)     )   ;   ssh_iau      (:,:)   = 0._wp
[2128]330#endif
[8030]331#if defined key_cice && defined key_asminc
[9213]332      ALLOCATE( ndaice_da    (jpi,jpj)     )   ;   ndaice_da    (:,:)   = 0._wp
[8030]333#endif
[9213]334      !
335      IF ( ln_trainc .OR. ln_dyninc .OR.   &       !--------------------------------------
336         & ln_sshinc .OR. ln_seaiceinc   ) THEN    ! Read the increments from file
337         !                                         !--------------------------------------
[2128]338         CALL iom_open( c_asminc, inum )
[9213]339         !
340         CALL iom_get( inum, 'time'       , zdate_inc   ) 
[2128]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
[9213]345         !
[2128]346         IF(lwp) THEN
347            WRITE(numout,*) 
[9213]348            WRITE(numout,*) 'asm_inc_init : Assimilation increments valid between dates ', z_inc_dateb,' and ', z_inc_datef
[2128]349            WRITE(numout,*) '~~~~~~~~~~~~'
350         ENDIF
[9213]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' )
[2128]356
[6140]357         IF ( ( ln_asmdin ).AND.( zdate_inc /= ditdin_date ) ) &
[2128]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_autoglo, 'bckint', t_bkginc, 1 )
363            CALL iom_get( inum, jpdom_autoglo, '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_autoglo, 'bckinu', u_bkginc, 1 )             
375            CALL iom_get( inum, jpdom_autoglo, '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_autoglo, '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
[3764]394         IF ( ln_seaiceinc ) THEN
395            CALL iom_get( inum, jpdom_autoglo, '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
[9213]402         !
[2128]403         CALL iom_close( inum )
[9213]404         !
[2128]405      ENDIF
[9213]406      !
407      !                                            !--------------------------------------
408      IF ( ln_dyninc .AND. nn_divdmp > 0 ) THEN    ! Apply divergence damping filter
409         !                                         !--------------------------------------
[9019]410         ALLOCATE( zhdiv(jpi,jpj) ) 
[5836]411         !
412         DO jt = 1, nn_divdmp
413            !
[9019]414            DO jk = 1, jpkm1           ! zhdiv = e1e1 * div
415               zhdiv(:,:) = 0._wp
[3764]416               DO jj = 2, jpjm1
417                  DO ji = fs_2, fs_jpim1   ! vector opt.
[9019]418                     zhdiv(ji,jj) = (  e2u(ji  ,jj) * e3u_n(ji  ,jj,jk) * u_bkginc(ji  ,jj,jk)    &
419                        &            - e2u(ji-1,jj) * e3u_n(ji-1,jj,jk) * u_bkginc(ji-1,jj,jk)    &
420                        &            + e1v(ji,jj  ) * e3v_n(ji,jj  ,jk) * v_bkginc(ji,jj  ,jk)    &
421                        &            - e1v(ji,jj-1) * e3v_n(ji,jj-1,jk) * v_bkginc(ji,jj-1,jk)  ) / e3t_n(ji,jj,jk)
[3764]422                  END DO
[3294]423               END DO
[10425]424               CALL lbc_lnk( 'asminc', zhdiv, 'T', 1. )   ! lateral boundary cond. (no sign change)
[5836]425               !
[3764]426               DO jj = 2, jpjm1
427                  DO ji = fs_2, fs_jpim1   ! vector opt.
[6140]428                     u_bkginc(ji,jj,jk) = u_bkginc(ji,jj,jk)                         &
[9019]429                        &               + 0.2_wp * ( zhdiv(ji+1,jj) - zhdiv(ji  ,jj) ) * r1_e1u(ji,jj) * umask(ji,jj,jk)
[6140]430                     v_bkginc(ji,jj,jk) = v_bkginc(ji,jj,jk)                         &
[9019]431                        &               + 0.2_wp * ( zhdiv(ji,jj+1) - zhdiv(ji,jj  ) ) * r1_e2v(ji,jj) * vmask(ji,jj,jk) 
[3764]432                  END DO
[3294]433               END DO
434            END DO
[5836]435            !
[3764]436         END DO
[5836]437         !
[9019]438         DEALLOCATE( zhdiv ) 
[5836]439         !
[3294]440      ENDIF
[9213]441      !
442      !                             !-----------------------------------------------------
443      IF ( ln_asmdin ) THEN         ! Allocate and initialize the background state arrays
444         !                          !-----------------------------------------------------
[5836]445         !
[9213]446         ALLOCATE( t_bkg  (jpi,jpj,jpk) )   ;   t_bkg  (:,:,:) = 0._wp
447         ALLOCATE( s_bkg  (jpi,jpj,jpk) )   ;   s_bkg  (:,:,:) = 0._wp
448         ALLOCATE( u_bkg  (jpi,jpj,jpk) )   ;   u_bkg  (:,:,:) = 0._wp
449         ALLOCATE( v_bkg  (jpi,jpj,jpk) )   ;   v_bkg  (:,:,:) = 0._wp
450         ALLOCATE( ssh_bkg(jpi,jpj)     )   ;   ssh_bkg(:,:)   = 0._wp
[5836]451         !
452         !
[2128]453         !--------------------------------------------------------------------
454         ! Read from file the background state at analysis time
455         !--------------------------------------------------------------------
[5836]456         !
[2128]457         CALL iom_open( c_asmdin, inum )
[5836]458         !
[3764]459         CALL iom_get( inum, 'rdastp', zdate_bkg ) 
[5836]460         !
[2128]461         IF(lwp) THEN
462            WRITE(numout,*) 
[9213]463            WRITE(numout,*) '   ==>>>  Assimilation background state valid at : ', zdate_bkg
464            WRITE(numout,*)
[2128]465         ENDIF
[5836]466         !
[9213]467         IF ( zdate_bkg /= ditdin_date )   &
[2128]468            & CALL ctl_warn( ' Validity time of assimilation background state does', &
469            &                ' not agree with Direct Initialization time' )
[5836]470         !
[2128]471         IF ( ln_trainc ) THEN   
472            CALL iom_get( inum, jpdom_autoglo, 'tn', t_bkg )
473            CALL iom_get( inum, jpdom_autoglo, 'sn', s_bkg )
474            t_bkg(:,:,:) = t_bkg(:,:,:) * tmask(:,:,:)
475            s_bkg(:,:,:) = s_bkg(:,:,:) * tmask(:,:,:)
476         ENDIF
[5836]477         !
[2128]478         IF ( ln_dyninc ) THEN   
479            CALL iom_get( inum, jpdom_autoglo, 'un', u_bkg )
480            CALL iom_get( inum, jpdom_autoglo, 'vn', v_bkg )
481            u_bkg(:,:,:) = u_bkg(:,:,:) * umask(:,:,:)
482            v_bkg(:,:,:) = v_bkg(:,:,:) * vmask(:,:,:)
483         ENDIF
[5836]484         !
[2128]485         IF ( ln_sshinc ) THEN
486            CALL iom_get( inum, jpdom_autoglo, 'sshn', ssh_bkg )
487            ssh_bkg(:,:) = ssh_bkg(:,:) * tmask(:,:,1)
488         ENDIF
[5836]489         !
[2128]490         CALL iom_close( inum )
[5836]491         !
[2128]492      ENDIF
[2392]493      !
[9213]494      IF(lwp) WRITE(numout,*) '   ==>>>   Euler time step switch is ', neuler
495      !
496      IF( lk_asminc ) THEN                            !==  data assimilation  ==!
497         IF( ln_bkgwri )   CALL asm_bkg_wri( nit000 - 1 )      ! Output background fields
498         IF( ln_asmdin ) THEN                                  ! Direct initialization
[10954]499            IF( ln_trainc )   CALL tra_asm_inc( nit000 - 1, Kmm )      ! Tracers
[9213]500            IF( ln_dyninc )   CALL dyn_asm_inc( nit000 - 1 )      ! Dynamics
501            IF( ln_sshinc )   CALL ssh_asm_inc( nit000 - 1 )      ! SSH
502         ENDIF
503      ENDIF
504      !
[2128]505   END SUBROUTINE asm_inc_init
[9213]506   
507   
[10954]508   SUBROUTINE tra_asm_inc( kt, Kmm )
[2128]509      !!----------------------------------------------------------------------
510      !!                    ***  ROUTINE tra_asm_inc  ***
511      !!         
512      !! ** Purpose : Apply the tracer (T and S) assimilation increments
513      !!
514      !! ** Method  : Direct initialization or Incremental Analysis Updating
515      !!
516      !! ** Action  :
517      !!----------------------------------------------------------------------
[5836]518      INTEGER, INTENT(IN) ::   kt   ! Current time step
[10954]519      INTEGER, INTENT(IN) ::   Kmm  ! Current time level index
[2392]520      !
[5836]521      INTEGER  :: ji, jj, jk
522      INTEGER  :: it
[2128]523      REAL(wp) :: zincwgt  ! IAU weight for current time step
[3764]524      REAL (wp), DIMENSION(jpi,jpj,jpk) :: fzptnz ! 3d freezing point values
[2392]525      !!----------------------------------------------------------------------
[5836]526      !
[3764]527      ! freezing point calculation taken from oc_fz_pt (but calculated for all depths)
528      ! used to prevent the applied increments taking the temperature below the local freezing point
[4990]529      DO jk = 1, jpkm1
[6140]530        CALL eos_fzp( tsn(:,:,jk,jp_sal), fzptnz(:,:,jk), gdept_n(:,:,jk) )
[4990]531      END DO
[5836]532         !
533         !                             !--------------------------------------
534      IF ( ln_asmiau ) THEN            ! Incremental Analysis Updating
535         !                             !--------------------------------------
536         !
[2128]537         IF ( ( kt >= nitiaustr_r ).AND.( kt <= nitiaufin_r ) ) THEN
[5836]538            !
[2128]539            it = kt - nit000 + 1
540            zincwgt = wgtiau(it) / rdt   ! IAU weight for the current time step
[5836]541            !
[2128]542            IF(lwp) THEN
543               WRITE(numout,*) 
[4990]544               WRITE(numout,*) 'tra_asm_inc : Tracer IAU at time step = ', kt,' with IAU weight = ', wgtiau(it)
[2128]545               WRITE(numout,*) '~~~~~~~~~~~~'
546            ENDIF
[5836]547            !
[2128]548            ! Update the tracer tendencies
549            DO jk = 1, jpkm1
[3764]550               IF (ln_temnofreeze) THEN
551                  ! Do not apply negative increments if the temperature will fall below freezing
552                  WHERE(t_bkginc(:,:,jk) > 0.0_wp .OR. &
553                     &   tsn(:,:,jk,jp_tem) + tsa(:,:,jk,jp_tem) + t_bkginc(:,:,jk) * wgtiau(it) > fzptnz(:,:,jk) ) 
554                     tsa(:,:,jk,jp_tem) = tsa(:,:,jk,jp_tem) + t_bkginc(:,:,jk) * zincwgt 
555                  END WHERE
556               ELSE
557                  tsa(:,:,jk,jp_tem) = tsa(:,:,jk,jp_tem) + t_bkginc(:,:,jk) * zincwgt 
558               ENDIF
559               IF (ln_salfix) THEN
560                  ! Do not apply negative increments if the salinity will fall below a specified
561                  ! minimum value salfixmin
562                  WHERE(s_bkginc(:,:,jk) > 0.0_wp .OR. &
563                     &   tsn(:,:,jk,jp_sal) + tsa(:,:,jk,jp_sal) + s_bkginc(:,:,jk) * wgtiau(it) > salfixmin ) 
564                     tsa(:,:,jk,jp_sal) = tsa(:,:,jk,jp_sal) + s_bkginc(:,:,jk) * zincwgt
565                  END WHERE
566               ELSE
567                  tsa(:,:,jk,jp_sal) = tsa(:,:,jk,jp_sal) + s_bkginc(:,:,jk) * zincwgt
568               ENDIF
[2128]569            END DO
[5836]570            !
[2128]571         ENDIF
[5836]572         !
[2128]573         IF ( kt == nitiaufin_r + 1  ) THEN   ! For bias crcn to work
574            DEALLOCATE( t_bkginc )
575            DEALLOCATE( s_bkginc )
576         ENDIF
[5836]577         !                             !--------------------------------------
578      ELSEIF ( ln_asmdin ) THEN        ! Direct Initialization
579         !                             !--------------------------------------
580         !           
[2128]581         IF ( kt == nitdin_r ) THEN
[5836]582            !
[2128]583            neuler = 0  ! Force Euler forward step
[5836]584            !
[2128]585            ! Initialize the now fields with the background + increment
[3764]586            IF (ln_temnofreeze) THEN
587               ! Do not apply negative increments if the temperature will fall below freezing
[4990]588               WHERE( t_bkginc(:,:,:) > 0.0_wp .OR. tsn(:,:,:,jp_tem) + t_bkginc(:,:,:) > fzptnz(:,:,:) ) 
[3764]589                  tsn(:,:,:,jp_tem) = t_bkg(:,:,:) + t_bkginc(:,:,:)   
590               END WHERE
591            ELSE
592               tsn(:,:,:,jp_tem) = t_bkg(:,:,:) + t_bkginc(:,:,:)   
593            ENDIF
[2128]594            IF (ln_salfix) THEN
[3764]595               ! Do not apply negative increments if the salinity will fall below a specified
596               ! minimum value salfixmin
[4990]597               WHERE( s_bkginc(:,:,:) > 0.0_wp .OR. tsn(:,:,:,jp_sal) + s_bkginc(:,:,:) > salfixmin ) 
[3764]598                  tsn(:,:,:,jp_sal) = s_bkg(:,:,:) + s_bkginc(:,:,:)   
599               END WHERE
600            ELSE
601               tsn(:,:,:,jp_sal) = s_bkg(:,:,:) + s_bkginc(:,:,:)   
[2128]602            ENDIF
603
[4990]604            tsb(:,:,:,:) = tsn(:,:,:,:)                 ! Update before fields
[2128]605
[4990]606            CALL eos( tsb, rhd, rhop, gdept_0(:,:,:) )  ! Before potential and in situ densities
607!!gm  fabien
608!            CALL eos( tsb, rhd, rhop )                ! Before potential and in situ densities
609!!gm
610
[10954]611            IF( ln_zps .AND. .NOT. lk_c1d .AND. .NOT. ln_isfcav)           &
612               &  CALL zps_hde    ( kt, Kmm, jpts, tsb, gtsu, gtsv,        &  ! Partial steps: before horizontal gradient
613               &                              rhd, gru , grv               )  ! of t, s, rd at the last ocean level
614            IF( ln_zps .AND. .NOT. lk_c1d .AND.       ln_isfcav)                       &
615               &  CALL zps_hde_isf( nit000, Kmm, jpts, tsb, gtsu, gtsv, gtui, gtvi,    &  ! Partial steps for top cell (ISF)
616               &                                  rhd, gru , grv , grui, grvi          )  ! of t, s, rd at the last ocean level
[6140]617
[2128]618            DEALLOCATE( t_bkginc )
619            DEALLOCATE( s_bkginc )
620            DEALLOCATE( t_bkg    )
621            DEALLOCATE( s_bkg    )
622         ENDIF
[2392]623         
[2128]624      ENDIF
[3764]625      ! Perhaps the following call should be in step
[5836]626      IF ( ln_seaiceinc  )   CALL seaice_asm_inc ( kt )   ! apply sea ice concentration increment
[2392]627      !
[2128]628   END SUBROUTINE tra_asm_inc
629
[2392]630
[2128]631   SUBROUTINE dyn_asm_inc( kt )
632      !!----------------------------------------------------------------------
633      !!                    ***  ROUTINE dyn_asm_inc  ***
634      !!         
635      !! ** Purpose : Apply the dynamics (u and v) assimilation increments.
636      !!
637      !! ** Method  : Direct initialization or Incremental Analysis Updating.
638      !!
639      !! ** Action  :
640      !!----------------------------------------------------------------------
641      INTEGER, INTENT(IN) :: kt   ! Current time step
[2392]642      !
[2128]643      INTEGER :: jk
644      INTEGER :: it
645      REAL(wp) :: zincwgt  ! IAU weight for current time step
[2392]646      !!----------------------------------------------------------------------
[5836]647      !
648      !                          !--------------------------------------------
649      IF ( ln_asmiau ) THEN      ! Incremental Analysis Updating
650         !                       !--------------------------------------------
651         !
[2128]652         IF ( ( kt >= nitiaustr_r ).AND.( kt <= nitiaufin_r ) ) THEN
[5836]653            !
[2128]654            it = kt - nit000 + 1
655            zincwgt = wgtiau(it) / rdt   ! IAU weight for the current time step
[5836]656            !
[2128]657            IF(lwp) THEN
658               WRITE(numout,*) 
[5836]659               WRITE(numout,*) 'dyn_asm_inc : Dynamics IAU at time step = ', kt,' with IAU weight = ', wgtiau(it)
[2128]660               WRITE(numout,*) '~~~~~~~~~~~~'
661            ENDIF
[5836]662            !
[2128]663            ! Update the dynamic tendencies
664            DO jk = 1, jpkm1
665               ua(:,:,jk) = ua(:,:,jk) + u_bkginc(:,:,jk) * zincwgt
666               va(:,:,jk) = va(:,:,jk) + v_bkginc(:,:,jk) * zincwgt
667            END DO
[5836]668            !
[2128]669            IF ( kt == nitiaufin_r ) THEN
670               DEALLOCATE( u_bkginc )
671               DEALLOCATE( v_bkginc )
672            ENDIF
[5836]673            !
[2128]674         ENDIF
[5836]675         !                          !-----------------------------------------
676      ELSEIF ( ln_asmdin ) THEN     ! Direct Initialization
677         !                          !-----------------------------------------
678         !         
[2128]679         IF ( kt == nitdin_r ) THEN
[5836]680            !
[2128]681            neuler = 0                    ! Force Euler forward step
[5836]682            !
[2128]683            ! Initialize the now fields with the background + increment
684            un(:,:,:) = u_bkg(:,:,:) + u_bkginc(:,:,:)
685            vn(:,:,:) = v_bkg(:,:,:) + v_bkginc(:,:,:) 
[5836]686            !
[2128]687            ub(:,:,:) = un(:,:,:)         ! Update before fields
688            vb(:,:,:) = vn(:,:,:)
[5836]689            !
[2128]690            DEALLOCATE( u_bkg    )
691            DEALLOCATE( v_bkg    )
692            DEALLOCATE( u_bkginc )
693            DEALLOCATE( v_bkginc )
694         ENDIF
[2392]695         !
[2128]696      ENDIF
[2392]697      !
[2128]698   END SUBROUTINE dyn_asm_inc
699
[2392]700
[2128]701   SUBROUTINE ssh_asm_inc( kt )
702      !!----------------------------------------------------------------------
703      !!                    ***  ROUTINE ssh_asm_inc  ***
704      !!         
705      !! ** Purpose : Apply the sea surface height assimilation increment.
706      !!
707      !! ** Method  : Direct initialization or Incremental Analysis Updating.
708      !!
709      !! ** Action  :
710      !!----------------------------------------------------------------------
711      INTEGER, INTENT(IN) :: kt   ! Current time step
[2392]712      !
[2128]713      INTEGER :: it
[3764]714      INTEGER :: jk
[2128]715      REAL(wp) :: zincwgt  ! IAU weight for current time step
[2392]716      !!----------------------------------------------------------------------
[5836]717      !
718      !                             !-----------------------------------------
719      IF ( ln_asmiau ) THEN         ! Incremental Analysis Updating
720         !                          !-----------------------------------------
721         !
[2128]722         IF ( ( kt >= nitiaustr_r ).AND.( kt <= nitiaufin_r ) ) THEN
[5836]723            !
[2128]724            it = kt - nit000 + 1
725            zincwgt = wgtiau(it) / rdt   ! IAU weight for the current time step
[5836]726            !
[2128]727            IF(lwp) THEN
728               WRITE(numout,*) 
729               WRITE(numout,*) 'ssh_asm_inc : SSH IAU at time step = ', &
730                  &  kt,' with IAU weight = ', wgtiau(it)
731               WRITE(numout,*) '~~~~~~~~~~~~'
732            ENDIF
[5836]733            !
[2128]734            ! Save the tendency associated with the IAU weighted SSH increment
735            ! (applied in dynspg.*)
736#if defined key_asminc
737            ssh_iau(:,:) = ssh_bkginc(:,:) * zincwgt
738#endif
[5836]739            !
[8547]740         ELSE IF( kt == nitiaufin_r+1 ) THEN
[8546]741            !
[9465]742            ! test on ssh_bkginc needed as ssh_asm_inc is called twice by time step
743            IF (ALLOCATED(ssh_bkginc)) DEALLOCATE( ssh_bkginc )
744            !
745#if defined key_asminc
[8546]746            ssh_iau(:,:) = 0._wp
[9465]747#endif
[8546]748            !
[2128]749         ENDIF
[5836]750         !                          !-----------------------------------------
751      ELSEIF ( ln_asmdin ) THEN     ! Direct Initialization
752         !                          !-----------------------------------------
753         !
[2128]754         IF ( kt == nitdin_r ) THEN
[5836]755            !
756            neuler = 0                                   ! Force Euler forward step
757            !
758            sshn(:,:) = ssh_bkg(:,:) + ssh_bkginc(:,:)   ! Initialize the now fields the background + increment
759            !
760            sshb(:,:) = sshn(:,:)                        ! Update before fields
[6140]761            e3t_b(:,:,:) = e3t_n(:,:,:)
762!!gm why not e3u_b, e3v_b, gdept_b ????
[5836]763            !
[2128]764            DEALLOCATE( ssh_bkg    )
765            DEALLOCATE( ssh_bkginc )
[5836]766            !
[2128]767         ENDIF
[2392]768         !
[2128]769      ENDIF
[2392]770      !
[2128]771   END SUBROUTINE ssh_asm_inc
772
[9213]773
[9023]774   SUBROUTINE ssh_asm_div( kt, phdivn )
775      !!----------------------------------------------------------------------
776      !!                  ***  ROUTINE ssh_asm_div  ***
777      !!
778      !! ** Purpose :   ssh increment with z* is incorporated via a correction of the local divergence         
779      !!                across all the water column
780      !!
781      !! ** Method  :
782      !!                CAUTION : sshiau is positive (inflow) decreasing the
783      !!                          divergence and expressed in m/s
784      !!
785      !! ** Action  :   phdivn   decreased by the ssh increment
786      !!----------------------------------------------------------------------
787      INTEGER, INTENT(IN) :: kt                               ! ocean time-step index
788      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   phdivn   ! horizontal divergence
789      !!
790      INTEGER  ::   jk                                        ! dummy loop index
791      REAL(wp), DIMENSION(:,:)  , POINTER       ::   ztim     ! local array
792      !!----------------------------------------------------------------------
793      !
794#if defined key_asminc
795      CALL ssh_asm_inc( kt ) !==   (calculate increments)
796      !
797      IF( ln_linssh ) THEN
798         phdivn(:,:,1) = phdivn(:,:,1) - ssh_iau(:,:) / e3t_n(:,:,1) * tmask(:,:,1)
799      ELSE
[9044]800         ALLOCATE( ztim(jpi,jpj) )
[9023]801         ztim(:,:) = ssh_iau(:,:) / ( ht_n(:,:) + 1.0 - ssmask(:,:) )
802         DO jk = 1, jpkm1                                 
803            phdivn(:,:,jk) = phdivn(:,:,jk) - ztim(:,:) * tmask(:,:,jk) 
804         END DO
805         !
[9044]806         DEALLOCATE(ztim)
[9023]807      ENDIF
808#endif
809      !
810   END SUBROUTINE ssh_asm_div
[3785]811
[9213]812
[3764]813   SUBROUTINE seaice_asm_inc( kt, kindic )
814      !!----------------------------------------------------------------------
815      !!                    ***  ROUTINE seaice_asm_inc  ***
816      !!         
817      !! ** Purpose : Apply the sea ice assimilation increment.
818      !!
819      !! ** Method  : Direct initialization or Incremental Analysis Updating.
820      !!
821      !! ** Action  :
822      !!
823      !!----------------------------------------------------------------------
[5836]824      INTEGER, INTENT(in)           ::   kt       ! Current time step
[3785]825      INTEGER, INTENT(in), OPTIONAL ::   kindic   ! flag for disabling the deallocation
826      !
827      INTEGER  ::   it
828      REAL(wp) ::   zincwgt   ! IAU weight for current time step
[9570]829#if defined key_si3
[9656]830      REAL(wp), DIMENSION(jpi,jpj) ::   zofrld, zohicif, zseaicendg, zhicifinc
[3785]831      REAL(wp) ::   zhicifmin = 0.5_wp      ! ice minimum depth in metres
[3764]832#endif
[3785]833      !!----------------------------------------------------------------------
[5836]834      !
835      !                             !-----------------------------------------
836      IF ( ln_asmiau ) THEN         ! Incremental Analysis Updating
837         !                          !-----------------------------------------
838         !
[3764]839         IF ( ( kt >= nitiaustr_r ).AND.( kt <= nitiaufin_r ) ) THEN
[5836]840            !
[3764]841            it = kt - nit000 + 1
842            zincwgt = wgtiau(it)      ! IAU weight for the current time step
843            ! note this is not a tendency so should not be divided by rdt (as with the tracer and other increments)
[5836]844            !
[3764]845            IF(lwp) THEN
846               WRITE(numout,*) 
[5836]847               WRITE(numout,*) 'seaice_asm_inc : sea ice conc IAU at time step = ', kt,' with IAU weight = ', wgtiau(it)
[3764]848               WRITE(numout,*) '~~~~~~~~~~~~'
849            ENDIF
[5836]850            !
[9656]851            ! Sea-ice : SI3 case
[5836]852            !
[9570]853#if defined key_si3
[9019]854            zofrld (:,:) = 1._wp - at_i(:,:)
855            zohicif(:,:) = hm_i(:,:)
[3785]856            !
[9019]857            at_i  (:,:) = 1. - MIN( MAX( 1.-at_i  (:,:) - seaice_bkginc(:,:) * zincwgt, 0.0_wp), 1.0_wp)
858            at_i_b(:,:) = 1. - MIN( MAX( 1.-at_i_b(:,:) - seaice_bkginc(:,:) * zincwgt, 0.0_wp), 1.0_wp)
859            fr_i(:,:) = at_i(:,:)        ! adjust ice fraction
[3785]860            !
[9019]861            zseaicendg(:,:) = zofrld(:,:) - (1. - at_i(:,:))   ! find out actual sea ice nudge applied
[3785]862            !
[3764]863            ! Nudge sea ice depth to bring it up to a required minimum depth
[9019]864            WHERE( zseaicendg(:,:) > 0.0_wp .AND. hm_i(:,:) < zhicifmin ) 
865               zhicifinc(:,:) = (zhicifmin - hm_i(:,:)) * zincwgt   
[3764]866            ELSEWHERE
867               zhicifinc(:,:) = 0.0_wp
868            END WHERE
[3785]869            !
870            ! nudge ice depth
[9019]871            hm_i (:,:) = hm_i (:,:) + zhicifinc(:,:)
[3785]872            !
873            ! seaice salinity balancing (to add)
[3764]874#endif
[9213]875            !
[4245]876#if defined key_cice && defined key_asminc
[3785]877            ! Sea-ice : CICE case. Pass ice increment tendency into CICE
[3764]878            ndaice_da(:,:) = seaice_bkginc(:,:) * zincwgt / rdt
879#endif
[9213]880            !
[3764]881            IF ( kt == nitiaufin_r ) THEN
882               DEALLOCATE( seaice_bkginc )
883            ENDIF
[9213]884            !
[3764]885         ELSE
[9213]886            !
[4245]887#if defined key_cice && defined key_asminc
[5836]888            ndaice_da(:,:) = 0._wp        ! Sea-ice : CICE case. Zero ice increment tendency into CICE
[3764]889#endif
[9213]890            !
[3764]891         ENDIF
[5836]892         !                          !-----------------------------------------
893      ELSEIF ( ln_asmdin ) THEN     ! Direct Initialization
894         !                          !-----------------------------------------
895         !
[3764]896         IF ( kt == nitdin_r ) THEN
[5836]897            !
[3764]898            neuler = 0                    ! Force Euler forward step
[5836]899            !
[9656]900            ! Sea-ice : SI3 case
[5836]901            !
[9570]902#if defined key_si3
[9019]903            zofrld (:,:) = 1._wp - at_i(:,:)
904            zohicif(:,:) = hm_i(:,:)
[3785]905            !
[3764]906            ! Initialize the now fields the background + increment
[9019]907            at_i(:,:) = 1. - MIN( MAX( 1.-at_i(:,:) - seaice_bkginc(:,:), 0.0_wp), 1.0_wp)
908            at_i_b(:,:) = at_i(:,:) 
909            fr_i(:,:) = at_i(:,:)        ! adjust ice fraction
[3785]910            !
[9019]911            zseaicendg(:,:) = zofrld(:,:) - (1. - at_i(:,:))   ! find out actual sea ice nudge applied
912            !
[3764]913            ! Nudge sea ice depth to bring it up to a required minimum depth
[9019]914            WHERE( zseaicendg(:,:) > 0.0_wp .AND. hm_i(:,:) < zhicifmin ) 
[10355]915               zhicifinc(:,:) = zhicifmin - hm_i(:,:)
[3764]916            ELSEWHERE
[9019]917               zhicifinc(:,:) = 0.0_wp
[3764]918            END WHERE
[3785]919            !
920            ! nudge ice depth
[9019]921            hm_i (:,:) = hm_i (:,:) + zhicifinc(:,:)
[3785]922            !
923            ! seaice salinity balancing (to add)
[3764]924#endif
[5836]925            !
[4245]926#if defined key_cice && defined key_asminc
927            ! Sea-ice : CICE case. Pass ice increment tendency into CICE
[3764]928           ndaice_da(:,:) = seaice_bkginc(:,:) / rdt
929#endif
[5836]930            IF ( .NOT. PRESENT(kindic) ) THEN
931               DEALLOCATE( seaice_bkginc )
932            END IF
933            !
[3764]934         ELSE
[5836]935            !
936#if defined key_cice && defined key_asminc
937            ndaice_da(:,:) = 0._wp     ! Sea-ice : CICE case. Zero ice increment tendency into CICE
[3764]938#endif
[5836]939            !
[3764]940         ENDIF
941
[9570]942!#if defined defined key_si3 || defined key_cice
[3764]943!
944!            IF (ln_seaicebal ) THEN       
945!             !! balancing salinity increments
946!             !! simple case from limflx.F90 (doesn't include a mass flux)
947!             !! assumption is that as ice concentration is reduced or increased
948!             !! the snow and ice depths remain constant
949!             !! note that snow is being created where ice concentration is being increased
950!             !! - could be more sophisticated and
951!             !! not do this (but would need to alter h_snow)
952!
953!             usave(:,:,:)=sb(:,:,:)   ! use array as a temporary store
954!
955!             DO jj = 1, jpj
956!               DO ji = 1, jpi
957!           ! calculate change in ice and snow mass per unit area
958!           ! positive values imply adding salt to the ocean (results from ice formation)
959!           ! fwf : ice formation and melting
960!
961!                 zfons = ( -nfresh_da(ji,jj)*soce + nfsalt_da(ji,jj) )*rdt
962!
963!           ! change salinity down to mixed layer depth
964!                 mld=hmld_kara(ji,jj)
965!
966!           ! prevent small mld
967!           ! less than 10m can cause salinity instability
968!                 IF (mld < 10) mld=10
969!
970!           ! set to bottom of a level
971!                 DO jk = jpk-1, 2, -1
972!                   IF ((mld > gdepw(ji,jj,jk)) .and. (mld < gdepw(ji,jj,jk+1))) THEN
973!                     mld=gdepw(ji,jj,jk+1)
974!                     jkmax=jk
975!                   ENDIF
976!                 ENDDO
977!
978!            ! avoid applying salinity balancing in shallow water or on land
979!            !
980!
981!            ! dsal_ocn (psu kg m^-2) / (kg m^-3 * m)
982!
983!                 dsal_ocn=0.0_wp
984!                 sal_thresh=5.0_wp        ! minimum salinity threshold for salinity balancing
985!
986!                 if (tmask(ji,jj,1) > 0 .AND. tmask(ji,jj,jkmax) > 0 ) &
987!                              dsal_ocn = zfons / (rhop(ji,jj,1) * mld)
988!
989!           ! put increments in for levels in the mixed layer
990!           ! but prevent salinity below a threshold value
991!
992!                   DO jk = 1, jkmax             
993!
994!                     IF (dsal_ocn > 0.0_wp .or. sb(ji,jj,jk)+dsal_ocn > sal_thresh) THEN
995!                           sb(ji,jj,jk) = sb(ji,jj,jk) + dsal_ocn
996!                           sn(ji,jj,jk) = sn(ji,jj,jk) + dsal_ocn
997!                     ENDIF
998!
999!                   ENDDO
1000!
1001!      !            !  salt exchanges at the ice/ocean interface
1002!      !            zpmess         = zfons / rdt_ice    ! rdt_ice is ice timestep
1003!      !
1004!      !! Adjust fsalt. A +ve fsalt means adding salt to ocean
1005!      !!           fsalt(ji,jj) =  fsalt(ji,jj) + zpmess     ! adjust fsalt 
1006!      !!               
1007!      !!           emps(ji,jj) = emps(ji,jj) + zpmess        ! or adjust emps (see icestp1d)
1008!      !!                                                     ! E-P (kg m-2 s-2)
1009!      !            emp(ji,jj) = emp(ji,jj) + zpmess          ! E-P (kg m-2 s-2)
1010!               ENDDO !ji
1011!             ENDDO !jj!
1012!
1013!            ENDIF !ln_seaicebal
1014!
1015!#endif
[5836]1016         !
[3764]1017      ENDIF
[5836]1018      !
[3764]1019   END SUBROUTINE seaice_asm_inc
[3785]1020   
[2392]1021   !!======================================================================
[2128]1022END MODULE asminc
Note: See TracBrowser for help on using the repository browser.