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.
sshwzv.F90 in trunk/NEMOGCM/NEMO/OPA_SRC/DYN – NEMO

source: trunk/NEMOGCM/NEMO/OPA_SRC/DYN/sshwzv.F90 @ 7753

Last change on this file since 7753 was 7753, checked in by mocavero, 7 years ago

Reverting trunk to remove OpenMP

  • Property svn:keywords set to Id
File size: 12.8 KB
RevLine 
[1565]1MODULE sshwzv   
[3]2   !!==============================================================================
[1438]3   !!                       ***  MODULE  sshwzv  ***
4   !! Ocean dynamics : sea surface height and vertical velocity
[3]5   !!==============================================================================
[1438]6   !! History :  3.1  !  2009-02  (G. Madec, M. Leclair)  Original code
[2528]7   !!            3.3  !  2010-04  (M. Leclair, G. Madec)  modified LF-RA
8   !!             -   !  2010-05  (K. Mogensen, A. Weaver, M. Martin, D. Lea) Assimilation interface
9   !!             -   !  2010-09  (D.Storkey and E.O'Dea) bug fixes for BDY module
[4292]10   !!            3.3  !  2011-10  (M. Leclair) split former ssh_wzv routine and remove all vvl related work
[3]11   !!----------------------------------------------------------------------
[1438]12
[3]13   !!----------------------------------------------------------------------
[6140]14   !!   ssh_nxt       : after ssh
15   !!   ssh_swp       : filter ans swap the ssh arrays
16   !!   wzv           : compute now vertical velocity
[1438]17   !!----------------------------------------------------------------------
[6140]18   USE oce            ! ocean dynamics and tracers variables
19   USE dom_oce        ! ocean space and time domain variables
20   USE sbc_oce        ! surface boundary condition: ocean
21   USE domvvl         ! Variable volume
22   USE divhor         ! horizontal divergence
23   USE phycst         ! physical constants
[7646]24   USE bdy_oce   , ONLY: ln_bdy, bdytmask
[6140]25   USE bdydyn2d       ! bdy_ssh routine
[2528]26#if defined key_agrif
[2486]27   USE agrif_opa_interp
[2528]28#endif
29#if defined key_asminc   
[6140]30   USE   asminc       ! Assimilation increment
[2528]31#endif
[6140]32   !
33   USE in_out_manager ! I/O manager
34   USE restart        ! only for lrst_oce
35   USE prtctl         ! Print control
36   USE lbclnk         ! ocean lateral boundary condition (or mpp link)
37   USE lib_mpp        ! MPP library
38   USE wrk_nemo       ! Memory Allocation
39   USE timing         ! Timing
[6152]40   USE wet_dry         ! Wetting/Drying flux limting
[592]41
[3]42   IMPLICIT NONE
43   PRIVATE
44
[1438]45   PUBLIC   ssh_nxt    ! called by step.F90
[4292]46   PUBLIC   wzv        ! called by step.F90
47   PUBLIC   ssh_swp    ! called by step.F90
[3]48
49   !! * Substitutions
[1438]50#  include "vectopt_loop_substitute.h90"
[3]51   !!----------------------------------------------------------------------
[2528]52   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
[888]53   !! $Id$
[2715]54   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
[592]55   !!----------------------------------------------------------------------
[3]56CONTAINS
57
[4292]58   SUBROUTINE ssh_nxt( kt )
[3]59      !!----------------------------------------------------------------------
[4292]60      !!                ***  ROUTINE ssh_nxt  ***
[1438]61      !!                   
[4292]62      !! ** Purpose :   compute the after ssh (ssha)
[3]63      !!
[4292]64      !! ** Method  : - Using the incompressibility hypothesis, the ssh increment
65      !!      is computed by integrating the horizontal divergence and multiply by
66      !!      by the time step.
[3]67      !!
[5836]68      !! ** action  :   ssha, after sea surface height
[2528]69      !!
70      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
[3]71      !!----------------------------------------------------------------------
[5836]72      INTEGER, INTENT(in) ::   kt   ! time step
[4292]73      !
[7753]74      INTEGER  ::   jk            ! dummy loop indice
[5836]75      REAL(wp) ::   z2dt, zcoef   ! local scalars
76      REAL(wp), POINTER, DIMENSION(:,:  ) ::   zhdiv   ! 2D workspace
[3]77      !!----------------------------------------------------------------------
[3294]78      !
[5836]79      IF( nn_timing == 1 )   CALL timing_start('ssh_nxt')
[3294]80      !
[5836]81      CALL wrk_alloc( jpi,jpj,   zhdiv ) 
[3294]82      !
[3]83      IF( kt == nit000 ) THEN
84         IF(lwp) WRITE(numout,*)
[4292]85         IF(lwp) WRITE(numout,*) 'ssh_nxt : after sea surface height'
[1438]86         IF(lwp) WRITE(numout,*) '~~~~~~~ '
[3]87      ENDIF
[2528]88      !
[7646]89      z2dt = 2._wp * rdt                          ! set time step size (Euler/Leapfrog)
[2715]90      IF( neuler == 0 .AND. kt == nit000 )   z2dt = rdt
[7646]91      zcoef = 0.5_wp * r1_rau0
[3]92
[1438]93      !                                           !------------------------------!
94      !                                           !   After Sea Surface Height   !
95      !                                           !------------------------------!
[7646]96      IF(ln_wd) THEN
[7753]97         CALL wad_lmt(sshb, zcoef * (emp_b(:,:) + emp(:,:)), z2dt)
98      ENDIF
[7646]99
[7753]100      CALL div_hor( kt )                               ! Horizontal divergence
[7646]101      !
[7753]102      zhdiv(:,:) = 0._wp
[1438]103      DO jk = 1, jpkm1                                 ! Horizontal divergence of barotropic transports
[7753]104        zhdiv(:,:) = zhdiv(:,:) + e3t_n(:,:,jk) * hdivn(:,:,jk)
[1438]105      END DO
106      !                                                ! Sea surface elevation time stepping
[4338]107      ! In time-split case we need a first guess of the ssh after (using the baroclinic timestep) in order to
108      ! compute the vertical velocity which can be used to compute the non-linear terms of the momentum equations.
109      !
[7753]110      ssha(:,:) = (  sshb(:,:) - z2dt * ( zcoef * ( emp_b(:,:) + emp(:,:) ) + zhdiv(:,:) )  ) * ssmask(:,:)
111
[5930]112      IF ( .NOT.ln_dynspg_ts ) THEN
113         ! These lines are not necessary with time splitting since
114         ! boundary condition on sea level is set during ts loop
[5836]115# if defined key_agrif
[5930]116         CALL agrif_ssh( kt )
[5836]117# endif
[7646]118         IF( ln_bdy ) THEN
[5930]119            CALL lbc_lnk( ssha, 'T', 1. )    ! Not sure that's necessary
120            CALL bdy_ssh( ssha )             ! Duplicate sea level across open boundaries
121         ENDIF
[4292]122      ENDIF
[4486]123
[3764]124#if defined key_asminc
[5836]125      IF( lk_asminc .AND. ln_sshinc .AND. ln_asmiau ) THEN     ! Include the IAU weighted SSH increment
[3764]126         CALL ssh_asm_inc( kt )
[7753]127         ssha(:,:) = ssha(:,:) + z2dt * ssh_iau(:,:)
[3764]128      ENDIF
129#endif
[4292]130      !                                           !------------------------------!
131      !                                           !           outputs            !
132      !                                           !------------------------------!
133      !
134      IF(ln_ctl)   CALL prt_ctl( tab2d_1=ssha, clinfo1=' ssha  - : ', mask1=tmask, ovlap=1 )
135      !
136      CALL wrk_dealloc( jpi, jpj, zhdiv ) 
137      !
138      IF( nn_timing == 1 )  CALL timing_stop('ssh_nxt')
139      !
140   END SUBROUTINE ssh_nxt
141
142   
143   SUBROUTINE wzv( kt )
144      !!----------------------------------------------------------------------
145      !!                ***  ROUTINE wzv  ***
146      !!                   
147      !! ** Purpose :   compute the now vertical velocity
148      !!
149      !! ** Method  : - Using the incompressibility hypothesis, the vertical
150      !!      velocity is computed by integrating the horizontal divergence 
151      !!      from the bottom to the surface minus the scale factor evolution.
152      !!        The boundary conditions are w=0 at the bottom (no flux) and.
153      !!
154      !! ** action  :   wn      : now vertical velocity
155      !!
156      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
157      !!----------------------------------------------------------------------
[5836]158      INTEGER, INTENT(in) ::   kt   ! time step
[4292]159      !
[5836]160      INTEGER  ::   ji, jj, jk   ! dummy loop indices
161      REAL(wp) ::   z1_2dt       ! local scalars
[4292]162      REAL(wp), POINTER, DIMENSION(:,:  ) ::  z2d
163      REAL(wp), POINTER, DIMENSION(:,:,:) ::  z3d, zhdiv
164      !!----------------------------------------------------------------------
165      !
[5836]166      IF( nn_timing == 1 )   CALL timing_start('wzv')
167      !
[4292]168      IF( kt == nit000 ) THEN
169         IF(lwp) WRITE(numout,*)
170         IF(lwp) WRITE(numout,*) 'wzv : now vertical velocity '
171         IF(lwp) WRITE(numout,*) '~~~~~ '
172         !
[7753]173         wn(:,:,jpk) = 0._wp                  ! bottom boundary condition: w=0 (set once for all)
[4292]174      ENDIF
175      !                                           !------------------------------!
176      !                                           !     Now Vertical Velocity    !
177      !                                           !------------------------------!
178      z1_2dt = 1. / ( 2. * rdt )                         ! set time step size (Euler/Leapfrog)
179      IF( neuler == 0 .AND. kt == nit000 )   z1_2dt = 1. / rdt
180      !
181      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN      ! z_tilde and layer cases
182         CALL wrk_alloc( jpi, jpj, jpk, zhdiv ) 
183         !
184         DO jk = 1, jpkm1
185            ! horizontal divergence of thickness diffusion transport ( velocity multiplied by e3t)
[4338]186            ! - ML - note: computation already done in dom_vvl_sf_nxt. Could be optimized (not critical and clearer this way)
[4292]187            DO jj = 2, jpjm1
188               DO ji = fs_2, fs_jpim1   ! vector opt.
[5836]189                  zhdiv(ji,jj,jk) = r1_e1e2t(ji,jj) * ( un_td(ji,jj,jk) - un_td(ji-1,jj,jk) + vn_td(ji,jj,jk) - vn_td(ji,jj-1,jk) )
[4292]190               END DO
[592]191            END DO
192         END DO
[4292]193         CALL lbc_lnk(zhdiv, 'T', 1.)  ! - ML - Perhaps not necessary: not used for horizontal "connexions"
194         !                             ! Is it problematic to have a wrong vertical velocity in boundary cells?
195         !                             ! Same question holds for hdivn. Perhaps just for security
196         DO jk = jpkm1, 1, -1                       ! integrate from the bottom the hor. divergence
197            ! computation of w
[7753]198            wn(:,:,jk) = wn(:,:,jk+1) - (  e3t_n(:,:,jk) * hdivn(:,:,jk) + zhdiv(:,:,jk)    &
199               &                         + z1_2dt * ( e3t_a(:,:,jk) - e3t_b(:,:,jk) )     ) * tmask(:,:,jk)
[4292]200         END DO
201         !          IF( ln_vvl_layer ) wn(:,:,:) = 0.e0
202         CALL wrk_dealloc( jpi, jpj, jpk, zhdiv ) 
203      ELSE   ! z_star and linear free surface cases
204         DO jk = jpkm1, 1, -1                       ! integrate from the bottom the hor. divergence
[7753]205            ! computation of w
206            wn(:,:,jk) = wn(:,:,jk+1) - (  e3t_n(:,:,jk) * hdivn(:,:,jk)                 &
207               &                         + z1_2dt * ( e3t_a(:,:,jk) - e3t_b(:,:,jk) )  ) * tmask(:,:,jk)
[4292]208         END DO
[1438]209      ENDIF
[592]210
[7646]211      IF( ln_bdy ) THEN
[4327]212         DO jk = 1, jpkm1
[7753]213            wn(:,:,jk) = wn(:,:,jk) * bdytmask(:,:)
[4327]214         END DO
215      ENDIF
[4292]216      !
217      IF( nn_timing == 1 )  CALL timing_stop('wzv')
[5836]218      !
219   END SUBROUTINE wzv
[592]220
221
[4292]222   SUBROUTINE ssh_swp( kt )
[1438]223      !!----------------------------------------------------------------------
224      !!                    ***  ROUTINE ssh_nxt  ***
225      !!
226      !! ** Purpose :   achieve the sea surface  height time stepping by
227      !!              applying Asselin time filter and swapping the arrays
[4292]228      !!              ssha  already computed in ssh_nxt 
[1438]229      !!
[2528]230      !! ** Method  : - apply Asselin time fiter to now ssh (excluding the forcing
231      !!              from the filter, see Leclair and Madec 2010) and swap :
232      !!                sshn = ssha + atfp * ( sshb -2 sshn + ssha )
233      !!                            - atfp * rdt * ( emp_b - emp ) / rau0
234      !!                sshn = ssha
[1438]235      !!
236      !! ** action  : - sshb, sshn   : before & now sea surface height
237      !!                               ready for the next time step
[2528]238      !!
239      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
[1438]240      !!----------------------------------------------------------------------
[2528]241      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
[6140]242      !
243      REAL(wp) ::   zcoef   ! local scalar
[1438]244      !!----------------------------------------------------------------------
[3294]245      !
[4292]246      IF( nn_timing == 1 )  CALL timing_start('ssh_swp')
[3294]247      !
[1438]248      IF( kt == nit000 ) THEN
249         IF(lwp) WRITE(numout,*)
[4292]250         IF(lwp) WRITE(numout,*) 'ssh_swp : Asselin time filter and swap of sea surface height'
[1438]251         IF(lwp) WRITE(numout,*) '~~~~~~~ '
252      ENDIF
[6140]253      !              !==  Euler time-stepping: no filter, just swap  ==!
254      IF(  ( neuler == 0 .AND. kt == nit000 ) .OR.    &
255         & ( ln_bt_fw    .AND. ln_dynspg_ts )      ) THEN
[7753]256         sshb(:,:) = sshn(:,:)                              ! before <-- now
257         sshn(:,:) = ssha(:,:)                              ! now    <-- after  (before already = now)
[5836]258         !
[6140]259      ELSE           !==  Leap-Frog time-stepping: Asselin filter + swap  ==!
260         !                                                  ! before <-- now filtered
[7753]261         sshb(:,:) = sshn(:,:) + atfp * ( sshb(:,:) - 2 * sshn(:,:) + ssha(:,:) )
[6140]262         IF( .NOT.ln_linssh ) THEN                          ! before <-- with forcing removed
263            zcoef = atfp * rdt * r1_rau0
[7753]264            sshb(:,:) = sshb(:,:) - zcoef * (     emp_b(:,:) - emp   (:,:)   &
265               &                             -    rnf_b(:,:) + rnf   (:,:)   &
266               &                             + fwfisf_b(:,:) - fwfisf(:,:)   ) * ssmask(:,:)
[6140]267         ENDIF
[7753]268         sshn(:,:) = ssha(:,:)                              ! now <-- after
[1438]269      ENDIF
270      !
[2528]271      IF(ln_ctl)   CALL prt_ctl( tab2d_1=sshb, clinfo1=' sshb  - : ', mask1=tmask, ovlap=1 )
272      !
[6140]273      IF( nn_timing == 1 )   CALL timing_stop('ssh_swp')
[3294]274      !
[4292]275   END SUBROUTINE ssh_swp
[3]276
277   !!======================================================================
[1565]278END MODULE sshwzv
Note: See TracBrowser for help on using the repository browser.