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

source: branches/2012/dev_MERGE_2012/NEMOGCM/NEMO/OPA_SRC/ASM/asminc.F90 @ 3764

Last change on this file since 3764 was 3764, checked in by smasson, 11 years ago

dev_MERGE_2012: report bugfixes done in the trunk from r3555 to r3763 into dev_MERGE_2012

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