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

source: branches/UKMO/AMM15_v3_6_STABLE_package_collate_amm7ps45/NEMOGCM/NEMO/OPA_SRC/ASM/asminc.F90 @ 15670

Last change on this file since 15670 was 15670, checked in by petesykes, 2 years ago

Adding PS45 AMM7 changes

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