1 | MODULE stprk3 |
---|
2 | !!====================================================================== |
---|
3 | !! *** MODULE stpRK3 *** |
---|
4 | !! Time-stepping : manager of the ocean, tracer and ice time stepping |
---|
5 | !! using a 3rd order Rung Kuta with fixed or quasi-eulerian coordinate |
---|
6 | !!====================================================================== |
---|
7 | !! History : 4.5 ! 2021-01 (S. Techene, G. Madec, N. Ducousso, F. Lemarie) Original code |
---|
8 | !! NEMO |
---|
9 | !!---------------------------------------------------------------------- |
---|
10 | #if defined key_qco || defined key_linssh |
---|
11 | !!---------------------------------------------------------------------- |
---|
12 | !! 'key_qco' Quasi-Eulerian vertical coordinate |
---|
13 | !! OR |
---|
14 | !! 'key_linssh Fixed in time vertical coordinate |
---|
15 | !!---------------------------------------------------------------------- |
---|
16 | |
---|
17 | !!---------------------------------------------------------------------- |
---|
18 | !! stp_RK3 : NEMO 3rd order Runge-Kutta time-stepping |
---|
19 | !!---------------------------------------------------------------------- |
---|
20 | USE step_oce ! time stepping used modules |
---|
21 | USE domqco ! quasi-eulerian coordinate (dom_qco_r3c routine) |
---|
22 | USE stprk3_stg ! RK3 stages |
---|
23 | USE stp2d ! external mode solver |
---|
24 | USE dynspg_ts, ONLY: un_adv, vn_adv ! updated Kmm barotropic transport |
---|
25 | |
---|
26 | IMPLICIT NONE |
---|
27 | PRIVATE |
---|
28 | |
---|
29 | PUBLIC stp_RK3 ! called by nemogcm.F90 |
---|
30 | |
---|
31 | ! !** time level indices **! |
---|
32 | INTEGER, PUBLIC :: Nbb, Nnn, Naa, Nrhs !: used by nemo_init |
---|
33 | |
---|
34 | !! * Substitutions |
---|
35 | # include "do_loop_substitute.h90" |
---|
36 | # include "domzgr_substitute.h90" |
---|
37 | !!---------------------------------------------------------------------- |
---|
38 | !! NEMO/OCE 4.0 , NEMO Consortium (2018) |
---|
39 | !! $Id: step.F90 12377 2020-02-12 14:39:06Z acc $ |
---|
40 | !! Software governed by the CeCILL license (see ./LICENSE) |
---|
41 | !!---------------------------------------------------------------------- |
---|
42 | CONTAINS |
---|
43 | |
---|
44 | #if defined key_agrif |
---|
45 | RECURSIVE SUBROUTINE stp_RK3( ) |
---|
46 | INTEGER :: kstp ! ocean time-step index |
---|
47 | #else |
---|
48 | SUBROUTINE stp_RK3( kstp ) |
---|
49 | INTEGER, INTENT(in) :: kstp ! ocean time-step index |
---|
50 | #endif |
---|
51 | !!---------------------------------------------------------------------- |
---|
52 | !! *** ROUTINE stp_RK3 *** |
---|
53 | !! |
---|
54 | !! ** Purpose : - Time stepping of OCE (momentum and active tracer Eqs.) (RK3) |
---|
55 | !! - Time stepping of SI3 (dynamic and thermodynamic Eqs.) (FBS) |
---|
56 | !! - Time stepping of TRC (passive tracer Eqs.) |
---|
57 | !! |
---|
58 | !! ** Method : -1- Update forcings and data |
---|
59 | !! -2- Update ocean physics |
---|
60 | !! -3- Compute the after (Naa) ssh and velocity |
---|
61 | !! -4- diagnostics and output at Now (Nnn) |
---|
62 | !! -4- Compute the after (Naa) T-S |
---|
63 | !! -5- Update now |
---|
64 | !! -6- Update the horizontal velocity |
---|
65 | !! -7- Compute the diagnostics variables (rd,N2, hdiv,w) |
---|
66 | !! -8- Outputs and diagnostics |
---|
67 | !!---------------------------------------------------------------------- |
---|
68 | INTEGER :: ji, jj, jk, jtile ! dummy loop indice |
---|
69 | REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) :: zgdept |
---|
70 | !! --------------------------------------------------------------------- |
---|
71 | #if defined key_agrif |
---|
72 | IF( nstop > 0 ) RETURN ! avoid to go further if an error was detected during previous time step (child grid) |
---|
73 | kstp = nit000 + Agrif_Nb_Step() |
---|
74 | Kbb_a = Nbb ; Kmm_a = Nnn ; Krhs_a = Nrhs ! agrif_oce module copies of time level indices |
---|
75 | IF( lk_agrif_debug ) THEN |
---|
76 | IF( Agrif_Root() .AND. lwp) WRITE(*,*) '---' |
---|
77 | IF(lwp) WRITE(*,*) 'Grid Number', Agrif_Fixed(),' time step ', kstp, 'int tstep', Agrif_NbStepint() |
---|
78 | ENDIF |
---|
79 | IF( kstp == nit000 + 1 ) lk_agrif_fstep = .FALSE. |
---|
80 | # if defined key_xios |
---|
81 | IF( Agrif_Nbstepint() == 0 ) CALL iom_swap( cxios_context ) |
---|
82 | # endif |
---|
83 | #endif |
---|
84 | ! |
---|
85 | IF( ln_timing ) CALL timing_start('stp_RK3') |
---|
86 | ! |
---|
87 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
88 | ! update I/O and calendar |
---|
89 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
90 | ! |
---|
91 | IF( kstp == nit000 ) THEN ! initialize IOM context (must be done after nemo_init for AGRIF+XIOS+OASIS) |
---|
92 | CALL iom_init( cxios_context, ld_closedef=.FALSE. ) ! for model grid (including possible AGRIF zoom) |
---|
93 | IF( lk_diamlr ) CALL dia_mlr_iom_init ! with additional setup for multiple-linear-regression analysis |
---|
94 | CALL iom_init_closedef |
---|
95 | IF( ln_crs ) CALL iom_init( TRIM(cxios_context)//"_crs" ) ! for coarse grid |
---|
96 | ENDIF |
---|
97 | IF( kstp == nitrst .AND. lwxios ) THEN |
---|
98 | CALL iom_swap( cw_ocerst_cxt ) |
---|
99 | CALL iom_init_closedef( cw_ocerst_cxt ) |
---|
100 | CALL iom_setkt( kstp - nit000 + 1, cw_ocerst_cxt ) |
---|
101 | #if defined key_top |
---|
102 | CALL iom_swap( cw_toprst_cxt ) |
---|
103 | CALL iom_init_closedef( cw_toprst_cxt ) |
---|
104 | CALL iom_setkt( kstp - nit000 + 1, cw_toprst_cxt ) |
---|
105 | #endif |
---|
106 | ENDIF |
---|
107 | #if defined key_si3 |
---|
108 | IF( kstp + nn_fsbc - 1 == nitrst .AND. lwxios ) THEN |
---|
109 | CALL iom_swap( cw_icerst_cxt ) |
---|
110 | CALL iom_init_closedef( cw_icerst_cxt ) |
---|
111 | CALL iom_setkt( kstp - nit000 + 1, cw_icerst_cxt ) |
---|
112 | ENDIF |
---|
113 | #endif |
---|
114 | IF( kstp /= nit000 ) CALL day( kstp ) ! Calendar (day was already called at nit000 in day_init) |
---|
115 | CALL iom_setkt( kstp - nit000 + 1, cxios_context ) ! tell IOM we are at time step kstp |
---|
116 | IF( ln_crs ) CALL iom_setkt( kstp - nit000 + 1, TRIM(cxios_context)//"_crs" ) ! tell IOM we are at time step kstp |
---|
117 | |
---|
118 | |
---|
119 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
120 | ! Update external forcing (tides, open boundaries, ice shelf interaction and surface boundary condition (including sea-ice) |
---|
121 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
122 | IF( ln_tide ) CALL tide_update( kstp ) ! update tide potential |
---|
123 | IF( ln_apr_dyn ) CALL sbc_apr ( kstp ) ! atmospheric pressure (NB: call before bdy_dta which needs ssh_ib) |
---|
124 | IF( ln_bdy ) CALL bdy_dta ( kstp, Nnn ) ! update dynamic & tracer data at open boundaries |
---|
125 | IF( ln_isf ) CALL isf_stp ( kstp, Nnn ) ! update iceshelf geometry |
---|
126 | CALL sbc ( kstp, Nbb, Nnn ) ! Sea Boundary Condition (including sea-ice) |
---|
127 | |
---|
128 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
129 | ! Update stochastic parameters and random T/S fluctuations |
---|
130 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
131 | IF( ln_sto_eos ) CALL sto_par( kstp ) ! Stochastic parameters |
---|
132 | IF( ln_sto_eos ) CALL sto_pts( ts(:,:,:,:,Nnn) ) ! Random T/S fluctuations |
---|
133 | |
---|
134 | !!gm ocean physic computed at stage 3 ? |
---|
135 | |
---|
136 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
137 | ! Ocean physics update |
---|
138 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
139 | ! THERMODYNAMICS |
---|
140 | !!gm only Before is needed for bn2 and rab (except at stage 3 for n2) |
---|
141 | !!gm issue with Nnn used in rad(Nbb) |
---|
142 | CALL eos_rab( ts(:,:,:,:,Nbb), rab_b, Nnn ) ! before local thermal/haline expension ratio at T-points |
---|
143 | !!st CALL eos_rab( ts(:,:,:,:,Nnn), rab_n, Nnn ) ! now local thermal/haline expension ratio at T-points |
---|
144 | CALL bn2 ( ts(:,:,:,:,Nbb), rab_b, rn2b, Nnn ) ! before Brunt-Vaisala frequency |
---|
145 | !!st CALL bn2 ( ts(:,:,:,:,Nnn), rab_n, rn2, Nnn ) ! now Brunt-Vaisala frequency |
---|
146 | !!gm |
---|
147 | rab_n = rab_b |
---|
148 | rn2 = rn2b |
---|
149 | !!gm sh2 computed at the end of the time-step |
---|
150 | !!gm or call zdf_phy at the end ! |
---|
151 | ! VERTICAL PHYSICS |
---|
152 | !!st CALL zdf_phy( kstp, Nbb, Nnn, Nrhs ) ! vertical physics update (top/bot drag, avt, avs, avm + MLD) |
---|
153 | CALL zdf_phy( kstp, Nbb, Nbb, Nrhs ) ! vertical physics update (top/bot drag, avt, avs, avm + MLD) |
---|
154 | !!gm |
---|
155 | ! LATERAL PHYSICS |
---|
156 | ! |
---|
157 | IF( l_ldfslp ) THEN ! slope of lateral mixing |
---|
158 | !!gm gdep |
---|
159 | CALL eos( ts(:,:,:,:,Nbb), rhd, gdept_0(:,:,:) ) ! before in situ density |
---|
160 | |
---|
161 | IF( ln_zps .AND. .NOT. ln_isfcav) & |
---|
162 | & CALL zps_hde ( kstp, Nnn, jpts, ts(:,:,:,:,Nbb), gtsu, gtsv, & ! Partial steps: before horizontal gradient |
---|
163 | & rhd, gru , grv ) ! of t, s, rd at the last ocean level |
---|
164 | |
---|
165 | IF( ln_zps .AND. ln_isfcav) & |
---|
166 | & CALL zps_hde_isf( kstp, Nnn, jpts, ts(:,:,:,:,Nbb), gtsu, gtsv, gtui, gtvi, & ! Partial steps for top cell (ISF) |
---|
167 | & rhd, gru , grv , grui, grvi ) ! of t, s, rd at the first ocean level |
---|
168 | IF( ln_traldf_triad ) THEN |
---|
169 | CALL ldf_slp_triad( kstp, Nbb, Nnn ) ! before slope for triad operator |
---|
170 | ELSE |
---|
171 | CALL ldf_slp ( kstp, rhd, rn2b, Nbb, Nnn ) ! before slope for standard operator |
---|
172 | ENDIF |
---|
173 | ENDIF |
---|
174 | ! ! eddy diffusivity coeff. |
---|
175 | IF( l_ldftra_time .OR. l_ldfeiv_time ) CALL ldf_tra( kstp, Nbb, Nnn ) ! and/or eiv coeff. |
---|
176 | IF( l_ldfdyn_time ) CALL ldf_dyn( kstp, Nbb ) ! eddy viscosity coeff. |
---|
177 | |
---|
178 | |
---|
179 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
180 | ! RK3 : single first external mode computation |
---|
181 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
182 | |
---|
183 | CALL stp_2D( kstp, Nbb, Nbb, Naa, Nrhs ) ! out: ssh, (uu_b,vv_b) and (un_adv,vn_adv) at Naa |
---|
184 | |
---|
185 | |
---|
186 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
187 | ! RK3 time integration |
---|
188 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
189 | |
---|
190 | ! Stage 1 : |
---|
191 | CALL stp_RK3_stg( 1, kstp, Nbb, Nbb, Nrhs, Naa ) |
---|
192 | ! |
---|
193 | Nrhs = Nnn ; Nnn = Naa ; Naa = Nrhs ! Swap: Nbb unchanged, Nnn <==> Naa |
---|
194 | ! |
---|
195 | ! Stage 2 : |
---|
196 | CALL stp_RK3_stg( 2, kstp, Nbb, Nnn, Nrhs, Naa ) |
---|
197 | ! |
---|
198 | Nrhs = Nnn ; Nnn = Naa ; Naa = Nrhs ! Swap: Nbb unchanged, Nnn <==> Naa |
---|
199 | ! |
---|
200 | ! Stage 3 : |
---|
201 | CALL stp_RK3_stg( 3, kstp, Nbb, Nnn, Nrhs, Naa ) |
---|
202 | ! |
---|
203 | Nrhs = Nbb ; Nbb = Naa ; Naa = Nrhs ! Swap: Nnn unchanged, Nbb <==> Naa |
---|
204 | |
---|
205 | ! Swap time levels |
---|
206 | |
---|
207 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
208 | ! diagnostics and outputs |
---|
209 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
210 | |
---|
211 | !==>>> at Nbb no more Nnn |
---|
212 | |
---|
213 | IF( ln_diacfl ) CALL dia_cfl ( kstp, Nbb ) ! Courant number diagnostics |
---|
214 | CALL dia_hth ( kstp, Nbb ) ! Thermocline depth (20 degres isotherm depth) |
---|
215 | IF( ln_diadct ) CALL dia_dct ( kstp, Nbb ) ! Transports |
---|
216 | !!st CALL dia_ar5 ( kstp, Nbb ) ! ar5 diag |
---|
217 | CALL dia_ptr ( kstp, Nbb ) ! Poleward adv/ldf TRansports diagnostics |
---|
218 | CALL dia_wri ( kstp, Nbb ) ! ocean model: outputs |
---|
219 | IF( ln_crs ) CALL crs_fld ( kstp, Nbb ) ! ocean model: online field coarsening & output |
---|
220 | IF( lk_diadetide ) CALL dia_detide( kstp ) ! Weights computation for daily detiding of model diagnostics |
---|
221 | IF( lk_diamlr ) CALL dia_mlr ! Update time used in multiple-linear-regression analysis |
---|
222 | |
---|
223 | !!====>>>> to be modified for RK3 |
---|
224 | ! IF( ln_floats ) CALL flo_stp ( kstp, Nbb, Nnn ) ! drifting Floats |
---|
225 | ! IF( ln_diahsb ) CALL dia_hsb ( kstp, Nbb, Nnn ) ! - ML - global conservation diagnostics |
---|
226 | !!st |
---|
227 | IF( ln_diahsb ) CALL dia_hsb ( kstp, Nbb, Nbb ) ! - ML - global conservation diagnostics |
---|
228 | !!st |
---|
229 | |
---|
230 | !!gm : This does not only concern the dynamics ==>>> add a new title |
---|
231 | !!gm2: why ouput restart before AGRIF update? |
---|
232 | !! |
---|
233 | !!jc: That would be better, but see comment above |
---|
234 | !! |
---|
235 | !!====>>>> to be modified for RK3 |
---|
236 | IF( lrst_oce ) CALL rst_write ( kstp, Nbb, Nnn ) ! write output ocean restart file |
---|
237 | IF( ln_sto_eos ) CALL sto_rst_write( kstp ) ! write restart file for stochastic parameters |
---|
238 | |
---|
239 | #if defined key_agrif |
---|
240 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
241 | ! AGRIF recursive integration |
---|
242 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
243 | Kbb_a = Nbb; Kmm_a = Nnn; Krhs_a = Nrhs ! agrif_oce module copies of time level indices |
---|
244 | CALL Agrif_Integrate_ChildGrids( stp_RK3 ) ! allows to finish all the Child Grids before updating |
---|
245 | |
---|
246 | #endif |
---|
247 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
248 | ! Control |
---|
249 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
250 | CALL stp_ctl ( kstp, Nbb ) |
---|
251 | #if defined key_agrif |
---|
252 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
253 | ! AGRIF update |
---|
254 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
255 | IF( Agrif_NbStepint() == 0 .AND. nstop == 0 ) & |
---|
256 | & CALL Agrif_update_all( ) ! Update all components |
---|
257 | |
---|
258 | #endif |
---|
259 | IF( ln_diaobs .AND. nstop == 0 ) & |
---|
260 | & CALL dia_obs( kstp, Nnn ) ! obs-minus-model (assimilation) diags (after dynamics update) |
---|
261 | |
---|
262 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
263 | ! File manipulation at the end of the first time step |
---|
264 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
265 | IF( kstp == nit000 ) THEN ! 1st time step only |
---|
266 | CALL iom_close( numror ) ! close input ocean restart file |
---|
267 | IF( lrxios ) CALL iom_context_finalize( cr_ocerst_cxt ) |
---|
268 | IF(lwm) CALL FLUSH ( numond ) ! flush output namelist oce |
---|
269 | IF(lwm .AND. numoni /= -1 ) CALL FLUSH ( numoni ) ! flush output namelist ice (if exist) |
---|
270 | ENDIF |
---|
271 | |
---|
272 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
273 | ! Coupled mode |
---|
274 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
275 | IF( lk_oasis .AND. nstop == 0 ) CALL sbc_cpl_snd( kstp, Nbb, Nnn ) ! coupled mode : field exchanges |
---|
276 | ! |
---|
277 | #if defined key_xios |
---|
278 | !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> |
---|
279 | ! Finalize contextes if end of simulation or error detected |
---|
280 | !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< |
---|
281 | IF( kstp == nitend .OR. nstop > 0 ) THEN |
---|
282 | CALL iom_context_finalize( cxios_context ) ! needed for XIOS+AGRIF |
---|
283 | IF( ln_crs ) CALL iom_context_finalize( trim(cxios_context)//"_crs" ) ! |
---|
284 | ENDIF |
---|
285 | #endif |
---|
286 | ! |
---|
287 | IF( ln_timing ) CALL timing_stop('stp_RK3') |
---|
288 | ! |
---|
289 | END SUBROUTINE stp_RK3 |
---|
290 | |
---|
291 | |
---|
292 | SUBROUTINE mlf_baro_corr( Kmm, Kaa, puu, pvv ) |
---|
293 | !!---------------------------------------------------------------------- |
---|
294 | !! *** ROUTINE mlf_baro_corr *** |
---|
295 | !! |
---|
296 | !! ** Purpose : Finalize after horizontal velocity. |
---|
297 | !! |
---|
298 | !! ** Method : * Ensure after velocities transport matches time splitting |
---|
299 | !! estimate (ln_dynspg_ts=T) |
---|
300 | !! |
---|
301 | !! ** Action : puu(Kmm),pvv(Kmm) updated now horizontal velocity (ln_bt_fw=F) |
---|
302 | !! puu(Kaa),pvv(Kaa) after horizontal velocity |
---|
303 | !!---------------------------------------------------------------------- |
---|
304 | USE dynspg_ts, ONLY : un_adv, vn_adv ! updated Kmm barotropic transport |
---|
305 | !! |
---|
306 | INTEGER , INTENT(in ) :: Kmm, Kaa ! before and after time level indices |
---|
307 | REAL(wp), DIMENSION(jpi,jpj,jpk,jpt), INTENT(inout) :: puu, pvv ! velocities |
---|
308 | ! |
---|
309 | INTEGER :: jk ! dummy loop indices |
---|
310 | REAL(wp), DIMENSION(jpi,jpj) :: zue, zve |
---|
311 | !!---------------------------------------------------------------------- |
---|
312 | |
---|
313 | ! Ensure below that barotropic velocities match time splitting estimate |
---|
314 | ! Compute actual transport and replace it with ts estimate at "after" time step |
---|
315 | zue(:,:) = e3u(:,:,1,Kaa) * puu(:,:,1,Kaa) * umask(:,:,1) |
---|
316 | zve(:,:) = e3v(:,:,1,Kaa) * pvv(:,:,1,Kaa) * vmask(:,:,1) |
---|
317 | DO jk = 2, jpkm1 |
---|
318 | zue(:,:) = zue(:,:) + e3u(:,:,jk,Kaa) * puu(:,:,jk,Kaa) * umask(:,:,jk) |
---|
319 | zve(:,:) = zve(:,:) + e3v(:,:,jk,Kaa) * pvv(:,:,jk,Kaa) * vmask(:,:,jk) |
---|
320 | END DO |
---|
321 | DO jk = 1, jpkm1 |
---|
322 | puu(:,:,jk,Kaa) = ( puu(:,:,jk,Kaa) - zue(:,:) * r1_hu(:,:,Kaa) + uu_b(:,:,Kaa) ) * umask(:,:,jk) |
---|
323 | pvv(:,:,jk,Kaa) = ( pvv(:,:,jk,Kaa) - zve(:,:) * r1_hv(:,:,Kaa) + vv_b(:,:,Kaa) ) * vmask(:,:,jk) |
---|
324 | END DO |
---|
325 | ! |
---|
326 | !!st IF( .NOT.ln_bt_fw ) THEN |
---|
327 | !!st ! Remove advective velocity from "now velocities" |
---|
328 | !!st ! prior to asselin filtering |
---|
329 | !!st ! In the forward case, this is done below after asselin filtering |
---|
330 | !!st ! so that asselin contribution is removed at the same time |
---|
331 | !!st DO jk = 1, jpkm1 |
---|
332 | !!st puu(:,:,jk,Kmm) = ( puu(:,:,jk,Kmm) - un_adv(:,:)*r1_hu(:,:,Kmm) + uu_b(:,:,Kmm) )*umask(:,:,jk) |
---|
333 | !!st pvv(:,:,jk,Kmm) = ( pvv(:,:,jk,Kmm) - vn_adv(:,:)*r1_hv(:,:,Kmm) + vv_b(:,:,Kmm) )*vmask(:,:,jk) |
---|
334 | !!st END DO |
---|
335 | !!st ENDIF |
---|
336 | ! |
---|
337 | END SUBROUTINE mlf_baro_corr |
---|
338 | |
---|
339 | |
---|
340 | SUBROUTINE finalize_lbc( kt, Kbb, Kaa, puu, pvv, pts ) |
---|
341 | !!---------------------------------------------------------------------- |
---|
342 | !! *** ROUTINE finalize_lbc *** |
---|
343 | !! |
---|
344 | !! ** Purpose : Apply the boundary condition on the after velocity |
---|
345 | !! |
---|
346 | !! ** Method : * Apply lateral boundary conditions on after velocity |
---|
347 | !! at the local domain boundaries through lbc_lnk call, |
---|
348 | !! at the one-way open boundaries (ln_bdy=T), |
---|
349 | !! at the AGRIF zoom boundaries (lk_agrif=T) |
---|
350 | !! |
---|
351 | !! ** Action : puu(Kaa),pvv(Kaa) after horizontal velocity and tracers |
---|
352 | !!---------------------------------------------------------------------- |
---|
353 | #if defined key_agrif |
---|
354 | USE agrif_oce_interp |
---|
355 | #endif |
---|
356 | USE bdydyn ! ocean open boundary conditions (define bdy_dyn) |
---|
357 | !! |
---|
358 | INTEGER , INTENT(in ) :: kt ! ocean time-step index |
---|
359 | INTEGER , INTENT(in ) :: Kbb, Kaa ! before and after time level indices |
---|
360 | REAL(wp), DIMENSION(jpi,jpj,jpk,jpt) , INTENT(inout) :: puu, pvv ! velocities to be time filtered |
---|
361 | REAL(wp), DIMENSION(jpi,jpj,jpk,jpts,jpt), INTENT(inout) :: pts ! active tracers |
---|
362 | !!---------------------------------------------------------------------- |
---|
363 | ! |
---|
364 | ! Update after tracer and velocity on domain lateral boundaries |
---|
365 | ! |
---|
366 | # if defined key_agrif |
---|
367 | CALL Agrif_tra !* AGRIF zoom boundaries |
---|
368 | CALL Agrif_dyn( kt ) |
---|
369 | # endif |
---|
370 | ! ! local domain boundaries (T-point, unchanged sign) |
---|
371 | CALL lbc_lnk_multi( 'finalize_lbc', puu(:,:,:, Kaa), 'U', -1., pvv(:,:,: ,Kaa), 'V', -1. & |
---|
372 | & , pts(:,:,:,jp_tem,Kaa), 'T', 1., pts(:,:,:,jp_sal,Kaa), 'T', 1. ) |
---|
373 | ! |
---|
374 | ! !* BDY open boundaries |
---|
375 | IF( ln_bdy ) THEN |
---|
376 | CALL bdy_tra( kt, Kbb, pts, Kaa ) |
---|
377 | IF( ln_dynspg_exp ) CALL bdy_dyn( kt, Kbb, puu, pvv, Kaa ) |
---|
378 | IF( ln_dynspg_ts ) CALL bdy_dyn( kt, Kbb, puu, pvv, Kaa, dyn3d_only=.true. ) |
---|
379 | ENDIF |
---|
380 | ! |
---|
381 | END SUBROUTINE finalize_lbc |
---|
382 | |
---|
383 | #else |
---|
384 | !!---------------------------------------------------------------------- |
---|
385 | !! default option EMPTY MODULE qco not activated |
---|
386 | !!---------------------------------------------------------------------- |
---|
387 | #endif |
---|
388 | |
---|
389 | !!====================================================================== |
---|
390 | END MODULE stprk3 |
---|