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

source: branches/2011/dev_r2802_NOCL_prjhpg/NEMOGCM/NEMO/OPA_SRC/DYN/dynhpg.F90 @ 3080

Last change on this file since 3080 was 3080, checked in by hliu, 12 years ago

remove some used local variables from OPA_SRC/DYN/dynhpg.F90

  • Property svn:keywords set to Id
File size: 72.6 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   !!       hpg_prj  : s-coordinate (Pressure Jacobian with Cubic polynomial)
30   !!----------------------------------------------------------------------
31   USE oce             ! ocean dynamics and tracers
32   USE dom_oce         ! ocean space and time domain
33   USE phycst          ! physical constants
34   USE trdmod          ! ocean dynamics trends
35   USE trdmod_oce      ! ocean variables trends
36   USE in_out_manager  ! I/O manager
37   USE prtctl          ! Print control
38   USE lbclnk          ! lateral boundary condition
39   USE lib_mpp         ! MPP library
40
41   IMPLICIT NONE
42   PRIVATE
43
44   PUBLIC   dyn_hpg        ! routine called by step module
45   PUBLIC   dyn_hpg_init   ! routine called by opa module
46
47   !                                              !!* Namelist namdyn_hpg : hydrostatic pressure gradient
48   LOGICAL , PUBLIC ::   ln_hpg_zco    = .TRUE.    !: z-coordinate - full steps
49   LOGICAL , PUBLIC ::   ln_hpg_zps    = .FALSE.   !: z-coordinate - partial steps (interpolation)
50   LOGICAL , PUBLIC ::   ln_hpg_sco    = .FALSE.   !: s-coordinate (standard jacobian formulation)
51   LOGICAL , PUBLIC ::   ln_hpg_hel    = .FALSE.   !: s-coordinate (helsinki modification)
52   LOGICAL , PUBLIC ::   ln_hpg_wdj    = .FALSE.   !: s-coordinate (weighted density jacobian)
53   LOGICAL , PUBLIC ::   ln_hpg_djc    = .FALSE.   !: s-coordinate (Density Jacobian with Cubic polynomial)
54   LOGICAL , PUBLIC ::   ln_hpg_rot    = .FALSE.   !: s-coordinate (ROTated axes scheme)
55   LOGICAL , PUBLIC ::   ln_hpg_prj    = .FALSE.   !: s-coordinate (Pressure Jacobian scheme)
56   REAL(wp), PUBLIC ::   rn_gamma      = 0._wp     !: weighting coefficient
57   LOGICAL , PUBLIC ::   ln_dynhpg_imp = .FALSE.   !: semi-implicite hpg flag
58
59   INTEGER  ::   nhpg  =  0   ! = 0 to 7, type of pressure gradient scheme used ! (deduced from ln_hpg_... flags)
60
61   !! * Substitutions
62#  include "domzgr_substitute.h90"
63#  include "vectopt_loop_substitute.h90"
64   !!----------------------------------------------------------------------
65   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
66   !! $Id$
67   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
68   !!----------------------------------------------------------------------
69CONTAINS
70
71   SUBROUTINE dyn_hpg( kt )
72      !!---------------------------------------------------------------------
73      !!                  ***  ROUTINE dyn_hpg  ***
74      !!
75      !! ** Method  :   Call the hydrostatic pressure gradient routine
76      !!              using the scheme defined in the namelist
77      !!   
78      !! ** Action : - Update (ua,va) with the now hydrastatic pressure trend
79      !!             - Save the trend (l_trddyn=T)
80      !!----------------------------------------------------------------------
81      USE wrk_nemo, ONLY:   wrk_in_use, wrk_not_released
82      USE wrk_nemo, ONLY:   ztrdu => wrk_3d_1 , ztrdv => wrk_3d_2   ! 3D workspace
83      !!
84      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
85      !!----------------------------------------------------------------------
86      !
87      IF( wrk_in_use(3, 1,2) ) THEN
88         CALL ctl_stop('dyn_hpg: requested workspace arrays are unavailable')   ;   RETURN
89      ENDIF
90      !
91      IF( l_trddyn ) THEN                    ! Temporary saving of ua and va trends (l_trddyn)
92         ztrdu(:,:,:) = ua(:,:,:) 
93         ztrdv(:,:,:) = va(:,:,:) 
94      ENDIF     
95      !
96      SELECT CASE ( nhpg )      ! Hydrastatic pressure gradient computation
97      CASE (  0 )   ;   CALL hpg_zco    ( kt )      ! z-coordinate
98      CASE (  1 )   ;   CALL hpg_zps    ( kt )      ! z-coordinate plus partial steps (interpolation)
99      CASE (  2 )   ;   CALL hpg_sco    ( kt )      ! s-coordinate (standard jacobian formulation)
100      CASE (  3 )   ;   CALL hpg_hel    ( kt )      ! s-coordinate (helsinki modification)
101      CASE (  4 )   ;   CALL hpg_wdj    ( kt )      ! s-coordinate (weighted density jacobian)
102      CASE (  5 )   ;   CALL hpg_djc    ( kt )      ! s-coordinate (Density Jacobian with Cubic polynomial)
103      CASE (  6 )   ;   CALL hpg_rot    ( kt )      ! s-coordinate (ROTated axes scheme)
104      CASE (  7 )   ;   CALL hpg_prj    ( kt )      ! s-coordinate (Pressure Jacobian scheme)
105      END SELECT
106      !
107      IF( l_trddyn ) THEN      ! save the hydrostatic pressure gradient trends for momentum trend diagnostics
108         ztrdu(:,:,:) = ua(:,:,:) - ztrdu(:,:,:)
109         ztrdv(:,:,:) = va(:,:,:) - ztrdv(:,:,:)
110         CALL trd_mod( ztrdu, ztrdv, jpdyn_trd_hpg, 'DYN', kt )
111      ENDIF         
112      !
113      IF(ln_ctl)   CALL prt_ctl( tab3d_1=ua, clinfo1=' hpg  - Ua: ', mask1=umask,   &
114         &                       tab3d_2=va, clinfo2=       ' Va: ', mask2=vmask, clinfo3='dyn' )
115      !
116      IF( wrk_not_released(3, 1,2) )   CALL ctl_stop('dyn_hpg: failed to release workspace arrays')
117      !
118   END SUBROUTINE dyn_hpg
119
120
121   SUBROUTINE dyn_hpg_init
122      !!----------------------------------------------------------------------
123      !!                 ***  ROUTINE dyn_hpg_init  ***
124      !!
125      !! ** Purpose :   initializations for the hydrostatic pressure gradient
126      !!              computation and consistency control
127      !!
128      !! ** Action  :   Read the namelist namdyn_hpg and check the consistency
129      !!      with the type of vertical coordinate used (zco, zps, sco)
130      !!----------------------------------------------------------------------
131      INTEGER ::   ioptio = 0      ! temporary integer
132      !!
133      NAMELIST/namdyn_hpg/ ln_hpg_zco, ln_hpg_zps, ln_hpg_sco, ln_hpg_hel,    &
134         &                 ln_hpg_wdj, ln_hpg_djc, ln_hpg_rot, ln_hpg_prj, rn_gamma  , ln_dynhpg_imp
135      !!----------------------------------------------------------------------
136      !
137      REWIND( numnam )               ! Read Namelist namdyn_hpg
138      READ  ( numnam, namdyn_hpg )
139      !
140      IF(lwp) THEN                   ! Control print
141         WRITE(numout,*)
142         WRITE(numout,*) 'dyn_hpg_init : hydrostatic pressure gradient initialisation'
143         WRITE(numout,*) '~~~~~~~~~~~~'
144         WRITE(numout,*) '   Namelist namdyn_hpg : choice of hpg scheme'
145         WRITE(numout,*) '      z-coord. - full steps                             ln_hpg_zco    = ', ln_hpg_zco
146         WRITE(numout,*) '      z-coord. - partial steps (interpolation)          ln_hpg_zps    = ', ln_hpg_zps
147         WRITE(numout,*) '      s-coord. (standard jacobian formulation)          ln_hpg_sco    = ', ln_hpg_sco
148         WRITE(numout,*) '      s-coord. (helsinki modification)                  ln_hpg_hel    = ', ln_hpg_hel
149         WRITE(numout,*) '      s-coord. (weighted density jacobian)              ln_hpg_wdj    = ', ln_hpg_wdj
150         WRITE(numout,*) '      s-coord. (Density Jacobian: Cubic polynomial)     ln_hpg_djc    = ', ln_hpg_djc
151         WRITE(numout,*) '      s-coord. (ROTated axes scheme)                    ln_hpg_rot    = ', ln_hpg_rot
152         WRITE(numout,*) '      s-coord. (Pressure Jacobian: Cubic polynomial)    ln_hpg_prj    = ', ln_hpg_prj
153         WRITE(numout,*) '      weighting coeff. (wdj scheme)                     rn_gamma      = ', rn_gamma
154         WRITE(numout,*) '      time stepping: centered (F) or semi-implicit (T)  ln_dynhpg_imp = ', ln_dynhpg_imp
155      ENDIF
156      !
157      IF( lk_vvl .AND. .NOT. (ln_hpg_sco.OR.ln_hpg_prj) )   &
158         &   CALL ctl_stop('dyn_hpg_init : variable volume key_vvl require:&
159                           & the standard jacobian formulation hpg_sco or &
160                           & the pressure jacobian formulation hpg_prj')
161      !
162      !                               ! Set nhpg from ln_hpg_... flags
163      IF( ln_hpg_zco )   nhpg = 0
164      IF( ln_hpg_zps )   nhpg = 1
165      IF( ln_hpg_sco )   nhpg = 2
166      IF( ln_hpg_hel )   nhpg = 3
167      IF( ln_hpg_wdj )   nhpg = 4
168      IF( ln_hpg_djc )   nhpg = 5
169      IF( ln_hpg_rot )   nhpg = 6
170      IF( ln_hpg_prj )   nhpg = 7
171      !
172      !                               ! Consitency check
173      ioptio = 0 
174      IF( ln_hpg_zco )   ioptio = ioptio + 1
175      IF( ln_hpg_zps )   ioptio = ioptio + 1
176      IF( ln_hpg_sco )   ioptio = ioptio + 1
177      IF( ln_hpg_hel )   ioptio = ioptio + 1
178      IF( ln_hpg_wdj )   ioptio = ioptio + 1
179      IF( ln_hpg_djc )   ioptio = ioptio + 1
180      IF( ln_hpg_rot )   ioptio = ioptio + 1
181      IF( ln_hpg_prj )   ioptio = ioptio + 1
182      IF( ioptio /= 1 )   CALL ctl_stop( 'NO or several hydrostatic pressure gradient options used' )
183      !
184   END SUBROUTINE dyn_hpg_init
185
186
187   SUBROUTINE hpg_zco( kt )
188      !!---------------------------------------------------------------------
189      !!                  ***  ROUTINE hpg_zco  ***
190      !!
191      !! ** Method  :   z-coordinate case, levels are horizontal surfaces.
192      !!      The now hydrostatic pressure gradient at a given level, jk,
193      !!      is computed by taking the vertical integral of the in-situ
194      !!      density gradient along the model level from the suface to that
195      !!      level:    zhpi = grav .....
196      !!                zhpj = grav .....
197      !!      add it to the general momentum trend (ua,va).
198      !!            ua = ua - 1/e1u * zhpi
199      !!            va = va - 1/e2v * zhpj
200      !!
201      !! ** Action : - Update (ua,va) with the now hydrastatic pressure trend
202      !!----------------------------------------------------------------------
203      USE oce, ONLY:   zhpi => ta , zhpj => sa   ! (ta,sa) used as 3D workspace
204      !!
205      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
206      !!
207      INTEGER  ::   ji, jj, jk       ! dummy loop indices
208      REAL(wp) ::   zcoef0, zcoef1   ! temporary scalars
209      !!----------------------------------------------------------------------
210     
211      IF( kt == nit000 ) THEN
212         IF(lwp) WRITE(numout,*)
213         IF(lwp) WRITE(numout,*) 'dyn:hpg_zco : hydrostatic pressure gradient trend'
214         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   z-coordinate case '
215      ENDIF
216     
217      zcoef0 = - grav * 0.5_wp      ! Local constant initialization
218
219      ! Surface value
220      DO jj = 2, jpjm1
221         DO ji = fs_2, fs_jpim1   ! vector opt.
222            zcoef1 = zcoef0 * fse3w(ji,jj,1)
223            ! hydrostatic pressure gradient
224            zhpi(ji,jj,1) = zcoef1 * ( rhd(ji+1,jj,1) - rhd(ji,jj,1) ) / e1u(ji,jj)
225            zhpj(ji,jj,1) = zcoef1 * ( rhd(ji,jj+1,1) - rhd(ji,jj,1) ) / e2v(ji,jj)
226            ! add to the general momentum trend
227            ua(ji,jj,1) = ua(ji,jj,1) + zhpi(ji,jj,1)
228            va(ji,jj,1) = va(ji,jj,1) + zhpj(ji,jj,1)
229         END DO
230      END DO
231      !
232      ! interior value (2=<jk=<jpkm1)
233      DO jk = 2, jpkm1
234         DO jj = 2, jpjm1
235            DO ji = fs_2, fs_jpim1   ! vector opt.
236               zcoef1 = zcoef0 * fse3w(ji,jj,jk)
237               ! hydrostatic pressure gradient
238               zhpi(ji,jj,jk) = zhpi(ji,jj,jk-1)   &
239                  &           + zcoef1 * (  ( rhd(ji+1,jj,jk)+rhd(ji+1,jj,jk-1) )   &
240                  &                       - ( rhd(ji  ,jj,jk)+rhd(ji  ,jj,jk-1) )  ) / e1u(ji,jj)
241
242               zhpj(ji,jj,jk) = zhpj(ji,jj,jk-1)   &
243                  &           + zcoef1 * (  ( rhd(ji,jj+1,jk)+rhd(ji,jj+1,jk-1) )   &
244                  &                       - ( rhd(ji,jj,  jk)+rhd(ji,jj  ,jk-1) )  ) / e2v(ji,jj)
245               ! add to the general momentum trend
246               ua(ji,jj,jk) = ua(ji,jj,jk) + zhpi(ji,jj,jk)
247               va(ji,jj,jk) = va(ji,jj,jk) + zhpj(ji,jj,jk)
248            END DO
249         END DO
250      END DO
251      !
252   END SUBROUTINE hpg_zco
253
254
255   SUBROUTINE hpg_zps( kt )
256      !!---------------------------------------------------------------------
257      !!                 ***  ROUTINE hpg_zps  ***
258      !!                   
259      !! ** Method  :   z-coordinate plus partial steps case.  blahblah...
260      !!
261      !! ** Action  : - Update (ua,va) with the now hydrastatic pressure trend
262      !!----------------------------------------------------------------------
263      USE oce, ONLY:   zhpi => ta , zhpj => sa   ! (ta,sa) used as 3D workspace
264      !!
265      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
266      !!
267      INTEGER  ::   ji, jj, jk                       ! dummy loop indices
268      INTEGER  ::   iku, ikv                         ! temporary integers
269      REAL(wp) ::   zcoef0, zcoef1, zcoef2, zcoef3   ! temporary scalars
270      !!----------------------------------------------------------------------
271
272      IF( kt == nit000 ) THEN
273         IF(lwp) WRITE(numout,*)
274         IF(lwp) WRITE(numout,*) 'dyn:hpg_zps : hydrostatic pressure gradient trend'
275         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   z-coordinate with partial steps - vector optimization'
276      ENDIF
277
278      ! Local constant initialization
279      zcoef0 = - grav * 0.5_wp
280
281      !  Surface value (also valid in partial step case)
282      DO jj = 2, jpjm1
283         DO ji = fs_2, fs_jpim1   ! vector opt.
284            zcoef1 = zcoef0 * fse3w(ji,jj,1)
285            ! hydrostatic pressure gradient
286            zhpi(ji,jj,1) = zcoef1 * ( rhd(ji+1,jj  ,1) - rhd(ji,jj,1) ) / e1u(ji,jj)
287            zhpj(ji,jj,1) = zcoef1 * ( rhd(ji  ,jj+1,1) - rhd(ji,jj,1) ) / e2v(ji,jj)
288            ! add to the general momentum trend
289            ua(ji,jj,1) = ua(ji,jj,1) + zhpi(ji,jj,1)
290            va(ji,jj,1) = va(ji,jj,1) + zhpj(ji,jj,1)
291         END DO
292      END DO
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      ! partial steps correction at the last level  (use gru & grv computed in zpshde.F90)
315# if defined key_vectopt_loop
316         jj = 1
317         DO ji = jpi+2, jpij-jpi-1   ! vector opt. (forced unrolling)
318# else
319      DO jj = 2, jpjm1
320         DO ji = 2, jpim1
321# endif
322            iku = mbku(ji,jj)
323            ikv = mbkv(ji,jj)
324            zcoef2 = zcoef0 * MIN( fse3w(ji,jj,iku), fse3w(ji+1,jj  ,iku) )
325            zcoef3 = zcoef0 * MIN( fse3w(ji,jj,ikv), fse3w(ji  ,jj+1,ikv) )
326            IF( iku > 1 ) THEN            ! on i-direction (level 2 or more)
327               ua  (ji,jj,iku) = ua(ji,jj,iku) - zhpi(ji,jj,iku)         ! subtract old value
328               zhpi(ji,jj,iku) = zhpi(ji,jj,iku-1)                   &   ! compute the new one
329                  &            + zcoef2 * ( rhd(ji+1,jj,iku-1) - rhd(ji,jj,iku-1) + gru(ji,jj) ) / e1u(ji,jj)
330               ua  (ji,jj,iku) = ua(ji,jj,iku) + zhpi(ji,jj,iku)         ! add the new one to the general momentum trend
331            ENDIF
332            IF( ikv > 1 ) THEN            ! on j-direction (level 2 or more)
333               va  (ji,jj,ikv) = va(ji,jj,ikv) - zhpj(ji,jj,ikv)         ! subtract old value
334               zhpj(ji,jj,ikv) = zhpj(ji,jj,ikv-1)                   &   ! compute the new one
335                  &            + zcoef3 * ( rhd(ji,jj+1,ikv-1) - rhd(ji,jj,ikv-1) + grv(ji,jj) ) / e2v(ji,jj)
336               va  (ji,jj,ikv) = va(ji,jj,ikv) + zhpj(ji,jj,ikv)         ! add the new one to the general momentum trend
337            ENDIF
338# if ! defined key_vectopt_loop
339         END DO
340# endif
341      END DO
342      !
343   END SUBROUTINE hpg_zps
344
345
346   SUBROUTINE hpg_sco( kt )
347      !!---------------------------------------------------------------------
348      !!                  ***  ROUTINE hpg_sco  ***
349      !!
350      !! ** Method  :   s-coordinate case. Jacobian scheme.
351      !!      The now hydrostatic pressure gradient at a given level, jk,
352      !!      is computed by taking the vertical integral of the in-situ
353      !!      density gradient along the model level from the suface to that
354      !!      level. s-coordinates (ln_sco): a corrective term is added
355      !!      to the horizontal pressure gradient :
356      !!         zhpi = grav .....  + 1/e1u mi(rhd) di[ grav dep3w ]
357      !!         zhpj = grav .....  + 1/e2v mj(rhd) dj[ grav dep3w ]
358      !!      add it to the general momentum trend (ua,va).
359      !!         ua = ua - 1/e1u * zhpi
360      !!         va = va - 1/e2v * zhpj
361      !!
362      !! ** Action : - Update (ua,va) with the now hydrastatic pressure trend
363      !!----------------------------------------------------------------------
364      USE oce, ONLY:   zhpi => ta , zhpj => sa   ! (ta,sa) used as 3D workspace
365      !!
366      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
367      !!
368      INTEGER  ::   ji, jj, jk                 ! dummy loop indices
369      REAL(wp) ::   zcoef0, zuap, zvap, znad   ! temporary scalars
370      !!----------------------------------------------------------------------
371
372      IF( kt == nit000 ) THEN
373         IF(lwp) WRITE(numout,*)
374         IF(lwp) WRITE(numout,*) 'dyn:hpg_sco : hydrostatic pressure gradient trend'
375         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   s-coordinate case, OPA original scheme used'
376      ENDIF
377
378      ! Local constant initialization
379      zcoef0 = - grav * 0.5_wp
380      ! To use density and not density anomaly
381      IF ( lk_vvl ) THEN   ;     znad = 1._wp          ! Variable volume
382      ELSE                 ;     znad = 0._wp         ! Fixed volume
383      ENDIF
384
385      ! Surface value
386      DO jj = 2, jpjm1
387         DO ji = fs_2, fs_jpim1   ! vector opt.   
388            ! hydrostatic pressure gradient along s-surfaces
389            zhpi(ji,jj,1) = zcoef0 / e1u(ji,jj) * ( fse3w(ji+1,jj  ,1) * ( znad + rhd(ji+1,jj  ,1) )   &
390               &                                  - fse3w(ji  ,jj  ,1) * ( znad + rhd(ji  ,jj  ,1) ) )
391            zhpj(ji,jj,1) = zcoef0 / e2v(ji,jj) * ( fse3w(ji  ,jj+1,1) * ( znad + rhd(ji  ,jj+1,1) )   &
392               &                                  - fse3w(ji  ,jj  ,1) * ( znad + rhd(ji  ,jj  ,1) ) )
393            ! s-coordinate pressure gradient correction
394            zuap = -zcoef0 * ( rhd   (ji+1,jj,1) + rhd   (ji,jj,1) + 2._wp * znad )   &
395               &           * ( fsde3w(ji+1,jj,1) - fsde3w(ji,jj,1) ) / e1u(ji,jj)
396            zvap = -zcoef0 * ( rhd   (ji,jj+1,1) + rhd   (ji,jj,1) + 2._wp * znad )   &
397               &           * ( fsde3w(ji,jj+1,1) - fsde3w(ji,jj,1) ) / e2v(ji,jj)
398            ! add to the general momentum trend
399            ua(ji,jj,1) = ua(ji,jj,1) + zhpi(ji,jj,1) + zuap
400            va(ji,jj,1) = va(ji,jj,1) + zhpj(ji,jj,1) + zvap
401         END DO 
402      END DO   
403           
404      ! interior value (2=<jk=<jpkm1)
405      DO jk = 2, jpkm1                                 
406         DO jj = 2, jpjm1     
407            DO ji = fs_2, fs_jpim1   ! vector opt.     
408               ! hydrostatic pressure gradient along s-surfaces
409               zhpi(ji,jj,jk) = zhpi(ji,jj,jk-1) + zcoef0 / e1u(ji,jj)   & 
410                  &           * (  fse3w(ji+1,jj,jk) * ( rhd(ji+1,jj,jk) + rhd(ji+1,jj,jk-1) + 2*znad )   & 
411                  &              - fse3w(ji  ,jj,jk) * ( rhd(ji  ,jj,jk) + rhd(ji  ,jj,jk-1) + 2*znad )  )
412               zhpj(ji,jj,jk) = zhpj(ji,jj,jk-1) + zcoef0 / e2v(ji,jj)   &
413                  &           * (  fse3w(ji,jj+1,jk) * ( rhd(ji,jj+1,jk) + rhd(ji,jj+1,jk-1) + 2*znad )   &
414                  &              - fse3w(ji,jj  ,jk) * ( rhd(ji,jj,  jk) + rhd(ji,jj  ,jk-1) + 2*znad )  )
415               ! s-coordinate pressure gradient correction
416               zuap = -zcoef0 * ( rhd   (ji+1,jj  ,jk) + rhd   (ji,jj,jk) + 2._wp * znad )   &
417                  &           * ( fsde3w(ji+1,jj  ,jk) - fsde3w(ji,jj,jk) ) / e1u(ji,jj)
418               zvap = -zcoef0 * ( rhd   (ji  ,jj+1,jk) + rhd   (ji,jj,jk) + 2._wp * znad )   &
419                  &           * ( fsde3w(ji  ,jj+1,jk) - fsde3w(ji,jj,jk) ) / e2v(ji,jj)
420               ! add to the general momentum trend
421               ua(ji,jj,jk) = ua(ji,jj,jk) + zhpi(ji,jj,jk) + zuap
422               va(ji,jj,jk) = va(ji,jj,jk) + zhpj(ji,jj,jk) + zvap
423            END DO
424         END DO
425      END DO
426      !
427   END SUBROUTINE hpg_sco
428
429
430   SUBROUTINE hpg_hel( kt )
431      !!---------------------------------------------------------------------
432      !!                  ***  ROUTINE hpg_hel  ***
433      !!
434      !! ** Method  :   s-coordinate case.
435      !!      The now hydrostatic pressure gradient at a given level
436      !!      jk is computed by taking the vertical integral of the in-situ
437      !!      density gradient along the model level from the suface to that
438      !!      level. s-coordinates (ln_sco): a corrective term is added
439      !!      to the horizontal pressure gradient :
440      !!         zhpi = grav .....  + 1/e1u mi(rhd) di[ grav dep3w ]
441      !!         zhpj = grav .....  + 1/e2v mj(rhd) dj[ grav dep3w ]
442      !!      add it to the general momentum trend (ua,va).
443      !!         ua = ua - 1/e1u * zhpi
444      !!         va = va - 1/e2v * zhpj
445      !!
446      !! ** Action : - Update (ua,va) with the now hydrastatic pressure trend
447      !!             - Save the trend (l_trddyn=T)
448      !!----------------------------------------------------------------------
449      USE oce, ONLY:   zhpi => ta , zhpj => sa   ! (ta,sa) used as 3D workspace
450      !!
451      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
452      !!
453      INTEGER  ::   ji, jj, jk           ! dummy loop indices
454      REAL(wp) ::   zcoef0, zuap, zvap   ! temporary scalars
455      !!----------------------------------------------------------------------
456
457      IF( kt == nit000 ) THEN
458         IF(lwp) WRITE(numout,*)
459         IF(lwp) WRITE(numout,*) 'dyn:hpg_hel : hydrostatic pressure gradient trend'
460         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   s-coordinate case, helsinki modified scheme'
461      ENDIF
462
463      ! Local constant initialization
464      zcoef0 = - grav * 0.5_wp
465 
466      ! Surface value
467      DO jj = 2, jpjm1
468         DO ji = fs_2, fs_jpim1   ! vector opt.
469            ! hydrostatic pressure gradient along s-surfaces
470            zhpi(ji,jj,1) = zcoef0 / e1u(ji,jj) * ( fse3t(ji+1,jj  ,1) * rhd(ji+1,jj  ,1)  &
471               &                                  - fse3t(ji  ,jj  ,1) * rhd(ji  ,jj  ,1) )
472            zhpj(ji,jj,1) = zcoef0 / e2v(ji,jj) * ( fse3t(ji  ,jj+1,1) * rhd(ji  ,jj+1,1)  &
473               &                                  - fse3t(ji  ,jj  ,1) * rhd(ji  ,jj  ,1) )
474            ! s-coordinate pressure gradient correction
475            zuap = -zcoef0 * ( rhd   (ji+1,jj,1) + rhd   (ji,jj,1) )   &
476               &           * ( fsdept(ji+1,jj,1) - fsdept(ji,jj,1) ) / e1u(ji,jj)
477            zvap = -zcoef0 * ( rhd   (ji,jj+1,1) + rhd   (ji,jj,1) )   &
478               &           * ( fsdept(ji,jj+1,1) - fsdept(ji,jj,1) ) / e2v(ji,jj)
479            ! add to the general momentum trend
480            ua(ji,jj,1) = ua(ji,jj,1) + zhpi(ji,jj,1) + zuap
481            va(ji,jj,1) = va(ji,jj,1) + zhpj(ji,jj,1) + zvap
482         END DO
483      END DO
484      !
485      ! interior value (2=<jk=<jpkm1)
486      DO jk = 2, jpkm1
487         DO jj = 2, jpjm1
488            DO ji = fs_2, fs_jpim1   ! vector opt.
489               ! hydrostatic pressure gradient along s-surfaces
490               zhpi(ji,jj,jk) = zhpi(ji,jj,jk-1) &
491                  &           +  zcoef0 / e1u(ji,jj) * ( fse3t(ji+1,jj,jk  ) * rhd(ji+1,jj,jk)     &
492                  &                                     -fse3t(ji  ,jj,jk  ) * rhd(ji  ,jj,jk)   ) &
493                  &           +  zcoef0 / e1u(ji,jj) * ( fse3t(ji+1,jj,jk-1) * rhd(ji+1,jj,jk-1)   &
494                  &                                     -fse3t(ji  ,jj,jk-1) * rhd(ji  ,jj,jk-1) )
495               zhpj(ji,jj,jk) = zhpj(ji,jj,jk-1) &
496                  &           +  zcoef0 / e2v(ji,jj) * ( fse3t(ji,jj+1,jk  ) * rhd(ji,jj+1,jk)     &
497                  &                                     -fse3t(ji,jj  ,jk  ) * rhd(ji,jj,  jk)   ) &
498                  &           +  zcoef0 / e2v(ji,jj) * ( fse3t(ji,jj+1,jk-1) * rhd(ji,jj+1,jk-1)   &
499                  &                                     -fse3t(ji,jj  ,jk-1) * rhd(ji,jj,  jk-1) )
500               ! s-coordinate pressure gradient correction
501               zuap = - zcoef0 * ( rhd   (ji+1,jj,jk) + rhd   (ji,jj,jk) )   &
502                  &            * ( fsdept(ji+1,jj,jk) - fsdept(ji,jj,jk) ) / e1u(ji,jj)
503               zvap = - zcoef0 * ( rhd   (ji,jj+1,jk) + rhd   (ji,jj,jk) )   &
504                  &            * ( fsdept(ji,jj+1,jk) - fsdept(ji,jj,jk) ) / e2v(ji,jj)
505               ! add to the general momentum trend
506               ua(ji,jj,jk) = ua(ji,jj,jk) + zhpi(ji,jj,jk) + zuap
507               va(ji,jj,jk) = va(ji,jj,jk) + zhpj(ji,jj,jk) + zvap
508            END DO
509         END DO
510      END DO
511      !
512   END SUBROUTINE hpg_hel
513
514
515   SUBROUTINE hpg_wdj( kt )
516      !!---------------------------------------------------------------------
517      !!                  ***  ROUTINE hpg_wdj  ***
518      !!
519      !! ** Method  :   Weighted Density Jacobian (wdj) scheme (song 1998)
520      !!      The weighting coefficients from the namelist parameter rn_gamma
521      !!      (alpha=0.5-rn_gamma ; beta=1-alpha=0.5+rn_gamma
522      !!
523      !! Reference : Song, Mon. Wea. Rev., 126, 3213-3230, 1998.
524      !!----------------------------------------------------------------------
525      USE oce, ONLY:   zhpi => ta , zhpj => sa   ! (ta,sa) used as 3D workspace
526      !!
527      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
528      !!
529      INTEGER  ::   ji, jj, jk           ! dummy loop indices
530      REAL(wp) ::   zcoef0, zuap, zvap   ! temporary scalars
531      REAL(wp) ::   zalph , zbeta        !    "         "
532      !!----------------------------------------------------------------------
533
534      IF( kt == nit000 ) THEN
535         IF(lwp) WRITE(numout,*)
536         IF(lwp) WRITE(numout,*) 'dyn:hpg_wdj : hydrostatic pressure gradient trend'
537         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   Weighted Density Jacobian'
538      ENDIF
539
540      ! Local constant initialization
541      zcoef0 = - grav * 0.5_wp
542      zalph  = 0.5_wp - rn_gamma    ! weighting coefficients (alpha=0.5-rn_gamma
543      zbeta  = 0.5_wp + rn_gamma    !                        (beta =1-alpha=0.5+rn_gamma
544
545      ! Surface value (no ponderation)
546      DO jj = 2, jpjm1
547         DO ji = fs_2, fs_jpim1   ! vector opt.
548            ! hydrostatic pressure gradient along s-surfaces
549            zhpi(ji,jj,1) = zcoef0 / e1u(ji,jj) * (  fse3w(ji+1,jj  ,1) * rhd(ji+1,jj  ,1)   &
550               &                                   - fse3w(ji  ,jj  ,1) * rhd(ji  ,jj  ,1)  )
551            zhpj(ji,jj,1) = zcoef0 / e2v(ji,jj) * (  fse3w(ji  ,jj+1,1) * rhd(ji  ,jj+1,1)   &
552               &                                   - fse3w(ji  ,jj  ,1) * rhd(ji,  jj  ,1)  )
553            ! s-coordinate pressure gradient correction
554            zuap = -zcoef0 * ( rhd   (ji+1,jj,1) + rhd   (ji,jj,1) )   &
555               &           * ( fsde3w(ji+1,jj,1) - fsde3w(ji,jj,1) ) / e1u(ji,jj)
556            zvap = -zcoef0 * ( rhd   (ji,jj+1,1) + rhd   (ji,jj,1) )   &
557               &           * ( fsde3w(ji,jj+1,1) - fsde3w(ji,jj,1) ) / e2v(ji,jj)
558            ! add to the general momentum trend
559            ua(ji,jj,1) = ua(ji,jj,1) + zhpi(ji,jj,1) + zuap
560            va(ji,jj,1) = va(ji,jj,1) + zhpj(ji,jj,1) + zvap
561         END DO
562      END DO
563
564      ! Interior value (2=<jk=<jpkm1) (weighted with zalph & zbeta)
565      DO jk = 2, jpkm1
566         DO jj = 2, jpjm1
567            DO ji = fs_2, fs_jpim1   ! vector opt.
568               zhpi(ji,jj,jk) = zhpi(ji,jj,jk-1) + zcoef0 / e1u(ji,jj)                            &
569                  &           * (   (            fsde3w(ji+1,jj,jk  ) + fsde3w(ji,jj,jk  )        &
570                  &                            - fsde3w(ji+1,jj,jk-1) - fsde3w(ji,jj,jk-1)    )   &
571                  &               * (  zalph * ( rhd   (ji+1,jj,jk-1) - rhd   (ji,jj,jk-1) )      &
572                  &                  + zbeta * ( rhd   (ji+1,jj,jk  ) - rhd   (ji,jj,jk  ) )  )   &
573                  &             -   (            rhd   (ji+1,jj,jk  ) + rhd   (ji,jj,jk  )        &
574                  &                           - rhd   (ji+1,jj,jk-1) - rhd   (ji,jj,jk-1)     )   &
575                  &               * (  zalph * ( fsde3w(ji+1,jj,jk-1) - fsde3w(ji,jj,jk-1) )      &
576                  &                  + zbeta * ( fsde3w(ji+1,jj,jk  ) - fsde3w(ji,jj,jk  ) )  )  )
577               zhpj(ji,jj,jk) = zhpj(ji,jj,jk-1) + zcoef0 / e2v(ji,jj)                            &
578                  &           * (   (           fsde3w(ji,jj+1,jk  ) + fsde3w(ji,jj,jk  )         &
579                  &                           - fsde3w(ji,jj+1,jk-1) - fsde3w(ji,jj,jk-1)     )   &
580                  &               * (  zalph * ( rhd   (ji,jj+1,jk-1) - rhd   (ji,jj,jk-1) )      &
581                  &                  + zbeta * ( rhd   (ji,jj+1,jk  ) - rhd   (ji,jj,jk  ) )  )   &
582                  &             -   (            rhd   (ji,jj+1,jk  ) + rhd   (ji,jj,jk  )        &
583                  &                            - rhd   (ji,jj+1,jk-1) - rhd   (ji,jj,jk-1)    )   &
584                  &               * (  zalph * ( fsde3w(ji,jj+1,jk-1) - fsde3w(ji,jj,jk-1) )      &
585                  &                  + zbeta * ( fsde3w(ji,jj+1,jk  ) - fsde3w(ji,jj,jk  ) )  )  )
586               ! add to the general momentum trend
587               ua(ji,jj,jk) = ua(ji,jj,jk) + zhpi(ji,jj,jk)
588               va(ji,jj,jk) = va(ji,jj,jk) + zhpj(ji,jj,jk)
589            END DO
590         END DO
591      END DO
592      !
593   END SUBROUTINE hpg_wdj
594
595
596   SUBROUTINE hpg_djc( kt )
597      !!---------------------------------------------------------------------
598      !!                  ***  ROUTINE hpg_djc  ***
599      !!
600      !! ** Method  :   Density Jacobian with Cubic polynomial scheme
601      !!
602      !! Reference: Shchepetkin and McWilliams, J. Geophys. Res., 108(C3), 3090, 2003
603      !!----------------------------------------------------------------------
604      USE wrk_nemo, ONLY:   wrk_in_use, wrk_not_released
605      USE oce     , ONLY:   zhpi  => ta        , zhpj => sa       ! (ta,sa) used as 3D workspace
606      USE wrk_nemo, ONLY:   drhox => wrk_3d_1  , dzx  => wrk_3d_2
607      USE wrk_nemo, ONLY:   drhou => wrk_3d_3  , dzu  => wrk_3d_4 , rho_i => wrk_3d_5
608      USE wrk_nemo, ONLY:   drhoy => wrk_3d_6  , dzy  => wrk_3d_7
609      USE wrk_nemo, ONLY:   drhov => wrk_3d_8  , dzv  => wrk_3d_9 , rho_j => wrk_3d_10
610      USE wrk_nemo, ONLY:   drhoz => wrk_3d_11 , dzz  => wrk_3d_12 
611      USE wrk_nemo, ONLY:   drhow => wrk_3d_13 , dzw  => wrk_3d_14
612      USE wrk_nemo, ONLY:   rho_k => wrk_3d_15
613      !!
614      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
615      !!
616      INTEGER  ::   ji, jj, jk          ! dummy loop indices
617      REAL(wp) ::   zcoef0, zep, cffw   ! temporary scalars
618      REAL(wp) ::   z1_10, cffu, cffx   !    "         "
619      REAL(wp) ::   z1_12, cffv, cffy   !    "         "
620      !!----------------------------------------------------------------------
621
622      IF( wrk_in_use(3, 1,2,3,4,5,6,7,8,9,10,11,12,13,14,15) ) THEN
623         CALL ctl_stop('dyn:hpg_djc: requested workspace arrays unavailable')   ;   RETURN
624      ENDIF
625
626      IF( kt == nit000 ) THEN
627         IF(lwp) WRITE(numout,*)
628         IF(lwp) WRITE(numout,*) 'dyn:hpg_djc : hydrostatic pressure gradient trend'
629         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   s-coordinate case, density Jacobian with cubic polynomial scheme'
630      ENDIF
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 wrk_nemo, ONLY:   wrk_in_use, wrk_not_released
836      USE oce     , ONLY:   zhpi    => ta       , zhpj    => sa       ! (ta,sa) used as 3D workspace
837      USE wrk_nemo, ONLY:   zdistr  => wrk_2d_1 , zsina   => wrk_2d_2 , zcosa  => wrk_2d_3
838      USE wrk_nemo, ONLY:   zhpiorg => wrk_3d_1 , zhpirot => wrk_3d_2
839      USE wrk_nemo, ONLY:   zhpitra => wrk_3d_3 , zhpine  => wrk_3d_4
840      USE wrk_nemo, ONLY:   zhpjorg => wrk_3d_5 , zhpjrot => wrk_3d_6
841      USE wrk_nemo, ONLY:   zhpjtra => wrk_3d_7 , zhpjne  => wrk_3d_8
842      !!
843      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
844      !!
845      INTEGER  ::   ji, jj, jk          ! dummy loop indices
846      REAL(wp) ::   zforg, zcoef0, zuap, zmskd1, zmskd1m   ! temporary scalar
847      REAL(wp) ::   zfrot        , zvap, zmskd2, zmskd2m   !    "         "
848      !!----------------------------------------------------------------------
849
850      IF( wrk_in_use(2, 1,2,3)             .OR.   &
851          wrk_in_use(3, 1,2,3,4,5,6,7,8) ) THEN
852         CALL ctl_stop('dyn:hpg_rot: requested workspace arrays unavailable')   ;   RETURN
853      ENDIF
854
855      IF( kt == nit000 ) THEN
856         IF(lwp) WRITE(numout,*)
857         IF(lwp) WRITE(numout,*) 'dyn:hpg_rot : hydrostatic pressure gradient trend'
858         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   s-coordinate case, rotated axes scheme used'
859      ENDIF
860
861      ! -------------------------------
862      !  Local constant initialization
863      ! -------------------------------
864      zcoef0 = - grav * 0.5_wp
865      zforg  = 0.95_wp
866      zfrot  = 1._wp - zforg
867
868      ! inverse of the distance between 2 diagonal T-points (defined at F-point) (here zcoef0/distance)
869      zdistr(:,:) = zcoef0 / SQRT( e1f(:,:)*e1f(:,:) + e2f(:,:)*e1f(:,:) )
870
871      ! sinus and cosinus of diagonal angle at F-point
872      zsina(:,:) = ATAN2( e2f(:,:), e1f(:,:) )
873      zcosa(:,:) = COS( zsina(:,:) )
874      zsina(:,:) = SIN( zsina(:,:) )
875
876      ! ---------------
877      !  Surface value
878      ! ---------------
879      ! compute and add to the general trend the pressure gradients along the axes
880      DO jj = 2, jpjm1
881         DO ji = fs_2, fs_jpim1   ! vector opt.
882            ! hydrostatic pressure gradient along s-surfaces
883            zhpiorg(ji,jj,1) = zcoef0 / e1u(ji,jj) * (  fse3t(ji+1,jj,1) * rhd(ji+1,jj,1)   &
884               &                                      - fse3t(ji  ,jj,1) * rhd(ji  ,jj,1)  )
885            zhpjorg(ji,jj,1) = zcoef0 / e2v(ji,jj) * (  fse3t(ji,jj+1,1) * rhd(ji,jj+1,1)   &
886               &                                      - fse3t(ji,jj  ,1) * rhd(ji,jj  ,1)  )
887            ! s-coordinate pressure gradient correction
888            zuap = -zcoef0 * ( rhd   (ji+1,jj  ,1) + rhd   (ji,jj,1) )   &
889               &           * ( fsdept(ji+1,jj  ,1) - fsdept(ji,jj,1) ) / e1u(ji,jj)
890            zvap = -zcoef0 * ( rhd   (ji  ,jj+1,1) + rhd   (ji,jj,1) )   &
891               &           * ( fsdept(ji  ,jj+1,1) - fsdept(ji,jj,1) ) / e2v(ji,jj)
892            ! add to the general momentum trend
893            ua(ji,jj,1) = ua(ji,jj,1) + zforg * ( zhpiorg(ji,jj,1) + zuap )
894            va(ji,jj,1) = va(ji,jj,1) + zforg * ( zhpjorg(ji,jj,1) + zvap )
895         END DO
896      END DO
897
898      ! compute the pressure gradients in the diagonal directions
899      DO jj = 1, jpjm1
900         DO ji = 1, fs_jpim1   ! vector opt.
901            zmskd1 = tmask(ji+1,jj+1,1) * tmask(ji  ,jj,1)      ! mask in the 1st diagnonal
902            zmskd2 = tmask(ji  ,jj+1,1) * tmask(ji+1,jj,1)      ! mask in the 2nd diagnonal
903            ! hydrostatic pressure gradient along s-surfaces
904            zhpitra(ji,jj,1) = zdistr(ji,jj) * zmskd1 * (  fse3t(ji+1,jj+1,1) * rhd(ji+1,jj+1,1)   &
905               &                                         - fse3t(ji  ,jj  ,1) * rhd(ji  ,jj  ,1)  )
906            zhpjtra(ji,jj,1) = zdistr(ji,jj) * zmskd2 * (  fse3t(ji  ,jj+1,1) * rhd(ji  ,jj+1,1)   &
907               &                                         - fse3t(ji+1,jj  ,1) * rhd(ji+1,jj  ,1)  )
908            ! s-coordinate pressure gradient correction
909            zuap = -zdistr(ji,jj) * zmskd1 * ( rhd   (ji+1,jj+1,1) + rhd   (ji  ,jj,1) )   &
910               &                           * ( fsdept(ji+1,jj+1,1) - fsdept(ji  ,jj,1) )
911            zvap = -zdistr(ji,jj) * zmskd2 * ( rhd   (ji  ,jj+1,1) + rhd   (ji+1,jj,1) )   &
912               &                           * ( fsdept(ji  ,jj+1,1) - fsdept(ji+1,jj,1) )
913            ! back rotation
914            zhpine(ji,jj,1) = zcosa(ji,jj) * ( zhpitra(ji,jj,1) + zuap )   &
915               &            - zsina(ji,jj) * ( zhpjtra(ji,jj,1) + zvap )
916            zhpjne(ji,jj,1) = zsina(ji,jj) * ( zhpitra(ji,jj,1) + zuap )   &
917               &            + zcosa(ji,jj) * ( zhpjtra(ji,jj,1) + zvap )
918         END DO
919      END DO
920
921      ! interpolate and add to the general trend the diagonal gradient
922      DO jj = 2, jpjm1
923         DO ji = fs_2, fs_jpim1   ! vector opt.
924            ! averaging
925            zhpirot(ji,jj,1) = 0.5 * ( zhpine(ji,jj,1) + zhpine(ji  ,jj-1,1) )
926            zhpjrot(ji,jj,1) = 0.5 * ( zhpjne(ji,jj,1) + zhpjne(ji-1,jj  ,1) )
927            ! add to the general momentum trend
928            ua(ji,jj,1) = ua(ji,jj,1) + zfrot * zhpirot(ji,jj,1) 
929            va(ji,jj,1) = va(ji,jj,1) + zfrot * zhpjrot(ji,jj,1) 
930         END DO
931      END DO
932
933      ! -----------------
934      ! 2. interior value (2=<jk=<jpkm1)
935      ! -----------------
936      ! compute and add to the general trend the pressure gradients along the axes
937      DO jk = 2, jpkm1
938         DO jj = 2, jpjm1
939            DO ji = fs_2, fs_jpim1   ! vector opt.
940               ! hydrostatic pressure gradient along s-surfaces
941               zhpiorg(ji,jj,jk) = zhpiorg(ji,jj,jk-1)                                                 &
942                  &              +  zcoef0 / e1u(ji,jj) * (  fse3t(ji+1,jj,jk  ) * rhd(ji+1,jj,jk  )   &
943                  &                                        - fse3t(ji  ,jj,jk  ) * rhd(ji  ,jj,jk  )   &
944                  &                                        + fse3t(ji+1,jj,jk-1) * rhd(ji+1,jj,jk-1)   &
945                  &                                        - fse3t(ji  ,jj,jk-1) * rhd(ji  ,jj,jk-1)  )
946               zhpjorg(ji,jj,jk) = zhpjorg(ji,jj,jk-1)                                                 &
947                  &              +  zcoef0 / e2v(ji,jj) * (  fse3t(ji,jj+1,jk  ) * rhd(ji,jj+1,jk  )   &
948                  &                                        - fse3t(ji,jj  ,jk  ) * rhd(ji,jj,  jk  )   &
949                  &                                        + fse3t(ji,jj+1,jk-1) * rhd(ji,jj+1,jk-1)   &
950                  &                                        - fse3t(ji,jj  ,jk-1) * rhd(ji,jj,  jk-1)  )
951               ! s-coordinate pressure gradient correction
952               zuap = - zcoef0 * ( rhd   (ji+1,jj  ,jk) + rhd   (ji,jj,jk) )   &
953                  &            * ( fsdept(ji+1,jj  ,jk) - fsdept(ji,jj,jk) ) / e1u(ji,jj)
954               zvap = - zcoef0 * ( rhd   (ji  ,jj+1,jk) + rhd   (ji,jj,jk) )   &
955                  &            * ( fsdept(ji  ,jj+1,jk) - fsdept(ji,jj,jk) ) / e2v(ji,jj)
956               ! add to the general momentum trend
957               ua(ji,jj,jk) = ua(ji,jj,jk) + zforg*( zhpiorg(ji,jj,jk) + zuap )
958               va(ji,jj,jk) = va(ji,jj,jk) + zforg*( zhpjorg(ji,jj,jk) + zvap )
959            END DO
960         END DO
961      END DO
962
963      ! compute the pressure gradients in the diagonal directions
964      DO jk = 2, jpkm1
965         DO jj = 1, jpjm1
966            DO ji = 1, fs_jpim1   ! vector opt.
967               zmskd1  = tmask(ji+1,jj+1,jk  ) * tmask(ji  ,jj,jk  )      ! level jk   mask in the 1st diagnonal
968               zmskd1m = tmask(ji+1,jj+1,jk-1) * tmask(ji  ,jj,jk-1)      ! level jk-1    "               "     
969               zmskd2  = tmask(ji  ,jj+1,jk  ) * tmask(ji+1,jj,jk  )      ! level jk   mask in the 2nd diagnonal
970               zmskd2m = tmask(ji  ,jj+1,jk-1) * tmask(ji+1,jj,jk-1)      ! level jk-1    "               "     
971               ! hydrostatic pressure gradient along s-surfaces
972               zhpitra(ji,jj,jk) = zhpitra(ji,jj,jk-1)                                                       &
973                  &              + zdistr(ji,jj) * zmskd1  * ( fse3t(ji+1,jj+1,jk  ) * rhd(ji+1,jj+1,jk)     &
974                  &                                           -fse3t(ji  ,jj  ,jk  ) * rhd(ji  ,jj  ,jk) )   &
975                  &              + zdistr(ji,jj) * zmskd1m * ( fse3t(ji+1,jj+1,jk-1) * rhd(ji+1,jj+1,jk-1)   &
976                  &                                           -fse3t(ji  ,jj  ,jk-1) * rhd(ji  ,jj  ,jk-1) )
977               zhpjtra(ji,jj,jk) = zhpjtra(ji,jj,jk-1)                                                       &
978                  &              + zdistr(ji,jj) * zmskd2  * ( fse3t(ji  ,jj+1,jk  ) * rhd(ji  ,jj+1,jk)     &
979                  &                                           -fse3t(ji+1,jj  ,jk  ) * rhd(ji+1,jj,  jk) )   &
980                  &              + zdistr(ji,jj) * zmskd2m * ( fse3t(ji  ,jj+1,jk-1) * rhd(ji  ,jj+1,jk-1)   &
981                  &                                           -fse3t(ji+1,jj  ,jk-1) * rhd(ji+1,jj,  jk-1) )
982               ! s-coordinate pressure gradient correction
983               zuap = - zdistr(ji,jj) * zmskd1 * ( rhd   (ji+1,jj+1,jk) + rhd   (ji  ,jj,jk) )   &
984                  &                            * ( fsdept(ji+1,jj+1,jk) - fsdept(ji  ,jj,jk) )
985               zvap = - zdistr(ji,jj) * zmskd2 * ( rhd   (ji  ,jj+1,jk) + rhd   (ji+1,jj,jk) )   &
986                  &                            * ( fsdept(ji  ,jj+1,jk) - fsdept(ji+1,jj,jk) )
987               ! back rotation
988               zhpine(ji,jj,jk) = zcosa(ji,jj) * ( zhpitra(ji,jj,jk) + zuap )   &
989                  &             - zsina(ji,jj) * ( zhpjtra(ji,jj,jk) + zvap )
990               zhpjne(ji,jj,jk) = zsina(ji,jj) * ( zhpitra(ji,jj,jk) + zuap )   &
991                  &             + zcosa(ji,jj) * ( zhpjtra(ji,jj,jk) + zvap )
992            END DO
993         END DO
994      END DO
995
996      ! interpolate and add to the general trend
997      DO jk = 2, jpkm1
998         DO jj = 2, jpjm1
999            DO ji = fs_2, fs_jpim1   ! vector opt.
1000               ! averaging
1001               zhpirot(ji,jj,jk) = 0.5 * ( zhpine(ji,jj,jk) + zhpine(ji  ,jj-1,jk) )
1002               zhpjrot(ji,jj,jk) = 0.5 * ( zhpjne(ji,jj,jk) + zhpjne(ji-1,jj  ,jk) )
1003               ! add to the general momentum trend
1004               ua(ji,jj,jk) = ua(ji,jj,jk) + zfrot * zhpirot(ji,jj,jk) 
1005               va(ji,jj,jk) = va(ji,jj,jk) + zfrot * zhpjrot(ji,jj,jk) 
1006            END DO
1007         END DO
1008      END DO
1009      !
1010      IF( wrk_not_released(2, 1,2,3)           .OR.   &
1011          wrk_not_released(3, 1,2,3,4,5,6,7,8) )   CALL ctl_stop('dyn:hpg_rot: failed to release workspace arrays')
1012      !
1013   END SUBROUTINE hpg_rot
1014
1015   
1016   SUBROUTINE hpg_prj( kt )
1017      !!---------------------------------------------------------------------
1018      !!                  ***  ROUTINE hpg_prj  ***
1019      !!
1020      !! ** Method  :   s-coordinate case.
1021      !!      Reformulate the horizontal hydrostatical pressure gradient
1022      !!      term using Pressure Jacobian. A new correction term
1023      !!      is developed to eliminate the sigma-coordinate error.
1024      !!
1025      !! ** Action : - Update (ua,va) with the now hydrastatic pressure trend
1026      !!             - Save the trend (l_trddyn=T)
1027      !!
1028      !!----------------------------------------------------------------------
1029
1030      USE wrk_nemo, ONLY:   wrk_in_use, wrk_not_released
1031      USE oce     , ONLY:   fsde3w_tmp => ta       ! (ta,sa) used as 3D workspace
1032      USE oce     , ONLY:      rhd_tmp => sa 
1033      USE wrk_nemo, ONLY:         zhpi => wrk_3d_3 
1034      USE wrk_nemo, ONLY:           zu => wrk_3d_4 
1035      USE wrk_nemo, ONLY:           zv => wrk_3d_5
1036      USE wrk_nemo, ONLY:          fsp => wrk_3d_6 
1037      USE wrk_nemo, ONLY:          xsp => wrk_3d_7
1038      USE wrk_nemo, ONLY:          asp => wrk_3d_8
1039      USE wrk_nemo, ONLY:          bsp => wrk_3d_9
1040      USE wrk_nemo, ONLY:          csp => wrk_3d_10
1041      USE wrk_nemo, ONLY:          dsp => wrk_3d_11
1042      !!
1043      !!----------------------------------------------------------------------
1044      !!
1045      INTEGER, PARAMETER  :: polynomial_type = 1    ! 1: cubic spline, 2: linear
1046      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
1047      !!
1048      INTEGER  ::   ji, jj, jk, jkk                ! dummy loop indices
1049      REAL(wp) ::   zcoef0, znad               ! temporary scalars
1050      !!
1051      !! The local varialbes for the correction term
1052      INTEGER  :: jk1, jis, jid, jjs, jjd
1053      REAL(wp) :: zuijk, zvijk, pwes, pwed, pnss, pnsd, deps
1054      REAL(wp) :: rhdt1 
1055      REAL(wp) :: dpdx1, dpdx2, dpdy1, dpdy2
1056      INTEGER  :: bhitwe, bhitns
1057      !!----------------------------------------------------------------------
1058
1059      IF( wrk_in_use(3, 3,4,5,6,7,8,9,10,11) ) THEN
1060         CALL ctl_stop('dyn:hpg_prj: requested workspace arrays unavailable')   ;   RETURN
1061      ENDIF
1062
1063      IF( kt == nit000 ) THEN
1064         IF(lwp) WRITE(numout,*)
1065         IF(lwp) WRITE(numout,*) 'dyn:hpg_prj : hydrostatic pressure gradient trend'
1066         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   s-coordinate case, cubic spline pressure Jacobian'
1067      ENDIF
1068
1069      !!----------------------------------------------------------------------
1070      ! Local constant initialization
1071      zcoef0 = - grav 
1072      znad = 0.0_wp
1073      IF(lk_vvl) znad = 1._wp
1074
1075      ! Save fsde3w and rhd
1076      fsde3w_tmp(:,:,:) = fsde3w(:,:,:) 
1077      rhd_tmp(:,:,:)    = rhd(:,:,:)   
1078
1079      ! Clean 3-D work arraies
1080      zhpi(:,:,:) = 0.
1081     
1082      !how to add vector opt.? N.B., jpi&jpi rather than jpim1&jpjm1 are needed here
1083
1084      ! Preparing vertical density profile for hybrid-sco coordinate
1085      DO jj = 1, jpj
1086        DO ji = 1, jpi   
1087          jk = mbathy(ji,jj)
1088          IF(jk <= 0) THEN; rhd(ji,jj,:) = 0._wp
1089            ELSE IF(jk == 1) THEN; rhd(ji,jj, jk+1:jpk) = rhd(ji,jj,jk)
1090            ELSE IF(jk < jpkm1) THEN
1091            DO jkk = jk+1, jpk
1092              rhd(ji,jj,jkk) = interp1(fsde3w(ji,jj,jkk), fsde3w(ji,jj,jkk-1),&
1093                                       fsde3w(ji,jj,jkk-2), rhd(ji,jj,jkk-1), rhd(ji,jj,jkk-2))
1094            END DO
1095          END IF
1096        END DO
1097      END DO
1098
1099      DO jj = 1, jpj
1100        DO ji = 1, jpi
1101          fsde3w(ji,jj,1) = 0.5_wp * fse3w(ji,jj,1)
1102          fsde3w(ji,jj,1) = fsde3w(ji,jj,1) - sshn(ji,jj) * znad
1103          DO jk = 2, jpk
1104            fsde3w(ji,jj,jk) = fsde3w(ji,jj,jk-1) + fse3w(ji,jj,jk)
1105          END DO
1106        END DO
1107      END DO
1108
1109      DO jk = 1, jpkm1
1110        DO jj = 1, jpj
1111          DO ji = 1, jpi
1112            fsp(ji,jj,jk) = rhd(ji,jj,jk)
1113            xsp(ji,jj,jk) = fsde3w(ji,jj,jk)
1114          END DO
1115        END DO
1116      END DO
1117
1118      !                 ! Construct the vertical density profile with the
1119      !                 !Constrained cubic spline interpolation
1120      CALL cspline(fsp,xsp,asp,bsp,csp,dsp,polynomial_type)     
1121
1122      ! Calculate the hydrostatic pressure at T(ji,jj,1)
1123      DO jj = 2, jpj
1124        DO ji = 2, jpi 
1125          rhdt1 = rhd(ji,jj,1) - interp3(fsde3w(ji,jj,1),asp(ji,jj,1), &
1126                                         bsp(ji,jj,1),csp(ji,jj,1),dsp(ji,jj,1)) &
1127                               * 0.5_wp * fsde3w(ji,jj,1)
1128          rhdt1 = max(rhdt1, 1000._wp - rau0)        ! no lighter than fresh water
1129
1130     !                  ! assuming linear profile across the top half surface layer
1131          zhpi(ji,jj,1) =  0.5_wp * fse3w(ji,jj,1) * rhdt1 
1132        END DO
1133      END DO
1134
1135      ! Calculate the pressure at T(ji,jj,2:jpkm1)
1136      DO jk = 2, jpkm1                                 
1137        DO jj = 2, jpj     
1138          DO ji = 2, jpi
1139            zhpi(ji,jj,jk) = zhpi(ji,jj,jk-1) + &
1140                             integ2(fsde3w(ji,jj,jk-1),fsde3w(ji,jj,jk),&
1141                                    asp(ji,jj,jk-1),bsp(ji,jj,jk-1),&
1142                                    csp(ji,jj,jk-1),dsp(ji,jj,jk-1))
1143          END DO
1144        END DO
1145      END DO
1146
1147      ! Z coordinate of U(ji,jj,1:jpkm1) and V(ji,jj,1:jpkm1)
1148
1149      DO jj = 2, jpjm1     
1150        DO ji = 2, jpim1 
1151          zu(ji,jj,1) = - ( fse3u(ji,jj,1) - sshu_n(ji,jj) * znad)
1152          zv(ji,jj,1) = - ( fse3v(ji,jj,1) - sshv_n(ji,jj) * znad)
1153        END DO
1154      END DO
1155
1156      DO jk = 2, jpkm1                                 
1157        DO jj = 2, jpjm1     
1158          DO ji = 2, jpim1 
1159            zu(ji,jj,jk) = zu(ji,jj,jk-1)- fse3u(ji,jj,jk)
1160            zv(ji,jj,jk) = zv(ji,jj,jk-1)- fse3v(ji,jj,jk)
1161          END DO
1162        END DO
1163      END DO
1164               
1165      DO jk = 1, jpkm1                                 
1166        DO jj = 2, jpjm1     
1167          DO ji = 2, jpim1 
1168            zu(ji,jj,jk) = zu(ji,jj,jk) + 0.5_wp * fse3u(ji,jj,jk)
1169            zv(ji,jj,jk) = zv(ji,jj,jk) + 0.5_wp * fse3v(ji,jj,jk)
1170          END DO
1171        END DO
1172      END DO
1173
1174      DO jk = 1, jpkm1                                 
1175        DO jj = 2, jpjm1     
1176          DO ji = 2, jpim1 
1177            pwes = 0._wp; pwed = 0._wp
1178            pnss = 0._wp; pnsd = 0._wp
1179            zuijk = zu(ji,jj,jk)
1180            zvijk = zv(ji,jj,jk)
1181
1182            !!!!!     for u equation
1183            IF(jk <= mbku(ji,jj)) THEN
1184            IF(-fsde3w(ji+1,jj,mbku(ji,jj)) >= -fsde3w(ji,jj,mbku(ji,jj))) THEN
1185              jis = ji + 1; jid = ji
1186            ELSE
1187              jis = ji;     jid = ji +1
1188            ENDIF
1189
1190            !     !integrate the pressure on the shallow side
1191            jk1 = jk 
1192            bhitwe = 0
1193            DO WHILE ( -fsde3w(jis,jj,jk1) > zuijk )
1194              IF(jk1 == mbku(ji,jj)) THEN
1195                bhitwe = 1
1196                EXIT
1197              ENDIF
1198              deps = min(fsde3w(jis,jj,jk1+1), -zuijk)
1199              pwes = pwes + & 
1200                   integ2(fsde3w(jis,jj,jk1),deps,&
1201                          asp(jis,jj,jk1),bsp(jis,jj,jk1),&
1202                          csp(jis,jj,jk1),dsp(jis,jj,jk1))
1203              jk1 = jk1 + 1
1204            END DO
1205           
1206            IF(bhitwe == 1) THEN
1207              zuijk = -fsde3w(jis,jj,jk1)
1208            ENDIF
1209
1210            !     !integrate the pressure on the deep side
1211            jk1 = jk 
1212            bhitwe = 0
1213            DO WHILE ( -fsde3w(jid,jj,jk1) < zuijk )
1214              IF(jk1 == 1) THEN
1215                bhitwe = 1
1216                EXIT
1217              END IF
1218              deps = max(fsde3w(jid,jj,jk1-1), -zuijk)
1219              pwed = pwed + & 
1220                   integ2(deps,fsde3w(jid,jj,jk1),&
1221                          asp(jid,jj,jk1-1),bsp(jid,jj,jk1-1),&
1222                          csp(jid,jj,jk1-1),dsp(jid,jj,jk1-1))
1223              jk1 = jk1 - 1
1224            END DO
1225           
1226            IF(bhitwe == 1) THEN
1227              deps = fsde3w(jid,jj,1) + min(zuijk, sshn(jid,jj)*znad)
1228              rhdt1 = rhd(jid,jj,1) - interp3(fsde3w(jid,jj,1),asp(jid,jj,1), &
1229                                      bsp(jid,jj,1),csp(jid,jj,1),dsp(jid,jj,1)) * deps
1230              rhdt1 = max(rhdt1, 1000._wp - rau0)        ! no lighter than fresh water
1231              pwed = pwed + 0.5_wp * (rhd(jid,jj,1) + rhdt1) * deps
1232            ENDIF
1233
1234
1235            dpdx1 = zcoef0 / e1u(ji,jj) * (zhpi(ji+1,jj,jk) - zhpi(ji,jj,jk))
1236            IF(lk_vvl) THEN
1237              dpdx2 = zcoef0 / e1u(ji,jj) * & 
1238                      (REAL(jis-jid, wp) * (pwes + pwed) + (sshn(ji+1,jj)-sshn(ji,jj))) 
1239             ELSE
1240              dpdx2 = zcoef0 / e1u(ji,jj) * REAL(jis-jid, wp) * (pwes + pwed) 
1241            ENDIF
1242
1243            ua(ji,jj,jk) = ua(ji,jj,jk) + (dpdx1 + dpdx2)*&
1244               &           umask(ji,jj,jk)*tmask(ji,jj,jk)*tmask(ji+1,jj,jk)
1245            ENDIF
1246 
1247            !!!!!     for v equation
1248            IF(jk <= mbkv(ji,jj)) THEN
1249            IF(-fsde3w(ji,jj+1,mbkv(ji,jj)) >= -fsde3w(ji,jj,mbkv(ji,jj))) THEN
1250              jjs = jj + 1; jjd = jj
1251            ELSE
1252              jjs = jj    ; jjd = jj + 1
1253            ENDIF
1254
1255            !     !integrate the pressure on the shallow side
1256            jk1 = jk 
1257            bhitns = 0
1258            DO WHILE ( -fsde3w(ji,jjs,jk1) > zvijk )
1259              IF(jk1 == mbkv(ji,jj)) THEN
1260                bhitns = 1
1261                EXIT
1262              ENDIF
1263              deps = min(fsde3w(ji,jjs,jk1+1), -zvijk)
1264              pnss = pnss + & 
1265                   integ2(fsde3w(ji,jjs,jk1),deps,&
1266                          asp(ji,jjs,jk1),bsp(ji,jjs,jk1),&
1267                          csp(ji,jjs,jk1),dsp(ji,jjs,jk1))
1268              jk1 = jk1 + 1
1269            END DO
1270           
1271            IF(bhitns == 1) THEN
1272              zvijk = -fsde3w(ji,jjs,jk1)
1273            ENDIF
1274
1275            !     !integrate the pressure on the deep side
1276            jk1 = jk 
1277            bhitns = 0
1278            DO WHILE ( -fsde3w(ji,jjd,jk1) < zvijk )
1279              IF(jk1 == 1) THEN
1280                bhitns = 1
1281                EXIT
1282              END IF
1283              deps = max(fsde3w(ji,jjd,jk1-1), -zvijk)
1284              pnsd = pnsd + & 
1285                   integ2(deps,fsde3w(ji,jjd,jk1),&
1286                          asp(ji,jjd,jk1-1),bsp(ji,jjd,jk1-1),&
1287                          csp(ji,jjd,jk1-1),dsp(ji,jjd,jk1-1))
1288              jk1 = jk1 - 1
1289            END DO
1290           
1291            IF(bhitns == 1) THEN
1292              deps = fsde3w(ji,jjd,1) + min(zvijk, sshn(ji,jjd)*znad)
1293              rhdt1 = rhd(ji,jjd,1) - interp3(fsde3w(ji,jjd,1),asp(ji,jjd,1), &
1294                                      bsp(ji,jjd,1),csp(ji,jjd,1),dsp(ji,jjd,1)) * deps
1295              rhdt1 = max(rhdt1, 1000._wp - rau0)        ! no lighter than fresh water
1296              pnsd = pnsd + 0.5_wp * (rhd(ji,jjd,1) + rhdt1) * deps
1297            ENDIF
1298
1299
1300            dpdy1 = zcoef0 / e2v(ji,jj) * (zhpi(ji,jj+1,jk) - zhpi(ji,jj,jk))
1301            if(lk_vvl) then
1302                dpdy2 = zcoef0 / e2v(ji,jj) * &
1303                        (REAL(jjs-jjd, wp) * (pnss + pnsd) + (sshn(ji,jj+1)-sshn(ji,jj))) 
1304              else
1305                dpdy2 = zcoef0 / e2v(ji,jj) * REAL(jjs-jjd, wp) * (pnss + pnsd ) 
1306            end if
1307
1308            va(ji,jj,jk) = va(ji,jj,jk) + (dpdy1 + dpdy2)*&
1309            &              vmask(ji,jj,jk)*tmask(ji,jj,jk)*tmask(ji,jj+1,jk)
1310            ENDIF
1311
1312                   
1313           END DO
1314        END DO
1315      END DO
1316
1317      ! Restore fsde3w and rhd
1318      fsde3w(:,:,:) = fsde3w_tmp(:,:,:)
1319      rhd(:,:,:)    = rhd_tmp(:,:,:)
1320
1321      !
1322      IF( wrk_not_released(3, 3,4,5,6,7,8,9,10,11) )   &
1323         CALL ctl_stop('dyn:hpg_prj: failed to release workspace arrays')
1324      !
1325   END SUBROUTINE hpg_prj
1326
1327   SUBROUTINE cspline(fsp, xsp, asp, bsp, csp, dsp, polynomial_type)
1328      !!----------------------------------------------------------------------
1329      !!                 ***  ROUTINE cspline  ***
1330      !!       
1331      !! ** Purpose :   constrained cubic spline interpolation
1332      !!         
1333      !! ** Method  :   f(x) = asp + bsp*x + csp*x^2 + dsp*x^3
1334      !! Reference: K.W. Brodlie, A review of mehtods for curve and function
1335      !!                          drawing, 1980
1336      !!
1337      !!----------------------------------------------------------------------
1338      IMPLICIT NONE
1339      REAL(wp), DIMENSION(:,:,:), INTENT(in)  :: fsp, xsp           ! value and coordinate
1340      REAL(wp), DIMENSION(:,:,:), INTENT(out) :: asp, bsp, csp, dsp ! coefficients of
1341                                                                    ! the interpoated function
1342      INTEGER, INTENT(in) :: polynomial_type                        ! 1: cubic spline
1343                                                                    ! 2: Linear
1344
1345      ! Local Variables     
1346      INTEGER  ::   ji, jj, jk                 ! dummy loop indices
1347      INTEGER  ::   jpi, jpj, jpkm1
1348      REAL(wp) :: df1, df2, ddf1, ddf2, tmp1, tmp2, dxtmp
1349      REAL(wp) :: dxtmp1, dxtmp2, alpha
1350      REAL(wp) :: df(size(fsp,3))
1351      !!----------------------------------------------------------------------
1352
1353      jpi   = size(fsp,1)
1354      jpj   = size(fsp,2)
1355      jpkm1 = size(fsp,3) - 1
1356
1357      ! Clean output arrays
1358      asp = 0.0_wp
1359      bsp = 0.0_wp
1360      csp = 0.0_wp
1361      dsp = 0.0_wp
1362     
1363     
1364      Do ji = 1, jpi
1365      Do jj = 1, jpj
1366
1367      If (polynomial_type == 1) Then     !Constrained Cubic Spline
1368          Do jk = 2, jpkm1-1
1369             dxtmp1 = xsp(ji,jj,jk)   - xsp(ji,jj,jk-1) 
1370             dxtmp2 = xsp(ji,jj,jk+1) - xsp(ji,jj,jk) 
1371             df1    = ( fsp(ji,jj,jk)   - fsp(ji,jj,jk-1) ) / dxtmp1
1372             df2    = ( fsp(ji,jj,jk+1) - fsp(ji,jj,jk)   ) / dxtmp2
1373 
1374             alpha = ( dxtmp1 + 2._wp * dxtmp2 ) / ( dxtmp1 + dxtmp2 ) / 3._wp
1375             
1376             If(df1 * df2 <= 0._wp) Then
1377               df(jk) = 0._wp
1378             Else
1379               df(jk) = df1 * df2 / ( ( 1._wp - alpha ) * df1 + alpha * df2 )
1380             End If
1381          End Do
1382 
1383          df(1)     = 1.5_wp * ( fsp(ji,jj,2) - fsp(ji,jj,1) ) / ( xsp(ji,jj,2) - xsp(ji,jj,1) ) - &
1384                    & 0.5_wp * df(2)
1385          df(jpkm1) = 1.5_wp * ( fsp(ji,jj,jpkm1) - fsp(ji,jj,jpkm1-1) ) / &
1386                    &          ( xsp(ji,jj,jpkm1) - xsp(ji,jj,jpkm1-1) ) - &
1387                    & 0.5_wp * df(jpkm1 - 1)
1388 
1389          Do jk = 1, jpkm1 - 1
1390             dxtmp = xsp(ji,jj,jk+1) - xsp(ji,jj,jk) 
1391             tmp1  = (df(jk+1) + 2._wp * df(jk)) / dxtmp
1392             tmp2  =  6._wp * (fsp(ji,jj,jk+1) - fsp(ji,jj,jk)) / dxtmp / dxtmp
1393             ddf1  = -2._wp * tmp1 + tmp2 
1394             tmp1  = (2._wp * df(jk+1) + df(jk)) / dxtmp
1395             ddf2  =  2._wp * tmp1 - tmp2 
1396 
1397             dsp(ji,jj,jk) = (ddf2 - ddf1) / 6._wp / dxtmp
1398             csp(ji,jj,jk) = ( xsp(ji,jj,jk+1) * ddf1 - xsp(ji,jj,jk)*ddf2 ) / 2._wp / dxtmp
1399             bsp(ji,jj,jk) = ( fsp(ji,jj,jk+1) - fsp(ji,jj,jk) ) / dxtmp - & 
1400                           & csp(ji,jj,jk) * ( xsp(ji,jj,jk+1) + xsp(ji,jj,jk) ) - &
1401                           & dsp(ji,jj,jk) * ( xsp(ji,jj,jk+1)**2 + &
1402                           &                   xsp(ji,jj,jk+1) * xsp(ji,jj,jk) + &
1403                           &                   xsp(ji,jj,jk)**2 )
1404             asp(ji,jj,jk) = fsp(ji,jj,jk) - bsp(ji,jj,jk) * xsp(ji,jj,jk) - &
1405                           &                 csp(ji,jj,jk) * xsp(ji,jj,jk)**2 - &
1406                           &                 dsp(ji,jj,jk) * xsp(ji,jj,jk)**3
1407          End Do
1408
1409       Else If (polynomial_type == 2) Then     !Linear
1410
1411        Do jk = 1, jpkm1-1
1412           dxtmp =xsp(ji,jj,jk+1) - xsp(ji,jj,jk) 
1413           tmp1 = fsp(ji,jj,jk+1) - fsp(ji,jj,jk)
1414 
1415           dsp(ji,jj,jk) = 0._wp
1416           csp(ji,jj,jk) = 0._wp
1417           bsp(ji,jj,jk) = tmp1 / dxtmp
1418           asp(ji,jj,jk) = fsp(ji,jj,jk) - bsp(ji,jj,jk) * xsp(ji,jj,jk)
1419        End Do
1420
1421       Else
1422        CALL ctl_stop( 'invalid polynomial type in cspline' )
1423      End If
1424
1425      End Do
1426      End Do
1427     
1428   End Subroutine cspline
1429
1430
1431   FUNCTION interp1(x, xl, xr, fl, fr)  RESULT(f) 
1432      !!----------------------------------------------------------------------
1433      !!                 ***  ROUTINE interp1  ***
1434      !!       
1435      !! ** Purpose :   1-d linear interpolation
1436      !!         
1437      !! ** Method  : 
1438      !!                interpolation is straigt forward
1439      !!                extrapolation is also permitted (no value limit)
1440      !!
1441      !! H.Liu, Jan 2009,  POL
1442      !!----------------------------------------------------------------------
1443      IMPLICIT NONE
1444      REAL(wp), INTENT(in) ::  x, xl, xr, fl, fr   
1445      REAL(wp)             ::  f ! result of the interpolation (extrapolation)
1446      REAL(wp)             ::  deltx
1447      !!----------------------------------------------------------------------
1448
1449      deltx = xr - xl
1450      IF(abs(deltx) <= 10._wp * EPSILON(x)) THEN
1451        f = 0.5_wp * (fl + fr)
1452       ELSE
1453        f = ( (x - xl ) * fr - ( x - xr ) * fl ) / deltx
1454      END IF
1455     
1456   END FUNCTION interp1
1457
1458   FUNCTION interp2(x, a, b, c, d)  RESULT(f) 
1459      !!----------------------------------------------------------------------
1460      !!                 ***  ROUTINE interp1  ***
1461      !!       
1462      !! ** Purpose :   1-d constrained cubic spline interpolation
1463      !!         
1464      !! ** Method  :  cubic spline interpolation
1465      !!
1466      !!----------------------------------------------------------------------
1467      IMPLICIT NONE
1468      REAL(wp), INTENT(in) ::  x, a, b, c, d   
1469      REAL(wp)             ::  f ! value from the interpolation
1470      !!----------------------------------------------------------------------
1471
1472      f = a + x* ( b + x * ( c + d * x ) ) 
1473
1474   END FUNCTION interp2
1475
1476
1477   FUNCTION interp3(x, a, b, c, d)  RESULT(f) 
1478      !!----------------------------------------------------------------------
1479      !!                 ***  ROUTINE interp1  ***
1480      !!       
1481      !! ** Purpose :   deriavtive of a cubic spline function
1482      !!         
1483      !! ** Method  : 
1484      !!
1485      !!----------------------------------------------------------------------
1486      IMPLICIT NONE
1487      REAL(wp), INTENT(in) ::  x, a, b, c, d   
1488      REAL(wp)             ::  f ! value from the interpolation
1489      !!----------------------------------------------------------------------
1490
1491      f = b + x * ( 2._wp * c + 3._wp * d * x)
1492
1493   END FUNCTION interp3
1494
1495   
1496   FUNCTION integ2(xl, xr, a, b, c, d)  RESULT(f) 
1497      !!----------------------------------------------------------------------
1498      !!                 ***  ROUTINE interp1  ***
1499      !!       
1500      !! ** Purpose :   1-d constrained cubic spline integration
1501      !!         
1502      !! ** Method  :  integrate polynomial a+bx+cx^2+dx^3 from xl to xr
1503      !!
1504      !!----------------------------------------------------------------------
1505      IMPLICIT NONE
1506      REAL(wp), INTENT(in) ::  xl, xr, a, b, c, d   
1507      REAL(wp)             ::  a1, a2,a3     
1508      REAL(wp)             ::  f                   ! integration result
1509      !!----------------------------------------------------------------------
1510
1511      a1 = 0.5_wp * b 
1512      a2 = c / 3.0_wp 
1513      a3 = 0.25_wp * d 
1514
1515      f  = xr * ( a + xr * ( a1 + xr * ( a2 + a3 * xr ) ) ) - &
1516         & xl * ( a + xl * ( a1 + xl * ( a2 + a3 * xl ) ) )
1517
1518   END FUNCTION integ2
1519
1520
1521   !!======================================================================
1522END MODULE dynhpg
Note: See TracBrowser for help on using the repository browser.