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

source: branches/UKMO/dev_r5518_AMM15_package/NEMOGCM/NEMO/OPA_SRC/ASM/asminc.F90 @ 10249

Last change on this file since 10249 was 10249, checked in by kingr, 5 years ago

Merged AMM15_v3_6_STABLE_package_collate@10237

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