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 @ 10946

Last change on this file since 10946 was 10946, checked in by acc, 5 years ago

2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps : Convert STO, TRD and USR modules and all knock on effects of these conversions. Note change to USR module may have implications for the TEST CASES (not tested yet). Standard SETTE tested only

  • Property svn:keywords set to Id
File size: 19.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, Kmm, Krhs )
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  : - tsb & tsn ready for the next time step
87      !!----------------------------------------------------------------------
88      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
89      INTEGER, INTENT(in) ::   Kmm, Krhs  ! 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', tsa(:,:,:,jp_tem), 'T', 1., tsa(:,:,:,jp_sal), 'T', 1. )
111      !
112      IF( ln_bdy )   CALL bdy_tra( kt )  ! BDY open boundaries
113 
114      ! set time step size (Euler/Leapfrog)
115      IF( neuler == 0 .AND. kt == nit000 ) THEN   ;   r2dt =        rdt   ! at nit000             (Euler)
116      ELSEIF( kt <= nit000 + 1 )           THEN   ;   r2dt = 2._wp* rdt   ! at nit000 or nit000+1 (Leapfrog)
117      ENDIF
118
119      ! trends computation initialisation
120      IF( l_trdtra )   THEN                   
121         ALLOCATE( ztrdt(jpi,jpj,jpk) , ztrds(jpi,jpj,jpk) )
122         ztrdt(:,:,jpk) = 0._wp
123         ztrds(:,:,jpk) = 0._wp
124         IF( ln_traldf_iso ) THEN              ! diagnose the "pure" Kz diffusive trend
125            CALL trd_tra( kt, Kmm, Krhs, 'TRA', jp_tem, jptra_zdfp, ztrdt )
126            CALL trd_tra( kt, Kmm, Krhs, 'TRA', jp_sal, jptra_zdfp, ztrds )
127         ENDIF
128         ! total trend for the non-time-filtered variables.
129         zfact = 1.0 / rdt
130         ! G Nurser 23 Mar 2017. Recalculate trend as Delta(e3t*T)/e3tn; e3tn cancel from tsn terms
131         DO jk = 1, jpkm1
132            ztrdt(:,:,jk) = ( tsa(:,:,jk,jp_tem)*e3t_a(:,:,jk) / e3t_n(:,:,jk) - tsn(:,:,jk,jp_tem)) * zfact
133            ztrds(:,:,jk) = ( tsa(:,:,jk,jp_sal)*e3t_a(:,:,jk) / e3t_n(:,:,jk) - tsn(:,:,jk,jp_sal)) * zfact
134         END DO
135         CALL trd_tra( kt, Kmm, Krhs, 'TRA', jp_tem, jptra_tot, ztrdt )
136         CALL trd_tra( kt, Kmm, Krhs, 'TRA', jp_sal, jptra_tot, ztrds )
137         IF( ln_linssh ) THEN       ! linear sea surface height only
138            ! Store now fields before applying the Asselin filter
139            ! in order to calculate Asselin filter trend later.
140            ztrdt(:,:,:) = tsn(:,:,:,jp_tem) 
141            ztrds(:,:,:) = tsn(:,:,:,jp_sal)
142         ENDIF
143      ENDIF
144
145      IF( neuler == 0 .AND. kt == nit000 ) THEN       ! Euler time-stepping at first time-step (only swap)
146         DO jn = 1, jpts
147            DO jk = 1, jpkm1
148               tsn(:,:,jk,jn) = tsa(:,:,jk,jn)   
149            END DO
150         END DO
151         IF (l_trdtra .AND. .NOT. ln_linssh ) THEN   ! Zero Asselin filter contribution must be explicitly written out since for vvl
152            !                                        ! Asselin filter is output by tra_nxt_vvl that is not called on this time step
153            ztrdt(:,:,:) = 0._wp
154            ztrds(:,:,:) = 0._wp
155            CALL trd_tra( kt, Kmm, Krhs, 'TRA', jp_tem, jptra_atf, ztrdt )
156            CALL trd_tra( kt, Kmm, Krhs, 'TRA', jp_sal, jptra_atf, ztrds )
157         END IF
158         !
159      ELSE                                            ! Leap-Frog + Asselin filter time stepping
160         !
161         IF( ln_linssh ) THEN   ;   CALL tra_nxt_fix( kt, nit000,      'TRA', tsb, tsn, tsa, jpts )  ! linear free surface
162         ELSE                   ;   CALL tra_nxt_vvl( kt, Kmm, Krhs, nit000, rdt, 'TRA', tsb, tsn, tsa,   &
163           &                                                                sbc_tsc, sbc_tsc_b, jpts )  ! non-linear free surface
164         ENDIF
165         !
166         CALL lbc_lnk_multi( 'tranxt', tsb(:,:,:,jp_tem), 'T', 1., tsb(:,:,:,jp_sal), 'T', 1., &
167                  &          tsn(:,:,:,jp_tem), 'T', 1., tsn(:,:,:,jp_sal), 'T', 1., &
168                  &          tsa(:,:,:,jp_tem), 'T', 1., tsa(:,:,:,jp_sal), 'T', 1.  )
169         !
170      ENDIF     
171      !
172      IF( l_trdtra .AND. ln_linssh ) THEN      ! trend of the Asselin filter (tb filtered - tb)/dt     
173         zfact = 1._wp / r2dt             
174         DO jk = 1, jpkm1
175            ztrdt(:,:,jk) = ( tsb(:,:,jk,jp_tem) - ztrdt(:,:,jk) ) * zfact
176            ztrds(:,:,jk) = ( tsb(:,:,jk,jp_sal) - ztrds(:,:,jk) ) * zfact
177         END DO
178         CALL trd_tra( kt, Kmm, Krhs, 'TRA', jp_tem, jptra_atf, ztrdt )
179         CALL trd_tra( kt, Kmm, Krhs, 'TRA', jp_sal, jptra_atf, ztrds )
180      END IF
181      IF( l_trdtra )   DEALLOCATE( ztrdt , ztrds )
182      !
183      !                        ! control print
184      IF(ln_ctl)   CALL prt_ctl( tab3d_1=tsn(:,:,:,jp_tem), clinfo1=' nxt  - Tn: ', mask1=tmask,   &
185         &                       tab3d_2=tsn(:,:,:,jp_sal), clinfo2=       ' Sn: ', mask2=tmask )
186      !
187      IF( ln_timing )   CALL timing_stop('tra_nxt')
188      !
189   END SUBROUTINE tra_nxt
190
191
192   SUBROUTINE tra_nxt_fix( kt, kit000, cdtype, ptb, ptn, pta, kjpt )
193      !!----------------------------------------------------------------------
194      !!                   ***  ROUTINE tra_nxt_fix  ***
195      !!
196      !! ** Purpose :   fixed volume: apply the Asselin time filter and
197      !!                swap the tracer fields.
198      !!
199      !! ** Method  : - Apply a Asselin time filter on now fields.
200      !!              - swap tracer fields to prepare the next time_step.
201      !!
202      !! ** Action  : - tsb & tsn ready for the next time step
203      !!----------------------------------------------------------------------
204      INTEGER                              , INTENT(in   ) ::  kt        ! ocean time-step index
205      INTEGER                              , INTENT(in   ) ::  kit000    ! first time step index
206      CHARACTER(len=3)                     , INTENT(in   ) ::  cdtype    ! =TRA or TRC (tracer indicator)
207      INTEGER                              , INTENT(in   ) ::  kjpt      ! number of tracers
208      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::  ptb       ! before tracer fields
209      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::  ptn       ! now tracer fields
210      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::  pta       ! tracer trend
211      !
212      INTEGER  ::   ji, jj, jk, jn   ! dummy loop indices
213      REAL(wp) ::   ztn, ztd         ! local scalars
214      !!----------------------------------------------------------------------
215      !
216      IF( kt == kit000 )  THEN
217         IF(lwp) WRITE(numout,*)
218         IF(lwp) WRITE(numout,*) 'tra_nxt_fix : time stepping', cdtype
219         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
220      ENDIF
221      !
222      DO jn = 1, kjpt
223         !
224         DO jk = 1, jpkm1
225            DO jj = 2, jpjm1
226               DO ji = fs_2, fs_jpim1
227                  ztn = ptn(ji,jj,jk,jn)                                   
228                  ztd = pta(ji,jj,jk,jn) - 2._wp * ztn + ptb(ji,jj,jk,jn)  ! time laplacian on tracers
229                  !
230                  ptb(ji,jj,jk,jn) = ztn + atfp * ztd                      ! ptb <-- filtered ptn
231                  ptn(ji,jj,jk,jn) = pta(ji,jj,jk,jn)                      ! ptn <-- pta
232               END DO
233           END DO
234         END DO
235         !
236      END DO
237      !
238   END SUBROUTINE tra_nxt_fix
239
240
241   SUBROUTINE tra_nxt_vvl( kt, Kmm, Krhs, kit000, p2dt, cdtype, ptb, ptn, pta, psbc_tc, psbc_tc_b, kjpt )
242      !!----------------------------------------------------------------------
243      !!                   ***  ROUTINE tra_nxt_vvl  ***
244      !!
245      !! ** Purpose :   Time varying volume: apply the Asselin time filter 
246      !!                and swap the tracer fields.
247      !!
248      !! ** Method  : - Apply a thickness weighted Asselin time filter on now fields.
249      !!              - swap tracer fields to prepare the next time_step.
250      !!             tb  = ( e3t_n*tn + atfp*[ e3t_b*tb - 2 e3t_n*tn + e3t_a*ta ] )
251      !!                  /( e3t_n    + atfp*[ e3t_b    - 2 e3t_n    + e3t_a    ] )
252      !!             tn  = ta
253      !!
254      !! ** Action  : - tsb & tsn ready for the next time step
255      !!----------------------------------------------------------------------
256      INTEGER                              , INTENT(in   ) ::  kt        ! ocean time-step index
257      INTEGER                              , INTENT(in   ) ::  Kmm, Krhs ! time level indices
258      INTEGER                              , INTENT(in   ) ::  kit000    ! first time step index
259      REAL(wp)                             , INTENT(in   ) ::  p2dt      ! time-step
260      CHARACTER(len=3)                     , INTENT(in   ) ::  cdtype    ! =TRA or TRC (tracer indicator)
261      INTEGER                              , INTENT(in   ) ::  kjpt      ! number of tracers
262      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::  ptb       ! before tracer fields
263      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::  ptn       ! now tracer fields
264      REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) ::  pta       ! tracer trend
265      REAL(wp), DIMENSION(jpi,jpj    ,kjpt), INTENT(in   ) ::  psbc_tc   ! surface tracer content
266      REAL(wp), DIMENSION(jpi,jpj    ,kjpt), INTENT(in   ) ::  psbc_tc_b ! before surface tracer content
267      !
268      LOGICAL  ::   ll_traqsr, ll_rnf, ll_isf   ! local logical
269      INTEGER  ::   ji, jj, jk, jn              ! dummy loop indices
270      REAL(wp) ::   zfact, zfact1, ztc_a , ztc_n , ztc_b , ztc_f , ztc_d    ! local scalar
271      REAL(wp) ::   zfact2, ze3t_b, ze3t_n, ze3t_a, ze3t_f, ze3t_d   !   -      -
272      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:) ::   ztrd_atf
273      !!----------------------------------------------------------------------
274      !
275      IF( kt == kit000 )  THEN
276         IF(lwp) WRITE(numout,*)
277         IF(lwp) WRITE(numout,*) 'tra_nxt_vvl : time stepping', cdtype
278         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~'
279      ENDIF
280      !
281      IF( cdtype == 'TRA' )  THEN   
282         ll_traqsr  = ln_traqsr        ! active  tracers case  and  solar penetration
283         ll_rnf     = ln_rnf           ! active  tracers case  and  river runoffs
284         ll_isf     = ln_isf           ! active  tracers case  and  ice shelf melting
285      ELSE                          ! passive tracers case
286         ll_traqsr  = .FALSE.          ! NO solar penetration
287         ll_rnf     = .FALSE.          ! NO river runoffs ????          !!gm BUG ? 
288         ll_isf     = .FALSE.          ! NO ice shelf melting/freezing  !!gm BUG ??
289      ENDIF
290      !
291      IF( ( l_trdtra .AND. cdtype == 'TRA' ) .OR. ( l_trdtrc .AND. cdtype == 'TRC' ) )   THEN
292         ALLOCATE( ztrd_atf(jpi,jpj,jpk,kjpt) )
293         ztrd_atf(:,:,:,:) = 0.0_wp
294      ENDIF
295      zfact = 1._wp / p2dt
296      zfact1 = atfp * p2dt
297      zfact2 = zfact1 * r1_rau0
298      DO jn = 1, kjpt     
299         DO jk = 1, jpkm1
300            DO jj = 2, jpjm1
301               DO ji = fs_2, fs_jpim1
302                  ze3t_b = e3t_b(ji,jj,jk)
303                  ze3t_n = e3t_n(ji,jj,jk)
304                  ze3t_a = e3t_a(ji,jj,jk)
305                  !                                         ! tracer content at Before, now and after
306                  ztc_b  = ptb(ji,jj,jk,jn) * ze3t_b
307                  ztc_n  = ptn(ji,jj,jk,jn) * ze3t_n
308                  ztc_a  = pta(ji,jj,jk,jn) * ze3t_a
309                  !
310                  ze3t_d = ze3t_a - 2. * ze3t_n + ze3t_b
311                  ztc_d  = ztc_a  - 2. * ztc_n  + ztc_b
312                  !
313                  ze3t_f = ze3t_n + atfp * ze3t_d
314                  ztc_f  = ztc_n  + atfp * ztc_d
315                  !
316                  IF( jk == mikt(ji,jj) ) THEN           ! first level
317                     ze3t_f = ze3t_f - zfact2 * ( (emp_b(ji,jj)    - emp(ji,jj)   )  &
318                            &                   + (fwfisf_b(ji,jj) - fwfisf(ji,jj))  )
319                     ztc_f  = ztc_f  - zfact1 * ( psbc_tc(ji,jj,jn) - psbc_tc_b(ji,jj,jn) )
320                  ENDIF
321                  IF( ln_rnf_depth ) THEN
322                     ! Rivers are not just at the surface must go down to nk_rnf(ji,jj)
323                     IF( mikt(ji,jj) <=jk .and. jk <= nk_rnf(ji,jj)  ) THEN
324                        ze3t_f = ze3t_f - zfact2 * ( - (rnf_b(ji,jj) - rnf(ji,jj)   )  ) &
325                    &                            * ( e3t_n(ji,jj,jk) / h_rnf(ji,jj) ) 
326                     ENDIF
327                  ELSE
328                     IF( jk == mikt(ji,jj) ) THEN           ! first level
329                        ze3t_f = ze3t_f - zfact2 * ( - (rnf_b(ji,jj)    - rnf(ji,jj)   ) ) 
330                     ENDIF
331                  ENDIF
332
333                  !
334                  ! solar penetration (temperature only)
335                  IF( ll_traqsr .AND. jn == jp_tem .AND. jk <= nksr )                            & 
336                     &     ztc_f  = ztc_f  - zfact1 * ( qsr_hc(ji,jj,jk) - qsr_hc_b(ji,jj,jk) ) 
337                     !
338                  ! river runoff
339                  IF( ll_rnf .AND. jk <= nk_rnf(ji,jj) )                                          &
340                     &     ztc_f  = ztc_f  - zfact1 * ( rnf_tsc(ji,jj,jn) - rnf_tsc_b(ji,jj,jn) ) & 
341                     &                              * e3t_n(ji,jj,jk) / h_rnf(ji,jj)
342                     !
343                  ! ice shelf
344                  IF( ll_isf ) THEN
345                     ! level fully include in the Losch_2008 ice shelf boundary layer
346                     IF ( jk >= misfkt(ji,jj) .AND. jk < misfkb(ji,jj) )                          &
347                        ztc_f  = ztc_f  - zfact1 * ( risf_tsc(ji,jj,jn) - risf_tsc_b(ji,jj,jn) )  &
348                               &                 * e3t_n(ji,jj,jk) * r1_hisf_tbl (ji,jj)
349                     ! level partially include in Losch_2008 ice shelf boundary layer
350                     IF ( jk == misfkb(ji,jj) )                                                   &
351                        ztc_f  = ztc_f  - zfact1 * ( risf_tsc(ji,jj,jn) - risf_tsc_b(ji,jj,jn) )  &
352                               &                 * e3t_n(ji,jj,jk) * r1_hisf_tbl (ji,jj) * ralpha(ji,jj)
353                  END IF
354                  !
355                  ze3t_f = 1.e0 / ze3t_f
356                  ptb(ji,jj,jk,jn) = ztc_f * ze3t_f       ! ptb <-- ptn filtered
357                  ptn(ji,jj,jk,jn) = pta(ji,jj,jk,jn)     ! ptn <-- pta
358                  !
359                  IF( ( l_trdtra .and. cdtype == 'TRA' ) .OR. ( l_trdtrc .and. cdtype == 'TRC' ) ) THEN
360                     ztrd_atf(ji,jj,jk,jn) = (ztc_f - ztc_n) * zfact/ze3t_n
361                  ENDIF
362                  !
363               END DO
364            END DO
365         END DO
366         !
367      END DO
368      !
369      IF( ( l_trdtra .AND. cdtype == 'TRA' ) .OR. ( l_trdtrc .AND. cdtype == 'TRC' ) )   THEN
370         IF( l_trdtra .AND. cdtype == 'TRA' ) THEN
371            CALL trd_tra( kt, Kmm, Krhs, cdtype, jp_tem, jptra_atf, ztrd_atf(:,:,:,jp_tem) )
372            CALL trd_tra( kt, Kmm, Krhs, cdtype, jp_sal, jptra_atf, ztrd_atf(:,:,:,jp_sal) )
373         ENDIF
374         IF( l_trdtrc .AND. cdtype == 'TRC' ) THEN
375            DO jn = 1, kjpt
376               CALL trd_tra( kt, Kmm, Krhs, cdtype, jn, jptra_atf, ztrd_atf(:,:,:,jn) )
377            END DO
378         ENDIF
379         DEALLOCATE( ztrd_atf )
380      ENDIF
381      !
382   END SUBROUTINE tra_nxt_vvl
383
384   !!======================================================================
385END MODULE tranxt
Note: See TracBrowser for help on using the repository browser.