Newer
Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
MODULE dynspg_ts
!! Includes ROMS wd scheme with diagnostic outputs ; puu(:,:,:,Kmm) and puu(:,:,:,Krhs) updates are commented out !
!!======================================================================
!! *** MODULE dynspg_ts ***
!! Ocean dynamics: surface pressure gradient trend, split-explicit scheme
!!======================================================================
!! History : 1.0 ! 2004-12 (L. Bessieres, G. Madec) Original code
!! - ! 2005-11 (V. Garnier, G. Madec) optimization
!! - ! 2006-08 (S. Masson) distributed restart using iom
!! 2.0 ! 2007-07 (D. Storkey) calls to BDY routines
!! - ! 2008-01 (R. Benshila) change averaging method
!! 3.2 ! 2009-07 (R. Benshila, G. Madec) Complete revisit associated to vvl reactivation
!! 3.3 ! 2010-09 (D. Storkey, E. O'Dea) update for BDY for Shelf configurations
!! 3.3 ! 2011-03 (R. Benshila, R. Hordoir, P. Oddo) update calculation of ub_b
!! 3.5 ! 2013-07 (J. Chanut) Switch to Forward-backward time stepping
!! 3.6 ! 2013-11 (A. Coward) Update for z-tilde compatibility
!! 3.7 ! 2015-11 (J. Chanut) free surface simplification
!! - ! 2016-12 (G. Madec, E. Clementi) update for Stoke-Drift divergence
!! 4.0 ! 2017-05 (G. Madec) drag coef. defined at t-point (zdfdrg.F90)
!!---------------------------------------------------------------------
!!----------------------------------------------------------------------
!! dyn_spg_ts : compute surface pressure gradient trend using a time-splitting scheme
!! dyn_spg_ts_init: initialisation of the time-splitting scheme
!! ts_wgt : set time-splitting weights for temporal averaging (or not)
!! ts_rst : read/write time-splitting fields in restart file
!!----------------------------------------------------------------------
USE oce ! ocean dynamics and tracers
USE dom_oce ! ocean space and time domain
USE sbc_oce ! surface boundary condition: ocean
USE isf_oce ! ice shelf variable (fwfisf)
USE zdf_oce ! vertical physics: variables
USE zdfdrg ! vertical physics: top/bottom drag coef.
USE sbcapr ! surface boundary condition: atmospheric pressure
USE dynadv , ONLY: ln_dynadv_vec
USE dynvor ! vortivity scheme indicators
USE phycst ! physical constants
USE dynvor ! vorticity term
USE wet_dry ! wetting/drying flux limter
USE bdy_oce ! open boundary
USE bdyvol ! open boundary volume conservation
USE bdytides ! open boundary condition data
USE bdydyn2d ! open boundary conditions on barotropic variables
USE tide_mod !
USE sbcwave ! surface wave
#if defined key_agrif
USE agrif_oce_interp ! agrif
USE agrif_oce
#endif
#if defined key_asminc
USE asminc ! Assimilation increment
#endif
!
USE in_out_manager ! I/O manager
USE lib_mpp ! distributed memory computing library
USE lbclnk ! ocean lateral boundary conditions (or mpp link)
USE prtctl ! Print control
USE iom ! IOM library
USE restart ! only for lrst_oce
USE iom ! to remove
IMPLICIT NONE
PRIVATE
PUBLIC dyn_spg_ts ! called by dyn_spg
PUBLIC dyn_spg_ts_init ! - - dyn_spg_init
!! Time filtered arrays at baroclinic time step:
REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: un_adv , vn_adv !: Advection vel. at "now" barocl. step
!
INTEGER, SAVE :: icycle ! Number of barotropic sub-steps for each internal step nn_e <= 2.5 nn_e
REAL(wp),SAVE :: rDt_e ! Barotropic time step
!
REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:) :: wgtbtp1, wgtbtp2 ! 1st & 2nd weights used in time filtering of barotropic fields
REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) :: zwz ! ff_f/h at F points
REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) :: ftnw, ftne ! triad of coriolis parameter
REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) :: ftsw, ftse ! (only used with een vorticity scheme)
REAL(wp) :: r1_12 = 1._wp / 12._wp ! local ratios
REAL(wp) :: r1_8 = 0.125_wp !
REAL(wp) :: r1_4 = 0.25_wp !
REAL(wp) :: r1_2 = 0.5_wp !
!! * Substitutions
# include "do_loop_substitute.h90"
# include "domzgr_substitute.h90"
!!----------------------------------------------------------------------
!! NEMO/OCE 4.0 , NEMO Consortium (2018)
!! $Id: dynspg_ts.F90 15489 2021-11-10 09:18:39Z jchanut $
!! Software governed by the CeCILL license (see ./LICENSE)
!!----------------------------------------------------------------------
CONTAINS
INTEGER FUNCTION dyn_spg_ts_alloc()
!!----------------------------------------------------------------------
!! *** routine dyn_spg_ts_alloc ***
!!----------------------------------------------------------------------
INTEGER :: ierr(3)
!!----------------------------------------------------------------------
ierr(:) = 0
!
ALLOCATE( wgtbtp1(3*nn_e), wgtbtp2(3*nn_e), zwz(jpi,jpj), STAT=ierr(1) )
IF( ln_dynvor_een .OR. ln_dynvor_eeT ) &
& ALLOCATE( ftnw(jpi,jpj) , ftne(jpi,jpj) , ftsw(jpi,jpj) , ftse(jpi,jpj), STAT=ierr(2) )
!
ALLOCATE( un_adv(jpi,jpj), vn_adv(jpi,jpj) , STAT=ierr(3) )
!
dyn_spg_ts_alloc = MAXVAL( ierr(:) )
!
CALL mpp_sum( 'dynspg_ts', dyn_spg_ts_alloc )
IF( dyn_spg_ts_alloc /= 0 ) CALL ctl_stop( 'STOP', 'dyn_spg_ts_alloc: failed to allocate arrays' )
!
END FUNCTION dyn_spg_ts_alloc
SUBROUTINE dyn_spg_ts( kt, Kbb, Kmm, Krhs, puu, pvv, pssh, puu_b, pvv_b, Kaa, k_only_ADV )
!!----------------------------------------------------------------------
!!
!! ** Purpose : - Compute the now trend due to the explicit time stepping
!! of the quasi-linear barotropic system, and add it to the
!! general momentum trend.
!!
!! ** Method : - split-explicit schem (time splitting) :
!! Barotropic variables are advanced from internal time steps
!! "n" to "n+1" if ln_bt_fw=T
!! or from
!! "n-1" to "n+1" if ln_bt_fw=F
!! thanks to a generalized forward-backward time stepping (see ref. below).
!!
!! ** Action :
!! -Update the filtered free surface at step "n+1" : pssh(:,:,Kaa)
!! -Update filtered barotropic velocities at step "n+1" : puu_b(:,:,:,Kaa), vv_b(:,:,:,Kaa)
!! -Compute barotropic advective fluxes at step "n" : un_adv, vn_adv
!! These are used to advect tracers and are compliant with discrete
!! continuity equation taken at the baroclinic time steps. This
!! ensures tracers conservation.
!! - (puu(:,:,:,Krhs), pvv(:,:,:,Krhs)) momentum trend updated with barotropic component.
!!
!! References : Shchepetkin and McWilliams, Ocean Modelling, 2005.
!!---------------------------------------------------------------------
INTEGER , INTENT( in ) :: kt ! ocean time-step index
INTEGER , INTENT( in ) :: Kbb, Kmm, Krhs, Kaa ! ocean time level indices
REAL(wp), DIMENSION(jpi,jpj,jpk,jpt), INTENT(inout) :: puu, pvv ! ocean velocities and RHS of momentum equation
REAL(wp), DIMENSION(jpi,jpj,jpt) , INTENT(inout) :: pssh, puu_b, pvv_b ! SSH and barotropic velocities at main time levels
INTEGER , OPTIONAL , INTENT( in ) :: k_only_ADV ! only Advection in the RHS
!
INTEGER :: ji, jj, jk, jn ! dummy loop indices
LOGICAL :: ll_fw_start ! =T : forward integration
LOGICAL :: ll_init ! =T : special startup of 2d equations
INTEGER :: noffset ! local integers : time offset for bdy update
REAL(wp) :: r1_Dt_b, z1_hu, z1_hv ! local scalars
REAL(wp) :: za0, za1, za2, za3 ! - -
REAL(wp) :: zztmp, zldg ! - -
REAL(wp) :: zhu_bck, zhv_bck, zhdiv ! - -
REAL(wp) :: zun_save, zvn_save ! - -
REAL(wp), DIMENSION(jpi,jpj) :: zu_trd, zu_frc, zu_spg
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
REAL(wp), DIMENSION(jpi,jpj) :: zv_trd, zv_frc, zv_spg
REAL(wp), DIMENSION(jpi,jpj) :: zsshu_a, zhup2_e, zhtp2_e
REAL(wp), DIMENSION(jpi,jpj) :: zsshv_a, zhvp2_e, zsshp2_e
REAL(wp), DIMENSION(jpi,jpj) :: zCdU_u, zCdU_v ! top/bottom stress at u- & v-points
REAL(wp), DIMENSION(jpi,jpj) :: zhU, zhV ! fluxes
!!st#if defined key_qco
!!st REAL(wp), DIMENSION(jpi, jpj, jpk) :: ze3u, ze3v
!!st#endif
!
REAL(wp) :: zwdramp ! local scalar - only used if ln_wd_dl = .True.
INTEGER :: iwdg, jwdg, kwdg ! short-hand values for the indices of the output point
REAL(wp) :: zepsilon, zgamma ! - -
REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: zcpx, zcpy ! Wetting/Dying gravity filter coef.
REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: ztwdmask, zuwdmask, zvwdmask ! ROMS wetting and drying masks at t,u,v points
REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: zuwdav2, zvwdav2 ! averages over the sub-steps of zuwdmask and zvwdmask
REAL(wp) :: zt0substep ! Time of day at the beginning of the time substep
!!----------------------------------------------------------------------
!
IF( ln_wd_il ) ALLOCATE( zcpx(jpi,jpj), zcpy(jpi,jpj) )
! !* Allocate temporary arrays
IF( ln_wd_dl ) ALLOCATE( ztwdmask(jpi,jpj), zuwdmask(jpi,jpj), zvwdmask(jpi,jpj), zuwdav2(jpi,jpj), zvwdav2(jpi,jpj))
!
zwdramp = r_rn_wdmin1 ! simplest ramp
! zwdramp = 1._wp / (rn_wdmin2 - rn_wdmin1) ! more general ramp
! ! inverse of baroclinic time step
r1_Dt_b = 1._wp / rDt
!
ll_init = ln_bt_av ! if no time averaging, then no specific restart
ll_fw_start = .FALSE.
! ! time offset in steps for bdy data update
IF( .NOT.ln_bt_fw ) THEN ; noffset = - nn_e
ELSE ; noffset = 0
ENDIF
!
IF( kt == nit000 ) THEN !* initialisation
!
IF(lwp) WRITE(numout,*)
IF(lwp) WRITE(numout,*) 'dyn_spg_ts : surface pressure gradient trend'
IF(lwp) WRITE(numout,*) '~~~~~~~~~~ free surface with time splitting'
Loading
Loading full blame...