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

source: branches/UKMO/CO5_package_branch/NEMOGCM/NEMO/OPA_SRC/ASM/asminc.F90 @ 7367

Last change on this file since 7367 was 7367, checked in by deazer, 7 years ago

Contains merged code for CO5 reference.

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