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

Last change on this file was 14794, checked in by jwhile, 3 years ago

Updates for Gulf18, see UTILS ticket 442

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