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/2021/dev_r14116_HPC-10_mcastril_Mixed_Precision_implementation/src/OCE/ASM – NEMO

source: NEMO/branches/2021/dev_r14116_HPC-10_mcastril_Mixed_Precision_implementation/src/OCE/ASM/asminc.F90 @ 15540

Last change on this file since 15540 was 15540, checked in by sparonuz, 3 years ago

Mixed precision version, tested up to 30 years on ORCA2.

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