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.
tranxt.F90 in NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/TRA – NEMO

source: NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/TRA/tranxt.F90 @ 10969

Last change on this file since 10969 was 10957, checked in by davestorkey, 5 years ago

2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps : BDY module (including some removal of redundant code). Tested in AMM12 and ORCA1 (not full SETTE test at this stage).

  • Property svn:keywords set to Id
File size: 20.6 KB
Line 
1MODULE tranxt
2   !!======================================================================
3   !!                       ***  MODULE  tranxt  ***
4   !! Ocean active tracers:  time stepping on 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   !!----------------------------------------------------------------------
20
21   !!----------------------------------------------------------------------
22   !!   tra_nxt       : time stepping on tracers
23   !!   tra_nxt_fix   : time stepping on tracers : fixed    volume case
24   !!   tra_nxt_vvl   : time stepping on tracers : variable volume case
25   !!----------------------------------------------------------------------
26   USE oce             ! ocean dynamics and tracers variables
27   USE dom_oce         ! ocean space and time domain variables
28   USE sbc_oce         ! surface boundary condition: ocean
29   USE sbcrnf          ! river runoffs
30   USE sbcisf          ! ice shelf melting
31   USE zdf_oce         ! ocean vertical mixing
32   USE domvvl          ! variable volume
33   USE trd_oce         ! trends: ocean variables
34   USE trdtra          ! trends manager: tracers
35   USE traqsr          ! penetrative solar radiation (needed for nksr)
36   USE phycst          ! physical constant
37   USE ldftra          ! lateral physics : tracers
38   USE ldfslp          ! lateral physics : slopes
39   USE bdy_oce  , ONLY : ln_bdy
40   USE bdytra          ! open boundary condition (bdy_tra routine)
41   !
42   USE in_out_manager  ! I/O manager
43   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
44   USE prtctl          ! Print control
45   USE timing          ! Timing
46#if defined key_agrif
47   USE agrif_oce_interp
48#endif
49
50   IMPLICIT NONE
51   PRIVATE
52
53   PUBLIC   tra_nxt       ! routine called by step.F90
54   PUBLIC   tra_nxt_fix   ! to be used in trcnxt
55   PUBLIC   tra_nxt_vvl   ! to be used in trcnxt
56
57   !! * Substitutions
58#  include "vectopt_loop_substitute.h90"
59   !!----------------------------------------------------------------------
60   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
61   !! $Id$
62   !! Software governed by the CeCILL license (see ./LICENSE)
63   !!----------------------------------------------------------------------
64CONTAINS
65
66   SUBROUTINE tra_nxt( kt, Kbb, Kmm, Krhs, Kaa )
67      !!----------------------------------------------------------------------
68      !!                   ***  ROUTINE tranxt  ***
69      !!
70      !! ** Purpose :   Apply the boundary condition on the after temperature 
71      !!             and salinity fields, achieved the time stepping by adding
72      !!             the Asselin filter on now fields and swapping the 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(Kbb) & ts(Kmm) ready for the next time step
87      !!----------------------------------------------------------------------
88      INTEGER, INTENT(in) ::   kt                   ! ocean time-step index
89      INTEGER, INTENT(in) ::   Kbb, Kmm, Krhs, Kaa  ! time level indices
90      !!
91      INTEGER  ::   ji, jj, jk, jn   ! dummy loop indices
92      REAL(wp) ::   zfact            ! local scalars
93      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) ::   ztrdt, ztrds
94      !!----------------------------------------------------------------------
95      !
96      IF( ln_timing )   CALL timing_start( 'tra_nxt')
97      !
98      IF( kt == nit000 ) THEN
99         IF(lwp) WRITE(numout,*)
100         IF(lwp) WRITE(numout,*) 'tra_nxt : achieve the time stepping by Asselin filter and array swap'
101         IF(lwp) WRITE(numout,*) '~~~~~~~'
102      ENDIF
103
104      ! Update after tracer on domain lateral boundaries
105      !
106#if defined key_agrif
107      CALL Agrif_tra                     ! AGRIF zoom boundaries
108#endif
109      !                                              ! local domain boundaries  (T-point, unchanged sign)
110      CALL lbc_lnk_multi( 'tranxt', ts(:,:,:,jp_tem,Krhs), 'T', 1., ts(:,:,:,jp_sal,Krhs), 'T', 1. )
111      !
112      !! IMMERSE development : Following the general pattern for the new code we want to pass in the
113      !!                       velocities to bdy_dyn as arguments so here we use "ts" even
114      !!                       though we haven't converted the tracer names in the rest of tranxt.F90
115      !!                       because it will be completely rewritten. DS.
116      IF( ln_bdy )   CALL bdy_tra( kt, Kbb, ts, Kaa )  ! BDY open boundaries
117 
118      ! set time step size (Euler/Leapfrog)
119      IF( neuler == 0 .AND. kt == nit000 ) THEN   ;   r2dt =        rdt   ! at nit000             (Euler)
120      ELSEIF( kt <= nit000 + 1 )           THEN   ;   r2dt = 2._wp* rdt   ! at nit000 or nit000+1 (Leapfrog)
121      ENDIF
122
123      ! trends computation initialisation
124      IF( l_trdtra )   THEN                   
125         ALLOCATE( ztrdt(jpi,jpj,jpk) , ztrds(jpi,jpj,jpk) )
126         ztrdt(:,:,jpk) = 0._wp
127         ztrds(:,:,jpk) = 0._wp
128         IF( ln_traldf_iso ) THEN              ! diagnose the "pure" Kz diffusive trend
129            CALL trd_tra( kt, Kmm, Krhs, 'TRA', jp_tem, jptra_zdfp, ztrdt )
130            CALL trd_tra( kt, Kmm, Krhs, 'TRA', jp_sal, jptra_zdfp, ztrds )
131         ENDIF
132         ! total trend for the non-time-filtered variables.
133         zfact = 1.0 / rdt
134         ! G Nurser 23 Mar 2017. Recalculate trend as Delta(e3t*T)/e3tn; e3tn cancel from ts(Kmm) terms
135         DO jk = 1, jpkm1
136            ztrdt(:,:,jk) = ( ts(:,:,jk,jp_tem,Krhs)*e3t(:,:,jk,Krhs) / e3t(:,:,jk,Kmm) - ts(:,:,jk,jp_tem,Kmm)) * zfact
137            ztrds(:,:,jk) = ( ts(:,:,jk,jp_sal,Krhs)*e3t(:,:,jk,Krhs) / e3t(:,:,jk,Kmm) - ts(:,:,jk,jp_sal,Kmm)) * zfact
138         END DO
139         CALL trd_tra( kt, Kmm, Krhs, 'TRA', jp_tem, jptra_tot, ztrdt )
140         CALL trd_tra( kt, Kmm, Krhs, 'TRA', jp_sal, jptra_tot, ztrds )
141         IF( ln_linssh ) THEN       ! linear sea surface height only
142            ! Store now fields before applying the Asselin filter
143            ! in order to calculate Asselin filter trend later.
144            ztrdt(:,:,:) = ts(:,:,:,jp_tem,Kmm) 
145            ztrds(:,:,:) = ts(:,:,:,jp_sal,Kmm)
146         ENDIF
147      ENDIF
148
149      IF( neuler == 0 .AND. kt == nit000 ) THEN       ! Euler time-stepping at first time-step (only swap)
150         DO jn = 1, jpts
151            DO jk = 1, jpkm1
152               ts(:,:,jk,jn,Kmm) = ts(:,:,jk,jn,Krhs)   
153            END DO
154         END DO
155         IF (l_trdtra .AND. .NOT. ln_linssh ) THEN   ! Zero Asselin filter contribution must be explicitly written out since for vvl
156            !                                        ! Asselin filter is output by tra_nxt_vvl that is not called on this time step
157            ztrdt(:,:,:) = 0._wp
158            ztrds(:,:,:) = 0._wp
159            CALL trd_tra( kt, Kmm, Krhs, 'TRA', jp_tem, jptra_atf, ztrdt )
160            CALL trd_tra( kt, Kmm, Krhs, 'TRA', jp_sal, jptra_atf, ztrds )
161         END IF
162         !
163      ELSE                                            ! Leap-Frog + Asselin filter time stepping
164         !
165         IF( ln_linssh ) THEN   ;   CALL tra_nxt_fix( kt,      Kmm,       nit000,      'TRA',                           &
166           &                                                   ts(:,:,:,:,Kbb), ts(:,:,:,:,Kmm), ts(:,:,:,:,Krhs), jpts )  ! linear free surface
167         ELSE                   ;   CALL tra_nxt_vvl( kt, Kbb, Kmm, Krhs, nit000, rdt, 'TRA',                           &
168           &                                                   ts(:,:,:,:,Kbb), ts(:,:,:,:,Kmm), ts(:,:,:,:,Krhs),      &
169           &                                                   sbc_tsc        , sbc_tsc_b                        , jpts )  ! non-linear free surface
170         ENDIF
171         !
172         CALL lbc_lnk_multi( 'tranxt', ts(:,:,:,jp_tem,Kbb) , 'T', 1., ts(:,:,:,jp_sal,Kbb) , 'T', 1., &
173                  &                    ts(:,:,:,jp_tem,Kmm) , 'T', 1., ts(:,:,:,jp_sal,Kmm) , 'T', 1., &
174                  &                    ts(:,:,:,jp_tem,Krhs), 'T', 1., ts(:,:,:,jp_sal,Krhs), 'T', 1.  )
175         !
176      ENDIF     
177      !
178      IF( l_trdtra .AND. ln_linssh ) THEN      ! trend of the Asselin filter (tb filtered - tb)/dt     
179         zfact = 1._wp / r2dt             
180         DO jk = 1, jpkm1
181            ztrdt(:,:,jk) = ( ts(:,:,jk,jp_tem,Kbb) - ztrdt(:,:,jk) ) * zfact
182            ztrds(:,:,jk) = ( ts(:,:,jk,jp_sal,Kbb) - ztrds(:,:,jk) ) * zfact
183         END DO
184         CALL trd_tra( kt, Kmm, Krhs, 'TRA', jp_tem, jptra_atf, ztrdt )
185         CALL trd_tra( kt, Kmm, Krhs, 'TRA', jp_sal, jptra_atf, ztrds )
186      END IF
187      IF( l_trdtra )   DEALLOCATE( ztrdt , ztrds )
188      !
189      !                        ! control print
190      IF(ln_ctl)   CALL prt_ctl( tab3d_1=ts(:,:,:,jp_tem,Kmm), clinfo1=' nxt  - Tn: ', mask1=tmask,   &
191         &                       tab3d_2=ts(:,:,:,jp_sal,Kmm), clinfo2=       ' Sn: ', mask2=tmask )
192      !
193      IF( ln_timing )   CALL timing_stop('tra_nxt')
194      !
195   END SUBROUTINE tra_nxt
196
197
198   SUBROUTINE tra_nxt_fix( kt, Kmm, kit000, cdtype, ptb, ptn, pta, kjpt )
199      !!----------------------------------------------------------------------
200      !!                   ***  ROUTINE tra_nxt_fix  ***
201      !!
202      !! ** Purpose :   fixed volume: apply the Asselin time filter and
203      !!                swap the tracer fields.
204      !!
205      !! ** Method  : - Apply a Asselin time filter on now fields.
206      !!              - swap tracer fields to prepare the next time_step.
207      !!
208      !! ** Action  : - ptb & ptn ready for the next time step
209      !!----------------------------------------------------------------------
210      INTEGER                              , INTENT(in   ) ::  kt        ! ocean time-step index
211      INTEGER                              , INTENT(in   ) ::  Kmm       ! time level index
212      INTEGER                              , INTENT(in   ) ::  kit000    ! first time step index
213      CHARACTER(len=3)                     , INTENT(in   ) ::  cdtype    ! =TRA or TRC (tracer indicator)
214      INTEGER                              , INTENT(in   ) ::  kjpt      ! number of tracers
215      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::  ptb       ! before tracer fields
216      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::  ptn       ! now tracer fields
217      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::  pta       ! tracer trend
218      !
219      INTEGER  ::   ji, jj, jk, jn   ! dummy loop indices
220      REAL(wp) ::   ztn, ztd         ! local scalars
221      !!----------------------------------------------------------------------
222      !
223      IF( kt == kit000 )  THEN
224         IF(lwp) WRITE(numout,*)
225         IF(lwp) WRITE(numout,*) 'tra_nxt_fix : time stepping', cdtype
226         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
227      ENDIF
228      !
229      DO jn = 1, kjpt
230         !
231         DO jk = 1, jpkm1
232            DO jj = 2, jpjm1
233               DO ji = fs_2, fs_jpim1
234                  ztn = ptn(ji,jj,jk,jn)                                   
235                  ztd = pta(ji,jj,jk,jn) - 2._wp * ztn + ptb(ji,jj,jk,jn)  ! time laplacian on tracers
236                  !
237                  ptb(ji,jj,jk,jn) = ztn + atfp * ztd                      ! ptb <-- filtered ptn
238                  ptn(ji,jj,jk,jn) = pta(ji,jj,jk,jn)                      ! ptn <-- pta
239               END DO
240           END DO
241         END DO
242         !
243      END DO
244      !
245   END SUBROUTINE tra_nxt_fix
246
247
248   SUBROUTINE tra_nxt_vvl( kt, Kbb, Kmm, Krhs, kit000, p2dt, cdtype, ptb, ptn, pta, psbc_tc, psbc_tc_b, kjpt )
249      !!----------------------------------------------------------------------
250      !!                   ***  ROUTINE tra_nxt_vvl  ***
251      !!
252      !! ** Purpose :   Time varying volume: apply the Asselin time filter 
253      !!                and swap the tracer fields.
254      !!
255      !! ** Method  : - Apply a thickness weighted Asselin time filter on now fields.
256      !!              - swap tracer fields to prepare the next time_step.
257      !!             tb  = ( e3t(Kmm)*tn + atfp*[ e3t(Kbb)*tb - 2 e3t(Kmm)*tn + e3t_a*ta ] )
258      !!                  /( e3t(Kmm)    + atfp*[ e3t(Kbb)    - 2 e3t(Kmm)    + e3t(Krhs)    ] )
259      !!             tn  = ta
260      !!
261      !! ** Action  : - ptb & ptn ready for the next time step
262      !!----------------------------------------------------------------------
263      INTEGER                              , INTENT(in   ) ::  kt        ! ocean time-step index
264      INTEGER                              , INTENT(in   ) ::  Kbb, Kmm, Krhs ! time level indices
265      INTEGER                              , INTENT(in   ) ::  kit000    ! first time step index
266      REAL(wp)                             , INTENT(in   ) ::  p2dt      ! time-step
267      CHARACTER(len=3)                     , INTENT(in   ) ::  cdtype    ! =TRA or TRC (tracer indicator)
268      INTEGER                              , INTENT(in   ) ::  kjpt      ! number of tracers
269      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::  ptb       ! before tracer fields
270      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::  ptn       ! now tracer fields
271      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::  pta       ! tracer trend
272      REAL(wp), DIMENSION(jpi,jpj    ,kjpt), INTENT(in   ) ::  psbc_tc   ! surface tracer content
273      REAL(wp), DIMENSION(jpi,jpj    ,kjpt), INTENT(in   ) ::  psbc_tc_b ! before surface tracer content
274      !
275      LOGICAL  ::   ll_traqsr, ll_rnf, ll_isf   ! local logical
276      INTEGER  ::   ji, jj, jk, jn              ! dummy loop indices
277      REAL(wp) ::   zfact, zfact1, ztc_a , ztc_n , ztc_b , ztc_f , ztc_d    ! local scalar
278      REAL(wp) ::   zfact2, ze3t_b, ze3t_n, ze3t_a, ze3t_f, ze3t_d   !   -      -
279      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:) ::   ztrd_atf
280      !!----------------------------------------------------------------------
281      !
282      IF( kt == kit000 )  THEN
283         IF(lwp) WRITE(numout,*)
284         IF(lwp) WRITE(numout,*) 'tra_nxt_vvl : time stepping', cdtype
285         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
286      ENDIF
287      !
288      IF( cdtype == 'TRA' )  THEN   
289         ll_traqsr  = ln_traqsr        ! active  tracers case  and  solar penetration
290         ll_rnf     = ln_rnf           ! active  tracers case  and  river runoffs
291         ll_isf     = ln_isf           ! active  tracers case  and  ice shelf melting
292      ELSE                          ! passive tracers case
293         ll_traqsr  = .FALSE.          ! NO solar penetration
294         ll_rnf     = .FALSE.          ! NO river runoffs ????          !!gm BUG ? 
295         ll_isf     = .FALSE.          ! NO ice shelf melting/freezing  !!gm BUG ??
296      ENDIF
297      !
298      IF( ( l_trdtra .AND. cdtype == 'TRA' ) .OR. ( l_trdtrc .AND. cdtype == 'TRC' ) )   THEN
299         ALLOCATE( ztrd_atf(jpi,jpj,jpk,kjpt) )
300         ztrd_atf(:,:,:,:) = 0.0_wp
301      ENDIF
302      zfact = 1._wp / p2dt
303      zfact1 = atfp * p2dt
304      zfact2 = zfact1 * r1_rau0
305      DO jn = 1, kjpt     
306         DO jk = 1, jpkm1
307            DO jj = 2, jpjm1
308               DO ji = fs_2, fs_jpim1
309                  ze3t_b = e3t(ji,jj,jk,Kbb)
310                  ze3t_n = e3t(ji,jj,jk,Kmm)
311                  ze3t_a = e3t(ji,jj,jk,Krhs)
312                  !                                         ! tracer content at Before, now and after
313                  ztc_b  = ptb(ji,jj,jk,jn) * ze3t_b
314                  ztc_n  = ptn(ji,jj,jk,jn) * ze3t_n
315                  ztc_a  = pta(ji,jj,jk,jn) * ze3t_a
316                  !
317                  ze3t_d = ze3t_a - 2. * ze3t_n + ze3t_b
318                  ztc_d  = ztc_a  - 2. * ztc_n  + ztc_b
319                  !
320                  ze3t_f = ze3t_n + atfp * ze3t_d
321                  ztc_f  = ztc_n  + atfp * ztc_d
322                  !
323                  IF( jk == mikt(ji,jj) ) THEN           ! first level
324                     ze3t_f = ze3t_f - zfact2 * ( (emp_b(ji,jj)    - emp(ji,jj)   )  &
325                            &                   + (fwfisf_b(ji,jj) - fwfisf(ji,jj))  )
326                     ztc_f  = ztc_f  - zfact1 * ( psbc_tc(ji,jj,jn) - psbc_tc_b(ji,jj,jn) )
327                  ENDIF
328                  IF( ln_rnf_depth ) THEN
329                     ! Rivers are not just at the surface must go down to nk_rnf(ji,jj)
330                     IF( mikt(ji,jj) <=jk .and. jk <= nk_rnf(ji,jj)  ) THEN
331                        ze3t_f = ze3t_f - zfact2 * ( - (rnf_b(ji,jj) - rnf(ji,jj)   )  ) &
332                    &                            * ( e3t(ji,jj,jk,Kmm) / h_rnf(ji,jj) ) 
333                     ENDIF
334                  ELSE
335                     IF( jk == mikt(ji,jj) ) THEN           ! first level
336                        ze3t_f = ze3t_f - zfact2 * ( - (rnf_b(ji,jj)    - rnf(ji,jj)   ) ) 
337                     ENDIF
338                  ENDIF
339
340                  !
341                  ! solar penetration (temperature only)
342                  IF( ll_traqsr .AND. jn == jp_tem .AND. jk <= nksr )                            & 
343                     &     ztc_f  = ztc_f  - zfact1 * ( qsr_hc(ji,jj,jk) - qsr_hc_b(ji,jj,jk) ) 
344                     !
345                  ! river runoff
346                  IF( ll_rnf .AND. jk <= nk_rnf(ji,jj) )                                          &
347                     &     ztc_f  = ztc_f  - zfact1 * ( rnf_tsc(ji,jj,jn) - rnf_tsc_b(ji,jj,jn) ) & 
348                     &                              * e3t(ji,jj,jk,Kmm) / h_rnf(ji,jj)
349                     !
350                  ! ice shelf
351                  IF( ll_isf ) THEN
352                     ! level fully include in the Losch_2008 ice shelf boundary layer
353                     IF ( jk >= misfkt(ji,jj) .AND. jk < misfkb(ji,jj) )                          &
354                        ztc_f  = ztc_f  - zfact1 * ( risf_tsc(ji,jj,jn) - risf_tsc_b(ji,jj,jn) )  &
355                               &                 * e3t(ji,jj,jk,Kmm) * r1_hisf_tbl (ji,jj)
356                     ! level partially include in Losch_2008 ice shelf boundary layer
357                     IF ( jk == misfkb(ji,jj) )                                                   &
358                        ztc_f  = ztc_f  - zfact1 * ( risf_tsc(ji,jj,jn) - risf_tsc_b(ji,jj,jn) )  &
359                               &                 * e3t(ji,jj,jk,Kmm) * r1_hisf_tbl (ji,jj) * ralpha(ji,jj)
360                  END IF
361                  !
362                  ze3t_f = 1.e0 / ze3t_f
363                  ptb(ji,jj,jk,jn) = ztc_f * ze3t_f       ! ptb <-- ptn filtered
364                  ptn(ji,jj,jk,jn) = pta(ji,jj,jk,jn)     ! ptn <-- pta
365                  !
366                  IF( ( l_trdtra .and. cdtype == 'TRA' ) .OR. ( l_trdtrc .and. cdtype == 'TRC' ) ) THEN
367                     ztrd_atf(ji,jj,jk,jn) = (ztc_f - ztc_n) * zfact/ze3t_n
368                  ENDIF
369                  !
370               END DO
371            END DO
372         END DO
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, Krhs, cdtype, jp_tem, jptra_atf, ztrd_atf(:,:,:,jp_tem) )
379            CALL trd_tra( kt, Kmm, Krhs, 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, Krhs, cdtype, jn, jptra_atf, ztrd_atf(:,:,:,jn) )
384            END DO
385         ENDIF
386         DEALLOCATE( ztrd_atf )
387      ENDIF
388      !
389   END SUBROUTINE tra_nxt_vvl
390
391   !!======================================================================
392END MODULE tranxt
Note: See TracBrowser for help on using the repository browser.