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 branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/DYN – NEMO

source: branches/2011/DEV_r2739_STFC_dCSE/NEMOGCM/NEMO/OPA_SRC/DYN/dynhpg.F90 @ 4457

Last change on this file since 4457 was 4457, checked in by trackstand2, 10 years ago

Add mbkmax and z-first opt to dynhpg (hpg_sco)

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