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

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

source: branches/2014/dev_r4650_UKMO7_STARTHOUR/NEMOGCM/NEMO/OPA_SRC/ASM/asminc.F90 @ 5980

Last change on this file since 5980 was 5980, checked in by timgraham, 8 years ago

Upgraded to v3.6 revision of trunk (r5518)

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