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 branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC/DYN – NEMO

source: branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/OPA_SRC/DYN/sshwzv.F90 @ 7525

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

changes after review

  • Property svn:keywords set to Id
File size: 14.8 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        !
25   USE bdy_par        !
26   USE bdydyn2d       ! bdy_ssh routine
27#if defined key_agrif
28   USE agrif_opa_interp
29#endif
30#if defined key_asminc   
31   USE   asminc       ! Assimilation increment
32#endif
33   !
34   USE in_out_manager ! I/O manager
35   USE restart        ! only for lrst_oce
36   USE prtctl         ! Print control
37   USE lbclnk         ! ocean lateral boundary condition (or mpp link)
38   USE lib_mpp        ! MPP library
39   USE wrk_nemo       ! Memory Allocation
40   USE timing         ! Timing
41   USE wet_dry         ! Wetting/Drying flux limting
42
43   IMPLICIT NONE
44   PRIVATE
45
46   PUBLIC   ssh_nxt    ! called by step.F90
47   PUBLIC   wzv        ! called by step.F90
48   PUBLIC   ssh_swp    ! called by step.F90
49
50   !! * Substitutions
51#  include "vectopt_loop_substitute.h90"
52   !!----------------------------------------------------------------------
53   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
54   !! $Id$
55   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
56   !!----------------------------------------------------------------------
57CONTAINS
58
59   SUBROUTINE ssh_nxt( kt )
60      !!----------------------------------------------------------------------
61      !!                ***  ROUTINE ssh_nxt  ***
62      !!                   
63      !! ** Purpose :   compute the after ssh (ssha)
64      !!
65      !! ** Method  : - Using the incompressibility hypothesis, the ssh increment
66      !!      is computed by integrating the horizontal divergence and multiply by
67      !!      by the time step.
68      !!
69      !! ** action  :   ssha, after sea surface height
70      !!
71      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
72      !!----------------------------------------------------------------------
73      INTEGER, INTENT(in) ::   kt   ! time step
74      !
75      INTEGER  ::   jk, jj, ji            ! dummy loop indice
76      REAL(wp) ::   z2dt, zcoef   ! local scalars
77      REAL(wp), POINTER, DIMENSION(:,:  ) ::   zhdiv   ! 2D workspace
78      !!----------------------------------------------------------------------
79      !
80      IF( nn_timing == 1 )   CALL timing_start('ssh_nxt')
81      !
82      CALL wrk_alloc( jpi,jpj,   zhdiv ) 
83      !
84      IF( kt == nit000 ) THEN
85         IF(lwp) WRITE(numout,*)
86         IF(lwp) WRITE(numout,*) 'ssh_nxt : after sea surface height'
87         IF(lwp) WRITE(numout,*) '~~~~~~~ '
88      ENDIF
89      !
90      CALL div_hor( kt )                              ! Horizontal divergence
91      !
92      z2dt = 2._wp * rdt                              ! set time step size (Euler/Leapfrog)
93      IF( neuler == 0 .AND. kt == nit000 )   z2dt = rdt
94
95      !                                           !------------------------------!
96      !                                           !   After Sea Surface Height   !
97      !                                           !------------------------------!
98!$OMP PARALLEL
99!$OMP DO schedule(static) private(jj, ji)
100      DO jj = 1, jpj
101         DO ji = 1, jpi
102            zhdiv(ji,jj) = 0._wp
103         END DO
104      END DO           
105      DO jk = 1, jpkm1                                 ! Horizontal divergence of barotropic transports
106!$OMP DO schedule(static) private(jj, ji)
107         DO jj = 1, jpj
108            DO ji = 1, jpi   ! vector opt.
109               zhdiv(ji,jj) = zhdiv(ji,jj) + e3t_n(ji,jj,jk) * hdivn(ji,jj,jk)
110            END DO
111         END DO           
112      END DO
113!$OMP END PARALLEL
114      !                                                ! Sea surface elevation time stepping
115      ! In time-split case we need a first guess of the ssh after (using the baroclinic timestep) in order to
116      ! compute the vertical velocity which can be used to compute the non-linear terms of the momentum equations.
117      !
118      zcoef = 0.5_wp * r1_rau0
119
120      IF(ln_wd) CALL wad_lmt(sshb, zcoef * (emp_b(:,:) + emp(:,:)), z2dt)
121!$OMP PARALLEL DO schedule(static) private(jj, ji)
122         DO jj = 1, jpj
123            DO ji = 1, jpi
124               ssha(ji,jj) = (  sshb(ji,jj) - z2dt * ( zcoef * ( emp_b(ji,jj) + emp(ji,jj) ) + zhdiv(ji,jj) )  ) * ssmask(ji,jj)
125            END DO
126         END DO           
127      IF ( .NOT.ln_dynspg_ts ) THEN
128         ! These lines are not necessary with time splitting since
129         ! boundary condition on sea level is set during ts loop
130# if defined key_agrif
131         CALL agrif_ssh( kt )
132# endif
133# if defined key_bdy
134         IF( lk_bdy ) THEN
135            CALL lbc_lnk( ssha, 'T', 1. )    ! Not sure that's necessary
136            CALL bdy_ssh( ssha )             ! Duplicate sea level across open boundaries
137         ENDIF
138# endif
139      ENDIF
140
141#if defined key_asminc
142      IF( lk_asminc .AND. ln_sshinc .AND. ln_asmiau ) THEN     ! Include the IAU weighted SSH increment
143         CALL ssh_asm_inc( kt )
144!$OMP PARALLEL DO schedule(static) private(jj, ji)
145         DO jj = 1, jpj
146            DO ji = 1, jpi
147               ssha(ji,jj) = ssha(ji,jj) + z2dt * ssh_iau(ji,jj)
148            END DO
149         END DO           
150      ENDIF
151#endif
152      !                                           !------------------------------!
153      !                                           !           outputs            !
154      !                                           !------------------------------!
155      !
156      IF(ln_ctl)   CALL prt_ctl( tab2d_1=ssha, clinfo1=' ssha  - : ', mask1=tmask, ovlap=1 )
157      !
158      CALL wrk_dealloc( jpi, jpj, zhdiv ) 
159      !
160      IF( nn_timing == 1 )  CALL timing_stop('ssh_nxt')
161      !
162   END SUBROUTINE ssh_nxt
163
164   
165   SUBROUTINE wzv( kt )
166      !!----------------------------------------------------------------------
167      !!                ***  ROUTINE wzv  ***
168      !!                   
169      !! ** Purpose :   compute the now vertical velocity
170      !!
171      !! ** Method  : - Using the incompressibility hypothesis, the vertical
172      !!      velocity is computed by integrating the horizontal divergence 
173      !!      from the bottom to the surface minus the scale factor evolution.
174      !!        The boundary conditions are w=0 at the bottom (no flux) and.
175      !!
176      !! ** action  :   wn      : now vertical velocity
177      !!
178      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
179      !!----------------------------------------------------------------------
180      INTEGER, INTENT(in) ::   kt   ! time step
181      !
182      INTEGER  ::   ji, jj, jk   ! dummy loop indices
183      REAL(wp) ::   z1_2dt       ! local scalars
184      REAL(wp), POINTER, DIMENSION(:,:  ) ::  z2d
185      REAL(wp), POINTER, DIMENSION(:,:,:) ::  z3d, zhdiv
186      !!----------------------------------------------------------------------
187      !
188      IF( nn_timing == 1 )   CALL timing_start('wzv')
189      !
190      IF( kt == nit000 ) THEN
191         IF(lwp) WRITE(numout,*)
192         IF(lwp) WRITE(numout,*) 'wzv : now vertical velocity '
193         IF(lwp) WRITE(numout,*) '~~~~~ '
194         !
195!$OMP PARALLEL DO schedule(static) private(jj, ji)
196         DO jj = 1, jpj
197            DO ji = 1, jpi
198               wn(ji,jj,jpk) = 0._wp                  ! bottom boundary condition: w=0 (set once for all)
199            END DO
200         END DO
201      ENDIF
202      !                                           !------------------------------!
203      !                                           !     Now Vertical Velocity    !
204      !                                           !------------------------------!
205      z1_2dt = 1. / ( 2. * rdt )                         ! set time step size (Euler/Leapfrog)
206      IF( neuler == 0 .AND. kt == nit000 )   z1_2dt = 1. / rdt
207      !
208      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN      ! z_tilde and layer cases
209         CALL wrk_alloc( jpi, jpj, jpk, zhdiv ) 
210!$OMP PARALLEL DO schedule(static) private(jk, jj, ji)
211         DO jk = 1, jpkm1
212            ! horizontal divergence of thickness diffusion transport ( velocity multiplied by e3t)
213            ! - ML - note: computation already done in dom_vvl_sf_nxt. Could be optimized (not critical and clearer this way)
214            DO jj = 2, jpjm1
215               DO ji = fs_2, fs_jpim1   ! vector opt.
216                  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) )
217               END DO
218            END DO
219         END DO
220         CALL lbc_lnk(zhdiv, 'T', 1.)  ! - ML - Perhaps not necessary: not used for horizontal "connexions"
221         !                             ! Is it problematic to have a wrong vertical velocity in boundary cells?
222         !                             ! Same question holds for hdivn. Perhaps just for security
223         DO jk = jpkm1, 1, -1                       ! integrate from the bottom the hor. divergence
224            ! computation of w
225!$OMP PARALLEL DO schedule(static) private(jj, ji)
226            DO jj = 1, jpj
227               DO ji = 1, jpi   ! vector opt.
228                  wn(ji,jj,jk) = wn(ji,jj,jk+1) - ( e3t_n(ji,jj,jk) * hdivn(ji,jj,jk) + zhdiv(ji,jj,jk)    &
229                  &                         + z1_2dt * ( e3t_a(ji,jj,jk) - e3t_b(ji,jj,jk) )     ) * tmask(ji,jj,jk)
230               END DO
231            END DO
232         END DO
233         !          IF( ln_vvl_layer ) wn(:,:,:) = 0.e0
234         CALL wrk_dealloc( jpi, jpj, jpk, zhdiv ) 
235      ELSE   ! z_star and linear free surface cases
236         DO jk = jpkm1, 1, -1                       ! integrate from the bottom the hor. divergence
237!$OMP PARALLEL DO schedule(static) private(jj, ji)
238            DO jj = 1, jpj
239               DO ji = 1, jpi   ! vector opt.
240                  ! computation of w
241                  wn(ji,jj,jk) = wn(ji,jj,jk+1) - (  e3t_n(ji,jj,jk) * hdivn(ji,jj,jk)                 &
242                  &                         + z1_2dt * ( e3t_a(ji,jj,jk) - e3t_b(ji,jj,jk) )  ) * tmask(ji,jj,jk)
243                END DO
244            END DO
245         END DO
246      ENDIF
247
248#if defined key_bdy
249      IF( lk_bdy ) THEN
250!$OMP PARALLEL DO schedule(static) private(jk)
251         DO jk = 1, jpkm1
252            wn(:,:,jk) = wn(:,:,jk) * bdytmask(:,:)
253         END DO
254      ENDIF
255#endif
256      !
257      IF( nn_timing == 1 )  CALL timing_stop('wzv')
258      !
259   END SUBROUTINE wzv
260
261
262   SUBROUTINE ssh_swp( kt )
263      !!----------------------------------------------------------------------
264      !!                    ***  ROUTINE ssh_nxt  ***
265      !!
266      !! ** Purpose :   achieve the sea surface  height time stepping by
267      !!              applying Asselin time filter and swapping the arrays
268      !!              ssha  already computed in ssh_nxt 
269      !!
270      !! ** Method  : - apply Asselin time fiter to now ssh (excluding the forcing
271      !!              from the filter, see Leclair and Madec 2010) and swap :
272      !!                sshn = ssha + atfp * ( sshb -2 sshn + ssha )
273      !!                            - atfp * rdt * ( emp_b - emp ) / rau0
274      !!                sshn = ssha
275      !!
276      !! ** action  : - sshb, sshn   : before & now sea surface height
277      !!                               ready for the next time step
278      !!
279      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
280      !!----------------------------------------------------------------------
281      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
282      !
283      INTEGER  ::   ji, jj, jk   ! dummy loop indices
284      REAL(wp) ::   zcoef   ! local scalar
285      !!----------------------------------------------------------------------
286      !
287      IF( nn_timing == 1 )  CALL timing_start('ssh_swp')
288      !
289      IF( kt == nit000 ) THEN
290         IF(lwp) WRITE(numout,*)
291         IF(lwp) WRITE(numout,*) 'ssh_swp : Asselin time filter and swap of sea surface height'
292         IF(lwp) WRITE(numout,*) '~~~~~~~ '
293      ENDIF
294      !              !==  Euler time-stepping: no filter, just swap  ==!
295      IF(  ( neuler == 0 .AND. kt == nit000 ) .OR.    &
296         & ( ln_bt_fw    .AND. ln_dynspg_ts )      ) THEN 
297!$OMP PARALLEL DO schedule(static) private(jj, ji)
298         DO jj = 1, jpj
299            DO ji = 1, jpi
300               sshb(ji,jj) = sshn(ji,jj)                              ! before <-- now
301               sshn(ji,jj) = ssha(ji,jj)                              ! now    <-- after  (before already = now)
302            END DO
303         END DO           
304         !
305      ELSE           !==  Leap-Frog time-stepping: Asselin filter + swap  ==!
306         !                                                  ! before <-- now filtered
307!$OMP PARALLEL DO schedule(static) private(jj, ji)
308         DO jj = 1, jpj
309            DO ji = 1, jpi
310               sshb(ji,jj) = sshn(ji,jj) + atfp * ( sshb(ji,jj) - 2 * sshn(ji,jj) + ssha(ji,jj) )
311            END DO
312         END DO           
313         IF( .NOT.ln_linssh ) THEN                          ! before <-- with forcing removed
314            zcoef = atfp * rdt * r1_rau0
315!$OMP PARALLEL DO schedule(static) private(jj, ji)
316            DO jj = 1, jpj
317               DO ji = 1, jpi
318                  sshb(ji,jj) = sshb(ji,jj) - zcoef * (     emp_b(ji,jj) - emp   (ji,jj)   &
319                  &                             -    rnf_b(ji,jj) + rnf   (ji,jj)   &
320                  &                             + fwfisf_b(ji,jj) - fwfisf(ji,jj)   ) * ssmask(ji,jj)
321               END DO
322            END DO           
323         ENDIF
324!$OMP PARALLEL DO schedule(static) private(jj, ji)
325         DO jj = 1, jpj
326            DO ji = 1, jpi
327               sshn(ji,jj) = ssha(ji,jj)                              ! now <-- after
328            END DO
329         END DO           
330      ENDIF
331      !
332      IF(ln_ctl)   CALL prt_ctl( tab2d_1=sshb, clinfo1=' sshb  - : ', mask1=tmask, ovlap=1 )
333      !
334      IF( nn_timing == 1 )   CALL timing_stop('ssh_swp')
335      !
336   END SUBROUTINE ssh_swp
337
338   !!======================================================================
339END MODULE sshwzv
Note: See TracBrowser for help on using the repository browser.