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.
domvvl.F90 in NEMO/branches/2020/dev_r14116_HPC-04_mcastril_Mixed_Precision_implementation_final/src/OCE/DOM – NEMO

source: NEMO/branches/2020/dev_r14116_HPC-04_mcastril_Mixed_Precision_implementation_final/src/OCE/DOM/domvvl.F90 @ 14219

Last change on this file since 14219 was 14219, checked in by mcastril, 4 years ago

Add Mixed Precision support by Oriol Tintó

  • Property svn:keywords set to Id
File size: 52.8 KB
Line 
1MODULE domvvl
2   !!======================================================================
3   !!                       ***  MODULE domvvl   ***
4   !! Ocean :
5   !!======================================================================
6   !! History :  2.0  !  2006-06  (B. Levier, L. Marie)  original code
7   !!            3.1  !  2009-02  (G. Madec, M. Leclair, R. Benshila)  pure z* coordinate
8   !!            3.3  !  2011-10  (M. Leclair) totally rewrote domvvl: vvl option includes z_star and z_tilde coordinates
9   !!            3.6  !  2014-11  (P. Mathiot) add ice shelf capability
10   !!            4.1  !  2019-08  (A. Coward, D. Storkey) rename dom_vvl_sf_swp -> dom_vvl_sf_update for new timestepping
11   !!             -   !  2020-02  (G. Madec, S. Techene) introduce ssh to h0 ratio
12   !!----------------------------------------------------------------------
13
14   USE oce             ! ocean dynamics and tracers
15   USE phycst          ! physical constant
16   USE dom_oce         ! ocean space and time domain
17   USE sbc_oce         ! ocean surface boundary condition
18   USE wet_dry         ! wetting and drying
19   USE usrdef_istate   ! user defined initial state (wad only)
20   USE restart         ! ocean restart
21   !
22   USE in_out_manager  ! I/O manager
23   USE iom             ! I/O manager library
24   USE lib_mpp         ! distributed memory computing library
25   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
26   USE timing          ! Timing
27
28   IMPLICIT NONE
29   PRIVATE
30
31   !                                                      !!* Namelist nam_vvl
32   LOGICAL , PUBLIC :: ln_vvl_zstar           = .FALSE.    ! zstar  vertical coordinate
33   LOGICAL , PUBLIC :: ln_vvl_ztilde          = .FALSE.    ! ztilde vertical coordinate
34   LOGICAL , PUBLIC :: ln_vvl_layer           = .FALSE.    ! level  vertical coordinate
35   LOGICAL , PUBLIC :: ln_vvl_ztilde_as_zstar = .FALSE.    ! ztilde vertical coordinate
36   LOGICAL , PUBLIC :: ln_vvl_zstar_at_eqtor  = .FALSE.    ! ztilde vertical coordinate
37   LOGICAL , PUBLIC :: ln_vvl_kepe            = .FALSE.    ! kinetic/potential energy transfer
38   !                                                       ! conservation: not used yet
39   REAL(wp)         :: rn_ahe3                             ! thickness diffusion coefficient
40   REAL(wp)         :: rn_rst_e3t                          ! ztilde to zstar restoration timescale [days]
41   REAL(wp)         :: rn_lf_cutoff                        ! cutoff frequency for low-pass filter  [days]
42   REAL(wp)         :: rn_zdef_max                         ! maximum fractional e3t deformation
43   LOGICAL , PUBLIC :: ln_vvl_dbg = .FALSE.                ! debug control prints
44
45   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: un_td, vn_td                ! thickness diffusion transport
46   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: hdiv_lf                     ! low frequency part of hz divergence
47   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: tilde_e3t_b, tilde_e3t_n    ! baroclinic scale factors
48   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: tilde_e3t_a, dtilde_e3t_a   ! baroclinic scale factors
49   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:)   :: frq_rst_e3t                 ! retoring period for scale factors
50   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:)   :: frq_rst_hdv                 ! retoring period for low freq. divergence
51
52#if defined key_qco   ||   defined key_linssh
53   !!----------------------------------------------------------------------
54   !!   'key_qco'                        Quasi-Eulerian vertical coordinate
55   !!       OR         EMPTY MODULE
56   !!   'key_linssh'                        Fix in time vertical coordinate
57   !!----------------------------------------------------------------------
58#else
59   !!----------------------------------------------------------------------
60   !!   Default key      Old management of time varying vertical coordinate
61   !!----------------------------------------------------------------------
62
63   !!----------------------------------------------------------------------
64   !!   dom_vvl_init     : define initial vertical scale factors, depths and column thickness
65   !!   dom_vvl_sf_nxt   : Compute next vertical scale factors
66   !!   dom_vvl_sf_update   : Swap vertical scale factors and update the vertical grid
67   !!   dom_vvl_interpol : Interpolate vertical scale factors from one grid point to another
68   !!   dom_vvl_rst      : read/write restart file
69   !!   dom_vvl_ctl      : Check the vvl options
70   !!----------------------------------------------------------------------
71
72   PUBLIC  dom_vvl_init       ! called by domain.F90
73   PUBLIC  dom_vvl_zgr        ! called by isfcpl.F90
74   PUBLIC  dom_vvl_sf_nxt     ! called by step.F90
75   PUBLIC  dom_vvl_sf_update  ! called by step.F90
76   PUBLIC  dom_vvl_interpol   ! called by dynnxt.F90
77
78   !! * Substitutions
79#  include "do_loop_substitute.h90"
80#  include "single_precision_substitute.h90"
81   !!----------------------------------------------------------------------
82   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
83   !! $Id$
84   !! Software governed by the CeCILL license (see ./LICENSE)
85   !!----------------------------------------------------------------------
86CONTAINS
87
88   INTEGER FUNCTION dom_vvl_alloc()
89      !!----------------------------------------------------------------------
90      !!                ***  FUNCTION dom_vvl_alloc  ***
91      !!----------------------------------------------------------------------
92      IF( ln_vvl_zstar )   dom_vvl_alloc = 0
93      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN
94         ALLOCATE( tilde_e3t_b(jpi,jpj,jpk)  , tilde_e3t_n(jpi,jpj,jpk) , tilde_e3t_a(jpi,jpj,jpk) ,   &
95            &      dtilde_e3t_a(jpi,jpj,jpk) , un_td  (jpi,jpj,jpk)     , vn_td  (jpi,jpj,jpk)     ,   &
96            &      STAT = dom_vvl_alloc        )
97         CALL mpp_sum ( 'domvvl', dom_vvl_alloc )
98         IF( dom_vvl_alloc /= 0 )   CALL ctl_stop( 'STOP', 'dom_vvl_alloc: failed to allocate arrays' )
99         un_td = 0._wp
100         vn_td = 0._wp
101      ENDIF
102      IF( ln_vvl_ztilde ) THEN
103         ALLOCATE( frq_rst_e3t(jpi,jpj) , frq_rst_hdv(jpi,jpj) , hdiv_lf(jpi,jpj,jpk) , STAT= dom_vvl_alloc )
104         CALL mpp_sum ( 'domvvl', dom_vvl_alloc )
105         IF( dom_vvl_alloc /= 0 )   CALL ctl_stop( 'STOP', 'dom_vvl_alloc: failed to allocate arrays' )
106      ENDIF
107      !
108   END FUNCTION dom_vvl_alloc
109
110
111   SUBROUTINE dom_vvl_init( Kbb, Kmm, Kaa )
112      !!----------------------------------------------------------------------
113      !!                ***  ROUTINE dom_vvl_init  ***
114      !!
115      !! ** Purpose :  Initialization of all scale factors, depths
116      !!               and water column heights
117      !!
118      !! ** Method  :  - use restart file and/or initialize
119      !!               - interpolate scale factors
120      !!
121      !! ** Action  : - e3t_(n/b) and tilde_e3t_(n/b)
122      !!              - Regrid: e3[u/v](:,:,:,Kmm)
123      !!                        e3[u/v](:,:,:,Kmm)
124      !!                        e3w(:,:,:,Kmm)
125      !!                        e3[u/v]w_b
126      !!                        e3[u/v]w_n
127      !!                        gdept(:,:,:,Kmm), gdepw(:,:,:,Kmm) and gde3w
128      !!              - h(t/u/v)_0
129      !!              - frq_rst_e3t and frq_rst_hdv
130      !!
131      !! Reference  : Leclair, M., and G. Madec, 2011, Ocean Modelling.
132      !!----------------------------------------------------------------------
133      INTEGER, INTENT(in) :: Kbb, Kmm, Kaa
134      !
135      IF(lwp) WRITE(numout,*)
136      IF(lwp) WRITE(numout,*) 'dom_vvl_init : Variable volume activated'
137      IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~'
138      !
139      CALL dom_vvl_ctl     ! choose vertical coordinate (z_star, z_tilde or layer)
140      !
141      !                    ! Allocate module arrays
142      IF( dom_vvl_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'dom_vvl_init : unable to allocate arrays' )
143      !
144      !                    ! Read or initialize e3t_(b/n), tilde_e3t_(b/n) and hdiv_lf
145      CALL dom_vvl_rst( nit000, Kbb, Kmm, 'READ' )
146      e3t(:,:,jpk,Kaa) = e3t_0(:,:,jpk)  ! last level always inside the sea floor set one for all
147      !
148      CALL dom_vvl_zgr(Kbb, Kmm, Kaa) ! interpolation scale factor, depth and water column
149      !
150   END SUBROUTINE dom_vvl_init
151
152
153   SUBROUTINE dom_vvl_zgr(Kbb, Kmm, Kaa)
154      !!----------------------------------------------------------------------
155      !!                ***  ROUTINE dom_vvl_init  ***
156      !!
157      !! ** Purpose :  Interpolation of all scale factors,
158      !!               depths and water column heights
159      !!
160      !! ** Method  :  - interpolate scale factors
161      !!
162      !! ** Action  : - e3t_(n/b) and tilde_e3t_(n/b)
163      !!              - Regrid: e3(u/v)_n
164      !!                        e3(u/v)_b
165      !!                        e3w_n
166      !!                        e3(u/v)w_b
167      !!                        e3(u/v)w_n
168      !!                        gdept_n, gdepw_n and gde3w_n
169      !!              - h(t/u/v)_0
170      !!              - frq_rst_e3t and frq_rst_hdv
171      !!
172      !! Reference  : Leclair, M., and G. Madec, 2011, Ocean Modelling.
173      !!----------------------------------------------------------------------
174      INTEGER, INTENT(in) :: Kbb, Kmm, Kaa
175      !!----------------------------------------------------------------------
176      INTEGER ::   ji, jj, jk
177      INTEGER ::   ii0, ii1, ij0, ij1
178      REAL(wp)::   zcoef
179      !!----------------------------------------------------------------------
180      !
181      !                    !== Set of all other vertical scale factors  ==!  (now and before)
182      !                                ! Horizontal interpolation of e3t
183      CALL dom_vvl_interpol( CASTWP(e3t(:,:,:,Kbb)), e3u(:,:,:,Kbb), 'U' )    ! from T to U
184      CALL dom_vvl_interpol( CASTWP(e3t(:,:,:,Kmm)), e3u(:,:,:,Kmm), 'U' )
185      CALL dom_vvl_interpol( CASTWP(e3t(:,:,:,Kbb)), e3v(:,:,:,Kbb), 'V' )    ! from T to V
186      CALL dom_vvl_interpol( CASTWP(e3t(:,:,:,Kmm)), e3v(:,:,:,Kmm), 'V' )
187      CALL dom_vvl_interpol( e3u(:,:,:,Kmm), e3f(:,:,:), 'F' )    ! from U to F
188      !                                ! Vertical interpolation of e3t,u,v
189      CALL dom_vvl_interpol( CASTWP(e3t(:,:,:,Kmm)), e3w (:,:,:,Kmm), 'W'  )  ! from T to W
190      CALL dom_vvl_interpol( CASTWP(e3t(:,:,:,Kbb)), e3w (:,:,:,Kbb), 'W'  )
191      CALL dom_vvl_interpol( e3u(:,:,:,Kmm), e3uw(:,:,:,Kmm), 'UW' )  ! from U to UW
192      CALL dom_vvl_interpol( e3u(:,:,:,Kbb), e3uw(:,:,:,Kbb), 'UW' )
193      CALL dom_vvl_interpol( e3v(:,:,:,Kmm), e3vw(:,:,:,Kmm), 'VW' )  ! from V to UW
194      CALL dom_vvl_interpol( e3v(:,:,:,Kbb), e3vw(:,:,:,Kbb), 'VW' )
195
196      ! We need to define e3[tuv]_a for AGRIF initialisation (should not be a problem for the restartability...)
197      e3t(:,:,:,Kaa) = e3t(:,:,:,Kmm)
198      e3u(:,:,:,Kaa) = e3u(:,:,:,Kmm)
199      e3v(:,:,:,Kaa) = e3v(:,:,:,Kmm)
200      !
201      !                    !==  depth of t and w-point  ==!   (set the isf depth as it is in the initial timestep)
202      gdept(:,:,1,Kmm) = 0.5_wp * e3w(:,:,1,Kmm)       ! reference to the ocean surface (used for MLD and light penetration)
203      gdepw(:,:,1,Kmm) = 0.0_wp
204      gde3w(:,:,1) = gdept(:,:,1,Kmm) - ssh(:,:,Kmm)  ! reference to a common level z=0 for hpg
205      gdept(:,:,1,Kbb) = 0.5_wp * e3w(:,:,1,Kbb)
206      gdepw(:,:,1,Kbb) = 0.0_wp
207      DO_3D( 1, 1, 1, 1, 2, jpk )                     ! vertical sum
208         !    zcoef = tmask - wmask    ! 0 everywhere tmask = wmask, ie everywhere expect at jk = mikt
209         !                             ! 1 everywhere from mbkt to mikt + 1 or 1 (if no isf)
210         !                             ! 0.5 where jk = mikt
211!!gm ???????   BUG ?  gdept(:,:,:,Kmm) as well as gde3w  does not include the thickness of ISF ??
212         zcoef = ( tmask(ji,jj,jk) - wmask(ji,jj,jk) )
213         gdepw(ji,jj,jk,Kmm) = gdepw(ji,jj,jk-1,Kmm) + e3t(ji,jj,jk-1,Kmm)
214         gdept(ji,jj,jk,Kmm) =      zcoef  * ( gdepw(ji,jj,jk  ,Kmm) + 0.5 * e3w(ji,jj,jk,Kmm))  &
215            &                + (1-zcoef) * ( gdept(ji,jj,jk-1,Kmm) +       e3w(ji,jj,jk,Kmm))
216         gde3w(ji,jj,jk) = gdept(ji,jj,jk,Kmm) - ssh(ji,jj,Kmm)
217         gdepw(ji,jj,jk,Kbb) = gdepw(ji,jj,jk-1,Kbb) + e3t(ji,jj,jk-1,Kbb)
218         gdept(ji,jj,jk,Kbb) =      zcoef  * ( gdepw(ji,jj,jk  ,Kbb) + 0.5 * e3w(ji,jj,jk,Kbb))  &
219            &                + (1-zcoef) * ( gdept(ji,jj,jk-1,Kbb) +       e3w(ji,jj,jk,Kbb))
220      END_3D
221      !
222      !                    !==  thickness of the water column  !!   (ocean portion only)
223      ht(:,:) = e3t(:,:,1,Kmm) * tmask(:,:,1)   !!gm  BUG  :  this should be 1/2 * e3w(k=1) ....
224      hu(:,:,Kbb) = e3u(:,:,1,Kbb) * umask(:,:,1)
225      hu(:,:,Kmm) = e3u(:,:,1,Kmm) * umask(:,:,1)
226      hv(:,:,Kbb) = e3v(:,:,1,Kbb) * vmask(:,:,1)
227      hv(:,:,Kmm) = e3v(:,:,1,Kmm) * vmask(:,:,1)
228      DO jk = 2, jpkm1
229         ht(:,:) = ht(:,:) + e3t(:,:,jk,Kmm) * tmask(:,:,jk)
230         hu(:,:,Kbb) = hu(:,:,Kbb) + e3u(:,:,jk,Kbb) * umask(:,:,jk)
231         hu(:,:,Kmm) = hu(:,:,Kmm) + e3u(:,:,jk,Kmm) * umask(:,:,jk)
232         hv(:,:,Kbb) = hv(:,:,Kbb) + e3v(:,:,jk,Kbb) * vmask(:,:,jk)
233         hv(:,:,Kmm) = hv(:,:,Kmm) + e3v(:,:,jk,Kmm) * vmask(:,:,jk)
234      END DO
235      !
236      !                    !==  inverse of water column thickness   ==!   (u- and v- points)
237      r1_hu(:,:,Kbb) = ssumask(:,:) / ( hu(:,:,Kbb) + 1._wp - ssumask(:,:) )    ! _i mask due to ISF
238      r1_hu(:,:,Kmm) = ssumask(:,:) / ( hu(:,:,Kmm) + 1._wp - ssumask(:,:) )
239      r1_hv(:,:,Kbb) = ssvmask(:,:) / ( hv(:,:,Kbb) + 1._wp - ssvmask(:,:) )
240      r1_hv(:,:,Kmm) = ssvmask(:,:) / ( hv(:,:,Kmm) + 1._wp - ssvmask(:,:) )
241
242      !                    !==   z_tilde coordinate case  ==!   (Restoring frequencies)
243      IF( ln_vvl_ztilde ) THEN
244!!gm : idea: add here a READ in a file of custumized restoring frequency
245         !                                   ! Values in days provided via the namelist
246         !                                   ! use rsmall to avoid possible division by zero errors with faulty settings
247         frq_rst_e3t(:,:) = 2._wp * rpi / ( MAX( rn_rst_e3t  , rsmall ) * 86400.0_wp )
248         frq_rst_hdv(:,:) = 2._wp * rpi / ( MAX( rn_lf_cutoff, rsmall ) * 86400.0_wp )
249         !
250         IF( ln_vvl_ztilde_as_zstar ) THEN   ! z-star emulation using z-tile
251            frq_rst_e3t(:,:) = 0._wp               !Ignore namelist settings
252            frq_rst_hdv(:,:) = 1._wp / rn_Dt
253         ENDIF
254         IF ( ln_vvl_zstar_at_eqtor ) THEN   ! use z-star in vicinity of the Equator
255            DO_2D( 1, 1, 1, 1 )
256!!gm  case |gphi| >= 6 degrees is useless   initialized just above by default
257               IF( ABS(gphit(ji,jj)) >= 6.) THEN
258                  ! values outside the equatorial band and transition zone (ztilde)
259                  frq_rst_e3t(ji,jj) =  2.0_wp * rpi / ( MAX( rn_rst_e3t  , rsmall ) * 86400.e0_wp )
260                  frq_rst_hdv(ji,jj) =  2.0_wp * rpi / ( MAX( rn_lf_cutoff, rsmall ) * 86400.e0_wp )
261               ELSEIF( ABS(gphit(ji,jj)) <= 2.5) THEN    ! Equator strip ==> z-star
262                  ! values inside the equatorial band (ztilde as zstar)
263                  frq_rst_e3t(ji,jj) =  0.0_wp
264                  frq_rst_hdv(ji,jj) =  1.0_wp / rn_Dt
265               ELSE                                      ! transition band (2.5 to 6 degrees N/S)
266                  !                                      ! (linearly transition from z-tilde to z-star)
267                  frq_rst_e3t(ji,jj) = 0.0_wp + (frq_rst_e3t(ji,jj)-0.0_wp)*0.5_wp   &
268                     &            * (  1.0_wp - COS( rad*(ABS(gphit(ji,jj))-2.5_wp)  &
269                     &                                          * 180._wp / 3.5_wp ) )
270                  frq_rst_hdv(ji,jj) = (1.0_wp / rn_Dt)                                &
271                     &            + (  frq_rst_hdv(ji,jj)-(1.e0_wp / rn_Dt) )*0.5_wp   &
272                     &            * (  1._wp  - COS( rad*(ABS(gphit(ji,jj))-2.5_wp)  &
273                     &                                          * 180._wp / 3.5_wp ) )
274               ENDIF
275            END_2D
276            IF( cn_cfg == "orca" .OR. cn_cfg == "ORCA" ) THEN
277               IF( nn_cfg == 3 ) THEN   ! ORCA2: Suppress ztilde in the Foxe Basin for ORCA2
278                  ii0 = 103 + nn_hls - 1   ;   ii1 = 111 + nn_hls - 1
279                  ij0 = 128 + nn_hls       ;   ij1 = 135 + nn_hls
280                  frq_rst_e3t( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) ) =  0.0_wp
281                  frq_rst_hdv( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) ) =  1.e0_wp / rn_Dt
282               ENDIF
283            ENDIF
284         ENDIF
285      ENDIF
286      !
287   END SUBROUTINE dom_vvl_zgr
288
289
290   SUBROUTINE dom_vvl_sf_nxt( kt, Kbb, Kmm, Kaa, kcall )
291      !!----------------------------------------------------------------------
292      !!                ***  ROUTINE dom_vvl_sf_nxt  ***
293      !!
294      !! ** Purpose :  - compute the after scale factors used in tra_zdf, dynnxt,
295      !!                 tranxt and dynspg routines
296      !!
297      !! ** Method  :  - z_star case:  Repartition of ssh INCREMENT proportionnaly to the level thickness.
298      !!               - z_tilde_case: after scale factor increment =
299      !!                                    high frequency part of horizontal divergence
300      !!                                  + retsoring towards the background grid
301      !!                                  + thickness difusion
302      !!                               Then repartition of ssh INCREMENT proportionnaly
303      !!                               to the "baroclinic" level thickness.
304      !!
305      !! ** Action  :  - hdiv_lf    : restoring towards full baroclinic divergence in z_tilde case
306      !!               - tilde_e3t_a: after increment of vertical scale factor
307      !!                              in z_tilde case
308      !!               - e3(t/u/v)_a
309      !!
310      !! Reference  : Leclair, M., and Madec, G. 2011, Ocean Modelling.
311      !!----------------------------------------------------------------------
312      INTEGER, INTENT( in )           ::   kt             ! time step
313      INTEGER, INTENT( in )           ::   Kbb, Kmm, Kaa  ! time step
314      INTEGER, INTENT( in ), OPTIONAL ::   kcall          ! optional argument indicating call sequence
315      !
316      INTEGER                ::   ji, jj, jk            ! dummy loop indices
317      INTEGER , DIMENSION(3) ::   ijk_max, ijk_min      ! temporary integers
318      REAL(wp)               ::   z_tmin, z_tmax        ! local scalars
319      LOGICAL                ::   ll_do_bclinic         ! local logical
320      REAL(dp), DIMENSION(jpi,jpj)     :: z_scale
321      REAL(wp), DIMENSION(jpi,jpj)     ::   zht, zwu, zwv, zhdiv
322      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::   ze3t
323      LOGICAL , DIMENSION(:,:,:), ALLOCATABLE ::   llmsk
324      !!----------------------------------------------------------------------
325      !
326      IF( ln_linssh )   RETURN      ! No calculation in linear free surface
327      !
328      IF( ln_timing )   CALL timing_start('dom_vvl_sf_nxt')
329      !
330      IF( kt == nit000 ) THEN
331         IF(lwp) WRITE(numout,*)
332         IF(lwp) WRITE(numout,*) 'dom_vvl_sf_nxt : compute after scale factors'
333         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~~'
334      ENDIF
335
336      ll_do_bclinic = .TRUE.
337      IF( PRESENT(kcall) ) THEN
338         IF( kcall == 2 .AND. ln_vvl_ztilde )   ll_do_bclinic = .FALSE.
339      ENDIF
340
341      ! ******************************* !
342      ! After acale factors at t-points !
343      ! ******************************* !
344      !                                                ! --------------------------------------------- !
345      !                                                ! z_star coordinate and barotropic z-tilde part !
346      !                                                ! --------------------------------------------- !
347      !
348      z_scale(:,:) = ( ssh(:,:,Kaa) - ssh(:,:,Kbb) ) * ssmask(:,:) / ( ht_0(:,:) + ssh(:,:,Kmm) + 1. - ssmask(:,:) )
349      DO jk = 1, jpkm1
350         ! formally this is the same as e3t(:,:,:,Kaa) = e3t_0*(1+ssha/ht_0)
351         e3t(:,:,jk,Kaa) = e3t(:,:,jk,Kbb) + e3t(:,:,jk,Kmm) * z_scale(:,:) * tmask(:,:,jk)
352      END DO
353      !
354      IF( (ln_vvl_ztilde .OR. ln_vvl_layer) .AND. ll_do_bclinic ) THEN   ! z_tilde or layer coordinate !
355         !                                                               ! ------baroclinic part------ !
356         ! I - initialization
357         ! ==================
358
359         ! 1 - barotropic divergence
360         ! -------------------------
361         zhdiv(:,:) = 0._wp
362         zht(:,:)   = 0._wp
363         DO jk = 1, jpkm1
364            zhdiv(:,:) = zhdiv(:,:) + e3t(:,:,jk,Kmm) * hdiv(:,:,jk)
365            zht  (:,:) = zht  (:,:) + e3t(:,:,jk,Kmm) * tmask(:,:,jk)
366         END DO
367         zhdiv(:,:) = zhdiv(:,:) / ( zht(:,:) + 1. - tmask_i(:,:) )
368
369         ! 2 - Low frequency baroclinic horizontal divergence  (z-tilde case only)
370         ! --------------------------------------------------
371         IF( ln_vvl_ztilde ) THEN
372            IF( kt > nit000 ) THEN
373               DO jk = 1, jpkm1
374                  hdiv_lf(:,:,jk) = hdiv_lf(:,:,jk) - rn_Dt * frq_rst_hdv(:,:)   &
375                     &          * ( hdiv_lf(:,:,jk) - e3t(:,:,jk,Kmm) * ( hdiv(:,:,jk) - zhdiv(:,:) ) )
376               END DO
377            ENDIF
378         ENDIF
379
380         ! II - after z_tilde increments of vertical scale factors
381         ! =======================================================
382         tilde_e3t_a(:,:,:) = 0._wp  ! tilde_e3t_a used to store tendency terms
383
384         ! 1 - High frequency divergence term
385         ! ----------------------------------
386         IF( ln_vvl_ztilde ) THEN     ! z_tilde case
387            DO jk = 1, jpkm1
388               tilde_e3t_a(:,:,jk) = tilde_e3t_a(:,:,jk) - ( e3t(:,:,jk,Kmm) * ( hdiv(:,:,jk) - zhdiv(:,:) ) - hdiv_lf(:,:,jk) )
389            END DO
390         ELSE                         ! layer case
391            DO jk = 1, jpkm1
392               tilde_e3t_a(:,:,jk) = tilde_e3t_a(:,:,jk) -   e3t(:,:,jk,Kmm) * ( hdiv(:,:,jk) - zhdiv(:,:) ) * tmask(:,:,jk)
393            END DO
394         ENDIF
395
396         ! 2 - Restoring term (z-tilde case only)
397         ! ------------------
398         IF( ln_vvl_ztilde ) THEN
399            DO jk = 1, jpk
400               tilde_e3t_a(:,:,jk) = tilde_e3t_a(:,:,jk) - frq_rst_e3t(:,:) * tilde_e3t_b(:,:,jk)
401            END DO
402         ENDIF
403
404         ! 3 - Thickness diffusion term
405         ! ----------------------------
406         zwu(:,:) = 0._wp
407         zwv(:,:) = 0._wp
408         DO_3D( 1, 0, 1, 0, 1, jpkm1 )   ! a - first derivative: diffusive fluxes
409            un_td(ji,jj,jk) = rn_ahe3 * umask(ji,jj,jk) * e2_e1u(ji,jj)           &
410               &            * ( tilde_e3t_b(ji,jj,jk) - tilde_e3t_b(ji+1,jj  ,jk) )
411            vn_td(ji,jj,jk) = rn_ahe3 * vmask(ji,jj,jk) * e1_e2v(ji,jj)           &
412               &            * ( tilde_e3t_b(ji,jj,jk) - tilde_e3t_b(ji  ,jj+1,jk) )
413            zwu(ji,jj) = zwu(ji,jj) + un_td(ji,jj,jk)
414            zwv(ji,jj) = zwv(ji,jj) + vn_td(ji,jj,jk)
415         END_3D
416         DO_2D( 1, 1, 1, 1 )             ! b - correction for last oceanic u-v points
417            un_td(ji,jj,mbku(ji,jj)) = un_td(ji,jj,mbku(ji,jj)) - zwu(ji,jj)
418            vn_td(ji,jj,mbkv(ji,jj)) = vn_td(ji,jj,mbkv(ji,jj)) - zwv(ji,jj)
419         END_2D
420         DO_3D( 0, 0, 0, 0, 1, jpkm1 )   ! c - second derivative: divergence of diffusive fluxes
421            tilde_e3t_a(ji,jj,jk) = tilde_e3t_a(ji,jj,jk) + (   un_td(ji-1,jj  ,jk) - un_td(ji,jj,jk)    &
422               &                                          +     vn_td(ji  ,jj-1,jk) - vn_td(ji,jj,jk)    &
423               &                                            ) * r1_e1e2t(ji,jj)
424         END_3D
425         !                               ! d - thickness diffusion transport: boundary conditions
426         !                             (stored for tracer advction and continuity equation)
427         CALL lbc_lnk_multi( 'domvvl', un_td , 'U' , -1._wp, vn_td , 'V' , -1._wp)
428         ! 4 - Time stepping of baroclinic scale factors
429         ! ---------------------------------------------
430         CALL lbc_lnk( 'domvvl', tilde_e3t_a(:,:,:), 'T', 1._wp )
431         tilde_e3t_a(:,:,:) = tilde_e3t_b(:,:,:) + rDt * tmask(:,:,:) * tilde_e3t_a(:,:,:)
432
433         ! Maximum deformation control
434         ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~
435         ALLOCATE( ze3t(jpi,jpj,jpk), llmsk(jpi,jpj,jpk) )
436         DO_3D( 0, 0, 0, 0, 1, jpkm1 )
437            ze3t(ji,jj,jk) = tilde_e3t_a(ji,jj,jk) / e3t_0(ji,jj,jk) * tmask(ji,jj,jk) * tmask_i(ji,jj)
438         END_3D
439         !
440         llmsk(   1:Nis1,:,:) = .FALSE.   ! exclude halos from the checked region
441         llmsk(Nie1: jpi,:,:) = .FALSE.
442         llmsk(:,   1:Njs1,:) = .FALSE.
443         llmsk(:,Nje1: jpj,:) = .FALSE.
444         !
445         llmsk(Nis0:Nie0,Njs0:Nje0,:) = tmask(Nis0:Nie0,Njs0:Nje0,:) == 1._wp                  ! define only the inner domain
446         z_tmax = MAXVAL( ze3t(:,:,:), mask = llmsk )   ;   CALL mpp_max( 'domvvl', z_tmax )   ! max over the global domain
447         z_tmin = MINVAL( ze3t(:,:,:), mask = llmsk )   ;   CALL mpp_min( 'domvvl', z_tmin )   ! min over the global domain
448         ! - ML - test: for the moment, stop simulation for too large e3_t variations
449         IF( ( z_tmax >  rn_zdef_max ) .OR. ( z_tmin < - rn_zdef_max ) ) THEN
450            CALL mpp_maxloc( 'domvvl', ze3t, llmsk, z_tmax, ijk_max )
451            CALL mpp_minloc( 'domvvl', ze3t, llmsk, z_tmin, ijk_min )
452            IF (lwp) THEN
453               WRITE(numout, *) 'MAX( tilde_e3t_a(:,:,:) / e3t_0(:,:,:) ) =', z_tmax
454               WRITE(numout, *) 'at i, j, k=', ijk_max
455               WRITE(numout, *) 'MIN( tilde_e3t_a(:,:,:) / e3t_0(:,:,:) ) =', z_tmin
456               WRITE(numout, *) 'at i, j, k=', ijk_min
457               CALL ctl_stop( 'STOP', 'MAX( ABS( tilde_e3t_a(:,:,: ) ) / e3t_0(:,:,:) ) too high')
458            ENDIF
459         ENDIF
460         DEALLOCATE( ze3t, llmsk )
461         ! - ML - end test
462         ! - ML - Imposing these limits will cause a baroclinicity error which is corrected for below
463         tilde_e3t_a(:,:,:) = MIN( tilde_e3t_a(:,:,:),   rn_zdef_max * e3t_0(:,:,:) )
464         tilde_e3t_a(:,:,:) = MAX( tilde_e3t_a(:,:,:), - rn_zdef_max * e3t_0(:,:,:) )
465
466         !
467         ! "tilda" change in the after scale factor
468         ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
469         DO jk = 1, jpkm1
470            dtilde_e3t_a(:,:,jk) = tilde_e3t_a(:,:,jk) - tilde_e3t_b(:,:,jk)
471         END DO
472         ! III - Barotropic repartition of the sea surface height over the baroclinic profile
473         ! ==================================================================================
474         ! add ( ssh increment + "baroclinicity error" ) proportionly to e3t(n)
475         ! - ML - baroclinicity error should be better treated in the future
476         !        i.e. locally and not spread over the water column.
477         !        (keep in mind that the idea is to reduce Eulerian velocity as much as possible)
478         zht(:,:) = 0.
479         DO jk = 1, jpkm1
480            zht(:,:)  = zht(:,:) + tilde_e3t_a(:,:,jk) * tmask(:,:,jk)
481         END DO
482         z_scale(:,:) =  - zht(:,:) / ( ht_0(:,:) + ssh(:,:,Kmm) + 1. - ssmask(:,:) )
483         DO jk = 1, jpkm1
484            dtilde_e3t_a(:,:,jk) = dtilde_e3t_a(:,:,jk) + e3t(:,:,jk,Kmm) * z_scale(:,:) * tmask(:,:,jk)
485         END DO
486
487      ENDIF
488
489      IF( ln_vvl_ztilde .OR. ln_vvl_layer )  THEN   ! z_tilde or layer coordinate !
490      !                                           ! ---baroclinic part--------- !
491         DO jk = 1, jpkm1
492            e3t(:,:,jk,Kaa) = e3t(:,:,jk,Kaa) + dtilde_e3t_a(:,:,jk) * tmask(:,:,jk)
493         END DO
494      ENDIF
495
496      IF( ln_vvl_dbg .AND. .NOT. ll_do_bclinic ) THEN   ! - ML - test: control prints for debuging
497         !
498         IF( lwp ) WRITE(numout, *) 'kt =', kt
499         IF ( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN
500            z_tmax = MAXVAL( tmask(:,:,1) * tmask_i(:,:) * ABS( zht(:,:) ) )
501            CALL mpp_max( 'domvvl', z_tmax )                             ! max over the global domain
502            IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(SUM(tilde_e3t_a))) =', z_tmax
503         END IF
504         !
505         zht(:,:) = 0.0_wp
506         DO jk = 1, jpkm1
507            zht(:,:) = zht(:,:) + e3t(:,:,jk,Kmm) * tmask(:,:,jk)
508         END DO
509         z_tmax = MAXVAL( tmask(:,:,1) * tmask_i(:,:) * ABS( ht_0(:,:) + ssh(:,:,Kmm) - zht(:,:) ) )
510         CALL mpp_max( 'domvvl', z_tmax )                                ! max over the global domain
511         IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(ht_0+sshn-SUM(e3t(:,:,:,Kmm)))) =', z_tmax
512         !
513         zht(:,:) = 0.0_wp
514         DO jk = 1, jpkm1
515            zht(:,:) = zht(:,:) + e3t(:,:,jk,Kaa) * tmask(:,:,jk)
516         END DO
517         z_tmax = MAXVAL( tmask(:,:,1) * tmask_i(:,:) * ABS( ht_0(:,:) + ssh(:,:,Kaa) - zht(:,:) ) )
518         CALL mpp_max( 'domvvl', z_tmax )                                ! max over the global domain
519         IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(ht_0+ssha-SUM(e3t(:,:,:,Kaa)))) =', z_tmax
520         !
521         zht(:,:) = 0.0_wp
522         DO jk = 1, jpkm1
523            zht(:,:) = zht(:,:) + e3t(:,:,jk,Kbb) * tmask(:,:,jk)
524         END DO
525         z_tmax = MAXVAL( tmask(:,:,1) * tmask_i(:,:) * ABS( ht_0(:,:) + ssh(:,:,Kbb) - zht(:,:) ) )
526         CALL mpp_max( 'domvvl', z_tmax )                                ! max over the global domain
527         IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(ht_0+sshb-SUM(e3t(:,:,:,Kbb)))) =', z_tmax
528         !
529         z_tmax = MAXVAL( tmask(:,:,1) *  ABS( ssh(:,:,Kbb) ) )
530         CALL mpp_max( 'domvvl', z_tmax )                                ! max over the global domain
531         IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(ssh(:,:,Kbb)))) =', z_tmax
532         !
533         z_tmax = MAXVAL( tmask(:,:,1) *  ABS( ssh(:,:,Kmm) ) )
534         CALL mpp_max( 'domvvl', z_tmax )                                ! max over the global domain
535         IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(ssh(:,:,Kmm)))) =', z_tmax
536         !
537         z_tmax = MAXVAL( tmask(:,:,1) *  ABS( ssh(:,:,Kaa) ) )
538         CALL mpp_max( 'domvvl', z_tmax )                                ! max over the global domain
539         IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(ssh(:,:,Kaa)))) =', z_tmax
540      END IF
541
542      ! *********************************** !
543      ! After scale factors at u- v- points !
544      ! *********************************** !
545
546      CALL dom_vvl_interpol( CASTWP(e3t(:,:,:,Kaa)), e3u(:,:,:,Kaa), 'U' )
547      CALL dom_vvl_interpol( CASTWP(e3t(:,:,:,Kaa)), e3v(:,:,:,Kaa), 'V' )
548
549      ! *********************************** !
550      ! After depths at u- v points         !
551      ! *********************************** !
552
553      hu(:,:,Kaa) = e3u(:,:,1,Kaa) * umask(:,:,1)
554      hv(:,:,Kaa) = e3v(:,:,1,Kaa) * vmask(:,:,1)
555      DO jk = 2, jpkm1
556         hu(:,:,Kaa) = hu(:,:,Kaa) + e3u(:,:,jk,Kaa) * umask(:,:,jk)
557         hv(:,:,Kaa) = hv(:,:,Kaa) + e3v(:,:,jk,Kaa) * vmask(:,:,jk)
558      END DO
559      !                                        ! Inverse of the local depth
560!!gm BUG ?  don't understand the use of umask_i here .....
561      r1_hu(:,:,Kaa) = ssumask(:,:) / ( hu(:,:,Kaa) + 1._wp - ssumask(:,:) )
562      r1_hv(:,:,Kaa) = ssvmask(:,:) / ( hv(:,:,Kaa) + 1._wp - ssvmask(:,:) )
563      !
564      IF( ln_timing )   CALL timing_stop('dom_vvl_sf_nxt')
565      !
566   END SUBROUTINE dom_vvl_sf_nxt
567
568
569   SUBROUTINE dom_vvl_sf_update( kt, Kbb, Kmm, Kaa )
570      !!----------------------------------------------------------------------
571      !!                ***  ROUTINE dom_vvl_sf_update  ***
572      !!
573      !! ** Purpose :  for z tilde case: compute time filter and swap of scale factors
574      !!               compute all depths and related variables for next time step
575      !!               write outputs and restart file
576      !!
577      !! ** Method  :  - swap of e3t with trick for volume/tracer conservation (ONLY FOR Z TILDE CASE)
578      !!               - reconstruct scale factor at other grid points (interpolate)
579      !!               - recompute depths and water height fields
580      !!
581      !! ** Action  :  - tilde_e3t_(b/n) ready for next time step
582      !!               - Recompute:
583      !!                    e3(u/v)_b
584      !!                    e3w(:,:,:,Kmm)
585      !!                    e3(u/v)w_b
586      !!                    e3(u/v)w_n
587      !!                    gdept(:,:,:,Kmm), gdepw(:,:,:,Kmm)  and gde3w
588      !!                    h(u/v) and h(u/v)r
589      !!
590      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
591      !!              Leclair, M., and G. Madec, 2011, Ocean Modelling.
592      !!----------------------------------------------------------------------
593      INTEGER, INTENT( in ) ::   kt              ! time step
594      INTEGER, INTENT( in ) ::   Kbb, Kmm, Kaa   ! time level indices
595      !
596      INTEGER  ::   ji, jj, jk   ! dummy loop indices
597      REAL(wp) ::   zcoef        ! local scalar
598      !!----------------------------------------------------------------------
599      !
600      IF( ln_linssh )   RETURN      ! No calculation in linear free surface
601      !
602      IF( ln_timing )   CALL timing_start('dom_vvl_sf_update')
603      !
604      IF( kt == nit000 )   THEN
605         IF(lwp) WRITE(numout,*)
606         IF(lwp) WRITE(numout,*) 'dom_vvl_sf_update : - interpolate scale factors and compute depths for next time step'
607         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~~~~~'
608      ENDIF
609      !
610      ! Time filter and swap of scale factors
611      ! =====================================
612      ! - ML - e3(t/u/v)_b are allready computed in dynnxt.
613      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN
614         IF( l_1st_euler ) THEN
615            tilde_e3t_b(:,:,:) = tilde_e3t_n(:,:,:)
616         ELSE
617            tilde_e3t_b(:,:,:) = tilde_e3t_n(:,:,:) &
618            &         + rn_atfp * ( tilde_e3t_b(:,:,:) - 2.0_wp * tilde_e3t_n(:,:,:) + tilde_e3t_a(:,:,:) )
619         ENDIF
620         tilde_e3t_n(:,:,:) = tilde_e3t_a(:,:,:)
621      ENDIF
622
623      ! Compute all missing vertical scale factor and depths
624      ! ====================================================
625      ! Horizontal scale factor interpolations
626      ! --------------------------------------
627      ! - ML - e3u(:,:,:,Kbb) and e3v(:,:,:,Kbb) are already computed in dynnxt
628      ! - JC - hu(:,:,:,Kbb), hv(:,:,:,:,Kbb), hur_b, hvr_b also
629
630      CALL dom_vvl_interpol( e3u(:,:,:,Kmm), e3f(:,:,:), 'F'  )
631
632      ! Vertical scale factor interpolations
633      CALL dom_vvl_interpol( CASTWP(e3t(:,:,:,Kmm)),  e3w(:,:,:,Kmm), 'W'  )
634      CALL dom_vvl_interpol( e3u(:,:,:,Kmm), e3uw(:,:,:,Kmm), 'UW' )
635      CALL dom_vvl_interpol( e3v(:,:,:,Kmm), e3vw(:,:,:,Kmm), 'VW' )
636      CALL dom_vvl_interpol( CASTWP(e3t(:,:,:,Kbb)),  e3w(:,:,:,Kbb), 'W'  )
637      CALL dom_vvl_interpol( e3u(:,:,:,Kbb), e3uw(:,:,:,Kbb), 'UW' )
638      CALL dom_vvl_interpol( e3v(:,:,:,Kbb), e3vw(:,:,:,Kbb), 'VW' )
639
640      ! t- and w- points depth (set the isf depth as it is in the initial step)
641      gdept(:,:,1,Kmm) = 0.5_wp * e3w(:,:,1,Kmm)
642      gdepw(:,:,1,Kmm) = 0.0_wp
643      gde3w(:,:,1) = gdept(:,:,1,Kmm) - ssh(:,:,Kmm)
644      DO_3D( 1, 1, 1, 1, 2, jpk )
645        !    zcoef = (tmask(ji,jj,jk) - wmask(ji,jj,jk))   ! 0 everywhere tmask = wmask, ie everywhere expect at jk = mikt
646                                                           ! 1 for jk = mikt
647         zcoef = (tmask(ji,jj,jk) - wmask(ji,jj,jk))
648         gdepw(ji,jj,jk,Kmm) = gdepw(ji,jj,jk-1,Kmm) + e3t(ji,jj,jk-1,Kmm)
649         gdept(ji,jj,jk,Kmm) =    zcoef  * ( gdepw(ji,jj,jk  ,Kmm) + 0.5 * e3w(ji,jj,jk,Kmm) )  &
650             &             + (1-zcoef) * ( gdept(ji,jj,jk-1,Kmm) +       e3w(ji,jj,jk,Kmm) )
651         gde3w(ji,jj,jk) = gdept(ji,jj,jk,Kmm) - ssh(ji,jj,Kmm)
652      END_3D
653
654      ! Local depth and Inverse of the local depth of the water
655      ! -------------------------------------------------------
656      !
657      ht(:,:) = e3t(:,:,1,Kmm) * tmask(:,:,1)
658      DO jk = 2, jpkm1
659         ht(:,:) = ht(:,:) + e3t(:,:,jk,Kmm) * tmask(:,:,jk)
660      END DO
661
662      ! write restart file
663      ! ==================
664      IF( lrst_oce  )   CALL dom_vvl_rst( kt, Kbb, Kmm, 'WRITE' )
665      !
666      IF( ln_timing )   CALL timing_stop('dom_vvl_sf_update')
667      !
668   END SUBROUTINE dom_vvl_sf_update
669
670
671   SUBROUTINE dom_vvl_interpol( pe3_in, pe3_out, pout )
672      !!---------------------------------------------------------------------
673      !!                  ***  ROUTINE dom_vvl__interpol  ***
674      !!
675      !! ** Purpose :   interpolate scale factors from one grid point to another
676      !!
677      !! ** Method  :   e3_out = e3_0 + interpolation(e3_in - e3_0)
678      !!                - horizontal interpolation: grid cell surface averaging
679      !!                - vertical interpolation: simple averaging
680      !!----------------------------------------------------------------------
681      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in   ) ::  pe3_in    ! input e3 to be interpolated
682      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::  pe3_out   ! output interpolated e3
683      CHARACTER(LEN=*)                , INTENT(in   ) ::  pout      ! grid point of out scale factors
684      !                                                             !   =  'U', 'V', 'W, 'F', 'UW' or 'VW'
685      !
686      INTEGER ::   ji, jj, jk                                       ! dummy loop indices
687      REAL(wp) ::  zlnwd                                            ! =1./0. when ln_wd_il = T/F
688      !!----------------------------------------------------------------------
689      !
690      IF(ln_wd_il) THEN
691        zlnwd = 1.0_wp
692      ELSE
693        zlnwd = 0.0_wp
694      END IF
695      !
696      SELECT CASE ( pout )    !==  type of interpolation  ==!
697         !
698      CASE( 'U' )                   !* from T- to U-point : hor. surface weighted mean
699         DO_3D( 1, 0, 1, 0, 1, jpk )
700            pe3_out(ji,jj,jk) = 0.5_wp * (  umask(ji,jj,jk) * (1.0_wp - zlnwd) + zlnwd ) * r1_e1e2u(ji,jj)   &
701               &                       * (   e1e2t(ji  ,jj) * ( pe3_in(ji  ,jj,jk) - e3t_0(ji  ,jj,jk) )     &
702               &                           + e1e2t(ji+1,jj) * ( pe3_in(ji+1,jj,jk) - e3t_0(ji+1,jj,jk) ) )
703         END_3D
704         CALL lbc_lnk( 'domvvl', pe3_out(:,:,:), 'U', 1._wp )
705         pe3_out(:,:,:) = pe3_out(:,:,:) + e3u_0(:,:,:)
706         !
707      CASE( 'V' )                   !* from T- to V-point : hor. surface weighted mean
708         DO_3D( 1, 0, 1, 0, 1, jpk )
709            pe3_out(ji,jj,jk) = 0.5_wp * ( vmask(ji,jj,jk)  * (1.0_wp - zlnwd) + zlnwd ) * r1_e1e2v(ji,jj)   &
710               &                       * (   e1e2t(ji,jj  ) * ( pe3_in(ji,jj  ,jk) - e3t_0(ji,jj  ,jk) )     &
711               &                           + e1e2t(ji,jj+1) * ( pe3_in(ji,jj+1,jk) - e3t_0(ji,jj+1,jk) ) )
712         END_3D
713         CALL lbc_lnk( 'domvvl', pe3_out(:,:,:), 'V', 1._wp )
714         pe3_out(:,:,:) = pe3_out(:,:,:) + e3v_0(:,:,:)
715         !
716      CASE( 'F' )                   !* from U-point to F-point : hor. surface weighted mean
717         DO_3D( 1, 0, 1, 0, 1, jpk )
718            pe3_out(ji,jj,jk) = 0.5_wp * (  umask(ji,jj,jk) * umask(ji,jj+1,jk) * (1.0_wp - zlnwd) + zlnwd ) &
719               &                       *    r1_e1e2f(ji,jj)                                                  &
720               &                       * (   e1e2u(ji,jj  ) * ( pe3_in(ji,jj  ,jk) - e3u_0(ji,jj  ,jk) )     &
721               &                           + e1e2u(ji,jj+1) * ( pe3_in(ji,jj+1,jk) - e3u_0(ji,jj+1,jk) ) )
722         END_3D
723         CALL lbc_lnk( 'domvvl', pe3_out(:,:,:), 'F', 1._wp )
724         pe3_out(:,:,:) = pe3_out(:,:,:) + e3f_0(:,:,:)
725         !
726      CASE( 'W' )                   !* from T- to W-point : vertical simple mean
727         !
728         pe3_out(:,:,1) = e3w_0(:,:,1) + pe3_in(:,:,1) - e3t_0(:,:,1)
729         ! - ML - The use of mask in this formulea enables the special treatment of the last w-point without indirect adressing
730!!gm BUG? use here wmask in case of ISF ?  to be checked
731         DO jk = 2, jpk
732            pe3_out(:,:,jk) = e3w_0(:,:,jk) + ( 1.0_wp - 0.5_wp * ( tmask(:,:,jk) * (1.0_wp - zlnwd) + zlnwd ) )   &
733               &                            * ( pe3_in(:,:,jk-1) - e3t_0(:,:,jk-1) )                               &
734               &                            +            0.5_wp * ( tmask(:,:,jk) * (1.0_wp - zlnwd) + zlnwd )     &
735               &                            * ( pe3_in(:,:,jk  ) - e3t_0(:,:,jk  ) )
736         END DO
737         !
738      CASE( 'UW' )                  !* from U- to UW-point : vertical simple mean
739         !
740         pe3_out(:,:,1) = e3uw_0(:,:,1) + pe3_in(:,:,1) - e3u_0(:,:,1)
741         ! - ML - The use of mask in this formaula enables the special treatment of the last w- point without indirect adressing
742!!gm BUG? use here wumask in case of ISF ?  to be checked
743         DO jk = 2, jpk
744            pe3_out(:,:,jk) = e3uw_0(:,:,jk) + ( 1.0_wp - 0.5_wp * ( umask(:,:,jk) * (1.0_wp - zlnwd) + zlnwd ) )  &
745               &                             * ( pe3_in(:,:,jk-1) - e3u_0(:,:,jk-1) )                              &
746               &                             +            0.5_wp * ( umask(:,:,jk) * (1.0_wp - zlnwd) + zlnwd )    &
747               &                             * ( pe3_in(:,:,jk  ) - e3u_0(:,:,jk  ) )
748         END DO
749         !
750      CASE( 'VW' )                  !* from V- to VW-point : vertical simple mean
751         !
752         pe3_out(:,:,1) = e3vw_0(:,:,1) + pe3_in(:,:,1) - e3v_0(:,:,1)
753         ! - ML - The use of mask in this formaula enables the special treatment of the last w- point without indirect adressing
754!!gm BUG? use here wvmask in case of ISF ?  to be checked
755         DO jk = 2, jpk
756            pe3_out(:,:,jk) = e3vw_0(:,:,jk) + ( 1.0_wp - 0.5_wp * ( vmask(:,:,jk) * (1.0_wp - zlnwd) + zlnwd ) )  &
757               &                             * ( pe3_in(:,:,jk-1) - e3v_0(:,:,jk-1) )                              &
758               &                             +            0.5_wp * ( vmask(:,:,jk) * (1.0_wp - zlnwd) + zlnwd )    &
759               &                             * ( pe3_in(:,:,jk  ) - e3v_0(:,:,jk  ) )
760         END DO
761      END SELECT
762      !
763   END SUBROUTINE dom_vvl_interpol
764
765
766   SUBROUTINE dom_vvl_rst( kt, Kbb, Kmm, cdrw )
767      !!---------------------------------------------------------------------
768      !!                   ***  ROUTINE dom_vvl_rst  ***
769      !!
770      !! ** Purpose :   Read or write VVL file in restart file
771      !!
772      !! ** Method  : * restart comes from a linear ssh simulation :
773      !!                   an attempt to read e3t_n stops simulation
774      !!              * restart comes from a z-star, z-tilde, or layer :
775      !!                   read e3t_n and e3t_b
776      !!              * restart comes from a z-star :
777      !!                   set tilde_e3t_n, tilde_e3t_n, and hdiv_lf to 0
778      !!              * restart comes from layer :
779      !!                   read tilde_e3t_n and tilde_e3t_b
780      !!                   set hdiv_lf to 0
781      !!              * restart comes from a z-tilde:
782      !!                   read tilde_e3t_n, tilde_e3t_b, and hdiv_lf
783      !!
784      !!              NB: if l_1st_euler = T (ln_1st_euler or ssh_b not found)
785      !!                   Kbb fields set to Kmm ones
786      !!----------------------------------------------------------------------
787      INTEGER         , INTENT(in) ::   kt        ! ocean time-step
788      INTEGER         , INTENT(in) ::   Kbb, Kmm  ! ocean time level indices
789      CHARACTER(len=*), INTENT(in) ::   cdrw      ! "READ"/"WRITE" flag
790      !
791      INTEGER ::   ji, jj, jk      ! dummy loop indices
792      INTEGER ::   id3, id4, id5   ! local integers
793      !!----------------------------------------------------------------------
794      !
795      !                                      !=====================!
796      IF( TRIM(cdrw) == 'READ' ) THEN        !  Read / initialise  !
797         !                                   !=====================!
798         !
799         IF( ln_rstart ) THEN                   !==  Read the restart file  ==!
800            !
801            CALL rst_read_open                                          !*  open the restart file if necessary
802            !                                         ! --------- !
803            !                                         ! all cases !
804            !                                         ! --------- !
805            !
806            id3 = iom_varid( numror, 'tilde_e3t_b', ldstop = .FALSE. )  !*  check presence
807            id4 = iom_varid( numror, 'tilde_e3t_n', ldstop = .FALSE. )
808            id5 = iom_varid( numror, 'hdiv_lf'    , ldstop = .FALSE. )
809            !
810            !                                                           !*  scale factors
811            IF(lwp) WRITE(numout,*)    '          Kmm scale factor read in the restart file'
812            CALL iom_get( numror, jpdom_auto, 'e3t_n', e3t(:,:,:,Kmm) )
813            WHERE ( tmask(:,:,:) == 0.0_wp )
814               e3t(:,:,:,Kmm) = e3t_0(:,:,:)
815            END WHERE
816            IF( l_1st_euler ) THEN                       ! euler
817               IF(lwp) WRITE(numout,*) '          Euler first time step : e3t(Kbb) = e3t(Kmm)'
818               e3t(:,:,:,Kbb) = e3t(:,:,:,Kmm)
819            ELSE                                         ! leap frog
820               IF(lwp) WRITE(numout,*) '          Kbb scale factor read in the restart file'
821               CALL iom_get( numror, jpdom_auto, 'e3t_b', e3t(:,:,:,Kbb) )
822               WHERE ( tmask(:,:,:) == 0.0_wp )
823                  e3t(:,:,:,Kbb) = e3t_0(:,:,:)
824               END WHERE
825            ENDIF
826            !                                         ! ------------ !
827            IF( ln_vvl_zstar ) THEN                   !  z_star case !
828               !                                      ! ------------ !
829               IF( MIN( id3, id4 ) > 0 ) THEN
830                  CALL ctl_stop( 'dom_vvl_rst: z_star cannot restart from a z_tilde or layer run' )
831               ENDIF
832               !                                      ! ------------------------ !
833            ELSE                                      !  z_tilde and layer cases !
834               !                                      ! ------------------------ !
835               !
836               IF( id4 > 0 ) THEN                                       !*  scale factor increments
837                  IF(lwp) WRITE(numout,*)    '          Kmm scale factor increments read in the restart file'
838                  CALL iom_get( numror, jpdom_auto, 'tilde_e3t_n', tilde_e3t_n(:,:,:) )
839                  IF( l_1st_euler ) THEN                 ! euler
840                     IF(lwp) WRITE(numout,*) '          Euler first time step : tilde_e3t(Kbb) = tilde_e3t(Kmm)'
841                     tilde_e3t_b(:,:,:) = tilde_e3t_n(:,:,:)
842                  ELSE                                   ! leap frog
843                     IF(lwp) WRITE(numout,*) '          Kbb scale factor increments read in the restart file'
844                     CALL iom_get( numror, jpdom_auto, 'tilde_e3t_b', tilde_e3t_b(:,:,:) )
845                  ENDIF
846               ELSE
847                  tilde_e3t_b(:,:,:) = 0.0_wp
848                  tilde_e3t_n(:,:,:) = 0.0_wp
849               ENDIF
850               !                                      ! ------------ !
851               IF( ln_vvl_ztilde ) THEN               ! z_tilde case !
852                  !                                   ! ------------ !
853                  IF( id5 > 0 ) THEN  ! required array exists
854                     CALL iom_get( numror, jpdom_auto, 'hdiv_lf', hdiv_lf(:,:,:) )
855                  ELSE                ! array is missing
856                     hdiv_lf(:,:,:) = 0.0_wp
857                  ENDIF
858               ENDIF
859            ENDIF
860            !
861         ELSE                                   !==  Initialize at "rest" with ssh  ==!
862            !
863            DO jk = 1, jpk
864               e3t(:,:,jk,Kmm) =  e3t_0(:,:,jk) * ( 1._wp + ssh(:,:,Kmm) * r1_ht_0(:,:) * tmask(:,:,jk) )
865            END DO
866            e3t(:,:,:,Kbb) = e3t(:,:,:,Kmm)
867            !
868            IF( ln_vvl_ztilde .OR. ln_vvl_layer) THEN
869               tilde_e3t_b(:,:,:) = 0._wp
870               tilde_e3t_n(:,:,:) = 0._wp
871               IF( ln_vvl_ztilde ) hdiv_lf(:,:,:) = 0._wp
872            ENDIF
873         ENDIF
874         !                                       !=======================!
875      ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN       !  Create restart file  !
876         !                                       !=======================!
877         !
878         IF(lwp) WRITE(numout,*) '---- dom_vvl_rst ----'
879         !                                           ! --------- !
880         !                                           ! all cases !
881         !                                           ! --------- !
882         CALL iom_rstput( kt, nitrst, numrow, 'e3t_b', e3t(:,:,:,Kbb) )
883         CALL iom_rstput( kt, nitrst, numrow, 'e3t_n', e3t(:,:,:,Kmm) )
884         !                                           ! ----------------------- !
885         IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN  ! z_tilde and layer cases !
886            !                                        ! ----------------------- !
887            CALL iom_rstput( kt, nitrst, numrow, 'tilde_e3t_b', tilde_e3t_b(:,:,:))
888            CALL iom_rstput( kt, nitrst, numrow, 'tilde_e3t_n', tilde_e3t_n(:,:,:))
889         END IF
890         !                                           ! -------------!
891         IF( ln_vvl_ztilde ) THEN                    ! z_tilde case !
892            !                                        ! ------------ !
893            CALL iom_rstput( kt, nitrst, numrow, 'hdiv_lf', hdiv_lf(:,:,:))
894         ENDIF
895         !
896      ENDIF
897      !
898   END SUBROUTINE dom_vvl_rst
899
900
901   SUBROUTINE dom_vvl_ctl
902      !!---------------------------------------------------------------------
903      !!                  ***  ROUTINE dom_vvl_ctl  ***
904      !!
905      !! ** Purpose :   Control the consistency between namelist options
906      !!                for vertical coordinate
907      !!----------------------------------------------------------------------
908      INTEGER ::   ioptio, ios
909      !!
910      NAMELIST/nam_vvl/ ln_vvl_zstar, ln_vvl_ztilde, ln_vvl_layer, ln_vvl_ztilde_as_zstar, &
911         &              ln_vvl_zstar_at_eqtor      , rn_ahe3     , rn_rst_e3t            , &
912         &              rn_lf_cutoff               , rn_zdef_max , ln_vvl_dbg                ! not yet implemented: ln_vvl_kepe
913      !!----------------------------------------------------------------------
914      !
915      READ  ( numnam_ref, nam_vvl, IOSTAT = ios, ERR = 901)
916901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'nam_vvl in reference namelist' )
917      READ  ( numnam_cfg, nam_vvl, IOSTAT = ios, ERR = 902 )
918902   IF( ios >  0 ) CALL ctl_nam ( ios , 'nam_vvl in configuration namelist' )
919      IF(lwm) WRITE ( numond, nam_vvl )
920      !
921      IF(lwp) THEN                    ! Namelist print
922         WRITE(numout,*)
923         WRITE(numout,*) 'dom_vvl_ctl : choice/control of the variable vertical coordinate'
924         WRITE(numout,*) '~~~~~~~~~~~'
925         WRITE(numout,*) '   Namelist nam_vvl : chose a vertical coordinate'
926         WRITE(numout,*) '      zstar                      ln_vvl_zstar   = ', ln_vvl_zstar
927         WRITE(numout,*) '      ztilde                     ln_vvl_ztilde  = ', ln_vvl_ztilde
928         WRITE(numout,*) '      layer                      ln_vvl_layer   = ', ln_vvl_layer
929         WRITE(numout,*) '      ztilde as zstar   ln_vvl_ztilde_as_zstar  = ', ln_vvl_ztilde_as_zstar
930         WRITE(numout,*) '      ztilde near the equator    ln_vvl_zstar_at_eqtor  = ', ln_vvl_zstar_at_eqtor
931         WRITE(numout,*) '      !'
932         WRITE(numout,*) '      thickness diffusion coefficient                      rn_ahe3      = ', rn_ahe3
933         WRITE(numout,*) '      maximum e3t deformation fractional change            rn_zdef_max  = ', rn_zdef_max
934         IF( ln_vvl_ztilde_as_zstar ) THEN
935            WRITE(numout,*) '      ztilde running in zstar emulation mode (ln_vvl_ztilde_as_zstar=T) '
936            WRITE(numout,*) '         ignoring namelist timescale parameters and using:'
937            WRITE(numout,*) '            hard-wired : z-tilde to zstar restoration timescale (days)'
938            WRITE(numout,*) '                         rn_rst_e3t     = 0.e0'
939            WRITE(numout,*) '            hard-wired : z-tilde cutoff frequency of low-pass filter (days)'
940            WRITE(numout,*) '                         rn_lf_cutoff   = 1.0/rn_Dt'
941         ELSE
942            WRITE(numout,*) '      z-tilde to zstar restoration timescale (days)        rn_rst_e3t   = ', rn_rst_e3t
943            WRITE(numout,*) '      z-tilde cutoff frequency of low-pass filter (days)   rn_lf_cutoff = ', rn_lf_cutoff
944         ENDIF
945         WRITE(numout,*) '         debug prints flag                                 ln_vvl_dbg   = ', ln_vvl_dbg
946      ENDIF
947      !
948      ioptio = 0                      ! Parameter control
949      IF( ln_vvl_ztilde_as_zstar )   ln_vvl_ztilde = .true.
950      IF( ln_vvl_zstar           )   ioptio = ioptio + 1
951      IF( ln_vvl_ztilde          )   ioptio = ioptio + 1
952      IF( ln_vvl_layer           )   ioptio = ioptio + 1
953      !
954      IF( ioptio /= 1 )   CALL ctl_stop( 'Choose ONE vertical coordinate in namelist nam_vvl' )
955      !
956      IF(lwp) THEN                   ! Print the choice
957         WRITE(numout,*)
958         IF( ln_vvl_zstar           ) WRITE(numout,*) '      ==>>>   zstar vertical coordinate is used'
959         IF( ln_vvl_ztilde          ) WRITE(numout,*) '      ==>>>   ztilde vertical coordinate is used'
960         IF( ln_vvl_layer           ) WRITE(numout,*) '      ==>>>   layer vertical coordinate is used'
961         IF( ln_vvl_ztilde_as_zstar ) WRITE(numout,*) '      ==>>>   to emulate a zstar coordinate'
962      ENDIF
963      !
964#if defined key_agrif
965      IF( (.NOT.Agrif_Root()).AND.(.NOT.ln_vvl_zstar) )   CALL ctl_stop( 'AGRIF is implemented with zstar coordinate only' )
966#endif
967      !
968   END SUBROUTINE dom_vvl_ctl
969
970#endif
971
972   !!======================================================================
973END MODULE domvvl
Note: See TracBrowser for help on using the repository browser.