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/UKMO/dev_r5518_GO6_package_inc_asm/NEMOGCM/NEMO/OPA_SRC/ASM – NEMO

source: branches/UKMO/dev_r5518_GO6_package_inc_asm/NEMOGCM/NEMO/OPA_SRC/ASM/asminc.F90 @ 7955

Last change on this file since 7955 was 7955, checked in by jwhile, 8 years ago

Minor bug fix

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