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

source: branches/UKMO/dev_r5518_GO6_package_inc_asm/NEMOGCM/NEMO/OPA_SRC/ASM/asminc.F90 @ 7979

Last change on this file since 7979 was 7979, checked in by jwhile, 7 years ago

Minor update to ASM for seaice

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