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 NEMO/branches/2020/dev_r14116_HPC-04_mcastril_Mixed_Precision_implementation_final/src/OCE/ASM – NEMO

source: NEMO/branches/2020/dev_r14116_HPC-04_mcastril_Mixed_Precision_implementation_final/src/OCE/ASM/asminc.F90 @ 14219

Last change on this file since 14219 was 14219, checked in by mcastril, 3 years ago

Add Mixed Precision support by Oriol Tintó

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