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

source: branches/dev_r2586_dynamic_mem/NEMOGCM/NEMO/OPA_SRC/DYN/dynhpg.F90 @ 2636

Last change on this file since 2636 was 2636, checked in by gm, 13 years ago

dynamic mem: #785 ; move ctl_stop & warn in lib_mpp to avoid a circular dependency + ctl_stop improvment

  • Property svn:keywords set to Id
File size: 53.3 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 in_out_manager  ! I/O manager
34   USE trdmod          ! ocean dynamics trends
35   USE trdmod_oce      ! ocean variables trends
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   ! use ta as 3D workspace
196      USE oce, ONLY :   zhpj => sa   ! use sa as 3D workspace
197      !!
198      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
199      !!
200      INTEGER  ::   ji, jj, jk       ! dummy loop indices
201      REAL(wp) ::   zcoef0, zcoef1   ! temporary scalars
202      !!----------------------------------------------------------------------
203     
204      IF( kt == nit000 ) THEN
205         IF(lwp) WRITE(numout,*)
206         IF(lwp) WRITE(numout,*) 'dyn:hpg_zco : hydrostatic pressure gradient trend'
207         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   z-coordinate case '
208      ENDIF
209     
210      ! Local constant initialization
211      zcoef0 = - grav * 0.5_wp
212
213      ! Surface value
214      DO jj = 2, jpjm1
215         DO ji = fs_2, fs_jpim1   ! vector opt.
216            zcoef1 = zcoef0 * fse3w(ji,jj,1)
217            ! hydrostatic pressure gradient
218            zhpi(ji,jj,1) = zcoef1 * ( rhd(ji+1,jj,1) - rhd(ji,jj,1) ) / e1u(ji,jj)
219            zhpj(ji,jj,1) = zcoef1 * ( rhd(ji,jj+1,1) - rhd(ji,jj,1) ) / e2v(ji,jj)
220            ! add to the general momentum trend
221            ua(ji,jj,1) = ua(ji,jj,1) + zhpi(ji,jj,1)
222            va(ji,jj,1) = va(ji,jj,1) + zhpj(ji,jj,1)
223         END DO
224      END DO
225      !
226      ! interior value (2=<jk=<jpkm1)
227      DO jk = 2, jpkm1
228         DO jj = 2, jpjm1
229            DO ji = fs_2, fs_jpim1   ! vector opt.
230               zcoef1 = zcoef0 * fse3w(ji,jj,jk)
231               ! hydrostatic pressure gradient
232               zhpi(ji,jj,jk) = zhpi(ji,jj,jk-1)   &
233                  &           + zcoef1 * (  ( rhd(ji+1,jj,jk)+rhd(ji+1,jj,jk-1) )   &
234                  &                       - ( rhd(ji  ,jj,jk)+rhd(ji  ,jj,jk-1) )  ) / e1u(ji,jj)
235
236               zhpj(ji,jj,jk) = zhpj(ji,jj,jk-1)   &
237                  &           + zcoef1 * (  ( rhd(ji,jj+1,jk)+rhd(ji,jj+1,jk-1) )   &
238                  &                       - ( rhd(ji,jj,  jk)+rhd(ji,jj  ,jk-1) )  ) / e2v(ji,jj)
239               ! add to the general momentum trend
240               ua(ji,jj,jk) = ua(ji,jj,jk) + zhpi(ji,jj,jk)
241               va(ji,jj,jk) = va(ji,jj,jk) + zhpj(ji,jj,jk)
242            END DO
243         END DO
244      END DO
245      !
246   END SUBROUTINE hpg_zco
247
248
249   SUBROUTINE hpg_zps( kt )
250      !!---------------------------------------------------------------------
251      !!                 ***  ROUTINE hpg_zps  ***
252      !!                   
253      !! ** Method  :   z-coordinate plus partial steps case.  blahblah...
254      !!
255      !! ** Action  : - Update (ua,va) with the now hydrastatic pressure trend
256      !!----------------------------------------------------------------------
257      USE oce, ONLY :   zhpi => ta   ! use ta as 3D workspace
258      USE oce, ONLY :   zhpj => sa   ! use sa as 3D workspace
259      !!
260      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
261      !!
262      INTEGER  ::   ji, jj, jk                       ! dummy loop indices
263      INTEGER  ::   iku, ikv                         ! temporary integers
264      REAL(wp) ::   zcoef0, zcoef1, zcoef2, zcoef3   ! temporary scalars
265      !!----------------------------------------------------------------------
266
267      IF( kt == nit000 ) THEN
268         IF(lwp) WRITE(numout,*)
269         IF(lwp) WRITE(numout,*) 'dyn:hpg_zps : hydrostatic pressure gradient trend'
270         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   z-coordinate with partial steps - vector optimization'
271      ENDIF
272
273      ! Local constant initialization
274      zcoef0 = - grav * 0.5_wp
275
276      !  Surface value (also valid in partial step case)
277      DO jj = 2, jpjm1
278         DO ji = fs_2, fs_jpim1   ! vector opt.
279            zcoef1 = zcoef0 * fse3w(ji,jj,1)
280            ! hydrostatic pressure gradient
281            zhpi(ji,jj,1) = zcoef1 * ( rhd(ji+1,jj  ,1) - rhd(ji,jj,1) ) / e1u(ji,jj)
282            zhpj(ji,jj,1) = zcoef1 * ( rhd(ji  ,jj+1,1) - rhd(ji,jj,1) ) / e2v(ji,jj)
283            ! add to the general momentum trend
284            ua(ji,jj,1) = ua(ji,jj,1) + zhpi(ji,jj,1)
285            va(ji,jj,1) = va(ji,jj,1) + zhpj(ji,jj,1)
286         END DO
287      END DO
288
289      ! interior value (2=<jk=<jpkm1)
290      DO jk = 2, jpkm1
291         DO jj = 2, jpjm1
292            DO ji = fs_2, fs_jpim1   ! vector opt.
293               zcoef1 = zcoef0 * fse3w(ji,jj,jk)
294               ! hydrostatic pressure gradient
295               zhpi(ji,jj,jk) = zhpi(ji,jj,jk-1)   &
296                  &           + zcoef1 * (  ( rhd(ji+1,jj,jk) + rhd(ji+1,jj,jk-1) )   &
297                  &                       - ( rhd(ji  ,jj,jk) + rhd(ji  ,jj,jk-1) )  ) / e1u(ji,jj)
298
299               zhpj(ji,jj,jk) = zhpj(ji,jj,jk-1)   &
300                  &           + zcoef1 * (  ( rhd(ji,jj+1,jk) + rhd(ji,jj+1,jk-1) )   &
301                  &                       - ( rhd(ji,jj,  jk) + rhd(ji,jj  ,jk-1) )  ) / e2v(ji,jj)
302               ! add to the general momentum trend
303               ua(ji,jj,jk) = ua(ji,jj,jk) + zhpi(ji,jj,jk)
304               va(ji,jj,jk) = va(ji,jj,jk) + zhpj(ji,jj,jk)
305            END DO
306         END DO
307      END DO
308
309      ! partial steps correction at the last level  (use gru & grv computed in zpshde.F90)
310# if defined key_vectopt_loop
311         jj = 1
312         DO ji = jpi+2, jpij-jpi-1   ! vector opt. (forced unrolling)
313# else
314      DO jj = 2, jpjm1
315         DO ji = 2, jpim1
316# endif
317            iku = mbku(ji,jj)
318            ikv = mbkv(ji,jj)
319            zcoef2 = zcoef0 * MIN( fse3w(ji,jj,iku), fse3w(ji+1,jj  ,iku) )
320            zcoef3 = zcoef0 * MIN( fse3w(ji,jj,ikv), fse3w(ji  ,jj+1,ikv) )
321            IF( iku > 1 ) THEN            ! on i-direction (level 2 or more)
322               ua  (ji,jj,iku) = ua(ji,jj,iku) - zhpi(ji,jj,iku)         ! subtract old value
323               zhpi(ji,jj,iku) = zhpi(ji,jj,iku-1)                   &   ! compute the new one
324                  &            + zcoef2 * ( rhd(ji+1,jj,iku-1) - rhd(ji,jj,iku-1) + gru(ji,jj) ) / e1u(ji,jj)
325               ua  (ji,jj,iku) = ua(ji,jj,iku) + zhpi(ji,jj,iku)         ! add the new one to the general momentum trend
326            ENDIF
327            IF( ikv > 1 ) THEN            ! on j-direction (level 2 or more)
328               va  (ji,jj,ikv) = va(ji,jj,ikv) - zhpj(ji,jj,ikv)         ! subtract old value
329               zhpj(ji,jj,ikv) = zhpj(ji,jj,ikv-1)                   &   ! compute the new one
330                  &            + zcoef3 * ( rhd(ji,jj+1,ikv-1) - rhd(ji,jj,ikv-1) + grv(ji,jj) ) / e2v(ji,jj)
331               va  (ji,jj,ikv) = va(ji,jj,ikv) + zhpj(ji,jj,ikv)         ! add the new one to the general momentum trend
332            ENDIF
333# if ! defined key_vectopt_loop
334         END DO
335# endif
336      END DO
337      !
338   END SUBROUTINE hpg_zps
339
340
341   SUBROUTINE hpg_sco( kt )
342      !!---------------------------------------------------------------------
343      !!                  ***  ROUTINE hpg_sco  ***
344      !!
345      !! ** Method  :   s-coordinate case. Jacobian scheme.
346      !!      The now hydrostatic pressure gradient at a given level, jk,
347      !!      is computed by taking the vertical integral of the in-situ
348      !!      density gradient along the model level from the suface to that
349      !!      level. s-coordinates (ln_sco): a corrective term is added
350      !!      to the horizontal pressure gradient :
351      !!         zhpi = grav .....  + 1/e1u mi(rhd) di[ grav dep3w ]
352      !!         zhpj = grav .....  + 1/e2v mj(rhd) dj[ grav dep3w ]
353      !!      add it to the general momentum trend (ua,va).
354      !!         ua = ua - 1/e1u * zhpi
355      !!         va = va - 1/e2v * zhpj
356      !!
357      !! ** Action : - Update (ua,va) with the now hydrastatic pressure trend
358      !!----------------------------------------------------------------------
359      USE oce, ONLY :   zhpi => ta   ! use ta as 3D workspace
360      USE oce, ONLY :   zhpj => sa   ! use sa as 3D workspace
361      !!
362      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
363      !!
364      INTEGER  ::   ji, jj, jk                 ! dummy loop indices
365      REAL(wp) ::   zcoef0, zuap, zvap, znad   ! temporary scalars
366      !!----------------------------------------------------------------------
367
368      IF( kt == nit000 ) THEN
369         IF(lwp) WRITE(numout,*)
370         IF(lwp) WRITE(numout,*) 'dyn:hpg_sco : hydrostatic pressure gradient trend'
371         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   s-coordinate case, OPA original scheme used'
372      ENDIF
373
374      ! Local constant initialization
375      zcoef0 = - grav * 0.5_wp
376      ! To use density and not density anomaly
377      IF ( lk_vvl ) THEN   ;     znad = 1._wp          ! Variable volume
378      ELSE                 ;     znad = 0._wp         ! Fixed volume
379      ENDIF
380
381      ! Surface value
382      DO jj = 2, jpjm1
383         DO ji = fs_2, fs_jpim1   ! vector opt.   
384            ! hydrostatic pressure gradient along s-surfaces
385            zhpi(ji,jj,1) = zcoef0 / e1u(ji,jj) * ( fse3w(ji+1,jj  ,1) * ( znad + rhd(ji+1,jj  ,1) )   &
386               &                                  - fse3w(ji  ,jj  ,1) * ( znad + rhd(ji  ,jj  ,1) ) )
387            zhpj(ji,jj,1) = zcoef0 / e2v(ji,jj) * ( fse3w(ji  ,jj+1,1) * ( znad + rhd(ji  ,jj+1,1) )   &
388               &                                  - fse3w(ji  ,jj  ,1) * ( znad + rhd(ji  ,jj  ,1) ) )
389            ! s-coordinate pressure gradient correction
390            zuap = -zcoef0 * ( rhd   (ji+1,jj,1) + rhd   (ji,jj,1) + 2._wp * znad )   &
391               &           * ( fsde3w(ji+1,jj,1) - fsde3w(ji,jj,1) ) / e1u(ji,jj)
392            zvap = -zcoef0 * ( rhd   (ji,jj+1,1) + rhd   (ji,jj,1) + 2._wp * znad )   &
393               &           * ( fsde3w(ji,jj+1,1) - fsde3w(ji,jj,1) ) / e2v(ji,jj)
394            ! add to the general momentum trend
395            ua(ji,jj,1) = ua(ji,jj,1) + zhpi(ji,jj,1) + zuap
396            va(ji,jj,1) = va(ji,jj,1) + zhpj(ji,jj,1) + zvap
397         END DO 
398      END DO   
399           
400      ! interior value (2=<jk=<jpkm1)
401      DO jk = 2, jpkm1                                 
402         DO jj = 2, jpjm1     
403            DO ji = fs_2, fs_jpim1   ! vector opt.     
404               ! hydrostatic pressure gradient along s-surfaces
405               zhpi(ji,jj,jk) = zhpi(ji,jj,jk-1) + zcoef0 / e1u(ji,jj)   & 
406                  &           * (  fse3w(ji+1,jj,jk) * ( rhd(ji+1,jj,jk) + rhd(ji+1,jj,jk-1) + 2*znad )   & 
407                  &              - fse3w(ji  ,jj,jk) * ( rhd(ji  ,jj,jk) + rhd(ji  ,jj,jk-1) + 2*znad )  )
408               zhpj(ji,jj,jk) = zhpj(ji,jj,jk-1) + zcoef0 / e2v(ji,jj)   &
409                  &           * (  fse3w(ji,jj+1,jk) * ( rhd(ji,jj+1,jk) + rhd(ji,jj+1,jk-1) + 2*znad )   &
410                  &              - fse3w(ji,jj  ,jk) * ( rhd(ji,jj,  jk) + rhd(ji,jj  ,jk-1) + 2*znad )  )
411               ! s-coordinate pressure gradient correction
412               zuap = -zcoef0 * ( rhd   (ji+1,jj  ,jk) + rhd   (ji,jj,jk) + 2._wp * znad )   &
413                  &           * ( fsde3w(ji+1,jj  ,jk) - fsde3w(ji,jj,jk) ) / e1u(ji,jj)
414               zvap = -zcoef0 * ( rhd   (ji  ,jj+1,jk) + rhd   (ji,jj,jk) + 2._wp * znad )   &
415                  &           * ( fsde3w(ji  ,jj+1,jk) - fsde3w(ji,jj,jk) ) / e2v(ji,jj)
416               ! add to the general momentum trend
417               ua(ji,jj,jk) = ua(ji,jj,jk) + zhpi(ji,jj,jk) + zuap
418               va(ji,jj,jk) = va(ji,jj,jk) + zhpj(ji,jj,jk) + zvap
419            END DO
420         END DO
421      END DO
422      !
423   END SUBROUTINE hpg_sco
424
425
426   SUBROUTINE hpg_hel( kt )
427      !!---------------------------------------------------------------------
428      !!                  ***  ROUTINE hpg_hel  ***
429      !!
430      !! ** Method  :   s-coordinate case.
431      !!      The now hydrostatic pressure gradient at a given level
432      !!      jk is computed by taking the vertical integral of the in-situ
433      !!      density gradient along the model level from the suface to that
434      !!      level. s-coordinates (ln_sco): a corrective term is added
435      !!      to the horizontal pressure gradient :
436      !!         zhpi = grav .....  + 1/e1u mi(rhd) di[ grav dep3w ]
437      !!         zhpj = grav .....  + 1/e2v mj(rhd) dj[ grav dep3w ]
438      !!      add it to the general momentum trend (ua,va).
439      !!         ua = ua - 1/e1u * zhpi
440      !!         va = va - 1/e2v * zhpj
441      !!
442      !! ** Action : - Update (ua,va) with the now hydrastatic pressure trend
443      !!             - Save the trend (l_trddyn=T)
444      !!----------------------------------------------------------------------
445      USE oce, ONLY :   zhpi => ta   ! use ta as 3D workspace
446      USE oce, ONLY :   zhpj => sa   ! use sa as 3D workspace
447      !!
448      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
449      !!
450      INTEGER  ::   ji, jj, jk           ! dummy loop indices
451      REAL(wp) ::   zcoef0, zuap, zvap   ! temporary scalars
452      !!----------------------------------------------------------------------
453
454      IF( kt == nit000 ) THEN
455         IF(lwp) WRITE(numout,*)
456         IF(lwp) WRITE(numout,*) 'dyn:hpg_hel : hydrostatic pressure gradient trend'
457         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   s-coordinate case, helsinki modified scheme'
458      ENDIF
459
460      ! Local constant initialization
461      zcoef0 = - grav * 0.5_wp
462 
463      ! Surface value
464      DO jj = 2, jpjm1
465         DO ji = fs_2, fs_jpim1   ! vector opt.
466            ! hydrostatic pressure gradient along s-surfaces
467            zhpi(ji,jj,1) = zcoef0 / e1u(ji,jj) * ( fse3t(ji+1,jj  ,1) * rhd(ji+1,jj  ,1)  &
468               &                                  - fse3t(ji  ,jj  ,1) * rhd(ji  ,jj  ,1) )
469            zhpj(ji,jj,1) = zcoef0 / e2v(ji,jj) * ( fse3t(ji  ,jj+1,1) * rhd(ji  ,jj+1,1)  &
470               &                                  - fse3t(ji  ,jj  ,1) * rhd(ji  ,jj  ,1) )
471            ! s-coordinate pressure gradient correction
472            zuap = -zcoef0 * ( rhd   (ji+1,jj,1) + rhd   (ji,jj,1) )   &
473               &           * ( fsdept(ji+1,jj,1) - fsdept(ji,jj,1) ) / e1u(ji,jj)
474            zvap = -zcoef0 * ( rhd   (ji,jj+1,1) + rhd   (ji,jj,1) )   &
475               &           * ( fsdept(ji,jj+1,1) - fsdept(ji,jj,1) ) / e2v(ji,jj)
476            ! add to the general momentum trend
477            ua(ji,jj,1) = ua(ji,jj,1) + zhpi(ji,jj,1) + zuap
478            va(ji,jj,1) = va(ji,jj,1) + zhpj(ji,jj,1) + zvap
479         END DO
480      END DO
481      !
482      ! interior value (2=<jk=<jpkm1)
483      DO jk = 2, jpkm1
484         DO jj = 2, jpjm1
485            DO ji = fs_2, fs_jpim1   ! vector opt.
486               ! hydrostatic pressure gradient along s-surfaces
487               zhpi(ji,jj,jk) = zhpi(ji,jj,jk-1) &
488                  &           +  zcoef0 / e1u(ji,jj) * ( fse3t(ji+1,jj,jk  ) * rhd(ji+1,jj,jk)     &
489                  &                                     -fse3t(ji  ,jj,jk  ) * rhd(ji  ,jj,jk)   ) &
490                  &           +  zcoef0 / e1u(ji,jj) * ( fse3t(ji+1,jj,jk-1) * rhd(ji+1,jj,jk-1)   &
491                  &                                     -fse3t(ji  ,jj,jk-1) * rhd(ji  ,jj,jk-1) )
492               zhpj(ji,jj,jk) = zhpj(ji,jj,jk-1) &
493                  &           +  zcoef0 / e2v(ji,jj) * ( fse3t(ji,jj+1,jk  ) * rhd(ji,jj+1,jk)     &
494                  &                                     -fse3t(ji,jj  ,jk  ) * rhd(ji,jj,  jk)   ) &
495                  &           +  zcoef0 / e2v(ji,jj) * ( fse3t(ji,jj+1,jk-1) * rhd(ji,jj+1,jk-1)   &
496                  &                                     -fse3t(ji,jj  ,jk-1) * rhd(ji,jj,  jk-1) )
497               ! s-coordinate pressure gradient correction
498               zuap = - zcoef0 * ( rhd   (ji+1,jj,jk) + rhd   (ji,jj,jk) )   &
499                  &            * ( fsdept(ji+1,jj,jk) - fsdept(ji,jj,jk) ) / e1u(ji,jj)
500               zvap = - zcoef0 * ( rhd   (ji,jj+1,jk) + rhd   (ji,jj,jk) )   &
501                  &            * ( fsdept(ji,jj+1,jk) - fsdept(ji,jj,jk) ) / e2v(ji,jj)
502               ! add to the general momentum trend
503               ua(ji,jj,jk) = ua(ji,jj,jk) + zhpi(ji,jj,jk) + zuap
504               va(ji,jj,jk) = va(ji,jj,jk) + zhpj(ji,jj,jk) + zvap
505            END DO
506         END DO
507      END DO
508      !
509   END SUBROUTINE hpg_hel
510
511
512   SUBROUTINE hpg_wdj( kt )
513      !!---------------------------------------------------------------------
514      !!                  ***  ROUTINE hpg_wdj  ***
515      !!
516      !! ** Method  :   Weighted Density Jacobian (wdj) scheme (song 1998)
517      !!      The weighting coefficients from the namelist parameter rn_gamma
518      !!      (alpha=0.5-rn_gamma ; beta=1-alpha=0.5+rn_gamma
519      !!
520      !! Reference : Song, Mon. Wea. Rev., 126, 3213-3230, 1998.
521      !!----------------------------------------------------------------------
522      USE oce, ONLY :   zhpi => ta   ! use ta as 3D workspace
523      USE oce, ONLY :   zhpj => sa   ! use sa as 3D workspace
524      !!
525      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
526      !!
527      INTEGER  ::   ji, jj, jk           ! dummy loop indices
528      REAL(wp) ::   zcoef0, zuap, zvap   ! temporary scalars
529      REAL(wp) ::   zalph , zbeta        !    "         "
530      !!----------------------------------------------------------------------
531
532      IF( kt == nit000 ) THEN
533         IF(lwp) WRITE(numout,*)
534         IF(lwp) WRITE(numout,*) 'dyn:hpg_wdj : hydrostatic pressure gradient trend'
535         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   Weighted Density Jacobian'
536      ENDIF
537
538      ! Local constant initialization
539      zcoef0 = - grav * 0.5_wp
540      zalph  = 0.5_wp - rn_gamma    ! weighting coefficients (alpha=0.5-rn_gamma
541      zbeta  = 0.5_wp + rn_gamma    !                        (beta =1-alpha=0.5+rn_gamma
542
543      ! Surface value (no ponderation)
544      DO jj = 2, jpjm1
545         DO ji = fs_2, fs_jpim1   ! vector opt.
546            ! hydrostatic pressure gradient along s-surfaces
547            zhpi(ji,jj,1) = zcoef0 / e1u(ji,jj) * (  fse3w(ji+1,jj  ,1) * rhd(ji+1,jj  ,1)   &
548               &                                   - fse3w(ji  ,jj  ,1) * rhd(ji  ,jj  ,1)  )
549            zhpj(ji,jj,1) = zcoef0 / e2v(ji,jj) * (  fse3w(ji  ,jj+1,1) * rhd(ji  ,jj+1,1)   &
550               &                                   - fse3w(ji  ,jj  ,1) * rhd(ji,  jj  ,1)  )
551            ! s-coordinate pressure gradient correction
552            zuap = -zcoef0 * ( rhd   (ji+1,jj,1) + rhd   (ji,jj,1) )   &
553               &           * ( fsde3w(ji+1,jj,1) - fsde3w(ji,jj,1) ) / e1u(ji,jj)
554            zvap = -zcoef0 * ( rhd   (ji,jj+1,1) + rhd   (ji,jj,1) )   &
555               &           * ( fsde3w(ji,jj+1,1) - fsde3w(ji,jj,1) ) / e2v(ji,jj)
556            ! add to the general momentum trend
557            ua(ji,jj,1) = ua(ji,jj,1) + zhpi(ji,jj,1) + zuap
558            va(ji,jj,1) = va(ji,jj,1) + zhpj(ji,jj,1) + zvap
559         END DO
560      END DO
561
562      ! Interior value (2=<jk=<jpkm1) (weighted with zalph & zbeta)
563      DO jk = 2, jpkm1
564         DO jj = 2, jpjm1
565            DO ji = fs_2, fs_jpim1   ! vector opt.
566               zhpi(ji,jj,jk) = zhpi(ji,jj,jk-1) + zcoef0 / e1u(ji,jj)                            &
567                  &           * (   (            fsde3w(ji+1,jj,jk  ) + fsde3w(ji,jj,jk  )        &
568                  &                            - fsde3w(ji+1,jj,jk-1) - fsde3w(ji,jj,jk-1)    )   &
569                  &               * (  zalph * ( rhd   (ji+1,jj,jk-1) - rhd   (ji,jj,jk-1) )      &
570                  &                  + zbeta * ( rhd   (ji+1,jj,jk  ) - rhd   (ji,jj,jk  ) )  )   &
571                  &             -   (            rhd   (ji+1,jj,jk  ) + rhd   (ji,jj,jk  )        &
572                  &                           - rhd   (ji+1,jj,jk-1) - rhd   (ji,jj,jk-1)     )   &
573                  &               * (  zalph * ( fsde3w(ji+1,jj,jk-1) - fsde3w(ji,jj,jk-1) )      &
574                  &                  + zbeta * ( fsde3w(ji+1,jj,jk  ) - fsde3w(ji,jj,jk  ) )  )  )
575               zhpj(ji,jj,jk) = zhpj(ji,jj,jk-1) + zcoef0 / e2v(ji,jj)                            &
576                  &           * (   (           fsde3w(ji,jj+1,jk  ) + fsde3w(ji,jj,jk  )         &
577                  &                           - fsde3w(ji,jj+1,jk-1) - fsde3w(ji,jj,jk-1)     )   &
578                  &               * (  zalph * ( rhd   (ji,jj+1,jk-1) - rhd   (ji,jj,jk-1) )      &
579                  &                  + zbeta * ( rhd   (ji,jj+1,jk  ) - rhd   (ji,jj,jk  ) )  )   &
580                  &             -   (            rhd   (ji,jj+1,jk  ) + rhd   (ji,jj,jk  )        &
581                  &                            - rhd   (ji,jj+1,jk-1) - rhd   (ji,jj,jk-1)    )   &
582                  &               * (  zalph * ( fsde3w(ji,jj+1,jk-1) - fsde3w(ji,jj,jk-1) )      &
583                  &                  + zbeta * ( fsde3w(ji,jj+1,jk  ) - fsde3w(ji,jj,jk  ) )  )  )
584               ! add to the general momentum trend
585               ua(ji,jj,jk) = ua(ji,jj,jk) + zhpi(ji,jj,jk)
586               va(ji,jj,jk) = va(ji,jj,jk) + zhpj(ji,jj,jk)
587            END DO
588         END DO
589      END DO
590      !
591   END SUBROUTINE hpg_wdj
592
593
594   SUBROUTINE hpg_djc( kt )
595      !!---------------------------------------------------------------------
596      !!                  ***  ROUTINE hpg_djc  ***
597      !!
598      !! ** Method  :   Density Jacobian with Cubic polynomial scheme
599      !!
600      !! Reference: Shchepetkin and McWilliams, J. Geophys. Res., 108(C3), 3090, 2003
601      !!----------------------------------------------------------------------
602      USE oce, ONLY :   zhpi => ta   ! use ta as 3D workspace
603      USE oce, ONLY :   zhpj => sa   ! use sa as 3D workspace
604      USE wrk_nemo, ONLY: wrk_in_use, wrk_not_released
605      USE wrk_nemo, ONLY: drhox => wrk_3d_1  , dzx => wrk_3d_2
606      USE wrk_nemo, ONLY: drhou => wrk_3d_3  , dzu => wrk_3d_4 , rho_i => wrk_3d_5
607      USE wrk_nemo, ONLY: drhoy => wrk_3d_6  , dzy => wrk_3d_7
608      USE wrk_nemo, ONLY: drhov => wrk_3d_8  , dzv => wrk_3d_9 , rho_j => wrk_3d_10
609      USE wrk_nemo, ONLY: drhoz => wrk_3d_11 , dzz => wrk_3d_12 
610      USE wrk_nemo, ONLY: drhow => wrk_3d_13 , dzw => wrk_3d_14
611      USE wrk_nemo, ONLY: rho_k => wrk_3d_15
612      !!
613      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
614      !!
615      INTEGER  ::   ji, jj, jk          ! dummy loop indices
616      REAL(wp) ::   zcoef0, zep, cffw   ! temporary scalars
617      REAL(wp) ::   z1_10, cffu, cffx   !    "         "
618      REAL(wp) ::   z1_12, cffv, cffy   !    "         "
619      !!----------------------------------------------------------------------
620
621      IF( wrk_in_use(3, 1,2,3,4,5,6,7,8,9,10,11,12,13,14,15) ) THEN
622         CALL ctl_stop('dyn:hpg_djc : requested workspace arrays unavailable')   ;   RETURN
623      ENDIF
624
625      IF( kt == nit000 ) THEN
626         IF(lwp) WRITE(numout,*)
627         IF(lwp) WRITE(numout,*) 'dyn:hpg_djc : hydrostatic pressure gradient trend'
628         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   s-coordinate case, density Jacobian with cubic polynomial scheme'
629      ENDIF
630
631
632      ! Local constant initialization
633      zcoef0 = - grav * 0.5_wp
634      z1_10  = 1._wp / 10._wp
635      z1_12  = 1._wp / 12._wp
636
637      !----------------------------------------------------------------------------------------
638      !  compute and store in provisional arrays elementary vertical and horizontal differences
639      !----------------------------------------------------------------------------------------
640
641!!bug gm   Not a true bug, but... dzz=e3w  for dzx, dzy verify what it is really
642
643      DO jk = 2, jpkm1
644         DO jj = 2, jpjm1
645            DO ji = fs_2, fs_jpim1   ! vector opt.
646               drhoz(ji,jj,jk) = rhd   (ji  ,jj  ,jk) - rhd   (ji,jj,jk-1)
647               dzz  (ji,jj,jk) = fsde3w(ji  ,jj  ,jk) - fsde3w(ji,jj,jk-1)
648               drhox(ji,jj,jk) = rhd   (ji+1,jj  ,jk) - rhd   (ji,jj,jk  )
649               dzx  (ji,jj,jk) = fsde3w(ji+1,jj  ,jk) - fsde3w(ji,jj,jk  )
650               drhoy(ji,jj,jk) = rhd   (ji  ,jj+1,jk) - rhd   (ji,jj,jk  )
651               dzy  (ji,jj,jk) = fsde3w(ji  ,jj+1,jk) - fsde3w(ji,jj,jk  )
652            END DO
653         END DO
654      END DO
655
656      !-------------------------------------------------------------------------
657      ! compute harmonic averages using eq. 5.18
658      !-------------------------------------------------------------------------
659      zep = 1.e-15
660
661!!bug  gm  drhoz not defined at level 1 and used (jk-1 with jk=2)
662!!bug  gm  idem for drhox, drhoy et ji=jpi and jj=jpj
663
664      DO jk = 2, jpkm1
665         DO jj = 2, jpjm1
666            DO ji = fs_2, fs_jpim1   ! vector opt.
667               cffw = 2._wp * drhoz(ji  ,jj  ,jk) * drhoz(ji,jj,jk-1)
668
669               cffu = 2._wp * drhox(ji+1,jj  ,jk) * drhox(ji,jj,jk  )
670               cffx = 2._wp * dzx  (ji+1,jj  ,jk) * dzx  (ji,jj,jk  )
671 
672               cffv = 2._wp * drhoy(ji  ,jj+1,jk) * drhoy(ji,jj,jk  )
673               cffy = 2._wp * dzy  (ji  ,jj+1,jk) * dzy  (ji,jj,jk  )
674
675               IF( cffw > zep) THEN
676                  drhow(ji,jj,jk) = 2._wp *   drhoz(ji,jj,jk) * drhoz(ji,jj,jk-1)   &
677                     &                    / ( drhoz(ji,jj,jk) + drhoz(ji,jj,jk-1) )
678               ELSE
679                  drhow(ji,jj,jk) = 0._wp
680               ENDIF
681
682               dzw(ji,jj,jk) = 2._wp *   dzz(ji,jj,jk) * dzz(ji,jj,jk-1)   &
683                  &                  / ( dzz(ji,jj,jk) + dzz(ji,jj,jk-1) )
684
685               IF( cffu > zep ) THEN
686                  drhou(ji,jj,jk) = 2._wp *   drhox(ji+1,jj,jk) * drhox(ji,jj,jk)   &
687                     &                    / ( drhox(ji+1,jj,jk) + drhox(ji,jj,jk) )
688               ELSE
689                  drhou(ji,jj,jk ) = 0._wp
690               ENDIF
691
692               IF( cffx > zep ) THEN
693                  dzu(ji,jj,jk) = 2._wp *   dzx(ji+1,jj,jk) * dzx(ji,jj,jk)   &
694                     &                  / ( dzx(ji+1,jj,jk) + dzx(ji,jj,jk) )
695               ELSE
696                  dzu(ji,jj,jk) = 0._wp
697               ENDIF
698
699               IF( cffv > zep ) THEN
700                  drhov(ji,jj,jk) = 2._wp *   drhoy(ji,jj+1,jk) * drhoy(ji,jj,jk)   &
701                     &                    / ( drhoy(ji,jj+1,jk) + drhoy(ji,jj,jk) )
702               ELSE
703                  drhov(ji,jj,jk) = 0._wp
704               ENDIF
705
706               IF( cffy > zep ) THEN
707                  dzv(ji,jj,jk) = 2._wp *   dzy(ji,jj+1,jk) * dzy(ji,jj,jk)   &
708                     &                  / ( dzy(ji,jj+1,jk) + dzy(ji,jj,jk) )
709               ELSE
710                  dzv(ji,jj,jk) = 0._wp
711               ENDIF
712
713            END DO
714         END DO
715      END DO
716
717      !----------------------------------------------------------------------------------
718      ! apply boundary conditions at top and bottom using 5.36-5.37
719      !----------------------------------------------------------------------------------
720      drhow(:,:, 1 ) = 1.5_wp * ( drhoz(:,:, 2 ) - drhoz(:,:,  1  ) ) - 0.5_wp * drhow(:,:,  2  )
721      drhou(:,:, 1 ) = 1.5_wp * ( drhox(:,:, 2 ) - drhox(:,:,  1  ) ) - 0.5_wp * drhou(:,:,  2  )
722      drhov(:,:, 1 ) = 1.5_wp * ( drhoy(:,:, 2 ) - drhoy(:,:,  1  ) ) - 0.5_wp * drhov(:,:,  2  )
723
724      drhow(:,:,jpk) = 1.5_wp * ( drhoz(:,:,jpk) - drhoz(:,:,jpkm1) ) - 0.5_wp * drhow(:,:,jpkm1)
725      drhou(:,:,jpk) = 1.5_wp * ( drhox(:,:,jpk) - drhox(:,:,jpkm1) ) - 0.5_wp * drhou(:,:,jpkm1)
726      drhov(:,:,jpk) = 1.5_wp * ( drhoy(:,:,jpk) - drhoy(:,:,jpkm1) ) - 0.5_wp * drhov(:,:,jpkm1)
727
728
729      !--------------------------------------------------------------
730      ! Upper half of top-most grid box, compute and store
731      !-------------------------------------------------------------
732
733!!bug gm   :  e3w-de3w = 0.5*e3w  ....  and de3w(2)-de3w(1)=e3w(2) ....   to be verified
734!          true if de3w is really defined as the sum of the e3w scale factors as, it seems to me, it should be
735
736      DO jj = 2, jpjm1
737         DO ji = fs_2, fs_jpim1   ! vector opt.
738            rho_k(ji,jj,1) = -grav * ( fse3w(ji,jj,1) - fsde3w(ji,jj,1) )               &
739               &                   * (  rhd(ji,jj,1)                                    &
740               &                     + 0.5_wp * ( rhd(ji,jj,2) - rhd(ji,jj,1) )         &
741               &                              * ( fse3w (ji,jj,1) - fsde3w(ji,jj,1) )   &
742               &                              / ( fsde3w(ji,jj,2) - fsde3w(ji,jj,1) )  ) 
743         END DO
744      END DO
745
746!!bug gm    : here also, simplification is possible
747!!bug gm    : optimisation: 1/10 and 1/12 the division should be done before the loop
748
749      DO jk = 2, jpkm1
750         DO jj = 2, jpjm1
751            DO ji = fs_2, fs_jpim1   ! vector opt.
752
753               rho_k(ji,jj,jk) = zcoef0 * ( rhd   (ji,jj,jk) + rhd   (ji,jj,jk-1) )                                   &
754                  &                     * ( fsde3w(ji,jj,jk) - fsde3w(ji,jj,jk-1) )                                   &
755                  &            - grav * z1_10 * (                                                                     &
756                  &     ( drhow (ji,jj,jk) - drhow (ji,jj,jk-1) )                                                     &
757                  &   * ( fsde3w(ji,jj,jk) - fsde3w(ji,jj,jk-1) - z1_12 * ( dzw  (ji,jj,jk) + dzw  (ji,jj,jk-1) ) )   &
758                  &   - ( dzw   (ji,jj,jk) - dzw   (ji,jj,jk-1) )                                                     &
759                  &   * ( rhd   (ji,jj,jk) - rhd   (ji,jj,jk-1) - z1_12 * ( drhow(ji,jj,jk) + drhow(ji,jj,jk-1) ) )   &
760                  &                             )
761
762               rho_i(ji,jj,jk) = zcoef0 * ( rhd   (ji+1,jj,jk) + rhd   (ji,jj,jk) )                                   &
763                  &                     * ( fsde3w(ji+1,jj,jk) - fsde3w(ji,jj,jk) )                                   &
764                  &            - grav* z1_10 * (                                                                      &
765                  &     ( drhou (ji+1,jj,jk) - drhou (ji,jj,jk) )                                                     &
766                  &   * ( fsde3w(ji+1,jj,jk) - fsde3w(ji,jj,jk) - z1_12 * ( dzu  (ji+1,jj,jk) + dzu  (ji,jj,jk) ) )   &
767                  &   - ( dzu   (ji+1,jj,jk) - dzu   (ji,jj,jk) )                                                     &
768                  &   * ( rhd   (ji+1,jj,jk) - rhd   (ji,jj,jk) - z1_12 * ( drhou(ji+1,jj,jk) + drhou(ji,jj,jk) ) )   &
769                  &                            )
770
771               rho_j(ji,jj,jk) = zcoef0 * ( rhd   (ji,jj+1,jk) + rhd   (ji,jj,jk) )                                   &
772                  &                     * ( fsde3w(ji,jj+1,jk) - fsde3w(ji,jj,jk) )                                   &
773                  &            - grav* z1_10 * (                                                                      &
774                  &     ( drhov (ji,jj+1,jk) - drhov (ji,jj,jk) )                                                     &
775                  &   * ( fsde3w(ji,jj+1,jk) - fsde3w(ji,jj,jk) - z1_12 * ( dzv  (ji,jj+1,jk) + dzv  (ji,jj,jk) ) )   &
776                  &   - ( dzv   (ji,jj+1,jk) - dzv   (ji,jj,jk) )                                                     &
777                  &   * ( rhd   (ji,jj+1,jk) - rhd   (ji,jj,jk) - z1_12 * ( drhov(ji,jj+1,jk) + drhov(ji,jj,jk) ) )   &
778                  &                            )
779
780            END DO
781         END DO
782      END DO
783      CALL lbc_lnk(rho_k,'W',1.)
784      CALL lbc_lnk(rho_i,'U',1.)
785      CALL lbc_lnk(rho_j,'V',1.)
786
787
788      ! ---------------
789      !  Surface value
790      ! ---------------
791      DO jj = 2, jpjm1
792         DO ji = fs_2, fs_jpim1   ! vector opt.
793            zhpi(ji,jj,1) = ( rho_k(ji+1,jj  ,1) - rho_k(ji,jj,1) - rho_i(ji,jj,1) ) / e1u(ji,jj)
794            zhpj(ji,jj,1) = ( rho_k(ji  ,jj+1,1) - rho_k(ji,jj,1) - rho_j(ji,jj,1) ) / e2v(ji,jj)
795            ! add to the general momentum trend
796            ua(ji,jj,1) = ua(ji,jj,1) + zhpi(ji,jj,1)
797            va(ji,jj,1) = va(ji,jj,1) + zhpj(ji,jj,1)
798         END DO
799      END DO
800
801      ! ----------------
802      !  interior value   (2=<jk=<jpkm1)
803      ! ----------------
804      DO jk = 2, jpkm1
805         DO jj = 2, jpjm1 
806            DO ji = fs_2, fs_jpim1   ! vector opt.
807               ! hydrostatic pressure gradient along s-surfaces
808               zhpi(ji,jj,jk) = zhpi(ji,jj,jk-1)                                &
809                  &           + (  ( rho_k(ji+1,jj,jk) - rho_k(ji,jj,jk  ) )    &
810                  &              - ( rho_i(ji  ,jj,jk) - rho_i(ji,jj,jk-1) )  ) / e1u(ji,jj)
811               zhpj(ji,jj,jk) = zhpj(ji,jj,jk-1)                                &
812                  &           + (  ( rho_k(ji,jj+1,jk) - rho_k(ji,jj,jk  ) )    &
813                  &               -( rho_j(ji,jj  ,jk) - rho_j(ji,jj,jk-1) )  ) / e2v(ji,jj)
814               ! add to the general momentum trend
815               ua(ji,jj,jk) = ua(ji,jj,jk) + zhpi(ji,jj,jk)
816               va(ji,jj,jk) = va(ji,jj,jk) + zhpj(ji,jj,jk)
817            END DO
818         END DO
819      END DO
820      !
821      IF( wrk_not_released(3, 1,2,3,4,5,6,7,8,9,10,11,12,13,14,15) )   &
822         CALL ctl_stop('dyn:hpg_djc : failed to release workspace arrays')
823      !
824   END SUBROUTINE hpg_djc
825
826
827   SUBROUTINE hpg_rot( kt )
828      !!---------------------------------------------------------------------
829      !!                  ***  ROUTINE hpg_rot  ***
830      !!
831      !! ** Method  :   rotated axes scheme (Thiem and Berntsen 2005)
832      !!
833      !! Reference: Thiem & Berntsen, Ocean Modelling, In press, 2005.
834      !!----------------------------------------------------------------------
835      USE oce, ONLY :   zhpi => ta   ! use ta as 3D workspace
836      USE oce, ONLY :   zhpj => sa   ! use sa as 3D workspace
837      USE wrk_nemo, ONLY: wrk_in_use, wrk_not_released
838      USE wrk_nemo, ONLY: zdistr  => wrk_2d_1 , zsina   => wrk_2d_2 , zcosa  => wrk_2d_3
839      USE wrk_nemo, ONLY: zhpiorg => wrk_3d_1 , zhpirot => wrk_3d_2
840      USE wrk_nemo, ONLY: zhpitra => wrk_3d_3 , zhpine  => wrk_3d_4
841      USE wrk_nemo, ONLY: zhpjorg => wrk_3d_5 , zhpjrot => wrk_3d_6
842      USE wrk_nemo, ONLY: zhpjtra => wrk_3d_7 , zhpjne  => wrk_3d_8
843      !!
844      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
845      !!
846      INTEGER  ::   ji, jj, jk          ! dummy loop indices
847      REAL(wp) ::   zforg, zcoef0, zuap, zmskd1, zmskd1m   ! temporary scalar
848      REAL(wp) ::   zfrot        , zvap, zmskd2, zmskd2m   !    "         "
849      !!----------------------------------------------------------------------
850
851      IF( wrk_in_use(2, 1,2,3) .OR.      &
852          wrk_in_use(3, 1,2,3,4,5,6,7,8) ) THEN
853         CALL ctl_stop('dyn:hpg_rot : requested workspace arrays unavailable')   ;   RETURN
854      END IF
855
856      IF( kt == nit000 ) THEN
857         IF(lwp) WRITE(numout,*)
858         IF(lwp) WRITE(numout,*) 'dyn:hpg_rot : hydrostatic pressure gradient trend'
859         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   s-coordinate case, rotated axes scheme used'
860      ENDIF
861
862      ! -------------------------------
863      !  Local constant initialization
864      ! -------------------------------
865      zcoef0 = - grav * 0.5_wp
866      zforg  = 0.95_wp
867      zfrot  = 1._wp - zforg
868
869      ! inverse of the distance between 2 diagonal T-points (defined at F-point) (here zcoef0/distance)
870      zdistr(:,:) = zcoef0 / SQRT( e1f(:,:)*e1f(:,:) + e2f(:,:)*e1f(:,:) )
871
872      ! sinus and cosinus of diagonal angle at F-point
873      zsina(:,:) = ATAN2( e2f(:,:), e1f(:,:) )
874      zcosa(:,:) = COS( zsina(:,:) )
875      zsina(:,:) = SIN( zsina(:,:) )
876
877      ! ---------------
878      !  Surface value
879      ! ---------------
880      ! compute and add to the general trend the pressure gradients along the axes
881      DO jj = 2, jpjm1
882         DO ji = fs_2, fs_jpim1   ! vector opt.
883            ! hydrostatic pressure gradient along s-surfaces
884            zhpiorg(ji,jj,1) = zcoef0 / e1u(ji,jj) * (  fse3t(ji+1,jj,1) * rhd(ji+1,jj,1)   &
885               &                                      - fse3t(ji  ,jj,1) * rhd(ji  ,jj,1)  )
886            zhpjorg(ji,jj,1) = zcoef0 / e2v(ji,jj) * (  fse3t(ji,jj+1,1) * rhd(ji,jj+1,1)   &
887               &                                      - fse3t(ji,jj  ,1) * rhd(ji,jj  ,1)  )
888            ! s-coordinate pressure gradient correction
889            zuap = -zcoef0 * ( rhd   (ji+1,jj  ,1) + rhd   (ji,jj,1) )   &
890               &           * ( fsdept(ji+1,jj  ,1) - fsdept(ji,jj,1) ) / e1u(ji,jj)
891            zvap = -zcoef0 * ( rhd   (ji  ,jj+1,1) + rhd   (ji,jj,1) )   &
892               &           * ( fsdept(ji  ,jj+1,1) - fsdept(ji,jj,1) ) / e2v(ji,jj)
893            ! add to the general momentum trend
894            ua(ji,jj,1) = ua(ji,jj,1) + zforg * ( zhpiorg(ji,jj,1) + zuap )
895            va(ji,jj,1) = va(ji,jj,1) + zforg * ( zhpjorg(ji,jj,1) + zvap )
896         END DO
897      END DO
898
899      ! compute the pressure gradients in the diagonal directions
900      DO jj = 1, jpjm1
901         DO ji = 1, fs_jpim1   ! vector opt.
902            zmskd1 = tmask(ji+1,jj+1,1) * tmask(ji  ,jj,1)      ! mask in the 1st diagnonal
903            zmskd2 = tmask(ji  ,jj+1,1) * tmask(ji+1,jj,1)      ! mask in the 2nd diagnonal
904            ! hydrostatic pressure gradient along s-surfaces
905            zhpitra(ji,jj,1) = zdistr(ji,jj) * zmskd1 * (  fse3t(ji+1,jj+1,1) * rhd(ji+1,jj+1,1)   &
906               &                                         - fse3t(ji  ,jj  ,1) * rhd(ji  ,jj  ,1)  )
907            zhpjtra(ji,jj,1) = zdistr(ji,jj) * zmskd2 * (  fse3t(ji  ,jj+1,1) * rhd(ji  ,jj+1,1)   &
908               &                                         - fse3t(ji+1,jj  ,1) * rhd(ji+1,jj  ,1)  )
909            ! s-coordinate pressure gradient correction
910            zuap = -zdistr(ji,jj) * zmskd1 * ( rhd   (ji+1,jj+1,1) + rhd   (ji  ,jj,1) )   &
911               &                           * ( fsdept(ji+1,jj+1,1) - fsdept(ji  ,jj,1) )
912            zvap = -zdistr(ji,jj) * zmskd2 * ( rhd   (ji  ,jj+1,1) + rhd   (ji+1,jj,1) )   &
913               &                           * ( fsdept(ji  ,jj+1,1) - fsdept(ji+1,jj,1) )
914            ! back rotation
915            zhpine(ji,jj,1) = zcosa(ji,jj) * ( zhpitra(ji,jj,1) + zuap )   &
916               &            - zsina(ji,jj) * ( zhpjtra(ji,jj,1) + zvap )
917            zhpjne(ji,jj,1) = zsina(ji,jj) * ( zhpitra(ji,jj,1) + zuap )   &
918               &            + zcosa(ji,jj) * ( zhpjtra(ji,jj,1) + zvap )
919         END DO
920      END DO
921
922      ! interpolate and add to the general trend the diagonal gradient
923      DO jj = 2, jpjm1
924         DO ji = fs_2, fs_jpim1   ! vector opt.
925            ! averaging
926            zhpirot(ji,jj,1) = 0.5 * ( zhpine(ji,jj,1) + zhpine(ji  ,jj-1,1) )
927            zhpjrot(ji,jj,1) = 0.5 * ( zhpjne(ji,jj,1) + zhpjne(ji-1,jj  ,1) )
928            ! add to the general momentum trend
929            ua(ji,jj,1) = ua(ji,jj,1) + zfrot * zhpirot(ji,jj,1) 
930            va(ji,jj,1) = va(ji,jj,1) + zfrot * zhpjrot(ji,jj,1) 
931         END DO
932      END DO
933
934      ! -----------------
935      ! 2. interior value (2=<jk=<jpkm1)
936      ! -----------------
937      ! compute and add to the general trend the pressure gradients along the axes
938      DO jk = 2, jpkm1
939         DO jj = 2, jpjm1
940            DO ji = fs_2, fs_jpim1   ! vector opt.
941               ! hydrostatic pressure gradient along s-surfaces
942               zhpiorg(ji,jj,jk) = zhpiorg(ji,jj,jk-1)                                                 &
943                  &              +  zcoef0 / e1u(ji,jj) * (  fse3t(ji+1,jj,jk  ) * rhd(ji+1,jj,jk  )   &
944                  &                                        - fse3t(ji  ,jj,jk  ) * rhd(ji  ,jj,jk  )   &
945                  &                                        + fse3t(ji+1,jj,jk-1) * rhd(ji+1,jj,jk-1)   &
946                  &                                        - fse3t(ji  ,jj,jk-1) * rhd(ji  ,jj,jk-1)  )
947               zhpjorg(ji,jj,jk) = zhpjorg(ji,jj,jk-1)                                                 &
948                  &              +  zcoef0 / e2v(ji,jj) * (  fse3t(ji,jj+1,jk  ) * rhd(ji,jj+1,jk  )   &
949                  &                                        - fse3t(ji,jj  ,jk  ) * rhd(ji,jj,  jk  )   &
950                  &                                        + fse3t(ji,jj+1,jk-1) * rhd(ji,jj+1,jk-1)   &
951                  &                                        - fse3t(ji,jj  ,jk-1) * rhd(ji,jj,  jk-1)  )
952               ! s-coordinate pressure gradient correction
953               zuap = - zcoef0 * ( rhd   (ji+1,jj  ,jk) + rhd   (ji,jj,jk) )   &
954                  &            * ( fsdept(ji+1,jj  ,jk) - fsdept(ji,jj,jk) ) / e1u(ji,jj)
955               zvap = - zcoef0 * ( rhd   (ji  ,jj+1,jk) + rhd   (ji,jj,jk) )   &
956                  &            * ( fsdept(ji  ,jj+1,jk) - fsdept(ji,jj,jk) ) / e2v(ji,jj)
957               ! add to the general momentum trend
958               ua(ji,jj,jk) = ua(ji,jj,jk) + zforg*( zhpiorg(ji,jj,jk) + zuap )
959               va(ji,jj,jk) = va(ji,jj,jk) + zforg*( zhpjorg(ji,jj,jk) + zvap )
960            END DO
961         END DO
962      END DO
963
964      ! compute the pressure gradients in the diagonal directions
965      DO jk = 2, jpkm1
966         DO jj = 1, jpjm1
967            DO ji = 1, fs_jpim1   ! vector opt.
968               zmskd1  = tmask(ji+1,jj+1,jk  ) * tmask(ji  ,jj,jk  )      ! level jk   mask in the 1st diagnonal
969               zmskd1m = tmask(ji+1,jj+1,jk-1) * tmask(ji  ,jj,jk-1)      ! level jk-1    "               "     
970               zmskd2  = tmask(ji  ,jj+1,jk  ) * tmask(ji+1,jj,jk  )      ! level jk   mask in the 2nd diagnonal
971               zmskd2m = tmask(ji  ,jj+1,jk-1) * tmask(ji+1,jj,jk-1)      ! level jk-1    "               "     
972               ! hydrostatic pressure gradient along s-surfaces
973               zhpitra(ji,jj,jk) = zhpitra(ji,jj,jk-1)                                                       &
974                  &              + zdistr(ji,jj) * zmskd1  * ( fse3t(ji+1,jj+1,jk  ) * rhd(ji+1,jj+1,jk)     &
975                  &                                           -fse3t(ji  ,jj  ,jk  ) * rhd(ji  ,jj  ,jk) )   &
976                  &              + zdistr(ji,jj) * zmskd1m * ( fse3t(ji+1,jj+1,jk-1) * rhd(ji+1,jj+1,jk-1)   &
977                  &                                           -fse3t(ji  ,jj  ,jk-1) * rhd(ji  ,jj  ,jk-1) )
978               zhpjtra(ji,jj,jk) = zhpjtra(ji,jj,jk-1)                                                       &
979                  &              + zdistr(ji,jj) * zmskd2  * ( fse3t(ji  ,jj+1,jk  ) * rhd(ji  ,jj+1,jk)     &
980                  &                                           -fse3t(ji+1,jj  ,jk  ) * rhd(ji+1,jj,  jk) )   &
981                  &              + zdistr(ji,jj) * zmskd2m * ( fse3t(ji  ,jj+1,jk-1) * rhd(ji  ,jj+1,jk-1)   &
982                  &                                           -fse3t(ji+1,jj  ,jk-1) * rhd(ji+1,jj,  jk-1) )
983               ! s-coordinate pressure gradient correction
984               zuap = - zdistr(ji,jj) * zmskd1 * ( rhd   (ji+1,jj+1,jk) + rhd   (ji  ,jj,jk) )   &
985                  &                            * ( fsdept(ji+1,jj+1,jk) - fsdept(ji  ,jj,jk) )
986               zvap = - zdistr(ji,jj) * zmskd2 * ( rhd   (ji  ,jj+1,jk) + rhd   (ji+1,jj,jk) )   &
987                  &                            * ( fsdept(ji  ,jj+1,jk) - fsdept(ji+1,jj,jk) )
988               ! back rotation
989               zhpine(ji,jj,jk) = zcosa(ji,jj) * ( zhpitra(ji,jj,jk) + zuap )   &
990                  &             - zsina(ji,jj) * ( zhpjtra(ji,jj,jk) + zvap )
991               zhpjne(ji,jj,jk) = zsina(ji,jj) * ( zhpitra(ji,jj,jk) + zuap )   &
992                  &             + zcosa(ji,jj) * ( zhpjtra(ji,jj,jk) + zvap )
993            END DO
994         END DO
995      END DO
996
997      ! interpolate and add to the general trend
998      DO jk = 2, jpkm1
999         DO jj = 2, jpjm1
1000            DO ji = fs_2, fs_jpim1   ! vector opt.
1001               ! averaging
1002               zhpirot(ji,jj,jk) = 0.5 * ( zhpine(ji,jj,jk) + zhpine(ji  ,jj-1,jk) )
1003               zhpjrot(ji,jj,jk) = 0.5 * ( zhpjne(ji,jj,jk) + zhpjne(ji-1,jj  ,jk) )
1004               ! add to the general momentum trend
1005               ua(ji,jj,jk) = ua(ji,jj,jk) + zfrot * zhpirot(ji,jj,jk) 
1006               va(ji,jj,jk) = va(ji,jj,jk) + zfrot * zhpjrot(ji,jj,jk) 
1007            END DO
1008         END DO
1009      END DO
1010      !
1011      IF( wrk_not_released(2, 1,2,3)  .OR.     &
1012          wrk_not_released(3, 1,2,3,4,5,6,7,8) )   CALL ctl_stop('dyn:hpg_rot : failed to release workspace arrays')
1013      !
1014   END SUBROUTINE hpg_rot
1015
1016   !!======================================================================
1017END MODULE dynhpg
Note: See TracBrowser for help on using the repository browser.