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_v3.4_asm_nemovar_community/NEMOGCM/NEMO/OPA_SRC/ASM – NEMO

source: branches/UKMO/dev_r5518_v3.4_asm_nemovar_community/NEMOGCM/NEMO/OPA_SRC/ASM/asminc.F90 @ 6630

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

Adpated changes required to apply 2D surft increments in AMM7 from NEMO3.4 branch /branches/dev/frwe/vn3.4_ASM_NEMOVAR_community.

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