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

source: branches/2014/dev_r4650_UKMO7_STARTHOUR/NEMOGCM/NEMO/OPA_SRC/ASM/asminc.F90 @ 4811

Last change on this file since 4811 was 4811, checked in by djlea, 10 years ago

Update reference namelists with nn_time0. Fix to daymod. Fix to offline domain.F90 to get nn_time0 from the namelist. Assimilation date time updates.

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