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

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

source: branches/UKMO/dev_r5518_obs_oper_update_icethick/NEMOGCM/NEMO/OPA_SRC/ASM/asminc.F90 @ 10181

Last change on this file since 10181 was 10181, checked in by emmafiedler, 6 years ago

Working version of ice thickness assimilation updates

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