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.
dynnxt.F90 in branches/2017/dev_merge_2017/NEMOGCM/NEMO/OPA_SRC/DYN – NEMO

source: branches/2017/dev_merge_2017/NEMOGCM/NEMO/OPA_SRC/DYN/dynnxt.F90 @ 9226

Last change on this file since 9226 was 9226, checked in by gm, 6 years ago

dev_merge_2017: #2002 : dynnxt.F90 and icbtrj.F90: minor bug correction

  • Property svn:keywords set to Id
File size: 19.6 KB
RevLine 
[3]1MODULE dynnxt
[1502]2   !!=========================================================================
[3]3   !!                       ***  MODULE  dynnxt  ***
4   !! Ocean dynamics: time stepping
[1502]5   !!=========================================================================
[1438]6   !! History :  OPA  !  1987-02  (P. Andrich, D. L Hostis)  Original code
7   !!                 !  1990-10  (C. Levy, G. Madec)
8   !!            7.0  !  1993-03  (M. Guyon)  symetrical conditions
9   !!            8.0  !  1997-02  (G. Madec & M. Imbard)  opa, release 8.0
10   !!            8.2  !  1997-04  (A. Weaver)  Euler forward step
11   !!             -   !  1997-06  (G. Madec)  lateral boudary cond., lbc routine
12   !!    NEMO    1.0  !  2002-08  (G. Madec)  F90: Free form and module
13   !!             -   !  2002-10  (C. Talandier, A-M. Treguier) Open boundary cond.
14   !!            2.0  !  2005-11  (V. Garnier) Surface pressure gradient organization
15   !!            2.3  !  2007-07  (D. Storkey) Calls to BDY routines.
[1502]16   !!            3.2  !  2009-06  (G. Madec, R.Benshila)  re-introduce the vvl option
[2528]17   !!            3.3  !  2010-09  (D. Storkey, E.O'Dea) Bug fix for BDY module
[2723]18   !!            3.3  !  2011-03  (P. Oddo) Bug fix for time-splitting+(BDY-OBC) and not VVL
[4292]19   !!            3.5  !  2013-07  (J. Chanut) Compliant with time splitting changes
[6140]20   !!            3.6  !  2014-04  (G. Madec) add the diagnostic of the time filter trends
[5930]21   !!            3.7  !  2015-11  (J. Chanut) Free surface simplification
[1502]22   !!-------------------------------------------------------------------------
[1438]23 
[1502]24   !!-------------------------------------------------------------------------
[6140]25   !!   dyn_nxt       : obtain the next (after) horizontal velocity
[1502]26   !!-------------------------------------------------------------------------
[6140]27   USE oce            ! ocean dynamics and tracers
28   USE dom_oce        ! ocean space and time domain
29   USE sbc_oce        ! Surface boundary condition: ocean fields
[9023]30   USE sbcrnf         ! river runoffs
[6140]31   USE phycst         ! physical constants
32   USE dynadv         ! dynamics: vector invariant versus flux form
33   USE dynspg_ts      ! surface pressure gradient: split-explicit scheme
34   USE domvvl         ! variable volume
[7646]35   USE bdy_oce   , ONLY: ln_bdy
[6140]36   USE bdydta         ! ocean open boundary conditions
37   USE bdydyn         ! ocean open boundary conditions
38   USE bdyvol         ! ocean open boundary condition (bdy_vol routines)
39   USE trd_oce        ! trends: ocean variables
40   USE trddyn         ! trend manager: dynamics
41   USE trdken         ! trend manager: kinetic energy
[4990]42   !
[6140]43   USE in_out_manager ! I/O manager
44   USE iom            ! I/O manager library
45   USE lbclnk         ! lateral boundary condition (or mpp link)
46   USE lib_mpp        ! MPP library
47   USE prtctl         ! Print control
48   USE timing         ! Timing
[2528]49#if defined key_agrif
50   USE agrif_opa_interp
51#endif
[3]52
53   IMPLICIT NONE
54   PRIVATE
55
[1438]56   PUBLIC    dyn_nxt   ! routine called by step.F90
57
[2715]58   !!----------------------------------------------------------------------
[9019]59   !! NEMO/OPA 4.0 , NEMO Consortium (2017)
[1438]60   !! $Id$
[2715]61   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
62   !!----------------------------------------------------------------------
[3]63CONTAINS
64
65   SUBROUTINE dyn_nxt ( kt )
66      !!----------------------------------------------------------------------
67      !!                  ***  ROUTINE dyn_nxt  ***
68      !!                   
[5930]69      !! ** Purpose :   Finalize after horizontal velocity. Apply the boundary
70      !!             condition on the after velocity, achieve the time stepping
[1502]71      !!             by applying the Asselin filter on now fields and swapping
72      !!             the fields.
[3]73      !!
[5930]74      !! ** Method  : * Ensure after velocities transport matches time splitting
75      !!             estimate (ln_dynspg_ts=T)
[3]76      !!
[1502]77      !!              * Apply lateral boundary conditions on after velocity
78      !!             at the local domain boundaries through lbc_lnk call,
[7646]79      !!             at the one-way open boundaries (ln_bdy=T),
[4990]80      !!             at the AGRIF zoom   boundaries (lk_agrif=T)
[3]81      !!
[1502]82      !!              * Apply the time filter applied and swap of the dynamics
83      !!             arrays to start the next time step:
84      !!                (ub,vb) = (un,vn) + atfp [ (ub,vb) + (ua,va) - 2 (un,vn) ]
85      !!                (un,vn) = (ua,va).
[6140]86      !!             Note that with flux form advection and non linear free surface,
87      !!             the time filter is applied on thickness weighted velocity.
88      !!             As a result, dyn_nxt MUST be called after tra_nxt.
[1502]89      !!
90      !! ** Action :   ub,vb   filtered before horizontal velocity of next time-step
91      !!               un,vn   now horizontal velocity of next time-step
[3]92      !!----------------------------------------------------------------------
93      INTEGER, INTENT( in ) ::   kt      ! ocean time-step index
[2715]94      !
[3]95      INTEGER  ::   ji, jj, jk   ! dummy loop indices
[6140]96      INTEGER  ::   ikt          ! local integers
97      REAL(wp) ::   zue3a, zue3n, zue3b, zuf, zcoef    ! local scalars
[4990]98      REAL(wp) ::   zve3a, zve3n, zve3b, zvf, z1_2dt   !   -      -
[9019]99      REAL(wp), ALLOCATABLE, DIMENSION(:,:)   ::   zue, zve
100      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) ::   ze3u_f, ze3v_f, zua, zva 
[1502]101      !!----------------------------------------------------------------------
[3294]102      !
[9019]103      IF( ln_timing    )   CALL timing_start('dyn_nxt')
104      IF( ln_dynspg_ts )   ALLOCATE( zue(jpi,jpj)     , zve(jpi,jpj)     )
105      IF( l_trddyn     )   ALLOCATE( zua(jpi,jpj,jpk) , zva(jpi,jpj,jpk) )
[3294]106      !
[3]107      IF( kt == nit000 ) THEN
108         IF(lwp) WRITE(numout,*)
109         IF(lwp) WRITE(numout,*) 'dyn_nxt : time stepping'
110         IF(lwp) WRITE(numout,*) '~~~~~~~'
111      ENDIF
112
[5930]113      IF ( ln_dynspg_ts ) THEN
114         ! Ensure below that barotropic velocities match time splitting estimate
115         ! Compute actual transport and replace it with ts estimate at "after" time step
[7753]116         zue(:,:) = e3u_a(:,:,1) * ua(:,:,1) * umask(:,:,1)
117         zve(:,:) = e3v_a(:,:,1) * va(:,:,1) * vmask(:,:,1)
[5930]118         DO jk = 2, jpkm1
[7753]119            zue(:,:) = zue(:,:) + e3u_a(:,:,jk) * ua(:,:,jk) * umask(:,:,jk)
120            zve(:,:) = zve(:,:) + e3v_a(:,:,jk) * va(:,:,jk) * vmask(:,:,jk)
[1502]121         END DO
122         DO jk = 1, jpkm1
[7753]123            ua(:,:,jk) = ( ua(:,:,jk) - zue(:,:) * r1_hu_a(:,:) + ua_b(:,:) ) * umask(:,:,jk)
124            va(:,:,jk) = ( va(:,:,jk) - zve(:,:) * r1_hv_a(:,:) + va_b(:,:) ) * vmask(:,:,jk)
[592]125         END DO
[6140]126         !
127         IF( .NOT.ln_bt_fw ) THEN
[5930]128            ! Remove advective velocity from "now velocities"
129            ! prior to asselin filtering     
130            ! In the forward case, this is done below after asselin filtering   
131            ! so that asselin contribution is removed at the same time
132            DO jk = 1, jpkm1
[9023]133               un(:,:,jk) = ( un(:,:,jk) - un_adv(:,:)*r1_hu_n(:,:) + un_b(:,:) )*umask(:,:,jk)
134               vn(:,:,jk) = ( vn(:,:,jk) - vn_adv(:,:)*r1_hv_n(:,:) + vn_b(:,:) )*vmask(:,:,jk)
[7753]135            END DO 
[5930]136         ENDIF
[4292]137      ENDIF
138
[1502]139      ! Update after velocity on domain lateral boundaries
140      ! --------------------------------------------------     
[5930]141# if defined key_agrif
142      CALL Agrif_dyn( kt )             !* AGRIF zoom boundaries
143# endif
144      !
[9090]145      CALL lbc_lnk_multi( ua, 'U', -1., va, 'V', -1. )     !* local domain boundaries
[1502]146      !
147      !                                !* BDY open boundaries
[7646]148      IF( ln_bdy .AND. ln_dynspg_exp )   CALL bdy_dyn( kt )
149      IF( ln_bdy .AND. ln_dynspg_ts  )   CALL bdy_dyn( kt, dyn3d_only=.true. )
[3294]150
151!!$   Do we need a call to bdy_vol here??
152      !
[4990]153      IF( l_trddyn ) THEN             ! prepare the atf trend computation + some diagnostics
154         z1_2dt = 1._wp / (2. * rdt)        ! Euler or leap-frog time step
155         IF( neuler == 0 .AND. kt == nit000 )   z1_2dt = 1._wp / rdt
156         !
157         !                                  ! Kinetic energy and Conversion
158         IF( ln_KE_trd  )   CALL trd_dyn( ua, va, jpdyn_ken, kt )
159         !
160         IF( ln_dyn_trd ) THEN              ! 3D output: total momentum trends
[7753]161            zua(:,:,:) = ( ua(:,:,:) - ub(:,:,:) ) * z1_2dt
162            zva(:,:,:) = ( va(:,:,:) - vb(:,:,:) ) * z1_2dt
[4990]163            CALL iom_put( "utrd_tot", zua )        ! total momentum trends, except the asselin time filter
164            CALL iom_put( "vtrd_tot", zva )
165         ENDIF
166         !
[7753]167         zua(:,:,:) = un(:,:,:)             ! save the now velocity before the asselin filter
168         zva(:,:,:) = vn(:,:,:)             ! (caution: there will be a shift by 1 timestep in the
169         !                                  !  computation of the asselin filter trends)
[4990]170      ENDIF
171
[1438]172      ! Time filter and swap of dynamics arrays
173      ! ------------------------------------------
[1502]174      IF( neuler == 0 .AND. kt == nit000 ) THEN        !* Euler at first time-step: only swap
175         DO jk = 1, jpkm1
[9226]176            un(:,:,jk) = ua(:,:,jk)                         ! un <-- ua
[7753]177            vn(:,:,jk) = va(:,:,jk)
[1438]178         END DO
[9226]179         IF( .NOT.ln_linssh ) THEN                          ! e3._b <-- e3._n
180!!gm BUG ????    I don't understand why it is not : e3._n <-- e3._a 
[4292]181            DO jk = 1, jpkm1
[9226]182!               e3t_b(:,:,jk) = e3t_n(:,:,jk)
183!               e3u_b(:,:,jk) = e3u_n(:,:,jk)
184!               e3v_b(:,:,jk) = e3v_n(:,:,jk)
185               !
186               e3t_n(:,:,jk) = e3t_a(:,:,jk)
187               e3u_n(:,:,jk) = e3u_a(:,:,jk)
188               e3v_n(:,:,jk) = e3v_a(:,:,jk)
[6140]189            END DO
[9226]190!!gm BUG end
[4292]191         ENDIF
[9226]192                                                            !
193         
[1502]194      ELSE                                             !* Leap-Frog : Asselin filter and swap
[2528]195         !                                ! =============!
[6140]196         IF( ln_linssh ) THEN             ! Fixed volume !
[2528]197            !                             ! =============!
[1502]198            DO jk = 1, jpkm1                             
[592]199               DO jj = 1, jpj
[1502]200                  DO ji = 1, jpi   
[4990]201                     zuf = un(ji,jj,jk) + atfp * ( ub(ji,jj,jk) - 2._wp * un(ji,jj,jk) + ua(ji,jj,jk) )
202                     zvf = vn(ji,jj,jk) + atfp * ( vb(ji,jj,jk) - 2._wp * vn(ji,jj,jk) + va(ji,jj,jk) )
[1502]203                     !
204                     ub(ji,jj,jk) = zuf                      ! ub <-- filtered velocity
205                     vb(ji,jj,jk) = zvf
206                     un(ji,jj,jk) = ua(ji,jj,jk)             ! un <-- ua
207                     vn(ji,jj,jk) = va(ji,jj,jk)
208                  END DO
209               END DO
210            END DO
[2528]211            !                             ! ================!
212         ELSE                             ! Variable volume !
213            !                             ! ================!
[4292]214            ! Before scale factor at t-points
215            ! (used as a now filtered scale factor until the swap)
216            ! ----------------------------------------------------
[9023]217            DO jk = 1, jpkm1
218               e3t_b(:,:,jk) = e3t_n(:,:,jk) + atfp * ( e3t_b(:,:,jk) - 2._wp * e3t_n(:,:,jk) + e3t_a(:,:,jk) )
219            END DO
220            ! Add volume filter correction: compatibility with tracer advection scheme
221            ! => time filter + conservation correction (only at the first level)
222            zcoef = atfp * rdt * r1_rau0
223            IF ( .NOT. ln_isf ) THEN   ! if no ice shelf melting
224               e3t_b(:,:,1) = e3t_b(:,:,1) - zcoef * ( emp_b(:,:) - emp(:,:) &
225                              &                      - rnf_b(:,:) + rnf(:,:) ) * tmask(:,:,1)
226            ELSE                     ! if ice shelf melting
227               DO jj = 1, jpj
228                  DO ji = 1, jpi
229                     ikt = mikt(ji,jj)
230                     e3t_b(ji,jj,ikt) = e3t_b(ji,jj,ikt) - zcoef * (  emp_b   (ji,jj) - emp   (ji,jj)  &
231                        &                                           - rnf_b   (ji,jj) + rnf   (ji,jj)  &
232                        &                                           + fwfisf_b(ji,jj) - fwfisf(ji,jj)  ) * tmask(ji,jj,ikt)
233                  END DO
[6140]234               END DO
[4292]235               ! Add volume filter correction: compatibility with tracer advection scheme
236               ! => time filter + conservation correction (only at the first level)
[6140]237               zcoef = atfp * rdt * r1_rau0
238               IF ( .NOT. ln_isf ) THEN   ! if no ice shelf melting
[7753]239                  e3t_b(:,:,1) = e3t_b(:,:,1) - zcoef * ( emp_b(:,:) - emp(:,:) &
[9023]240                                 &                      ) * tmask(:,:,1)
241                  IF( ln_rnf_depth ) THEN
242                     DO jk = 1, jpkm1 ! Deal with Rivers separetely, as can be through depth too, not sure for ice shelf case yet
243                        DO jj = 1, jpj
244                           DO ji = 1, jpi
245                              IF( mikt(ji,jj) <= jk .and. jk <=  nk_rnf(ji,jj)  ) THEN
[9119]246                                 e3t_b(ji,jj,jk) =   e3t_b(ji,jj,jk) - zcoef *  ( - rnf_b(ji,jj) + rnf(ji,jj) ) &
247                                      &          * ( e3t_n(ji,jj,jk) / h_rnf(ji,jj) ) * tmask(ji,jj,jk)
[9023]248                              ENDIF
249                           ENDDO
250                        ENDDO
251                     ENDDO
252                  ELSE
253                      e3t_b(:,:,1) = e3t_b(:,:,1) - zcoef *  ( -rnf_b(:,:) + rnf(:,:))*tmask(:,:,1)
254                  ENDIF
[5643]255               ELSE                     ! if ice shelf melting
[6140]256                  DO jj = 1, jpj
257                     DO ji = 1, jpi
258                        ikt = mikt(ji,jj)
259                        e3t_b(ji,jj,ikt) = e3t_b(ji,jj,ikt) - zcoef * (  emp_b   (ji,jj) - emp   (ji,jj)  &
260                           &                                           - rnf_b   (ji,jj) + rnf   (ji,jj)  &
261                           &                                           + fwfisf_b(ji,jj) - fwfisf(ji,jj)  ) * tmask(ji,jj,ikt)
[5643]262                     END DO
263                  END DO
264               END IF
[9023]265            END IF
[2528]266            !
[6140]267            IF( ln_dynadv_vec ) THEN      ! Asselin filter applied on velocity
268               ! Before filtered scale factor at (u/v)-points
269               CALL dom_vvl_interpol( e3t_b(:,:,:), e3u_b(:,:,:), 'U' )
270               CALL dom_vvl_interpol( e3t_b(:,:,:), e3v_b(:,:,:), 'V' )
[4292]271               DO jk = 1, jpkm1
272                  DO jj = 1, jpj
[2528]273                     DO ji = 1, jpi
[4292]274                        zuf = un(ji,jj,jk) + atfp * ( ub(ji,jj,jk) - 2._wp * un(ji,jj,jk) + ua(ji,jj,jk) )
275                        zvf = vn(ji,jj,jk) + atfp * ( vb(ji,jj,jk) - 2._wp * vn(ji,jj,jk) + va(ji,jj,jk) )
[2528]276                        !
277                        ub(ji,jj,jk) = zuf                      ! ub <-- filtered velocity
278                        vb(ji,jj,jk) = zvf
279                        un(ji,jj,jk) = ua(ji,jj,jk)             ! un <-- ua
280                        vn(ji,jj,jk) = va(ji,jj,jk)
281                     END DO
282                  END DO
283               END DO
284               !
[6140]285            ELSE                          ! Asselin filter applied on thickness weighted velocity
286               !
[9019]287               ALLOCATE( ze3u_f(jpi,jpj,jpk) , ze3v_f(jpi,jpj,jpk) )
[6140]288               ! Before filtered scale factor at (u/v)-points stored in ze3u_f, ze3v_f
289               CALL dom_vvl_interpol( e3t_b(:,:,:), ze3u_f, 'U' )
290               CALL dom_vvl_interpol( e3t_b(:,:,:), ze3v_f, 'V' )
[4292]291               DO jk = 1, jpkm1
292                  DO jj = 1, jpj
[4312]293                     DO ji = 1, jpi                 
[6140]294                        zue3a = e3u_a(ji,jj,jk) * ua(ji,jj,jk)
295                        zve3a = e3v_a(ji,jj,jk) * va(ji,jj,jk)
296                        zue3n = e3u_n(ji,jj,jk) * un(ji,jj,jk)
297                        zve3n = e3v_n(ji,jj,jk) * vn(ji,jj,jk)
298                        zue3b = e3u_b(ji,jj,jk) * ub(ji,jj,jk)
299                        zve3b = e3v_b(ji,jj,jk) * vb(ji,jj,jk)
[2528]300                        !
[3294]301                        zuf = ( zue3n + atfp * ( zue3b - 2._wp * zue3n  + zue3a ) ) / ze3u_f(ji,jj,jk)
302                        zvf = ( zve3n + atfp * ( zve3b - 2._wp * zve3n  + zve3a ) ) / ze3v_f(ji,jj,jk)
[2528]303                        !
[3294]304                        ub(ji,jj,jk) = zuf                     ! ub <-- filtered velocity
[2528]305                        vb(ji,jj,jk) = zvf
[3294]306                        un(ji,jj,jk) = ua(ji,jj,jk)            ! un <-- ua
[2528]307                        vn(ji,jj,jk) = va(ji,jj,jk)
308                     END DO
309                  END DO
310               END DO
[7753]311               e3u_b(:,:,1:jpkm1) = ze3u_f(:,:,1:jpkm1)        ! e3u_b <-- filtered scale factor
312               e3v_b(:,:,1:jpkm1) = ze3v_f(:,:,1:jpkm1)
[6140]313               !
[9019]314               DEALLOCATE( ze3u_f , ze3v_f )
[2528]315            ENDIF
316            !
[3]317         ENDIF
[2528]318         !
[6140]319         IF( ln_dynspg_ts .AND. ln_bt_fw ) THEN
[4312]320            ! Revert "before" velocities to time split estimate
321            ! Doing it here also means that asselin filter contribution is removed 
[7753]322            zue(:,:) = e3u_b(:,:,1) * ub(:,:,1) * umask(:,:,1)
323            zve(:,:) = e3v_b(:,:,1) * vb(:,:,1) * vmask(:,:,1)   
[4990]324            DO jk = 2, jpkm1
[7753]325               zue(:,:) = zue(:,:) + e3u_b(:,:,jk) * ub(:,:,jk) * umask(:,:,jk)
326               zve(:,:) = zve(:,:) + e3v_b(:,:,jk) * vb(:,:,jk) * vmask(:,:,jk)   
[4370]327            END DO
328            DO jk = 1, jpkm1
[7753]329               ub(:,:,jk) = ub(:,:,jk) - (zue(:,:) * r1_hu_n(:,:) - un_b(:,:)) * umask(:,:,jk)
330               vb(:,:,jk) = vb(:,:,jk) - (zve(:,:) * r1_hv_n(:,:) - vn_b(:,:)) * vmask(:,:,jk)
[4292]331            END DO
332         ENDIF
333         !
334      ENDIF ! neuler =/0
[4354]335      !
336      ! Set "now" and "before" barotropic velocities for next time step:
337      ! JC: Would be more clever to swap variables than to make a full vertical
338      ! integration
339      !
[4370]340      !
[6140]341      IF(.NOT.ln_linssh ) THEN
[7753]342         hu_b(:,:) = e3u_b(:,:,1) * umask(:,:,1)
343         hv_b(:,:) = e3v_b(:,:,1) * vmask(:,:,1)
[6140]344         DO jk = 2, jpkm1
[7753]345            hu_b(:,:) = hu_b(:,:) + e3u_b(:,:,jk) * umask(:,:,jk)
346            hv_b(:,:) = hv_b(:,:) + e3v_b(:,:,jk) * vmask(:,:,jk)
[4354]347         END DO
[7753]348         r1_hu_b(:,:) = ssumask(:,:) / ( hu_b(:,:) + 1._wp - ssumask(:,:) )
349         r1_hv_b(:,:) = ssvmask(:,:) / ( hv_b(:,:) + 1._wp - ssvmask(:,:) )
[4354]350      ENDIF
351      !
[7753]352      un_b(:,:) = e3u_a(:,:,1) * un(:,:,1) * umask(:,:,1)
353      ub_b(:,:) = e3u_b(:,:,1) * ub(:,:,1) * umask(:,:,1)
354      vn_b(:,:) = e3v_a(:,:,1) * vn(:,:,1) * vmask(:,:,1)
355      vb_b(:,:) = e3v_b(:,:,1) * vb(:,:,1) * vmask(:,:,1)
[6140]356      DO jk = 2, jpkm1
[7753]357         un_b(:,:) = un_b(:,:) + e3u_a(:,:,jk) * un(:,:,jk) * umask(:,:,jk)
358         ub_b(:,:) = ub_b(:,:) + e3u_b(:,:,jk) * ub(:,:,jk) * umask(:,:,jk)
359         vn_b(:,:) = vn_b(:,:) + e3v_a(:,:,jk) * vn(:,:,jk) * vmask(:,:,jk)
360         vb_b(:,:) = vb_b(:,:) + e3v_b(:,:,jk) * vb(:,:,jk) * vmask(:,:,jk)
[4354]361      END DO
[7753]362      un_b(:,:) = un_b(:,:) * r1_hu_a(:,:)
363      vn_b(:,:) = vn_b(:,:) * r1_hv_a(:,:)
364      ub_b(:,:) = ub_b(:,:) * r1_hu_b(:,:)
365      vb_b(:,:) = vb_b(:,:) * r1_hv_b(:,:)
[4354]366      !
[6140]367      IF( .NOT.ln_dynspg_ts ) THEN        ! output the barotropic currents
368         CALL iom_put(  "ubar", un_b(:,:) )
369         CALL iom_put(  "vbar", vn_b(:,:) )
370      ENDIF
[4990]371      IF( l_trddyn ) THEN                ! 3D output: asselin filter trends on momentum
[7753]372         zua(:,:,:) = ( ub(:,:,:) - zua(:,:,:) ) * z1_2dt
373         zva(:,:,:) = ( vb(:,:,:) - zva(:,:,:) ) * z1_2dt
[4990]374         CALL trd_dyn( zua, zva, jpdyn_atf, kt )
375      ENDIF
376      !
[1438]377      IF(ln_ctl)   CALL prt_ctl( tab3d_1=un, clinfo1=' nxt  - Un: ', mask1=umask,   &
378         &                       tab3d_2=vn, clinfo2=' Vn: '       , mask2=vmask )
[6140]379      !
[9019]380      IF( ln_dynspg_ts )   DEALLOCATE( zue, zve )
381      IF( l_trddyn     )   DEALLOCATE( zua, zva )
382      IF( ln_timing    )   CALL timing_stop('dyn_nxt')
[2715]383      !
[3]384   END SUBROUTINE dyn_nxt
385
[1502]386   !!=========================================================================
[3]387END MODULE dynnxt
Note: See TracBrowser for help on using the repository browser.