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

source: branches/2011/dev_r2787_LOCEAN3_TRA_TRP/NEMOGCM/NEMO/OPA_SRC/DYN/dynhpg.F90 @ 2789

Last change on this file since 2789 was 2789, checked in by cetlod, 13 years ago

Implementation of the merge of TRA/TRP : first guess, see ticket #842

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