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 @ 7037

Last change on this file since 7037 was 7037, checked in by mocavero, 8 years ago

ORCA2_LIM_PISCES hybrid version update

  • Property svn:keywords set to Id
File size: 13.9 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 WORKSHARE
100      zhdiv(:,:) = 0._wp
101!$OMP END WORKSHARE
102      DO jk = 1, jpkm1                                 ! Horizontal divergence of barotropic transports
103!$OMP DO schedule(static) private(jj, ji)
104            DO jj = 1, jpj
105               DO ji = 1, jpi   ! vector opt.
106        zhdiv(ji,jj) = zhdiv(ji,jj) + e3t_n(ji,jj,jk) * hdivn(ji,jj,jk)
107               END DO
108            END DO           
109      END DO
110!$OMP END PARALLEL
111      !                                                ! Sea surface elevation time stepping
112      ! In time-split case we need a first guess of the ssh after (using the baroclinic timestep) in order to
113      ! compute the vertical velocity which can be used to compute the non-linear terms of the momentum equations.
114      !
115      zcoef = 0.5_wp * r1_rau0
116
117      IF(ln_wd) CALL wad_lmt(sshb, zcoef * (emp_b(:,:) + emp(:,:)), z2dt)
118!$OMP PARALLEL WORKSHARE
119      ssha(:,:) = (  sshb(:,:) - z2dt * ( zcoef * ( emp_b(:,:) + emp(:,:) ) + zhdiv(:,:) )  ) * ssmask(:,:)
120!$OMP END PARALLEL WORKSHARE
121      IF ( .NOT.ln_dynspg_ts ) THEN
122         ! These lines are not necessary with time splitting since
123         ! boundary condition on sea level is set during ts loop
124# if defined key_agrif
125         CALL agrif_ssh( kt )
126# endif
127# if defined key_bdy
128         IF( lk_bdy ) THEN
129            CALL lbc_lnk( ssha, 'T', 1. )    ! Not sure that's necessary
130            CALL bdy_ssh( ssha )             ! Duplicate sea level across open boundaries
131         ENDIF
132# endif
133      ENDIF
134
135#if defined key_asminc
136      IF( lk_asminc .AND. ln_sshinc .AND. ln_asmiau ) THEN     ! Include the IAU weighted SSH increment
137         CALL ssh_asm_inc( kt )
138!$OMP PARALLEL WORKSHARE
139         ssha(:,:) = ssha(:,:) + z2dt * ssh_iau(:,:)
140!$OMP END PARALLEL WORKSHARE
141      ENDIF
142#endif
143      !                                           !------------------------------!
144      !                                           !           outputs            !
145      !                                           !------------------------------!
146      !
147      IF(ln_ctl)   CALL prt_ctl( tab2d_1=ssha, clinfo1=' ssha  - : ', mask1=tmask, ovlap=1 )
148      !
149      CALL wrk_dealloc( jpi, jpj, zhdiv ) 
150      !
151      IF( nn_timing == 1 )  CALL timing_stop('ssh_nxt')
152      !
153   END SUBROUTINE ssh_nxt
154
155   
156   SUBROUTINE wzv( kt )
157      !!----------------------------------------------------------------------
158      !!                ***  ROUTINE wzv  ***
159      !!                   
160      !! ** Purpose :   compute the now vertical velocity
161      !!
162      !! ** Method  : - Using the incompressibility hypothesis, the vertical
163      !!      velocity is computed by integrating the horizontal divergence 
164      !!      from the bottom to the surface minus the scale factor evolution.
165      !!        The boundary conditions are w=0 at the bottom (no flux) and.
166      !!
167      !! ** action  :   wn      : now vertical velocity
168      !!
169      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
170      !!----------------------------------------------------------------------
171      INTEGER, INTENT(in) ::   kt   ! time step
172      !
173      INTEGER  ::   ji, jj, jk   ! dummy loop indices
174      REAL(wp) ::   z1_2dt       ! local scalars
175      REAL(wp), POINTER, DIMENSION(:,:  ) ::  z2d
176      REAL(wp), POINTER, DIMENSION(:,:,:) ::  z3d, zhdiv
177      !!----------------------------------------------------------------------
178      !
179      IF( nn_timing == 1 )   CALL timing_start('wzv')
180      !
181      IF( kt == nit000 ) THEN
182         IF(lwp) WRITE(numout,*)
183         IF(lwp) WRITE(numout,*) 'wzv : now vertical velocity '
184         IF(lwp) WRITE(numout,*) '~~~~~ '
185         !
186         wn(:,:,jpk) = 0._wp                  ! bottom boundary condition: w=0 (set once for all)
187      ENDIF
188      !                                           !------------------------------!
189      !                                           !     Now Vertical Velocity    !
190      !                                           !------------------------------!
191      z1_2dt = 1. / ( 2. * rdt )                         ! set time step size (Euler/Leapfrog)
192      IF( neuler == 0 .AND. kt == nit000 )   z1_2dt = 1. / rdt
193      !
194      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN      ! z_tilde and layer cases
195         CALL wrk_alloc( jpi, jpj, jpk, zhdiv ) 
196!$OMP PARALLEL DO schedule(static) private(jk, jj, ji)
197         DO jk = 1, jpkm1
198            ! horizontal divergence of thickness diffusion transport ( velocity multiplied by e3t)
199            ! - ML - note: computation already done in dom_vvl_sf_nxt. Could be optimized (not critical and clearer this way)
200            DO jj = 2, jpjm1
201               DO ji = fs_2, fs_jpim1   ! vector opt.
202                  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) )
203               END DO
204            END DO
205         END DO
206         CALL lbc_lnk(zhdiv, 'T', 1.)  ! - ML - Perhaps not necessary: not used for horizontal "connexions"
207         !                             ! Is it problematic to have a wrong vertical velocity in boundary cells?
208         !                             ! Same question holds for hdivn. Perhaps just for security
209         DO jk = jpkm1, 1, -1                       ! integrate from the bottom the hor. divergence
210            ! computation of w
211!$OMP PARALLEL DO schedule(static) private(jj, ji)
212            DO jj = 1, jpj
213               DO ji = 1, jpi   ! vector opt.
214            wn(ji,jj,jk) = wn(ji,jj,jk+1) - ( e3t_n(ji,jj,jk) * hdivn(ji,jj,jk) + zhdiv(ji,jj,jk)    &
215               &                         + z1_2dt * ( e3t_a(ji,jj,jk) - e3t_b(ji,jj,jk) )     ) * tmask(ji,jj,jk)
216               END DO
217            END DO
218         END DO
219         !          IF( ln_vvl_layer ) wn(:,:,:) = 0.e0
220         CALL wrk_dealloc( jpi, jpj, jpk, zhdiv ) 
221      ELSE   ! z_star and linear free surface cases
222         DO jk = jpkm1, 1, -1                       ! integrate from the bottom the hor. divergence
223!$OMP PARALLEL DO schedule(static) private(jj, ji)
224            DO jj = 1, jpj
225               DO ji = 1, jpi   ! vector opt.
226            ! computation of w
227            wn(ji,jj,jk) = wn(ji,jj,jk+1) - (  e3t_n(ji,jj,jk) * hdivn(ji,jj,jk)                 &
228               &                         + z1_2dt * ( e3t_a(ji,jj,jk) - e3t_b(ji,jj,jk) )  ) * tmask(ji,jj,jk)
229                END DO
230            END DO
231         END DO
232      ENDIF
233
234#if defined key_bdy
235      IF( lk_bdy ) THEN
236!$OMP PARALLEL DO schedule(static) private(jk)
237         DO jk = 1, jpkm1
238            wn(:,:,jk) = wn(:,:,jk) * bdytmask(:,:)
239         END DO
240      ENDIF
241#endif
242      !
243      IF( nn_timing == 1 )  CALL timing_stop('wzv')
244      !
245   END SUBROUTINE wzv
246
247
248   SUBROUTINE ssh_swp( kt )
249      !!----------------------------------------------------------------------
250      !!                    ***  ROUTINE ssh_nxt  ***
251      !!
252      !! ** Purpose :   achieve the sea surface  height time stepping by
253      !!              applying Asselin time filter and swapping the arrays
254      !!              ssha  already computed in ssh_nxt 
255      !!
256      !! ** Method  : - apply Asselin time fiter to now ssh (excluding the forcing
257      !!              from the filter, see Leclair and Madec 2010) and swap :
258      !!                sshn = ssha + atfp * ( sshb -2 sshn + ssha )
259      !!                            - atfp * rdt * ( emp_b - emp ) / rau0
260      !!                sshn = ssha
261      !!
262      !! ** action  : - sshb, sshn   : before & now sea surface height
263      !!                               ready for the next time step
264      !!
265      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
266      !!----------------------------------------------------------------------
267      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
268      !
269      REAL(wp) ::   zcoef   ! local scalar
270      !!----------------------------------------------------------------------
271      !
272      IF( nn_timing == 1 )  CALL timing_start('ssh_swp')
273      !
274      IF( kt == nit000 ) THEN
275         IF(lwp) WRITE(numout,*)
276         IF(lwp) WRITE(numout,*) 'ssh_swp : Asselin time filter and swap of sea surface height'
277         IF(lwp) WRITE(numout,*) '~~~~~~~ '
278      ENDIF
279      !              !==  Euler time-stepping: no filter, just swap  ==!
280      IF(  ( neuler == 0 .AND. kt == nit000 ) .OR.    &
281         & ( ln_bt_fw    .AND. ln_dynspg_ts )      ) THEN 
282!$OMP PARALLEL WORKSHARE
283         sshb(:,:) = sshn(:,:)                              ! before <-- now
284         sshn(:,:) = ssha(:,:)                              ! now    <-- after  (before already = now)
285!$OMP END PARALLEL WORKSHARE
286         !
287      ELSE           !==  Leap-Frog time-stepping: Asselin filter + swap  ==!
288         !                                                  ! before <-- now filtered
289!$OMP PARALLEL WORKSHARE
290         sshb(:,:) = sshn(:,:) + atfp * ( sshb(:,:) - 2 * sshn(:,:) + ssha(:,:) )
291!$OMP END PARALLEL WORKSHARE
292         IF( .NOT.ln_linssh ) THEN                          ! before <-- with forcing removed
293            zcoef = atfp * rdt * r1_rau0
294!$OMP PARALLEL WORKSHARE
295            sshb(:,:) = sshb(:,:) - zcoef * (     emp_b(:,:) - emp   (:,:)   &
296               &                             -    rnf_b(:,:) + rnf   (:,:)   &
297               &                             + fwfisf_b(:,:) - fwfisf(:,:)   ) * ssmask(:,:)
298!$OMP END PARALLEL WORKSHARE
299         ENDIF
300!$OMP PARALLEL WORKSHARE
301         sshn(:,:) = ssha(:,:)                              ! now <-- after
302!$OMP END PARALLEL WORKSHARE
303      ENDIF
304      !
305      IF(ln_ctl)   CALL prt_ctl( tab2d_1=sshb, clinfo1=' sshb  - : ', mask1=tmask, ovlap=1 )
306      !
307      IF( nn_timing == 1 )   CALL timing_stop('ssh_swp')
308      !
309   END SUBROUTINE ssh_swp
310
311   !!======================================================================
312END MODULE sshwzv
Note: See TracBrowser for help on using the repository browser.