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 NEMO/branches/UKMO/NEMO_4.0.1_FKOSM_m11715/src/OCE/DYN – NEMO

source: NEMO/branches/UKMO/NEMO_4.0.1_FKOSM_m11715/src/OCE/DYN/dynnxt.F90 @ 13454

Last change on this file since 13454 was 13454, checked in by dancopsey, 4 years ago

Merge in Catherine Guiavarch's change to implement implicit ice drag to improve stability. Her changes were implemented in:

https://code.metoffice.gov.uk/trac/roses-u/changeset?reponame=&new=161415%40b%2Fv%2F1%2F1%2F4&old=161163%40b%2Fv%2F1%2F1%2F4

File size: 20.1 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
[9361]31   USE sbcisf         ! ice shelf
[6140]32   USE phycst         ! physical constants
33   USE dynadv         ! dynamics: vector invariant versus flux form
34   USE dynspg_ts      ! surface pressure gradient: split-explicit scheme
35   USE domvvl         ! variable volume
[7646]36   USE bdy_oce   , ONLY: ln_bdy
[6140]37   USE bdydta         ! ocean open boundary conditions
38   USE bdydyn         ! ocean open boundary conditions
39   USE bdyvol         ! ocean open boundary condition (bdy_vol routines)
40   USE trd_oce        ! trends: ocean variables
41   USE trddyn         ! trend manager: dynamics
42   USE trdken         ! trend manager: kinetic energy
[4990]43   !
[6140]44   USE in_out_manager ! I/O manager
45   USE iom            ! I/O manager library
46   USE lbclnk         ! lateral boundary condition (or mpp link)
47   USE lib_mpp        ! MPP library
48   USE prtctl         ! Print control
49   USE timing         ! Timing
[13454]50   USE zdfdrg, ONLY: ln_drgice_imp, rCdU_top
[2528]51#if defined key_agrif
[9570]52   USE agrif_oce_interp
[2528]53#endif
[3]54
55   IMPLICIT NONE
56   PRIVATE
57
[1438]58   PUBLIC    dyn_nxt   ! routine called by step.F90
59
[2715]60   !!----------------------------------------------------------------------
[9598]61   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
[11715]62   !! $Id$
[10068]63   !! Software governed by the CeCILL license (see ./LICENSE)
[2715]64   !!----------------------------------------------------------------------
[3]65CONTAINS
66
67   SUBROUTINE dyn_nxt ( kt )
68      !!----------------------------------------------------------------------
69      !!                  ***  ROUTINE dyn_nxt  ***
70      !!                   
[5930]71      !! ** Purpose :   Finalize after horizontal velocity. Apply the boundary
72      !!             condition on the after velocity, achieve the time stepping
[1502]73      !!             by applying the Asselin filter on now fields and swapping
74      !!             the fields.
[3]75      !!
[5930]76      !! ** Method  : * Ensure after velocities transport matches time splitting
77      !!             estimate (ln_dynspg_ts=T)
[3]78      !!
[1502]79      !!              * Apply lateral boundary conditions on after velocity
80      !!             at the local domain boundaries through lbc_lnk call,
[7646]81      !!             at the one-way open boundaries (ln_bdy=T),
[4990]82      !!             at the AGRIF zoom   boundaries (lk_agrif=T)
[3]83      !!
[1502]84      !!              * Apply the time filter applied and swap of the dynamics
85      !!             arrays to start the next time step:
86      !!                (ub,vb) = (un,vn) + atfp [ (ub,vb) + (ua,va) - 2 (un,vn) ]
87      !!                (un,vn) = (ua,va).
[6140]88      !!             Note that with flux form advection and non linear free surface,
89      !!             the time filter is applied on thickness weighted velocity.
90      !!             As a result, dyn_nxt MUST be called after tra_nxt.
[1502]91      !!
92      !! ** Action :   ub,vb   filtered before horizontal velocity of next time-step
93      !!               un,vn   now horizontal velocity of next time-step
[3]94      !!----------------------------------------------------------------------
95      INTEGER, INTENT( in ) ::   kt      ! ocean time-step index
[2715]96      !
[3]97      INTEGER  ::   ji, jj, jk   ! dummy loop indices
[6140]98      INTEGER  ::   ikt          ! local integers
99      REAL(wp) ::   zue3a, zue3n, zue3b, zuf, zcoef    ! local scalars
[4990]100      REAL(wp) ::   zve3a, zve3n, zve3b, zvf, z1_2dt   !   -      -
[9019]101      REAL(wp), ALLOCATABLE, DIMENSION(:,:)   ::   zue, zve
[13454]102      REAL(wp), ALLOCATABLE, DIMENSION(:,:)   ::   zutau, zvtau
[9019]103      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) ::   ze3u_f, ze3v_f, zua, zva 
[1502]104      !!----------------------------------------------------------------------
[3294]105      !
[9019]106      IF( ln_timing    )   CALL timing_start('dyn_nxt')
107      IF( ln_dynspg_ts )   ALLOCATE( zue(jpi,jpj)     , zve(jpi,jpj)     )
108      IF( l_trddyn     )   ALLOCATE( zua(jpi,jpj,jpk) , zva(jpi,jpj,jpk) )
[3294]109      !
[3]110      IF( kt == nit000 ) THEN
111         IF(lwp) WRITE(numout,*)
112         IF(lwp) WRITE(numout,*) 'dyn_nxt : time stepping'
113         IF(lwp) WRITE(numout,*) '~~~~~~~'
114      ENDIF
115
[5930]116      IF ( ln_dynspg_ts ) THEN
117         ! Ensure below that barotropic velocities match time splitting estimate
118         ! Compute actual transport and replace it with ts estimate at "after" time step
[7753]119         zue(:,:) = e3u_a(:,:,1) * ua(:,:,1) * umask(:,:,1)
120         zve(:,:) = e3v_a(:,:,1) * va(:,:,1) * vmask(:,:,1)
[5930]121         DO jk = 2, jpkm1
[7753]122            zue(:,:) = zue(:,:) + e3u_a(:,:,jk) * ua(:,:,jk) * umask(:,:,jk)
123            zve(:,:) = zve(:,:) + e3v_a(:,:,jk) * va(:,:,jk) * vmask(:,:,jk)
[1502]124         END DO
125         DO jk = 1, jpkm1
[7753]126            ua(:,:,jk) = ( ua(:,:,jk) - zue(:,:) * r1_hu_a(:,:) + ua_b(:,:) ) * umask(:,:,jk)
127            va(:,:,jk) = ( va(:,:,jk) - zve(:,:) * r1_hv_a(:,:) + va_b(:,:) ) * vmask(:,:,jk)
[592]128         END DO
[6140]129         !
130         IF( .NOT.ln_bt_fw ) THEN
[5930]131            ! Remove advective velocity from "now velocities"
132            ! prior to asselin filtering     
133            ! In the forward case, this is done below after asselin filtering   
134            ! so that asselin contribution is removed at the same time
135            DO jk = 1, jpkm1
[9023]136               un(:,:,jk) = ( un(:,:,jk) - un_adv(:,:)*r1_hu_n(:,:) + un_b(:,:) )*umask(:,:,jk)
137               vn(:,:,jk) = ( vn(:,:,jk) - vn_adv(:,:)*r1_hv_n(:,:) + vn_b(:,:) )*vmask(:,:,jk)
[7753]138            END DO 
[5930]139         ENDIF
[4292]140      ENDIF
141
[1502]142      ! Update after velocity on domain lateral boundaries
143      ! --------------------------------------------------     
[5930]144# if defined key_agrif
145      CALL Agrif_dyn( kt )             !* AGRIF zoom boundaries
146# endif
147      !
[10425]148      CALL lbc_lnk_multi( 'dynnxt', ua, 'U', -1., va, 'V', -1. )     !* local domain boundaries
[1502]149      !
150      !                                !* BDY open boundaries
[7646]151      IF( ln_bdy .AND. ln_dynspg_exp )   CALL bdy_dyn( kt )
152      IF( ln_bdy .AND. ln_dynspg_ts  )   CALL bdy_dyn( kt, dyn3d_only=.true. )
[3294]153
154!!$   Do we need a call to bdy_vol here??
155      !
[4990]156      IF( l_trddyn ) THEN             ! prepare the atf trend computation + some diagnostics
157         z1_2dt = 1._wp / (2. * rdt)        ! Euler or leap-frog time step
158         IF( neuler == 0 .AND. kt == nit000 )   z1_2dt = 1._wp / rdt
159         !
160         !                                  ! Kinetic energy and Conversion
161         IF( ln_KE_trd  )   CALL trd_dyn( ua, va, jpdyn_ken, kt )
162         !
163         IF( ln_dyn_trd ) THEN              ! 3D output: total momentum trends
[7753]164            zua(:,:,:) = ( ua(:,:,:) - ub(:,:,:) ) * z1_2dt
165            zva(:,:,:) = ( va(:,:,:) - vb(:,:,:) ) * z1_2dt
[4990]166            CALL iom_put( "utrd_tot", zua )        ! total momentum trends, except the asselin time filter
167            CALL iom_put( "vtrd_tot", zva )
168         ENDIF
169         !
[7753]170         zua(:,:,:) = un(:,:,:)             ! save the now velocity before the asselin filter
171         zva(:,:,:) = vn(:,:,:)             ! (caution: there will be a shift by 1 timestep in the
172         !                                  !  computation of the asselin filter trends)
[4990]173      ENDIF
174
[1438]175      ! Time filter and swap of dynamics arrays
176      ! ------------------------------------------
[1502]177      IF( neuler == 0 .AND. kt == nit000 ) THEN        !* Euler at first time-step: only swap
178         DO jk = 1, jpkm1
[9226]179            un(:,:,jk) = ua(:,:,jk)                         ! un <-- ua
[7753]180            vn(:,:,jk) = va(:,:,jk)
[1438]181         END DO
[9226]182         IF( .NOT.ln_linssh ) THEN                          ! e3._b <-- e3._n
183!!gm BUG ????    I don't understand why it is not : e3._n <-- e3._a 
[4292]184            DO jk = 1, jpkm1
[9226]185!               e3t_b(:,:,jk) = e3t_n(:,:,jk)
186!               e3u_b(:,:,jk) = e3u_n(:,:,jk)
187!               e3v_b(:,:,jk) = e3v_n(:,:,jk)
188               !
189               e3t_n(:,:,jk) = e3t_a(:,:,jk)
190               e3u_n(:,:,jk) = e3u_a(:,:,jk)
191               e3v_n(:,:,jk) = e3v_a(:,:,jk)
[6140]192            END DO
[9226]193!!gm BUG end
[4292]194         ENDIF
[9226]195                                                            !
196         
[1502]197      ELSE                                             !* Leap-Frog : Asselin filter and swap
[2528]198         !                                ! =============!
[6140]199         IF( ln_linssh ) THEN             ! Fixed volume !
[2528]200            !                             ! =============!
[1502]201            DO jk = 1, jpkm1                             
[592]202               DO jj = 1, jpj
[1502]203                  DO ji = 1, jpi   
[4990]204                     zuf = un(ji,jj,jk) + atfp * ( ub(ji,jj,jk) - 2._wp * un(ji,jj,jk) + ua(ji,jj,jk) )
205                     zvf = vn(ji,jj,jk) + atfp * ( vb(ji,jj,jk) - 2._wp * vn(ji,jj,jk) + va(ji,jj,jk) )
[1502]206                     !
207                     ub(ji,jj,jk) = zuf                      ! ub <-- filtered velocity
208                     vb(ji,jj,jk) = zvf
209                     un(ji,jj,jk) = ua(ji,jj,jk)             ! un <-- ua
210                     vn(ji,jj,jk) = va(ji,jj,jk)
211                  END DO
212               END DO
213            END DO
[2528]214            !                             ! ================!
215         ELSE                             ! Variable volume !
216            !                             ! ================!
[4292]217            ! Before scale factor at t-points
218            ! (used as a now filtered scale factor until the swap)
219            ! ----------------------------------------------------
[9023]220            DO jk = 1, jpkm1
221               e3t_b(:,:,jk) = e3t_n(:,:,jk) + atfp * ( e3t_b(:,:,jk) - 2._wp * e3t_n(:,:,jk) + e3t_a(:,:,jk) )
222            END DO
223            ! Add volume filter correction: compatibility with tracer advection scheme
224            ! => time filter + conservation correction (only at the first level)
225            zcoef = atfp * rdt * r1_rau0
[9361]226
227            e3t_b(:,:,1) = e3t_b(:,:,1) - zcoef * ( emp_b(:,:) - emp(:,:) ) * tmask(:,:,1)
228
229            IF ( ln_rnf ) THEN
230               IF( ln_rnf_depth ) THEN
231                  DO jk = 1, jpkm1 ! Deal with Rivers separetely, as can be through depth too
232                     DO jj = 1, jpj
233                        DO ji = 1, jpi
234                           IF( jk <=  nk_rnf(ji,jj)  ) THEN
235                               e3t_b(ji,jj,jk) =   e3t_b(ji,jj,jk) - zcoef *  ( - rnf_b(ji,jj) + rnf(ji,jj) ) &
[9119]236                                      &          * ( e3t_n(ji,jj,jk) / h_rnf(ji,jj) ) * tmask(ji,jj,jk)
[9361]237                           ENDIF
[9023]238                        ENDDO
239                     ENDDO
[9361]240                  ENDDO
241               ELSE
242                  e3t_b(:,:,1) = e3t_b(:,:,1) - zcoef *  ( -rnf_b(:,:) + rnf(:,:))*tmask(:,:,1)
243               ENDIF
244            END IF
245
246            IF ( ln_isf ) THEN   ! if ice shelf melting
247               DO jk = 1, jpkm1 ! Deal with isf separetely, as can be through depth too
[6140]248                  DO jj = 1, jpj
249                     DO ji = 1, jpi
[10349]250                        IF( misfkt(ji,jj) <=jk .and. jk < misfkb(ji,jj)  ) THEN
251                           e3t_b(ji,jj,jk) = e3t_b(ji,jj,jk) - zcoef * ( fwfisf_b(ji,jj) - fwfisf(ji,jj) ) &
[9361]252                                &          * ( e3t_n(ji,jj,jk) * r1_hisf_tbl(ji,jj) ) * tmask(ji,jj,jk)
[10349]253                        ELSEIF ( jk==misfkb(ji,jj) ) THEN
254                           e3t_b(ji,jj,jk) = e3t_b(ji,jj,jk) - zcoef * ( fwfisf_b(ji,jj) - fwfisf(ji,jj) ) &
255                                &          * ( e3t_n(ji,jj,jk) * r1_hisf_tbl(ji,jj) ) * ralpha(ji,jj) * tmask(ji,jj,jk)
[9361]256                        ENDIF
[5643]257                     END DO
258                  END DO
[9361]259               END DO
[9023]260            END IF
[2528]261            !
[6140]262            IF( ln_dynadv_vec ) THEN      ! Asselin filter applied on velocity
263               ! Before filtered scale factor at (u/v)-points
264               CALL dom_vvl_interpol( e3t_b(:,:,:), e3u_b(:,:,:), 'U' )
265               CALL dom_vvl_interpol( e3t_b(:,:,:), e3v_b(:,:,:), 'V' )
[4292]266               DO jk = 1, jpkm1
267                  DO jj = 1, jpj
[2528]268                     DO ji = 1, jpi
[4292]269                        zuf = un(ji,jj,jk) + atfp * ( ub(ji,jj,jk) - 2._wp * un(ji,jj,jk) + ua(ji,jj,jk) )
270                        zvf = vn(ji,jj,jk) + atfp * ( vb(ji,jj,jk) - 2._wp * vn(ji,jj,jk) + va(ji,jj,jk) )
[2528]271                        !
272                        ub(ji,jj,jk) = zuf                      ! ub <-- filtered velocity
273                        vb(ji,jj,jk) = zvf
274                        un(ji,jj,jk) = ua(ji,jj,jk)             ! un <-- ua
275                        vn(ji,jj,jk) = va(ji,jj,jk)
276                     END DO
277                  END DO
278               END DO
279               !
[6140]280            ELSE                          ! Asselin filter applied on thickness weighted velocity
281               !
[9019]282               ALLOCATE( ze3u_f(jpi,jpj,jpk) , ze3v_f(jpi,jpj,jpk) )
[6140]283               ! Before filtered scale factor at (u/v)-points stored in ze3u_f, ze3v_f
284               CALL dom_vvl_interpol( e3t_b(:,:,:), ze3u_f, 'U' )
285               CALL dom_vvl_interpol( e3t_b(:,:,:), ze3v_f, 'V' )
[4292]286               DO jk = 1, jpkm1
287                  DO jj = 1, jpj
[4312]288                     DO ji = 1, jpi                 
[6140]289                        zue3a = e3u_a(ji,jj,jk) * ua(ji,jj,jk)
290                        zve3a = e3v_a(ji,jj,jk) * va(ji,jj,jk)
291                        zue3n = e3u_n(ji,jj,jk) * un(ji,jj,jk)
292                        zve3n = e3v_n(ji,jj,jk) * vn(ji,jj,jk)
293                        zue3b = e3u_b(ji,jj,jk) * ub(ji,jj,jk)
294                        zve3b = e3v_b(ji,jj,jk) * vb(ji,jj,jk)
[2528]295                        !
[3294]296                        zuf = ( zue3n + atfp * ( zue3b - 2._wp * zue3n  + zue3a ) ) / ze3u_f(ji,jj,jk)
297                        zvf = ( zve3n + atfp * ( zve3b - 2._wp * zve3n  + zve3a ) ) / ze3v_f(ji,jj,jk)
[2528]298                        !
[3294]299                        ub(ji,jj,jk) = zuf                     ! ub <-- filtered velocity
[2528]300                        vb(ji,jj,jk) = zvf
[3294]301                        un(ji,jj,jk) = ua(ji,jj,jk)            ! un <-- ua
[2528]302                        vn(ji,jj,jk) = va(ji,jj,jk)
303                     END DO
304                  END DO
305               END DO
[7753]306               e3u_b(:,:,1:jpkm1) = ze3u_f(:,:,1:jpkm1)        ! e3u_b <-- filtered scale factor
307               e3v_b(:,:,1:jpkm1) = ze3v_f(:,:,1:jpkm1)
[6140]308               !
[9019]309               DEALLOCATE( ze3u_f , ze3v_f )
[2528]310            ENDIF
311            !
[3]312         ENDIF
[2528]313         !
[6140]314         IF( ln_dynspg_ts .AND. ln_bt_fw ) THEN
[4312]315            ! Revert "before" velocities to time split estimate
316            ! Doing it here also means that asselin filter contribution is removed 
[7753]317            zue(:,:) = e3u_b(:,:,1) * ub(:,:,1) * umask(:,:,1)
318            zve(:,:) = e3v_b(:,:,1) * vb(:,:,1) * vmask(:,:,1)   
[4990]319            DO jk = 2, jpkm1
[7753]320               zue(:,:) = zue(:,:) + e3u_b(:,:,jk) * ub(:,:,jk) * umask(:,:,jk)
321               zve(:,:) = zve(:,:) + e3v_b(:,:,jk) * vb(:,:,jk) * vmask(:,:,jk)   
[4370]322            END DO
323            DO jk = 1, jpkm1
[7753]324               ub(:,:,jk) = ub(:,:,jk) - (zue(:,:) * r1_hu_n(:,:) - un_b(:,:)) * umask(:,:,jk)
325               vb(:,:,jk) = vb(:,:,jk) - (zve(:,:) * r1_hv_n(:,:) - vn_b(:,:)) * vmask(:,:,jk)
[4292]326            END DO
327         ENDIF
328         !
329      ENDIF ! neuler =/0
[4354]330      !
331      ! Set "now" and "before" barotropic velocities for next time step:
332      ! JC: Would be more clever to swap variables than to make a full vertical
333      ! integration
334      !
[4370]335      !
[6140]336      IF(.NOT.ln_linssh ) THEN
[7753]337         hu_b(:,:) = e3u_b(:,:,1) * umask(:,:,1)
338         hv_b(:,:) = e3v_b(:,:,1) * vmask(:,:,1)
[6140]339         DO jk = 2, jpkm1
[7753]340            hu_b(:,:) = hu_b(:,:) + e3u_b(:,:,jk) * umask(:,:,jk)
341            hv_b(:,:) = hv_b(:,:) + e3v_b(:,:,jk) * vmask(:,:,jk)
[4354]342         END DO
[7753]343         r1_hu_b(:,:) = ssumask(:,:) / ( hu_b(:,:) + 1._wp - ssumask(:,:) )
344         r1_hv_b(:,:) = ssvmask(:,:) / ( hv_b(:,:) + 1._wp - ssvmask(:,:) )
[4354]345      ENDIF
346      !
[7753]347      un_b(:,:) = e3u_a(:,:,1) * un(:,:,1) * umask(:,:,1)
348      ub_b(:,:) = e3u_b(:,:,1) * ub(:,:,1) * umask(:,:,1)
349      vn_b(:,:) = e3v_a(:,:,1) * vn(:,:,1) * vmask(:,:,1)
350      vb_b(:,:) = e3v_b(:,:,1) * vb(:,:,1) * vmask(:,:,1)
[6140]351      DO jk = 2, jpkm1
[7753]352         un_b(:,:) = un_b(:,:) + e3u_a(:,:,jk) * un(:,:,jk) * umask(:,:,jk)
353         ub_b(:,:) = ub_b(:,:) + e3u_b(:,:,jk) * ub(:,:,jk) * umask(:,:,jk)
354         vn_b(:,:) = vn_b(:,:) + e3v_a(:,:,jk) * vn(:,:,jk) * vmask(:,:,jk)
355         vb_b(:,:) = vb_b(:,:) + e3v_b(:,:,jk) * vb(:,:,jk) * vmask(:,:,jk)
[4354]356      END DO
[7753]357      un_b(:,:) = un_b(:,:) * r1_hu_a(:,:)
358      vn_b(:,:) = vn_b(:,:) * r1_hv_a(:,:)
359      ub_b(:,:) = ub_b(:,:) * r1_hu_b(:,:)
360      vb_b(:,:) = vb_b(:,:) * r1_hv_b(:,:)
[4354]361      !
[6140]362      IF( .NOT.ln_dynspg_ts ) THEN        ! output the barotropic currents
363         CALL iom_put(  "ubar", un_b(:,:) )
364         CALL iom_put(  "vbar", vn_b(:,:) )
365      ENDIF
[4990]366      IF( l_trddyn ) THEN                ! 3D output: asselin filter trends on momentum
[7753]367         zua(:,:,:) = ( ub(:,:,:) - zua(:,:,:) ) * z1_2dt
368         zva(:,:,:) = ( vb(:,:,:) - zva(:,:,:) ) * z1_2dt
[4990]369         CALL trd_dyn( zua, zva, jpdyn_atf, kt )
370      ENDIF
371      !
[13454]372      IF ( iom_use("utau") ) THEN
373         IF ( ln_drgice_imp.OR.ln_isfcav ) THEN
374            ALLOCATE(zutau(jpi,jpj)) 
375            DO jj = 2, jpjm1
376               DO ji = 2, jpim1
377                  jk = miku(ji,jj) 
378                  zutau(ji,jj) = utau(ji,jj) & 
379                  &  + 0.5_wp * rau0 * (rCdU_top(ji+1,jj)+rCdU_top(ji,jj)) * ua(ji,jj,jk) 
380               END DO
381            END DO
382            CALL lbc_lnk( 'dynnxt' , zutau, 'U', -1.)
383            CALL iom_put(  "utau", zutau(:,:) )
384            DEALLOCATE(zutau)
385         ELSE
386            CALL iom_put(  "utau", utau(:,:) )
387         ENDIF
388      ENDIF
389      !
390      IF ( iom_use("vtau") ) THEN
391         IF ( ln_drgice_imp.OR.ln_isfcav ) THEN
392            ALLOCATE(zvtau(jpi,jpj))
393            DO jj = 2, jpjm1
394               DO ji = 2, jpim1
395                  jk = mikv(ji,jj)
396                  zvtau(ji,jj) = vtau(ji,jj) &
397                  &  + 0.5_wp * rau0 * (rCdU_top(ji,jj+1)+rCdU_top(ji,jj)) * va(ji,jj,jk)
398               END DO
399            END DO
400            CALL lbc_lnk( 'dynnxt' , zvtau, 'V', -1.)
401            CALL iom_put(  "vtau", zvtau(:,:) )
402            DEALLOCATE(zvtau)
403         ELSE
404            CALL iom_put(  "vtau", vtau(:,:) )
405         ENDIF
406      ENDIF
407      !
[1438]408      IF(ln_ctl)   CALL prt_ctl( tab3d_1=un, clinfo1=' nxt  - Un: ', mask1=umask,   &
409         &                       tab3d_2=vn, clinfo2=' Vn: '       , mask2=vmask )
[6140]410      !
[9019]411      IF( ln_dynspg_ts )   DEALLOCATE( zue, zve )
412      IF( l_trddyn     )   DEALLOCATE( zua, zva )
413      IF( ln_timing    )   CALL timing_stop('dyn_nxt')
[2715]414      !
[3]415   END SUBROUTINE dyn_nxt
416
[1502]417   !!=========================================================================
[3]418END MODULE dynnxt
Note: See TracBrowser for help on using the repository browser.