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

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

update trunk with OpenMP parallelization

  • Property svn:keywords set to Id
File size: 14.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   , 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 wrk_nemo       ! Memory Allocation
39   USE timing         ! Timing
40   USE wet_dry         ! Wetting/Drying flux limting
41
42   IMPLICIT NONE
43   PRIVATE
44
45   PUBLIC   ssh_nxt    ! called by step.F90
46   PUBLIC   wzv        ! called by step.F90
47   PUBLIC   ssh_swp    ! called by step.F90
48
49   !! * Substitutions
50#  include "vectopt_loop_substitute.h90"
51   !!----------------------------------------------------------------------
52   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
53   !! $Id$
54   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
55   !!----------------------------------------------------------------------
56CONTAINS
57
58   SUBROUTINE ssh_nxt( kt )
59      !!----------------------------------------------------------------------
60      !!                ***  ROUTINE ssh_nxt  ***
61      !!                   
62      !! ** Purpose :   compute the after ssh (ssha)
63      !!
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.
67      !!
68      !! ** action  :   ssha, after sea surface height
69      !!
70      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
71      !!----------------------------------------------------------------------
72      INTEGER, INTENT(in) ::   kt   ! time step
73      !
74      INTEGER  ::   jk, jj, ji            ! dummy loop indice
75      REAL(wp) ::   z2dt, zcoef   ! local scalars
76      REAL(wp), POINTER, DIMENSION(:,:  ) ::   zhdiv   ! 2D workspace
77      !!----------------------------------------------------------------------
78      !
79      IF( nn_timing == 1 )   CALL timing_start('ssh_nxt')
80      !
81      CALL wrk_alloc( jpi,jpj,   zhdiv ) 
82      !
83      IF( kt == nit000 ) THEN
84         IF(lwp) WRITE(numout,*)
85         IF(lwp) WRITE(numout,*) 'ssh_nxt : after sea surface height'
86         IF(lwp) WRITE(numout,*) '~~~~~~~ '
87      ENDIF
88      !
89      z2dt = 2._wp * rdt                          ! set time step size (Euler/Leapfrog)
90      IF( neuler == 0 .AND. kt == nit000 )   z2dt = rdt
91      zcoef = 0.5_wp * r1_rau0
92
93      !                                           !------------------------------!
94      !                                           !   After Sea Surface Height   !
95      !                                           !------------------------------!
96      IF(ln_wd) THEN
97        CALL wad_lmt(sshb, zcoef * (emp_b(:,:) + emp(:,:)), z2dt)
98      END IF
99
100      CALL div_hor( kt )                              ! Horizontal divergence
101      !
102!$OMP PARALLEL
103!$OMP DO schedule(static) private(jj, ji)
104      DO jj = 1, jpj
105         DO ji = 1, jpi
106            zhdiv(ji,jj) = 0._wp
107         END DO
108      END DO           
109      DO jk = 1, jpkm1                                 ! Horizontal divergence of barotropic transports
110!$OMP DO schedule(static) private(jj, ji)
111         DO jj = 1, jpj
112            DO ji = 1, jpi
113               zhdiv(ji,jj) = zhdiv(ji,jj) + e3t_n(ji,jj,jk) * hdivn(ji,jj,jk)
114            END DO
115         END DO           
116      END DO
117      !                                                ! Sea surface elevation time stepping
118      ! In time-split case we need a first guess of the ssh after (using the baroclinic timestep) in order to
119      ! compute the vertical velocity which can be used to compute the non-linear terms of the momentum equations.
120      !
121!$OMP 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!$OMP END PARALLEL
128      IF ( .NOT.ln_dynspg_ts ) THEN
129         ! These lines are not necessary with time splitting since
130         ! boundary condition on sea level is set during ts loop
131# if defined key_agrif
132         CALL agrif_ssh( kt )
133# endif
134         IF( ln_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
140#if defined key_asminc
141      IF( lk_asminc .AND. ln_sshinc .AND. ln_asmiau ) THEN     ! Include the IAU weighted SSH increment
142         CALL ssh_asm_inc( kt )
143!$OMP PARALLEL DO schedule(static) private(jj, ji)
144         DO jj = 1, jpj
145            DO ji = 1, jpi
146               ssha(ji,jj) = ssha(ji,jj) + z2dt * ssh_iau(ji,jj)
147            END DO
148         END DO           
149      ENDIF
150#endif
151      !                                           !------------------------------!
152      !                                           !           outputs            !
153      !                                           !------------------------------!
154      !
155      IF(ln_ctl)   CALL prt_ctl( tab2d_1=ssha, clinfo1=' ssha  - : ', mask1=tmask, ovlap=1 )
156      !
157      CALL wrk_dealloc( jpi, jpj, zhdiv ) 
158      !
159      IF( nn_timing == 1 )  CALL timing_stop('ssh_nxt')
160      !
161   END SUBROUTINE ssh_nxt
162
163   
164   SUBROUTINE wzv( kt )
165      !!----------------------------------------------------------------------
166      !!                ***  ROUTINE wzv  ***
167      !!                   
168      !! ** Purpose :   compute the now vertical velocity
169      !!
170      !! ** Method  : - Using the incompressibility hypothesis, the vertical
171      !!      velocity is computed by integrating the horizontal divergence 
172      !!      from the bottom to the surface minus the scale factor evolution.
173      !!        The boundary conditions are w=0 at the bottom (no flux) and.
174      !!
175      !! ** action  :   wn      : now vertical velocity
176      !!
177      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
178      !!----------------------------------------------------------------------
179      INTEGER, INTENT(in) ::   kt   ! time step
180      !
181      INTEGER  ::   ji, jj, jk   ! dummy loop indices
182      REAL(wp) ::   z1_2dt       ! local scalars
183      REAL(wp), POINTER, DIMENSION(:,:  ) ::  z2d
184      REAL(wp), POINTER, DIMENSION(:,:,:) ::  z3d, zhdiv
185      !!----------------------------------------------------------------------
186      !
187      IF( nn_timing == 1 )   CALL timing_start('wzv')
188      !
189      IF( kt == nit000 ) THEN
190         IF(lwp) WRITE(numout,*)
191         IF(lwp) WRITE(numout,*) 'wzv : now vertical velocity '
192         IF(lwp) WRITE(numout,*) '~~~~~ '
193         !
194!$OMP PARALLEL DO schedule(static) private(jj, ji)
195         DO jj = 1, jpj
196            DO ji = 1, jpi
197               wn(ji,jj,jpk) = 0._wp                  ! bottom boundary condition: w=0 (set once for all)
198            END DO
199         END DO           
200      ENDIF
201      !                                           !------------------------------!
202      !                                           !     Now Vertical Velocity    !
203      !                                           !------------------------------!
204      z1_2dt = 1. / ( 2. * rdt )                         ! set time step size (Euler/Leapfrog)
205      IF( neuler == 0 .AND. kt == nit000 )   z1_2dt = 1. / rdt
206      !
207      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN      ! z_tilde and layer cases
208         CALL wrk_alloc( jpi, jpj, jpk, zhdiv ) 
209!$OMP PARALLEL DO schedule(static) private(jk, jj, ji)
210         !
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( ln_bdy ) THEN
249!$OMP PARALLEL DO schedule(static) private(jk,jj,ji)
250         DO jk = 1, jpkm1
251            DO jj = 1, jpj
252               DO ji = 1, jpi
253                  wn(ji,jj,jk) = wn(ji,jj,jk) * bdytmask(ji,jj)
254               END DO
255            END DO
256         END DO
257      ENDIF
258      !
259      IF( nn_timing == 1 )  CALL timing_stop('wzv')
260      !
261   END SUBROUTINE wzv
262
263
264   SUBROUTINE ssh_swp( kt )
265      !!----------------------------------------------------------------------
266      !!                    ***  ROUTINE ssh_nxt  ***
267      !!
268      !! ** Purpose :   achieve the sea surface  height time stepping by
269      !!              applying Asselin time filter and swapping the arrays
270      !!              ssha  already computed in ssh_nxt 
271      !!
272      !! ** Method  : - apply Asselin time fiter to now ssh (excluding the forcing
273      !!              from the filter, see Leclair and Madec 2010) and swap :
274      !!                sshn = ssha + atfp * ( sshb -2 sshn + ssha )
275      !!                            - atfp * rdt * ( emp_b - emp ) / rau0
276      !!                sshn = ssha
277      !!
278      !! ** action  : - sshb, sshn   : before & now sea surface height
279      !!                               ready for the next time step
280      !!
281      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
282      !!----------------------------------------------------------------------
283      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
284      !
285      INTEGER  ::   ji, jj, jk   ! dummy loop indices
286      REAL(wp) ::   zcoef   ! local scalar
287      !!----------------------------------------------------------------------
288      !
289      IF( nn_timing == 1 )  CALL timing_start('ssh_swp')
290      !
291      IF( kt == nit000 ) THEN
292         IF(lwp) WRITE(numout,*)
293         IF(lwp) WRITE(numout,*) 'ssh_swp : Asselin time filter and swap of sea surface height'
294         IF(lwp) WRITE(numout,*) '~~~~~~~ '
295      ENDIF
296      !              !==  Euler time-stepping: no filter, just swap  ==!
297      IF(  ( neuler == 0 .AND. kt == nit000 ) .OR.    &
298         & ( ln_bt_fw    .AND. ln_dynspg_ts )      ) THEN 
299!$OMP PARALLEL DO schedule(static) private(jj, ji)
300         DO jj = 1, jpj
301            DO ji = 1, jpi
302               sshb(ji,jj) = sshn(ji,jj)                              ! before <-- now
303               sshn(ji,jj) = ssha(ji,jj)                              ! now    <-- after  (before already = now)
304            END DO
305         END DO           
306         !
307      ELSE           !==  Leap-Frog time-stepping: Asselin filter + swap  ==!
308         !                                                  ! before <-- now filtered
309!$OMP PARALLEL DO schedule(static) private(jj, ji)
310         DO jj = 1, jpj
311            DO ji = 1, jpi
312               sshb(ji,jj) = sshn(ji,jj) + atfp * ( sshb(ji,jj) - 2 * sshn(ji,jj) + ssha(ji,jj) )
313            END DO
314         END DO           
315         IF( .NOT.ln_linssh ) THEN                          ! before <-- with forcing removed
316            zcoef = atfp * rdt * r1_rau0
317!$OMP PARALLEL DO schedule(static) private(jj, ji)
318            DO jj = 1, jpj
319               DO ji = 1, jpi
320                  sshb(ji,jj) = sshb(ji,jj) - zcoef * (     emp_b(ji,jj) - emp   (ji,jj)   &
321                  &                             -    rnf_b(ji,jj) + rnf   (ji,jj)   &
322                  &                             + fwfisf_b(ji,jj) - fwfisf(ji,jj)   ) * ssmask(ji,jj)
323               END DO
324            END DO           
325         ENDIF
326!$OMP PARALLEL DO schedule(static) private(jj, ji)
327         DO jj = 1, jpj
328            DO ji = 1, jpi
329               sshn(ji,jj) = ssha(ji,jj)                              ! now <-- after
330            END DO
331         END DO           
332      ENDIF
333      !
334      IF(ln_ctl)   CALL prt_ctl( tab2d_1=sshb, clinfo1=' sshb  - : ', mask1=tmask, ovlap=1 )
335      !
336      IF( nn_timing == 1 )   CALL timing_stop('ssh_swp')
337      !
338   END SUBROUTINE ssh_swp
339
340   !!======================================================================
341END MODULE sshwzv
Note: See TracBrowser for help on using the repository browser.