1 | MODULE tranxt |
---|
2 | !!====================================================================== |
---|
3 | !! *** MODULE tranxt *** |
---|
4 | !! Ocean active tracers: time stepping on temperature and salinity |
---|
5 | !!====================================================================== |
---|
6 | !! History : OPA ! 1991-11 (G. Madec) Original code |
---|
7 | !! 7.0 ! 1993-03 (M. Guyon) symetrical conditions |
---|
8 | !! 8.0 ! 1996-02 (G. Madec & M. Imbard) opa release 8.0 |
---|
9 | !! - ! 1996-04 (A. Weaver) Euler forward step |
---|
10 | !! 8.2 ! 1999-02 (G. Madec, N. Grima) semi-implicit pressure grad. |
---|
11 | !! NEMO 1.0 ! 2002-08 (G. Madec) F90: Free form and module |
---|
12 | !! - ! 2002-11 (C. Talandier, A-M Treguier) Open boundaries |
---|
13 | !! - ! 2005-04 (C. Deltel) Add Asselin trend in the ML budget |
---|
14 | !! 2.0 ! 2006-02 (L. Debreu, C. Mazauric) Agrif implementation |
---|
15 | !! 3.0 ! 2008-06 (G. Madec) time stepping always done in trazdf |
---|
16 | !! 3.1 ! 2009-02 (G. Madec, R. Benshila) re-introduce the vvl option |
---|
17 | !! 3.3 ! 2010-04 (M. Leclair, G. Madec) semi-implicit hpg with asselin filter + modified LF-RA |
---|
18 | !! - ! 2010-05 (C. Ethe, G. Madec) merge TRC-TRA |
---|
19 | !!---------------------------------------------------------------------- |
---|
20 | |
---|
21 | !!---------------------------------------------------------------------- |
---|
22 | !! tra_nxt : time stepping on tracers |
---|
23 | !! tra_nxt_fix : time stepping on tracers : fixed volume case |
---|
24 | !! tra_nxt_vvl : time stepping on tracers : variable volume case |
---|
25 | !!---------------------------------------------------------------------- |
---|
26 | USE oce ! ocean dynamics and tracers variables |
---|
27 | USE dom_oce ! ocean space and time domain variables |
---|
28 | USE sbc_oce ! surface boundary condition: ocean |
---|
29 | USE sbcrnf ! river runoffs |
---|
30 | USE sbcisf ! ice shelf melting |
---|
31 | USE zdf_oce ! ocean vertical mixing |
---|
32 | USE domvvl ! variable volume |
---|
33 | USE trd_oce ! trends: ocean variables |
---|
34 | USE trdtra ! trends manager: tracers |
---|
35 | USE traqsr ! penetrative solar radiation (needed for nksr) |
---|
36 | USE phycst ! physical constant |
---|
37 | USE ldftra ! lateral physics : tracers |
---|
38 | USE ldfslp ! lateral physics : slopes |
---|
39 | USE bdy_oce , ONLY : ln_bdy |
---|
40 | USE bdytra ! open boundary condition (bdy_tra routine) |
---|
41 | ! |
---|
42 | USE in_out_manager ! I/O manager |
---|
43 | USE lbclnk ! ocean lateral boundary conditions (or mpp link) |
---|
44 | USE prtctl ! Print control |
---|
45 | USE timing ! Timing |
---|
46 | #if defined key_agrif |
---|
47 | USE agrif_oce_interp |
---|
48 | #endif |
---|
49 | |
---|
50 | IMPLICIT NONE |
---|
51 | PRIVATE |
---|
52 | |
---|
53 | PUBLIC tra_nxt ! routine called by step.F90 |
---|
54 | PUBLIC tra_nxt_fix ! to be used in trcnxt |
---|
55 | PUBLIC tra_nxt_vvl ! to be used in trcnxt |
---|
56 | |
---|
57 | !! * Substitutions |
---|
58 | # include "vectopt_loop_substitute.h90" |
---|
59 | !!---------------------------------------------------------------------- |
---|
60 | !! NEMO/OCE 4.0 , NEMO Consortium (2018) |
---|
61 | !! $Id$ |
---|
62 | !! Software governed by the CeCILL license (see ./LICENSE) |
---|
63 | !!---------------------------------------------------------------------- |
---|
64 | CONTAINS |
---|
65 | |
---|
66 | SUBROUTINE tra_nxt( kt, Kbb, Kmm, Krhs ) |
---|
67 | !!---------------------------------------------------------------------- |
---|
68 | !! *** ROUTINE tranxt *** |
---|
69 | !! |
---|
70 | !! ** Purpose : Apply the boundary condition on the after temperature |
---|
71 | !! and salinity fields, achieved the time stepping by adding |
---|
72 | !! the Asselin filter on now fields and swapping the fields. |
---|
73 | !! |
---|
74 | !! ** Method : At this stage of the computation, ta and sa are the |
---|
75 | !! after temperature and salinity as the time stepping has |
---|
76 | !! been performed in trazdf_imp or trazdf_exp module. |
---|
77 | !! |
---|
78 | !! - Apply lateral boundary conditions on (ta,sa) |
---|
79 | !! at the local domain boundaries through lbc_lnk call, |
---|
80 | !! at the one-way open boundaries (ln_bdy=T), |
---|
81 | !! at the AGRIF zoom boundaries (lk_agrif=T) |
---|
82 | !! |
---|
83 | !! - Update lateral boundary conditions on AGRIF children |
---|
84 | !! domains (lk_agrif=T) |
---|
85 | !! |
---|
86 | !! ** Action : - ts(Kbb) & ts(Kmm) ready for the next time step |
---|
87 | !!---------------------------------------------------------------------- |
---|
88 | INTEGER, INTENT(in) :: kt ! ocean time-step index |
---|
89 | INTEGER, INTENT(in) :: Kbb, Kmm, Krhs ! time level indices |
---|
90 | !! |
---|
91 | INTEGER :: ji, jj, jk, jn ! dummy loop indices |
---|
92 | REAL(wp) :: zfact ! local scalars |
---|
93 | REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) :: ztrdt, ztrds |
---|
94 | !!---------------------------------------------------------------------- |
---|
95 | ! |
---|
96 | IF( ln_timing ) CALL timing_start( 'tra_nxt') |
---|
97 | ! |
---|
98 | IF( kt == nit000 ) THEN |
---|
99 | IF(lwp) WRITE(numout,*) |
---|
100 | IF(lwp) WRITE(numout,*) 'tra_nxt : achieve the time stepping by Asselin filter and array swap' |
---|
101 | IF(lwp) WRITE(numout,*) '~~~~~~~' |
---|
102 | ENDIF |
---|
103 | |
---|
104 | ! Update after tracer on domain lateral boundaries |
---|
105 | ! |
---|
106 | #if defined key_agrif |
---|
107 | CALL Agrif_tra ! AGRIF zoom boundaries |
---|
108 | #endif |
---|
109 | ! ! local domain boundaries (T-point, unchanged sign) |
---|
110 | CALL lbc_lnk_multi( 'tranxt', ts(:,:,:,jp_tem,Krhs), 'T', 1., ts(:,:,:,jp_sal,Krhs), 'T', 1. ) |
---|
111 | ! |
---|
112 | IF( ln_bdy ) CALL bdy_tra( kt ) ! BDY open boundaries |
---|
113 | |
---|
114 | ! set time step size (Euler/Leapfrog) |
---|
115 | IF( neuler == 0 .AND. kt == nit000 ) THEN ; r2dt = rdt ! at nit000 (Euler) |
---|
116 | ELSEIF( kt <= nit000 + 1 ) THEN ; r2dt = 2._wp* rdt ! at nit000 or nit000+1 (Leapfrog) |
---|
117 | ENDIF |
---|
118 | |
---|
119 | ! trends computation initialisation |
---|
120 | IF( l_trdtra ) THEN |
---|
121 | ALLOCATE( ztrdt(jpi,jpj,jpk) , ztrds(jpi,jpj,jpk) ) |
---|
122 | ztrdt(:,:,jpk) = 0._wp |
---|
123 | ztrds(:,:,jpk) = 0._wp |
---|
124 | IF( ln_traldf_iso ) THEN ! diagnose the "pure" Kz diffusive trend |
---|
125 | CALL trd_tra( kt, Kmm, Krhs, 'TRA', jp_tem, jptra_zdfp, ztrdt ) |
---|
126 | CALL trd_tra( kt, Kmm, Krhs, 'TRA', jp_sal, jptra_zdfp, ztrds ) |
---|
127 | ENDIF |
---|
128 | ! total trend for the non-time-filtered variables. |
---|
129 | zfact = 1.0 / rdt |
---|
130 | ! G Nurser 23 Mar 2017. Recalculate trend as Delta(e3t*T)/e3tn; e3tn cancel from ts(Kmm) terms |
---|
131 | DO jk = 1, jpkm1 |
---|
132 | ztrdt(:,:,jk) = ( ts(:,:,jk,jp_tem,Krhs)*e3t(:,:,jk,Krhs) / e3t(:,:,jk,Kmm) - ts(:,:,jk,jp_tem,Kmm)) * zfact |
---|
133 | ztrds(:,:,jk) = ( ts(:,:,jk,jp_sal,Krhs)*e3t(:,:,jk,Krhs) / e3t(:,:,jk,Kmm) - ts(:,:,jk,jp_sal,Kmm)) * zfact |
---|
134 | END DO |
---|
135 | CALL trd_tra( kt, Kmm, Krhs, 'TRA', jp_tem, jptra_tot, ztrdt ) |
---|
136 | CALL trd_tra( kt, Kmm, Krhs, 'TRA', jp_sal, jptra_tot, ztrds ) |
---|
137 | IF( ln_linssh ) THEN ! linear sea surface height only |
---|
138 | ! Store now fields before applying the Asselin filter |
---|
139 | ! in order to calculate Asselin filter trend later. |
---|
140 | ztrdt(:,:,:) = ts(:,:,:,jp_tem,Kmm) |
---|
141 | ztrds(:,:,:) = ts(:,:,:,jp_sal,Kmm) |
---|
142 | ENDIF |
---|
143 | ENDIF |
---|
144 | |
---|
145 | IF( neuler == 0 .AND. kt == nit000 ) THEN ! Euler time-stepping at first time-step (only swap) |
---|
146 | DO jn = 1, jpts |
---|
147 | DO jk = 1, jpkm1 |
---|
148 | ts(:,:,jk,jn,Kmm) = ts(:,:,jk,jn,Krhs) |
---|
149 | END DO |
---|
150 | END DO |
---|
151 | IF (l_trdtra .AND. .NOT. ln_linssh ) THEN ! Zero Asselin filter contribution must be explicitly written out since for vvl |
---|
152 | ! ! Asselin filter is output by tra_nxt_vvl that is not called on this time step |
---|
153 | ztrdt(:,:,:) = 0._wp |
---|
154 | ztrds(:,:,:) = 0._wp |
---|
155 | CALL trd_tra( kt, Kmm, Krhs, 'TRA', jp_tem, jptra_atf, ztrdt ) |
---|
156 | CALL trd_tra( kt, Kmm, Krhs, 'TRA', jp_sal, jptra_atf, ztrds ) |
---|
157 | END IF |
---|
158 | ! |
---|
159 | ELSE ! Leap-Frog + Asselin filter time stepping |
---|
160 | ! |
---|
161 | IF( ln_linssh ) THEN ; CALL tra_nxt_fix( kt, Kmm, nit000, 'TRA', & |
---|
162 | & ts(:,:,:,:,Kbb), ts(:,:,:,:,Kmm), ts(:,:,:,:,Krhs), jpts ) ! linear free surface |
---|
163 | ELSE ; CALL tra_nxt_vvl( kt, Kbb, Kmm, Krhs, nit000, rdt, 'TRA', & |
---|
164 | & ts(:,:,:,:,Kbb), ts(:,:,:,:,Kmm), ts(:,:,:,:,Krhs), & |
---|
165 | & sbc_tsc , sbc_tsc_b , jpts ) ! non-linear free surface |
---|
166 | ENDIF |
---|
167 | ! |
---|
168 | CALL lbc_lnk_multi( 'tranxt', ts(:,:,:,jp_tem,Kbb) , 'T', 1., ts(:,:,:,jp_sal,Kbb) , 'T', 1., & |
---|
169 | & ts(:,:,:,jp_tem,Kmm) , 'T', 1., ts(:,:,:,jp_sal,Kmm) , 'T', 1., & |
---|
170 | & ts(:,:,:,jp_tem,Krhs), 'T', 1., ts(:,:,:,jp_sal,Krhs), 'T', 1. ) |
---|
171 | ! |
---|
172 | ENDIF |
---|
173 | ! |
---|
174 | IF( l_trdtra .AND. ln_linssh ) THEN ! trend of the Asselin filter (tb filtered - tb)/dt |
---|
175 | zfact = 1._wp / r2dt |
---|
176 | DO jk = 1, jpkm1 |
---|
177 | ztrdt(:,:,jk) = ( ts(:,:,jk,jp_tem,Kbb) - ztrdt(:,:,jk) ) * zfact |
---|
178 | ztrds(:,:,jk) = ( ts(:,:,jk,jp_sal,Kbb) - ztrds(:,:,jk) ) * zfact |
---|
179 | END DO |
---|
180 | CALL trd_tra( kt, Kmm, Krhs, 'TRA', jp_tem, jptra_atf, ztrdt ) |
---|
181 | CALL trd_tra( kt, Kmm, Krhs, 'TRA', jp_sal, jptra_atf, ztrds ) |
---|
182 | END IF |
---|
183 | IF( l_trdtra ) DEALLOCATE( ztrdt , ztrds ) |
---|
184 | ! |
---|
185 | ! ! control print |
---|
186 | IF(ln_ctl) CALL prt_ctl( tab3d_1=ts(:,:,:,jp_tem,Kmm), clinfo1=' nxt - Tn: ', mask1=tmask, & |
---|
187 | & tab3d_2=ts(:,:,:,jp_sal,Kmm), clinfo2= ' Sn: ', mask2=tmask ) |
---|
188 | ! |
---|
189 | IF( ln_timing ) CALL timing_stop('tra_nxt') |
---|
190 | ! |
---|
191 | END SUBROUTINE tra_nxt |
---|
192 | |
---|
193 | |
---|
194 | SUBROUTINE tra_nxt_fix( kt, Kmm, kit000, cdtype, ptb, ptn, pta, kjpt ) |
---|
195 | !!---------------------------------------------------------------------- |
---|
196 | !! *** ROUTINE tra_nxt_fix *** |
---|
197 | !! |
---|
198 | !! ** Purpose : fixed volume: apply the Asselin time filter and |
---|
199 | !! swap the tracer fields. |
---|
200 | !! |
---|
201 | !! ** Method : - Apply a Asselin time filter on now fields. |
---|
202 | !! - swap tracer fields to prepare the next time_step. |
---|
203 | !! |
---|
204 | !! ** Action : - ptb & ptn ready for the next time step |
---|
205 | !!---------------------------------------------------------------------- |
---|
206 | INTEGER , INTENT(in ) :: kt ! ocean time-step index |
---|
207 | INTEGER , INTENT(in ) :: Kmm ! time level index |
---|
208 | INTEGER , INTENT(in ) :: kit000 ! first time step index |
---|
209 | CHARACTER(len=3) , INTENT(in ) :: cdtype ! =TRA or TRC (tracer indicator) |
---|
210 | INTEGER , INTENT(in ) :: kjpt ! number of tracers |
---|
211 | REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) :: ptb ! before tracer fields |
---|
212 | REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) :: ptn ! now tracer fields |
---|
213 | REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) :: pta ! tracer trend |
---|
214 | ! |
---|
215 | INTEGER :: ji, jj, jk, jn ! dummy loop indices |
---|
216 | REAL(wp) :: ztn, ztd ! local scalars |
---|
217 | !!---------------------------------------------------------------------- |
---|
218 | ! |
---|
219 | IF( kt == kit000 ) THEN |
---|
220 | IF(lwp) WRITE(numout,*) |
---|
221 | IF(lwp) WRITE(numout,*) 'tra_nxt_fix : time stepping', cdtype |
---|
222 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~~' |
---|
223 | ENDIF |
---|
224 | ! |
---|
225 | DO jn = 1, kjpt |
---|
226 | ! |
---|
227 | DO jk = 1, jpkm1 |
---|
228 | DO jj = 2, jpjm1 |
---|
229 | DO ji = fs_2, fs_jpim1 |
---|
230 | ztn = ptn(ji,jj,jk,jn) |
---|
231 | ztd = pta(ji,jj,jk,jn) - 2._wp * ztn + ptb(ji,jj,jk,jn) ! time laplacian on tracers |
---|
232 | ! |
---|
233 | ptb(ji,jj,jk,jn) = ztn + atfp * ztd ! ptb <-- filtered ptn |
---|
234 | ptn(ji,jj,jk,jn) = pta(ji,jj,jk,jn) ! ptn <-- pta |
---|
235 | END DO |
---|
236 | END DO |
---|
237 | END DO |
---|
238 | ! |
---|
239 | END DO |
---|
240 | ! |
---|
241 | END SUBROUTINE tra_nxt_fix |
---|
242 | |
---|
243 | |
---|
244 | SUBROUTINE tra_nxt_vvl( kt, Kbb, Kmm, Krhs, kit000, p2dt, cdtype, ptb, ptn, pta, psbc_tc, psbc_tc_b, kjpt ) |
---|
245 | !!---------------------------------------------------------------------- |
---|
246 | !! *** ROUTINE tra_nxt_vvl *** |
---|
247 | !! |
---|
248 | !! ** Purpose : Time varying volume: apply the Asselin time filter |
---|
249 | !! and swap the tracer fields. |
---|
250 | !! |
---|
251 | !! ** Method : - Apply a thickness weighted Asselin time filter on now fields. |
---|
252 | !! - swap tracer fields to prepare the next time_step. |
---|
253 | !! tb = ( e3t(Kmm)*tn + atfp*[ e3t(Kbb)*tb - 2 e3t(Kmm)*tn + e3t_a*ta ] ) |
---|
254 | !! /( e3t(Kmm) + atfp*[ e3t(Kbb) - 2 e3t(Kmm) + e3t(Krhs) ] ) |
---|
255 | !! tn = ta |
---|
256 | !! |
---|
257 | !! ** Action : - ptb & ptn ready for the next time step |
---|
258 | !!---------------------------------------------------------------------- |
---|
259 | INTEGER , INTENT(in ) :: kt ! ocean time-step index |
---|
260 | INTEGER , INTENT(in ) :: Kbb, Kmm, Krhs ! time level indices |
---|
261 | INTEGER , INTENT(in ) :: kit000 ! first time step index |
---|
262 | REAL(wp) , INTENT(in ) :: p2dt ! time-step |
---|
263 | CHARACTER(len=3) , INTENT(in ) :: cdtype ! =TRA or TRC (tracer indicator) |
---|
264 | INTEGER , INTENT(in ) :: kjpt ! number of tracers |
---|
265 | REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) :: ptb ! before tracer fields |
---|
266 | REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) :: ptn ! now tracer fields |
---|
267 | REAL(wp), DIMENSION(jpi,jpj,jpk,kjpt), INTENT(inout) :: pta ! tracer trend |
---|
268 | REAL(wp), DIMENSION(jpi,jpj ,kjpt), INTENT(in ) :: psbc_tc ! surface tracer content |
---|
269 | REAL(wp), DIMENSION(jpi,jpj ,kjpt), INTENT(in ) :: psbc_tc_b ! before surface tracer content |
---|
270 | ! |
---|
271 | LOGICAL :: ll_traqsr, ll_rnf, ll_isf ! local logical |
---|
272 | INTEGER :: ji, jj, jk, jn ! dummy loop indices |
---|
273 | REAL(wp) :: zfact, zfact1, ztc_a , ztc_n , ztc_b , ztc_f , ztc_d ! local scalar |
---|
274 | REAL(wp) :: zfact2, ze3t_b, ze3t_n, ze3t_a, ze3t_f, ze3t_d ! - - |
---|
275 | REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:) :: ztrd_atf |
---|
276 | !!---------------------------------------------------------------------- |
---|
277 | ! |
---|
278 | IF( kt == kit000 ) THEN |
---|
279 | IF(lwp) WRITE(numout,*) |
---|
280 | IF(lwp) WRITE(numout,*) 'tra_nxt_vvl : time stepping', cdtype |
---|
281 | IF(lwp) WRITE(numout,*) '~~~~~~~~~~~' |
---|
282 | ENDIF |
---|
283 | ! |
---|
284 | IF( cdtype == 'TRA' ) THEN |
---|
285 | ll_traqsr = ln_traqsr ! active tracers case and solar penetration |
---|
286 | ll_rnf = ln_rnf ! active tracers case and river runoffs |
---|
287 | ll_isf = ln_isf ! active tracers case and ice shelf melting |
---|
288 | ELSE ! passive tracers case |
---|
289 | ll_traqsr = .FALSE. ! NO solar penetration |
---|
290 | ll_rnf = .FALSE. ! NO river runoffs ???? !!gm BUG ? |
---|
291 | ll_isf = .FALSE. ! NO ice shelf melting/freezing !!gm BUG ?? |
---|
292 | ENDIF |
---|
293 | ! |
---|
294 | IF( ( l_trdtra .AND. cdtype == 'TRA' ) .OR. ( l_trdtrc .AND. cdtype == 'TRC' ) ) THEN |
---|
295 | ALLOCATE( ztrd_atf(jpi,jpj,jpk,kjpt) ) |
---|
296 | ztrd_atf(:,:,:,:) = 0.0_wp |
---|
297 | ENDIF |
---|
298 | zfact = 1._wp / p2dt |
---|
299 | zfact1 = atfp * p2dt |
---|
300 | zfact2 = zfact1 * r1_rau0 |
---|
301 | DO jn = 1, kjpt |
---|
302 | DO jk = 1, jpkm1 |
---|
303 | DO jj = 2, jpjm1 |
---|
304 | DO ji = fs_2, fs_jpim1 |
---|
305 | ze3t_b = e3t(ji,jj,jk,Kbb) |
---|
306 | ze3t_n = e3t(ji,jj,jk,Kmm) |
---|
307 | ze3t_a = e3t(ji,jj,jk,Krhs) |
---|
308 | ! ! tracer content at Before, now and after |
---|
309 | ztc_b = ptb(ji,jj,jk,jn) * ze3t_b |
---|
310 | ztc_n = ptn(ji,jj,jk,jn) * ze3t_n |
---|
311 | ztc_a = pta(ji,jj,jk,jn) * ze3t_a |
---|
312 | ! |
---|
313 | ze3t_d = ze3t_a - 2. * ze3t_n + ze3t_b |
---|
314 | ztc_d = ztc_a - 2. * ztc_n + ztc_b |
---|
315 | ! |
---|
316 | ze3t_f = ze3t_n + atfp * ze3t_d |
---|
317 | ztc_f = ztc_n + atfp * ztc_d |
---|
318 | ! |
---|
319 | IF( jk == mikt(ji,jj) ) THEN ! first level |
---|
320 | ze3t_f = ze3t_f - zfact2 * ( (emp_b(ji,jj) - emp(ji,jj) ) & |
---|
321 | & + (fwfisf_b(ji,jj) - fwfisf(ji,jj)) ) |
---|
322 | ztc_f = ztc_f - zfact1 * ( psbc_tc(ji,jj,jn) - psbc_tc_b(ji,jj,jn) ) |
---|
323 | ENDIF |
---|
324 | IF( ln_rnf_depth ) THEN |
---|
325 | ! Rivers are not just at the surface must go down to nk_rnf(ji,jj) |
---|
326 | IF( mikt(ji,jj) <=jk .and. jk <= nk_rnf(ji,jj) ) THEN |
---|
327 | ze3t_f = ze3t_f - zfact2 * ( - (rnf_b(ji,jj) - rnf(ji,jj) ) ) & |
---|
328 | & * ( e3t(ji,jj,jk,Kmm) / h_rnf(ji,jj) ) |
---|
329 | ENDIF |
---|
330 | ELSE |
---|
331 | IF( jk == mikt(ji,jj) ) THEN ! first level |
---|
332 | ze3t_f = ze3t_f - zfact2 * ( - (rnf_b(ji,jj) - rnf(ji,jj) ) ) |
---|
333 | ENDIF |
---|
334 | ENDIF |
---|
335 | |
---|
336 | ! |
---|
337 | ! solar penetration (temperature only) |
---|
338 | IF( ll_traqsr .AND. jn == jp_tem .AND. jk <= nksr ) & |
---|
339 | & ztc_f = ztc_f - zfact1 * ( qsr_hc(ji,jj,jk) - qsr_hc_b(ji,jj,jk) ) |
---|
340 | ! |
---|
341 | ! river runoff |
---|
342 | IF( ll_rnf .AND. jk <= nk_rnf(ji,jj) ) & |
---|
343 | & ztc_f = ztc_f - zfact1 * ( rnf_tsc(ji,jj,jn) - rnf_tsc_b(ji,jj,jn) ) & |
---|
344 | & * e3t(ji,jj,jk,Kmm) / h_rnf(ji,jj) |
---|
345 | ! |
---|
346 | ! ice shelf |
---|
347 | IF( ll_isf ) THEN |
---|
348 | ! level fully include in the Losch_2008 ice shelf boundary layer |
---|
349 | IF ( jk >= misfkt(ji,jj) .AND. jk < misfkb(ji,jj) ) & |
---|
350 | ztc_f = ztc_f - zfact1 * ( risf_tsc(ji,jj,jn) - risf_tsc_b(ji,jj,jn) ) & |
---|
351 | & * e3t(ji,jj,jk,Kmm) * r1_hisf_tbl (ji,jj) |
---|
352 | ! level partially include in Losch_2008 ice shelf boundary layer |
---|
353 | IF ( jk == misfkb(ji,jj) ) & |
---|
354 | ztc_f = ztc_f - zfact1 * ( risf_tsc(ji,jj,jn) - risf_tsc_b(ji,jj,jn) ) & |
---|
355 | & * e3t(ji,jj,jk,Kmm) * r1_hisf_tbl (ji,jj) * ralpha(ji,jj) |
---|
356 | END IF |
---|
357 | ! |
---|
358 | ze3t_f = 1.e0 / ze3t_f |
---|
359 | ptb(ji,jj,jk,jn) = ztc_f * ze3t_f ! ptb <-- ptn filtered |
---|
360 | ptn(ji,jj,jk,jn) = pta(ji,jj,jk,jn) ! ptn <-- pta |
---|
361 | ! |
---|
362 | IF( ( l_trdtra .and. cdtype == 'TRA' ) .OR. ( l_trdtrc .and. cdtype == 'TRC' ) ) THEN |
---|
363 | ztrd_atf(ji,jj,jk,jn) = (ztc_f - ztc_n) * zfact/ze3t_n |
---|
364 | ENDIF |
---|
365 | ! |
---|
366 | END DO |
---|
367 | END DO |
---|
368 | END DO |
---|
369 | ! |
---|
370 | END DO |
---|
371 | ! |
---|
372 | IF( ( l_trdtra .AND. cdtype == 'TRA' ) .OR. ( l_trdtrc .AND. cdtype == 'TRC' ) ) THEN |
---|
373 | IF( l_trdtra .AND. cdtype == 'TRA' ) THEN |
---|
374 | CALL trd_tra( kt, Kmm, Krhs, cdtype, jp_tem, jptra_atf, ztrd_atf(:,:,:,jp_tem) ) |
---|
375 | CALL trd_tra( kt, Kmm, Krhs, cdtype, jp_sal, jptra_atf, ztrd_atf(:,:,:,jp_sal) ) |
---|
376 | ENDIF |
---|
377 | IF( l_trdtrc .AND. cdtype == 'TRC' ) THEN |
---|
378 | DO jn = 1, kjpt |
---|
379 | CALL trd_tra( kt, Kmm, Krhs, cdtype, jn, jptra_atf, ztrd_atf(:,:,:,jn) ) |
---|
380 | END DO |
---|
381 | ENDIF |
---|
382 | DEALLOCATE( ztrd_atf ) |
---|
383 | ENDIF |
---|
384 | ! |
---|
385 | END SUBROUTINE tra_nxt_vvl |
---|
386 | |
---|
387 | !!====================================================================== |
---|
388 | END MODULE tranxt |
---|