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

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

source: branches/UKMO/r6232_obs_oper_update/NEMOGCM/NEMO/OPA_SRC/ASM/asminc.F90 @ 11202

Last change on this file since 11202 was 11202, checked in by jcastill, 5 years ago

Copy of branch branches/UKMO/dev_r5518_obs_oper_update@11130 without namelist_ref changes to allow merging with coupling and biogeochemistry branches

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