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/UKMO/NEMO_4.0.4_geo-adj_ovf/src/OCE/TRA – NEMO

source: NEMO/branches/UKMO/NEMO_4.0.4_geo-adj_ovf/src/OCE/TRA/tranxt.F90 @ 15559

Last change on this file since 15559 was 15559, checked in by dbruciaferri, 3 years ago

introducing key_geo_adj

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