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.
dynhpg.F90 in trunk/NEMOGCM/NEMO/OPA_SRC/DYN – NEMO

source: trunk/NEMOGCM/NEMO/OPA_SRC/DYN/dynhpg.F90 @ 2715

Last change on this file since 2715 was 2715, checked in by rblod, 13 years ago

First attempt to put dynamic allocation on the trunk

  • Property svn:keywords set to Id
File size: 53.1 KB
Line 
1MODULE dynhpg
2   !!======================================================================
3   !!                       ***  MODULE  dynhpg  ***
4   !! Ocean dynamics:  hydrostatic pressure gradient trend
5   !!======================================================================
6   !! History :  OPA  !  1987-09  (P. Andrich, M.-A. Foujols)  hpg_zco: Original code
7   !!            5.0  !  1991-11  (G. Madec)
8   !!            7.0  !  1996-01  (G. Madec)  hpg_sco: Original code for s-coordinates
9   !!            8.0  !  1997-05  (G. Madec)  split dynber into dynkeg and dynhpg
10   !!            8.5  !  2002-07  (G. Madec)  F90: Free form and module
11   !!            8.5  !  2002-08  (A. Bozec)  hpg_zps: Original code
12   !!   NEMO     1.0  !  2005-10  (A. Beckmann, B.W. An)  various s-coordinate options
13   !!                 !         Original code for hpg_ctl, hpg_hel hpg_wdj, hpg_djc, hpg_rot
14   !!             -   !  2005-11  (G. Madec) style & small optimisation
15   !!            3.3  !  2010-10  (C. Ethe, G. Madec) reorganisation of initialisation phase
16   !!----------------------------------------------------------------------
17
18   !!----------------------------------------------------------------------
19   !!   dyn_hpg      : update the momentum trend with the now horizontal
20   !!                  gradient of the hydrostatic pressure
21   !!   dyn_hpg_init : initialisation and control of options
22   !!       hpg_zco  : z-coordinate scheme
23   !!       hpg_zps  : z-coordinate plus partial steps (interpolation)
24   !!       hpg_sco  : s-coordinate (standard jacobian formulation)
25   !!       hpg_hel  : s-coordinate (helsinki modification)
26   !!       hpg_wdj  : s-coordinate (weighted density jacobian)
27   !!       hpg_djc  : s-coordinate (Density Jacobian with Cubic polynomial)
28   !!       hpg_rot  : s-coordinate (ROTated axes scheme)
29   !!----------------------------------------------------------------------
30   USE oce             ! ocean dynamics and tracers
31   USE dom_oce         ! ocean space and time domain
32   USE phycst          ! physical constants
33   USE trdmod          ! ocean dynamics trends
34   USE trdmod_oce      ! ocean variables trends
35   USE in_out_manager  ! I/O manager
36   USE prtctl          ! Print control
37   USE lbclnk          ! lateral boundary condition
38   USE lib_mpp         ! MPP library
39
40   IMPLICIT NONE
41   PRIVATE
42
43   PUBLIC   dyn_hpg        ! routine called by step module
44   PUBLIC   dyn_hpg_init   ! routine called by opa module
45
46   !                                              !!* Namelist namdyn_hpg : hydrostatic pressure gradient
47   LOGICAL , PUBLIC ::   ln_hpg_zco    = .TRUE.    !: z-coordinate - full steps
48   LOGICAL , PUBLIC ::   ln_hpg_zps    = .FALSE.   !: z-coordinate - partial steps (interpolation)
49   LOGICAL , PUBLIC ::   ln_hpg_sco    = .FALSE.   !: s-coordinate (standard jacobian formulation)
50   LOGICAL , PUBLIC ::   ln_hpg_hel    = .FALSE.   !: s-coordinate (helsinki modification)
51   LOGICAL , PUBLIC ::   ln_hpg_wdj    = .FALSE.   !: s-coordinate (weighted density jacobian)
52   LOGICAL , PUBLIC ::   ln_hpg_djc    = .FALSE.   !: s-coordinate (Density Jacobian with Cubic polynomial)
53   LOGICAL , PUBLIC ::   ln_hpg_rot    = .FALSE.   !: s-coordinate (ROTated axes scheme)
54   REAL(wp), PUBLIC ::   rn_gamma      = 0._wp     !: weighting coefficient
55   LOGICAL , PUBLIC ::   ln_dynhpg_imp = .FALSE.   !: semi-implicite hpg flag
56
57   INTEGER  ::   nhpg  =  0   ! = 0 to 6, type of pressure gradient scheme used ! (deduced from ln_hpg_... flags)
58
59   !! * Substitutions
60#  include "domzgr_substitute.h90"
61#  include "vectopt_loop_substitute.h90"
62   !!----------------------------------------------------------------------
63   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
64   !! $Id$
65   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
66   !!----------------------------------------------------------------------
67CONTAINS
68
69   SUBROUTINE dyn_hpg( kt )
70      !!---------------------------------------------------------------------
71      !!                  ***  ROUTINE dyn_hpg  ***
72      !!
73      !! ** Method  :   Call the hydrostatic pressure gradient routine
74      !!              using the scheme defined in the namelist
75      !!   
76      !! ** Action : - Update (ua,va) with the now hydrastatic pressure trend
77      !!             - Save the trend (l_trddyn=T)
78      !!----------------------------------------------------------------------
79      USE wrk_nemo, ONLY:   wrk_in_use, wrk_not_released
80      USE wrk_nemo, ONLY:   ztrdu => wrk_3d_1 , ztrdv => wrk_3d_2   ! 3D workspace
81      !!
82      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
83      !!----------------------------------------------------------------------
84      !
85      IF( wrk_in_use(3, 1,2) ) THEN
86         CALL ctl_stop('dyn_hpg: requested workspace arrays are unavailable')   ;   RETURN
87      ENDIF
88      !
89      IF( l_trddyn ) THEN                    ! Temporary saving of ua and va trends (l_trddyn)
90         ztrdu(:,:,:) = ua(:,:,:) 
91         ztrdv(:,:,:) = va(:,:,:) 
92      ENDIF     
93      !
94      SELECT CASE ( nhpg )      ! Hydrastatic pressure gradient computation
95      CASE (  0 )   ;   CALL hpg_zco    ( kt )      ! z-coordinate
96      CASE (  1 )   ;   CALL hpg_zps    ( kt )      ! z-coordinate plus partial steps (interpolation)
97      CASE (  2 )   ;   CALL hpg_sco    ( kt )      ! s-coordinate (standard jacobian formulation)
98      CASE (  3 )   ;   CALL hpg_hel    ( kt )      ! s-coordinate (helsinki modification)
99      CASE (  4 )   ;   CALL hpg_wdj    ( kt )      ! s-coordinate (weighted density jacobian)
100      CASE (  5 )   ;   CALL hpg_djc    ( kt )      ! s-coordinate (Density Jacobian with Cubic polynomial)
101      CASE (  6 )   ;   CALL hpg_rot    ( kt )      ! s-coordinate (ROTated axes scheme)
102      END SELECT
103      !
104      IF( l_trddyn ) THEN      ! save the hydrostatic pressure gradient trends for momentum trend diagnostics
105         ztrdu(:,:,:) = ua(:,:,:) - ztrdu(:,:,:)
106         ztrdv(:,:,:) = va(:,:,:) - ztrdv(:,:,:)
107         CALL trd_mod( ztrdu, ztrdv, jpdyn_trd_hpg, 'DYN', kt )
108      ENDIF         
109      !
110      IF(ln_ctl)   CALL prt_ctl( tab3d_1=ua, clinfo1=' hpg  - Ua: ', mask1=umask,   &
111         &                       tab3d_2=va, clinfo2=       ' Va: ', mask2=vmask, clinfo3='dyn' )
112      !
113      IF( wrk_not_released(3, 1,2) )   CALL ctl_stop('dyn_hpg: failed to release workspace arrays')
114      !
115   END SUBROUTINE dyn_hpg
116
117
118   SUBROUTINE dyn_hpg_init
119      !!----------------------------------------------------------------------
120      !!                 ***  ROUTINE dyn_hpg_init  ***
121      !!
122      !! ** Purpose :   initializations for the hydrostatic pressure gradient
123      !!              computation and consistency control
124      !!
125      !! ** Action  :   Read the namelist namdyn_hpg and check the consistency
126      !!      with the type of vertical coordinate used (zco, zps, sco)
127      !!----------------------------------------------------------------------
128      INTEGER ::   ioptio = 0      ! temporary integer
129      !!
130      NAMELIST/namdyn_hpg/ ln_hpg_zco, ln_hpg_zps, ln_hpg_sco, ln_hpg_hel,    &
131         &                 ln_hpg_wdj, ln_hpg_djc, ln_hpg_rot, rn_gamma  , ln_dynhpg_imp
132      !!----------------------------------------------------------------------
133      !
134      REWIND( numnam )               ! Read Namelist namdyn_hpg
135      READ  ( numnam, namdyn_hpg )
136      !
137      IF(lwp) THEN                   ! Control print
138         WRITE(numout,*)
139         WRITE(numout,*) 'dyn_hpg_init : hydrostatic pressure gradient initialisation'
140         WRITE(numout,*) '~~~~~~~~~~~~'
141         WRITE(numout,*) '   Namelist namdyn_hpg : choice of hpg scheme'
142         WRITE(numout,*) '      z-coord. - full steps                             ln_hpg_zco    = ', ln_hpg_zco
143         WRITE(numout,*) '      z-coord. - partial steps (interpolation)          ln_hpg_zps    = ', ln_hpg_zps
144         WRITE(numout,*) '      s-coord. (standard jacobian formulation)          ln_hpg_sco    = ', ln_hpg_sco
145         WRITE(numout,*) '      s-coord. (helsinki modification)                  ln_hpg_hel    = ', ln_hpg_hel
146         WRITE(numout,*) '      s-coord. (weighted density jacobian)              ln_hpg_wdj    = ', ln_hpg_wdj
147         WRITE(numout,*) '      s-coord. (Density Jacobian: Cubic polynomial)     ln_hpg_djc    = ', ln_hpg_djc
148         WRITE(numout,*) '      s-coord. (ROTated axes scheme)                    ln_hpg_rot    = ', ln_hpg_rot
149         WRITE(numout,*) '      weighting coeff. (wdj scheme)                     rn_gamma      = ', rn_gamma
150         WRITE(numout,*) '      time stepping: centered (F) or semi-implicit (T)  ln_dynhpg_imp = ', ln_dynhpg_imp
151      ENDIF
152      !
153      IF( lk_vvl .AND. .NOT. ln_hpg_sco )   &
154         &   CALL ctl_stop('dyn_hpg_init : variable volume key_vvl require the standard jacobian formulation hpg_sco')
155      !
156      !                               ! Set nhpg from ln_hpg_... flags
157      IF( ln_hpg_zco )   nhpg = 0
158      IF( ln_hpg_zps )   nhpg = 1
159      IF( ln_hpg_sco )   nhpg = 2
160      IF( ln_hpg_hel )   nhpg = 3
161      IF( ln_hpg_wdj )   nhpg = 4
162      IF( ln_hpg_djc )   nhpg = 5
163      IF( ln_hpg_rot )   nhpg = 6
164      !
165      !                               ! Consitency check
166      ioptio = 0 
167      IF( ln_hpg_zco )   ioptio = ioptio + 1
168      IF( ln_hpg_zps )   ioptio = ioptio + 1
169      IF( ln_hpg_sco )   ioptio = ioptio + 1
170      IF( ln_hpg_hel )   ioptio = ioptio + 1
171      IF( ln_hpg_wdj )   ioptio = ioptio + 1
172      IF( ln_hpg_djc )   ioptio = ioptio + 1
173      IF( ln_hpg_rot )   ioptio = ioptio + 1
174      IF( ioptio /= 1 )   CALL ctl_stop( 'NO or several hydrostatic pressure gradient options used' )
175      !
176   END SUBROUTINE dyn_hpg_init
177
178
179   SUBROUTINE hpg_zco( kt )
180      !!---------------------------------------------------------------------
181      !!                  ***  ROUTINE hpg_zco  ***
182      !!
183      !! ** Method  :   z-coordinate case, levels are horizontal surfaces.
184      !!      The now hydrostatic pressure gradient at a given level, jk,
185      !!      is computed by taking the vertical integral of the in-situ
186      !!      density gradient along the model level from the suface to that
187      !!      level:    zhpi = grav .....
188      !!                zhpj = grav .....
189      !!      add it to the general momentum trend (ua,va).
190      !!            ua = ua - 1/e1u * zhpi
191      !!            va = va - 1/e2v * zhpj
192      !!
193      !! ** Action : - Update (ua,va) with the now hydrastatic pressure trend
194      !!----------------------------------------------------------------------
195      USE oce, ONLY:   zhpi => ta , zhpj => sa   ! (ta,sa) used as 3D workspace
196      !!
197      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
198      !!
199      INTEGER  ::   ji, jj, jk       ! dummy loop indices
200      REAL(wp) ::   zcoef0, zcoef1   ! temporary scalars
201      !!----------------------------------------------------------------------
202     
203      IF( kt == nit000 ) THEN
204         IF(lwp) WRITE(numout,*)
205         IF(lwp) WRITE(numout,*) 'dyn:hpg_zco : hydrostatic pressure gradient trend'
206         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   z-coordinate case '
207      ENDIF
208     
209      zcoef0 = - grav * 0.5_wp      ! Local constant initialization
210
211      ! Surface value
212      DO jj = 2, jpjm1
213         DO ji = fs_2, fs_jpim1   ! vector opt.
214            zcoef1 = zcoef0 * fse3w(ji,jj,1)
215            ! hydrostatic pressure gradient
216            zhpi(ji,jj,1) = zcoef1 * ( rhd(ji+1,jj,1) - rhd(ji,jj,1) ) / e1u(ji,jj)
217            zhpj(ji,jj,1) = zcoef1 * ( rhd(ji,jj+1,1) - rhd(ji,jj,1) ) / e2v(ji,jj)
218            ! add to the general momentum trend
219            ua(ji,jj,1) = ua(ji,jj,1) + zhpi(ji,jj,1)
220            va(ji,jj,1) = va(ji,jj,1) + zhpj(ji,jj,1)
221         END DO
222      END DO
223      !
224      ! interior value (2=<jk=<jpkm1)
225      DO jk = 2, jpkm1
226         DO jj = 2, jpjm1
227            DO ji = fs_2, fs_jpim1   ! vector opt.
228               zcoef1 = zcoef0 * fse3w(ji,jj,jk)
229               ! hydrostatic pressure gradient
230               zhpi(ji,jj,jk) = zhpi(ji,jj,jk-1)   &
231                  &           + zcoef1 * (  ( rhd(ji+1,jj,jk)+rhd(ji+1,jj,jk-1) )   &
232                  &                       - ( rhd(ji  ,jj,jk)+rhd(ji  ,jj,jk-1) )  ) / e1u(ji,jj)
233
234               zhpj(ji,jj,jk) = zhpj(ji,jj,jk-1)   &
235                  &           + zcoef1 * (  ( rhd(ji,jj+1,jk)+rhd(ji,jj+1,jk-1) )   &
236                  &                       - ( rhd(ji,jj,  jk)+rhd(ji,jj  ,jk-1) )  ) / e2v(ji,jj)
237               ! add to the general momentum trend
238               ua(ji,jj,jk) = ua(ji,jj,jk) + zhpi(ji,jj,jk)
239               va(ji,jj,jk) = va(ji,jj,jk) + zhpj(ji,jj,jk)
240            END DO
241         END DO
242      END DO
243      !
244   END SUBROUTINE hpg_zco
245
246
247   SUBROUTINE hpg_zps( kt )
248      !!---------------------------------------------------------------------
249      !!                 ***  ROUTINE hpg_zps  ***
250      !!                   
251      !! ** Method  :   z-coordinate plus partial steps case.  blahblah...
252      !!
253      !! ** Action  : - Update (ua,va) with the now hydrastatic pressure trend
254      !!----------------------------------------------------------------------
255      USE oce, ONLY:   zhpi => ta , zhpj => sa   ! (ta,sa) used as 3D workspace
256      !!
257      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
258      !!
259      INTEGER  ::   ji, jj, jk                       ! dummy loop indices
260      INTEGER  ::   iku, ikv                         ! temporary integers
261      REAL(wp) ::   zcoef0, zcoef1, zcoef2, zcoef3   ! temporary scalars
262      !!----------------------------------------------------------------------
263
264      IF( kt == nit000 ) THEN
265         IF(lwp) WRITE(numout,*)
266         IF(lwp) WRITE(numout,*) 'dyn:hpg_zps : hydrostatic pressure gradient trend'
267         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   z-coordinate with partial steps - vector optimization'
268      ENDIF
269
270      ! Local constant initialization
271      zcoef0 = - grav * 0.5_wp
272
273      !  Surface value (also valid in partial step case)
274      DO jj = 2, jpjm1
275         DO ji = fs_2, fs_jpim1   ! vector opt.
276            zcoef1 = zcoef0 * fse3w(ji,jj,1)
277            ! hydrostatic pressure gradient
278            zhpi(ji,jj,1) = zcoef1 * ( rhd(ji+1,jj  ,1) - rhd(ji,jj,1) ) / e1u(ji,jj)
279            zhpj(ji,jj,1) = zcoef1 * ( rhd(ji  ,jj+1,1) - rhd(ji,jj,1) ) / e2v(ji,jj)
280            ! add to the general momentum trend
281            ua(ji,jj,1) = ua(ji,jj,1) + zhpi(ji,jj,1)
282            va(ji,jj,1) = va(ji,jj,1) + zhpj(ji,jj,1)
283         END DO
284      END DO
285
286      ! interior value (2=<jk=<jpkm1)
287      DO jk = 2, jpkm1
288         DO jj = 2, jpjm1
289            DO ji = fs_2, fs_jpim1   ! vector opt.
290               zcoef1 = zcoef0 * fse3w(ji,jj,jk)
291               ! hydrostatic pressure gradient
292               zhpi(ji,jj,jk) = zhpi(ji,jj,jk-1)   &
293                  &           + zcoef1 * (  ( rhd(ji+1,jj,jk) + rhd(ji+1,jj,jk-1) )   &
294                  &                       - ( rhd(ji  ,jj,jk) + rhd(ji  ,jj,jk-1) )  ) / e1u(ji,jj)
295
296               zhpj(ji,jj,jk) = zhpj(ji,jj,jk-1)   &
297                  &           + zcoef1 * (  ( rhd(ji,jj+1,jk) + rhd(ji,jj+1,jk-1) )   &
298                  &                       - ( rhd(ji,jj,  jk) + rhd(ji,jj  ,jk-1) )  ) / e2v(ji,jj)
299               ! add to the general momentum trend
300               ua(ji,jj,jk) = ua(ji,jj,jk) + zhpi(ji,jj,jk)
301               va(ji,jj,jk) = va(ji,jj,jk) + zhpj(ji,jj,jk)
302            END DO
303         END DO
304      END DO
305
306      ! partial steps correction at the last level  (use gru & grv computed in zpshde.F90)
307# if defined key_vectopt_loop
308         jj = 1
309         DO ji = jpi+2, jpij-jpi-1   ! vector opt. (forced unrolling)
310# else
311      DO jj = 2, jpjm1
312         DO ji = 2, jpim1
313# endif
314            iku = mbku(ji,jj)
315            ikv = mbkv(ji,jj)
316            zcoef2 = zcoef0 * MIN( fse3w(ji,jj,iku), fse3w(ji+1,jj  ,iku) )
317            zcoef3 = zcoef0 * MIN( fse3w(ji,jj,ikv), fse3w(ji  ,jj+1,ikv) )
318            IF( iku > 1 ) THEN            ! on i-direction (level 2 or more)
319               ua  (ji,jj,iku) = ua(ji,jj,iku) - zhpi(ji,jj,iku)         ! subtract old value
320               zhpi(ji,jj,iku) = zhpi(ji,jj,iku-1)                   &   ! compute the new one
321                  &            + zcoef2 * ( rhd(ji+1,jj,iku-1) - rhd(ji,jj,iku-1) + gru(ji,jj) ) / e1u(ji,jj)
322               ua  (ji,jj,iku) = ua(ji,jj,iku) + zhpi(ji,jj,iku)         ! add the new one to the general momentum trend
323            ENDIF
324            IF( ikv > 1 ) THEN            ! on j-direction (level 2 or more)
325               va  (ji,jj,ikv) = va(ji,jj,ikv) - zhpj(ji,jj,ikv)         ! subtract old value
326               zhpj(ji,jj,ikv) = zhpj(ji,jj,ikv-1)                   &   ! compute the new one
327                  &            + zcoef3 * ( rhd(ji,jj+1,ikv-1) - rhd(ji,jj,ikv-1) + grv(ji,jj) ) / e2v(ji,jj)
328               va  (ji,jj,ikv) = va(ji,jj,ikv) + zhpj(ji,jj,ikv)         ! add the new one to the general momentum trend
329            ENDIF
330# if ! defined key_vectopt_loop
331         END DO
332# endif
333      END DO
334      !
335   END SUBROUTINE hpg_zps
336
337
338   SUBROUTINE hpg_sco( kt )
339      !!---------------------------------------------------------------------
340      !!                  ***  ROUTINE hpg_sco  ***
341      !!
342      !! ** Method  :   s-coordinate case. Jacobian scheme.
343      !!      The now hydrostatic pressure gradient at a given level, jk,
344      !!      is computed by taking the vertical integral of the in-situ
345      !!      density gradient along the model level from the suface to that
346      !!      level. s-coordinates (ln_sco): a corrective term is added
347      !!      to the horizontal pressure gradient :
348      !!         zhpi = grav .....  + 1/e1u mi(rhd) di[ grav dep3w ]
349      !!         zhpj = grav .....  + 1/e2v mj(rhd) dj[ grav dep3w ]
350      !!      add it to the general momentum trend (ua,va).
351      !!         ua = ua - 1/e1u * zhpi
352      !!         va = va - 1/e2v * zhpj
353      !!
354      !! ** Action : - Update (ua,va) with the now hydrastatic pressure trend
355      !!----------------------------------------------------------------------
356      USE oce, ONLY:   zhpi => ta , zhpj => sa   ! (ta,sa) used as 3D workspace
357      !!
358      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
359      !!
360      INTEGER  ::   ji, jj, jk                 ! dummy loop indices
361      REAL(wp) ::   zcoef0, zuap, zvap, znad   ! temporary scalars
362      !!----------------------------------------------------------------------
363
364      IF( kt == nit000 ) THEN
365         IF(lwp) WRITE(numout,*)
366         IF(lwp) WRITE(numout,*) 'dyn:hpg_sco : hydrostatic pressure gradient trend'
367         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   s-coordinate case, OPA original scheme used'
368      ENDIF
369
370      ! Local constant initialization
371      zcoef0 = - grav * 0.5_wp
372      ! To use density and not density anomaly
373      IF ( lk_vvl ) THEN   ;     znad = 1._wp          ! Variable volume
374      ELSE                 ;     znad = 0._wp         ! Fixed volume
375      ENDIF
376
377      ! Surface value
378      DO jj = 2, jpjm1
379         DO ji = fs_2, fs_jpim1   ! vector opt.   
380            ! hydrostatic pressure gradient along s-surfaces
381            zhpi(ji,jj,1) = zcoef0 / e1u(ji,jj) * ( fse3w(ji+1,jj  ,1) * ( znad + rhd(ji+1,jj  ,1) )   &
382               &                                  - fse3w(ji  ,jj  ,1) * ( znad + rhd(ji  ,jj  ,1) ) )
383            zhpj(ji,jj,1) = zcoef0 / e2v(ji,jj) * ( fse3w(ji  ,jj+1,1) * ( znad + rhd(ji  ,jj+1,1) )   &
384               &                                  - fse3w(ji  ,jj  ,1) * ( znad + rhd(ji  ,jj  ,1) ) )
385            ! s-coordinate pressure gradient correction
386            zuap = -zcoef0 * ( rhd   (ji+1,jj,1) + rhd   (ji,jj,1) + 2._wp * znad )   &
387               &           * ( fsde3w(ji+1,jj,1) - fsde3w(ji,jj,1) ) / e1u(ji,jj)
388            zvap = -zcoef0 * ( rhd   (ji,jj+1,1) + rhd   (ji,jj,1) + 2._wp * znad )   &
389               &           * ( fsde3w(ji,jj+1,1) - fsde3w(ji,jj,1) ) / e2v(ji,jj)
390            ! add to the general momentum trend
391            ua(ji,jj,1) = ua(ji,jj,1) + zhpi(ji,jj,1) + zuap
392            va(ji,jj,1) = va(ji,jj,1) + zhpj(ji,jj,1) + zvap
393         END DO 
394      END DO   
395           
396      ! interior value (2=<jk=<jpkm1)
397      DO jk = 2, jpkm1                                 
398         DO jj = 2, jpjm1     
399            DO ji = fs_2, fs_jpim1   ! vector opt.     
400               ! hydrostatic pressure gradient along s-surfaces
401               zhpi(ji,jj,jk) = zhpi(ji,jj,jk-1) + zcoef0 / e1u(ji,jj)   & 
402                  &           * (  fse3w(ji+1,jj,jk) * ( rhd(ji+1,jj,jk) + rhd(ji+1,jj,jk-1) + 2*znad )   & 
403                  &              - fse3w(ji  ,jj,jk) * ( rhd(ji  ,jj,jk) + rhd(ji  ,jj,jk-1) + 2*znad )  )
404               zhpj(ji,jj,jk) = zhpj(ji,jj,jk-1) + zcoef0 / e2v(ji,jj)   &
405                  &           * (  fse3w(ji,jj+1,jk) * ( rhd(ji,jj+1,jk) + rhd(ji,jj+1,jk-1) + 2*znad )   &
406                  &              - fse3w(ji,jj  ,jk) * ( rhd(ji,jj,  jk) + rhd(ji,jj  ,jk-1) + 2*znad )  )
407               ! s-coordinate pressure gradient correction
408               zuap = -zcoef0 * ( rhd   (ji+1,jj  ,jk) + rhd   (ji,jj,jk) + 2._wp * znad )   &
409                  &           * ( fsde3w(ji+1,jj  ,jk) - fsde3w(ji,jj,jk) ) / e1u(ji,jj)
410               zvap = -zcoef0 * ( rhd   (ji  ,jj+1,jk) + rhd   (ji,jj,jk) + 2._wp * znad )   &
411                  &           * ( fsde3w(ji  ,jj+1,jk) - fsde3w(ji,jj,jk) ) / e2v(ji,jj)
412               ! add to the general momentum trend
413               ua(ji,jj,jk) = ua(ji,jj,jk) + zhpi(ji,jj,jk) + zuap
414               va(ji,jj,jk) = va(ji,jj,jk) + zhpj(ji,jj,jk) + zvap
415            END DO
416         END DO
417      END DO
418      !
419   END SUBROUTINE hpg_sco
420
421
422   SUBROUTINE hpg_hel( kt )
423      !!---------------------------------------------------------------------
424      !!                  ***  ROUTINE hpg_hel  ***
425      !!
426      !! ** Method  :   s-coordinate case.
427      !!      The now hydrostatic pressure gradient at a given level
428      !!      jk is computed by taking the vertical integral of the in-situ
429      !!      density gradient along the model level from the suface to that
430      !!      level. s-coordinates (ln_sco): a corrective term is added
431      !!      to the horizontal pressure gradient :
432      !!         zhpi = grav .....  + 1/e1u mi(rhd) di[ grav dep3w ]
433      !!         zhpj = grav .....  + 1/e2v mj(rhd) dj[ grav dep3w ]
434      !!      add it to the general momentum trend (ua,va).
435      !!         ua = ua - 1/e1u * zhpi
436      !!         va = va - 1/e2v * zhpj
437      !!
438      !! ** Action : - Update (ua,va) with the now hydrastatic pressure trend
439      !!             - Save the trend (l_trddyn=T)
440      !!----------------------------------------------------------------------
441      USE oce, ONLY:   zhpi => ta , zhpj => sa   ! (ta,sa) used as 3D workspace
442      !!
443      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
444      !!
445      INTEGER  ::   ji, jj, jk           ! dummy loop indices
446      REAL(wp) ::   zcoef0, zuap, zvap   ! temporary scalars
447      !!----------------------------------------------------------------------
448
449      IF( kt == nit000 ) THEN
450         IF(lwp) WRITE(numout,*)
451         IF(lwp) WRITE(numout,*) 'dyn:hpg_hel : hydrostatic pressure gradient trend'
452         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   s-coordinate case, helsinki modified scheme'
453      ENDIF
454
455      ! Local constant initialization
456      zcoef0 = - grav * 0.5_wp
457 
458      ! Surface value
459      DO jj = 2, jpjm1
460         DO ji = fs_2, fs_jpim1   ! vector opt.
461            ! hydrostatic pressure gradient along s-surfaces
462            zhpi(ji,jj,1) = zcoef0 / e1u(ji,jj) * ( fse3t(ji+1,jj  ,1) * rhd(ji+1,jj  ,1)  &
463               &                                  - fse3t(ji  ,jj  ,1) * rhd(ji  ,jj  ,1) )
464            zhpj(ji,jj,1) = zcoef0 / e2v(ji,jj) * ( fse3t(ji  ,jj+1,1) * rhd(ji  ,jj+1,1)  &
465               &                                  - fse3t(ji  ,jj  ,1) * rhd(ji  ,jj  ,1) )
466            ! s-coordinate pressure gradient correction
467            zuap = -zcoef0 * ( rhd   (ji+1,jj,1) + rhd   (ji,jj,1) )   &
468               &           * ( fsdept(ji+1,jj,1) - fsdept(ji,jj,1) ) / e1u(ji,jj)
469            zvap = -zcoef0 * ( rhd   (ji,jj+1,1) + rhd   (ji,jj,1) )   &
470               &           * ( fsdept(ji,jj+1,1) - fsdept(ji,jj,1) ) / e2v(ji,jj)
471            ! add to the general momentum trend
472            ua(ji,jj,1) = ua(ji,jj,1) + zhpi(ji,jj,1) + zuap
473            va(ji,jj,1) = va(ji,jj,1) + zhpj(ji,jj,1) + zvap
474         END DO
475      END DO
476      !
477      ! interior value (2=<jk=<jpkm1)
478      DO jk = 2, jpkm1
479         DO jj = 2, jpjm1
480            DO ji = fs_2, fs_jpim1   ! vector opt.
481               ! hydrostatic pressure gradient along s-surfaces
482               zhpi(ji,jj,jk) = zhpi(ji,jj,jk-1) &
483                  &           +  zcoef0 / e1u(ji,jj) * ( fse3t(ji+1,jj,jk  ) * rhd(ji+1,jj,jk)     &
484                  &                                     -fse3t(ji  ,jj,jk  ) * rhd(ji  ,jj,jk)   ) &
485                  &           +  zcoef0 / e1u(ji,jj) * ( fse3t(ji+1,jj,jk-1) * rhd(ji+1,jj,jk-1)   &
486                  &                                     -fse3t(ji  ,jj,jk-1) * rhd(ji  ,jj,jk-1) )
487               zhpj(ji,jj,jk) = zhpj(ji,jj,jk-1) &
488                  &           +  zcoef0 / e2v(ji,jj) * ( fse3t(ji,jj+1,jk  ) * rhd(ji,jj+1,jk)     &
489                  &                                     -fse3t(ji,jj  ,jk  ) * rhd(ji,jj,  jk)   ) &
490                  &           +  zcoef0 / e2v(ji,jj) * ( fse3t(ji,jj+1,jk-1) * rhd(ji,jj+1,jk-1)   &
491                  &                                     -fse3t(ji,jj  ,jk-1) * rhd(ji,jj,  jk-1) )
492               ! s-coordinate pressure gradient correction
493               zuap = - zcoef0 * ( rhd   (ji+1,jj,jk) + rhd   (ji,jj,jk) )   &
494                  &            * ( fsdept(ji+1,jj,jk) - fsdept(ji,jj,jk) ) / e1u(ji,jj)
495               zvap = - zcoef0 * ( rhd   (ji,jj+1,jk) + rhd   (ji,jj,jk) )   &
496                  &            * ( fsdept(ji,jj+1,jk) - fsdept(ji,jj,jk) ) / e2v(ji,jj)
497               ! add to the general momentum trend
498               ua(ji,jj,jk) = ua(ji,jj,jk) + zhpi(ji,jj,jk) + zuap
499               va(ji,jj,jk) = va(ji,jj,jk) + zhpj(ji,jj,jk) + zvap
500            END DO
501         END DO
502      END DO
503      !
504   END SUBROUTINE hpg_hel
505
506
507   SUBROUTINE hpg_wdj( kt )
508      !!---------------------------------------------------------------------
509      !!                  ***  ROUTINE hpg_wdj  ***
510      !!
511      !! ** Method  :   Weighted Density Jacobian (wdj) scheme (song 1998)
512      !!      The weighting coefficients from the namelist parameter rn_gamma
513      !!      (alpha=0.5-rn_gamma ; beta=1-alpha=0.5+rn_gamma
514      !!
515      !! Reference : Song, Mon. Wea. Rev., 126, 3213-3230, 1998.
516      !!----------------------------------------------------------------------
517      USE oce, ONLY:   zhpi => ta , zhpj => sa   ! (ta,sa) used as 3D workspace
518      !!
519      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
520      !!
521      INTEGER  ::   ji, jj, jk           ! dummy loop indices
522      REAL(wp) ::   zcoef0, zuap, zvap   ! temporary scalars
523      REAL(wp) ::   zalph , zbeta        !    "         "
524      !!----------------------------------------------------------------------
525
526      IF( kt == nit000 ) THEN
527         IF(lwp) WRITE(numout,*)
528         IF(lwp) WRITE(numout,*) 'dyn:hpg_wdj : hydrostatic pressure gradient trend'
529         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   Weighted Density Jacobian'
530      ENDIF
531
532      ! Local constant initialization
533      zcoef0 = - grav * 0.5_wp
534      zalph  = 0.5_wp - rn_gamma    ! weighting coefficients (alpha=0.5-rn_gamma
535      zbeta  = 0.5_wp + rn_gamma    !                        (beta =1-alpha=0.5+rn_gamma
536
537      ! Surface value (no ponderation)
538      DO jj = 2, jpjm1
539         DO ji = fs_2, fs_jpim1   ! vector opt.
540            ! hydrostatic pressure gradient along s-surfaces
541            zhpi(ji,jj,1) = zcoef0 / e1u(ji,jj) * (  fse3w(ji+1,jj  ,1) * rhd(ji+1,jj  ,1)   &
542               &                                   - fse3w(ji  ,jj  ,1) * rhd(ji  ,jj  ,1)  )
543            zhpj(ji,jj,1) = zcoef0 / e2v(ji,jj) * (  fse3w(ji  ,jj+1,1) * rhd(ji  ,jj+1,1)   &
544               &                                   - fse3w(ji  ,jj  ,1) * rhd(ji,  jj  ,1)  )
545            ! s-coordinate pressure gradient correction
546            zuap = -zcoef0 * ( rhd   (ji+1,jj,1) + rhd   (ji,jj,1) )   &
547               &           * ( fsde3w(ji+1,jj,1) - fsde3w(ji,jj,1) ) / e1u(ji,jj)
548            zvap = -zcoef0 * ( rhd   (ji,jj+1,1) + rhd   (ji,jj,1) )   &
549               &           * ( fsde3w(ji,jj+1,1) - fsde3w(ji,jj,1) ) / e2v(ji,jj)
550            ! add to the general momentum trend
551            ua(ji,jj,1) = ua(ji,jj,1) + zhpi(ji,jj,1) + zuap
552            va(ji,jj,1) = va(ji,jj,1) + zhpj(ji,jj,1) + zvap
553         END DO
554      END DO
555
556      ! Interior value (2=<jk=<jpkm1) (weighted with zalph & zbeta)
557      DO jk = 2, jpkm1
558         DO jj = 2, jpjm1
559            DO ji = fs_2, fs_jpim1   ! vector opt.
560               zhpi(ji,jj,jk) = zhpi(ji,jj,jk-1) + zcoef0 / e1u(ji,jj)                            &
561                  &           * (   (            fsde3w(ji+1,jj,jk  ) + fsde3w(ji,jj,jk  )        &
562                  &                            - fsde3w(ji+1,jj,jk-1) - fsde3w(ji,jj,jk-1)    )   &
563                  &               * (  zalph * ( rhd   (ji+1,jj,jk-1) - rhd   (ji,jj,jk-1) )      &
564                  &                  + zbeta * ( rhd   (ji+1,jj,jk  ) - rhd   (ji,jj,jk  ) )  )   &
565                  &             -   (            rhd   (ji+1,jj,jk  ) + rhd   (ji,jj,jk  )        &
566                  &                           - rhd   (ji+1,jj,jk-1) - rhd   (ji,jj,jk-1)     )   &
567                  &               * (  zalph * ( fsde3w(ji+1,jj,jk-1) - fsde3w(ji,jj,jk-1) )      &
568                  &                  + zbeta * ( fsde3w(ji+1,jj,jk  ) - fsde3w(ji,jj,jk  ) )  )  )
569               zhpj(ji,jj,jk) = zhpj(ji,jj,jk-1) + zcoef0 / e2v(ji,jj)                            &
570                  &           * (   (           fsde3w(ji,jj+1,jk  ) + fsde3w(ji,jj,jk  )         &
571                  &                           - fsde3w(ji,jj+1,jk-1) - fsde3w(ji,jj,jk-1)     )   &
572                  &               * (  zalph * ( rhd   (ji,jj+1,jk-1) - rhd   (ji,jj,jk-1) )      &
573                  &                  + zbeta * ( rhd   (ji,jj+1,jk  ) - rhd   (ji,jj,jk  ) )  )   &
574                  &             -   (            rhd   (ji,jj+1,jk  ) + rhd   (ji,jj,jk  )        &
575                  &                            - rhd   (ji,jj+1,jk-1) - rhd   (ji,jj,jk-1)    )   &
576                  &               * (  zalph * ( fsde3w(ji,jj+1,jk-1) - fsde3w(ji,jj,jk-1) )      &
577                  &                  + zbeta * ( fsde3w(ji,jj+1,jk  ) - fsde3w(ji,jj,jk  ) )  )  )
578               ! add to the general momentum trend
579               ua(ji,jj,jk) = ua(ji,jj,jk) + zhpi(ji,jj,jk)
580               va(ji,jj,jk) = va(ji,jj,jk) + zhpj(ji,jj,jk)
581            END DO
582         END DO
583      END DO
584      !
585   END SUBROUTINE hpg_wdj
586
587
588   SUBROUTINE hpg_djc( kt )
589      !!---------------------------------------------------------------------
590      !!                  ***  ROUTINE hpg_djc  ***
591      !!
592      !! ** Method  :   Density Jacobian with Cubic polynomial scheme
593      !!
594      !! Reference: Shchepetkin and McWilliams, J. Geophys. Res., 108(C3), 3090, 2003
595      !!----------------------------------------------------------------------
596      USE wrk_nemo, ONLY:   wrk_in_use, wrk_not_released
597      USE oce     , ONLY:   zhpi  => ta        , zhpj => sa       ! (ta,sa) used as 3D workspace
598      USE wrk_nemo, ONLY:   drhox => wrk_3d_1  , dzx  => wrk_3d_2
599      USE wrk_nemo, ONLY:   drhou => wrk_3d_3  , dzu  => wrk_3d_4 , rho_i => wrk_3d_5
600      USE wrk_nemo, ONLY:   drhoy => wrk_3d_6  , dzy  => wrk_3d_7
601      USE wrk_nemo, ONLY:   drhov => wrk_3d_8  , dzv  => wrk_3d_9 , rho_j => wrk_3d_10
602      USE wrk_nemo, ONLY:   drhoz => wrk_3d_11 , dzz  => wrk_3d_12 
603      USE wrk_nemo, ONLY:   drhow => wrk_3d_13 , dzw  => wrk_3d_14
604      USE wrk_nemo, ONLY:   rho_k => wrk_3d_15
605      !!
606      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
607      !!
608      INTEGER  ::   ji, jj, jk          ! dummy loop indices
609      REAL(wp) ::   zcoef0, zep, cffw   ! temporary scalars
610      REAL(wp) ::   z1_10, cffu, cffx   !    "         "
611      REAL(wp) ::   z1_12, cffv, cffy   !    "         "
612      !!----------------------------------------------------------------------
613
614      IF( wrk_in_use(3, 1,2,3,4,5,6,7,8,9,10,11,12,13,14,15) ) THEN
615         CALL ctl_stop('dyn:hpg_djc: requested workspace arrays unavailable')   ;   RETURN
616      ENDIF
617
618      IF( kt == nit000 ) THEN
619         IF(lwp) WRITE(numout,*)
620         IF(lwp) WRITE(numout,*) 'dyn:hpg_djc : hydrostatic pressure gradient trend'
621         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   s-coordinate case, density Jacobian with cubic polynomial scheme'
622      ENDIF
623
624      ! Local constant initialization
625      zcoef0 = - grav * 0.5_wp
626      z1_10  = 1._wp / 10._wp
627      z1_12  = 1._wp / 12._wp
628
629      !----------------------------------------------------------------------------------------
630      !  compute and store in provisional arrays elementary vertical and horizontal differences
631      !----------------------------------------------------------------------------------------
632
633!!bug gm   Not a true bug, but... dzz=e3w  for dzx, dzy verify what it is really
634
635      DO jk = 2, jpkm1
636         DO jj = 2, jpjm1
637            DO ji = fs_2, fs_jpim1   ! vector opt.
638               drhoz(ji,jj,jk) = rhd   (ji  ,jj  ,jk) - rhd   (ji,jj,jk-1)
639               dzz  (ji,jj,jk) = fsde3w(ji  ,jj  ,jk) - fsde3w(ji,jj,jk-1)
640               drhox(ji,jj,jk) = rhd   (ji+1,jj  ,jk) - rhd   (ji,jj,jk  )
641               dzx  (ji,jj,jk) = fsde3w(ji+1,jj  ,jk) - fsde3w(ji,jj,jk  )
642               drhoy(ji,jj,jk) = rhd   (ji  ,jj+1,jk) - rhd   (ji,jj,jk  )
643               dzy  (ji,jj,jk) = fsde3w(ji  ,jj+1,jk) - fsde3w(ji,jj,jk  )
644            END DO
645         END DO
646      END DO
647
648      !-------------------------------------------------------------------------
649      ! compute harmonic averages using eq. 5.18
650      !-------------------------------------------------------------------------
651      zep = 1.e-15
652
653!!bug  gm  drhoz not defined at level 1 and used (jk-1 with jk=2)
654!!bug  gm  idem for drhox, drhoy et ji=jpi and jj=jpj
655
656      DO jk = 2, jpkm1
657         DO jj = 2, jpjm1
658            DO ji = fs_2, fs_jpim1   ! vector opt.
659               cffw = 2._wp * drhoz(ji  ,jj  ,jk) * drhoz(ji,jj,jk-1)
660
661               cffu = 2._wp * drhox(ji+1,jj  ,jk) * drhox(ji,jj,jk  )
662               cffx = 2._wp * dzx  (ji+1,jj  ,jk) * dzx  (ji,jj,jk  )
663 
664               cffv = 2._wp * drhoy(ji  ,jj+1,jk) * drhoy(ji,jj,jk  )
665               cffy = 2._wp * dzy  (ji  ,jj+1,jk) * dzy  (ji,jj,jk  )
666
667               IF( cffw > zep) THEN
668                  drhow(ji,jj,jk) = 2._wp *   drhoz(ji,jj,jk) * drhoz(ji,jj,jk-1)   &
669                     &                    / ( drhoz(ji,jj,jk) + drhoz(ji,jj,jk-1) )
670               ELSE
671                  drhow(ji,jj,jk) = 0._wp
672               ENDIF
673
674               dzw(ji,jj,jk) = 2._wp *   dzz(ji,jj,jk) * dzz(ji,jj,jk-1)   &
675                  &                  / ( dzz(ji,jj,jk) + dzz(ji,jj,jk-1) )
676
677               IF( cffu > zep ) THEN
678                  drhou(ji,jj,jk) = 2._wp *   drhox(ji+1,jj,jk) * drhox(ji,jj,jk)   &
679                     &                    / ( drhox(ji+1,jj,jk) + drhox(ji,jj,jk) )
680               ELSE
681                  drhou(ji,jj,jk ) = 0._wp
682               ENDIF
683
684               IF( cffx > zep ) THEN
685                  dzu(ji,jj,jk) = 2._wp *   dzx(ji+1,jj,jk) * dzx(ji,jj,jk)   &
686                     &                  / ( dzx(ji+1,jj,jk) + dzx(ji,jj,jk) )
687               ELSE
688                  dzu(ji,jj,jk) = 0._wp
689               ENDIF
690
691               IF( cffv > zep ) THEN
692                  drhov(ji,jj,jk) = 2._wp *   drhoy(ji,jj+1,jk) * drhoy(ji,jj,jk)   &
693                     &                    / ( drhoy(ji,jj+1,jk) + drhoy(ji,jj,jk) )
694               ELSE
695                  drhov(ji,jj,jk) = 0._wp
696               ENDIF
697
698               IF( cffy > zep ) THEN
699                  dzv(ji,jj,jk) = 2._wp *   dzy(ji,jj+1,jk) * dzy(ji,jj,jk)   &
700                     &                  / ( dzy(ji,jj+1,jk) + dzy(ji,jj,jk) )
701               ELSE
702                  dzv(ji,jj,jk) = 0._wp
703               ENDIF
704
705            END DO
706         END DO
707      END DO
708
709      !----------------------------------------------------------------------------------
710      ! apply boundary conditions at top and bottom using 5.36-5.37
711      !----------------------------------------------------------------------------------
712      drhow(:,:, 1 ) = 1.5_wp * ( drhoz(:,:, 2 ) - drhoz(:,:,  1  ) ) - 0.5_wp * drhow(:,:,  2  )
713      drhou(:,:, 1 ) = 1.5_wp * ( drhox(:,:, 2 ) - drhox(:,:,  1  ) ) - 0.5_wp * drhou(:,:,  2  )
714      drhov(:,:, 1 ) = 1.5_wp * ( drhoy(:,:, 2 ) - drhoy(:,:,  1  ) ) - 0.5_wp * drhov(:,:,  2  )
715
716      drhow(:,:,jpk) = 1.5_wp * ( drhoz(:,:,jpk) - drhoz(:,:,jpkm1) ) - 0.5_wp * drhow(:,:,jpkm1)
717      drhou(:,:,jpk) = 1.5_wp * ( drhox(:,:,jpk) - drhox(:,:,jpkm1) ) - 0.5_wp * drhou(:,:,jpkm1)
718      drhov(:,:,jpk) = 1.5_wp * ( drhoy(:,:,jpk) - drhoy(:,:,jpkm1) ) - 0.5_wp * drhov(:,:,jpkm1)
719
720
721      !--------------------------------------------------------------
722      ! Upper half of top-most grid box, compute and store
723      !-------------------------------------------------------------
724
725!!bug gm   :  e3w-de3w = 0.5*e3w  ....  and de3w(2)-de3w(1)=e3w(2) ....   to be verified
726!          true if de3w is really defined as the sum of the e3w scale factors as, it seems to me, it should be
727
728      DO jj = 2, jpjm1
729         DO ji = fs_2, fs_jpim1   ! vector opt.
730            rho_k(ji,jj,1) = -grav * ( fse3w(ji,jj,1) - fsde3w(ji,jj,1) )               &
731               &                   * (  rhd(ji,jj,1)                                    &
732               &                     + 0.5_wp * ( rhd(ji,jj,2) - rhd(ji,jj,1) )         &
733               &                              * ( fse3w (ji,jj,1) - fsde3w(ji,jj,1) )   &
734               &                              / ( fsde3w(ji,jj,2) - fsde3w(ji,jj,1) )  ) 
735         END DO
736      END DO
737
738!!bug gm    : here also, simplification is possible
739!!bug gm    : optimisation: 1/10 and 1/12 the division should be done before the loop
740
741      DO jk = 2, jpkm1
742         DO jj = 2, jpjm1
743            DO ji = fs_2, fs_jpim1   ! vector opt.
744
745               rho_k(ji,jj,jk) = zcoef0 * ( rhd   (ji,jj,jk) + rhd   (ji,jj,jk-1) )                                   &
746                  &                     * ( fsde3w(ji,jj,jk) - fsde3w(ji,jj,jk-1) )                                   &
747                  &            - grav * z1_10 * (                                                                     &
748                  &     ( drhow (ji,jj,jk) - drhow (ji,jj,jk-1) )                                                     &
749                  &   * ( fsde3w(ji,jj,jk) - fsde3w(ji,jj,jk-1) - z1_12 * ( dzw  (ji,jj,jk) + dzw  (ji,jj,jk-1) ) )   &
750                  &   - ( dzw   (ji,jj,jk) - dzw   (ji,jj,jk-1) )                                                     &
751                  &   * ( rhd   (ji,jj,jk) - rhd   (ji,jj,jk-1) - z1_12 * ( drhow(ji,jj,jk) + drhow(ji,jj,jk-1) ) )   &
752                  &                             )
753
754               rho_i(ji,jj,jk) = zcoef0 * ( rhd   (ji+1,jj,jk) + rhd   (ji,jj,jk) )                                   &
755                  &                     * ( fsde3w(ji+1,jj,jk) - fsde3w(ji,jj,jk) )                                   &
756                  &            - grav* z1_10 * (                                                                      &
757                  &     ( drhou (ji+1,jj,jk) - drhou (ji,jj,jk) )                                                     &
758                  &   * ( fsde3w(ji+1,jj,jk) - fsde3w(ji,jj,jk) - z1_12 * ( dzu  (ji+1,jj,jk) + dzu  (ji,jj,jk) ) )   &
759                  &   - ( dzu   (ji+1,jj,jk) - dzu   (ji,jj,jk) )                                                     &
760                  &   * ( rhd   (ji+1,jj,jk) - rhd   (ji,jj,jk) - z1_12 * ( drhou(ji+1,jj,jk) + drhou(ji,jj,jk) ) )   &
761                  &                            )
762
763               rho_j(ji,jj,jk) = zcoef0 * ( rhd   (ji,jj+1,jk) + rhd   (ji,jj,jk) )                                   &
764                  &                     * ( fsde3w(ji,jj+1,jk) - fsde3w(ji,jj,jk) )                                   &
765                  &            - grav* z1_10 * (                                                                      &
766                  &     ( drhov (ji,jj+1,jk) - drhov (ji,jj,jk) )                                                     &
767                  &   * ( fsde3w(ji,jj+1,jk) - fsde3w(ji,jj,jk) - z1_12 * ( dzv  (ji,jj+1,jk) + dzv  (ji,jj,jk) ) )   &
768                  &   - ( dzv   (ji,jj+1,jk) - dzv   (ji,jj,jk) )                                                     &
769                  &   * ( rhd   (ji,jj+1,jk) - rhd   (ji,jj,jk) - z1_12 * ( drhov(ji,jj+1,jk) + drhov(ji,jj,jk) ) )   &
770                  &                            )
771
772            END DO
773         END DO
774      END DO
775      CALL lbc_lnk(rho_k,'W',1.)
776      CALL lbc_lnk(rho_i,'U',1.)
777      CALL lbc_lnk(rho_j,'V',1.)
778
779
780      ! ---------------
781      !  Surface value
782      ! ---------------
783      DO jj = 2, jpjm1
784         DO ji = fs_2, fs_jpim1   ! vector opt.
785            zhpi(ji,jj,1) = ( rho_k(ji+1,jj  ,1) - rho_k(ji,jj,1) - rho_i(ji,jj,1) ) / e1u(ji,jj)
786            zhpj(ji,jj,1) = ( rho_k(ji  ,jj+1,1) - rho_k(ji,jj,1) - rho_j(ji,jj,1) ) / e2v(ji,jj)
787            ! add to the general momentum trend
788            ua(ji,jj,1) = ua(ji,jj,1) + zhpi(ji,jj,1)
789            va(ji,jj,1) = va(ji,jj,1) + zhpj(ji,jj,1)
790         END DO
791      END DO
792
793      ! ----------------
794      !  interior value   (2=<jk=<jpkm1)
795      ! ----------------
796      DO jk = 2, jpkm1
797         DO jj = 2, jpjm1 
798            DO ji = fs_2, fs_jpim1   ! vector opt.
799               ! hydrostatic pressure gradient along s-surfaces
800               zhpi(ji,jj,jk) = zhpi(ji,jj,jk-1)                                &
801                  &           + (  ( rho_k(ji+1,jj,jk) - rho_k(ji,jj,jk  ) )    &
802                  &              - ( rho_i(ji  ,jj,jk) - rho_i(ji,jj,jk-1) )  ) / e1u(ji,jj)
803               zhpj(ji,jj,jk) = zhpj(ji,jj,jk-1)                                &
804                  &           + (  ( rho_k(ji,jj+1,jk) - rho_k(ji,jj,jk  ) )    &
805                  &               -( rho_j(ji,jj  ,jk) - rho_j(ji,jj,jk-1) )  ) / e2v(ji,jj)
806               ! add to the general momentum trend
807               ua(ji,jj,jk) = ua(ji,jj,jk) + zhpi(ji,jj,jk)
808               va(ji,jj,jk) = va(ji,jj,jk) + zhpj(ji,jj,jk)
809            END DO
810         END DO
811      END DO
812      !
813      IF( wrk_not_released(3, 1,2,3,4,5,6,7,8,9,10,11,12,13,14,15) )   &
814         CALL ctl_stop('dyn:hpg_djc: failed to release workspace arrays')
815      !
816   END SUBROUTINE hpg_djc
817
818
819   SUBROUTINE hpg_rot( kt )
820      !!---------------------------------------------------------------------
821      !!                  ***  ROUTINE hpg_rot  ***
822      !!
823      !! ** Method  :   rotated axes scheme (Thiem and Berntsen 2005)
824      !!
825      !! Reference: Thiem & Berntsen, Ocean Modelling, In press, 2005.
826      !!----------------------------------------------------------------------
827      USE wrk_nemo, ONLY:   wrk_in_use, wrk_not_released
828      USE oce     , ONLY:   zhpi    => ta       , zhpj    => sa       ! (ta,sa) used as 3D workspace
829      USE wrk_nemo, ONLY:   zdistr  => wrk_2d_1 , zsina   => wrk_2d_2 , zcosa  => wrk_2d_3
830      USE wrk_nemo, ONLY:   zhpiorg => wrk_3d_1 , zhpirot => wrk_3d_2
831      USE wrk_nemo, ONLY:   zhpitra => wrk_3d_3 , zhpine  => wrk_3d_4
832      USE wrk_nemo, ONLY:   zhpjorg => wrk_3d_5 , zhpjrot => wrk_3d_6
833      USE wrk_nemo, ONLY:   zhpjtra => wrk_3d_7 , zhpjne  => wrk_3d_8
834      !!
835      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
836      !!
837      INTEGER  ::   ji, jj, jk          ! dummy loop indices
838      REAL(wp) ::   zforg, zcoef0, zuap, zmskd1, zmskd1m   ! temporary scalar
839      REAL(wp) ::   zfrot        , zvap, zmskd2, zmskd2m   !    "         "
840      !!----------------------------------------------------------------------
841
842      IF( wrk_in_use(2, 1,2,3)             .OR.   &
843          wrk_in_use(3, 1,2,3,4,5,6,7,8) ) THEN
844         CALL ctl_stop('dyn:hpg_rot: requested workspace arrays unavailable')   ;   RETURN
845      ENDIF
846
847      IF( kt == nit000 ) THEN
848         IF(lwp) WRITE(numout,*)
849         IF(lwp) WRITE(numout,*) 'dyn:hpg_rot : hydrostatic pressure gradient trend'
850         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   s-coordinate case, rotated axes scheme used'
851      ENDIF
852
853      ! -------------------------------
854      !  Local constant initialization
855      ! -------------------------------
856      zcoef0 = - grav * 0.5_wp
857      zforg  = 0.95_wp
858      zfrot  = 1._wp - zforg
859
860      ! inverse of the distance between 2 diagonal T-points (defined at F-point) (here zcoef0/distance)
861      zdistr(:,:) = zcoef0 / SQRT( e1f(:,:)*e1f(:,:) + e2f(:,:)*e1f(:,:) )
862
863      ! sinus and cosinus of diagonal angle at F-point
864      zsina(:,:) = ATAN2( e2f(:,:), e1f(:,:) )
865      zcosa(:,:) = COS( zsina(:,:) )
866      zsina(:,:) = SIN( zsina(:,:) )
867
868      ! ---------------
869      !  Surface value
870      ! ---------------
871      ! compute and add to the general trend the pressure gradients along the axes
872      DO jj = 2, jpjm1
873         DO ji = fs_2, fs_jpim1   ! vector opt.
874            ! hydrostatic pressure gradient along s-surfaces
875            zhpiorg(ji,jj,1) = zcoef0 / e1u(ji,jj) * (  fse3t(ji+1,jj,1) * rhd(ji+1,jj,1)   &
876               &                                      - fse3t(ji  ,jj,1) * rhd(ji  ,jj,1)  )
877            zhpjorg(ji,jj,1) = zcoef0 / e2v(ji,jj) * (  fse3t(ji,jj+1,1) * rhd(ji,jj+1,1)   &
878               &                                      - fse3t(ji,jj  ,1) * rhd(ji,jj  ,1)  )
879            ! s-coordinate pressure gradient correction
880            zuap = -zcoef0 * ( rhd   (ji+1,jj  ,1) + rhd   (ji,jj,1) )   &
881               &           * ( fsdept(ji+1,jj  ,1) - fsdept(ji,jj,1) ) / e1u(ji,jj)
882            zvap = -zcoef0 * ( rhd   (ji  ,jj+1,1) + rhd   (ji,jj,1) )   &
883               &           * ( fsdept(ji  ,jj+1,1) - fsdept(ji,jj,1) ) / e2v(ji,jj)
884            ! add to the general momentum trend
885            ua(ji,jj,1) = ua(ji,jj,1) + zforg * ( zhpiorg(ji,jj,1) + zuap )
886            va(ji,jj,1) = va(ji,jj,1) + zforg * ( zhpjorg(ji,jj,1) + zvap )
887         END DO
888      END DO
889
890      ! compute the pressure gradients in the diagonal directions
891      DO jj = 1, jpjm1
892         DO ji = 1, fs_jpim1   ! vector opt.
893            zmskd1 = tmask(ji+1,jj+1,1) * tmask(ji  ,jj,1)      ! mask in the 1st diagnonal
894            zmskd2 = tmask(ji  ,jj+1,1) * tmask(ji+1,jj,1)      ! mask in the 2nd diagnonal
895            ! hydrostatic pressure gradient along s-surfaces
896            zhpitra(ji,jj,1) = zdistr(ji,jj) * zmskd1 * (  fse3t(ji+1,jj+1,1) * rhd(ji+1,jj+1,1)   &
897               &                                         - fse3t(ji  ,jj  ,1) * rhd(ji  ,jj  ,1)  )
898            zhpjtra(ji,jj,1) = zdistr(ji,jj) * zmskd2 * (  fse3t(ji  ,jj+1,1) * rhd(ji  ,jj+1,1)   &
899               &                                         - fse3t(ji+1,jj  ,1) * rhd(ji+1,jj  ,1)  )
900            ! s-coordinate pressure gradient correction
901            zuap = -zdistr(ji,jj) * zmskd1 * ( rhd   (ji+1,jj+1,1) + rhd   (ji  ,jj,1) )   &
902               &                           * ( fsdept(ji+1,jj+1,1) - fsdept(ji  ,jj,1) )
903            zvap = -zdistr(ji,jj) * zmskd2 * ( rhd   (ji  ,jj+1,1) + rhd   (ji+1,jj,1) )   &
904               &                           * ( fsdept(ji  ,jj+1,1) - fsdept(ji+1,jj,1) )
905            ! back rotation
906            zhpine(ji,jj,1) = zcosa(ji,jj) * ( zhpitra(ji,jj,1) + zuap )   &
907               &            - zsina(ji,jj) * ( zhpjtra(ji,jj,1) + zvap )
908            zhpjne(ji,jj,1) = zsina(ji,jj) * ( zhpitra(ji,jj,1) + zuap )   &
909               &            + zcosa(ji,jj) * ( zhpjtra(ji,jj,1) + zvap )
910         END DO
911      END DO
912
913      ! interpolate and add to the general trend the diagonal gradient
914      DO jj = 2, jpjm1
915         DO ji = fs_2, fs_jpim1   ! vector opt.
916            ! averaging
917            zhpirot(ji,jj,1) = 0.5 * ( zhpine(ji,jj,1) + zhpine(ji  ,jj-1,1) )
918            zhpjrot(ji,jj,1) = 0.5 * ( zhpjne(ji,jj,1) + zhpjne(ji-1,jj  ,1) )
919            ! add to the general momentum trend
920            ua(ji,jj,1) = ua(ji,jj,1) + zfrot * zhpirot(ji,jj,1) 
921            va(ji,jj,1) = va(ji,jj,1) + zfrot * zhpjrot(ji,jj,1) 
922         END DO
923      END DO
924
925      ! -----------------
926      ! 2. interior value (2=<jk=<jpkm1)
927      ! -----------------
928      ! compute and add to the general trend the pressure gradients along the axes
929      DO jk = 2, jpkm1
930         DO jj = 2, jpjm1
931            DO ji = fs_2, fs_jpim1   ! vector opt.
932               ! hydrostatic pressure gradient along s-surfaces
933               zhpiorg(ji,jj,jk) = zhpiorg(ji,jj,jk-1)                                                 &
934                  &              +  zcoef0 / e1u(ji,jj) * (  fse3t(ji+1,jj,jk  ) * rhd(ji+1,jj,jk  )   &
935                  &                                        - fse3t(ji  ,jj,jk  ) * rhd(ji  ,jj,jk  )   &
936                  &                                        + fse3t(ji+1,jj,jk-1) * rhd(ji+1,jj,jk-1)   &
937                  &                                        - fse3t(ji  ,jj,jk-1) * rhd(ji  ,jj,jk-1)  )
938               zhpjorg(ji,jj,jk) = zhpjorg(ji,jj,jk-1)                                                 &
939                  &              +  zcoef0 / e2v(ji,jj) * (  fse3t(ji,jj+1,jk  ) * rhd(ji,jj+1,jk  )   &
940                  &                                        - fse3t(ji,jj  ,jk  ) * rhd(ji,jj,  jk  )   &
941                  &                                        + fse3t(ji,jj+1,jk-1) * rhd(ji,jj+1,jk-1)   &
942                  &                                        - fse3t(ji,jj  ,jk-1) * rhd(ji,jj,  jk-1)  )
943               ! s-coordinate pressure gradient correction
944               zuap = - zcoef0 * ( rhd   (ji+1,jj  ,jk) + rhd   (ji,jj,jk) )   &
945                  &            * ( fsdept(ji+1,jj  ,jk) - fsdept(ji,jj,jk) ) / e1u(ji,jj)
946               zvap = - zcoef0 * ( rhd   (ji  ,jj+1,jk) + rhd   (ji,jj,jk) )   &
947                  &            * ( fsdept(ji  ,jj+1,jk) - fsdept(ji,jj,jk) ) / e2v(ji,jj)
948               ! add to the general momentum trend
949               ua(ji,jj,jk) = ua(ji,jj,jk) + zforg*( zhpiorg(ji,jj,jk) + zuap )
950               va(ji,jj,jk) = va(ji,jj,jk) + zforg*( zhpjorg(ji,jj,jk) + zvap )
951            END DO
952         END DO
953      END DO
954
955      ! compute the pressure gradients in the diagonal directions
956      DO jk = 2, jpkm1
957         DO jj = 1, jpjm1
958            DO ji = 1, fs_jpim1   ! vector opt.
959               zmskd1  = tmask(ji+1,jj+1,jk  ) * tmask(ji  ,jj,jk  )      ! level jk   mask in the 1st diagnonal
960               zmskd1m = tmask(ji+1,jj+1,jk-1) * tmask(ji  ,jj,jk-1)      ! level jk-1    "               "     
961               zmskd2  = tmask(ji  ,jj+1,jk  ) * tmask(ji+1,jj,jk  )      ! level jk   mask in the 2nd diagnonal
962               zmskd2m = tmask(ji  ,jj+1,jk-1) * tmask(ji+1,jj,jk-1)      ! level jk-1    "               "     
963               ! hydrostatic pressure gradient along s-surfaces
964               zhpitra(ji,jj,jk) = zhpitra(ji,jj,jk-1)                                                       &
965                  &              + zdistr(ji,jj) * zmskd1  * ( fse3t(ji+1,jj+1,jk  ) * rhd(ji+1,jj+1,jk)     &
966                  &                                           -fse3t(ji  ,jj  ,jk  ) * rhd(ji  ,jj  ,jk) )   &
967                  &              + zdistr(ji,jj) * zmskd1m * ( fse3t(ji+1,jj+1,jk-1) * rhd(ji+1,jj+1,jk-1)   &
968                  &                                           -fse3t(ji  ,jj  ,jk-1) * rhd(ji  ,jj  ,jk-1) )
969               zhpjtra(ji,jj,jk) = zhpjtra(ji,jj,jk-1)                                                       &
970                  &              + zdistr(ji,jj) * zmskd2  * ( fse3t(ji  ,jj+1,jk  ) * rhd(ji  ,jj+1,jk)     &
971                  &                                           -fse3t(ji+1,jj  ,jk  ) * rhd(ji+1,jj,  jk) )   &
972                  &              + zdistr(ji,jj) * zmskd2m * ( fse3t(ji  ,jj+1,jk-1) * rhd(ji  ,jj+1,jk-1)   &
973                  &                                           -fse3t(ji+1,jj  ,jk-1) * rhd(ji+1,jj,  jk-1) )
974               ! s-coordinate pressure gradient correction
975               zuap = - zdistr(ji,jj) * zmskd1 * ( rhd   (ji+1,jj+1,jk) + rhd   (ji  ,jj,jk) )   &
976                  &                            * ( fsdept(ji+1,jj+1,jk) - fsdept(ji  ,jj,jk) )
977               zvap = - zdistr(ji,jj) * zmskd2 * ( rhd   (ji  ,jj+1,jk) + rhd   (ji+1,jj,jk) )   &
978                  &                            * ( fsdept(ji  ,jj+1,jk) - fsdept(ji+1,jj,jk) )
979               ! back rotation
980               zhpine(ji,jj,jk) = zcosa(ji,jj) * ( zhpitra(ji,jj,jk) + zuap )   &
981                  &             - zsina(ji,jj) * ( zhpjtra(ji,jj,jk) + zvap )
982               zhpjne(ji,jj,jk) = zsina(ji,jj) * ( zhpitra(ji,jj,jk) + zuap )   &
983                  &             + zcosa(ji,jj) * ( zhpjtra(ji,jj,jk) + zvap )
984            END DO
985         END DO
986      END DO
987
988      ! interpolate and add to the general trend
989      DO jk = 2, jpkm1
990         DO jj = 2, jpjm1
991            DO ji = fs_2, fs_jpim1   ! vector opt.
992               ! averaging
993               zhpirot(ji,jj,jk) = 0.5 * ( zhpine(ji,jj,jk) + zhpine(ji  ,jj-1,jk) )
994               zhpjrot(ji,jj,jk) = 0.5 * ( zhpjne(ji,jj,jk) + zhpjne(ji-1,jj  ,jk) )
995               ! add to the general momentum trend
996               ua(ji,jj,jk) = ua(ji,jj,jk) + zfrot * zhpirot(ji,jj,jk) 
997               va(ji,jj,jk) = va(ji,jj,jk) + zfrot * zhpjrot(ji,jj,jk) 
998            END DO
999         END DO
1000      END DO
1001      !
1002      IF( wrk_not_released(2, 1,2,3)           .OR.   &
1003          wrk_not_released(3, 1,2,3,4,5,6,7,8) )   CALL ctl_stop('dyn:hpg_rot: failed to release workspace arrays')
1004      !
1005   END SUBROUTINE hpg_rot
1006
1007   !!======================================================================
1008END MODULE dynhpg
Note: See TracBrowser for help on using the repository browser.