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.
traatf.F90 in NEMO/trunk/src/OCE/TRA – NEMO

source: NEMO/trunk/src/OCE/TRA/traatf.F90 @ 12377

Last change on this file since 12377 was 12377, checked in by acc, 4 years ago

The big one. Merging all 2019 developments from the option 1 branch back onto the trunk.

This changeset reproduces 2019/dev_r11943_MERGE_2019 on the trunk using a 2-URL merge
onto a working copy of the trunk. I.e.:

svn merge --ignore-ancestry \

svn+ssh://acc@forge.ipsl.jussieu.fr/ipsl/forge/projets/nemo/svn/NEMO/trunk \
svn+ssh://acc@forge.ipsl.jussieu.fr/ipsl/forge/projets/nemo/svn/NEMO/branches/2019/dev_r11943_MERGE_2019 ./

The --ignore-ancestry flag avoids problems that may otherwise arise from the fact that
the merge history been trunk and branch may have been applied in a different order but
care has been taken before this step to ensure that all applicable fixes and updates
are present in the merge branch.

The trunk state just before this step has been branched to releases/release-4.0-HEAD
and that branch has been immediately tagged as releases/release-4.0.2. Any fixes
or additions in response to tickets on 4.0, 4.0.1 or 4.0.2 should be done on
releases/release-4.0-HEAD. From now on future 'point' releases (e.g. 4.0.2) will
remain unchanged with periodic releases as needs demand. Note release-4.0-HEAD is a
transitional naming convention. Future full releases, say 4.2, will have a release-4.2
branch which fulfills this role and the first point release (e.g. 4.2.0) will be made
immediately following the release branch creation.

2020 developments can be started from any trunk revision later than this one.

  • Property svn:keywords set to Id
