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

source: branches/2015/dev_MetOffice_merge_2015/NEMOGCM/NEMO/OPA_SRC/ASM/asminc.F90 @ 5997

Last change on this file since 5997 was 5997, checked in by timgraham, 8 years ago

Merged branches/2014/dev_r4650_UKMO7_STARTHOUR into branch

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