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/2017/dev_r7881_no_wrk_alloc/NEMOGCM/NEMO/OPA_SRC/DYN – NEMO

source: branches/2017/dev_r7881_no_wrk_alloc/NEMOGCM/NEMO/OPA_SRC/DYN/sshwzv.F90 @ 7910

Last change on this file since 7910 was 7910, checked in by timgraham, 7 years ago

All wrk_alloc removed

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