New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
asminc.F90 in branches/UKMO/test_moci_test_suite_namelist_read/NEMOGCM/NEMO/OPA_SRC/ASM – NEMO

source: branches/UKMO/test_moci_test_suite_namelist_read/NEMOGCM/NEMO/OPA_SRC/ASM/asminc.F90 @ 9383

Last change on this file since 9383 was 9383, checked in by andmirek, 6 years ago

#2050 fixes and changes

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