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/UKMO/2014_Surge_Modelling/NEMOGCM/NEMO/OPA_SRC/DYN – NEMO

source: branches/UKMO/2014_Surge_Modelling/NEMOGCM/NEMO/OPA_SRC/DYN/dynhpg.F90 @ 5746

Last change on this file since 5746 was 5727, checked in by rfurner, 9 years ago

some bug fixes for wetting and drying elements...still not working though

  • Property svn:keywords set to Id
File size: 62.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   !!            3.4  !  2011-11  (H. Liu) hpg_prj: Original code for s-coordinates
17   !!                 !           (A. Coward) suppression of hel, wdj and rot options
18   !!            3.6? !  2014-09  (H. Liu) add Wetting/Drying pressure filter
19   !!----------------------------------------------------------------------
20
21   !!----------------------------------------------------------------------
22   !!   dyn_hpg      : update the momentum trend with the now horizontal
23   !!                  gradient of the hydrostatic pressure
24   !!   dyn_hpg_init : initialisation and control of options
25   !!       hpg_zco  : z-coordinate scheme
26   !!       hpg_zps  : z-coordinate plus partial steps (interpolation)
27   !!       hpg_sco  : s-coordinate (standard jacobian formulation)
28   !!       hpg_djc  : s-coordinate (Density Jacobian with Cubic polynomial)
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   USE wrk_nemo        ! Memory Allocation
41   USE timing          ! Timing
42
43   IMPLICIT NONE
44   PRIVATE
45
46   PUBLIC   dyn_hpg        ! routine called by step module
47   PUBLIC   dyn_hpg_init   ! routine called by opa module
48
49   !                                    !!* Namelist namdyn_hpg : hydrostatic pressure gradient
50   LOGICAL , PUBLIC ::   ln_hpg_zco      !: z-coordinate - full steps
51   LOGICAL , PUBLIC ::   ln_hpg_zps      !: z-coordinate - partial steps (interpolation)
52   LOGICAL , PUBLIC ::   ln_hpg_sco      !: s-coordinate (standard jacobian formulation)
53   LOGICAL , PUBLIC ::   ln_hpg_djc      !: s-coordinate (Density Jacobian with Cubic polynomial)
54   LOGICAL , PUBLIC ::   ln_hpg_prj      !: s-coordinate (Pressure Jacobian scheme)
55   LOGICAL , PUBLIC ::   ln_dynhpg_imp   !: semi-implicite hpg flag
56
57   INTEGER , PUBLIC ::   nhpg  =  0   ! = 0 to 7, type of pressure gradient scheme used ! (deduced from ln_hpg_... flags) (PUBLIC for TAM)
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      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
80      REAL(wp), POINTER, DIMENSION(:,:,:) ::  ztrdu, ztrdv
81      !!----------------------------------------------------------------------
82      !
83      IF( nn_timing == 1 )  CALL timing_start('dyn_hpg')
84      !
85      IF( l_trddyn ) THEN                    ! Temporary saving of ua and va trends (l_trddyn)
86         CALL wrk_alloc( jpi,jpj,jpk, ztrdu, ztrdv )
87         ztrdu(:,:,:) = ua(:,:,:)
88         ztrdv(:,:,:) = va(:,:,:)
89      ENDIF
90      !
91      SELECT CASE ( nhpg )      ! Hydrostatic pressure gradient computation
92      CASE (  0 )   ;   CALL hpg_zco    ( kt )      ! z-coordinate
93      CASE (  1 )   ;   CALL hpg_zps    ( kt )      ! z-coordinate plus partial steps (interpolation)
94      CASE (  2 )   ;   CALL hpg_sco    ( kt )      ! s-coordinate (standard jacobian formulation)
95      CASE (  3 )   ;   CALL hpg_djc    ( kt )      ! s-coordinate (Density Jacobian with Cubic polynomial)
96      CASE (  4 )   ;   CALL hpg_prj    ( kt )      ! s-coordinate (Pressure Jacobian scheme)
97      END SELECT
98      !
99      IF( l_trddyn ) THEN      ! save the hydrostatic pressure gradient trends for momentum trend diagnostics
100         ztrdu(:,:,:) = ua(:,:,:) - ztrdu(:,:,:)
101         ztrdv(:,:,:) = va(:,:,:) - ztrdv(:,:,:)
102         CALL trd_mod( ztrdu, ztrdv, jpdyn_trd_hpg, 'DYN', kt )
103         CALL wrk_dealloc( jpi,jpj,jpk, ztrdu, ztrdv )
104      ENDIF
105      !
106      IF(ln_ctl)   CALL prt_ctl( tab3d_1=ua, clinfo1=' hpg  - Ua: ', mask1=umask,   &
107         &                       tab3d_2=va, clinfo2=       ' Va: ', mask2=vmask, clinfo3='dyn' )
108      !
109      IF( nn_timing == 1 )  CALL timing_stop('dyn_hpg')
110      !
111   END SUBROUTINE dyn_hpg
112
113
114   SUBROUTINE dyn_hpg_init
115      !!----------------------------------------------------------------------
116      !!                 ***  ROUTINE dyn_hpg_init  ***
117      !!
118      !! ** Purpose :   initializations for the hydrostatic pressure gradient
119      !!              computation and consistency control
120      !!
121      !! ** Action  :   Read the namelist namdyn_hpg and check the consistency
122      !!      with the type of vertical coordinate used (zco, zps, sco)
123      !!----------------------------------------------------------------------
124      INTEGER ::   ioptio = 0      ! temporary integer
125      INTEGER ::   ios             ! Local integer output status for namelist read
126      !!
127      NAMELIST/namdyn_hpg/ ln_hpg_zco, ln_hpg_zps, ln_hpg_sco,     &
128         &                 ln_hpg_djc, ln_hpg_prj, ln_dynhpg_imp
129      !!----------------------------------------------------------------------
130      !
131      REWIND( numnam_ref )              ! Namelist namdyn_hpg in reference namelist : Hydrostatic pressure gradient
132      READ  ( numnam_ref, namdyn_hpg, IOSTAT = ios, ERR = 901)
133901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namdyn_hpg in reference namelist', lwp )
134
135      REWIND( numnam_cfg )              ! Namelist namdyn_hpg in configuration namelist : Hydrostatic pressure gradient
136      READ  ( numnam_cfg, namdyn_hpg, IOSTAT = ios, ERR = 902 )
137902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namdyn_hpg in configuration namelist', lwp )
138      IF(lwm) WRITE ( numond, 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. (Density Jacobian: Cubic polynomial)     ln_hpg_djc    = ', ln_hpg_djc
149         WRITE(numout,*) '      s-coord. (Pressure Jacobian: Cubic polynomial)    ln_hpg_prj    = ', ln_hpg_prj
150         WRITE(numout,*) '      time stepping: centered (F) or semi-implicit (T)  ln_dynhpg_imp = ', ln_dynhpg_imp
151      ENDIF
152      !
153      IF( ln_hpg_djc )   &
154         &   CALL ctl_stop('dyn_hpg_init : Density Jacobian: Cubic polynominal method &
155                           & currently disabled (bugs under investigation). Please select &
156                           & either  ln_hpg_sco or  ln_hpg_prj instead')
157      !
158      IF( lk_vvl .AND. .NOT. (ln_hpg_sco.OR.ln_hpg_prj) )   &
159         &   CALL ctl_stop('dyn_hpg_init : variable volume key_vvl requires:&
160                           & the standard jacobian formulation hpg_sco or &
161                           & the pressure jacobian formulation hpg_prj')
162      !
163      !                               ! Set nhpg from ln_hpg_... flags
164      IF( ln_hpg_zco )   nhpg = 0
165      IF( ln_hpg_zps )   nhpg = 1
166      IF( ln_hpg_sco )   nhpg = 2
167      IF( ln_hpg_djc )   nhpg = 3
168      IF( ln_hpg_prj )   nhpg = 4
169      !
170      !                               ! Consistency check
171      ioptio = 0
172      IF( ln_hpg_zco )   ioptio = ioptio + 1
173      IF( ln_hpg_zps )   ioptio = ioptio + 1
174      IF( ln_hpg_sco )   ioptio = ioptio + 1
175      IF( ln_hpg_djc )   ioptio = ioptio + 1
176      IF( ln_hpg_prj )   ioptio = ioptio + 1
177      IF( ioptio /= 1 )   CALL ctl_stop( 'NO or several hydrostatic pressure gradient options used' )
178      !
179   END SUBROUTINE dyn_hpg_init
180
181
182   SUBROUTINE hpg_zco( kt )
183      !!---------------------------------------------------------------------
184      !!                  ***  ROUTINE hpg_zco  ***
185      !!
186      !! ** Method  :   z-coordinate case, levels are horizontal surfaces.
187      !!      The now hydrostatic pressure gradient at a given level, jk,
188      !!      is computed by taking the vertical integral of the in-situ
189      !!      density gradient along the model level from the suface to that
190      !!      level:    zhpi = grav .....
191      !!                zhpj = grav .....
192      !!      add it to the general momentum trend (ua,va).
193      !!            ua = ua - 1/e1u * zhpi
194      !!            va = va - 1/e2v * zhpj
195      !!
196      !! ** Action : - Update (ua,va) with the now hydrastatic pressure trend
197      !!----------------------------------------------------------------------
198      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
199      !!
200      INTEGER  ::   ji, jj, jk       ! dummy loop indices
201      REAL(wp) ::   zcoef0, zcoef1   ! temporary scalars
202      REAL(wp), POINTER, DIMENSION(:,:,:) ::  zhpi, zhpj
203      !!----------------------------------------------------------------------
204      !
205      CALL wrk_alloc( jpi,jpj,jpk, zhpi, zhpj )
206      !
207      IF( kt == nit000 ) THEN
208         IF(lwp) WRITE(numout,*)
209         IF(lwp) WRITE(numout,*) 'dyn:hpg_zco : hydrostatic pressure gradient trend'
210         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   z-coordinate case '
211      ENDIF
212
213      zcoef0 = - grav * 0.5_wp      ! Local constant initialization
214
215      ! Surface value
216      DO jj = 2, jpjm1
217         DO ji = fs_2, fs_jpim1   ! vector opt.
218            zcoef1 = zcoef0 * fse3w(ji,jj,1)
219            ! hydrostatic pressure gradient
220            zhpi(ji,jj,1) = zcoef1 * ( rhd(ji+1,jj,1) - rhd(ji,jj,1) ) / e1u(ji,jj)
221            zhpj(ji,jj,1) = zcoef1 * ( rhd(ji,jj+1,1) - rhd(ji,jj,1) ) / e2v(ji,jj)
222            ! add to the general momentum trend
223            ua(ji,jj,1) = ua(ji,jj,1) + zhpi(ji,jj,1)
224            va(ji,jj,1) = va(ji,jj,1) + zhpj(ji,jj,1)
225         END DO
226      END DO
227
228      !
229      ! interior value (2=<jk=<jpkm1)
230      DO jk = 2, jpkm1
231         DO jj = 2, jpjm1
232            DO ji = fs_2, fs_jpim1   ! vector opt.
233               zcoef1 = zcoef0 * fse3w(ji,jj,jk)
234               ! hydrostatic pressure gradient
235               zhpi(ji,jj,jk) = zhpi(ji,jj,jk-1)   &
236                  &           + zcoef1 * (  ( rhd(ji+1,jj,jk)+rhd(ji+1,jj,jk-1) )   &
237                  &                       - ( rhd(ji  ,jj,jk)+rhd(ji  ,jj,jk-1) )  ) / e1u(ji,jj)
238
239               zhpj(ji,jj,jk) = zhpj(ji,jj,jk-1)   &
240                  &           + zcoef1 * (  ( rhd(ji,jj+1,jk)+rhd(ji,jj+1,jk-1) )   &
241                  &                       - ( rhd(ji,jj,  jk)+rhd(ji,jj  ,jk-1) )  ) / e2v(ji,jj)
242               ! add to the general momentum trend
243               ua(ji,jj,jk) = ua(ji,jj,jk) + zhpi(ji,jj,jk)
244               va(ji,jj,jk) = va(ji,jj,jk) + zhpj(ji,jj,jk)
245            END DO
246         END DO
247      END DO
248      !
249      CALL wrk_dealloc( jpi,jpj,jpk, zhpi, zhpj )
250      !
251   END SUBROUTINE hpg_zco
252
253
254   SUBROUTINE hpg_zps( kt )
255      !!---------------------------------------------------------------------
256      !!                 ***  ROUTINE hpg_zps  ***
257      !!
258      !! ** Method  :   z-coordinate plus partial steps case.  blahblah...
259      !!
260      !! ** Action  : - Update (ua,va) with the now hydrastatic pressure trend
261      !!----------------------------------------------------------------------
262      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
263      !!
264      INTEGER  ::   ji, jj, jk                       ! dummy loop indices
265      INTEGER  ::   iku, ikv                         ! temporary integers
266      REAL(wp) ::   zcoef0, zcoef1, zcoef2, zcoef3   ! temporary scalars
267      REAL(wp), POINTER, DIMENSION(:,:,:) ::  zhpi, zhpj
268      !!----------------------------------------------------------------------
269      !
270      CALL wrk_alloc( jpi,jpj,jpk, zhpi, zhpj )
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
279      ! Local constant initialization
280      zcoef0 = - grav * 0.5_wp
281
282      !  Surface value (also valid in partial step case)
283      DO jj = 2, jpjm1
284         DO ji = fs_2, fs_jpim1   ! vector opt.
285            zcoef1 = zcoef0 * fse3w(ji,jj,1)
286            ! hydrostatic pressure gradient
287            zhpi(ji,jj,1) = zcoef1 * ( rhd(ji+1,jj  ,1) - rhd(ji,jj,1) ) / e1u(ji,jj)
288            zhpj(ji,jj,1) = zcoef1 * ( rhd(ji  ,jj+1,1) - rhd(ji,jj,1) ) / e2v(ji,jj)
289            ! add to the general momentum trend
290            ua(ji,jj,1) = ua(ji,jj,1) + zhpi(ji,jj,1)
291            va(ji,jj,1) = va(ji,jj,1) + zhpj(ji,jj,1)
292         END DO
293      END DO
294
295
296      ! interior value (2=<jk=<jpkm1)
297      DO jk = 2, jpkm1
298         DO jj = 2, jpjm1
299            DO ji = fs_2, fs_jpim1   ! vector opt.
300               zcoef1 = zcoef0 * fse3w(ji,jj,jk)
301               ! hydrostatic pressure gradient
302               zhpi(ji,jj,jk) = zhpi(ji,jj,jk-1)   &
303                  &           + zcoef1 * (  ( rhd(ji+1,jj,jk) + rhd(ji+1,jj,jk-1) )   &
304                  &                       - ( rhd(ji  ,jj,jk) + rhd(ji  ,jj,jk-1) )  ) / e1u(ji,jj)
305
306               zhpj(ji,jj,jk) = zhpj(ji,jj,jk-1)   &
307                  &           + zcoef1 * (  ( rhd(ji,jj+1,jk) + rhd(ji,jj+1,jk-1) )   &
308                  &                       - ( rhd(ji,jj,  jk) + rhd(ji,jj  ,jk-1) )  ) / e2v(ji,jj)
309               ! add to the general momentum trend
310               ua(ji,jj,jk) = ua(ji,jj,jk) + zhpi(ji,jj,jk)
311               va(ji,jj,jk) = va(ji,jj,jk) + zhpj(ji,jj,jk)
312            END DO
313         END DO
314      END DO
315
316
317      ! partial steps correction at the last level  (use gru & grv computed in zpshde.F90)
318# if defined key_vectopt_loop
319         jj = 1
320         DO ji = jpi+2, jpij-jpi-1   ! vector opt. (forced unrolling)
321# else
322      DO jj = 2, jpjm1
323         DO ji = 2, jpim1
324# endif
325            iku = mbku(ji,jj)
326            ikv = mbkv(ji,jj)
327            zcoef2 = zcoef0 * MIN( fse3w(ji,jj,iku), fse3w(ji+1,jj  ,iku) )
328            zcoef3 = zcoef0 * MIN( fse3w(ji,jj,ikv), fse3w(ji  ,jj+1,ikv) )
329            IF( iku > 1 ) THEN            ! on i-direction (level 2 or more)
330               ua  (ji,jj,iku) = ua(ji,jj,iku) - zhpi(ji,jj,iku)         ! subtract old value
331               zhpi(ji,jj,iku) = zhpi(ji,jj,iku-1)                   &   ! compute the new one
332                  &            + zcoef2 * ( rhd(ji+1,jj,iku-1) - rhd(ji,jj,iku-1) + gru(ji,jj) ) / e1u(ji,jj)
333               ua  (ji,jj,iku) = ua(ji,jj,iku) + zhpi(ji,jj,iku)         ! add the new one to the general momentum trend
334            ENDIF
335            IF( ikv > 1 ) THEN            ! on j-direction (level 2 or more)
336               va  (ji,jj,ikv) = va(ji,jj,ikv) - zhpj(ji,jj,ikv)         ! subtract old value
337               zhpj(ji,jj,ikv) = zhpj(ji,jj,ikv-1)                   &   ! compute the new one
338                  &            + zcoef3 * ( rhd(ji,jj+1,ikv-1) - rhd(ji,jj,ikv-1) + grv(ji,jj) ) / e2v(ji,jj)
339               va  (ji,jj,ikv) = va(ji,jj,ikv) + zhpj(ji,jj,ikv)         ! add the new one to the general momentum trend
340            ENDIF
341# if ! defined key_vectopt_loop
342         END DO
343# endif
344      END DO
345      !
346      CALL wrk_dealloc( jpi,jpj,jpk, zhpi, zhpj )
347      !
348   END SUBROUTINE hpg_zps
349
350
351   SUBROUTINE hpg_sco( kt )
352      !!---------------------------------------------------------------------
353      !!                  ***  ROUTINE hpg_sco  ***
354      !!
355      !! ** Method  :   s-coordinate case. Jacobian scheme.
356      !!      The now hydrostatic pressure gradient at a given level, jk,
357      !!      is computed by taking the vertical integral of the in-situ
358      !!      density gradient along the model level from the suface to that
359      !!      level. s-coordinates (ln_sco): a corrective term is added
360      !!      to the horizontal pressure gradient :
361      !!         zhpi = grav .....  + 1/e1u mi(rhd) di[ grav dep3w ]
362      !!         zhpj = grav .....  + 1/e2v mj(rhd) dj[ grav dep3w ]
363      !!      add it to the general momentum trend (ua,va).
364      !!         ua = ua - 1/e1u * zhpi
365      !!         va = va - 1/e2v * zhpj
366      !!
367      !! ** Action : - Update (ua,va) with the now hydrastatic pressure trend
368      !!----------------------------------------------------------------------
369      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
370      !!
371      INTEGER  ::   ji, jj, jk                 ! dummy loop indices
372      REAL(wp) ::   zcoef0, zuap, zvap, znad, ztmp   ! temporary scalars
373      LOGICAL  ::   ll_tmp1, ll_tmp2           ! local logical variables
374      REAL(wp), POINTER, DIMENSION(:,:,:) ::  zhpi, zhpj
375      REAL(wp), POINTER, DIMENSION(:,:)   ::  zcpx, zcpy    !W/D pressure filter
376      !!----------------------------------------------------------------------
377      !
378      CALL wrk_alloc( jpi,jpj,jpk, zhpi, zhpj )
379      IF(ln_wd) CALL wrk_alloc( jpi,jpj, zcpx, zcpy )
380      !
381      IF( kt == nit000 ) THEN
382         IF(lwp) WRITE(numout,*)
383         IF(lwp) WRITE(numout,*) 'dyn:hpg_sco : hydrostatic pressure gradient trend'
384         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   s-coordinate case, OPA original scheme used'
385      ENDIF
386
387      ! Local constant initialization
388      zcoef0 = - grav * 0.5_wp
389      ! To use density and not density anomaly
390      IF ( lk_vvl ) THEN   ;     znad = 1._wp          ! Variable volume
391      ELSE                 ;     znad = 0._wp         ! Fixed volume
392      ENDIF
393
394      IF(ln_wd) THEN
395        DO jj = 2, jpjm1
396           DO ji = 2, jpim1
397             IF ( tmask(ji+1,jj,1) == 0._wp) THEN
398                zcpx(ji,jj) = 1.0_wp
399             ELSE
400                ll_tmp1 = MIN(sshn(ji,jj), sshn(ji+1,jj)) > MAX(-bathy(ji,jj), -bathy(ji+1,jj))
401                ll_tmp2 = MAX(sshn(ji,jj), sshn(ji+1,jj)) > MAX(-bathy(ji,jj), -bathy(ji+1,jj)) +&
402                                                          & rn_wdmin1 + rn_wdmin2
403   
404                IF(ll_tmp1) THEN
405                  zcpx(ji,jj) = 1.0_wp
406                ELSE IF(ll_tmp2) THEN
407                   ! no worries about sshn(ji+1,jj)-sshn(ji,jj) = 0, it won't happen ! here
408                   zcpx(ji,jj) = ABS((sshn(ji+1,jj) + bathy(ji+1,jj) - sshn(ji,jj) - bathy(ji,jj)) /&
409                               &     (sshn(ji+1,jj) - sshn(ji,jj)))
410                ELSE
411                  zcpx(ji,jj) = 0._wp
412                END IF
413             ENDIF
414       
415             IF ( tmask(ji,jj+1,1) == 0._wp) THEN
416                zcpy(ji,jj) = 1.0_wp
417             ELSE
418                ll_tmp1 = MIN(sshn(ji,jj), sshn(ji,jj+1)) > MAX(-bathy(ji,jj), -bathy(ji,jj+1))
419                ll_tmp2 = MAX(sshn(ji,jj), sshn(ji,jj+1)) > MAX(-bathy(ji,jj), -bathy(ji,jj+1)) +&
420                                                       & rn_wdmin1 + rn_wdmin2
421
422                IF(ll_tmp1) THEN
423                  zcpy(ji,jj) = 1.0_wp
424                ELSE IF(ll_tmp2) THEN
425                  ! no worries about sshn(ji,jj+1)-sshn(ji,jj) = 0, it won't happen ! here
426                  zcpy(ji,jj) = ABS((sshn(ji,jj+1) + bathy(ji,jj+1) - sshn(ji,jj) - bathy(ji,jj)) /&
427                              &     (sshn(ji,jj+1) - sshn(ji,jj)))
428                ELSE
429                  zcpy(ji,jj) = 0._wp
430                END IF
431             ENDIF
432           END DO
433        END DO
434        CALL lbc_lnk( zcpx, 'U', 1._wp )    ;   CALL lbc_lnk( zcpy, 'V', 1._wp )
435      ENDIF
436
437#if defined key_surge
438      ! Surface value
439      DO jj = 2, jpjm1
440         DO ji = fs_2, fs_jpim1   ! vector opt.
441            ! hydrostatic pressure gradient along s-surfaces
442            zhpi(ji,jj,1) = zcoef0 / e1u(ji,jj) * ( fse3w(ji+1,jj  ,1) * ( znad )   &
443               &                                  - fse3w(ji  ,jj  ,1) * ( znad ) )
444            zhpj(ji,jj,1) = zcoef0 / e2v(ji,jj) * ( fse3w(ji  ,jj+1,1) * ( znad )   &
445               &                                  - fse3w(ji  ,jj  ,1) * ( znad ) )
446            ! s-coordinate pressure gradient correction
447            zuap = -zcoef0 * ( 2._wp * znad )   &
448               &           * ( fsde3w(ji+1,jj,1) - fsde3w(ji,jj,1) ) / e1u(ji,jj)
449            zvap = -zcoef0 * ( 2._wp * znad )   &
450               &           * ( fsde3w(ji,jj+1,1) - fsde3w(ji,jj,1) ) / e2v(ji,jj)
451
452            IF(ln_wd) THEN
453              zhpi(ji,jj,1) = zhpi(ji,jj,1) * zcpx(ji,jj)
454              zhpj(ji,jj,1) = zhpj(ji,jj,1) * zcpy(ji,jj)
455              zuap = zuap * zcpx(ji,jj)
456              zvap = zvap * zcpy(ji,jj)
457            ENDIF
458
459            ! add to the general momentum trend
460            ua(ji,jj,1) = ua(ji,jj,1) + ( zhpi(ji,jj,1) + zuap ) * umask(ji,jj,1)
461            va(ji,jj,1) = va(ji,jj,1) + ( zhpj(ji,jj,1) + zvap ) * vmask(ji,jj,1)
462         END DO
463      END DO
464
465      ! interior value (2=<jk=<jpkm1)
466      DO jk = 2, jpkm1
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,jk) = zhpi(ji,jj,jk-1) + zcoef0 / e1u(ji,jj)   &
471                  &           * (  fse3w(ji+1,jj,jk) * ( 2*znad )   &
472                  &              - fse3w(ji  ,jj,jk) * ( 2*znad )  )
473               zhpj(ji,jj,jk) = zhpj(ji,jj,jk-1) + zcoef0 / e2v(ji,jj)   &
474                  &           * (  fse3w(ji,jj+1,jk) * ( 2*znad )   &
475                  &              - fse3w(ji,jj  ,jk) * ( 2*znad )  )
476               ! s-coordinate pressure gradient correction
477               zuap = -zcoef0 * ( 2._wp * znad )   &
478                  &           * ( fsde3w(ji+1,jj  ,jk) - fsde3w(ji,jj,jk) ) / e1u(ji,jj)
479               zvap = -zcoef0 * ( 2._wp * znad )   &
480                  &           * ( fsde3w(ji  ,jj+1,jk) - fsde3w(ji,jj,jk) ) / e2v(ji,jj)
481
482               IF(ln_wd) THEN
483                 zhpi(ji,jj,jk) = zhpi(ji,jj,jk) * zcpx(ji,jj)
484                 zhpj(ji,jj,jk) = zhpj(ji,jj,jk) * zcpy(ji,jj)
485                 zuap = zuap * zcpx(ji,jj)
486                 zvap = zvap * zcpy(ji,jj)
487               ENDIF
488
489               ! add to the general momentum trend
490               ua(ji,jj,jk) = ua(ji,jj,jk) + ( zhpi(ji,jj,jk) + zuap ) * umask(ji,jj,jk)
491               va(ji,jj,jk) = va(ji,jj,jk) + ( zhpj(ji,jj,jk) + zvap ) * vmask(ji,jj,jk)
492            END DO
493         END DO
494      END DO
495      !
496#else
497      ! Surface value
498      DO jj = 2, jpjm1
499         DO ji = fs_2, fs_jpim1   ! vector opt.
500            ! hydrostatic pressure gradient along s-surfaces
501            zhpi(ji,jj,1) = zcoef0 / e1u(ji,jj) * ( fse3w(ji+1,jj  ,1) * ( znad + rhd(ji+1,jj  ,1) )   &
502               &                                  - fse3w(ji  ,jj  ,1) * ( znad + rhd(ji  ,jj  ,1) ) )
503            zhpj(ji,jj,1) = zcoef0 / e2v(ji,jj) * ( fse3w(ji  ,jj+1,1) * ( znad + rhd(ji  ,jj+1,1) )   &
504               &                                  - fse3w(ji  ,jj  ,1) * ( znad + rhd(ji  ,jj  ,1) ) )
505            ! s-coordinate pressure gradient correction
506            zuap = -zcoef0 * ( rhd   (ji+1,jj,1) + rhd   (ji,jj,1) + 2._wp * znad )   &
507               &           * ( fsde3w(ji+1,jj,1) - fsde3w(ji,jj,1) ) / e1u(ji,jj)
508            zvap = -zcoef0 * ( rhd   (ji,jj+1,1) + rhd   (ji,jj,1) + 2._wp * znad )   &
509               &           * ( fsde3w(ji,jj+1,1) - fsde3w(ji,jj,1) ) / e2v(ji,jj)
510
511            IF(ln_wd) THEN
512              zhpi(ji,jj,1) = zhpi(ji,jj,1) * zcpx(ji,jj)
513              zhpj(ji,jj,1) = zhpj(ji,jj,1) * zcpy(ji,jj)
514              zuap = zuap * zcpx(ji,jj)
515              zvap = zvap * zcpy(ji,jj)
516            ENDIF
517
518            ! add to the general momentum trend
519            ua(ji,jj,1) = ua(ji,jj,1) + zhpi(ji,jj,1) + zuap
520            va(ji,jj,1) = va(ji,jj,1) + zhpj(ji,jj,1) + zvap
521         END DO
522      END DO
523
524      ! interior value (2=<jk=<jpkm1)
525      DO jk = 2, jpkm1
526         DO jj = 2, jpjm1
527            DO ji = fs_2, fs_jpim1   ! vector opt.
528               ! hydrostatic pressure gradient along s-surfaces
529               zhpi(ji,jj,jk) = zhpi(ji,jj,jk-1) + zcoef0 / e1u(ji,jj)   &
530                  &           * (  fse3w(ji+1,jj,jk) * ( rhd(ji+1,jj,jk) + rhd(ji+1,jj,jk-1) + 2*znad )   &
531                  &              - fse3w(ji  ,jj,jk) * ( rhd(ji  ,jj,jk) + rhd(ji  ,jj,jk-1) + 2*znad )  )
532               zhpj(ji,jj,jk) = zhpj(ji,jj,jk-1) + zcoef0 / e2v(ji,jj)   &
533                  &           * (  fse3w(ji,jj+1,jk) * ( rhd(ji,jj+1,jk) + rhd(ji,jj+1,jk-1) + 2*znad )   &
534                  &              - fse3w(ji,jj  ,jk) * ( rhd(ji,jj,  jk) + rhd(ji,jj  ,jk-1) + 2*znad )  )
535               ! s-coordinate pressure gradient correction
536               zuap = -zcoef0 * ( rhd   (ji+1,jj  ,jk) + rhd   (ji,jj,jk) + 2._wp * znad )   &
537                  &           * ( fsde3w(ji+1,jj  ,jk) - fsde3w(ji,jj,jk) ) / e1u(ji,jj)
538               zvap = -zcoef0 * ( rhd   (ji  ,jj+1,jk) + rhd   (ji,jj,jk) + 2._wp * znad )   &
539                  &           * ( fsde3w(ji  ,jj+1,jk) - fsde3w(ji,jj,jk) ) / e2v(ji,jj)
540
541               IF(ln_wd) THEN
542                 zhpi(ji,jj,jk) = zhpi(ji,jj,jk) * zcpx(ji,jj)
543                 zhpj(ji,jj,jk) = zhpj(ji,jj,jk) * zcpy(ji,jj)
544                 zuap = zuap * zcpx(ji,jj)
545                 zvap = zvap * zcpy(ji,jj)
546               ENDIF
547
548               ! add to the general momentum trend
549               ua(ji,jj,jk) = ua(ji,jj,jk) + zhpi(ji,jj,jk) + zuap
550               va(ji,jj,jk) = va(ji,jj,jk) + zhpj(ji,jj,jk) + zvap
551            END DO
552         END DO
553      END DO
554#endif
555      !
556
557      CALL wrk_dealloc( jpi,jpj,jpk, zhpi, zhpj )
558      IF(ln_wd) CALL wrk_dealloc( jpi,jpj, zcpx, zcpy )
559      !
560   END SUBROUTINE hpg_sco
561
562   SUBROUTINE hpg_djc( kt )
563      !!---------------------------------------------------------------------
564      !!                  ***  ROUTINE hpg_djc  ***
565      !!
566      !! ** Method  :   Density Jacobian with Cubic polynomial scheme
567      !!
568      !! Reference: Shchepetkin and McWilliams, J. Geophys. Res., 108(C3), 3090, 2003
569      !!----------------------------------------------------------------------
570      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
571      !!
572      INTEGER  ::   ji, jj, jk          ! dummy loop indices
573      REAL(wp) ::   zcoef0, zep, cffw   ! temporary scalars
574      REAL(wp) ::   z1_10, cffu, cffx   !    "         "
575      REAL(wp) ::   z1_12, cffv, cffy   !    "         "
576      LOGICAL  ::   ll_tmp1, ll_tmp2    ! local logical variables
577      REAL(wp), POINTER, DIMENSION(:,:,:) ::  zhpi, zhpj
578      REAL(wp), POINTER, DIMENSION(:,:,:) ::  dzx, dzy, dzz, dzu, dzv, dzw
579      REAL(wp), POINTER, DIMENSION(:,:,:) ::  drhox, drhoy, drhoz, drhou, drhov, drhow
580      REAL(wp), POINTER, DIMENSION(:,:,:) ::  rho_i, rho_j, rho_k
581      REAL(wp), POINTER, DIMENSION(:,:)   ::  zcpx, zcpy    !W/D pressure filter
582      !!----------------------------------------------------------------------
583      !
584      CALL wrk_alloc( jpi, jpj, jpk, dzx  , dzy  , dzz  , dzu  , dzv  , dzw   )
585      CALL wrk_alloc( jpi, jpj, jpk, drhox, drhoy, drhoz, drhou, drhov, drhow )
586      CALL wrk_alloc( jpi, jpj, jpk, rho_i, rho_j, rho_k,  zhpi,  zhpj        )
587      IF(ln_wd) CALL wrk_alloc( jpi,jpj, zcpx, zcpy )
588      !
589
590      IF( kt == nit000 ) THEN
591         IF(lwp) WRITE(numout,*)
592         IF(lwp) WRITE(numout,*) 'dyn:hpg_djc : hydrostatic pressure gradient trend'
593         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   s-coordinate case, density Jacobian with cubic polynomial scheme'
594      ENDIF
595
596      IF(ln_wd) THEN
597        DO jj = 2, jpjm1
598           DO ji = 2, jpim1
599             ll_tmp1 = MIN(sshn(ji,jj), sshn(ji+1,jj)) > MAX(-bathy(ji,jj), -bathy(ji+1,jj))
600             ll_tmp2 = MAX(sshn(ji,jj), sshn(ji+1,jj)) > MAX(-bathy(ji,jj), -bathy(ji+1,jj)) +&
601                                                       & rn_wdmin1 + rn_wdmin2
602
603             IF(ll_tmp1) THEN
604               zcpx(ji,jj) = 1.0_wp
605             ELSE IF(ll_tmp2) THEN
606               ! no worries about sshn(ji+1,jj)-sshn(ji,jj) = 0, it won't happen ! here
607               zcpx(ji,jj) = ABS((sshn(ji+1,jj) + bathy(ji+1,jj) - sshn(ji,jj) - bathy(ji,jj)) /&
608                           &     (sshn(ji+1,jj) - sshn(ji,jj)))
609             ELSE
610               zcpx(ji,jj) = 0._wp
611             END IF
612     
613             ll_tmp1 = MIN(sshn(ji,jj), sshn(ji,jj+1)) > MAX(-bathy(ji,jj), -bathy(ji,jj+1))
614             ll_tmp2 = MAX(sshn(ji,jj), sshn(ji,jj+1)) > MAX(-bathy(ji,jj), -bathy(ji,jj+1)) +&
615                                                       & rn_wdmin1 + rn_wdmin2
616
617             IF(ll_tmp1) THEN
618               zcpy(ji,jj) = 1.0_wp
619             ELSE IF(ll_tmp2) THEN
620               ! no worries about sshn(ji,jj+1)-sshn(ji,jj) = 0, it won't happen ! here
621               zcpy(ji,jj) = ABS((sshn(ji,jj+1) + bathy(ji,jj+1) - sshn(ji,jj) - bathy(ji,jj)) /&
622                           &     (sshn(ji,jj+1) - sshn(ji,jj)))
623             ELSE
624               zcpy(ji,jj) = 0._wp
625             END IF
626           END DO
627        END DO
628        CALL lbc_lnk( zcpx, 'U', 1._wp )    ;   CALL lbc_lnk( zcpy, 'V', 1._wp )
629      ENDIF
630
631      ! Local constant initialization
632      zcoef0 = - grav * 0.5_wp
633      z1_10  = 1._wp / 10._wp
634      z1_12  = 1._wp / 12._wp
635
636      !----------------------------------------------------------------------------------------
637      !  compute and store in provisional arrays elementary vertical and horizontal differences
638      !----------------------------------------------------------------------------------------
639
640!!bug gm   Not a true bug, but... dzz=e3w  for dzx, dzy verify what it is really
641
642      DO jk = 2, jpkm1
643         DO jj = 2, jpjm1
644            DO ji = fs_2, fs_jpim1   ! vector opt.
645               drhoz(ji,jj,jk) = rhd   (ji  ,jj  ,jk) - rhd   (ji,jj,jk-1)
646               dzz  (ji,jj,jk) = fsde3w(ji  ,jj  ,jk) - fsde3w(ji,jj,jk-1)
647               drhox(ji,jj,jk) = rhd   (ji+1,jj  ,jk) - rhd   (ji,jj,jk  )
648               dzx  (ji,jj,jk) = fsde3w(ji+1,jj  ,jk) - fsde3w(ji,jj,jk  )
649               drhoy(ji,jj,jk) = rhd   (ji  ,jj+1,jk) - rhd   (ji,jj,jk  )
650               dzy  (ji,jj,jk) = fsde3w(ji  ,jj+1,jk) - fsde3w(ji,jj,jk  )
651            END DO
652         END DO
653      END DO
654
655      !-------------------------------------------------------------------------
656      ! compute harmonic averages using eq. 5.18
657      !-------------------------------------------------------------------------
658      zep = 1.e-15
659
660!!bug  gm  drhoz not defined at level 1 and used (jk-1 with jk=2)
661!!bug  gm  idem for drhox, drhoy et ji=jpi and jj=jpj
662
663      DO jk = 2, jpkm1
664         DO jj = 2, jpjm1
665            DO ji = fs_2, fs_jpim1   ! vector opt.
666               cffw = 2._wp * drhoz(ji  ,jj  ,jk) * drhoz(ji,jj,jk-1)
667
668               cffu = 2._wp * drhox(ji+1,jj  ,jk) * drhox(ji,jj,jk  )
669               cffx = 2._wp * dzx  (ji+1,jj  ,jk) * dzx  (ji,jj,jk  )
670
671               cffv = 2._wp * drhoy(ji  ,jj+1,jk) * drhoy(ji,jj,jk  )
672               cffy = 2._wp * dzy  (ji  ,jj+1,jk) * dzy  (ji,jj,jk  )
673
674               IF( cffw > zep) THEN
675                  drhow(ji,jj,jk) = 2._wp *   drhoz(ji,jj,jk) * drhoz(ji,jj,jk-1)   &
676                     &                    / ( drhoz(ji,jj,jk) + drhoz(ji,jj,jk-1) )
677               ELSE
678                  drhow(ji,jj,jk) = 0._wp
679               ENDIF
680
681               dzw(ji,jj,jk) = 2._wp *   dzz(ji,jj,jk) * dzz(ji,jj,jk-1)   &
682                  &                  / ( dzz(ji,jj,jk) + dzz(ji,jj,jk-1) )
683
684               IF( cffu > zep ) THEN
685                  drhou(ji,jj,jk) = 2._wp *   drhox(ji+1,jj,jk) * drhox(ji,jj,jk)   &
686                     &                    / ( drhox(ji+1,jj,jk) + drhox(ji,jj,jk) )
687               ELSE
688                  drhou(ji,jj,jk ) = 0._wp
689               ENDIF
690
691               IF( cffx > zep ) THEN
692                  dzu(ji,jj,jk) = 2._wp *   dzx(ji+1,jj,jk) * dzx(ji,jj,jk)   &
693                     &                  / ( dzx(ji+1,jj,jk) + dzx(ji,jj,jk) )
694               ELSE
695                  dzu(ji,jj,jk) = 0._wp
696               ENDIF
697
698               IF( cffv > zep ) THEN
699                  drhov(ji,jj,jk) = 2._wp *   drhoy(ji,jj+1,jk) * drhoy(ji,jj,jk)   &
700                     &                    / ( drhoy(ji,jj+1,jk) + drhoy(ji,jj,jk) )
701               ELSE
702                  drhov(ji,jj,jk) = 0._wp
703               ENDIF
704
705               IF( cffy > zep ) THEN
706                  dzv(ji,jj,jk) = 2._wp *   dzy(ji,jj+1,jk) * dzy(ji,jj,jk)   &
707                     &                  / ( dzy(ji,jj+1,jk) + dzy(ji,jj,jk) )
708               ELSE
709                  dzv(ji,jj,jk) = 0._wp
710               ENDIF
711
712            END DO
713         END DO
714      END DO
715
716      !----------------------------------------------------------------------------------
717      ! apply boundary conditions at top and bottom using 5.36-5.37
718      !----------------------------------------------------------------------------------
719      drhow(:,:, 1 ) = 1.5_wp * ( drhoz(:,:, 2 ) - drhoz(:,:,  1  ) ) - 0.5_wp * drhow(:,:,  2  )
720      drhou(:,:, 1 ) = 1.5_wp * ( drhox(:,:, 2 ) - drhox(:,:,  1  ) ) - 0.5_wp * drhou(:,:,  2  )
721      drhov(:,:, 1 ) = 1.5_wp * ( drhoy(:,:, 2 ) - drhoy(:,:,  1  ) ) - 0.5_wp * drhov(:,:,  2  )
722
723      drhow(:,:,jpk) = 1.5_wp * ( drhoz(:,:,jpk) - drhoz(:,:,jpkm1) ) - 0.5_wp * drhow(:,:,jpkm1)
724      drhou(:,:,jpk) = 1.5_wp * ( drhox(:,:,jpk) - drhox(:,:,jpkm1) ) - 0.5_wp * drhou(:,:,jpkm1)
725      drhov(:,:,jpk) = 1.5_wp * ( drhoy(:,:,jpk) - drhoy(:,:,jpkm1) ) - 0.5_wp * drhov(:,:,jpkm1)
726
727
728      !--------------------------------------------------------------
729      ! Upper half of top-most grid box, compute and store
730      !-------------------------------------------------------------
731
732!!bug gm   :  e3w-de3w = 0.5*e3w  ....  and de3w(2)-de3w(1)=e3w(2) ....   to be verified
733!          true if de3w is really defined as the sum of the e3w scale factors as, it seems to me, it should be
734
735      DO jj = 2, jpjm1
736         DO ji = fs_2, fs_jpim1   ! vector opt.
737            rho_k(ji,jj,1) = -grav * ( fse3w(ji,jj,1) - fsde3w(ji,jj,1) )               &
738               &                   * (  rhd(ji,jj,1)                                    &
739               &                     + 0.5_wp * ( rhd(ji,jj,2) - rhd(ji,jj,1) )         &
740               &                              * ( fse3w (ji,jj,1) - fsde3w(ji,jj,1) )   &
741               &                              / ( fsde3w(ji,jj,2) - fsde3w(ji,jj,1) )  )
742         END DO
743      END DO
744
745!!bug gm    : here also, simplification is possible
746!!bug gm    : optimisation: 1/10 and 1/12 the division should be done before the loop
747
748      DO jk = 2, jpkm1
749         DO jj = 2, jpjm1
750            DO ji = fs_2, fs_jpim1   ! vector opt.
751
752               rho_k(ji,jj,jk) = zcoef0 * ( rhd   (ji,jj,jk) + rhd   (ji,jj,jk-1) )                                   &
753                  &                     * ( fsde3w(ji,jj,jk) - fsde3w(ji,jj,jk-1) )                                   &
754                  &            - grav * z1_10 * (                                                                     &
755                  &     ( drhow (ji,jj,jk) - drhow (ji,jj,jk-1) )                                                     &
756                  &   * ( fsde3w(ji,jj,jk) - fsde3w(ji,jj,jk-1) - z1_12 * ( dzw  (ji,jj,jk) + dzw  (ji,jj,jk-1) ) )   &
757                  &   - ( dzw   (ji,jj,jk) - dzw   (ji,jj,jk-1) )                                                     &
758                  &   * ( rhd   (ji,jj,jk) - rhd   (ji,jj,jk-1) - z1_12 * ( drhow(ji,jj,jk) + drhow(ji,jj,jk-1) ) )   &
759                  &                             )
760
761               rho_i(ji,jj,jk) = zcoef0 * ( rhd   (ji+1,jj,jk) + rhd   (ji,jj,jk) )                                   &
762                  &                     * ( fsde3w(ji+1,jj,jk) - fsde3w(ji,jj,jk) )                                   &
763                  &            - grav* z1_10 * (                                                                      &
764                  &     ( drhou (ji+1,jj,jk) - drhou (ji,jj,jk) )                                                     &
765                  &   * ( fsde3w(ji+1,jj,jk) - fsde3w(ji,jj,jk) - z1_12 * ( dzu  (ji+1,jj,jk) + dzu  (ji,jj,jk) ) )   &
766                  &   - ( dzu   (ji+1,jj,jk) - dzu   (ji,jj,jk) )                                                     &
767                  &   * ( rhd   (ji+1,jj,jk) - rhd   (ji,jj,jk) - z1_12 * ( drhou(ji+1,jj,jk) + drhou(ji,jj,jk) ) )   &
768                  &                            )
769
770               rho_j(ji,jj,jk) = zcoef0 * ( rhd   (ji,jj+1,jk) + rhd   (ji,jj,jk) )                                   &
771                  &                     * ( fsde3w(ji,jj+1,jk) - fsde3w(ji,jj,jk) )                                   &
772                  &            - grav* z1_10 * (                                                                      &
773                  &     ( drhov (ji,jj+1,jk) - drhov (ji,jj,jk) )                                                     &
774                  &   * ( fsde3w(ji,jj+1,jk) - fsde3w(ji,jj,jk) - z1_12 * ( dzv  (ji,jj+1,jk) + dzv  (ji,jj,jk) ) )   &
775                  &   - ( dzv   (ji,jj+1,jk) - dzv   (ji,jj,jk) )                                                     &
776                  &   * ( rhd   (ji,jj+1,jk) - rhd   (ji,jj,jk) - z1_12 * ( drhov(ji,jj+1,jk) + drhov(ji,jj,jk) ) )   &
777                  &                            )
778
779            END DO
780         END DO
781      END DO
782      CALL lbc_lnk(rho_k,'W',1.)
783      CALL lbc_lnk(rho_i,'U',1.)
784      CALL lbc_lnk(rho_j,'V',1.)
785
786
787      ! ---------------
788      !  Surface value
789      ! ---------------
790      DO jj = 2, jpjm1
791         DO ji = fs_2, fs_jpim1   ! vector opt.
792            zhpi(ji,jj,1) = ( rho_k(ji+1,jj  ,1) - rho_k(ji,jj,1) - rho_i(ji,jj,1) ) / e1u(ji,jj)
793            zhpj(ji,jj,1) = ( rho_k(ji  ,jj+1,1) - rho_k(ji,jj,1) - rho_j(ji,jj,1) ) / e2v(ji,jj)
794            IF(ln_wd) THEN
795              zhpi(ji,jj,1) = zhpi(ji,jj,1) * zcpx(ji,jj)
796              zhpj(ji,jj,1) = zhpj(ji,jj,1) * zcpy(ji,jj)
797            ENDIF
798            ! add to the general momentum trend
799            ua(ji,jj,1) = ua(ji,jj,1) + zhpi(ji,jj,1)
800            va(ji,jj,1) = va(ji,jj,1) + zhpj(ji,jj,1)
801         END DO
802      END DO
803
804      ! ----------------
805      !  interior value   (2=<jk=<jpkm1)
806      ! ----------------
807      DO jk = 2, jpkm1
808         DO jj = 2, jpjm1
809            DO ji = fs_2, fs_jpim1   ! vector opt.
810               ! hydrostatic pressure gradient along s-surfaces
811               zhpi(ji,jj,jk) = zhpi(ji,jj,jk-1)                                &
812                  &           + (  ( rho_k(ji+1,jj,jk) - rho_k(ji,jj,jk  ) )    &
813                  &              - ( rho_i(ji  ,jj,jk) - rho_i(ji,jj,jk-1) )  ) / e1u(ji,jj)
814               zhpj(ji,jj,jk) = zhpj(ji,jj,jk-1)                                &
815                  &           + (  ( rho_k(ji,jj+1,jk) - rho_k(ji,jj,jk  ) )    &
816                  &               -( rho_j(ji,jj  ,jk) - rho_j(ji,jj,jk-1) )  ) / e2v(ji,jj)
817               IF(ln_wd) THEN
818                 zhpi(ji,jj,jk) = zhpi(ji,jj,jk) * zcpx(ji,jj)
819                 zhpj(ji,jj,jk) = zhpj(ji,jj,jk) * zcpy(ji,jj)
820               ENDIF
821               ! add to the general momentum trend
822               ua(ji,jj,jk) = ua(ji,jj,jk) + zhpi(ji,jj,jk)
823               va(ji,jj,jk) = va(ji,jj,jk) + zhpj(ji,jj,jk)
824            END DO
825         END DO
826      END DO
827      !
828      CALL wrk_dealloc( jpi, jpj, jpk, dzx  , dzy  , dzz  , dzu  , dzv  , dzw   )
829      CALL wrk_dealloc( jpi, jpj, jpk, drhox, drhoy, drhoz, drhou, drhov, drhow )
830      CALL wrk_dealloc( jpi, jpj, jpk, rho_i, rho_j, rho_k,  zhpi,  zhpj        )
831      IF(ln_wd) CALL wrk_dealloc( jpi,jpj, zcpx, zcpy )
832      !
833   END SUBROUTINE hpg_djc
834
835
836   SUBROUTINE hpg_prj( kt )
837      !!---------------------------------------------------------------------
838      !!                  ***  ROUTINE hpg_prj  ***
839      !!
840      !! ** Method  :   s-coordinate case.
841      !!      A Pressure-Jacobian horizontal pressure gradient method
842      !!      based on the constrained cubic-spline interpolation for
843      !!      all vertical coordinate systems
844      !!
845      !! ** Action : - Update (ua,va) with the now hydrastatic pressure trend
846      !!             - Save the trend (l_trddyn=T)
847      !!
848      !!----------------------------------------------------------------------
849      INTEGER, PARAMETER  :: polynomial_type = 1    ! 1: cubic spline, 2: linear
850      INTEGER, INTENT(in) ::   kt                   ! ocean time-step index
851      !!
852      INTEGER  ::   ji, jj, jk, jkk                 ! dummy loop indices
853      REAL(wp) ::   zcoef0, znad                    ! temporary scalars
854      !!
855      !! The local variables for the correction term
856      INTEGER  :: jk1, jis, jid, jjs, jjd
857      LOGICAL  :: ll_tmp1, ll_tmp2                  ! local logical variables
858      REAL(wp) :: zuijk, zvijk, zpwes, zpwed, zpnss, zpnsd, zdeps
859      REAL(wp) :: zrhdt1
860      REAL(wp) :: zdpdx1, zdpdx2, zdpdy1, zdpdy2
861      REAL(wp), POINTER, DIMENSION(:,:,:) ::   zdept, zrhh
862      REAL(wp), POINTER, DIMENSION(:,:,:) ::   zhpi, zu, zv, fsp, xsp, asp, bsp, csp, dsp
863      REAL(wp), POINTER, DIMENSION(:,:)   ::   sshu_n, sshv_n
864      REAL(wp), POINTER, DIMENSION(:,:)   ::   zcpx, zcpy    !W/D pressure filter
865      !!----------------------------------------------------------------------
866      !
867      CALL wrk_alloc( jpi,jpj,jpk, zhpi, zu, zv, fsp, xsp, asp, bsp, csp, dsp )
868      CALL wrk_alloc( jpi,jpj,jpk, zdept, zrhh )
869      CALL wrk_alloc( jpi,jpj, sshu_n, sshv_n )
870      IF(ln_wd) CALL wrk_alloc( jpi,jpj, zcpx, zcpy )
871      !
872      IF( kt == nit000 ) THEN
873         IF(lwp) WRITE(numout,*)
874         IF(lwp) WRITE(numout,*) 'dyn:hpg_prj : hydrostatic pressure gradient trend'
875         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   s-coordinate case, cubic spline pressure Jacobian'
876      ENDIF
877
878      !!----------------------------------------------------------------------
879      ! Local constant initialization
880      zcoef0 = - grav
881      znad = 0.0_wp
882      IF( lk_vvl ) znad = 1._wp
883
884      IF(ln_wd) THEN
885        DO jj = 2, jpjm1
886           DO ji = 2, jpim1
887             ll_tmp1 = MIN(sshn(ji,jj), sshn(ji+1,jj)) > MAX(-bathy(ji,jj), -bathy(ji+1,jj))
888             ll_tmp2 = MAX(sshn(ji,jj), sshn(ji+1,jj)) > MAX(-bathy(ji,jj), -bathy(ji+1,jj)) +&
889                                                       & rn_wdmin1 + rn_wdmin2
890
891             IF(ll_tmp1) THEN
892               zcpx(ji,jj) = 1.0_wp
893             ELSE IF(ll_tmp2) THEN
894               ! no worries about sshn(ji+1,jj)-sshn(ji,jj) = 0, it won't happen ! here
895               zcpx(ji,jj) = ABS((sshn(ji+1,jj) + bathy(ji+1,jj) - sshn(ji,jj) - bathy(ji,jj)) /&
896                           &     (sshn(ji+1,jj) - sshn(ji,jj)))
897             ELSE
898               zcpx(ji,jj) = 0._wp
899             END IF
900     
901             ll_tmp1 = MIN(sshn(ji,jj), sshn(ji,jj+1)) > MAX(-bathy(ji,jj), -bathy(ji,jj+1))
902             ll_tmp2 = MAX(sshn(ji,jj), sshn(ji,jj+1)) > MAX(-bathy(ji,jj), -bathy(ji,jj+1)) +&
903                                                       & rn_wdmin1 + rn_wdmin2
904
905             IF(ll_tmp1.OR.ll_tmp2) THEN
906               zcpy(ji,jj) = 1.0_wp
907             ELSE IF(ll_tmp2) THEN
908               ! no worries about sshn(ji,jj+1)-sshn(ji,jj) = 0, it won't happen ! here
909               zcpy(ji,jj) = ABS((sshn(ji,jj+1) + bathy(ji,jj+1) - sshn(ji,jj) - bathy(ji,jj)) /&
910                           &     (sshn(ji,jj+1) - sshn(ji,jj)))
911             ELSE
912               zcpy(ji,jj) = 0._wp
913             END IF
914           END DO
915        END DO
916        CALL lbc_lnk( zcpx, 'U', 1._wp )    ;   CALL lbc_lnk( zcpy, 'V', 1._wp )
917      ENDIF
918
919      ! Clean 3-D work arrays
920      zhpi(:,:,:) = 0._wp
921      zrhh(:,:,:) = rhd(:,:,:)
922
923      ! Preparing vertical density profile "zrhh(:,:,:)" for hybrid-sco coordinate
924      DO jj = 1, jpj
925        DO ji = 1, jpi
926          jk = mbathy(ji,jj)
927          IF( jk <= 0 ) THEN; zrhh(ji,jj,:) = 0._wp
928          ELSE IF(jk == 1) THEN; zrhh(ji,jj, jk+1:jpk) = rhd(ji,jj,jk)
929          ELSE IF(jk < jpkm1) THEN
930             DO jkk = jk+1, jpk
931                zrhh(ji,jj,jkk) = interp1(fsde3w(ji,jj,jkk),   fsde3w(ji,jj,jkk-1), &
932                                         fsde3w(ji,jj,jkk-2), rhd(ji,jj,jkk-1), rhd(ji,jj,jkk-2))
933             END DO
934          ENDIF
935        END DO
936      END DO
937
938      ! Transfer the depth of "T(:,:,:)" to vertical coordinate "zdept(:,:,:)"
939      DO jj = 1, jpj;   DO ji = 1, jpi
940          zdept(ji,jj,1) = 0.5_wp * fse3w(ji,jj,1) - sshn(ji,jj) * znad
941      END DO        ;   END DO
942
943      DO jk = 2, jpk;   DO jj = 1, jpj;   DO ji = 1, jpi
944          zdept(ji,jj,jk) = zdept(ji,jj,jk-1) + fse3w(ji,jj,jk)
945      END DO        ;   END DO        ;   END DO
946
947      fsp(:,:,:) = zrhh(:,:,:)
948      xsp(:,:,:) = zdept(:,:,:)
949
950      ! Construct the vertical density profile with the
951      ! constrained cubic spline interpolation
952      ! rho(z) = asp + bsp*z + csp*z^2 + dsp*z^3
953      CALL cspline(fsp,xsp,asp,bsp,csp,dsp,polynomial_type)
954
955      ! Integrate the hydrostatic pressure "zhpi(:,:,:)" at "T(ji,jj,1)"
956      DO jj = 2, jpj
957        DO ji = 2, jpi
958          zrhdt1 = zrhh(ji,jj,1) - interp3(zdept(ji,jj,1),asp(ji,jj,1), &
959                                         bsp(ji,jj,1),   csp(ji,jj,1), &
960                                         dsp(ji,jj,1) ) * 0.25_wp * fse3w(ji,jj,1)
961
962          ! assuming linear profile across the top half surface layer
963          zhpi(ji,jj,1) =  0.5_wp * fse3w(ji,jj,1) * zrhdt1
964        END DO
965      END DO
966
967      ! Calculate the pressure "zhpi(:,:,:)" at "T(ji,jj,2:jpkm1)"
968      DO jk = 2, jpkm1
969        DO jj = 2, jpj
970          DO ji = 2, jpi
971            zhpi(ji,jj,jk) = zhpi(ji,jj,jk-1) +                          &
972                             integ_spline(zdept(ji,jj,jk-1), zdept(ji,jj,jk),&
973                                    asp(ji,jj,jk-1),    bsp(ji,jj,jk-1), &
974                                    csp(ji,jj,jk-1),    dsp(ji,jj,jk-1))
975          END DO
976        END DO
977      END DO
978
979      ! Z coordinate of U(ji,jj,1:jpkm1) and V(ji,jj,1:jpkm1)
980
981      ! The following modification "sshu_n -> sshn" is a big mistake, this
982      ! should never happen here, remember to correct this in NEMO v3.6
983      ! trunk.  H.L.
984
985      !prepare sshu_n and sshv_n
986      DO jj = 1, jpjm1
987        DO ji = 1, jpim1
988          sshu_n(ji,jj) = (e12u(ji,jj) * sshn(ji,jj) + e12u(ji+1, jj) * sshn(ji+1,jj)) * &
989                        & r1_e12u(ji,jj) * umask(ji,jj,1) * 0.5_wp
990          sshv_n(ji,jj) = (e12v(ji,jj) * sshn(ji,jj) + e12v(ji+1, jj) * sshn(ji,jj+1)) * &
991                        & r1_e12v(ji,jj) * vmask(ji,jj,1) * 0.5_wp
992        END DO
993      END DO
994
995      CALL lbc_lnk (sshu_n, 'U', 1.)
996      CALL lbc_lnk (sshv_n, 'V', 1.)
997
998      DO jj = 2, jpjm1
999        DO ji = 2, jpim1
1000!          zu(ji,jj,1) = - ( fse3u(ji,jj,1) - sshn(ji,jj) * znad)    ! probable bug: changed from sshu_n for ztilde compilation
1001!          zv(ji,jj,1) = - ( fse3v(ji,jj,1) - sshn(ji,jj) * znad)    ! probable bug: changed from sshv_n for ztilde compilation
1002          zu(ji,jj,1) = - ( fse3u(ji,jj,1) - sshu_n(ji,jj) * znad)
1003          zv(ji,jj,1) = - ( fse3v(ji,jj,1) - sshv_n(ji,jj) * znad)
1004        END DO
1005      END DO
1006
1007      DO jk = 2, jpkm1
1008        DO jj = 2, jpjm1
1009          DO ji = 2, jpim1
1010            zu(ji,jj,jk) = zu(ji,jj,jk-1)- fse3u(ji,jj,jk)
1011            zv(ji,jj,jk) = zv(ji,jj,jk-1)- fse3v(ji,jj,jk)
1012          END DO
1013        END DO
1014      END DO
1015
1016      DO jk = 1, jpkm1
1017        DO jj = 2, jpjm1
1018          DO ji = 2, jpim1
1019            zu(ji,jj,jk) = zu(ji,jj,jk) + 0.5_wp * fse3u(ji,jj,jk)
1020            zv(ji,jj,jk) = zv(ji,jj,jk) + 0.5_wp * fse3v(ji,jj,jk)
1021          END DO
1022        END DO
1023      END DO
1024
1025      DO jk = 1, jpkm1
1026        DO jj = 2, jpjm1
1027          DO ji = 2, jpim1
1028            zu(ji,jj,jk) = min(zu(ji,jj,jk), max(-zdept(ji,jj,jk), -zdept(ji+1,jj,jk)))
1029            zu(ji,jj,jk) = max(zu(ji,jj,jk), min(-zdept(ji,jj,jk), -zdept(ji+1,jj,jk)))
1030            zv(ji,jj,jk) = min(zv(ji,jj,jk), max(-zdept(ji,jj,jk), -zdept(ji,jj+1,jk)))
1031            zv(ji,jj,jk) = max(zv(ji,jj,jk), min(-zdept(ji,jj,jk), -zdept(ji,jj+1,jk)))
1032          END DO
1033        END DO
1034      END DO
1035
1036
1037      DO jk = 1, jpkm1
1038        DO jj = 2, jpjm1
1039          DO ji = 2, jpim1
1040            zpwes = 0._wp; zpwed = 0._wp
1041            zpnss = 0._wp; zpnsd = 0._wp
1042            zuijk = zu(ji,jj,jk)
1043            zvijk = zv(ji,jj,jk)
1044
1045            !!!!!     for u equation
1046            IF( jk <= mbku(ji,jj) ) THEN
1047               IF( -zdept(ji+1,jj,jk) >= -zdept(ji,jj,jk) ) THEN
1048                 jis = ji + 1; jid = ji
1049               ELSE
1050                 jis = ji;     jid = ji +1
1051               ENDIF
1052
1053               ! integrate the pressure on the shallow side
1054               jk1 = jk
1055               DO WHILE ( -zdept(jis,jj,jk1) > zuijk )
1056                 IF( jk1 == mbku(ji,jj) ) THEN
1057                   zuijk = -zdept(jis,jj,jk1)
1058                   EXIT
1059                 ENDIF
1060                 zdeps = MIN(zdept(jis,jj,jk1+1), -zuijk)
1061                 zpwes = zpwes +                                    &
1062                      integ_spline(zdept(jis,jj,jk1), zdeps,            &
1063                             asp(jis,jj,jk1),    bsp(jis,jj,jk1), &
1064                             csp(jis,jj,jk1),    dsp(jis,jj,jk1))
1065                 jk1 = jk1 + 1
1066               END DO
1067
1068               ! integrate the pressure on the deep side
1069               jk1 = jk
1070               DO WHILE ( -zdept(jid,jj,jk1) < zuijk )
1071                 IF( jk1 == 1 ) THEN
1072                   zdeps = zdept(jid,jj,1) + MIN(zuijk, sshn(jid,jj)*znad)
1073                   zrhdt1 = zrhh(jid,jj,1) - interp3(zdept(jid,jj,1), asp(jid,jj,1), &
1074                                                     bsp(jid,jj,1),   csp(jid,jj,1), &
1075                                                     dsp(jid,jj,1)) * zdeps
1076                   zpwed  = zpwed + 0.5_wp * (zrhh(jid,jj,1) + zrhdt1) * zdeps
1077                   EXIT
1078                 ENDIF
1079                 zdeps = MAX(zdept(jid,jj,jk1-1), -zuijk)
1080                 zpwed = zpwed +                                        &
1081                        integ_spline(zdeps,              zdept(jid,jj,jk1), &
1082                               asp(jid,jj,jk1-1), bsp(jid,jj,jk1-1),  &
1083                               csp(jid,jj,jk1-1), dsp(jid,jj,jk1-1) )
1084                 jk1 = jk1 - 1
1085               END DO
1086
1087               ! update the momentum trends in u direction
1088
1089               zdpdx1 = zcoef0 / e1u(ji,jj) * (zhpi(ji+1,jj,jk) - zhpi(ji,jj,jk))
1090               IF( lk_vvl ) THEN
1091                 zdpdx2 = zcoef0 / e1u(ji,jj) * &
1092                         ( REAL(jis-jid, wp) * (zpwes + zpwed) + (sshn(ji+1,jj)-sshn(ji,jj)) )
1093                ELSE
1094                 zdpdx2 = zcoef0 / e1u(ji,jj) * REAL(jis-jid, wp) * (zpwes + zpwed)
1095               ENDIF
1096
1097               IF(ln_wd) THEN
1098                  zdpdx1 = zdpdx1 * zcpx(ji,jj)
1099                  zdpdx2 = zdpdx2 * zcpx(ji,jj)
1100                ENDIF
1101                ua(ji,jj,jk) = ua(ji,jj,jk) + (zdpdx1 + zdpdx2) * umask(ji,jj,jk) 
1102            ENDIF
1103
1104            !!!!!     for v equation
1105            IF( jk <= mbkv(ji,jj) ) THEN
1106               IF( -zdept(ji,jj+1,jk) >= -zdept(ji,jj,jk) ) THEN
1107                 jjs = jj + 1; jjd = jj
1108               ELSE
1109                 jjs = jj    ; jjd = jj + 1
1110               ENDIF
1111
1112               ! integrate the pressure on the shallow side
1113               jk1 = jk
1114               DO WHILE ( -zdept(ji,jjs,jk1) > zvijk )
1115                 IF( jk1 == mbkv(ji,jj) ) THEN
1116                   zvijk = -zdept(ji,jjs,jk1)
1117                   EXIT
1118                 ENDIF
1119                 zdeps = MIN(zdept(ji,jjs,jk1+1), -zvijk)
1120                 zpnss = zpnss +                                      &
1121                        integ_spline(zdept(ji,jjs,jk1), zdeps,            &
1122                               asp(ji,jjs,jk1),    bsp(ji,jjs,jk1), &
1123                               csp(ji,jjs,jk1),    dsp(ji,jjs,jk1) )
1124                 jk1 = jk1 + 1
1125               END DO
1126
1127               ! integrate the pressure on the deep side
1128               jk1 = jk
1129               DO WHILE ( -zdept(ji,jjd,jk1) < zvijk )
1130                 IF( jk1 == 1 ) THEN
1131                   zdeps = zdept(ji,jjd,1) + MIN(zvijk, sshn(ji,jjd)*znad)
1132                   zrhdt1 = zrhh(ji,jjd,1) - interp3(zdept(ji,jjd,1), asp(ji,jjd,1), &
1133                                                     bsp(ji,jjd,1),   csp(ji,jjd,1), &
1134                                                     dsp(ji,jjd,1) ) * zdeps
1135                   zpnsd  = zpnsd + 0.5_wp * (zrhh(ji,jjd,1) + zrhdt1) * zdeps
1136                   EXIT
1137                 ENDIF
1138                 zdeps = MAX(zdept(ji,jjd,jk1-1), -zvijk)
1139                 zpnsd = zpnsd +                                        &
1140                        integ_spline(zdeps,              zdept(ji,jjd,jk1), &
1141                               asp(ji,jjd,jk1-1), bsp(ji,jjd,jk1-1), &
1142                               csp(ji,jjd,jk1-1), dsp(ji,jjd,jk1-1) )
1143                 jk1 = jk1 - 1
1144               END DO
1145
1146
1147               ! update the momentum trends in v direction
1148
1149               zdpdy1 = zcoef0 / e2v(ji,jj) * (zhpi(ji,jj+1,jk) - zhpi(ji,jj,jk))
1150               IF( lk_vvl ) THEN
1151                   zdpdy2 = zcoef0 / e2v(ji,jj) * &
1152                           ( REAL(jjs-jjd, wp) * (zpnss + zpnsd) + (sshn(ji,jj+1)-sshn(ji,jj)) )
1153               ELSE
1154                   zdpdy2 = zcoef0 / e2v(ji,jj) * REAL(jjs-jjd, wp) * (zpnss + zpnsd )
1155               ENDIF
1156
1157               IF(ln_wd) THEN
1158                  zdpdy1 = zdpdy1 * zcpy(ji,jj)
1159                  zdpdy2 = zdpdy2 * zcpy(ji,jj)
1160               ENDIF
1161
1162               va(ji,jj,jk) = va(ji,jj,jk) + (zdpdy1 + zdpdy2) * vmask(ji,jj,jk)
1163            ENDIF
1164
1165
1166           END DO
1167        END DO
1168      END DO
1169      !
1170      CALL wrk_dealloc( jpi,jpj,jpk, zhpi, zu, zv, fsp, xsp, asp, bsp, csp, dsp )
1171      CALL wrk_dealloc( jpi,jpj,jpk, zdept, zrhh )
1172      CALL wrk_dealloc( jpi,jpj, sshu_n, sshv_n )
1173      IF(ln_wd) CALL wrk_dealloc( jpi,jpj, zcpx, zcpy )
1174      !
1175   END SUBROUTINE hpg_prj
1176
1177   SUBROUTINE cspline(fsp, xsp, asp, bsp, csp, dsp, polynomial_type)
1178      !!----------------------------------------------------------------------
1179      !!                 ***  ROUTINE cspline  ***
1180      !!
1181      !! ** Purpose :   constrained cubic spline interpolation
1182      !!
1183      !! ** Method  :   f(x) = asp + bsp*x + csp*x^2 + dsp*x^3
1184      !! Reference: CJC Kruger, Constrained Cubic Spline Interpoltation
1185      !!
1186      !!----------------------------------------------------------------------
1187      IMPLICIT NONE
1188      REAL(wp), DIMENSION(:,:,:), INTENT(in)  :: fsp, xsp           ! value and coordinate
1189      REAL(wp), DIMENSION(:,:,:), INTENT(out) :: asp, bsp, csp, dsp ! coefficients of
1190                                                                    ! the interpoated function
1191      INTEGER, INTENT(in) :: polynomial_type                        ! 1: cubic spline
1192                                                                    ! 2: Linear
1193
1194      ! Local Variables
1195      INTEGER  ::   ji, jj, jk                 ! dummy loop indices
1196      INTEGER  ::   jpi, jpj, jpkm1
1197      REAL(wp) ::   zdf1, zdf2, zddf1, zddf2, ztmp1, ztmp2, zdxtmp
1198      REAL(wp) ::   zdxtmp1, zdxtmp2, zalpha
1199      REAL(wp) ::   zdf(size(fsp,3))
1200      !!----------------------------------------------------------------------
1201
1202      jpi   = size(fsp,1)
1203      jpj   = size(fsp,2)
1204      jpkm1 = size(fsp,3) - 1
1205
1206
1207      IF (polynomial_type == 1) THEN     ! Constrained Cubic Spline
1208         DO ji = 1, jpi
1209            DO jj = 1, jpj
1210           !!Fritsch&Butland's method, 1984 (preferred, but more computation)
1211           !    DO jk = 2, jpkm1-1
1212           !       zdxtmp1 = xsp(ji,jj,jk)   - xsp(ji,jj,jk-1)
1213           !       zdxtmp2 = xsp(ji,jj,jk+1) - xsp(ji,jj,jk)
1214           !       zdf1    = ( fsp(ji,jj,jk)   - fsp(ji,jj,jk-1) ) / zdxtmp1
1215           !       zdf2    = ( fsp(ji,jj,jk+1) - fsp(ji,jj,jk)   ) / zdxtmp2
1216           !
1217           !       zalpha = ( zdxtmp1 + 2._wp * zdxtmp2 ) / ( zdxtmp1 + zdxtmp2 ) / 3._wp
1218           !
1219           !       IF(zdf1 * zdf2 <= 0._wp) THEN
1220           !           zdf(jk) = 0._wp
1221           !       ELSE
1222           !         zdf(jk) = zdf1 * zdf2 / ( ( 1._wp - zalpha ) * zdf1 + zalpha * zdf2 )
1223           !       ENDIF
1224           !    END DO
1225
1226           !!Simply geometric average
1227               DO jk = 2, jpkm1-1
1228                  zdf1 = (fsp(ji,jj,jk) - fsp(ji,jj,jk-1)) / (xsp(ji,jj,jk) - xsp(ji,jj,jk-1))
1229                  zdf2 = (fsp(ji,jj,jk+1) - fsp(ji,jj,jk)) / (xsp(ji,jj,jk+1) - xsp(ji,jj,jk))
1230
1231                  IF(zdf1 * zdf2 <= 0._wp) THEN
1232                     zdf(jk) = 0._wp
1233                  ELSE
1234                     zdf(jk) = 2._wp * zdf1 * zdf2 / (zdf1 + zdf2)
1235                  ENDIF
1236               END DO
1237
1238               zdf(1)     = 1.5_wp * ( fsp(ji,jj,2) - fsp(ji,jj,1) ) / &
1239                          &          ( xsp(ji,jj,2) - xsp(ji,jj,1) ) -  0.5_wp * zdf(2)
1240               zdf(jpkm1) = 1.5_wp * ( fsp(ji,jj,jpkm1) - fsp(ji,jj,jpkm1-1) ) / &
1241                          &          ( xsp(ji,jj,jpkm1) - xsp(ji,jj,jpkm1-1) ) - &
1242                          & 0.5_wp * zdf(jpkm1 - 1)
1243
1244               DO jk = 1, jpkm1 - 1
1245                 zdxtmp = xsp(ji,jj,jk+1) - xsp(ji,jj,jk)
1246                 ztmp1  = (zdf(jk+1) + 2._wp * zdf(jk)) / zdxtmp
1247                 ztmp2  =  6._wp * (fsp(ji,jj,jk+1) - fsp(ji,jj,jk)) / zdxtmp / zdxtmp
1248                 zddf1  = -2._wp * ztmp1 + ztmp2
1249                 ztmp1  = (2._wp * zdf(jk+1) + zdf(jk)) / zdxtmp
1250                 zddf2  =  2._wp * ztmp1 - ztmp2
1251
1252                 dsp(ji,jj,jk) = (zddf2 - zddf1) / 6._wp / zdxtmp
1253                 csp(ji,jj,jk) = ( xsp(ji,jj,jk+1) * zddf1 - xsp(ji,jj,jk)*zddf2 ) / 2._wp / zdxtmp
1254                 bsp(ji,jj,jk) = ( fsp(ji,jj,jk+1) - fsp(ji,jj,jk) ) / zdxtmp - &
1255                               & csp(ji,jj,jk) * ( xsp(ji,jj,jk+1) + xsp(ji,jj,jk) ) - &
1256                               & dsp(ji,jj,jk) * ((xsp(ji,jj,jk+1) + xsp(ji,jj,jk))**2 - &
1257                               &                   xsp(ji,jj,jk+1) * xsp(ji,jj,jk))
1258                 asp(ji,jj,jk) = fsp(ji,jj,jk) - xsp(ji,jj,jk) * (bsp(ji,jj,jk) + &
1259                               &                (xsp(ji,jj,jk) * (csp(ji,jj,jk) + &
1260                               &                 dsp(ji,jj,jk) * xsp(ji,jj,jk))))
1261               END DO
1262            END DO
1263         END DO
1264
1265      ELSE IF (polynomial_type == 2) THEN     ! Linear
1266         DO ji = 1, jpi
1267            DO jj = 1, jpj
1268               DO jk = 1, jpkm1-1
1269                  zdxtmp =xsp(ji,jj,jk+1) - xsp(ji,jj,jk)
1270                  ztmp1 = fsp(ji,jj,jk+1) - fsp(ji,jj,jk)
1271
1272                  dsp(ji,jj,jk) = 0._wp
1273                  csp(ji,jj,jk) = 0._wp
1274                  bsp(ji,jj,jk) = ztmp1 / zdxtmp
1275                  asp(ji,jj,jk) = fsp(ji,jj,jk) - bsp(ji,jj,jk) * xsp(ji,jj,jk)
1276               END DO
1277            END DO
1278         END DO
1279
1280      ELSE
1281           CALL ctl_stop( 'invalid polynomial type in cspline' )
1282      ENDIF
1283
1284
1285   END SUBROUTINE cspline
1286
1287
1288   FUNCTION interp1(x, xl, xr, fl, fr)  RESULT(f)
1289      !!----------------------------------------------------------------------
1290      !!                 ***  ROUTINE interp1  ***
1291      !!
1292      !! ** Purpose :   1-d linear interpolation
1293      !!
1294      !! ** Method  :
1295      !!                interpolation is straight forward
1296      !!                extrapolation is also permitted (no value limit)
1297      !!
1298      !!----------------------------------------------------------------------
1299      IMPLICIT NONE
1300      REAL(wp), INTENT(in) ::  x, xl, xr, fl, fr
1301      REAL(wp)             ::  f ! result of the interpolation (extrapolation)
1302      REAL(wp)             ::  zdeltx
1303      !!----------------------------------------------------------------------
1304
1305      zdeltx = xr - xl
1306      IF(abs(zdeltx) <= 10._wp * EPSILON(x)) THEN
1307        f = 0.5_wp * (fl + fr)
1308      ELSE
1309        f = ( (x - xl ) * fr - ( x - xr ) * fl ) / zdeltx
1310      ENDIF
1311
1312   END FUNCTION interp1
1313
1314   FUNCTION interp2(x, a, b, c, d)  RESULT(f)
1315      !!----------------------------------------------------------------------
1316      !!                 ***  ROUTINE interp1  ***
1317      !!
1318      !! ** Purpose :   1-d constrained cubic spline interpolation
1319      !!
1320      !! ** Method  :  cubic spline interpolation
1321      !!
1322      !!----------------------------------------------------------------------
1323      IMPLICIT NONE
1324      REAL(wp), INTENT(in) ::  x, a, b, c, d
1325      REAL(wp)             ::  f ! value from the interpolation
1326      !!----------------------------------------------------------------------
1327
1328      f = a + x* ( b + x * ( c + d * x ) )
1329
1330   END FUNCTION interp2
1331
1332
1333   FUNCTION interp3(x, a, b, c, d)  RESULT(f)
1334      !!----------------------------------------------------------------------
1335      !!                 ***  ROUTINE interp1  ***
1336      !!
1337      !! ** Purpose :   Calculate the first order of deriavtive of
1338      !!                a cubic spline function y=a+b*x+c*x^2+d*x^3
1339      !!
1340      !! ** Method  :   f=dy/dx=b+2*c*x+3*d*x^2
1341      !!
1342      !!----------------------------------------------------------------------
1343      IMPLICIT NONE
1344      REAL(wp), INTENT(in) ::  x, a, b, c, d
1345      REAL(wp)             ::  f ! value from the interpolation
1346      !!----------------------------------------------------------------------
1347
1348      f = b + x * ( 2._wp * c + 3._wp * d * x)
1349
1350   END FUNCTION interp3
1351
1352
1353   FUNCTION integ_spline(xl, xr, a, b, c, d)  RESULT(f)
1354      !!----------------------------------------------------------------------
1355      !!                 ***  ROUTINE interp1  ***
1356      !!
1357      !! ** Purpose :   1-d constrained cubic spline integration
1358      !!
1359      !! ** Method  :  integrate polynomial a+bx+cx^2+dx^3 from xl to xr
1360      !!
1361      !!----------------------------------------------------------------------
1362      IMPLICIT NONE
1363      REAL(wp), INTENT(in) ::  xl, xr, a, b, c, d
1364      REAL(wp)             ::  za1, za2, za3
1365      REAL(wp)             ::  f                   ! integration result
1366      !!----------------------------------------------------------------------
1367
1368      za1 = 0.5_wp * b
1369      za2 = c / 3.0_wp
1370      za3 = 0.25_wp * d
1371
1372      f  = xr * ( a + xr * ( za1 + xr * ( za2 + za3 * xr ) ) ) - &
1373         & xl * ( a + xl * ( za1 + xl * ( za2 + za3 * xl ) ) )
1374
1375   END FUNCTION integ_spline
1376
1377
1378   !!======================================================================
1379END MODULE dynhpg
1380
Note: See TracBrowser for help on using the repository browser.