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

source: branches/UKMO/dev5518_gulf18_changes/NEMOGCM/NEMO/OPA_SRC/ASM/asminc.F90 @ 14034

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

Updates to asminc

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