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

source: branches/2014/dev_CNRS0_NOC1_LDF/NEMOGCM/NEMO/OPA_SRC/ASM/asminc.F90 @ 4616

Last change on this file since 4616 was 4616, checked in by gm, 10 years ago

#1260 : see the associated wiki page for explanation

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