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 NEMO/branches/UKMO/dev_r10037_GPU/src/OCE/DYN – NEMO

source: NEMO/branches/UKMO/dev_r10037_GPU/src/OCE/DYN/sshwzv.F90 @ 11467

Last change on this file since 11467 was 11467, checked in by andmirek, 5 years ago

Ticket #2197 allocate arrays at the beggining of the run

File size: 12.3 KB
Line 
1MODULE sshwzv   
2   !!==============================================================================
3   !!                       ***  MODULE  sshwzv  ***
4   !! Ocean dynamics : sea surface height and vertical velocity
5   !!==============================================================================
6   !! History :  3.1  !  2009-02  (G. Madec, M. Leclair)  Original code
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
10   !!            3.3  !  2011-10  (M. Leclair) split former ssh_wzv routine and remove all vvl related work
11   !!----------------------------------------------------------------------
12
13   !!----------------------------------------------------------------------
14   !!   ssh_nxt       : after ssh
15   !!   ssh_swp       : filter ans swap the ssh arrays
16   !!   wzv           : compute now vertical velocity
17   !!----------------------------------------------------------------------
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
24   USE bdy_oce , ONLY : ln_bdy, bdytmask   ! Open BounDarY
25   USE bdydyn2d       ! bdy_ssh routine
26#if defined key_agrif
27   USE agrif_oce_interp
28#endif
29   !
30   USE in_out_manager ! I/O manager
31   USE restart        ! only for lrst_oce
32   USE prtctl         ! Print control
33   USE lbclnk         ! ocean lateral boundary condition (or mpp link)
34   USE lib_mpp        ! MPP library
35   USE timing         ! Timing
36   USE wet_dry        ! Wetting/Drying flux limiting
37
38   IMPLICIT NONE
39   PRIVATE
40
41   PUBLIC   ssh_nxt    ! called by step.F90
42   PUBLIC   wzv        ! called by step.F90
43   PUBLIC   ssh_swp    ! called by step.F90
44
45   !! * Substitutions
46#  include "vectopt_loop_substitute.h90"
47   !!----------------------------------------------------------------------
48   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
49   !! $Id$
50   !! Software governed by the CeCILL license (see ./LICENSE)
51   !!----------------------------------------------------------------------
52CONTAINS
53
54   SUBROUTINE ssh_nxt( kt )
55      !!----------------------------------------------------------------------
56      !!                ***  ROUTINE ssh_nxt  ***
57      !!                   
58      !! ** Purpose :   compute the after ssh (ssha)
59      !!
60      !! ** Method  : - Using the incompressibility hypothesis, the ssh increment
61      !!      is computed by integrating the horizontal divergence and multiply by
62      !!      by the time step.
63      !!
64      !! ** action  :   ssha, after sea surface height
65      !!
66      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
67      !!----------------------------------------------------------------------
68      use scoce, ONLY : zhdiv => scr2D1
69      INTEGER, INTENT(in) ::   kt   ! time step
70      !
71      INTEGER  ::   jk            ! dummy loop indice
72      REAL(wp) ::   z2dt, zcoef   ! local scalars
73      !!----------------------------------------------------------------------
74      !
75      IF( ln_timing )   CALL timing_start('ssh_nxt')
76      !
77      IF( kt == nit000 ) THEN
78         IF(lwp) WRITE(numout,*)
79         IF(lwp) WRITE(numout,*) 'ssh_nxt : after sea surface height'
80         IF(lwp) WRITE(numout,*) '~~~~~~~ '
81      ENDIF
82      !
83      z2dt = 2._wp * rdt                          ! set time step size (Euler/Leapfrog)
84      IF( neuler == 0 .AND. kt == nit000 )   z2dt = rdt
85      zcoef = 0.5_wp * r1_rau0
86
87      !                                           !------------------------------!
88      !                                           !   After Sea Surface Height   !
89      !                                           !------------------------------!
90      IF(ln_wd_il) THEN
91         CALL wad_lmt(sshb, zcoef * (emp_b(:,:) + emp(:,:)), z2dt)
92      ENDIF
93
94      CALL div_hor( kt )                               ! Horizontal divergence
95      !
96      zhdiv(:,:) = 0._wp
97      DO jk = 1, jpkm1                                 ! Horizontal divergence of barotropic transports
98        zhdiv(:,:) = zhdiv(:,:) + e3t_n(:,:,jk) * hdivn(:,:,jk)
99      END DO
100      !                                                ! Sea surface elevation time stepping
101      ! In time-split case we need a first guess of the ssh after (using the baroclinic timestep) in order to
102      ! compute the vertical velocity which can be used to compute the non-linear terms of the momentum equations.
103      !
104      ssha(:,:) = (  sshb(:,:) - z2dt * ( zcoef * ( emp_b(:,:) + emp(:,:) ) + zhdiv(:,:) )  ) * ssmask(:,:)
105      !
106#if defined key_agrif
107      CALL agrif_ssh( kt )
108#endif
109      !
110      IF ( .NOT.ln_dynspg_ts ) THEN
111         IF( ln_bdy ) THEN
112            CALL lbc_lnk( ssha, 'T', 1. )    ! Not sure that's necessary
113            CALL bdy_ssh( ssha )             ! Duplicate sea level across open boundaries
114         ENDIF
115      ENDIF
116      !                                           !------------------------------!
117      !                                           !           outputs            !
118      !                                           !------------------------------!
119      !
120      IF(ln_ctl)   CALL prt_ctl( tab2d_1=ssha, clinfo1=' ssha  - : ', mask1=tmask )
121      !
122      IF( ln_timing )   CALL timing_stop('ssh_nxt')
123      !
124   END SUBROUTINE ssh_nxt
125
126   
127   SUBROUTINE wzv( kt )
128      !!----------------------------------------------------------------------
129      !!                ***  ROUTINE wzv  ***
130      !!                   
131      !! ** Purpose :   compute the now vertical velocity
132      !!
133      !! ** Method  : - Using the incompressibility hypothesis, the vertical
134      !!      velocity is computed by integrating the horizontal divergence 
135      !!      from the bottom to the surface minus the scale factor evolution.
136      !!        The boundary conditions are w=0 at the bottom (no flux) and.
137      !!
138      !! ** action  :   wn      : now vertical velocity
139      !!
140      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
141      !!----------------------------------------------------------------------
142      USE scoce, ONLY : zhdiv => scr1
143      INTEGER, INTENT(in) ::   kt   ! time step
144      !
145      INTEGER  ::   ji, jj, jk   ! dummy loop indices
146      REAL(wp) ::   z1_2dt       ! local scalars
147      !!----------------------------------------------------------------------
148      !
149      IF( ln_timing )   CALL timing_start('wzv')
150      !
151      IF( kt == nit000 ) THEN
152         IF(lwp) WRITE(numout,*)
153         IF(lwp) WRITE(numout,*) 'wzv : now vertical velocity '
154         IF(lwp) WRITE(numout,*) '~~~~~ '
155         !
156         wn(:,:,jpk) = 0._wp                  ! bottom boundary condition: w=0 (set once for all)
157      ENDIF
158      !                                           !------------------------------!
159      !                                           !     Now Vertical Velocity    !
160      !                                           !------------------------------!
161      z1_2dt = 1. / ( 2. * rdt )                         ! set time step size (Euler/Leapfrog)
162      IF( neuler == 0 .AND. kt == nit000 )   z1_2dt = 1. / rdt
163      !
164      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN      ! z_tilde and layer cases
165         !
166         DO jk = 1, jpkm1
167            ! horizontal divergence of thickness diffusion transport ( velocity multiplied by e3t)
168            ! - ML - note: computation already done in dom_vvl_sf_nxt. Could be optimized (not critical and clearer this way)
169            DO jj = 2, jpjm1
170               DO ji = fs_2, fs_jpim1   ! vector opt.
171                  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) )
172               END DO
173            END DO
174         END DO
175         CALL lbc_lnk(zhdiv, 'T', 1.)  ! - ML - Perhaps not necessary: not used for horizontal "connexions"
176         !                             ! Is it problematic to have a wrong vertical velocity in boundary cells?
177         !                             ! Same question holds for hdivn. Perhaps just for security
178         DO jk = jpkm1, 1, -1                       ! integrate from the bottom the hor. divergence
179            ! computation of w
180            wn(:,:,jk) = wn(:,:,jk+1) - (  e3t_n(:,:,jk) * hdivn(:,:,jk) + zhdiv(:,:,jk)    &
181               &                         + z1_2dt * ( e3t_a(:,:,jk) - e3t_b(:,:,jk) )     ) * tmask(:,:,jk)
182         END DO
183         !          IF( ln_vvl_layer ) wn(:,:,:) = 0.e0
184      ELSE   ! z_star and linear free surface cases
185         DO jk = jpkm1, 1, -1                       ! integrate from the bottom the hor. divergence
186            ! computation of w
187            wn(:,:,jk) = wn(:,:,jk+1) - (  e3t_n(:,:,jk) * hdivn(:,:,jk)                 &
188               &                         + z1_2dt * ( e3t_a(:,:,jk) - e3t_b(:,:,jk) )  ) * tmask(:,:,jk)
189         END DO
190      ENDIF
191
192      IF( ln_bdy ) THEN
193         DO jk = 1, jpkm1
194            wn(:,:,jk) = wn(:,:,jk) * bdytmask(:,:)
195         END DO
196      ENDIF
197      !
198#if defined key_agrif 
199      IF( .NOT. AGRIF_Root() ) THEN
200         IF ((nbondi ==  1).OR.(nbondi == 2)) wn(nlci-1 , :     ,:) = 0.e0      ! east
201         IF ((nbondi == -1).OR.(nbondi == 2)) wn(2      , :     ,:) = 0.e0      ! west
202         IF ((nbondj ==  1).OR.(nbondj == 2)) wn(:      ,nlcj-1 ,:) = 0.e0      ! north
203         IF ((nbondj == -1).OR.(nbondj == 2)) wn(:      ,2      ,:) = 0.e0      ! south
204      ENDIF 
205#endif 
206      !
207      IF( ln_timing )   CALL timing_stop('wzv')
208      !
209   END SUBROUTINE wzv
210
211
212   SUBROUTINE ssh_swp( kt )
213      !!----------------------------------------------------------------------
214      !!                    ***  ROUTINE ssh_nxt  ***
215      !!
216      !! ** Purpose :   achieve the sea surface  height time stepping by
217      !!              applying Asselin time filter and swapping the arrays
218      !!              ssha  already computed in ssh_nxt 
219      !!
220      !! ** Method  : - apply Asselin time fiter to now ssh (excluding the forcing
221      !!              from the filter, see Leclair and Madec 2010) and swap :
222      !!                sshn = ssha + atfp * ( sshb -2 sshn + ssha )
223      !!                            - atfp * rdt * ( emp_b - emp ) / rau0
224      !!                sshn = ssha
225      !!
226      !! ** action  : - sshb, sshn   : before & now sea surface height
227      !!                               ready for the next time step
228      !!
229      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
230      !!----------------------------------------------------------------------
231      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
232      !
233      REAL(wp) ::   zcoef   ! local scalar
234      !!----------------------------------------------------------------------
235      !
236      IF( ln_timing )   CALL timing_start('ssh_swp')
237      !
238      IF( kt == nit000 ) THEN
239         IF(lwp) WRITE(numout,*)
240         IF(lwp) WRITE(numout,*) 'ssh_swp : Asselin time filter and swap of sea surface height'
241         IF(lwp) WRITE(numout,*) '~~~~~~~ '
242      ENDIF
243      !              !==  Euler time-stepping: no filter, just swap  ==!
244      IF ( neuler == 0 .AND. kt == nit000 ) THEN
245         sshn(:,:) = ssha(:,:)                              ! now    <-- after  (before already = now)
246         !
247      ELSE           !==  Leap-Frog time-stepping: Asselin filter + swap  ==!
248         !                                                  ! before <-- now filtered
249         sshb(:,:) = sshn(:,:) + atfp * ( sshb(:,:) - 2 * sshn(:,:) + ssha(:,:) )
250         IF( .NOT.ln_linssh ) THEN                          ! before <-- with forcing removed
251            zcoef = atfp * rdt * r1_rau0
252            sshb(:,:) = sshb(:,:) - zcoef * (     emp_b(:,:) - emp   (:,:)   &
253               &                             -    rnf_b(:,:) + rnf   (:,:)   &
254               &                             + fwfisf_b(:,:) - fwfisf(:,:)   ) * ssmask(:,:)
255         ENDIF
256         sshn(:,:) = ssha(:,:)                              ! now <-- after
257      ENDIF
258      !
259      IF(ln_ctl)   CALL prt_ctl( tab2d_1=sshb, clinfo1=' sshb  - : ', mask1=tmask )
260      !
261      IF( ln_timing )   CALL timing_stop('ssh_swp')
262      !
263   END SUBROUTINE ssh_swp
264
265   !!======================================================================
266END MODULE sshwzv
Note: See TracBrowser for help on using the repository browser.