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

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

source: branches/UKMO/AMM15_v3_6_STABLE_package_collate/NEMOGCM/NEMO/OPA_SRC/ASM/asminc.F90 @ 14190

Last change on this file since 14190 was 14190, checked in by kingr, 4 years ago

Allow heat/salt conservation when assimilating SSH to be namelist controlled.

File size: 65.7 KB
Line 
1MODULE asminc
2   !!======================================================================
3   !!                       ***  MODULE asminc  ***
4   !! Assimilation increment : Apply an increment generated by data
5   !!                          assimilation
6   !!======================================================================
7   !! History :       ! 2007-03  (M. Martin)  Met Office version
8   !!                 ! 2007-04  (A. Weaver)  calc_date original code
9   !!                 ! 2007-04  (A. Weaver)  Merge with OPAVAR/NEMOVAR
10   !!   NEMO     3.3  ! 2010-05  (D. Lea)  Update to work with NEMO v3.2
11   !!             -   ! 2010-05  (D. Lea)  add calc_month_len routine based on day_init
12   !!            3.4  ! 2012-10  (A. Weaver and K. Mogensen) Fix for direct initialization
13   !!----------------------------------------------------------------------
14
15   !!----------------------------------------------------------------------
16   !!   'key_asminc'   : Switch on the assimilation increment interface
17   !!----------------------------------------------------------------------
18   !!   asm_inc_init   : Initialize the increment arrays and IAU weights
19   !!   calc_date      : Compute the calendar date YYYYMMDD on a given step
20   !!   tra_asm_inc    : Apply the tracer (T and S) increments
21   !!   dyn_asm_inc    : Apply the dynamic (u and v) increments
22   !!   ssh_asm_inc    : Apply the SSH increment
23   !!   seaice_asm_inc : Apply the sea ice concentration increment
24   !!   sit_asm_inc    : Apply the sea ice thickness increment
25   !!----------------------------------------------------------------------
26   USE wrk_nemo         ! Memory Allocation
27   USE par_oce          ! Ocean space and time domain variables
28   USE dom_oce          ! Ocean space and time domain
29   USE domvvl           ! domain: variable volume level
30   USE oce              ! Dynamics and active tracers defined in memory
31   USE ldfdyn_oce       ! ocean dynamics: lateral physics
32   USE eosbn2           ! Equation of state - in situ and potential density
33   USE zpshde           ! Partial step : Horizontal Derivative
34   USE iom              ! Library to read input files
35   USE asmpar           ! Parameters for the assmilation interface
36   USE c1d              ! 1D initialization
37   USE in_out_manager   ! I/O manager
38   USE lib_mpp          ! MPP library
39#if defined key_lim2
40   USE ice_2            ! LIM2
41#endif
42   USE sbc_oce          ! Surface boundary condition variables.
43   USE zdfmxl, ONLY :  & 
44   &  hmld_tref,       &   
45#if defined key_karaml
46   &  hmld_kara,       &
47   &  ln_kara,         &
48#endif   
49   &  hmld,            & 
50   &  hmlp,            &
51   &  hmlpt
52#if defined key_bdy 
53   USE bdy_oce, ONLY: bdytmask 
54#endif 
55   USE asmbgc           ! Biogeochemistry assimilation
56
57   IMPLICIT NONE
58   PRIVATE
59   
60   PUBLIC   asm_inc_init   !: Initialize the increment arrays and IAU weights
61   PUBLIC   calc_date      !: Compute the calendar date YYYYMMDD on a given step
62   PUBLIC   tra_asm_inc    !: Apply the tracer (T and S) increments
63   PUBLIC   dyn_asm_inc    !: Apply the dynamic (u and v) increments
64   PUBLIC   ssh_asm_inc    !: Apply the SSH increment
65   PUBLIC   seaice_asm_inc !: Apply the seaice concentration increment
66   PUBLIC   sit_asm_inc    !: Apply the seaice thickness increment
67   PUBLIC   bgc_asm_inc    !: Apply the biogeochemistry increments
68
69#if defined key_asminc
70    LOGICAL, PUBLIC, PARAMETER :: lk_asminc = .TRUE.   !: Logical switch for assimilation increment interface
71#else
72    LOGICAL, PUBLIC, PARAMETER :: lk_asminc = .FALSE.  !: No assimilation increments
73#endif
74   LOGICAL, PUBLIC :: ln_bkgwri = .FALSE.      !: No output of the background state fields
75   LOGICAL, PUBLIC :: ln_avgbkg = .FALSE.      !: No output of the mean background state fields
76   LOGICAL, PUBLIC :: ln_asmiau = .FALSE.      !: No applying forcing with an assimilation increment
77   LOGICAL, PUBLIC :: ln_asmdin = .FALSE.      !: No direct initialization
78   LOGICAL, PUBLIC :: ln_trainc = .FALSE.      !: No tracer (T and S) assimilation increments
79   LOGICAL, PUBLIC :: ln_dyninc = .FALSE.      !: No dynamics (u and v) assimilation increments
80   LOGICAL, PUBLIC :: ln_ssh_hs_cons = .FALSE. !: Conserve heat and salt when adding SSH increment
81   LOGICAL, PUBLIC :: ln_sshinc = .FALSE.      !: No sea surface height assimilation increment
82   LOGICAL, PUBLIC :: ln_seaiceinc = .FALSE.   !: No sea ice concentration increment
83   LOGICAL, PUBLIC :: ln_sitinc = .FALSE.      !: No sea ice thickness increment
84   LOGICAL, PUBLIC :: lk_bgcinc = .FALSE.      !: No biogeochemistry increments
85   LOGICAL, PUBLIC :: ln_salfix = .FALSE.      !: Apply minimum salinity check
86   LOGICAL, PUBLIC :: ln_temnofreeze = .FALSE. !: Don't allow the temperature to drop below freezing
87   INTEGER, PUBLIC :: nn_divdmp                !: Apply divergence damping filter nn_divdmp times
88
89   REAL(wp), PUBLIC, DIMENSION(:,:,:), ALLOCATABLE ::   t_bkg   , s_bkg      !: Background temperature and salinity
90   REAL(wp), PUBLIC, DIMENSION(:,:,:), ALLOCATABLE ::   u_bkg   , v_bkg      !: Background u- & v- velocity components
91   REAL(wp), PUBLIC, DIMENSION(:,:,:), ALLOCATABLE ::   t_bkginc, s_bkginc   !: Increment to the background T & S
92   REAL(wp), PUBLIC, DIMENSION(:,:,:), ALLOCATABLE ::   u_bkginc, v_bkginc   !: Increment to the u- & v-components
93   REAL(wp), PUBLIC, DIMENSION(:)    , ALLOCATABLE ::   wgtiau               !: IAU weights for each time step
94#if defined key_asminc
95   REAL(wp), PUBLIC, DIMENSION(:,:), ALLOCATABLE ::   ssh_iau           !: IAU-weighted sea surface height increment
96#endif
97   !                                !!! time steps relative to the cycle interval [0,nitend-nit000-1]
98   INTEGER , PUBLIC ::   nitbkg      !: Time step of the background state used in the Jb term
99   INTEGER , PUBLIC ::   nitdin      !: Time step of the background state for direct initialization
100   INTEGER , PUBLIC ::   nitiaustr   !: Time step of the start of the IAU interval
101   INTEGER , PUBLIC ::   nitiaufin   !: Time step of the end of the IAU interval
102   INTEGER , PUBLIC ::   nitavgbkg   !: Number of timesteps to average assim bkg [0,nitavgbkg]
103   !
104   INTEGER , PUBLIC ::   niaufn      !: Type of IAU weighing function: = 0   Constant weighting
105   !                                 !: = 1   Linear hat-like, centred in middle of IAU interval
106   REAL(wp), PUBLIC ::   salfixmin   !: Ensure that the salinity is larger than this  value if (ln_salfix)
107
108   REAL(wp), DIMENSION(:,:), ALLOCATABLE ::   ssh_bkg, ssh_bkginc   ! Background sea surface height and its increment
109   REAL(wp), DIMENSION(:,:), ALLOCATABLE ::   seaice_bkginc         ! Increment to the background sea ice conc
110
111   INTEGER :: mld_choice        = 4   !: choice of mld criteria to use for physics assimilation
112                                      !: 1) hmld      - Turbocline/mixing depth                           [W points]
113                                      !: 2) hmlp      - Density criterion (0.01 kg/m^3 change from 10m)   [W points]
114                                      !: 3) hmld_kara - Kara MLD                                          [Interpolated]
115                                      !: 4) hmld_tref - Temperature criterion (0.2 K change from surface) [T points]
116
117
118   !! * Substitutions
119#  include "domzgr_substitute.h90"
120#  include "ldfdyn_substitute.h90"
121#  include "vectopt_loop_substitute.h90"
122   !!----------------------------------------------------------------------
123   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
124   !! $Id$
125   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
126   !!----------------------------------------------------------------------
127CONTAINS
128
129   SUBROUTINE asm_inc_init
130      !!----------------------------------------------------------------------
131      !!                    ***  ROUTINE asm_inc_init  ***
132      !!         
133      !! ** Purpose : Initialize the assimilation increment and IAU weights.
134      !!
135      !! ** Method  : Initialize the assimilation increment and IAU weights.
136      !!
137      !! ** Action  :
138      !!----------------------------------------------------------------------
139      INTEGER :: ji, jj, jk, jt  ! dummy loop indices
140      INTEGER :: imid, inum      ! local integers
141      INTEGER :: ios             ! Local integer output status for namelist read
142      INTEGER :: iiauper         ! Number of time steps in the IAU period
143      INTEGER :: icycper         ! Number of time steps in the cycle
144      INTEGER :: iitend_date     ! Date YYYYMMDD of final time step
145      INTEGER :: iitbkg_date     ! Date YYYYMMDD of background time step for Jb term
146      INTEGER :: iitdin_date     ! Date YYYYMMDD of background time step for DI
147      INTEGER :: iitiaustr_date  ! Date YYYYMMDD of IAU interval start time step
148      INTEGER :: iitiaufin_date  ! Date YYYYMMDD of IAU interval final time step
149      INTEGER :: isurfstat       ! Local integer for status of reading surft variable
150      INTEGER :: iitavgbkg_date  ! Date YYYYMMDD of end of assim bkg averaging period
151      !
152      REAL(wp) :: znorm        ! Normalization factor for IAU weights
153      REAL(wp) :: ztotwgt      ! Value of time-integrated IAU weights (should be equal to one)
154      REAL(wp) :: z_inc_dateb  ! Start date of interval on which increment is valid
155      REAL(wp) :: z_inc_datef  ! End date of interval on which increment is valid
156      REAL(wp) :: zdate_bkg    ! Date in background state file for DI
157      REAL(wp) :: zdate_inc    ! Time axis in increments file
158      !
159      REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: & 
160          &       t_bkginc_2d  ! file for reading in 2D 
161      !                        ! temperature increments
162      REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: & 
163          &       z_mld     ! Mixed layer depth
164         
165      REAL(wp), POINTER, DIMENSION(:,:) ::   hdiv   ! 2D workspace
166      !
167      LOGICAL :: lk_surft      ! Logical: T => Increments file contains surft variable
168                               !               so only apply surft increments.
169      !!
170      NAMELIST/nam_asminc/ ln_bkgwri, ln_avgbkg, ln_balwri,                &
171         &                 ln_trainc, ln_dyninc, ln_sshinc,                &
172         &                 ln_phytobal, ln_slchltotinc, ln_slchldiainc,    &
173         &                 ln_slchlnaninc, ln_slchlpicinc, ln_slchldininc, &
174         &                 ln_slchlnoninc, ln_schltotinc, ln_slphytotinc,  &
175         &                 ln_slphydiainc, ln_slphynoninc, ln_spco2inc,    &
176         &                 ln_sfco2inc, ln_plchltotinc, ln_pchltotinc,     &
177         &                 ln_plchldiainc, ln_plchlnaninc, ln_plchlpicinc, &
178         &                 ln_plchldininc,                                 & 
179         &                 ln_pno3inc, ln_psi4inc, ln_pdicinc, ln_palkinc, &
180         &                 ln_pphinc, ln_po2inc, ln_ppo4inc,               &
181         &                 ln_asmdin, ln_asmiau, ln_ssh_hs_cons,           &
182         &                 nitbkg, nitdin, nitiaustr, nitiaufin, niaufn,   &
183         &                 ln_salfix, salfixmin, nn_divdmp, nitavgbkg,     &
184         &                 mld_choice, mld_choice_bgc, rn_maxchlinc
185      !!----------------------------------------------------------------------
186
187      !-----------------------------------------------------------------------
188      ! Read Namelist nam_asminc : assimilation increment interface
189      !-----------------------------------------------------------------------
190
191      ! Set default values
192      ln_bkgwri = .FALSE.
193      ln_avgbkg = .FALSE.
194      ln_trainc = .FALSE.
195      ln_dyninc = .FALSE.
196      ln_sshinc = .FALSE.
197      ln_asmdin = .FALSE.
198      ln_asmiau = .TRUE.
199      ln_salfix = .FALSE.
200
201      ln_ssh_hs_cons = .FALSE.
202      ln_seaiceinc = .FALSE.
203      ln_sitinc = .FALSE.
204      ln_temnofreeze = .FALSE.
205
206      salfixmin = -9999
207      nitbkg    = 0
208      nitdin    = 0     
209      nitiaustr = 1
210      nitiaufin = 150
211      niaufn    = 0
212      nitavgbkg = 1
213
214      REWIND( numnam_ref )              ! Namelist nam_asminc in reference namelist : Assimilation increment
215      READ  ( numnam_ref, nam_asminc, IOSTAT = ios, ERR = 901)
216901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'nam_asminc in reference namelist', lwp )
217
218      REWIND( numnam_cfg )              ! Namelist nam_asminc in configuration namelist : Assimilation increment
219      READ  ( numnam_cfg, nam_asminc, IOSTAT = ios, ERR = 902 )
220902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'nam_asminc in configuration namelist', lwp )
221      IF(lwm) WRITE ( numond, nam_asminc )
222
223      ! Control print
224      IF(lwp) THEN
225         WRITE(numout,*)
226         WRITE(numout,*) 'asm_inc_init : Assimilation increment initialization :'
227         WRITE(numout,*) '~~~~~~~~~~~~'
228         WRITE(numout,*) '   Namelist nam_asminc : set assimilation increment parameters'
229         WRITE(numout,*) '      Logical switch for writing out background state          ln_bkgwri = ', ln_bkgwri
230         WRITE(numout,*) '      Logical switch for writing mean background state         ln_avgbkg = ', ln_avgbkg
231         WRITE(numout,*) '      Logical switch for writing out balancing increments      ln_balwri = ', ln_balwri
232         WRITE(numout,*) '      Logical switch for applying tracer increments            ln_trainc = ', ln_trainc
233         WRITE(numout,*) '      Logical switch for applying velocity increments          ln_dyninc = ', ln_dyninc
234         WRITE(numout,*) '      Logical switch for conserving heat/salt when applying SSH increments ln_ssh_hs_cons = ', ln_ssh_hs_cons
235         WRITE(numout,*) '      Logical switch for applying SSH increments               ln_sshinc = ', ln_sshinc
236         WRITE(numout,*) '      Logical switch for applying SIC increments               ln_seaiceinc = ', ln_seaiceinc
237         WRITE(numout,*) '      Logical switch for applying SIT increments               ln_sitinc = ', ln_sitinc
238         WRITE(numout,*) '      Logical switch for Direct Initialization (DI)            ln_asmdin = ', ln_asmdin
239         WRITE(numout,*) '      Logical switch for Incremental Analysis Updating (IAU)   ln_asmiau = ', ln_asmiau
240         WRITE(numout,*) '      Timestep of background in [0,nitend-nit000-1]            nitbkg    = ', nitbkg
241         WRITE(numout,*) '      Timestep of background for DI in [0,nitend-nit000-1]     nitdin    = ', nitdin
242         WRITE(numout,*) '      Timestep of start of IAU interval in [0,nitend-nit000-1] nitiaustr = ', nitiaustr
243         WRITE(numout,*) '      Timestep of end of IAU interval in [0,nitend-nit000-1]   nitiaufin = ', nitiaufin
244         WRITE(numout,*) '      Number of timesteps to average assim bkg [0,nitavgbkg]   nitavgbkg = ', nitavgbkg
245         WRITE(numout,*) '      Type of IAU weighting function                           niaufn    = ', niaufn
246         WRITE(numout,*) '      Logical switch for ensuring that the sa > salfixmin      ln_salfix = ', ln_salfix
247         WRITE(numout,*) '      Minimum salinity after applying the increments           salfixmin = ', salfixmin
248         WRITE(numout,*) '      Choice of MLD for physics assimilation                  mld_choice = ', mld_choice
249         WRITE(numout,*) '      Choice of MLD for BGC assimilation                  mld_choice_bgc = ', mld_choice_bgc
250         WRITE(numout,*) '      Maximum absolute chlorophyll increment (<=0 = off)    rn_maxchlinc = ', rn_maxchlinc
251         WRITE(numout,*) '      Logical switch for phytoplankton balancing               ln_phytobal = ', ln_phytobal
252         WRITE(numout,*) '      Logical switch for applying slchltot increments          ln_slchltotinc = ', ln_slchltotinc
253         WRITE(numout,*) '      Logical switch for applying slchldia increments          ln_slchldiainc = ', ln_slchldiainc
254         WRITE(numout,*) '      Logical switch for applying slchlnon increments          ln_slchlnoninc = ', ln_slchlnoninc
255         WRITE(numout,*) '      Logical switch for applying slchlnan increments          ln_slchlnaninc = ', ln_slchlnaninc
256         WRITE(numout,*) '      Logical switch for applying slchlpic increments          ln_slchlpicinc = ', ln_slchlpicinc
257         WRITE(numout,*) '      Logical switch for applying slchldin increments          ln_slchldininc = ', ln_slchldininc
258         WRITE(numout,*) '      Logical switch for applying schltot increments           ln_schltotinc = ', ln_schltotinc
259         WRITE(numout,*) '      Logical switch for applying slphytot increments          ln_slphytotinc = ', ln_slphytotinc
260         WRITE(numout,*) '      Logical switch for applying slphydia increments          ln_slphydiainc = ', ln_slphydiainc
261         WRITE(numout,*) '      Logical switch for applying slphynon increments          ln_slphynoninc = ', ln_slphynoninc
262         WRITE(numout,*) '      Logical switch for applying spco2 increments             ln_spco2inc = ', ln_spco2inc
263         WRITE(numout,*) '      Logical switch for applying sfco2 increments             ln_sfco2inc = ', ln_sfco2inc
264         WRITE(numout,*) '      Logical switch for applying plchltot increments          ln_plchltotinc = ', ln_plchltotinc
265         WRITE(numout,*) '      Logical switch for applying pchltot increments           ln_pchltotinc = ', ln_pchltotinc
266         WRITE(numout,*) '      Logical switch for applying plchldia increments          ln_plchldiainc = ', ln_plchldiainc
267         WRITE(numout,*) '      Logical switch for applying plchlnan increments          ln_plchlnaninc = ', ln_plchlnaninc
268         WRITE(numout,*) '      Logical switch for applying plchlpic increments          ln_plchlpicinc = ', ln_plchlpicinc
269         WRITE(numout,*) '      Logical switch for applying plchldin increments          ln_plchldininc = ', ln_plchldininc
270         WRITE(numout,*) '      Logical switch for applying pno3 increments              ln_pno3inc = ', ln_pno3inc
271         WRITE(numout,*) '      Logical switch for applying psi4 increments              ln_psi4inc = ', ln_psi4inc
272         WRITE(numout,*) '      Logical switch for applying ppo4 increments              ln_ppo4inc = ', ln_ppo4inc
273         WRITE(numout,*) '      Logical switch for applying pdic increments              ln_pdicinc = ', ln_pdicinc
274         WRITE(numout,*) '      Logical switch for applying palk increments              ln_palkinc = ', ln_palkinc
275         WRITE(numout,*) '      Logical switch for applying pph increments               ln_pphinc = ', ln_pphinc
276         WRITE(numout,*) '      Logical switch for applying po2 increments               ln_po2inc = ', ln_po2inc
277      ENDIF
278
279      nitbkg_r    = nitbkg    + nit000 - 1  ! Background time referenced to nit000
280      nitdin_r    = nitdin    + nit000 - 1  ! Background time for DI referenced to nit000
281      nitiaustr_r = nitiaustr + nit000 - 1  ! Start of IAU interval referenced to nit000
282      nitiaufin_r = nitiaufin + nit000 - 1  ! End of IAU interval referenced to nit000
283      nitavgbkg_r = nitavgbkg + nit000 - 1  ! Averaging period referenced to nit000
284
285      iiauper = nitiaufin_r - nitiaustr_r + 1  ! IAU interval length
286      icycper = nitend      - nit000      + 1  ! Cycle interval length
287
288      CALL calc_date( nit000, nitend     , ndate0, iitend_date    )     ! Date of final time step
289      CALL calc_date( nit000, nitbkg_r   , ndate0, iitbkg_date    )     ! Background time for Jb referenced to ndate0
290      CALL calc_date( nit000, nitdin_r   , ndate0, iitdin_date    )     ! Background time for DI referenced to ndate0
291      CALL calc_date( nit000, nitiaustr_r, ndate0, iitiaustr_date )     ! IAU start time referenced to ndate0
292      CALL calc_date( nit000, nitiaufin_r, ndate0, iitiaufin_date )     ! IAU end time referenced to ndate0
293      CALL calc_date( nit000, nitavgbkg_r, ndate0, iitavgbkg_date )     ! End of assim bkg averaging period referenced to ndate0
294      !
295      IF(lwp) THEN
296         WRITE(numout,*)
297         WRITE(numout,*) '   Time steps referenced to current cycle:'
298         WRITE(numout,*) '       iitrst      = ', nit000 - 1
299         WRITE(numout,*) '       nit000      = ', nit000
300         WRITE(numout,*) '       nitend      = ', nitend
301         WRITE(numout,*) '       nitbkg_r    = ', nitbkg_r
302         WRITE(numout,*) '       nitdin_r    = ', nitdin_r
303         WRITE(numout,*) '       nitiaustr_r = ', nitiaustr_r
304         WRITE(numout,*) '       nitiaufin_r = ', nitiaufin_r
305         WRITE(numout,*) '       nitavgbkg_r = ', nitavgbkg_r
306         WRITE(numout,*)
307         WRITE(numout,*) '   Dates referenced to current cycle:'
308         WRITE(numout,*) '       ndastp         = ', ndastp
309         WRITE(numout,*) '       ndate0         = ', ndate0
310         WRITE(numout,*) '       iitend_date    = ', iitend_date
311         WRITE(numout,*) '       iitbkg_date    = ', iitbkg_date
312         WRITE(numout,*) '       iitdin_date    = ', iitdin_date
313         WRITE(numout,*) '       iitiaustr_date = ', iitiaustr_date
314         WRITE(numout,*) '       iitiaufin_date = ', iitiaufin_date
315         WRITE(numout,*) '       iitavgbkg_date = ', iitavgbkg_date
316      ENDIF
317      IF ( ln_slchltotinc .OR. ln_slchldiainc .OR. ln_slchlnoninc .OR. &
318         & ln_slchlnaninc .OR. ln_slchlpicinc .OR. ln_slchldininc .OR. &
319         & ln_schltotinc  .OR. ln_slphytotinc .OR. ln_slphydiainc .OR. &
320         & ln_slphynoninc .OR. ln_spco2inc    .OR. ln_sfco2inc    .OR. &
321         & ln_plchltotinc .OR. ln_pchltotinc  .OR. ln_plchldiainc .OR. &
322         & ln_plchlnaninc .OR. ln_plchlpicinc .OR. ln_plchldininc .OR. &
323         & ln_pno3inc     .OR.                                         &
324         & ln_psi4inc     .OR. ln_pdicinc     .OR. ln_palkinc     .OR. &
325         & ln_pphinc      .OR. ln_po2inc      .OR. ln_ppo4inc ) THEN
326         lk_bgcinc = .TRUE.
327      ENDIF
328
329      IF ( nacc /= 0 ) &
330         & CALL ctl_stop( ' nacc /= 0 and key_asminc :',  &
331         &                ' Assimilation increments have only been implemented', &
332         &                ' for synchronous time stepping' )
333
334      IF ( ( ln_asmdin ).AND.( ln_asmiau ) )   &
335         & CALL ctl_stop( ' ln_asmdin and ln_asmiau :', &
336         &                ' Choose Direct Initialization OR Incremental Analysis Updating')
337
338      IF (      ( ( .NOT. ln_asmdin ).AND.( .NOT. ln_asmiau ) ) &
339         & .AND.( ( ln_trainc ).OR.( ln_dyninc ).OR.( ln_sshinc ).OR.( ln_seaiceinc ).OR. &
340         &        ( ln_sitinc ).OR.( lk_bgcinc ) )) &
341         & CALL ctl_stop( ' One or more of ln_trainc, ln_dyninc, ln_sshinc, ln_seaiceinc,', &
342         &                ' ln_sitinc and ln_(bgc-variable)inc is set to .true.', &
343         &                ' but ln_asmdin and ln_asmiau are both set to .false. :', &
344         &                ' Inconsistent options')
345
346      IF ( ( niaufn /= 0 ).AND.( niaufn /= 1 ) ) &
347         & CALL ctl_stop( ' niaufn /= 0 or niaufn /=1 :',  &
348         &                ' Type IAU weighting function is invalid')
349
350      IF ( ( .NOT. ln_trainc ).AND.( .NOT. ln_dyninc ).AND.( .NOT. ln_sshinc ).AND.( .NOT. ln_seaiceinc ) &
351         & .AND.( .NOT. ln_sitinc ).AND.( .NOT. lk_bgcinc ) )  &
352         & CALL ctl_warn( ' ln_trainc, ln_dyninc, ln_sshinc, ln_seaiceinc,', &
353         &                ' ln_sitinc and ln_(bgc-variable)inc are set to .false. :', &
354         &                ' The assimilation increments are not applied')
355
356      IF ( ( ln_asmiau ).AND.( nitiaustr == nitiaufin ) ) &
357         & CALL ctl_stop( ' nitiaustr = nitiaufin :',  &
358         &                ' IAU interval is of zero length')
359
360      IF ( ( ln_asmiau ).AND.( ( nitiaustr_r < nit000 ).OR.( nitiaufin_r > nitend ) ) ) &
361         & CALL ctl_stop( ' nitiaustr or nitiaufin :',  &
362         &                ' IAU starting or final time step is outside the cycle interval', &
363         &                 ' Valid range nit000 to nitend')
364
365      IF ( ( nitbkg_r < nit000 - 1 ).OR.( nitbkg_r > nitend ) ) &
366         & CALL ctl_stop( ' nitbkg :',  &
367         &                ' Background time step is outside the cycle interval')
368
369      IF ( ( nitdin_r < nit000 - 1 ).OR.( nitdin_r > nitend ) ) &
370         & CALL ctl_stop( ' nitdin :',  &
371         &                ' Background time step for Direct Initialization is outside', &
372         &                ' the cycle interval')
373
374      IF ( nitavgbkg_r > nitend ) &
375         & CALL ctl_stop( ' nitavgbkg_r :',  &
376         &                ' Assim bkg averaging period is outside', &
377         &                ' the cycle interval')
378     
379      IF ( lk_bgcinc ) CALL asm_bgc_check_options
380
381      IF ( nstop > 0 ) RETURN       ! if there are any errors then go no further
382
383      !--------------------------------------------------------------------
384      ! Initialize the Incremental Analysis Updating weighting function
385      !--------------------------------------------------------------------
386
387      IF ( ln_asmiau ) THEN
388
389         ALLOCATE( wgtiau( icycper ) )
390
391         wgtiau(:) = 0.0
392
393         IF ( niaufn == 0 ) THEN
394
395            !---------------------------------------------------------
396            ! Constant IAU forcing
397            !---------------------------------------------------------
398
399            DO jt = 1, iiauper
400               wgtiau(jt+nitiaustr-1) = 1.0 / REAL( iiauper )
401            END DO
402
403         ELSEIF ( niaufn == 1 ) THEN
404
405            !---------------------------------------------------------
406            ! Linear hat-like, centred in middle of IAU interval
407            !---------------------------------------------------------
408
409            ! Compute the normalization factor
410            znorm = 0.0
411            IF ( MOD( iiauper, 2 ) == 0 ) THEN  ! Even number of time steps in IAU interval
412               imid = iiauper / 2 
413               DO jt = 1, imid
414                  znorm = znorm + REAL( jt )
415               END DO
416               znorm = 2.0 * znorm
417            ELSE                               ! Odd number of time steps in IAU interval
418               imid = ( iiauper + 1 ) / 2       
419               DO jt = 1, imid - 1
420                  znorm = znorm + REAL( jt )
421               END DO
422               znorm = 2.0 * znorm + REAL( imid )
423            ENDIF
424            znorm = 1.0 / znorm
425
426            DO jt = 1, imid - 1
427               wgtiau(jt+nitiaustr-1) = REAL( jt ) * znorm
428            END DO
429            DO jt = imid, iiauper
430               wgtiau(jt+nitiaustr-1) = REAL( iiauper - jt + 1 ) * znorm
431            END DO
432
433         ENDIF
434
435         ! Test that the integral of the weights over the weighting interval equals 1
436          IF(lwp) THEN
437             WRITE(numout,*)
438             WRITE(numout,*) 'asm_inc_init : IAU weights'
439             WRITE(numout,*) '~~~~~~~~~~~~'
440             WRITE(numout,*) '             time step         IAU  weight'
441             WRITE(numout,*) '             =========     ====================='
442             ztotwgt = 0.0
443             DO jt = 1, icycper
444                ztotwgt = ztotwgt + wgtiau(jt)
445                WRITE(numout,*) '         ', jt, '       ', wgtiau(jt) 
446             END DO   
447             WRITE(numout,*) '         ==================================='
448             WRITE(numout,*) '         Time-integrated weight = ', ztotwgt
449             WRITE(numout,*) '         ==================================='
450          ENDIF
451         
452      ENDIF
453
454      !--------------------------------------------------------------------
455      ! Allocate and initialize the increment arrays
456      !--------------------------------------------------------------------
457
458      IF ( ln_trainc ) THEN
459         ALLOCATE( t_bkginc(jpi,jpj,jpk) )
460         ALLOCATE( s_bkginc(jpi,jpj,jpk) )
461         t_bkginc(:,:,:) = 0.0
462         s_bkginc(:,:,:) = 0.0
463      ENDIF
464      IF ( ln_dyninc ) THEN
465         ALLOCATE( u_bkginc(jpi,jpj,jpk) )
466         ALLOCATE( v_bkginc(jpi,jpj,jpk) )
467         u_bkginc(:,:,:) = 0.0
468         v_bkginc(:,:,:) = 0.0
469      ENDIF
470      IF ( ln_sshinc ) THEN
471         ALLOCATE( ssh_bkginc(jpi,jpj)   )
472         ssh_bkginc(:,:) = 0.0
473      ENDIF
474      IF ( ln_seaiceinc ) THEN
475         ALLOCATE( seaice_bkginc(jpi,jpj))
476         seaice_bkginc(:,:) = 0.0
477      ENDIF
478      IF ( ln_sitinc ) THEN
479         ALLOCATE( sit_bkginc(jpi,jpj))
480         sit_bkginc(:,:) = 0.0
481      ENDIF
482#if defined key_asminc
483      ALLOCATE( ssh_iau(jpi,jpj)      )
484      ssh_iau(:,:)    = 0.0
485#endif
486      IF ( ( ln_trainc ).OR.( ln_dyninc ).OR.( ln_sshinc ).OR.( ln_seaiceinc ) &
487         &  .OR.( ln_sitinc ).OR.( lk_bgcinc ) ) THEN
488
489         !--------------------------------------------------------------------
490         ! Read the increments from file
491         !--------------------------------------------------------------------
492
493         CALL iom_open( c_asminc, inum )
494
495         CALL iom_get( inum, 'time', zdate_inc ) 
496
497         CALL iom_get( inum, 'z_inc_dateb', z_inc_dateb )
498         CALL iom_get( inum, 'z_inc_datef', z_inc_datef )
499         z_inc_dateb = zdate_inc
500         z_inc_datef = zdate_inc
501
502         IF(lwp) THEN
503            WRITE(numout,*) 
504            WRITE(numout,*) 'asm_inc_init : Assimilation increments valid ', &
505               &            ' between dates ', NINT( z_inc_dateb ),' and ',  &
506               &            NINT( z_inc_datef )
507            WRITE(numout,*) '~~~~~~~~~~~~'
508         ENDIF
509
510         IF (     ( NINT( z_inc_dateb ) < ndastp      ) &
511            & .OR.( NINT( z_inc_datef ) > iitend_date ) ) &
512            & CALL ctl_warn( ' Validity time of assimilation increments is ', &
513            &                ' outside the assimilation interval' )
514
515         IF ( ( ln_asmdin ).AND.( NINT( zdate_inc ) /= iitdin_date ) ) &
516            & CALL ctl_warn( ' Validity time of assimilation increments does ', &
517            &                ' not agree with Direct Initialization time' )
518
519         IF ( ln_trainc ) THEN   
520
521            !Test if the increments file contains the surft variable.
522            isurfstat = iom_varid( inum, 'bckinsurft', ldstop = .FALSE. )
523            IF ( isurfstat == -1 ) THEN
524               lk_surft = .FALSE.
525            ELSE
526               lk_surft = .TRUE.
527               CALL ctl_warn( ' Applying 2D temperature increment to bottom of ML: ', &
528            &                 ' bckinsurft found in increments file.' )
529            ENDIF             
530
531            IF (lk_surft) THEN
532               
533               ALLOCATE(z_mld(jpi,jpj)) 
534               SELECT CASE(mld_choice) 
535               CASE(1) 
536                  z_mld = hmld 
537               CASE(2) 
538                  z_mld = hmlp 
539               CASE(3) 
540#if defined key_karaml
541                  IF ( ln_kara ) THEN
542                     z_mld = hmld_kara
543                  ELSE
544                     CALL ctl_stop("Kara mixed layer not calculated as ln_kara=.false.")
545                  ENDIF
546#else
547                  CALL ctl_stop("Kara mixed layer not defined in current version of NEMO")  ! JW: Safety feature, should be removed
548                                                                                            ! once the Kara mixed layer is available
549#endif
550               CASE(4) 
551                  z_mld = hmld_tref 
552               END SELECT       
553                     
554               ALLOCATE( t_bkginc_2d(jpi,jpj) ) 
555               CALL iom_get( inum, jpdom_autoglo, 'bckinsurft', t_bkginc_2d, 1) 
556#if defined key_bdy               
557               DO jk = 1,jpkm1 
558                  WHERE( z_mld(:,:) > fsdepw(:,:,jk) ) 
559                     t_bkginc(:,:,jk) = t_bkginc_2d(:,:) * 0.5 * &
560                     &       ( 1 + cos( (fsdept(:,:,jk)/z_mld(:,:) ) * rpi ) )
561                     
562                     t_bkginc(:,:,jk) = t_bkginc(:,:,jk) * bdytmask(:,:) 
563                  ELSEWHERE 
564                     t_bkginc(:,:,jk) = 0. 
565                  ENDWHERE 
566               ENDDO 
567#else
568               t_bkginc(:,:,:) = 0. 
569#endif               
570               s_bkginc(:,:,:) = 0. 
571               
572               DEALLOCATE(z_mld, t_bkginc_2d) 
573           
574            ELSE
575               
576               CALL iom_get( inum, jpdom_autoglo, 'bckint', t_bkginc, 1 )
577               CALL iom_get( inum, jpdom_autoglo, 'bckins', s_bkginc, 1 )
578               ! Apply the masks
579               t_bkginc(:,:,:) = t_bkginc(:,:,:) * tmask(:,:,:)
580               s_bkginc(:,:,:) = s_bkginc(:,:,:) * tmask(:,:,:)
581               ! Set missing increments to 0.0 rather than 1e+20
582               ! to allow for differences in masks
583               WHERE( ABS( t_bkginc(:,:,:) ) > 1.0e+10 ) t_bkginc(:,:,:) = 0.0
584               WHERE( ABS( s_bkginc(:,:,:) ) > 1.0e+10 ) s_bkginc(:,:,:) = 0.0
585         
586            ENDIF
587
588         ENDIF
589
590         IF ( ln_dyninc ) THEN   
591            CALL iom_get( inum, jpdom_autoglo, 'bckinu', u_bkginc, 1 )             
592            CALL iom_get( inum, jpdom_autoglo, 'bckinv', v_bkginc, 1 )             
593            ! Apply the masks
594            u_bkginc(:,:,:) = u_bkginc(:,:,:) * umask(:,:,:)
595            v_bkginc(:,:,:) = v_bkginc(:,:,:) * vmask(:,:,:)
596            ! Set missing increments to 0.0 rather than 1e+20
597            ! to allow for differences in masks
598            WHERE( ABS( u_bkginc(:,:,:) ) > 1.0e+10 ) u_bkginc(:,:,:) = 0.0
599            WHERE( ABS( v_bkginc(:,:,:) ) > 1.0e+10 ) v_bkginc(:,:,:) = 0.0
600         ENDIF
601       
602         IF ( ln_sshinc ) THEN
603            CALL iom_get( inum, jpdom_autoglo, 'bckineta', ssh_bkginc, 1 )
604            ! Apply the masks
605            ssh_bkginc(:,:) = ssh_bkginc(:,:) * tmask(:,:,1)
606            ! Set missing increments to 0.0 rather than 1e+20
607            ! to allow for differences in masks
608            WHERE( ABS( ssh_bkginc(:,:) ) > 1.0e+10 ) ssh_bkginc(:,:) = 0.0
609         ENDIF
610
611         IF ( ln_seaiceinc ) THEN
612            CALL iom_get( inum, jpdom_autoglo, 'bckinseaice', seaice_bkginc, 1 )
613            ! Apply the masks
614            seaice_bkginc(:,:) = seaice_bkginc(:,:) * tmask(:,:,1)
615            ! Set missing increments to 0.0 rather than 1e+20
616            ! to allow for differences in masks
617            WHERE( ABS( seaice_bkginc(:,:) ) > 1.0e+10 ) seaice_bkginc(:,:) = 0.0
618         ENDIF
619
620         IF ( lk_bgcinc ) THEN
621            CALL asm_bgc_init_incs( inum )
622         ENDIF
623
624         CALL iom_close( inum )
625 
626      ENDIF
627
628      !-----------------------------------------------------------------------
629      ! Apply divergence damping filter
630      !-----------------------------------------------------------------------
631
632      IF ( ln_dyninc .AND. nn_divdmp > 0 ) THEN
633
634         CALL wrk_alloc(jpi,jpj,hdiv) 
635
636         DO  jt = 1, nn_divdmp
637
638            DO jk = 1, jpkm1
639
640               hdiv(:,:) = 0._wp
641
642               DO jj = 2, jpjm1
643                  DO ji = fs_2, fs_jpim1   ! vector opt.
644                     hdiv(ji,jj) =   &
645                        (  e2u(ji  ,jj  ) * fse3u(ji  ,jj  ,jk) * u_bkginc(ji  ,jj  ,jk)     &
646                         - e2u(ji-1,jj  ) * fse3u(ji-1,jj  ,jk) * u_bkginc(ji-1,jj  ,jk)     &
647                         + e1v(ji  ,jj  ) * fse3v(ji  ,jj  ,jk) * v_bkginc(ji  ,jj  ,jk)     &
648                         - e1v(ji  ,jj-1) * fse3v(ji  ,jj-1,jk) * v_bkginc(ji  ,jj-1,jk)  )  &
649                         / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) )
650                  END DO
651               END DO
652
653               CALL lbc_lnk( hdiv, 'T', 1. )   ! lateral boundary cond. (no sign change)
654
655               DO jj = 2, jpjm1
656                  DO ji = fs_2, fs_jpim1   ! vector opt.
657                     u_bkginc(ji,jj,jk) = u_bkginc(ji,jj,jk) + 0.2_wp * ( e1t(ji+1,jj)*e2t(ji+1,jj) * hdiv(ji+1,jj)   &
658                                                                        - e1t(ji  ,jj)*e2t(ji  ,jj) * hdiv(ji  ,jj) ) &
659                                                                      / e1u(ji,jj) * umask(ji,jj,jk) 
660                     v_bkginc(ji,jj,jk) = v_bkginc(ji,jj,jk) + 0.2_wp * ( e1t(ji,jj+1)*e2t(ji,jj+1) * hdiv(ji,jj+1)   &
661                                                                        - e1t(ji,jj  )*e2t(ji,jj  ) * hdiv(ji,jj  ) ) &
662                                                                      / e2v(ji,jj) * vmask(ji,jj,jk) 
663                  END DO
664               END DO
665
666            END DO
667
668         END DO
669
670         CALL wrk_dealloc(jpi,jpj,hdiv) 
671
672      ENDIF
673
674
675
676      !-----------------------------------------------------------------------
677      ! Allocate and initialize the background state arrays
678      !-----------------------------------------------------------------------
679
680      IF ( ln_asmdin ) THEN
681
682         ALLOCATE( t_bkg(jpi,jpj,jpk) )
683         ALLOCATE( s_bkg(jpi,jpj,jpk) )
684         ALLOCATE( u_bkg(jpi,jpj,jpk) )
685         ALLOCATE( v_bkg(jpi,jpj,jpk) )
686         ALLOCATE( ssh_bkg(jpi,jpj)   )
687
688         t_bkg(:,:,:) = 0.0
689         s_bkg(:,:,:) = 0.0
690         u_bkg(:,:,:) = 0.0
691         v_bkg(:,:,:) = 0.0
692         ssh_bkg(:,:) = 0.0
693
694         !--------------------------------------------------------------------
695         ! Read from file the background state at analysis time
696         !--------------------------------------------------------------------
697
698         CALL iom_open( c_asmdin, inum )
699
700         CALL iom_get( inum, 'rdastp', zdate_bkg ) 
701       
702         IF(lwp) THEN
703            WRITE(numout,*) 
704            WRITE(numout,*) 'asm_inc_init : Assimilation background state valid at : ', &
705               &  NINT( zdate_bkg )
706            WRITE(numout,*) '~~~~~~~~~~~~'
707         ENDIF
708
709         IF ( NINT( zdate_bkg ) /= iitdin_date ) &
710            & CALL ctl_warn( ' Validity time of assimilation background state does', &
711            &                ' not agree with Direct Initialization time' )
712
713         IF ( ln_trainc ) THEN   
714            CALL iom_get( inum, jpdom_autoglo, 'tn', t_bkg )
715            CALL iom_get( inum, jpdom_autoglo, 'sn', s_bkg )
716            t_bkg(:,:,:) = t_bkg(:,:,:) * tmask(:,:,:)
717            s_bkg(:,:,:) = s_bkg(:,:,:) * tmask(:,:,:)
718         ENDIF
719
720         IF ( ln_dyninc ) THEN   
721            CALL iom_get( inum, jpdom_autoglo, 'un', u_bkg )
722            CALL iom_get( inum, jpdom_autoglo, 'vn', v_bkg )
723            u_bkg(:,:,:) = u_bkg(:,:,:) * umask(:,:,:)
724            v_bkg(:,:,:) = v_bkg(:,:,:) * vmask(:,:,:)
725         ENDIF
726       
727         IF ( ln_sshinc ) THEN
728            CALL iom_get( inum, jpdom_autoglo, 'sshn', ssh_bkg )
729            ssh_bkg(:,:) = ssh_bkg(:,:) * tmask(:,:,1)
730         ENDIF
731
732         CALL iom_close( inum )
733
734      ENDIF
735         
736      IF ( lk_bgcinc ) THEN
737         CALL asm_bgc_init_bkg
738      ENDIF
739      !
740   END SUBROUTINE asm_inc_init
741
742
743   SUBROUTINE calc_date( kit000, kt, kdate0, kdate )
744      !!----------------------------------------------------------------------
745      !!                    ***  ROUTINE calc_date  ***
746      !!         
747      !! ** Purpose : Compute the calendar date YYYYMMDD at a given time step.
748      !!
749      !! ** Method  : Compute the calendar date YYYYMMDD at a given time step.
750      !!
751      !! ** Action  :
752      !!----------------------------------------------------------------------
753      INTEGER, INTENT(IN) :: kit000  ! Initial time step
754      INTEGER, INTENT(IN) :: kt      ! Current time step referenced to kit000
755      INTEGER, INTENT(IN) :: kdate0  ! Initial date
756      INTEGER, INTENT(OUT) :: kdate  ! Current date reference to kdate0
757      !
758      INTEGER :: iyea0    ! Initial year
759      INTEGER :: imon0    ! Initial month
760      INTEGER :: iday0    ! Initial day
761      INTEGER :: iyea     ! Current year
762      INTEGER :: imon     ! Current month
763      INTEGER :: iday     ! Current day
764      INTEGER :: idaystp  ! Number of days between initial and current date
765      INTEGER :: idaycnt  ! Day counter
766
767      INTEGER, DIMENSION(12) ::   imonth_len    !: length in days of the months of the current year
768
769      !-----------------------------------------------------------------------
770      ! Compute the calendar date YYYYMMDD
771      !-----------------------------------------------------------------------
772
773      ! Initial date
774      iyea0 =   kdate0 / 10000
775      imon0 = ( kdate0 - ( iyea0 * 10000 ) ) / 100
776      iday0 =   kdate0 - ( iyea0 * 10000 ) - ( imon0 * 100 ) 
777
778      ! Check that kt >= kit000 - 1
779      IF ( kt < kit000 - 1 ) CALL ctl_stop( ' kt must be >= kit000 - 1')
780
781      ! If kt = kit000 - 1 then set the date to the restart date
782      IF ( kt == kit000 - 1 ) THEN
783
784         kdate = ndastp
785         RETURN
786
787      ENDIF
788
789      ! Compute the number of days from the initial date
790      idaystp = INT( REAL( kt - kit000 ) * rdt / 86400. )
791   
792      iday    = iday0
793      imon    = imon0
794      iyea    = iyea0
795      idaycnt = 0
796
797      CALL calc_month_len( iyea, imonth_len )
798
799      DO WHILE ( idaycnt < idaystp )
800         iday = iday + 1
801         IF ( iday > imonth_len(imon) )  THEN
802            iday = 1
803            imon = imon + 1
804         ENDIF
805         IF ( imon > 12 ) THEN
806            imon = 1
807            iyea = iyea + 1
808            CALL calc_month_len( iyea, imonth_len )  ! update month lengths
809         ENDIF                 
810         idaycnt = idaycnt + 1
811      END DO
812      !
813      kdate = iyea * 10000 + imon * 100 + iday
814      !
815   END SUBROUTINE
816
817
818   SUBROUTINE calc_month_len( iyear, imonth_len )
819      !!----------------------------------------------------------------------
820      !!                    ***  ROUTINE calc_month_len  ***
821      !!         
822      !! ** Purpose : Compute the number of days in a months given a year.
823      !!
824      !! ** Method  :
825      !!----------------------------------------------------------------------
826      INTEGER, DIMENSION(12) ::   imonth_len    !: length in days of the months of the current year
827      INTEGER :: iyear         !: year
828      !!----------------------------------------------------------------------
829      !
830      ! length of the month of the current year (from nleapy, read in namelist)
831      IF ( nleapy < 2 ) THEN
832         imonth_len(:) = (/ 31, 28, 31, 30, 31, 30, 31, 31, 30, 31, 30, 31 /)
833         IF ( nleapy == 1 ) THEN   ! we are using calendar with leap years
834            IF ( MOD(iyear, 4) == 0 .AND. ( MOD(iyear, 400) == 0 .OR. MOD(iyear, 100) /= 0 ) ) THEN
835               imonth_len(2) = 29
836            ENDIF
837         ENDIF
838      ELSE
839         imonth_len(:) = nleapy   ! all months with nleapy days per year
840      ENDIF
841      !
842   END SUBROUTINE
843
844
845   SUBROUTINE tra_asm_inc( kt )
846      !!----------------------------------------------------------------------
847      !!                    ***  ROUTINE tra_asm_inc  ***
848      !!         
849      !! ** Purpose : Apply the tracer (T and S) assimilation increments
850      !!
851      !! ** Method  : Direct initialization or Incremental Analysis Updating
852      !!
853      !! ** Action  :
854      !!----------------------------------------------------------------------
855      INTEGER, INTENT(IN) :: kt               ! Current time step
856      !
857      INTEGER :: ji,jj,jk
858      INTEGER :: it
859      REAL(wp) :: zincwgt  ! IAU weight for current time step
860      REAL (wp), DIMENSION(jpi,jpj,jpk) :: fzptnz ! 3d freezing point values
861      !!----------------------------------------------------------------------
862
863      ! freezing point calculation taken from oc_fz_pt (but calculated for all depths)
864      ! used to prevent the applied increments taking the temperature below the local freezing point
865
866      DO jk = 1, jpkm1
867        CALL eos_fzp( tsn(:,:,jk,jp_sal), fzptnz(:,:,jk), fsdept(:,:,jk) )
868      END DO
869
870      IF ( ln_asmiau ) THEN
871
872         !--------------------------------------------------------------------
873         ! Incremental Analysis Updating
874         !--------------------------------------------------------------------
875
876         IF ( ( kt >= nitiaustr_r ).AND.( kt <= nitiaufin_r ) ) THEN
877
878            it = kt - nit000 + 1
879            zincwgt = wgtiau(it) / rdt   ! IAU weight for the current time step
880
881            IF(lwp) THEN
882               WRITE(numout,*) 
883               WRITE(numout,*) 'tra_asm_inc : Tracer IAU at time step = ', kt,' with IAU weight = ', wgtiau(it)
884               WRITE(numout,*) '~~~~~~~~~~~~'
885            ENDIF
886
887            ! Update the tracer tendencies
888            DO jk = 1, jpkm1
889               IF (ln_temnofreeze) THEN
890                  ! Do not apply negative increments if the temperature will fall below freezing
891                  WHERE(t_bkginc(:,:,jk) > 0.0_wp .OR. &
892                     &   tsn(:,:,jk,jp_tem) + tsa(:,:,jk,jp_tem) + t_bkginc(:,:,jk) * wgtiau(it) > fzptnz(:,:,jk) ) 
893                     tsa(:,:,jk,jp_tem) = tsa(:,:,jk,jp_tem) + t_bkginc(:,:,jk) * zincwgt 
894                  END WHERE
895               ELSE
896                  tsa(:,:,jk,jp_tem) = tsa(:,:,jk,jp_tem) + t_bkginc(:,:,jk) * zincwgt 
897               ENDIF
898               IF (ln_salfix) THEN
899                  ! Do not apply negative increments if the salinity will fall below a specified
900                  ! minimum value salfixmin
901                  WHERE(s_bkginc(:,:,jk) > 0.0_wp .OR. &
902                     &   tsn(:,:,jk,jp_sal) + tsa(:,:,jk,jp_sal) + s_bkginc(:,:,jk) * wgtiau(it) > salfixmin ) 
903                     tsa(:,:,jk,jp_sal) = tsa(:,:,jk,jp_sal) + s_bkginc(:,:,jk) * zincwgt
904                  END WHERE
905               ELSE
906                  tsa(:,:,jk,jp_sal) = tsa(:,:,jk,jp_sal) + s_bkginc(:,:,jk) * zincwgt
907               ENDIF
908            END DO
909
910         ENDIF
911
912         IF ( kt == nitiaufin_r + 1  ) THEN   ! For bias crcn to work
913            DEALLOCATE( t_bkginc )
914            DEALLOCATE( s_bkginc )
915         ENDIF
916
917
918      ELSEIF ( ln_asmdin ) THEN
919
920         !--------------------------------------------------------------------
921         ! Direct Initialization
922         !--------------------------------------------------------------------
923           
924         IF ( kt == nitdin_r ) THEN
925
926            neuler = 0  ! Force Euler forward step
927
928            ! Initialize the now fields with the background + increment
929            IF (ln_temnofreeze) THEN
930               ! Do not apply negative increments if the temperature will fall below freezing
931               WHERE( t_bkginc(:,:,:) > 0.0_wp .OR. tsn(:,:,:,jp_tem) + t_bkginc(:,:,:) > fzptnz(:,:,:) ) 
932                  tsn(:,:,:,jp_tem) = t_bkg(:,:,:) + t_bkginc(:,:,:)   
933               END WHERE
934            ELSE
935               tsn(:,:,:,jp_tem) = t_bkg(:,:,:) + t_bkginc(:,:,:)   
936            ENDIF
937            IF (ln_salfix) THEN
938               ! Do not apply negative increments if the salinity will fall below a specified
939               ! minimum value salfixmin
940               WHERE( s_bkginc(:,:,:) > 0.0_wp .OR. tsn(:,:,:,jp_sal) + s_bkginc(:,:,:) > salfixmin ) 
941                  tsn(:,:,:,jp_sal) = s_bkg(:,:,:) + s_bkginc(:,:,:)   
942               END WHERE
943            ELSE
944               tsn(:,:,:,jp_sal) = s_bkg(:,:,:) + s_bkginc(:,:,:)   
945            ENDIF
946
947            tsb(:,:,:,:) = tsn(:,:,:,:)                 ! Update before fields
948
949            CALL eos( tsb, rhd, rhop, gdept_0(:,:,:) )  ! Before potential and in situ densities
950!!gm  fabien
951!            CALL eos( tsb, rhd, rhop )                ! Before potential and in situ densities
952!!gm
953
954
955            IF( ln_zps .AND. .NOT. lk_c1d .AND. .NOT. ln_isfcav)      &
956               &  CALL zps_hde    ( kt, jpts, tsb, gtsu, gtsv,        &  ! Partial steps: before horizontal gradient
957               &                              rhd, gru , grv          )  ! of t, s, rd at the last ocean level
958            IF( ln_zps .AND. .NOT. lk_c1d .AND.       ln_isfcav)      &
959               &  CALL zps_hde_isf( nit000, jpts, tsb, gtsu, gtsv,    &    ! Partial steps for top cell (ISF)
960               &                                  rhd, gru , grv , aru , arv , gzu , gzv , ge3ru , ge3rv ,   &
961               &                           gtui, gtvi, grui, grvi, arui, arvi, gzui, gzvi, ge3rui, ge3rvi    ) ! of t, s, rd at the last ocean level
962
963#if defined key_zdfkpp
964            CALL eos( tsn, rhd, fsdept_n(:,:,:) )                      ! Compute rhd
965!!gm fabien            CALL eos( tsn, rhd )                      ! Compute rhd
966#endif
967
968            DEALLOCATE( t_bkginc )
969            DEALLOCATE( s_bkginc )
970            DEALLOCATE( t_bkg    )
971            DEALLOCATE( s_bkg    )
972         ENDIF
973         
974      ENDIF
975      ! Perhaps the following call should be in step
976      IF   ( ln_seaiceinc  )   CALL seaice_asm_inc ( kt )   ! apply sea ice concentration increment
977      IF   ( ln_sitinc  )      CALL sit_asm_inc ( kt )      ! apply sea ice thickness increment
978      !
979   END SUBROUTINE tra_asm_inc
980
981
982   SUBROUTINE dyn_asm_inc( kt )
983      !!----------------------------------------------------------------------
984      !!                    ***  ROUTINE dyn_asm_inc  ***
985      !!         
986      !! ** Purpose : Apply the dynamics (u and v) assimilation increments.
987      !!
988      !! ** Method  : Direct initialization or Incremental Analysis Updating.
989      !!
990      !! ** Action  :
991      !!----------------------------------------------------------------------
992      INTEGER, INTENT(IN) :: kt   ! Current time step
993      !
994      INTEGER :: jk
995      INTEGER :: it
996      REAL(wp) :: zincwgt  ! IAU weight for current time step
997      !!----------------------------------------------------------------------
998
999      IF ( ln_asmiau ) THEN
1000
1001         !--------------------------------------------------------------------
1002         ! Incremental Analysis Updating
1003         !--------------------------------------------------------------------
1004
1005         IF ( ( kt >= nitiaustr_r ).AND.( kt <= nitiaufin_r ) ) THEN
1006
1007            it = kt - nit000 + 1
1008            zincwgt = wgtiau(it) / rdt   ! IAU weight for the current time step
1009
1010            IF(lwp) THEN
1011               WRITE(numout,*) 
1012               WRITE(numout,*) 'dyn_asm_inc : Dynamics IAU at time step = ', &
1013                  &  kt,' with IAU weight = ', wgtiau(it)
1014               WRITE(numout,*) '~~~~~~~~~~~~'
1015            ENDIF
1016
1017            ! Update the dynamic tendencies
1018            DO jk = 1, jpkm1
1019               ua(:,:,jk) = ua(:,:,jk) + u_bkginc(:,:,jk) * zincwgt
1020               va(:,:,jk) = va(:,:,jk) + v_bkginc(:,:,jk) * zincwgt
1021            END DO
1022           
1023            IF ( kt == nitiaufin_r ) THEN
1024               DEALLOCATE( u_bkginc )
1025               DEALLOCATE( v_bkginc )
1026            ENDIF
1027
1028         ENDIF
1029
1030      ELSEIF ( ln_asmdin ) THEN 
1031
1032         !--------------------------------------------------------------------
1033         ! Direct Initialization
1034         !--------------------------------------------------------------------
1035         
1036         IF ( kt == nitdin_r ) THEN
1037
1038            neuler = 0                    ! Force Euler forward step
1039
1040            ! Initialize the now fields with the background + increment
1041            un(:,:,:) = u_bkg(:,:,:) + u_bkginc(:,:,:)
1042            vn(:,:,:) = v_bkg(:,:,:) + v_bkginc(:,:,:) 
1043
1044            ub(:,:,:) = un(:,:,:)         ! Update before fields
1045            vb(:,:,:) = vn(:,:,:)
1046 
1047            DEALLOCATE( u_bkg    )
1048            DEALLOCATE( v_bkg    )
1049            DEALLOCATE( u_bkginc )
1050            DEALLOCATE( v_bkginc )
1051         ENDIF
1052         !
1053      ENDIF
1054      !
1055   END SUBROUTINE dyn_asm_inc
1056
1057
1058   SUBROUTINE ssh_asm_inc( kt )
1059      !!----------------------------------------------------------------------
1060      !!                    ***  ROUTINE ssh_asm_inc  ***
1061      !!         
1062      !! ** Purpose : Apply the sea surface height assimilation increment.
1063      !!
1064      !! ** Method  : Direct initialization or Incremental Analysis Updating.
1065      !!
1066      !! ** Action  :
1067      !!----------------------------------------------------------------------
1068      INTEGER, INTENT(IN) :: kt   ! Current time step
1069      !
1070      INTEGER :: it
1071      INTEGER :: jk
1072      REAL(wp) :: zincwgt  ! IAU weight for current time step
1073      !!----------------------------------------------------------------------
1074
1075      IF ( ln_asmiau ) THEN
1076
1077         !--------------------------------------------------------------------
1078         ! Incremental Analysis Updating
1079         !--------------------------------------------------------------------
1080
1081         IF ( ( kt >= nitiaustr_r ).AND.( kt <= nitiaufin_r ) ) THEN
1082
1083            it = kt - nit000 + 1
1084            zincwgt = wgtiau(it) / rdt   ! IAU weight for the current time step
1085
1086            IF(lwp) THEN
1087               WRITE(numout,*) 
1088               WRITE(numout,*) 'ssh_asm_inc : SSH IAU at time step = ', &
1089                  &  kt,' with IAU weight = ', wgtiau(it)
1090               WRITE(numout,*) '~~~~~~~~~~~~'
1091            ENDIF
1092
1093            ! Save the tendency associated with the IAU weighted SSH increment
1094            ! (applied in dynspg.*)
1095#if defined key_asminc
1096            ssh_iau(:,:) = ssh_bkginc(:,:) * zincwgt
1097#endif
1098            IF ( kt == nitiaufin_r ) THEN
1099               DEALLOCATE( ssh_bkginc )
1100            ENDIF
1101
1102#if defined key_asminc
1103         ELSE
1104            ssh_iau(:,:) = 0._wp
1105#endif
1106         ENDIF
1107
1108      ELSEIF ( ln_asmdin ) THEN
1109
1110         !--------------------------------------------------------------------
1111         ! Direct Initialization
1112         !--------------------------------------------------------------------
1113
1114         IF ( kt == nitdin_r ) THEN
1115
1116            neuler = 0                    ! Force Euler forward step
1117
1118            ! Initialize the now fields the background + increment
1119            sshn(:,:) = ssh_bkg(:,:) + ssh_bkginc(:,:) 
1120
1121            ! Update before fields
1122            sshb(:,:) = sshn(:,:)         
1123
1124            IF( lk_vvl ) THEN
1125               DO jk = 1, jpk
1126                  fse3t_b(:,:,jk) = fse3t_n(:,:,jk)
1127               END DO
1128            ENDIF
1129
1130            DEALLOCATE( ssh_bkg    )
1131            DEALLOCATE( ssh_bkginc )
1132
1133         ENDIF
1134         !
1135      ENDIF
1136      !
1137   END SUBROUTINE ssh_asm_inc
1138
1139   SUBROUTINE sit_asm_inc( kt, kindic )
1140      !!----------------------------------------------------------------------
1141      !!                    ***  ROUTINE sit_asm_inc  ***
1142      !!         
1143      !! ** Purpose : Apply the sea ice thickness assimilation increment.
1144      !!
1145      !! ** Method  : Direct initialization or Incremental Analysis Updating.
1146      !!
1147      !! ** Action  :
1148      !!
1149      !!----------------------------------------------------------------------
1150      IMPLICIT NONE
1151      !
1152      INTEGER, INTENT(in)           ::   kt   ! Current time step
1153      INTEGER, INTENT(in), OPTIONAL ::   kindic   ! flag for disabling the deallocation
1154      !
1155      INTEGER  ::   it
1156      REAL(wp) ::   zincwgt   ! IAU weight for current time step
1157      !!----------------------------------------------------------------------
1158
1159      IF ( ln_asmiau ) THEN
1160
1161         !--------------------------------------------------------------------
1162         ! Incremental Analysis Updating
1163         !--------------------------------------------------------------------
1164
1165         IF ( ( kt >= nitiaustr_r ).AND.( kt <= nitiaufin_r ) ) THEN
1166
1167            it = kt - nit000 + 1
1168            zincwgt = wgtiau(it)      ! IAU weight for the current time step
1169            ! note this is not a tendency so should not be divided by rdt (as with the tracer and other increments)
1170            ! EF: Actually CICE is expecting a tendency so is divided by rdt below
1171
1172            IF(lwp) THEN
1173               WRITE(numout,*) 
1174               WRITE(numout,*) 'sit_asm_inc : sea ice thick IAU at time step = ', &
1175                  &  kt,' with IAU weight = ', wgtiau(it)
1176               WRITE(numout,*) '~~~~~~~~~~~~'
1177            ENDIF
1178
1179#if defined key_cice && defined key_asminc
1180            ! Sea-ice thickness : CICE case. Pass ice thickness increment tendency into CICE
1181            ndsit_da(:,:) = sit_bkginc(:,:) * zincwgt / rdt
1182#endif
1183
1184            IF ( kt == nitiaufin_r ) THEN
1185               DEALLOCATE( sit_bkginc )
1186            ENDIF
1187
1188         ELSE
1189
1190#if defined key_cice && defined key_asminc
1191            ! Sea-ice thickness : CICE case. Zero ice increment tendency into CICE
1192            ndsit_da(:,:) = 0.0_wp
1193#endif
1194
1195         ENDIF
1196
1197      ELSEIF ( ln_asmdin ) THEN
1198
1199         !--------------------------------------------------------------------
1200         ! Direct Initialization
1201         !--------------------------------------------------------------------
1202
1203         IF ( kt == nitdin_r ) THEN
1204
1205            neuler = 0                    ! Force Euler forward step
1206
1207#if defined key_cice && defined key_asminc
1208            ! Sea-ice thickness : CICE case. Pass ice thickness increment tendency into CICE
1209           ndsit_da(:,:) = sit_bkginc(:,:) / rdt
1210#endif
1211           IF ( .NOT. PRESENT(kindic) ) THEN
1212              DEALLOCATE( sit_bkginc )
1213           END IF
1214
1215         ELSE
1216
1217#if defined key_cice && defined key_asminc
1218            ! Sea-ice thicnkness : CICE case. Zero ice thickness increment tendency into CICE
1219            ndsit_da(:,:) = 0.0_wp
1220#endif
1221         
1222         ENDIF
1223
1224      ENDIF
1225
1226   END SUBROUTINE sit_asm_inc
1227
1228   SUBROUTINE seaice_asm_inc( kt, kindic )
1229      !!----------------------------------------------------------------------
1230      !!                    ***  ROUTINE seaice_asm_inc  ***
1231      !!         
1232      !! ** Purpose : Apply the sea ice assimilation increment.
1233      !!
1234      !! ** Method  : Direct initialization or Incremental Analysis Updating.
1235      !!
1236      !! ** Action  :
1237      !!
1238      !!----------------------------------------------------------------------
1239      IMPLICIT NONE
1240      !
1241      INTEGER, INTENT(in)           ::   kt   ! Current time step
1242      INTEGER, INTENT(in), OPTIONAL ::   kindic   ! flag for disabling the deallocation
1243      !
1244      INTEGER  ::   it
1245      REAL(wp) ::   zincwgt   ! IAU weight for current time step
1246#if defined key_lim2
1247      REAL(wp), DIMENSION(jpi,jpj) ::   zofrld, zohicif, zseaicendg, zhicifinc  ! LIM
1248      REAL(wp) ::   zhicifmin = 0.5_wp      ! ice minimum depth in metres
1249#endif
1250      !!----------------------------------------------------------------------
1251
1252      IF ( ln_asmiau ) THEN
1253
1254         !--------------------------------------------------------------------
1255         ! Incremental Analysis Updating
1256         !--------------------------------------------------------------------
1257
1258         IF ( ( kt >= nitiaustr_r ).AND.( kt <= nitiaufin_r ) ) THEN
1259
1260            it = kt - nit000 + 1
1261            zincwgt = wgtiau(it)      ! IAU weight for the current time step
1262            ! note this is not a tendency so should not be divided by rdt (as with the tracer and other increments)
1263
1264            IF(lwp) THEN
1265               WRITE(numout,*) 
1266               WRITE(numout,*) 'seaice_asm_inc : sea ice conc IAU at time step = ', &
1267                  &  kt,' with IAU weight = ', wgtiau(it)
1268               WRITE(numout,*) '~~~~~~~~~~~~'
1269            ENDIF
1270
1271            ! Sea-ice : LIM-3 case (to add)
1272
1273#if defined key_lim2
1274            ! Sea-ice : LIM-2 case
1275            zofrld (:,:) = frld(:,:)
1276            zohicif(:,:) = hicif(:,:)
1277            !
1278            frld  = MIN( MAX( frld (:,:) - seaice_bkginc(:,:) * zincwgt, 0.0_wp), 1.0_wp)
1279            pfrld = MIN( MAX( pfrld(:,:) - seaice_bkginc(:,:) * zincwgt, 0.0_wp), 1.0_wp)
1280            fr_i(:,:) = 1.0_wp - frld(:,:)        ! adjust ice fraction
1281            !
1282            zseaicendg(:,:) = zofrld(:,:) - frld(:,:)   ! find out actual sea ice nudge applied
1283            !
1284            ! Nudge sea ice depth to bring it up to a required minimum depth
1285            WHERE( zseaicendg(:,:) > 0.0_wp .AND. hicif(:,:) < zhicifmin ) 
1286               zhicifinc(:,:) = (zhicifmin - hicif(:,:)) * zincwgt   
1287            ELSEWHERE
1288               zhicifinc(:,:) = 0.0_wp
1289            END WHERE
1290            !
1291            ! nudge ice depth
1292            hicif (:,:) = hicif (:,:) + zhicifinc(:,:)
1293            phicif(:,:) = phicif(:,:) + zhicifinc(:,:)       
1294            !
1295            ! seaice salinity balancing (to add)
1296#endif
1297
1298#if defined key_cice && defined key_asminc
1299            ! Sea-ice : CICE case. Pass ice increment tendency into CICE
1300            ndaice_da(:,:) = seaice_bkginc(:,:) * zincwgt / rdt
1301#endif
1302
1303            IF ( kt == nitiaufin_r ) THEN
1304               DEALLOCATE( seaice_bkginc )
1305            ENDIF
1306
1307         ELSE
1308
1309#if defined key_cice && defined key_asminc
1310            ! Sea-ice : CICE case. Zero ice increment tendency into CICE
1311            ndaice_da(:,:) = 0.0_wp
1312#endif
1313
1314         ENDIF
1315
1316      ELSEIF ( ln_asmdin ) THEN
1317
1318         !--------------------------------------------------------------------
1319         ! Direct Initialization
1320         !--------------------------------------------------------------------
1321
1322         IF ( kt == nitdin_r ) THEN
1323
1324            neuler = 0                    ! Force Euler forward step
1325
1326            ! Sea-ice : LIM-3 case (to add)
1327
1328#if defined key_lim2
1329            ! Sea-ice : LIM-2 case.
1330            zofrld(:,:)=frld(:,:)
1331            zohicif(:,:)=hicif(:,:)
1332            !
1333            ! Initialize the now fields the background + increment
1334            frld (:,:) = MIN( MAX( frld(:,:) - seaice_bkginc(:,:), 0.0_wp), 1.0_wp)
1335            pfrld(:,:) = frld(:,:) 
1336            fr_i (:,:) = 1.0_wp - frld(:,:)                ! adjust ice fraction
1337            zseaicendg(:,:) = zofrld(:,:) - frld(:,:)      ! find out actual sea ice nudge applied
1338            !
1339            ! Nudge sea ice depth to bring it up to a required minimum depth
1340            WHERE( zseaicendg(:,:) > 0.0_wp .AND. hicif(:,:) < zhicifmin ) 
1341               zhicifinc(:,:) = (zhicifmin - hicif(:,:)) * zincwgt   
1342            ELSEWHERE
1343               zhicifinc(:,:) = 0.0_wp
1344            END WHERE
1345            !
1346            ! nudge ice depth
1347            hicif (:,:) = hicif (:,:) + zhicifinc(:,:)
1348            phicif(:,:) = phicif(:,:)       
1349            !
1350            ! seaice salinity balancing (to add)
1351#endif
1352 
1353#if defined key_cice && defined key_asminc
1354            ! Sea-ice : CICE case. Pass ice increment tendency into CICE
1355           ndaice_da(:,:) = seaice_bkginc(:,:) / rdt
1356#endif
1357           IF ( .NOT. PRESENT(kindic) ) THEN
1358              DEALLOCATE( seaice_bkginc )
1359           END IF
1360
1361         ELSE
1362
1363#if defined key_cice && defined key_asminc
1364            ! Sea-ice : CICE case. Zero ice increment tendency into CICE
1365            ndaice_da(:,:) = 0.0_wp
1366#endif
1367         
1368         ENDIF
1369
1370!#if defined defined key_lim2 || defined key_cice
1371!
1372!            IF (ln_seaicebal ) THEN       
1373!             !! balancing salinity increments
1374!             !! simple case from limflx.F90 (doesn't include a mass flux)
1375!             !! assumption is that as ice concentration is reduced or increased
1376!             !! the snow and ice depths remain constant
1377!             !! note that snow is being created where ice concentration is being increased
1378!             !! - could be more sophisticated and
1379!             !! not do this (but would need to alter h_snow)
1380!
1381!             usave(:,:,:)=sb(:,:,:)   ! use array as a temporary store
1382!
1383!             DO jj = 1, jpj
1384!               DO ji = 1, jpi
1385!           ! calculate change in ice and snow mass per unit area
1386!           ! positive values imply adding salt to the ocean (results from ice formation)
1387!           ! fwf : ice formation and melting
1388!
1389!                 zfons = ( -nfresh_da(ji,jj)*soce + nfsalt_da(ji,jj) )*rdt
1390!
1391!           ! change salinity down to mixed layer depth
1392!                 mld=hmld_kara(ji,jj)
1393!
1394!           ! prevent small mld
1395!           ! less than 10m can cause salinity instability
1396!                 IF (mld < 10) mld=10
1397!
1398!           ! set to bottom of a level
1399!                 DO jk = jpk-1, 2, -1
1400!                   IF ((mld > gdepw(ji,jj,jk)) .and. (mld < gdepw(ji,jj,jk+1))) THEN
1401!                     mld=gdepw(ji,jj,jk+1)
1402!                     jkmax=jk
1403!                   ENDIF
1404!                 ENDDO
1405!
1406!            ! avoid applying salinity balancing in shallow water or on land
1407!            !
1408!
1409!            ! dsal_ocn (psu kg m^-2) / (kg m^-3 * m)
1410!
1411!                 dsal_ocn=0.0_wp
1412!                 sal_thresh=5.0_wp        ! minimum salinity threshold for salinity balancing
1413!
1414!                 if (tmask(ji,jj,1) > 0 .AND. tmask(ji,jj,jkmax) > 0 ) &
1415!                              dsal_ocn = zfons / (rhop(ji,jj,1) * mld)
1416!
1417!           ! put increments in for levels in the mixed layer
1418!           ! but prevent salinity below a threshold value
1419!
1420!                   DO jk = 1, jkmax             
1421!
1422!                     IF (dsal_ocn > 0.0_wp .or. sb(ji,jj,jk)+dsal_ocn > sal_thresh) THEN
1423!                           sb(ji,jj,jk) = sb(ji,jj,jk) + dsal_ocn
1424!                           sn(ji,jj,jk) = sn(ji,jj,jk) + dsal_ocn
1425!                     ENDIF
1426!
1427!                   ENDDO
1428!
1429!      !            !  salt exchanges at the ice/ocean interface
1430!      !            zpmess         = zfons / rdt_ice    ! rdt_ice is ice timestep
1431!      !
1432!      !! Adjust fsalt. A +ve fsalt means adding salt to ocean
1433!      !!           fsalt(ji,jj) =  fsalt(ji,jj) + zpmess     ! adjust fsalt 
1434!      !!               
1435!      !!           emps(ji,jj) = emps(ji,jj) + zpmess        ! or adjust emps (see icestp1d)
1436!      !!                                                     ! E-P (kg m-2 s-2)
1437!      !            emp(ji,jj) = emp(ji,jj) + zpmess          ! E-P (kg m-2 s-2)
1438!               ENDDO !ji
1439!             ENDDO !jj!
1440!
1441!            ENDIF !ln_seaicebal
1442!
1443!#endif
1444
1445      ENDIF
1446
1447   END SUBROUTINE seaice_asm_inc
1448
1449
1450   SUBROUTINE bgc_asm_inc( kt )
1451      !!----------------------------------------------------------------------
1452      !!                    ***  ROUTINE bgc_asm_inc  ***
1453      !!         
1454      !! ** Purpose : Apply the biogeochemistry assimilation increments
1455      !!
1456      !! ** Method  : Call relevant routines in asmbgc
1457      !!
1458      !! ** Action  : Call relevant routines in asmbgc
1459      !!
1460      !!----------------------------------------------------------------------
1461      !!
1462      INTEGER, INTENT(in   ) :: kt        ! Current time step
1463      !
1464      INTEGER                :: icycper   ! Dimension of wgtiau
1465      !!
1466      !!----------------------------------------------------------------------
1467     
1468      icycper = SIZE( wgtiau )
1469     
1470      ! Ocean colour variables first
1471      IF ( ln_slchltotinc .OR. ln_slchldiainc .OR. ln_slchlnoninc .OR. &
1472         & ln_slchlnaninc .OR. ln_slchlpicinc .OR. ln_slchldininc .OR. &
1473         & ln_schltotinc  .OR. ln_slphytotinc .OR. ln_slphydiainc .OR. &
1474         & ln_slphynoninc ) THEN
1475         CALL phyto2d_asm_inc( kt, ln_asmdin, ln_asmiau, icycper, wgtiau )
1476      ENDIF
1477     
1478      ! Surface pCO2/fCO2 next
1479      IF ( ln_sfco2inc .OR. ln_spco2inc ) THEN
1480         CALL pco2_asm_inc( kt, ln_asmdin, ln_asmiau, icycper, wgtiau, &
1481            &               ln_trainc, t_bkginc, s_bkginc )
1482      ENDIF
1483     
1484      ! Profile pH next
1485      IF ( ln_pphinc ) THEN
1486         CALL ph_asm_inc( kt, ln_asmdin, ln_asmiau, icycper, wgtiau, &
1487            &             ln_trainc, t_bkginc, s_bkginc )
1488      ENDIF
1489     
1490      ! Then chlorophyll profiles
1491      IF ( ln_plchltotinc .OR. ln_pchltotinc  .OR. ln_plchldiainc .OR. &
1492         & ln_plchlnaninc .OR. ln_plchlpicinc .OR. ln_plchldininc ) THEN
1493         CALL phyto3d_asm_inc( kt, ln_asmdin, ln_asmiau, icycper, wgtiau )
1494      ENDIF
1495     
1496      ! Remaining bgc profile variables
1497      IF ( ln_pno3inc .OR. ln_psi4inc .OR. ln_pdicinc .OR. &
1498         & ln_palkinc .OR. ln_po2inc  .OR. ln_ppo4inc ) THEN
1499         CALL bgc3d_asm_inc( kt, ln_asmdin, ln_asmiau, icycper, wgtiau )
1500      ENDIF
1501
1502   END SUBROUTINE bgc_asm_inc
1503   
1504   !!======================================================================
1505END MODULE asminc
Note: See TracBrowser for help on using the repository browser.