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

source: branches/2014/dev_r4650_UKMO2_ice_shelves/NEMOGCM/NEMO/OPA_SRC/DYN/dynhpg.F90 @ 4726

Last change on this file since 4726 was 4726, checked in by mathiot, 10 years ago

ISF branch: change name of 2 variables (icedep => risfdep and lmask => ssmask), cosmetic changes and add ldfslp key

  • Property svn:keywords set to Id
File size: 63.2 KB
RevLine 
[3]1MODULE dynhpg
2   !!======================================================================
3   !!                       ***  MODULE  dynhpg  ***
4   !! Ocean dynamics:  hydrostatic pressure gradient trend
5   !!======================================================================
[2528]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
[3764]13   !!                 !         Original code for hpg_ctl, hpg_hel hpg_wdj, hpg_djc, hpg_rot
[2528]14   !!             -   !  2005-11  (G. Madec) style & small optimisation
15   !!            3.3  !  2010-10  (C. Ethe, G. Madec) reorganisation of initialisation phase
[3294]16   !!            3.4  !  2011-11  (H. Liu) hpg_prj: Original code for s-coordinates
17   !!                 !           (A. Coward) suppression of hel, wdj and rot options
[503]18   !!----------------------------------------------------------------------
[3]19
20   !!----------------------------------------------------------------------
[455]21   !!   dyn_hpg      : update the momentum trend with the now horizontal
[3]22   !!                  gradient of the hydrostatic pressure
[2528]23   !!   dyn_hpg_init : initialisation and control of options
[455]24   !!       hpg_zco  : z-coordinate scheme
25   !!       hpg_zps  : z-coordinate plus partial steps (interpolation)
26   !!       hpg_sco  : s-coordinate (standard jacobian formulation)
27   !!       hpg_djc  : s-coordinate (Density Jacobian with Cubic polynomial)
[3294]28   !!       hpg_prj  : s-coordinate (Pressure Jacobian with Cubic polynomial)
[3]29   !!----------------------------------------------------------------------
30   USE oce             ! ocean dynamics and tracers
[4724]31   USE sbc_oce         ! surface variable (only for the flag with ice shelf)
[3]32   USE dom_oce         ! ocean space and time domain
33   USE phycst          ! physical constants
[3764]34   USE trdmod          ! ocean dynamics trends
[216]35   USE trdmod_oce      ! ocean variables trends
[2715]36   USE in_out_manager  ! I/O manager
[258]37   USE prtctl          ! Print control
[3764]38   USE lbclnk          ! lateral boundary condition
[2715]39   USE lib_mpp         ! MPP library
[4666]40   USE eosbn2          ! compute density
[3294]41   USE wrk_nemo        ! Memory Allocation
42   USE timing          ! Timing
[3]43
44   IMPLICIT NONE
45   PRIVATE
46
[2528]47   PUBLIC   dyn_hpg        ! routine called by step module
48   PUBLIC   dyn_hpg_init   ! routine called by opa module
[3]49
[4147]50   !                                    !!* Namelist namdyn_hpg : hydrostatic pressure gradient
51   LOGICAL , PUBLIC ::   ln_hpg_zco      !: z-coordinate - full steps
52   LOGICAL , PUBLIC ::   ln_hpg_zps      !: z-coordinate - partial steps (interpolation)
53   LOGICAL , PUBLIC ::   ln_hpg_sco      !: s-coordinate (standard jacobian formulation)
54   LOGICAL , PUBLIC ::   ln_hpg_djc      !: s-coordinate (Density Jacobian with Cubic polynomial)
55   LOGICAL , PUBLIC ::   ln_hpg_prj      !: s-coordinate (Pressure Jacobian scheme)
56   LOGICAL , PUBLIC ::   ln_dynhpg_imp   !: semi-implicite hpg flag
[455]57
[3764]58   INTEGER , PUBLIC ::   nhpg  =  0   ! = 0 to 7, type of pressure gradient scheme used ! (deduced from ln_hpg_... flags) (PUBLIC for TAM)
[455]59
[3]60   !! * Substitutions
61#  include "domzgr_substitute.h90"
62#  include "vectopt_loop_substitute.h90"
63   !!----------------------------------------------------------------------
[2528]64   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
65   !! $Id$
66   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
[3]67   !!----------------------------------------------------------------------
68CONTAINS
69
70   SUBROUTINE dyn_hpg( kt )
71      !!---------------------------------------------------------------------
72      !!                  ***  ROUTINE dyn_hpg  ***
73      !!
[3764]74      !! ** Method  :   Call the hydrostatic pressure gradient routine
[503]75      !!              using the scheme defined in the namelist
[3764]76      !!
[455]77      !! ** Action : - Update (ua,va) with the now hydrastatic pressure trend
78      !!             - Save the trend (l_trddyn=T)
[503]79      !!----------------------------------------------------------------------
80      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
[3294]81      REAL(wp), POINTER, DIMENSION(:,:,:) ::  ztrdu, ztrdv
[455]82      !!----------------------------------------------------------------------
[2528]83      !
[3294]84      IF( nn_timing == 1 )  CALL timing_start('dyn_hpg')
[2715]85      !
[503]86      IF( l_trddyn ) THEN                    ! Temporary saving of ua and va trends (l_trddyn)
[3294]87         CALL wrk_alloc( jpi,jpj,jpk, ztrdu, ztrdv )
[3764]88         ztrdu(:,:,:) = ua(:,:,:)
89         ztrdv(:,:,:) = va(:,:,:)
90      ENDIF
[2528]91      !
[3294]92      SELECT CASE ( nhpg )      ! Hydrostatic pressure gradient computation
[503]93      CASE (  0 )   ;   CALL hpg_zco    ( kt )      ! z-coordinate
94      CASE (  1 )   ;   CALL hpg_zps    ( kt )      ! z-coordinate plus partial steps (interpolation)
95      CASE (  2 )   ;   CALL hpg_sco    ( kt )      ! s-coordinate (standard jacobian formulation)
[3294]96      CASE (  3 )   ;   CALL hpg_djc    ( kt )      ! s-coordinate (Density Jacobian with Cubic polynomial)
97      CASE (  4 )   ;   CALL hpg_prj    ( kt )      ! s-coordinate (Pressure Jacobian scheme)
[455]98      END SELECT
[2528]99      !
[503]100      IF( l_trddyn ) THEN      ! save the hydrostatic pressure gradient trends for momentum trend diagnostics
[455]101         ztrdu(:,:,:) = ua(:,:,:) - ztrdu(:,:,:)
102         ztrdv(:,:,:) = va(:,:,:) - ztrdv(:,:,:)
[503]103         CALL trd_mod( ztrdu, ztrdv, jpdyn_trd_hpg, 'DYN', kt )
[3294]104         CALL wrk_dealloc( jpi,jpj,jpk, ztrdu, ztrdv )
[3764]105      ENDIF
[503]106      !
107      IF(ln_ctl)   CALL prt_ctl( tab3d_1=ua, clinfo1=' hpg  - Ua: ', mask1=umask,   &
108         &                       tab3d_2=va, clinfo2=       ' Va: ', mask2=vmask, clinfo3='dyn' )
109      !
[3294]110      IF( nn_timing == 1 )  CALL timing_stop('dyn_hpg')
[2715]111      !
[455]112   END SUBROUTINE dyn_hpg
113
114
[2528]115   SUBROUTINE dyn_hpg_init
[455]116      !!----------------------------------------------------------------------
[2528]117      !!                 ***  ROUTINE dyn_hpg_init  ***
[455]118      !!
119      !! ** Purpose :   initializations for the hydrostatic pressure gradient
120      !!              computation and consistency control
121      !!
[1601]122      !! ** Action  :   Read the namelist namdyn_hpg and check the consistency
[455]123      !!      with the type of vertical coordinate used (zco, zps, sco)
124      !!----------------------------------------------------------------------
125      INTEGER ::   ioptio = 0      ! temporary integer
[4147]126      INTEGER ::   ios             ! Local integer output status for namelist read
[1601]127      !!
[3294]128      NAMELIST/namdyn_hpg/ ln_hpg_zco, ln_hpg_zps, ln_hpg_sco,     &
129         &                 ln_hpg_djc, ln_hpg_prj, ln_dynhpg_imp
[455]130      !!----------------------------------------------------------------------
[2528]131      !
[4147]132      REWIND( numnam_ref )              ! Namelist namdyn_hpg in reference namelist : Hydrostatic pressure gradient
133      READ  ( numnam_ref, namdyn_hpg, IOSTAT = ios, ERR = 901)
134901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namdyn_hpg in reference namelist', lwp )
135
136      REWIND( numnam_cfg )              ! Namelist namdyn_hpg in configuration namelist : Hydrostatic pressure gradient
137      READ  ( numnam_cfg, namdyn_hpg, IOSTAT = ios, ERR = 902 )
138902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namdyn_hpg in configuration namelist', lwp )
[4666]139      WRITE ( numond, namdyn_hpg )
[2528]140      !
141      IF(lwp) THEN                   ! Control print
[455]142         WRITE(numout,*)
[2528]143         WRITE(numout,*) 'dyn_hpg_init : hydrostatic pressure gradient initialisation'
144         WRITE(numout,*) '~~~~~~~~~~~~'
[1601]145         WRITE(numout,*) '   Namelist namdyn_hpg : choice of hpg scheme'
146         WRITE(numout,*) '      z-coord. - full steps                             ln_hpg_zco    = ', ln_hpg_zco
147         WRITE(numout,*) '      z-coord. - partial steps (interpolation)          ln_hpg_zps    = ', ln_hpg_zps
148         WRITE(numout,*) '      s-coord. (standard jacobian formulation)          ln_hpg_sco    = ', ln_hpg_sco
149         WRITE(numout,*) '      s-coord. (Density Jacobian: Cubic polynomial)     ln_hpg_djc    = ', ln_hpg_djc
[3294]150         WRITE(numout,*) '      s-coord. (Pressure Jacobian: Cubic polynomial)    ln_hpg_prj    = ', ln_hpg_prj
[1601]151         WRITE(numout,*) '      time stepping: centered (F) or semi-implicit (T)  ln_dynhpg_imp = ', ln_dynhpg_imp
[455]152      ENDIF
[2528]153      !
[3294]154      IF( ln_hpg_djc )   &
155         &   CALL ctl_stop('dyn_hpg_init : Density Jacobian: Cubic polynominal method &
156                           & currently disabled (bugs under investigation). Please select &
157                           & either  ln_hpg_sco or  ln_hpg_prj instead')
[2528]158      !
[3294]159      IF( lk_vvl .AND. .NOT. (ln_hpg_sco.OR.ln_hpg_prj) )   &
160         &   CALL ctl_stop('dyn_hpg_init : variable volume key_vvl requires:&
161                           & the standard jacobian formulation hpg_sco or &
162                           & the pressure jacobian formulation hpg_prj')
163      !
[503]164      !                               ! Set nhpg from ln_hpg_... flags
[455]165      IF( ln_hpg_zco )   nhpg = 0
166      IF( ln_hpg_zps )   nhpg = 1
167      IF( ln_hpg_sco )   nhpg = 2
[3294]168      IF( ln_hpg_djc )   nhpg = 3
169      IF( ln_hpg_prj )   nhpg = 4
[2528]170      !
[3294]171      !                               ! Consistency check
[3764]172      ioptio = 0
[455]173      IF( ln_hpg_zco )   ioptio = ioptio + 1
174      IF( ln_hpg_zps )   ioptio = ioptio + 1
175      IF( ln_hpg_sco )   ioptio = ioptio + 1
176      IF( ln_hpg_djc )   ioptio = ioptio + 1
[3294]177      IF( ln_hpg_prj )   ioptio = ioptio + 1
[2715]178      IF( ioptio /= 1 )   CALL ctl_stop( 'NO or several hydrostatic pressure gradient options used' )
[4724]179      IF( (ln_hpg_zco .OR. ln_hpg_zps .OR. ln_hpg_djc .OR. ln_hpg_prj ) .AND. nn_isf .NE. 0 )   CALL ctl_stop( 'Only hpg_sco has been corrected to work with ice shelf cavity. Comparison in a GYRE simulation with bump in the middle show similar result than hpg_zps' )
[503]180      !
[2528]181   END SUBROUTINE dyn_hpg_init
[455]182
183
184   SUBROUTINE hpg_zco( kt )
185      !!---------------------------------------------------------------------
186      !!                  ***  ROUTINE hpg_zco  ***
187      !!
188      !! ** Method  :   z-coordinate case, levels are horizontal surfaces.
189      !!      The now hydrostatic pressure gradient at a given level, jk,
190      !!      is computed by taking the vertical integral of the in-situ
191      !!      density gradient along the model level from the suface to that
192      !!      level:    zhpi = grav .....
193      !!                zhpj = grav .....
[3]194      !!      add it to the general momentum trend (ua,va).
[455]195      !!            ua = ua - 1/e1u * zhpi
196      !!            va = va - 1/e2v * zhpj
[3764]197      !!
[3]198      !! ** Action : - Update (ua,va) with the now hydrastatic pressure trend
[503]199      !!----------------------------------------------------------------------
200      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
201      !!
202      INTEGER  ::   ji, jj, jk       ! dummy loop indices
203      REAL(wp) ::   zcoef0, zcoef1   ! temporary scalars
[3764]204      REAL(wp), POINTER, DIMENSION(:,:,:) ::  zhpi, zhpj
[3]205      !!----------------------------------------------------------------------
[3764]206      !
[3294]207      CALL wrk_alloc( jpi,jpj,jpk, zhpi, zhpj )
208      !
[3]209      IF( kt == nit000 ) THEN
210         IF(lwp) WRITE(numout,*)
[455]211         IF(lwp) WRITE(numout,*) 'dyn:hpg_zco : hydrostatic pressure gradient trend'
212         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   z-coordinate case '
[3]213      ENDIF
214
[3764]215      zcoef0 = - grav * 0.5_wp      ! Local constant initialization
216
[455]217      ! Surface value
[3]218      DO jj = 2, jpjm1
219         DO ji = fs_2, fs_jpim1   ! vector opt.
[455]220            zcoef1 = zcoef0 * fse3w(ji,jj,1)
221            ! hydrostatic pressure gradient
222            zhpi(ji,jj,1) = zcoef1 * ( rhd(ji+1,jj,1) - rhd(ji,jj,1) ) / e1u(ji,jj)
223            zhpj(ji,jj,1) = zcoef1 * ( rhd(ji,jj+1,1) - rhd(ji,jj,1) ) / e2v(ji,jj)
[3]224            ! add to the general momentum trend
[455]225            ua(ji,jj,1) = ua(ji,jj,1) + zhpi(ji,jj,1)
226            va(ji,jj,1) = va(ji,jj,1) + zhpj(ji,jj,1)
227         END DO
228      END DO
[3294]229
[503]230      !
[455]231      ! interior value (2=<jk=<jpkm1)
[3]232      DO jk = 2, jpkm1
[455]233         DO jj = 2, jpjm1
[3]234            DO ji = fs_2, fs_jpim1   ! vector opt.
[455]235               zcoef1 = zcoef0 * fse3w(ji,jj,jk)
236               ! hydrostatic pressure gradient
237               zhpi(ji,jj,jk) = zhpi(ji,jj,jk-1)   &
238                  &           + zcoef1 * (  ( rhd(ji+1,jj,jk)+rhd(ji+1,jj,jk-1) )   &
239                  &                       - ( rhd(ji  ,jj,jk)+rhd(ji  ,jj,jk-1) )  ) / e1u(ji,jj)
240
241               zhpj(ji,jj,jk) = zhpj(ji,jj,jk-1)   &
242                  &           + zcoef1 * (  ( rhd(ji,jj+1,jk)+rhd(ji,jj+1,jk-1) )   &
243                  &                       - ( rhd(ji,jj,  jk)+rhd(ji,jj  ,jk-1) )  ) / e2v(ji,jj)
[3]244               ! add to the general momentum trend
[455]245               ua(ji,jj,jk) = ua(ji,jj,jk) + zhpi(ji,jj,jk)
246               va(ji,jj,jk) = va(ji,jj,jk) + zhpj(ji,jj,jk)
[3]247            END DO
248         END DO
249      END DO
[503]250      !
[3294]251      CALL wrk_dealloc( jpi,jpj,jpk, zhpi, zhpj )
252      !
[455]253   END SUBROUTINE hpg_zco
[216]254
[3]255
[455]256   SUBROUTINE hpg_zps( kt )
[3]257      !!---------------------------------------------------------------------
[455]258      !!                 ***  ROUTINE hpg_zps  ***
[3764]259      !!
[455]260      !! ** Method  :   z-coordinate plus partial steps case.  blahblah...
[3764]261      !!
[3]262      !! ** Action  : - Update (ua,va) with the now hydrastatic pressure trend
[3764]263      !!----------------------------------------------------------------------
[503]264      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
265      !!
266      INTEGER  ::   ji, jj, jk                       ! dummy loop indices
267      INTEGER  ::   iku, ikv                         ! temporary integers
268      REAL(wp) ::   zcoef0, zcoef1, zcoef2, zcoef3   ! temporary scalars
[3764]269      REAL(wp), POINTER, DIMENSION(:,:,:) ::  zhpi, zhpj
[3]270      !!----------------------------------------------------------------------
[3294]271      !
272      CALL wrk_alloc( jpi,jpj,jpk, zhpi, zhpj )
273      !
[3]274      IF( kt == nit000 ) THEN
275         IF(lwp) WRITE(numout,*)
[455]276         IF(lwp) WRITE(numout,*) 'dyn:hpg_zps : hydrostatic pressure gradient trend'
[503]277         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   z-coordinate with partial steps - vector optimization'
[3]278      ENDIF
279
[3294]280
[503]281      ! Local constant initialization
[2528]282      zcoef0 = - grav * 0.5_wp
[3]283
[2528]284      !  Surface value (also valid in partial step case)
[3]285      DO jj = 2, jpjm1
286         DO ji = fs_2, fs_jpim1   ! vector opt.
[170]287            zcoef1 = zcoef0 * fse3w(ji,jj,1)
[3]288            ! hydrostatic pressure gradient
[455]289            zhpi(ji,jj,1) = zcoef1 * ( rhd(ji+1,jj  ,1) - rhd(ji,jj,1) ) / e1u(ji,jj)
290            zhpj(ji,jj,1) = zcoef1 * ( rhd(ji  ,jj+1,1) - rhd(ji,jj,1) ) / e2v(ji,jj)
[3]291            ! add to the general momentum trend
292            ua(ji,jj,1) = ua(ji,jj,1) + zhpi(ji,jj,1)
293            va(ji,jj,1) = va(ji,jj,1) + zhpj(ji,jj,1)
294         END DO
295      END DO
296
[3294]297
[503]298      ! interior value (2=<jk=<jpkm1)
[3]299      DO jk = 2, jpkm1
300         DO jj = 2, jpjm1
301            DO ji = fs_2, fs_jpim1   ! vector opt.
[170]302               zcoef1 = zcoef0 * fse3w(ji,jj,jk)
[3]303               ! hydrostatic pressure gradient
304               zhpi(ji,jj,jk) = zhpi(ji,jj,jk-1)   &
[455]305                  &           + zcoef1 * (  ( rhd(ji+1,jj,jk) + rhd(ji+1,jj,jk-1) )   &
306                  &                       - ( rhd(ji  ,jj,jk) + rhd(ji  ,jj,jk-1) )  ) / e1u(ji,jj)
[3]307
308               zhpj(ji,jj,jk) = zhpj(ji,jj,jk-1)   &
[455]309                  &           + zcoef1 * (  ( rhd(ji,jj+1,jk) + rhd(ji,jj+1,jk-1) )   &
310                  &                       - ( rhd(ji,jj,  jk) + rhd(ji,jj  ,jk-1) )  ) / e2v(ji,jj)
[3]311               ! add to the general momentum trend
312               ua(ji,jj,jk) = ua(ji,jj,jk) + zhpi(ji,jj,jk)
313               va(ji,jj,jk) = va(ji,jj,jk) + zhpj(ji,jj,jk)
[455]314            END DO
[3]315         END DO
316      END DO
317
[3294]318
[2528]319      ! partial steps correction at the last level  (use gru & grv computed in zpshde.F90)
[3]320# if defined key_vectopt_loop
321         jj = 1
322         DO ji = jpi+2, jpij-jpi-1   ! vector opt. (forced unrolling)
323# else
324      DO jj = 2, jpjm1
325         DO ji = 2, jpim1
326# endif
[2528]327            iku = mbku(ji,jj)
328            ikv = mbkv(ji,jj)
[3]329            zcoef2 = zcoef0 * MIN( fse3w(ji,jj,iku), fse3w(ji+1,jj  ,iku) )
330            zcoef3 = zcoef0 * MIN( fse3w(ji,jj,ikv), fse3w(ji  ,jj+1,ikv) )
[2528]331            IF( iku > 1 ) THEN            ! on i-direction (level 2 or more)
332               ua  (ji,jj,iku) = ua(ji,jj,iku) - zhpi(ji,jj,iku)         ! subtract old value
333               zhpi(ji,jj,iku) = zhpi(ji,jj,iku-1)                   &   ! compute the new one
334                  &            + zcoef2 * ( rhd(ji+1,jj,iku-1) - rhd(ji,jj,iku-1) + gru(ji,jj) ) / e1u(ji,jj)
335               ua  (ji,jj,iku) = ua(ji,jj,iku) + zhpi(ji,jj,iku)         ! add the new one to the general momentum trend
[3]336            ENDIF
[2528]337            IF( ikv > 1 ) THEN            ! on j-direction (level 2 or more)
338               va  (ji,jj,ikv) = va(ji,jj,ikv) - zhpj(ji,jj,ikv)         ! subtract old value
339               zhpj(ji,jj,ikv) = zhpj(ji,jj,ikv-1)                   &   ! compute the new one
340                  &            + zcoef3 * ( rhd(ji,jj+1,ikv-1) - rhd(ji,jj,ikv-1) + grv(ji,jj) ) / e2v(ji,jj)
341               va  (ji,jj,ikv) = va(ji,jj,ikv) + zhpj(ji,jj,ikv)         ! add the new one to the general momentum trend
[3]342            ENDIF
343# if ! defined key_vectopt_loop
344         END DO
345# endif
346      END DO
[503]347      !
[3294]348      CALL wrk_dealloc( jpi,jpj,jpk, zhpi, zhpj )
349      !
[455]350   END SUBROUTINE hpg_zps
[216]351
[3]352
[455]353   SUBROUTINE hpg_sco( kt )
[3]354      !!---------------------------------------------------------------------
[455]355      !!                  ***  ROUTINE hpg_sco  ***
[3]356      !!
[455]357      !! ** Method  :   s-coordinate case. Jacobian scheme.
358      !!      The now hydrostatic pressure gradient at a given level, jk,
359      !!      is computed by taking the vertical integral of the in-situ
[3]360      !!      density gradient along the model level from the suface to that
[455]361      !!      level. s-coordinates (ln_sco): a corrective term is added
362      !!      to the horizontal pressure gradient :
363      !!         zhpi = grav .....  + 1/e1u mi(rhd) di[ grav dep3w ]
364      !!         zhpj = grav .....  + 1/e2v mj(rhd) dj[ grav dep3w ]
[3]365      !!      add it to the general momentum trend (ua,va).
[455]366      !!         ua = ua - 1/e1u * zhpi
367      !!         va = va - 1/e2v * zhpj
[3]368      !!
369      !! ** Action : - Update (ua,va) with the now hydrastatic pressure trend
[503]370      !!----------------------------------------------------------------------
371      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
372      !!
[4666]373      INTEGER  ::   ji, jj, jk, iku, ikv, ikt, iktp1i, iktp1j                 ! dummy loop indices
374      REAL(wp) ::   zcoef0, zuap, zvap, znad, ze3wu, ze3wv, zuapint, zvapint, zhpjint, zhpiint, zdzwt, zdzwtjp1, zdzwtip1             ! temporary scalars
375      REAL(wp), POINTER, DIMENSION(:,:,:)   ::  zhpi, zhpj, zrhd
376      REAL(wp), POINTER, DIMENSION(:,:,:)   ::  ztstop
377      REAL(wp), POINTER, DIMENSION(:,:)     ::  ze3w, zp, zrhdtop_isf, zrhdtop_oce, ziceload, zdept, zpshpi, zpshpj
[3]378      !!----------------------------------------------------------------------
[3294]379      !
[4666]380      CALL wrk_alloc( jpi,jpj, 2, ztstop) 
381      CALL wrk_alloc( jpi,jpj,jpk, zhpi, zhpj, zrhd)
382      CALL wrk_alloc( jpi,jpj, ze3w, zp, zrhdtop_isf, zrhdtop_oce, ziceload, zdept, zpshpi, zpshpj) 
[3294]383      !
[4666]384     IF( kt == nit000 ) THEN
[3]385         IF(lwp) WRITE(numout,*)
[455]386         IF(lwp) WRITE(numout,*) 'dyn:hpg_sco : hydrostatic pressure gradient trend'
387         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   s-coordinate case, OPA original scheme used'
[3]388      ENDIF
389
[503]390      ! Local constant initialization
[2528]391      zcoef0 = - grav * 0.5_wp
[592]392      ! To use density and not density anomaly
[4666]393!      IF ( lk_vvl ) THEN   ;     znad = 1._wp          ! Variable volume
394!      ELSE                 ;     znad = 0._wp         ! Fixed volume
395!      ENDIF
396      znad=1._wp
397      ! iniitialised to 0. zhpi zhpi
398      zhpi(:,:,:)=0._wp ; zhpj(:,:,:)=0._wp
[4724]399     
400      ! assume density of water displaced by the ice shelf is at T=-1.9 and S=34.4 (rude)
[4666]401      ztstop(:,:,1)=-1.9_wp ; ztstop(:,:,2)=34.4_wp
[4724]402      zrhd = rhd ! save rhd
[4666]403      DO jk = 1, jpk
404           zdept(:,:)=gdept_1d(jk)
405           CALL eos(ztstop(:,:,:),zdept(:,:),rhd(:,:,jk))
406      END DO
407      WHERE ( tmask(:,:,:) == 1._wp)
[4724]408        rhd(:,:,:) = zrhd(:,:,:) ! replace wet cell by the saved rhd
[4666]409      END WHERE
410     
[4724]411      ! compute rhd at the ice/oce interface (ice shelf side)
[4726]412      CALL eos(ztstop,risfdep,zrhdtop_isf)
[4666]413
[4724]414      ! compute rhd at the ice/oce interface (ocean side)
[4666]415      DO ji=1,jpi
416        DO jj=1,jpj
417          ikt=mikt(ji,jj)
418          ztstop(ji,jj,1)=tsn(ji,jj,ikt,1)
419          ztstop(ji,jj,2)=tsn(ji,jj,ikt,2)
420        END DO
421      END DO
[4726]422      CALL eos(ztstop,risfdep,zrhdtop_oce)
[4666]423      !
424      ! Surface value + ice shelf gradient
[4724]425      ! compute pressure due to ice shelf load (used to compute hpgi/j for all the level from 1 to miku/v)
[4666]426      ziceload = 0._wp
427      DO jj = 1, jpj
428         DO ji = 1, jpi   ! vector opt.
429            ikt=mikt(ji,jj)
430            ziceload(ji,jj) = ziceload(ji,jj) + (znad + rhd(ji,jj,1) ) * fse3w(ji,jj,1) * (1._wp - tmask(ji,jj,1))
431            DO jk=2,ikt-1
432               ziceload(ji,jj) = ziceload(ji,jj) + (2._wp * znad + rhd(ji,jj,jk-1) + rhd(ji,jj,jk)) * fse3w(ji,jj,jk) &
433                  &                              * (1._wp - tmask(ji,jj,jk))
434            END DO
435            IF (ikt .GE. 2) ziceload(ji,jj) = ziceload(ji,jj) + (2._wp * znad + zrhdtop_isf(ji,jj) + rhd(ji,jj,ikt-1)) &
436                               &                              * ( fsdepw(ji,jj,ikt) - gdept_1d(ikt-1) )
437         END DO
438      END DO
[4724]439      riceload(:,:) = 0.0_wp ; riceload(:,:)=ziceload(:,:)  ! need to be saved for diaar5
440      ! compute zp from z=0 to first T wet point (correction due to zps not yet applied)
[455]441      DO jj = 2, jpjm1
[3764]442         DO ji = fs_2, fs_jpim1   ! vector opt.
[4666]443            ikt=mikt(ji,jj) ; iktp1i=mikt(ji+1,jj); iktp1j=mikt(ji,jj+1)
444            ! hydrostatic pressure gradient along s-surfaces and ice shelf pressure
[4724]445            ! we assume ISF is in isostatic equilibrium
[4666]446            zhpi(ji,jj,1) = zcoef0 / e1u(ji,jj) * ( 0.5_wp * fse3w(ji+1,jj  ,iktp1i)                                    &
447               &                                  * ( 2._wp * znad + rhd(ji+1,jj  ,iktp1i) + zrhdtop_oce(ji+1,jj  ) )   &
448               &                                  - 0.5_wp * fse3w(ji  ,jj  ,ikt   )                                    &
449               &                                  * ( 2._wp * znad + rhd(ji  ,jj  ,ikt   ) + zrhdtop_oce(ji  ,jj  ) )   &
450               &                                  + ( ziceload(ji+1,jj) - ziceload(ji,jj))                              ) 
451            zhpj(ji,jj,1) = zcoef0 / e2v(ji,jj) * ( 0.5_wp * fse3w(ji  ,jj+1,iktp1j)                                    &
452               &                                  * ( 2._wp * znad + rhd(ji  ,jj+1,iktp1j) + zrhdtop_oce(ji  ,jj+1) )   &
453               &                                  - 0.5_wp * fse3w(ji  ,jj  ,ikt   )                                    & 
454               &                                  * ( 2._wp * znad + rhd(ji  ,jj  ,ikt   ) + zrhdtop_oce(ji  ,jj  ) )   &
455               &                                  + ( ziceload(ji,jj+1) - ziceload(ji,jj) )                             ) 
[4724]456            ! s-coordinate pressure gradient correction (=0 if z coordinate)
[2528]457            zuap = -zcoef0 * ( rhd   (ji+1,jj,1) + rhd   (ji,jj,1) + 2._wp * znad )   &
[455]458               &           * ( fsde3w(ji+1,jj,1) - fsde3w(ji,jj,1) ) / e1u(ji,jj)
[2528]459            zvap = -zcoef0 * ( rhd   (ji,jj+1,1) + rhd   (ji,jj,1) + 2._wp * znad )   &
[455]460               &           * ( fsde3w(ji,jj+1,1) - fsde3w(ji,jj,1) ) / e2v(ji,jj)
461            ! add to the general momentum trend
[4666]462            ua(ji,jj,1) = ua(ji,jj,1) + (zhpi(ji,jj,1) + zuap) * umask(ji,jj,1)
463            va(ji,jj,1) = va(ji,jj,1) + (zhpj(ji,jj,1) + zvap) * vmask(ji,jj,1)
[3764]464         END DO
465      END DO
[4666]466      !
467      DO jj = 2, jpjm1
468         DO ji = fs_2, fs_jpim1   ! vector opt.
469            iku = miku(ji,jj) ; 
470            zpshpi(ji,jj)=0.0_wp ; zpshpj(ji,jj)=0.0_wp
471            ze3wu  = fse3w(ji+1,jj,iku+1)-fse3w(ji,jj,iku+1) 
472            ! u direction
473            IF ( iku .GT. 1 ) THEN
474               ! case iku
475               zhpi(ji,jj,iku)   =  zcoef0 / e1u(ji,jj) * ze3wu                                            &
476                      &                                 * ( rhd    (ji+1,jj,iku) + rhd   (ji,jj,iku)       &
477                      &                                   + SIGN(1._wp,ze3wu) * grui(ji,jj) + 2._wp * znad )
[4724]478               ! corrective term ( = 0 if z coordinate )
[4666]479               zuap              = -zcoef0 * ( arui(ji,jj) + 2._wp * znad ) * gzui(ji,jj) / e1u(ji,jj)
480               ! zhpi will be added in interior loop and zuapint will be removed in the interior loop
481               ua(ji,jj,iku)     = ua(ji,jj,iku) + zuap
[4724]482               ! in case of 2 cell water column, need to save the pressure gradient to compute the bottom pressure 
483               IF (mbku(ji,jj) == iku + 1) zpshpi(ji,jj)  = zhpi(ji,jj,iku)
[3764]484
[4666]485               ! zhpi will be added in interior loop and zuapint will be removed in the interior loop
486               ! case iku + 1 (remove the zphi term added in the interior loop)
487               zhpiint        =  zcoef0 / e1u(ji,jj)                                                             &   
488                  &           * (  fse3w(ji+1,jj  ,iku+1) * ( (rhd(ji+1,jj,iku+1) + znad)                          &
489                  &                                         + (rhd(ji+1,jj,iku  ) + znad) ) * tmask(ji+1,jj,iku)   &
490                  &              - fse3w(ji  ,jj  ,iku+1) * ( (rhd(ji  ,jj,iku+1) + znad)                          &
491                  &                                         + (rhd(ji  ,jj,iku  ) + znad) ) * tmask(ji  ,jj,iku)   )
492               zhpi(ji,jj,iku+1) =  zcoef0 / e1u(ji,jj) * ge3rui(ji,jj) - zhpiint 
493            END IF
494            ! v direction
495            ikv = mikv(ji,jj)
496            ze3wv  = fse3w(ji,jj+1,ikv+1)-fse3w(ji,jj,ikv+1) 
497            IF ( ikv .GT. 1 ) THEN
498               ! case ikv
499               zhpj(ji,jj,ikv)   =  zcoef0 / e2v(ji,jj) * ze3wv                                            &
500                     &                                  * ( rhd(ji,jj+1,ikv) + rhd   (ji,jj,ikv)           &
501                     &                                    + SIGN(1._wp,ze3wv) * grvi(ji,jj) + 2._wp * znad )
[4724]502               ! corrective term ( = 0 if z coordinate )
[4666]503               zvap              = -zcoef0 * ( arvi(ji,jj) + 2._wp * znad ) * gzvi(ji,jj) / e2v(ji,jj)
504               ! zhpi will be added in interior loop and zvapint will be removed in the interior loop
505               va(ji,jj,ikv)      = va(ji,jj,ikv) + zvap
[4724]506               ! in case of 2 cell water column, need to save the pressure gradient to compute the bottom pressure 
[4666]507               IF (mbkv(ji,jj) == ikv + 1)  zpshpj(ji,jj)  =  zhpj(ji,jj,ikv) 
508               ! zhpi will be added in interior loop and zvapint will be removed in the interior loop
509               ! case ikv + 1 (remove the zphj term added in the interior loop)
510               zhpjint        =  zcoef0 / e2v(ji,jj)                                                            &
511                  &           * (  fse3w(ji  ,jj+1,ikv+1) * ( (rhd(ji,jj+1,ikv+1) + znad)                         &
512                  &                                       + (rhd(ji,jj+1,ikv  ) + znad) ) * tmask(ji,jj+1,ikv)  &
513                  &              - fse3w(ji  ,jj  ,ikv+1) * ( (rhd(ji,jj  ,ikv+1) + znad)                         &
514                  &                                       + (rhd(ji,jj  ,ikv  ) + znad) ) * tmask(ji,jj  ,ikv)  )
515               zhpj(ji,jj,ikv+1) =  zcoef0 / e2v(ji,jj) * ge3rvi(ji,jj) - zhpjint
516            END IF
517         END DO
518      END DO
519
[503]520      ! interior value (2=<jk=<jpkm1)
[4666]521      DO jj = 2, jpjm1
522         DO ji = fs_2, fs_jpim1   ! vector opt.
523            iku=miku(ji,jj); ikv=mikv(ji,jj)
524            DO jk = 2, jpkm1
[455]525               ! hydrostatic pressure gradient along s-surfaces
526               ! s-coordinate pressure gradient correction
[4666]527               zhpi(ji,jj,jk) = zhpi(ji,jj,jk) + zhpi(ji,jj,jk-1)                                                              &
528                  &                            + zcoef0 / e1u(ji,jj)                                                           &
529                  &                                     * ( fse3w(ji+1,jj  ,jk) * ( (rhd(ji+1,jj,jk  ) + znad)                 &
530                  &                                                     + (rhd(ji+1,jj,jk-1) + znad) ) * tmask(ji+1,jj,jk-1)   &
531                  &                                       - fse3w(ji  ,jj  ,jk) * ( (rhd(ji  ,jj,jk  ) + znad)                 &
532                  &                                                     + (rhd(ji  ,jj,jk-1) + znad) ) * tmask(ji  ,jj,jk-1)   ) 
[4724]533             !corrective term, we mask this term for the first wet level beneath the ice shelf (contribution done in the upper bloc)
[4666]534               zuap = - zcoef0 * ( rhd   (ji+1,jj  ,jk) + rhd   (ji,jj,jk) + 2._wp * znad )   &
535                  &            * ( fsde3w(ji+1,jj  ,jk) - fsde3w(ji,jj,jk) ) / e1u(ji,jj) * umask(ji,jj,jk-1)
536               ua(ji,jj,jk) = ua(ji,jj,jk) + ( zhpi(ji,jj,jk) + zuap) * umask(ji,jj,jk)
537
538               zhpj(ji,jj,jk) = zhpj(ji,jj,jk) + zhpj(ji,jj,jk-1)                                                              &
539                  &                            + zcoef0 / e2v(ji,jj)                                                           &
540                  &                                     * ( fse3w(ji  ,jj+1,jk) * ( (rhd(ji,jj+1,jk  ) + znad)                 &
541                  &                                                     + (rhd(ji,jj+1,jk-1) + znad) ) * tmask(ji,jj+1,jk-1)   &
542                  &                                       - fse3w(ji  ,jj  ,jk) * ( (rhd(ji,jj  ,jk  ) + znad)                 &
543                  &                                                     + (rhd(ji,jj  ,jk-1) + znad) ) * tmask(ji,jj  ,jk-1)   )
[4724]544             !corrective term, we mask this term for the first wet level beneath the ice shelf (contribution done in the upper bloc)
[4666]545               zvap = - zcoef0 * ( rhd   (ji  ,jj+1,jk) + rhd   (ji,jj,jk) + 2._wp * znad )   &
546                  &            * ( fsde3w(ji  ,jj+1,jk) - fsde3w(ji,jj,jk) ) / e2v(ji,jj) * vmask(ji,jj,jk-1)
[455]547               ! add to the general momentum trend
[4666]548               va(ji,jj,jk) = va(ji,jj,jk) + ( zhpj(ji,jj,jk) + zvap ) * vmask(ji,jj,jk)
[455]549            END DO
550         END DO
551      END DO
[4666]552     
553      ! partial steps correction at the last level  (use gru & grv computed in zpshde.F90)
554# if defined key_vectopt_loop
555         jj = 1
556         DO ji = jpi+2, jpij-jpi-1   ! vector opt. (forced unrolling)
557# else
558      DO jj = 2, jpjm1
559         DO ji = 2, jpim1
560# endif
561            iku = mbku(ji,jj)
562            ikv = mbkv(ji,jj)
563
564            IF (iku .GT. 1) THEN
565               ! remove old value (interior case)
566               zuap            = -zcoef0 * ( rhd   (ji+1,jj  ,iku) + rhd   (ji,jj,iku) + 2._wp * znad )   &
567                     &                   * ( fsde3w(ji+1,jj  ,iku) - fsde3w(ji,jj,iku) ) / e1u(ji,jj)
568               ua(ji,jj,iku)   = ua(ji,jj,iku) - zhpi(ji,jj,iku) - zuap
569               ! put new value
[4724]570               ! -zpshpi to avoid double contribution of the partial step in the top layer
[4666]571               zuap            = -zcoef0 * ( aru(ji,jj) + 2._wp * znad ) * gzu(ji,jj)  / e1u(ji,jj)
572               zhpi(ji,jj,iku) =  zhpi(ji,jj,iku-1) + zcoef0 / e1u(ji,jj) * ge3ru(ji,jj) - zpshpi(ji,jj) 
573               ua(ji,jj,iku)   =  ua(ji,jj,iku) + zhpi(ji,jj,iku) + zuap
574            END IF
575            ! v direction
576            IF (ikv .GT. 1) THEN
577               ! remove old value (interior case)
578               zvap            = -zcoef0 * ( rhd   (ji  ,jj+1,ikv) + rhd   (ji,jj,ikv) + 2._wp * znad )   &
579                     &                   * ( fsde3w(ji  ,jj+1,ikv) - fsde3w(ji,jj,ikv) )   / e2v(ji,jj)
580               va(ji,jj,ikv)   = va(ji,jj,ikv) - zhpj(ji,jj,ikv) - zvap
581               ! put new value
[4724]582               ! -zpshpj to avoid double contribution of the partial step in the top layer
[4666]583               zvap            = -zcoef0 * ( arv(ji,jj) + 2._wp * znad ) * gzv(ji,jj)     / e2v(ji,jj)
584               zhpj(ji,jj,ikv) =  zhpj(ji,jj,ikv-1) + zcoef0 / e2v(ji,jj) * ge3rv(ji,jj) - zpshpj(ji,jj)   
585               va(ji,jj,ikv)   =  va(ji,jj,ikv) + zhpj(ji,jj,ikv) + zvap
586            END IF
587# if ! defined key_vectopt_loop
588         END DO
589# endif
590      END DO
591
592
593      rhd = zrhd
[503]594      !
[4666]595      CALL wrk_dealloc( jpi,jpj,2, ztstop)
596      CALL wrk_dealloc( jpi,jpj,jpk, zhpi, zhpj, zrhd)
597      CALL wrk_dealloc( jpi,jpj, ze3w, zp, zrhdtop_isf, zrhdtop_oce, ziceload, zdept, zpshpi, zpshpj)
[3294]598      !
[455]599   END SUBROUTINE hpg_sco
600
601   SUBROUTINE hpg_djc( kt )
602      !!---------------------------------------------------------------------
603      !!                  ***  ROUTINE hpg_djc  ***
604      !!
605      !! ** Method  :   Density Jacobian with Cubic polynomial scheme
[3764]606      !!
[503]607      !! Reference: Shchepetkin and McWilliams, J. Geophys. Res., 108(C3), 3090, 2003
[455]608      !!----------------------------------------------------------------------
[503]609      INTEGER, INTENT(in) ::   kt    ! ocean time-step index
610      !!
611      INTEGER  ::   ji, jj, jk          ! dummy loop indices
612      REAL(wp) ::   zcoef0, zep, cffw   ! temporary scalars
613      REAL(wp) ::   z1_10, cffu, cffx   !    "         "
614      REAL(wp) ::   z1_12, cffv, cffy   !    "         "
[3764]615      REAL(wp), POINTER, DIMENSION(:,:,:) ::  zhpi, zhpj
[3294]616      REAL(wp), POINTER, DIMENSION(:,:,:) ::  dzx, dzy, dzz, dzu, dzv, dzw
617      REAL(wp), POINTER, DIMENSION(:,:,:) ::  drhox, drhoy, drhoz, drhou, drhov, drhow
618      REAL(wp), POINTER, DIMENSION(:,:,:) ::  rho_i, rho_j, rho_k
[455]619      !!----------------------------------------------------------------------
[3294]620      !
[3764]621      CALL wrk_alloc( jpi, jpj, jpk, dzx  , dzy  , dzz  , dzu  , dzv  , dzw   )
622      CALL wrk_alloc( jpi, jpj, jpk, drhox, drhoy, drhoz, drhou, drhov, drhow )
623      CALL wrk_alloc( jpi, jpj, jpk, rho_i, rho_j, rho_k,  zhpi,  zhpj        )
[3294]624      !
[455]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'
[216]630      ENDIF
631
[503]632      ! Local constant initialization
[2528]633      zcoef0 = - grav * 0.5_wp
634      z1_10  = 1._wp / 10._wp
635      z1_12  = 1._wp / 12._wp
[455]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
[503]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
[455]663
664      DO jk = 2, jpkm1
665         DO jj = 2, jpjm1
666            DO ji = fs_2, fs_jpim1   ! vector opt.
[2528]667               cffw = 2._wp * drhoz(ji  ,jj  ,jk) * drhoz(ji,jj,jk-1)
[455]668
[2528]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  )
[3764]671
[2528]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  )
[455]674
675               IF( cffw > zep) THEN
[2528]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) )
[455]678               ELSE
[2528]679                  drhow(ji,jj,jk) = 0._wp
[455]680               ENDIF
681
[2528]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) )
[455]684
685               IF( cffu > zep ) THEN
[2528]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) )
[455]688               ELSE
[2528]689                  drhou(ji,jj,jk ) = 0._wp
[455]690               ENDIF
691
692               IF( cffx > zep ) THEN
[2528]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) )
[455]695               ELSE
[2528]696                  dzu(ji,jj,jk) = 0._wp
[455]697               ENDIF
698
699               IF( cffv > zep ) THEN
[2528]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) )
[455]702               ELSE
[2528]703                  drhov(ji,jj,jk) = 0._wp
[455]704               ENDIF
705
706               IF( cffy > zep ) THEN
[2528]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) )
[455]709               ELSE
[2528]710                  dzv(ji,jj,jk) = 0._wp
[455]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      !----------------------------------------------------------------------------------
[2528]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  )
[455]723
[2528]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)
[455]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.
[2528]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) )   &
[3764]742               &                              / ( fsde3w(ji,jj,2) - fsde3w(ji,jj,1) )  )
[455]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
[3764]805         DO jj = 2, jpjm1
[455]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
[503]820      !
[3764]821      CALL wrk_dealloc( jpi, jpj, jpk, dzx  , dzy  , dzz  , dzu  , dzv  , dzw   )
822      CALL wrk_dealloc( jpi, jpj, jpk, drhox, drhoy, drhoz, drhou, drhov, drhow )
823      CALL wrk_dealloc( jpi, jpj, jpk, rho_i, rho_j, rho_k,  zhpi,  zhpj        )
[2715]824      !
[455]825   END SUBROUTINE hpg_djc
826
827
[3294]828   SUBROUTINE hpg_prj( kt )
[455]829      !!---------------------------------------------------------------------
[3294]830      !!                  ***  ROUTINE hpg_prj  ***
[455]831      !!
[3294]832      !! ** Method  :   s-coordinate case.
833      !!      A Pressure-Jacobian horizontal pressure gradient method
834      !!      based on the constrained cubic-spline interpolation for
835      !!      all vertical coordinate systems
[455]836      !!
[3294]837      !! ** Action : - Update (ua,va) with the now hydrastatic pressure trend
838      !!             - Save the trend (l_trddyn=T)
839      !!
[455]840      !!----------------------------------------------------------------------
[3294]841      INTEGER, PARAMETER  :: polynomial_type = 1    ! 1: cubic spline, 2: linear
842      INTEGER, INTENT(in) ::   kt                   ! ocean time-step index
[503]843      !!
[3294]844      INTEGER  ::   ji, jj, jk, jkk                 ! dummy loop indices
845      REAL(wp) ::   zcoef0, znad                    ! temporary scalars
[503]846      !!
[3294]847      !! The local variables for the correction term
848      INTEGER  :: jk1, jis, jid, jjs, jjd
849      REAL(wp) :: zuijk, zvijk, zpwes, zpwed, zpnss, zpnsd, zdeps
[3764]850      REAL(wp) :: zrhdt1
[3294]851      REAL(wp) :: zdpdx1, zdpdx2, zdpdy1, zdpdy2
[3764]852      REAL(wp), POINTER, DIMENSION(:,:,:) ::   zdept, zrhh
[3294]853      REAL(wp), POINTER, DIMENSION(:,:,:) ::   zhpi, zu, zv, fsp, xsp, asp, bsp, csp, dsp
[455]854      !!----------------------------------------------------------------------
[3294]855      !
[3764]856      CALL wrk_alloc( jpi,jpj,jpk, zhpi, zu, zv, fsp, xsp, asp, bsp, csp, dsp )
857      CALL wrk_alloc( jpi,jpj,jpk, zdept, zrhh )
[3294]858      !
[455]859      IF( kt == nit000 ) THEN
860         IF(lwp) WRITE(numout,*)
[3294]861         IF(lwp) WRITE(numout,*) 'dyn:hpg_prj : hydrostatic pressure gradient trend'
862         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   s-coordinate case, cubic spline pressure Jacobian'
[3]863      ENDIF
864
[3294]865      !!----------------------------------------------------------------------
866      ! Local constant initialization
[3764]867      zcoef0 = - grav
[3294]868      znad = 0.0_wp
869      IF( lk_vvl ) znad = 1._wp
[3]870
[3294]871      ! Clean 3-D work arrays
872      zhpi(:,:,:) = 0._wp
873      zrhh(:,:,:) = rhd(:,:,:)
[3764]874
[3294]875      ! Preparing vertical density profile "zrhh(:,:,:)" for hybrid-sco coordinate
876      DO jj = 1, jpj
[3764]877        DO ji = 1, jpi
[3294]878          jk = mbathy(ji,jj)
879          IF( jk <= 0 ) THEN; zrhh(ji,jj,:) = 0._wp
880          ELSE IF(jk == 1) THEN; zrhh(ji,jj, jk+1:jpk) = rhd(ji,jj,jk)
881          ELSE IF(jk < jpkm1) THEN
882             DO jkk = jk+1, jpk
883                zrhh(ji,jj,jkk) = interp1(fsde3w(ji,jj,jkk),   fsde3w(ji,jj,jkk-1), &
884                                         fsde3w(ji,jj,jkk-2), rhd(ji,jj,jkk-1), rhd(ji,jj,jkk-2))
[3764]885             END DO
[3294]886          ENDIF
887        END DO
888      END DO
[3]889
[3632]890      ! Transfer the depth of "T(:,:,:)" to vertical coordinate "zdept(:,:,:)"
891      DO jj = 1, jpj;   DO ji = 1, jpi
892          zdept(ji,jj,1) = 0.5_wp * fse3w(ji,jj,1) - sshn(ji,jj) * znad
893      END DO        ;   END DO
[455]894
[3632]895      DO jk = 2, jpk;   DO jj = 1, jpj;   DO ji = 1, jpi
896          zdept(ji,jj,jk) = zdept(ji,jj,jk-1) + fse3w(ji,jj,jk)
897      END DO        ;   END DO        ;   END DO
[455]898
[3632]899      fsp(:,:,:) = zrhh(:,:,:)
900      xsp(:,:,:) = zdept(:,:,:)
901
[3764]902      ! Construct the vertical density profile with the
[3294]903      ! constrained cubic spline interpolation
904      ! rho(z) = asp + bsp*z + csp*z^2 + dsp*z^3
[3764]905      CALL cspline(fsp,xsp,asp,bsp,csp,dsp,polynomial_type)
[3294]906
907      ! Integrate the hydrostatic pressure "zhpi(:,:,:)" at "T(ji,jj,1)"
908      DO jj = 2, jpj
[3764]909        DO ji = 2, jpi
[3632]910          zrhdt1 = zrhh(ji,jj,1) - interp3(zdept(ji,jj,1),asp(ji,jj,1), &
[3294]911                                         bsp(ji,jj,1),   csp(ji,jj,1), &
[3632]912                                         dsp(ji,jj,1) ) * 0.25_wp * fse3w(ji,jj,1)
[3294]913
914          ! assuming linear profile across the top half surface layer
[3764]915          zhpi(ji,jj,1) =  0.5_wp * fse3w(ji,jj,1) * zrhdt1
[3294]916        END DO
[455]917      END DO
918
[3294]919      ! Calculate the pressure "zhpi(:,:,:)" at "T(ji,jj,2:jpkm1)"
[3764]920      DO jk = 2, jpkm1
921        DO jj = 2, jpj
[3294]922          DO ji = 2, jpi
923            zhpi(ji,jj,jk) = zhpi(ji,jj,jk-1) +                          &
[3632]924                             integ_spline(zdept(ji,jj,jk-1), zdept(ji,jj,jk),&
[3294]925                                    asp(ji,jj,jk-1),    bsp(ji,jj,jk-1), &
926                                    csp(ji,jj,jk-1),    dsp(ji,jj,jk-1))
927          END DO
928        END DO
[455]929      END DO
930
[3294]931      ! Z coordinate of U(ji,jj,1:jpkm1) and V(ji,jj,1:jpkm1)
[3764]932      DO jj = 2, jpjm1
933        DO ji = 2, jpim1
[4292]934          zu(ji,jj,1) = - ( fse3u(ji,jj,1) - sshn(ji,jj) * znad)    ! probable bug: changed from sshu_n for ztilde compilation
935          zv(ji,jj,1) = - ( fse3v(ji,jj,1) - sshn(ji,jj) * znad)    ! probable bug: changed from sshv_n for ztilde compilation
[3294]936        END DO
[455]937      END DO
938
[3764]939      DO jk = 2, jpkm1
940        DO jj = 2, jpjm1
941          DO ji = 2, jpim1
[3294]942            zu(ji,jj,jk) = zu(ji,jj,jk-1)- fse3u(ji,jj,jk)
943            zv(ji,jj,jk) = zv(ji,jj,jk-1)- fse3v(ji,jj,jk)
944          END DO
945        END DO
[455]946      END DO
[3764]947
948      DO jk = 1, jpkm1
949        DO jj = 2, jpjm1
950          DO ji = 2, jpim1
[3294]951            zu(ji,jj,jk) = zu(ji,jj,jk) + 0.5_wp * fse3u(ji,jj,jk)
952            zv(ji,jj,jk) = zv(ji,jj,jk) + 0.5_wp * fse3v(ji,jj,jk)
953          END DO
954        END DO
955      END DO
[455]956
[3632]957      DO jk = 1, jpkm1
958        DO jj = 2, jpjm1
959          DO ji = 2, jpim1
960            zu(ji,jj,jk) = min(zu(ji,jj,jk), max(-zdept(ji,jj,jk), -zdept(ji+1,jj,jk)))
961            zu(ji,jj,jk) = max(zu(ji,jj,jk), min(-zdept(ji,jj,jk), -zdept(ji+1,jj,jk)))
962            zv(ji,jj,jk) = min(zv(ji,jj,jk), max(-zdept(ji,jj,jk), -zdept(ji,jj+1,jk)))
963            zv(ji,jj,jk) = max(zv(ji,jj,jk), min(-zdept(ji,jj,jk), -zdept(ji,jj+1,jk)))
964          END DO
965        END DO
966      END DO
967
968
[3764]969      DO jk = 1, jpkm1
970        DO jj = 2, jpjm1
971          DO ji = 2, jpim1
[3294]972            zpwes = 0._wp; zpwed = 0._wp
973            zpnss = 0._wp; zpnsd = 0._wp
974            zuijk = zu(ji,jj,jk)
975            zvijk = zv(ji,jj,jk)
976
977            !!!!!     for u equation
978            IF( jk <= mbku(ji,jj) ) THEN
[3632]979               IF( -zdept(ji+1,jj,jk) >= -zdept(ji,jj,jk) ) THEN
[3294]980                 jis = ji + 1; jid = ji
981               ELSE
982                 jis = ji;     jid = ji +1
983               ENDIF
984
985               ! integrate the pressure on the shallow side
[3764]986               jk1 = jk
[3632]987               DO WHILE ( -zdept(jis,jj,jk1) > zuijk )
[3294]988                 IF( jk1 == mbku(ji,jj) ) THEN
[3632]989                   zuijk = -zdept(jis,jj,jk1)
[3294]990                   EXIT
991                 ENDIF
[3632]992                 zdeps = MIN(zdept(jis,jj,jk1+1), -zuijk)
[3764]993                 zpwes = zpwes +                                    &
[3632]994                      integ_spline(zdept(jis,jj,jk1), zdeps,            &
[3294]995                             asp(jis,jj,jk1),    bsp(jis,jj,jk1), &
996                             csp(jis,jj,jk1),    dsp(jis,jj,jk1))
997                 jk1 = jk1 + 1
998               END DO
[3764]999
[3294]1000               ! integrate the pressure on the deep side
[3764]1001               jk1 = jk
[3632]1002               DO WHILE ( -zdept(jid,jj,jk1) < zuijk )
[3294]1003                 IF( jk1 == 1 ) THEN
[3632]1004                   zdeps = zdept(jid,jj,1) + MIN(zuijk, sshn(jid,jj)*znad)
1005                   zrhdt1 = zrhh(jid,jj,1) - interp3(zdept(jid,jj,1), asp(jid,jj,1), &
1006                                                     bsp(jid,jj,1),   csp(jid,jj,1), &
1007                                                     dsp(jid,jj,1)) * zdeps
1008                   zpwed  = zpwed + 0.5_wp * (zrhh(jid,jj,1) + zrhdt1) * zdeps
[3294]1009                   EXIT
1010                 ENDIF
[3632]1011                 zdeps = MAX(zdept(jid,jj,jk1-1), -zuijk)
[3764]1012                 zpwed = zpwed +                                        &
[3632]1013                        integ_spline(zdeps,              zdept(jid,jj,jk1), &
[3294]1014                               asp(jid,jj,jk1-1), bsp(jid,jj,jk1-1),  &
1015                               csp(jid,jj,jk1-1), dsp(jid,jj,jk1-1) )
1016                 jk1 = jk1 - 1
1017               END DO
[3764]1018
[3294]1019               ! update the momentum trends in u direction
1020
1021               zdpdx1 = zcoef0 / e1u(ji,jj) * (zhpi(ji+1,jj,jk) - zhpi(ji,jj,jk))
1022               IF( lk_vvl ) THEN
[3764]1023                 zdpdx2 = zcoef0 / e1u(ji,jj) * &
1024                         ( REAL(jis-jid, wp) * (zpwes + zpwed) + (sshn(ji+1,jj)-sshn(ji,jj)) )
[3294]1025                ELSE
[3764]1026                 zdpdx2 = zcoef0 / e1u(ji,jj) * REAL(jis-jid, wp) * (zpwes + zpwed)
[3294]1027               ENDIF
1028
1029               ua(ji,jj,jk) = ua(ji,jj,jk) + (zdpdx1 + zdpdx2) * &
1030               &           umask(ji,jj,jk) * tmask(ji,jj,jk) * tmask(ji+1,jj,jk)
1031            ENDIF
[3764]1032
[3294]1033            !!!!!     for v equation
1034            IF( jk <= mbkv(ji,jj) ) THEN
[3632]1035               IF( -zdept(ji,jj+1,jk) >= -zdept(ji,jj,jk) ) THEN
[3294]1036                 jjs = jj + 1; jjd = jj
1037               ELSE
1038                 jjs = jj    ; jjd = jj + 1
1039               ENDIF
1040
1041               ! integrate the pressure on the shallow side
[3764]1042               jk1 = jk
[3632]1043               DO WHILE ( -zdept(ji,jjs,jk1) > zvijk )
[3294]1044                 IF( jk1 == mbkv(ji,jj) ) THEN
[3632]1045                   zvijk = -zdept(ji,jjs,jk1)
[3294]1046                   EXIT
1047                 ENDIF
[3632]1048                 zdeps = MIN(zdept(ji,jjs,jk1+1), -zvijk)
[3764]1049                 zpnss = zpnss +                                      &
[3632]1050                        integ_spline(zdept(ji,jjs,jk1), zdeps,            &
[3294]1051                               asp(ji,jjs,jk1),    bsp(ji,jjs,jk1), &
1052                               csp(ji,jjs,jk1),    dsp(ji,jjs,jk1) )
1053                 jk1 = jk1 + 1
1054               END DO
[3764]1055
[3294]1056               ! integrate the pressure on the deep side
[3764]1057               jk1 = jk
[3632]1058               DO WHILE ( -zdept(ji,jjd,jk1) < zvijk )
[3294]1059                 IF( jk1 == 1 ) THEN
[3632]1060                   zdeps = zdept(ji,jjd,1) + MIN(zvijk, sshn(ji,jjd)*znad)
1061                   zrhdt1 = zrhh(ji,jjd,1) - interp3(zdept(ji,jjd,1), asp(ji,jjd,1), &
1062                                                     bsp(ji,jjd,1),   csp(ji,jjd,1), &
1063                                                     dsp(ji,jjd,1) ) * zdeps
1064                   zpnsd  = zpnsd + 0.5_wp * (zrhh(ji,jjd,1) + zrhdt1) * zdeps
[3294]1065                   EXIT
1066                 ENDIF
[3632]1067                 zdeps = MAX(zdept(ji,jjd,jk1-1), -zvijk)
[3764]1068                 zpnsd = zpnsd +                                        &
[3632]1069                        integ_spline(zdeps,              zdept(ji,jjd,jk1), &
[3294]1070                               asp(ji,jjd,jk1-1), bsp(ji,jjd,jk1-1), &
1071                               csp(ji,jjd,jk1-1), dsp(ji,jjd,jk1-1) )
1072                 jk1 = jk1 - 1
1073               END DO
1074
[3764]1075
[3294]1076               ! update the momentum trends in v direction
1077
1078               zdpdy1 = zcoef0 / e2v(ji,jj) * (zhpi(ji,jj+1,jk) - zhpi(ji,jj,jk))
1079               IF( lk_vvl ) THEN
1080                   zdpdy2 = zcoef0 / e2v(ji,jj) * &
[3764]1081                           ( REAL(jjs-jjd, wp) * (zpnss + zpnsd) + (sshn(ji,jj+1)-sshn(ji,jj)) )
[3294]1082               ELSE
[3764]1083                   zdpdy2 = zcoef0 / e2v(ji,jj) * REAL(jjs-jjd, wp) * (zpnss + zpnsd )
[3294]1084               ENDIF
1085
1086               va(ji,jj,jk) = va(ji,jj,jk) + (zdpdy1 + zdpdy2)*&
1087               &              vmask(ji,jj,jk)*tmask(ji,jj,jk)*tmask(ji,jj+1,jk)
1088            ENDIF
1089
[3764]1090
[3294]1091           END DO
1092        END DO
[455]1093      END DO
[503]1094      !
[3764]1095      CALL wrk_dealloc( jpi,jpj,jpk, zhpi, zu, zv, fsp, xsp, asp, bsp, csp, dsp )
1096      CALL wrk_dealloc( jpi,jpj,jpk, zdept, zrhh )
[2715]1097      !
[3294]1098   END SUBROUTINE hpg_prj
[455]1099
[3294]1100   SUBROUTINE cspline(fsp, xsp, asp, bsp, csp, dsp, polynomial_type)
1101      !!----------------------------------------------------------------------
1102      !!                 ***  ROUTINE cspline  ***
[3764]1103      !!
[3294]1104      !! ** Purpose :   constrained cubic spline interpolation
[3764]1105      !!
1106      !! ** Method  :   f(x) = asp + bsp*x + csp*x^2 + dsp*x^3
[3294]1107      !! Reference: CJC Kruger, Constrained Cubic Spline Interpoltation
1108      !!
1109      !!----------------------------------------------------------------------
1110      IMPLICIT NONE
1111      REAL(wp), DIMENSION(:,:,:), INTENT(in)  :: fsp, xsp           ! value and coordinate
[3764]1112      REAL(wp), DIMENSION(:,:,:), INTENT(out) :: asp, bsp, csp, dsp ! coefficients of
[3294]1113                                                                    ! the interpoated function
[3764]1114      INTEGER, INTENT(in) :: polynomial_type                        ! 1: cubic spline
[3294]1115                                                                    ! 2: Linear
1116
[3764]1117      ! Local Variables
[3294]1118      INTEGER  ::   ji, jj, jk                 ! dummy loop indices
1119      INTEGER  ::   jpi, jpj, jpkm1
1120      REAL(wp) ::   zdf1, zdf2, zddf1, zddf2, ztmp1, ztmp2, zdxtmp
1121      REAL(wp) ::   zdxtmp1, zdxtmp2, zalpha
1122      REAL(wp) ::   zdf(size(fsp,3))
1123      !!----------------------------------------------------------------------
1124
1125      jpi   = size(fsp,1)
1126      jpj   = size(fsp,2)
1127      jpkm1 = size(fsp,3) - 1
1128
[3764]1129
[3294]1130      IF (polynomial_type == 1) THEN     ! Constrained Cubic Spline
1131         DO ji = 1, jpi
1132            DO jj = 1, jpj
[3764]1133           !!Fritsch&Butland's method, 1984 (preferred, but more computation)
[3294]1134           !    DO jk = 2, jpkm1-1
[3764]1135           !       zdxtmp1 = xsp(ji,jj,jk)   - xsp(ji,jj,jk-1)
1136           !       zdxtmp2 = xsp(ji,jj,jk+1) - xsp(ji,jj,jk)
[3294]1137           !       zdf1    = ( fsp(ji,jj,jk)   - fsp(ji,jj,jk-1) ) / zdxtmp1
1138           !       zdf2    = ( fsp(ji,jj,jk+1) - fsp(ji,jj,jk)   ) / zdxtmp2
1139           !
1140           !       zalpha = ( zdxtmp1 + 2._wp * zdxtmp2 ) / ( zdxtmp1 + zdxtmp2 ) / 3._wp
[3764]1141           !
[3294]1142           !       IF(zdf1 * zdf2 <= 0._wp) THEN
1143           !           zdf(jk) = 0._wp
1144           !       ELSE
1145           !         zdf(jk) = zdf1 * zdf2 / ( ( 1._wp - zalpha ) * zdf1 + zalpha * zdf2 )
1146           !       ENDIF
1147           !    END DO
[3764]1148
[3294]1149           !!Simply geometric average
1150               DO jk = 2, jpkm1-1
1151                  zdf1 = (fsp(ji,jj,jk) - fsp(ji,jj,jk-1)) / (xsp(ji,jj,jk) - xsp(ji,jj,jk-1))
1152                  zdf2 = (fsp(ji,jj,jk+1) - fsp(ji,jj,jk)) / (xsp(ji,jj,jk+1) - xsp(ji,jj,jk))
[3764]1153
[3294]1154                  IF(zdf1 * zdf2 <= 0._wp) THEN
1155                     zdf(jk) = 0._wp
1156                  ELSE
1157                     zdf(jk) = 2._wp * zdf1 * zdf2 / (zdf1 + zdf2)
1158                  ENDIF
1159               END DO
[3764]1160
[3294]1161               zdf(1)     = 1.5_wp * ( fsp(ji,jj,2) - fsp(ji,jj,1) ) / &
1162                          &          ( xsp(ji,jj,2) - xsp(ji,jj,1) ) -  0.5_wp * zdf(2)
1163               zdf(jpkm1) = 1.5_wp * ( fsp(ji,jj,jpkm1) - fsp(ji,jj,jpkm1-1) ) / &
1164                          &          ( xsp(ji,jj,jpkm1) - xsp(ji,jj,jpkm1-1) ) - &
1165                          & 0.5_wp * zdf(jpkm1 - 1)
[3764]1166
[3294]1167               DO jk = 1, jpkm1 - 1
[3764]1168                 zdxtmp = xsp(ji,jj,jk+1) - xsp(ji,jj,jk)
[3294]1169                 ztmp1  = (zdf(jk+1) + 2._wp * zdf(jk)) / zdxtmp
1170                 ztmp2  =  6._wp * (fsp(ji,jj,jk+1) - fsp(ji,jj,jk)) / zdxtmp / zdxtmp
[3764]1171                 zddf1  = -2._wp * ztmp1 + ztmp2
[3294]1172                 ztmp1  = (2._wp * zdf(jk+1) + zdf(jk)) / zdxtmp
[3764]1173                 zddf2  =  2._wp * ztmp1 - ztmp2
1174
[3294]1175                 dsp(ji,jj,jk) = (zddf2 - zddf1) / 6._wp / zdxtmp
1176                 csp(ji,jj,jk) = ( xsp(ji,jj,jk+1) * zddf1 - xsp(ji,jj,jk)*zddf2 ) / 2._wp / zdxtmp
[3764]1177                 bsp(ji,jj,jk) = ( fsp(ji,jj,jk+1) - fsp(ji,jj,jk) ) / zdxtmp - &
[3294]1178                               & csp(ji,jj,jk) * ( xsp(ji,jj,jk+1) + xsp(ji,jj,jk) ) - &
1179                               & dsp(ji,jj,jk) * ((xsp(ji,jj,jk+1) + xsp(ji,jj,jk))**2 - &
1180                               &                   xsp(ji,jj,jk+1) * xsp(ji,jj,jk))
1181                 asp(ji,jj,jk) = fsp(ji,jj,jk) - xsp(ji,jj,jk) * (bsp(ji,jj,jk) + &
1182                               &                (xsp(ji,jj,jk) * (csp(ji,jj,jk) + &
1183                               &                 dsp(ji,jj,jk) * xsp(ji,jj,jk))))
1184               END DO
1185            END DO
1186         END DO
[3764]1187
[3294]1188      ELSE IF (polynomial_type == 2) THEN     ! Linear
1189         DO ji = 1, jpi
1190            DO jj = 1, jpj
1191               DO jk = 1, jpkm1-1
[3764]1192                  zdxtmp =xsp(ji,jj,jk+1) - xsp(ji,jj,jk)
[3294]1193                  ztmp1 = fsp(ji,jj,jk+1) - fsp(ji,jj,jk)
[3764]1194
[3294]1195                  dsp(ji,jj,jk) = 0._wp
1196                  csp(ji,jj,jk) = 0._wp
1197                  bsp(ji,jj,jk) = ztmp1 / zdxtmp
1198                  asp(ji,jj,jk) = fsp(ji,jj,jk) - bsp(ji,jj,jk) * xsp(ji,jj,jk)
1199               END DO
1200            END DO
1201         END DO
1202
1203      ELSE
1204           CALL ctl_stop( 'invalid polynomial type in cspline' )
1205      ENDIF
1206
[3764]1207
[3294]1208   END SUBROUTINE cspline
1209
1210
[3764]1211   FUNCTION interp1(x, xl, xr, fl, fr)  RESULT(f)
[3294]1212      !!----------------------------------------------------------------------
1213      !!                 ***  ROUTINE interp1  ***
[3764]1214      !!
[3294]1215      !! ** Purpose :   1-d linear interpolation
[3764]1216      !!
1217      !! ** Method  :
[3294]1218      !!                interpolation is straight forward
[3764]1219      !!                extrapolation is also permitted (no value limit)
[3294]1220      !!
1221      !!----------------------------------------------------------------------
1222      IMPLICIT NONE
[3764]1223      REAL(wp), INTENT(in) ::  x, xl, xr, fl, fr
[3294]1224      REAL(wp)             ::  f ! result of the interpolation (extrapolation)
1225      REAL(wp)             ::  zdeltx
1226      !!----------------------------------------------------------------------
1227
1228      zdeltx = xr - xl
1229      IF(abs(zdeltx) <= 10._wp * EPSILON(x)) THEN
1230        f = 0.5_wp * (fl + fr)
1231      ELSE
1232        f = ( (x - xl ) * fr - ( x - xr ) * fl ) / zdeltx
1233      ENDIF
[3764]1234
[3294]1235   END FUNCTION interp1
1236
[3764]1237   FUNCTION interp2(x, a, b, c, d)  RESULT(f)
[3294]1238      !!----------------------------------------------------------------------
1239      !!                 ***  ROUTINE interp1  ***
[3764]1240      !!
[3294]1241      !! ** Purpose :   1-d constrained cubic spline interpolation
[3764]1242      !!
[3294]1243      !! ** Method  :  cubic spline interpolation
1244      !!
1245      !!----------------------------------------------------------------------
1246      IMPLICIT NONE
[3764]1247      REAL(wp), INTENT(in) ::  x, a, b, c, d
[3294]1248      REAL(wp)             ::  f ! value from the interpolation
1249      !!----------------------------------------------------------------------
1250
[3764]1251      f = a + x* ( b + x * ( c + d * x ) )
[3294]1252
1253   END FUNCTION interp2
1254
1255
[3764]1256   FUNCTION interp3(x, a, b, c, d)  RESULT(f)
[3294]1257      !!----------------------------------------------------------------------
1258      !!                 ***  ROUTINE interp1  ***
[3764]1259      !!
[3294]1260      !! ** Purpose :   Calculate the first order of deriavtive of
1261      !!                a cubic spline function y=a+b*x+c*x^2+d*x^3
[3764]1262      !!
[3294]1263      !! ** Method  :   f=dy/dx=b+2*c*x+3*d*x^2
1264      !!
1265      !!----------------------------------------------------------------------
1266      IMPLICIT NONE
[3764]1267      REAL(wp), INTENT(in) ::  x, a, b, c, d
[3294]1268      REAL(wp)             ::  f ! value from the interpolation
1269      !!----------------------------------------------------------------------
1270
1271      f = b + x * ( 2._wp * c + 3._wp * d * x)
1272
1273   END FUNCTION interp3
1274
[3764]1275
1276   FUNCTION integ_spline(xl, xr, a, b, c, d)  RESULT(f)
[3294]1277      !!----------------------------------------------------------------------
1278      !!                 ***  ROUTINE interp1  ***
[3764]1279      !!
[3294]1280      !! ** Purpose :   1-d constrained cubic spline integration
1281      !!
[3764]1282      !! ** Method  :  integrate polynomial a+bx+cx^2+dx^3 from xl to xr
1283      !!
[3294]1284      !!----------------------------------------------------------------------
1285      IMPLICIT NONE
[3764]1286      REAL(wp), INTENT(in) ::  xl, xr, a, b, c, d
1287      REAL(wp)             ::  za1, za2, za3
[3294]1288      REAL(wp)             ::  f                   ! integration result
1289      !!----------------------------------------------------------------------
1290
[3764]1291      za1 = 0.5_wp * b
1292      za2 = c / 3.0_wp
1293      za3 = 0.25_wp * d
[3294]1294
1295      f  = xr * ( a + xr * ( za1 + xr * ( za2 + za3 * xr ) ) ) - &
1296         & xl * ( a + xl * ( za1 + xl * ( za2 + za3 * xl ) ) )
1297
[3632]1298   END FUNCTION integ_spline
[3294]1299
1300
[3]1301   !!======================================================================
1302END MODULE dynhpg
[3632]1303
Note: See TracBrowser for help on using the repository browser.