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

source: branches/UKMO/dev_r5518_25hr_mean_assim_bkg/NEMOGCM/NEMO/OPA_SRC/ASM/asminc.F90 @ 6761

Last change on this file since 6761 was 6761, checked in by kingr, 8 years ago

Merged branches/UKMO/dev_r5518_v3.4_asm_nemovar_community@6754

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