File size: 20.2 KB
Line 
1MODULE traatf
2   !!======================================================================
3   !!                       ***  MODULE  traatf  ***
4   !! Ocean active tracers:  Asselin time filtering for temperature and salinity
5   !!======================================================================
6   !! History :  OPA  !  1991-11  (G. Madec)  Original code
7   !!            7.0  !  1993-03  (M. Guyon)  symetrical conditions
8   !!            8.0  !  1996-02  (G. Madec & M. Imbard)  opa release 8.0
9   !!             -   !  1996-04  (A. Weaver)  Euler forward step
10   !!            8.2  !  1999-02  (G. Madec, N. Grima)  semi-implicit pressure grad.
11   !!  NEMO      1.0  !  2002-08  (G. Madec)  F90: Free form and module
12   !!             -   !  2002-11  (C. Talandier, A-M Treguier) Open boundaries
13   !!             -   !  2005-04  (C. Deltel) Add Asselin trend in the ML budget
14   !!            2.0  !  2006-02  (L. Debreu, C. Mazauric) Agrif implementation
15   !!            3.0  !  2008-06  (G. Madec)  time stepping always done in trazdf
16   !!            3.1  !  2009-02  (G. Madec, R. Benshila)  re-introduce the vvl option
17   !!            3.3  !  2010-04  (M. Leclair, G. Madec)  semi-implicit hpg with asselin filter + modified LF-RA
18   !!             -   !  2010-05  (C. Ethe, G. Madec)  merge TRC-TRA
19   !!            4.1  !  2019-08  (A. Coward, D. Storkey) rename tranxt.F90 -> traatf.F90. Now only does time filtering.
20   !!----------------------------------------------------------------------
21
22   !!----------------------------------------------------------------------
23   !!   tra_atf       : time filtering on tracers
24   !!   tra_atf_fix   : time filtering on tracers : fixed    volume case
25   !!   tra_atf_vvl   : time filtering on tracers : variable volume case
26   !!----------------------------------------------------------------------
27   USE oce             ! ocean dynamics and tracers variables
28   USE dom_oce         ! ocean space and time domain variables
29   USE sbc_oce         ! surface boundary condition: ocean
30   USE sbcrnf          ! river runoffs
31   USE isf_oce         ! ice shelf melting
32   USE zdf_oce         ! ocean vertical mixing
33   USE domvvl          ! variable volume
34   USE trd_oce         ! trends: ocean variables
35   USE trdtra          ! trends manager: tracers
36   USE traqsr          ! penetrative solar radiation (needed for nksr)
37   USE phycst          ! physical constant
38   USE ldftra          ! lateral physics : tracers
39   USE ldfslp          ! lateral physics : slopes
40   USE bdy_oce  , ONLY : ln_bdy
41   USE bdytra          ! open boundary condition (bdy_tra routine)
42   !
43   USE in_out_manager  ! I/O manager
44   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
45   USE prtctl          ! Print control
46   USE timing          ! Timing
47#if defined key_agrif
48   USE agrif_oce_interp
49#endif
50
51   IMPLICIT NONE
52   PRIVATE
53
54   PUBLIC   tra_atf       ! routine called by step.F90
55   PUBLIC   tra_atf_fix   ! to be used in trcnxt
56   PUBLIC   tra_atf_vvl   ! to be used in trcnxt
57
58   !! * Substitutions
59#  include "do_loop_substitute.h90"
60   !!----------------------------------------------------------------------
61   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
62   !! $Id$
63   !! Software governed by the CeCILL license (see ./LICENSE)
64   !!----------------------------------------------------------------------
65CONTAINS
66
67   SUBROUTINE tra_atf( kt, Kbb, Kmm, Kaa, pts )
68      !!----------------------------------------------------------------------
69      !!                   ***  ROUTINE traatf  ***
70      !!
71      !! ** Purpose :   Apply the boundary condition on the after temperature 
72      !!             and salinity fields and add the Asselin time filter on now fields.
73      !!
74      !! ** Method  :   At this stage of the computation, ta and sa are the
75      !!             after temperature and salinity as the time stepping has
76      !!             been performed in trazdf_imp or trazdf_exp module.
77      !!
78      !!              - Apply lateral boundary conditions on (ta,sa)
79      !!             at the local domain   boundaries through lbc_lnk call,
80      !!             at the one-way open boundaries (ln_bdy=T),
81      !!             at the AGRIF zoom   boundaries (lk_agrif=T)
82      !!
83      !!              - Update lateral boundary conditions on AGRIF children
84      !!             domains (lk_agrif=T)
85      !!
86      !! ** Action  : - ts(Kmm) time filtered
87      !!----------------------------------------------------------------------
88      INTEGER                                  , INTENT(in   ) :: kt             ! ocean time-step index
89      INTEGER                                  , INTENT(in   ) :: Kbb, Kmm, Kaa  ! time level indices
90      REAL(wp), DIMENSION(jpi,jpj,jpk,jpts,jpt), INTENT(inout) :: pts            ! active tracers
91      !!
92      INTEGER  ::   ji, jj, jk, jn   ! dummy loop indices
93      REAL(wp) ::   zfact            ! local scalars
94      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) ::   ztrdt, ztrds
95      !!----------------------------------------------------------------------
96      !
97      IF( ln_timing )   CALL timing_start( 'tra_atf')
98      !
99      IF( kt == nit000 ) THEN
100         IF(lwp) WRITE(numout,*)
101         IF(lwp) WRITE(numout,*) 'tra_atf : apply Asselin time filter to "now" fields'
102         IF(lwp) WRITE(numout,*) '~~~~~~~'
103      ENDIF
104
105      ! Update after tracer on domain lateral boundaries
106      !
107#if defined key_agrif
108      CALL Agrif_tra                     ! AGRIF zoom boundaries
109#endif
110      !                                              ! local domain boundaries  (T-point, unchanged sign)
111      CALL lbc_lnk_multi( 'traatf', pts(:,:,:,jp_tem,Kaa), 'T', 1., pts(:,:,:,jp_sal,Kaa), 'T', 1. )
112      !
113      IF( ln_bdy )   CALL bdy_tra( kt, Kbb, pts, Kaa )  ! BDY open boundaries
114 
115      ! set time step size (Euler/Leapfrog)
116      IF( neuler == 0 .AND. kt == nit000 ) THEN   ;   r2dt =        rdt   ! at nit000             (Euler)
117      ELSEIF( kt <= nit000 + 1 )           THEN   ;   r2dt = 2._wp* rdt   ! at nit000 or nit000+1 (Leapfrog)
118      ENDIF
119
120      ! trends computation initialisation
121      IF( l_trdtra )   THEN                   
122         ALLOCATE( ztrdt(jpi,jpj,jpk) , ztrds(jpi,jpj,jpk) )
123         ztrdt(:,:,jpk) = 0._wp
124         ztrds(:,:,jpk) = 0._wp
125         IF( ln_traldf_iso ) THEN              ! diagnose the "pure" Kz diffusive trend
126            CALL trd_tra( kt, Kmm, Kaa, 'TRA', jp_tem, jptra_zdfp, ztrdt )
127            CALL trd_tra( kt, Kmm, Kaa, 'TRA', jp_sal, jptra_zdfp, ztrds )
128         ENDIF
129         ! total trend for the non-time-filtered variables.
130         zfact = 1.0 / rdt
131         ! G Nurser 23 Mar 2017. Recalculate trend as Delta(e3t*T)/e3tn; e3tn cancel from pts(Kmm) terms
132         DO jk = 1, jpkm1
133            ztrdt(:,:,jk) = ( pts(:,:,jk,jp_tem,Kaa)*e3t(:,:,jk,Kaa) / e3t(:,:,jk,Kmm) - pts(:,:,jk,jp_tem,Kmm)) * zfact
134            ztrds(:,:,jk) = ( pts(:,:,jk,jp_sal,Kaa)*e3t(:,:,jk,Kaa) / e3t(:,:,jk,Kmm) - pts(:,:,jk,jp_sal,Kmm)) * zfact
135         END DO
136         CALL trd_tra( kt, Kmm, Kaa, 'TRA', jp_tem, jptra_tot, ztrdt )
137         CALL trd_tra( kt, Kmm, Kaa, 'TRA', jp_sal, jptra_tot, ztrds )
138         IF( ln_linssh ) THEN       ! linear sea surface height only
139            ! Store now fields before applying the Asselin filter
140            ! in order to calculate Asselin filter trend later.
141            ztrdt(:,:,:) = pts(:,:,:,jp_tem,Kmm) 
142            ztrds(:,:,:) = pts(:,:,:,jp_sal,Kmm)
143         ENDIF
144      ENDIF
145
146      IF( neuler == 0 .AND. kt == nit000 ) THEN       ! Euler time-stepping
147         !
148         IF (l_trdtra .AND. .NOT. ln_linssh ) THEN   ! Zero Asselin filter contribution must be explicitly written out since for vvl
149            !                                        ! Asselin filter is output by tra_atf_vvl that is not called on this time step
150            ztrdt(:,:,:) = 0._wp
151            ztrds(:,:,:) = 0._wp
152            CALL trd_tra( kt, Kmm, Kaa, 'TRA', jp_tem, jptra_atf, ztrdt )
153            CALL trd_tra( kt, Kmm, Kaa, 'TRA', jp_sal, jptra_atf, ztrds )
154         END IF
155         !
156      ELSE                                            ! Leap-Frog + Asselin filter time stepping
157         !
158         IF( ln_linssh ) THEN   ;   CALL tra_atf_fix( kt, Kbb, Kmm, Kaa, nit000,      'TRA', pts, jpts )  ! linear free surface
159         ELSE                   ;   CALL tra_atf_vvl( kt, Kbb, Kmm, Kaa, nit000, rdt, 'TRA', pts, sbc_tsc, sbc_tsc_b, jpts )  ! non-linear free surface
160         ENDIF
161         !
162         CALL lbc_lnk_multi( 'traatf', pts(:,:,:,jp_tem,Kbb) , 'T', 1., pts(:,:,:,jp_sal,Kbb) , 'T', 1., &
163                  &                    pts(:,:,:,jp_tem,Kmm) , 'T', 1., pts(:,:,:,jp_sal,Kmm) , 'T', 1., &
164                  &                    pts(:,:,:,jp_tem,Kaa), 'T', 1., pts(:,:,:,jp_sal,Kaa), 'T', 1.  )
165         !
166      ENDIF     
167      !
168      IF( l_trdtra .AND. ln_linssh ) THEN      ! trend of the Asselin filter (tb filtered - tb)/dt     
169         zfact = 1._wp / r2dt             
170         DO jk = 1, jpkm1
171            ztrdt(:,:,jk) = ( pts(:,:,jk,jp_tem,Kmm) - ztrdt(:,:,jk) ) * zfact
172            ztrds(:,:,jk) = ( pts(:,:,jk,jp_sal,Kmm) - ztrds(:,:,jk) ) * zfact
173         END DO
174         CALL trd_tra( kt, Kmm, Kaa, 'TRA', jp_tem, jptra_atf, ztrdt )
175         CALL trd_tra( kt, Kmm, Kaa, 'TRA', jp_sal, jptra_atf, ztrds )
176      END IF
177      IF( l_trdtra )   DEALLOCATE( ztrdt , ztrds )
178      !
179      !                        ! control print
180      IF(sn_cfctl%l_prtctl)   CALL prt_ctl( tab3d_1=pts(:,:,:,jp_tem,Kmm), clinfo1=' nxt  - Tn: ', mask1=tmask,   &
181         &                                  tab3d_2=pts(:,:,:,jp_sal,Kmm), clinfo2=       ' Sn: ', mask2=tmask )
182      !
183      IF( ln_timing )   CALL timing_stop('tra_atf')
184      !
185   END SUBROUTINE tra_atf
186
187
188   SUBROUTINE tra_atf_fix( kt, Kbb, Kmm, Kaa, kit000, cdtype, pt, kjpt )
189      !!----------------------------------------------------------------------
190      !!                   ***  ROUTINE tra_atf_fix  ***
191      !!
192      !! ** Purpose :   fixed volume: apply the Asselin time filter to the "now" field
193      !!
194      !! ** Method  : - Apply a Asselin time filter on now fields.
195      !!
196      !! ** Action  : - pt(Kmm) ready for the next time step
197      !!----------------------------------------------------------------------
198      INTEGER                                  , INTENT(in   ) ::  kt            ! ocean time-step index
199      INTEGER                                  , INTENT(in   ) ::  Kbb, Kmm, Kaa ! time level indices
200      INTEGER                                  , INTENT(in   ) ::  kit000        ! first time step index
201      CHARACTER(len=3)                         , INTENT(in   ) ::  cdtype        ! =TRA or TRC (tracer indicator)
202      INTEGER                                  , INTENT(in   ) ::  kjpt          ! number of tracers
203      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt,jpt), INTENT(inout) ::  pt            ! tracer fields
204      !
205      INTEGER  ::   ji, jj, jk, jn   ! dummy loop indices
206      REAL(wp) ::   ztn, ztd         ! local scalars
207      !!----------------------------------------------------------------------
208      !
209      IF( kt == kit000 )  THEN
210         IF(lwp) WRITE(numout,*)
211         IF(lwp) WRITE(numout,*) 'tra_atf_fix : time filtering', cdtype
212         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
213      ENDIF
214      !
215      DO jn = 1, kjpt
216         !
217         DO_3D_00_00( 1, jpkm1 )
218            ztn = pt(ji,jj,jk,jn,Kmm)                                   
219            ztd = pt(ji,jj,jk,jn,Kaa) - 2._wp * ztn + pt(ji,jj,jk,jn,Kbb)  ! time laplacian on tracers
220            !
221            pt(ji,jj,jk,jn,Kmm) = ztn + atfp * ztd                      ! pt <-- filtered pt
222         END_3D
223         !
224      END DO
225      !
226   END SUBROUTINE tra_atf_fix
227
228
229   SUBROUTINE tra_atf_vvl( kt, Kbb, Kmm, Kaa, kit000, p2dt, cdtype, pt, psbc_tc, psbc_tc_b, kjpt )
230      !!----------------------------------------------------------------------
231      !!                   ***  ROUTINE tra_atf_vvl  ***
232      !!
233      !! ** Purpose :   Time varying volume: apply the Asselin time filter 
234      !!
235      !! ** Method  : - Apply a thickness weighted Asselin time filter on now fields.
236      !!             pt(Kmm)  = ( e3t(Kmm)*pt(Kmm) + atfp*[ e3t(Kbb)*pt(Kbb) - 2 e3t(Kmm)*pt(Kmm) + e3t_a*pt(Kaa) ] )
237      !!                       /( e3t(Kmm)         + atfp*[ e3t(Kbb)         - 2 e3t(Kmm)         + e3t(Kaa)    ] )
238      !!
239      !! ** Action  : - pt(Kmm) ready for the next time step
240      !!----------------------------------------------------------------------
241      INTEGER                                  , INTENT(in   ) ::  kt        ! ocean time-step index
242      INTEGER                                  , INTENT(in   ) ::  Kbb, Kmm, Kaa ! time level indices
243      INTEGER                                  , INTENT(in   ) ::  kit000    ! first time step index
244      REAL(wp)                                 , INTENT(in   ) ::  p2dt      ! time-step
245      CHARACTER(len=3)                         , INTENT(in   ) ::  cdtype    ! =TRA or TRC (tracer indicator)
246      INTEGER                                  , INTENT(in   ) ::  kjpt      ! number of tracers
247      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt,jpt), INTENT(inout) ::  pt        ! tracer fields
248      REAL(wp), DIMENSION(jpi,jpj    ,kjpt)    , INTENT(in   ) ::  psbc_tc   ! surface tracer content
249      REAL(wp), DIMENSION(jpi,jpj    ,kjpt)    , INTENT(in   ) ::  psbc_tc_b ! before surface tracer content
250      !
251      LOGICAL  ::   ll_traqsr, ll_rnf, ll_isf   ! local logical
252      INTEGER  ::   ji, jj, jk, jn              ! dummy loop indices
253      REAL(wp) ::   zfact, zfact1, ztc_a , ztc_n , ztc_b , ztc_f , ztc_d    ! local scalar
254      REAL(wp) ::   zfact2, ze3t_b, ze3t_n, ze3t_a, ze3t_f, ze3t_d, zscale  !   -      -
255      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:) ::   ztrd_atf
256      !!----------------------------------------------------------------------
257      !
258      IF( kt == kit000 )  THEN
259         IF(lwp) WRITE(numout,*)
260         IF(lwp) WRITE(numout,*) 'tra_atf_vvl : time filtering', cdtype
261         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
262      ENDIF
263      !
264      IF( cdtype == 'TRA' )  THEN   
265         ll_traqsr  = ln_traqsr        ! active  tracers case  and  solar penetration
266         ll_rnf     = ln_rnf           ! active  tracers case  and  river runoffs
267         ll_isf     = ln_isf           ! active  tracers case  and  ice shelf melting
268      ELSE                          ! passive tracers case
269         ll_traqsr  = .FALSE.          ! NO solar penetration
270         ll_rnf     = .FALSE.          ! NO river runoffs ????          !!gm BUG ? 
271         ll_isf     = .FALSE.          ! NO ice shelf melting/freezing  !!gm BUG ??
272      ENDIF
273      !
274      IF( ( l_trdtra .AND. cdtype == 'TRA' ) .OR. ( l_trdtrc .AND. cdtype == 'TRC' ) )   THEN
275         ALLOCATE( ztrd_atf(jpi,jpj,jpk,kjpt) )
276         ztrd_atf(:,:,:,:) = 0.0_wp
277      ENDIF
278      zfact = 1._wp / p2dt
279      zfact1 = atfp * p2dt
280      zfact2 = zfact1 * r1_rau0
281      DO jn = 1, kjpt     
282         DO_3D_00_00( 1, jpkm1 )
283            ze3t_b = e3t(ji,jj,jk,Kbb)
284            ze3t_n = e3t(ji,jj,jk,Kmm)
285            ze3t_a = e3t(ji,jj,jk,Kaa)
286            !                                         ! tracer content at Before, now and after
287            ztc_b  = pt(ji,jj,jk,jn,Kbb) * ze3t_b
288            ztc_n  = pt(ji,jj,jk,jn,Kmm) * ze3t_n
289            ztc_a  = pt(ji,jj,jk,jn,Kaa) * ze3t_a
290            !
291            ze3t_d = ze3t_a - 2. * ze3t_n + ze3t_b
292            ztc_d  = ztc_a  - 2. * ztc_n  + ztc_b
293            !
294            ze3t_f = ze3t_n + atfp * ze3t_d
295            ztc_f  = ztc_n  + atfp * ztc_d
296            !
297            ! Add asselin correction on scale factors:
298            zscale = tmask(ji,jj,jk) * e3t(ji,jj,jk,Kmm) / ( ht(ji,jj) + 1._wp - ssmask(ji,jj) ) 
299            ze3t_f = ze3t_f - zfact2 * zscale * ( emp_b(ji,jj) - emp(ji,jj) ) 
300            IF ( ll_rnf ) ze3t_f = ze3t_f + zfact2 * zscale * (    rnf_b(ji,jj) -    rnf(ji,jj) ) 
301            IF ( ll_isf ) THEN
302               IF ( ln_isfcav_mlt ) ze3t_f = ze3t_f - zfact2 * zscale * ( fwfisf_cav_b(ji,jj) - fwfisf_cav(ji,jj) )
303               IF ( ln_isfpar_mlt ) ze3t_f = ze3t_f - zfact2 * zscale * ( fwfisf_par_b(ji,jj) - fwfisf_par(ji,jj) )
304            ENDIF
305            !
306            IF( jk == mikt(ji,jj) ) THEN           ! first level
307               ztc_f  = ztc_f  - zfact1 * ( psbc_tc(ji,jj,jn) - psbc_tc_b(ji,jj,jn) )
308            ENDIF
309            !
310            ! solar penetration (temperature only)
311            IF( ll_traqsr .AND. jn == jp_tem .AND. jk <= nksr )                            & 
312               &     ztc_f  = ztc_f  - zfact1 * ( qsr_hc(ji,jj,jk) - qsr_hc_b(ji,jj,jk) ) 
313               !
314            !
315            IF( ll_rnf .AND. jk <= nk_rnf(ji,jj) )                                          &
316               &     ztc_f  = ztc_f  - zfact1 * ( rnf_tsc(ji,jj,jn) - rnf_tsc_b(ji,jj,jn) ) & 
317               &                              * e3t(ji,jj,jk,Kmm) / h_rnf(ji,jj)
318
319            !
320            ! ice shelf
321            IF( ll_isf ) THEN
322               !
323               ! melt in the cavity
324               IF ( ln_isfcav_mlt ) THEN
325                  ! level fully include in the Losch_2008 ice shelf boundary layer
326                  IF ( jk >= misfkt_cav(ji,jj) .AND. jk < misfkb_cav(ji,jj) ) THEN
327                     ztc_f  = ztc_f  - zfact1 * ( risf_cav_tsc(ji,jj,jn) - risf_cav_tsc_b(ji,jj,jn) ) &
328                        &                     * e3t(ji,jj,jk,Kmm) / rhisf_tbl_cav(ji,jj)
329                  END IF
330                  ! level partially include in Losch_2008 ice shelf boundary layer
331                  IF ( jk == misfkb_cav(ji,jj) ) THEN
332                     ztc_f  = ztc_f  - zfact1 * ( risf_cav_tsc(ji,jj,jn) - risf_cav_tsc_b(ji,jj,jn) )  &
333                            &                 * e3t(ji,jj,jk,Kmm) / rhisf_tbl_cav(ji,jj) * rfrac_tbl_cav(ji,jj)
334                  END IF
335               END IF
336               !
337               ! parametrised melt (cavity closed)
338               IF ( ln_isfpar_mlt ) THEN
339                  ! level fully include in the Losch_2008 ice shelf boundary layer
340                  IF ( jk >= misfkt_par(ji,jj) .AND. jk < misfkb_par(ji,jj) ) THEN
341                     ztc_f  = ztc_f  - zfact1 * ( risf_par_tsc(ji,jj,jn) - risf_par_tsc_b(ji,jj,jn) )  &
342                            &                 * e3t(ji,jj,jk,Kmm) / rhisf_tbl_par(ji,jj)
343                  END IF
344                  ! level partially include in Losch_2008 ice shelf boundary layer
345                  IF ( jk == misfkb_par(ji,jj) ) THEN
346                     ztc_f  = ztc_f  - zfact1 * ( risf_par_tsc(ji,jj,jn) - risf_par_tsc_b(ji,jj,jn) )  &
347                            &                 * e3t(ji,jj,jk,Kmm) / rhisf_tbl_par(ji,jj) * rfrac_tbl_par(ji,jj)
348                  END IF
349               END IF
350               !
351               ! ice sheet coupling correction
352               IF ( ln_isfcpl ) THEN
353                  !
354                  ! at kt = nit000,  risfcpl_vol_n = 0 and risfcpl_vol_b = risfcpl_vol so contribution nul
355                  IF ( ln_rstart .AND. kt == nit000+1 ) THEN
356                     ztc_f  = ztc_f  + zfact1 * risfcpl_tsc(ji,jj,jk,jn) * r1_e1e2t(ji,jj)
357                     ! Shouldn't volume increment be spread according thanks to zscale  ?
358                     ze3t_f = ze3t_f - zfact1 * risfcpl_vol(ji,jj,jk   ) * r1_e1e2t(ji,jj)
359                  END IF
360                  !
361               END IF
362               !
363            END IF
364            !
365            ze3t_f = 1.e0 / ze3t_f
366            pt(ji,jj,jk,jn,Kmm) = ztc_f * ze3t_f    ! time filtered "now" field
367            !
368            IF( ( l_trdtra .and. cdtype == 'TRA' ) .OR. ( l_trdtrc .and. cdtype == 'TRC' ) ) THEN
369               ztrd_atf(ji,jj,jk,jn) = (ztc_f - ztc_n) * zfact/ze3t_n
370            ENDIF
371            !
372         END_3D
373         !
374      END DO
375      !
376      IF( ( l_trdtra .AND. cdtype == 'TRA' ) .OR. ( l_trdtrc .AND. cdtype == 'TRC' ) )   THEN
377         IF( l_trdtra .AND. cdtype == 'TRA' ) THEN
378            CALL trd_tra( kt, Kmm, Kaa, cdtype, jp_tem, jptra_atf, ztrd_atf(:,:,:,jp_tem) )
379            CALL trd_tra( kt, Kmm, Kaa, cdtype, jp_sal, jptra_atf, ztrd_atf(:,:,:,jp_sal) )
380         ENDIF
381         IF( l_trdtrc .AND. cdtype == 'TRC' ) THEN
382            DO jn = 1, kjpt
383               CALL trd_tra( kt, Kmm, Kaa, cdtype, jn, jptra_atf, ztrd_atf(:,:,:,jn) )
384            END DO
385         ENDIF
386         DEALLOCATE( ztrd_atf )
387      ENDIF
388      !
389   END SUBROUTINE tra_atf_vvl
390
391   !!======================================================================
392END MODULE traatf
Note: See TracBrowser for help on using the repository browser.