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

source: branches/2012/dev_MERGE_2012/NEMOGCM/NEMO/OPA_SRC/ASM/asminc.F90 @ 3785

Last change on this file since 3785 was 3785, checked in by gm, 11 years ago

dev_MERGE_2012: #1051 : bug correction to compile LIM3 with ASM and with ICB

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