1 | MODULE dtadyn |
---|
2 | !!====================================================================== |
---|
3 | !! *** MODULE dtadyn *** |
---|
4 | !! OFFLINE : interpolation of the physical fields |
---|
5 | !!===================================================================== |
---|
6 | |
---|
7 | !!---------------------------------------------------------------------- |
---|
8 | !! dta_dyn_init : initialization, namelist read, and parameters control |
---|
9 | !! dta_dyn : Interpolation of the fields |
---|
10 | !!---------------------------------------------------------------------- |
---|
11 | !! * Modules used |
---|
12 | USE oce ! ocean dynamics and tracers variables |
---|
13 | USE dom_oce ! ocean space and time domain variables |
---|
14 | USE zdf_oce ! ocean vertical physics |
---|
15 | USE in_out_manager ! I/O manager |
---|
16 | USE phycst ! physical constants |
---|
17 | USE sbc_oce |
---|
18 | USE ldfslp |
---|
19 | USE ldfeiv ! eddy induced velocity coef. (ldf_eiv routine) |
---|
20 | USE ldftra_oce ! ocean tracer lateral physics |
---|
21 | USE zdfmxl |
---|
22 | USE trabbl ! tracers: bottom boundary layer |
---|
23 | USE zdfddm ! vertical physics: double diffusion |
---|
24 | USE lbclnk ! ocean lateral boundary conditions (or mpp link) |
---|
25 | USE lib_mpp ! distributed memory computing library |
---|
26 | |
---|
27 | IMPLICIT NONE |
---|
28 | PRIVATE |
---|
29 | |
---|
30 | !! * Routine accessibility |
---|
31 | PUBLIC dta_dyn_init ! called by opa.F90 |
---|
32 | PUBLIC dta_dyn ! called by step.F90 |
---|
33 | |
---|
34 | !! * Module variables |
---|
35 | INTEGER , PUBLIC, PARAMETER :: jpflx = 7 |
---|
36 | INTEGER , PUBLIC, PARAMETER :: & |
---|
37 | jptaux = 1 , & ! indice in flux for taux |
---|
38 | jptauy = 2 , & ! indice in flux for tauy |
---|
39 | jpwind = 3 , & ! indice in flux for wind speed |
---|
40 | jpemp = 4 , & ! indice in flux for E-P |
---|
41 | jpice = 5 , & ! indice in flux for ice concentration |
---|
42 | jpqsr = 6 ! indice in flux for shortwave heat flux |
---|
43 | |
---|
44 | LOGICAL , PUBLIC :: & |
---|
45 | lperdyn = .TRUE. , & ! boolean for periodic fields or not |
---|
46 | lfirdyn = .TRUE. ! boolean for the first call or not |
---|
47 | |
---|
48 | INTEGER , PUBLIC :: & |
---|
49 | ndtadyn = 12 , & ! Number of dat in one year |
---|
50 | ndtatot = 12 , & ! Number of data in the input field |
---|
51 | nsptint = 1 , & ! type of spatial interpolation |
---|
52 | nficdyn = 2 ! number of dynamical fields |
---|
53 | |
---|
54 | INTEGER :: & |
---|
55 | ndyn1, ndyn2 , & |
---|
56 | nlecoff = 0 , & ! switch for the first read |
---|
57 | numfl_t, numfl_u, & |
---|
58 | numfl_v, numfl_w |
---|
59 | |
---|
60 | |
---|
61 | REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: & |
---|
62 | tdta , & ! temperature at two consecutive times |
---|
63 | sdta , & ! salinity at two consecutive times |
---|
64 | udta , & ! zonal velocity at two consecutive times |
---|
65 | vdta , & ! meridional velocity at two consecutive times |
---|
66 | wdta , & ! vertical velocity at two consecutive times |
---|
67 | hdivdta, & ! horizontal divergence |
---|
68 | avtdta ! vertical diffusivity coefficient |
---|
69 | |
---|
70 | |
---|
71 | #if defined key_ldfslp |
---|
72 | REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: & |
---|
73 | uslpdta , & ! zonal isopycnal slopes |
---|
74 | vslpdta , & ! meridional isopycnal slopes |
---|
75 | wslpidta , & ! zonal diapycnal slopes |
---|
76 | wslpjdta ! meridional diapycnal slopes |
---|
77 | #endif |
---|
78 | |
---|
79 | #if ! defined key_off_degrad |
---|
80 | |
---|
81 | # if defined key_traldf_c2d |
---|
82 | REAL(wp), DIMENSION(jpi,jpj,2) :: & |
---|
83 | ahtwdta ! Lateral diffusivity |
---|
84 | # if defined key_trcldf_eiv |
---|
85 | REAL(wp), DIMENSION(jpi,jpj,2) :: & |
---|
86 | aeiwdta ! G&M coefficient |
---|
87 | # endif |
---|
88 | # endif |
---|
89 | |
---|
90 | #else |
---|
91 | |
---|
92 | REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: & |
---|
93 | ahtudta, ahtvdta, ahtwdta ! Lateral diffusivity |
---|
94 | # if defined key_trcldf_eiv |
---|
95 | REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: & |
---|
96 | aeiudta, aeivdta, aeiwdta ! G&M coefficient |
---|
97 | # endif |
---|
98 | |
---|
99 | #endif |
---|
100 | # if defined key_diaeiv |
---|
101 | !! GM Velocity : to be used latter |
---|
102 | REAL(wp), DIMENSION(jpi,jpj,jpk,2) :: & |
---|
103 | eivudta, eivvdta, eivwdta |
---|
104 | # endif |
---|
105 | |
---|
106 | REAL(wp), DIMENSION(jpi,jpj,jpflx,2) :: & |
---|
107 | flxdta ! auxiliary 2-D forcing fields at two consecutive times |
---|
108 | REAL(wp), DIMENSION(jpi,jpj,2) :: & |
---|
109 | zmxldta ! mixed layer depth at two consecutive times |
---|
110 | |
---|
111 | #if defined key_trcbbl_dif || defined key_trcbbl_adv |
---|
112 | REAL(wp), DIMENSION(jpi,jpj,2) :: & |
---|
113 | bblxdta , & ! frequency of bbl in the x direction at 2 consecutive times |
---|
114 | bblydta ! frequency of bbl in the y direction at 2 consecutive times |
---|
115 | #endif |
---|
116 | |
---|
117 | !! * Substitutions |
---|
118 | # include "domzgr_substitute.h90" |
---|
119 | # include "vectopt_loop_substitute.h90" |
---|
120 | !!---------------------------------------------------------------------- |
---|
121 | !! OPA 9.0 , LOCEAN-IPSL (2005) |
---|
122 | !! $Id$ |
---|
123 | !! This software is governed by the CeCILL licence see modipsl/doc/NEMO_CeCILL.txt |
---|
124 | !!---------------------------------------------------------------------- |
---|
125 | |
---|
126 | CONTAINS |
---|
127 | |
---|
128 | SUBROUTINE dta_dyn_init |
---|
129 | !!---------------------------------------------------------------------- |
---|
130 | !! *** ROUTINE dta_dyn_init *** |
---|
131 | !! |
---|
132 | !! ** Purpose : initializations of parameters for the interpolation |
---|
133 | !! |
---|
134 | !! ** Method : |
---|
135 | !! |
---|
136 | !! History : |
---|
137 | !! ! original : 92-01 (M. Imbard: sub domain) |
---|
138 | !! ! 98-04 (L.Bopp MA Foujols: slopes for isopyc.) |
---|
139 | !! ! 98-05 (L. Bopp read output of coupled run) |
---|
140 | !! ! 05-03 (O. Aumont and A. El Moussaoui) F90 |
---|
141 | !!---------------------------------------------------------------------- |
---|
142 | !! * Modules used |
---|
143 | |
---|
144 | !! * Local declarations |
---|
145 | |
---|
146 | |
---|
147 | NAMELIST/namdyn/ ndtadyn, ndtatot, nsptint, & |
---|
148 | & nficdyn, lperdyn |
---|
149 | !!---------------------------------------------------------------------- |
---|
150 | |
---|
151 | ! Define the dynamical input parameters |
---|
152 | ! ====================================== |
---|
153 | |
---|
154 | ! Read Namelist namdyn : Lateral physics on tracers |
---|
155 | REWIND( numnam ) |
---|
156 | READ ( numnam, namdyn ) |
---|
157 | |
---|
158 | IF(lwp) THEN |
---|
159 | WRITE(numout,*) |
---|
160 | WRITE(numout,*) 'namdyn : offline dynamical selection' |
---|
161 | WRITE(numout,*) '~~~~~~~' |
---|
162 | WRITE(numout,*) ' Namelist namdyn : set parameters for the lecture of the dynamical fields' |
---|
163 | WRITE(numout,*) |
---|
164 | WRITE(numout,*) ' number of elements in the FILE for a year ndtadyn = ' , ndtadyn |
---|
165 | WRITE(numout,*) ' total number of elements in the FILE ndtatot = ' , ndtatot |
---|
166 | WRITE(numout,*) ' type of interpolation nsptint = ' , nsptint |
---|
167 | WRITE(numout,*) ' number of dynamics FILE nficdyn = ' , nficdyn |
---|
168 | WRITE(numout,*) ' loop on the same FILE lperdyn = ' , lperdyn |
---|
169 | WRITE(numout,*) ' ' |
---|
170 | ENDIF |
---|
171 | |
---|
172 | END SUBROUTINE dta_dyn_init |
---|
173 | |
---|
174 | SUBROUTINE dta_dyn(kt) |
---|
175 | !!---------------------------------------------------------------------- |
---|
176 | !! *** ROUTINE dta_dyn *** |
---|
177 | !! |
---|
178 | !! ** Purpose : Prepares dynamics and physics fields from an |
---|
179 | !! OPA9 simulation for an off-line simulation |
---|
180 | !! for passive tracer |
---|
181 | !! |
---|
182 | !! ** Method : calculates the position of DATA to read READ DATA |
---|
183 | !! (example month changement) computes slopes IF needed |
---|
184 | !! interpolates DATA IF needed |
---|
185 | !! |
---|
186 | !! ** History : |
---|
187 | !! ! original : 92-01 (M. Imbard: sub domain) |
---|
188 | !! ! addition : 98-04 (L.Bopp MA Foujols: slopes for isopyc.) |
---|
189 | !! ! addition : 98-05 (L. Bopp read output of coupled run) |
---|
190 | !! ! addition : 05-03 (O. Aumont and A. El Moussaoui) F90 |
---|
191 | !! ! addition : 05-12 (C. Ethe) Adapted for DEGINT |
---|
192 | !!---------------------------------------------------------------------- |
---|
193 | !! * Modules used |
---|
194 | USE eosbn2 |
---|
195 | |
---|
196 | !! * Arguments |
---|
197 | INTEGER, INTENT( in ) :: kt ! ocean time-step index |
---|
198 | |
---|
199 | !! * Local declarations |
---|
200 | INTEGER :: iper, iperm1, iswap |
---|
201 | |
---|
202 | REAL(wp) :: zpdtan, zpdtpe, zdemi, zt |
---|
203 | REAL(wp) :: zweigh, zweighm1 |
---|
204 | |
---|
205 | REAL(wp), DIMENSION(jpi,jpj,jpflx) :: & |
---|
206 | flx ! auxiliary field for 2-D surface boundary conditions |
---|
207 | |
---|
208 | |
---|
209 | ! 0. Initialization |
---|
210 | ! ----------------- |
---|
211 | |
---|
212 | IF (lfirdyn) THEN |
---|
213 | ! |
---|
214 | ! time step MUST BE nint000 |
---|
215 | ! |
---|
216 | IF (kt.ne.nit000) THEN |
---|
217 | IF (lwp) THEN |
---|
218 | WRITE (numout,*) ' kt MUST BE EQUAL to nit000. kt=',kt & |
---|
219 | ,' nit000=',nit000 |
---|
220 | END IF |
---|
221 | STOP 'dtadyn' |
---|
222 | END if |
---|
223 | ! Initialize the parameters of the interpolation |
---|
224 | CALL dta_dyn_init |
---|
225 | ENDIF |
---|
226 | |
---|
227 | |
---|
228 | zpdtan = raass / rdt |
---|
229 | zpdtpe = ((zpdtan / FLOAT (ndtadyn))) |
---|
230 | zdemi = zpdtpe * 0.5 |
---|
231 | zt = (FLOAT (kt) + zdemi ) / zpdtpe |
---|
232 | |
---|
233 | zweigh = zt - FLOAT(INT(zt)) |
---|
234 | zweighm1 = 1. - zweigh |
---|
235 | |
---|
236 | IF (lperdyn) THEN |
---|
237 | iperm1 = MOD(INT(zt),ndtadyn) |
---|
238 | ELSE |
---|
239 | iperm1 = MOD(INT(zt),(ndtatot-1)) + 1 |
---|
240 | ENDIF |
---|
241 | iper = iperm1 + 1 |
---|
242 | IF (iperm1 == 0) THEN |
---|
243 | IF (lperdyn) THEN |
---|
244 | iperm1 = ndtadyn |
---|
245 | ELSE |
---|
246 | IF (lfirdyn) THEN |
---|
247 | IF (lwp) THEN |
---|
248 | WRITE (numout,*) ' dynamic file is not periodic ' |
---|
249 | WRITE (numout,*) ' with or without interpolation, ' |
---|
250 | WRITE (numout,*) ' we take the first value' |
---|
251 | WRITE (numout,*) ' for the previous period ' |
---|
252 | WRITE (numout,*) ' iperm1 = 0 ' |
---|
253 | END IF |
---|
254 | END IF |
---|
255 | END IF |
---|
256 | END IF |
---|
257 | |
---|
258 | iswap = 0 |
---|
259 | |
---|
260 | ! 1. First call lfirdyn = true |
---|
261 | ! ---------------------------- |
---|
262 | |
---|
263 | IF (lfirdyn) THEN |
---|
264 | ! |
---|
265 | ! store the information of the period read |
---|
266 | ! |
---|
267 | ndyn1 = iperm1 |
---|
268 | ndyn2 = iper |
---|
269 | |
---|
270 | IF (lwp) THEN |
---|
271 | WRITE (numout,*) & |
---|
272 | ' dynamics DATA READ for the period ndyn1 =',ndyn1, & |
---|
273 | & ' and for the period ndyn2 = ',ndyn2 |
---|
274 | WRITE (numout,*) ' time step is :',kt |
---|
275 | WRITE (numout,*) ' we have ndtadyn = ',ndtadyn,& |
---|
276 | & ' records in the dynamic FILE for one year' |
---|
277 | END IF |
---|
278 | ! |
---|
279 | ! DATA READ for the iperm1 period |
---|
280 | ! |
---|
281 | IF( iperm1 /= 0 ) THEN |
---|
282 | CALL dynrea( kt, iperm1 ) |
---|
283 | ELSE |
---|
284 | CALL dynrea( kt, 1 ) |
---|
285 | ENDIF |
---|
286 | ! |
---|
287 | ! Computes dynamical fields |
---|
288 | ! |
---|
289 | tn(:,:,:)=tdta(:,:,:,2) |
---|
290 | sn(:,:,:)=sdta(:,:,:,2) |
---|
291 | avt(:,:,:)=avtdta(:,:,:,2) |
---|
292 | |
---|
293 | |
---|
294 | IF(lwp) THEN |
---|
295 | WRITE(numout,*)' temperature ' |
---|
296 | WRITE(numout,*) |
---|
297 | CALL prihre(tn(1,1,1),jpi,jpj,1,jpi,20,1,jpj,20,1.,numout) |
---|
298 | WRITE(numout,*) ' level = ',jpk/2 |
---|
299 | CALL prihre(tn(1,1,jpk/2),jpi,jpj,1,jpi,20,1,jpj,20,1.,numout) |
---|
300 | WRITE(numout,*) ' level = ',jpkm1 |
---|
301 | CALL prihre(tn(1,1,jpkm1),jpi,jpj,1,jpi,20,1,jpj,20,1.,numout) |
---|
302 | ENDIF |
---|
303 | |
---|
304 | #if defined key_ldfslp |
---|
305 | CALL eos( tn, sn, rhd, rhop ) ! Time-filtered in situ density |
---|
306 | CALL bn2( tn, sn, rn2 ) ! before Brunt-Vaisala frequency |
---|
307 | CALL zdf_mxl( kt ) ! mixed layer depth |
---|
308 | CALL ldf_slp( kt, rhd, rn2 ) |
---|
309 | |
---|
310 | uslpdta(:,:,:,2)=uslp(:,:,:) |
---|
311 | vslpdta(:,:,:,2)=vslp(:,:,:) |
---|
312 | wslpidta(:,:,:,2)=wslpi(:,:,:) |
---|
313 | wslpjdta(:,:,:,2)=wslpj(:,:,:) |
---|
314 | #endif |
---|
315 | ! |
---|
316 | ! swap from record 2 to 1 |
---|
317 | ! |
---|
318 | udta(:,:,:,1)=udta(:,:,:,2) |
---|
319 | vdta(:,:,:,1)=vdta(:,:,:,2) |
---|
320 | wdta(:,:,:,1)=wdta(:,:,:,2) |
---|
321 | #if defined key_trc_diatrd |
---|
322 | hdivdta(:,:,:,1)=hdivdta(:,:,:,2) |
---|
323 | #endif |
---|
324 | avtdta(:,:,:,1)=avtdta(:,:,:,2) |
---|
325 | tdta(:,:,:,1)=tdta(:,:,:,2) |
---|
326 | sdta(:,:,:,1)=sdta(:,:,:,2) |
---|
327 | #if defined key_ldfslp |
---|
328 | uslpdta(:,:,:,1)=uslpdta(:,:,:,2) |
---|
329 | vslpdta(:,:,:,1)=vslpdta(:,:,:,2) |
---|
330 | wslpidta(:,:,:,1)=wslpidta(:,:,:,2) |
---|
331 | wslpjdta(:,:,:,1)=wslpjdta(:,:,:,2) |
---|
332 | #endif |
---|
333 | flxdta(:,:,:,1) = flxdta(:,:,:,2) |
---|
334 | zmxldta(:,:,1)=zmxldta(:,:,2) |
---|
335 | #if ! defined key_off_degrad |
---|
336 | |
---|
337 | # if defined key_traldf_c2d |
---|
338 | ahtwdta(:,:,1)= ahtwdta(:,:,2) |
---|
339 | # if defined key_trcldf_eiv |
---|
340 | aeiwdta(:,:,1)= aeiwdta(:,:,2) |
---|
341 | # endif |
---|
342 | # endif |
---|
343 | |
---|
344 | #else |
---|
345 | ahtudta(:,:,:,1) = ahtudta(:,:,:,2) |
---|
346 | ahtvdta(:,:,:,1) = ahtvdta(:,:,:,2) |
---|
347 | ahtwdta(:,:,:,1) = ahtwdta(:,:,:,2) |
---|
348 | # if defined key_trcldf_eiv |
---|
349 | aeiudta(:,:,:,1) = aeiudta(:,:,:,2) |
---|
350 | aeivdta(:,:,:,1) = aeivdta(:,:,:,2) |
---|
351 | aeiwdta(:,:,:,1) = aeiwdta(:,:,:,2) |
---|
352 | # endif |
---|
353 | |
---|
354 | #endif |
---|
355 | |
---|
356 | #if defined key_trcbbl_dif || defined key_trcbbl_adv |
---|
357 | bblxdta(:,:,1)=bblxdta(:,:,2) |
---|
358 | bblydta(:,:,1)=bblydta(:,:,2) |
---|
359 | #endif |
---|
360 | ! |
---|
361 | ! indicates a swap |
---|
362 | ! |
---|
363 | iswap = 1 |
---|
364 | ! |
---|
365 | ! DATA READ for the iper period |
---|
366 | ! |
---|
367 | CALL dynrea( kt, iper ) |
---|
368 | ! |
---|
369 | ! Computes wdta (and slopes if key_trahdfiso) |
---|
370 | ! |
---|
371 | tn(:,:,:)=tdta(:,:,:,2) |
---|
372 | sn(:,:,:)=sdta(:,:,:,2) |
---|
373 | avt(:,:,:)=avtdta(:,:,:,2) |
---|
374 | |
---|
375 | |
---|
376 | #if defined key_ldfslp |
---|
377 | CALL eos( tn, sn, rhd, rhop ) ! Time-filtered in situ density |
---|
378 | CALL bn2( tn, sn, rn2 ) ! before Brunt-Vaisala frequency |
---|
379 | CALL zdf_mxl( kt ) ! mixed layer depth |
---|
380 | CALL ldf_slp( kt, rhd, rn2 ) |
---|
381 | |
---|
382 | uslpdta(:,:,:,2)=uslp(:,:,:) |
---|
383 | vslpdta(:,:,:,2)=vslp(:,:,:) |
---|
384 | wslpidta(:,:,:,2)=wslpi(:,:,:) |
---|
385 | wslpjdta(:,:,:,2)=wslpj(:,:,:) |
---|
386 | #endif |
---|
387 | ! |
---|
388 | ! trace the first CALL |
---|
389 | ! |
---|
390 | lfirdyn=.FALSE. |
---|
391 | ENDIF |
---|
392 | ! |
---|
393 | ! and now what we have to DO at every time step |
---|
394 | ! |
---|
395 | ! check the validity of the period in memory |
---|
396 | ! |
---|
397 | IF (iperm1.NE.ndyn1) THEN |
---|
398 | IF (iperm1.EQ.0.) THEN |
---|
399 | IF (lwp) THEN |
---|
400 | WRITE (numout,*) ' dynamic file is not periodic ' |
---|
401 | WRITE (numout,*) ' with or without interpolation, ' |
---|
402 | WRITE (numout,*) ' we take the last value' |
---|
403 | WRITE (numout,*) ' for the last period ' |
---|
404 | WRITE (numout,*) ' iperm1 = 12 ' |
---|
405 | WRITE (numout,*) ' iper = 13' |
---|
406 | ENDIF |
---|
407 | iperm1 = 12 |
---|
408 | iper =13 |
---|
409 | ENDIF |
---|
410 | ! |
---|
411 | ! we have to prepare a NEW READ of DATA |
---|
412 | ! |
---|
413 | ! swap from record 2 to 1 |
---|
414 | ! |
---|
415 | udta(:,:,:,1) = udta(:,:,:,2) |
---|
416 | vdta(:,:,:,1) = vdta(:,:,:,2) |
---|
417 | wdta(:,:,:,1) = wdta(:,:,:,2) |
---|
418 | #if defined key_trc_diatrd |
---|
419 | hdivdta(:,:,:,1)=hdivdta(:,:,:,2) |
---|
420 | #endif |
---|
421 | avtdta(:,:,:,1) = avtdta(:,:,:,2) |
---|
422 | tdta(:,:,:,1) = tdta(:,:,:,2) |
---|
423 | sdta(:,:,:,1) = sdta(:,:,:,2) |
---|
424 | #if defined key_ldfslp |
---|
425 | uslpdta(:,:,:,1) = uslpdta(:,:,:,2) |
---|
426 | vslpdta(:,:,:,1) = vslpdta(:,:,:,2) |
---|
427 | wslpidta(:,:,:,1) = wslpidta(:,:,:,2) |
---|
428 | wslpjdta(:,:,:,1) = wslpjdta(:,:,:,2) |
---|
429 | #endif |
---|
430 | flxdta(:,:,:,1) = flxdta(:,:,:,2) |
---|
431 | zmxldta(:,:,1) = zmxldta(:,:,2) |
---|
432 | |
---|
433 | #if ! defined key_off_degrad |
---|
434 | |
---|
435 | # if defined key_traldf_c2d |
---|
436 | ahtwdta(:,:,1)= ahtwdta(:,:,2) |
---|
437 | # if defined key_trcldf_eiv |
---|
438 | aeiwdta(:,:,1)= aeiwdta(:,:,2) |
---|
439 | # endif |
---|
440 | # endif |
---|
441 | |
---|
442 | #else |
---|
443 | ahtudta(:,:,:,1) = ahtudta(:,:,:,2) |
---|
444 | ahtvdta(:,:,:,1) = ahtvdta(:,:,:,2) |
---|
445 | ahtwdta(:,:,:,1) = ahtwdta(:,:,:,2) |
---|
446 | # if defined key_trcldf_eiv |
---|
447 | aeiudta(:,:,:,1) = aeiudta(:,:,:,2) |
---|
448 | aeivdta(:,:,:,1) = aeivdta(:,:,:,2) |
---|
449 | aeiwdta(:,:,:,1) = aeiwdta(:,:,:,2) |
---|
450 | # endif |
---|
451 | |
---|
452 | #endif |
---|
453 | |
---|
454 | #if defined key_trcbbl_dif || defined key_trcbbl_adv |
---|
455 | bblxdta(:,:,1) = bblxdta(:,:,2) |
---|
456 | bblydta(:,:,1) = bblydta(:,:,2) |
---|
457 | #endif |
---|
458 | ! |
---|
459 | ! indicates a swap |
---|
460 | ! |
---|
461 | iswap = 1 |
---|
462 | ! |
---|
463 | ! READ DATA for the iper period |
---|
464 | ! |
---|
465 | CALL dynrea( kt, iper ) |
---|
466 | ! |
---|
467 | ! Computes wdta (and slopes if key_trahdfiso) |
---|
468 | ! |
---|
469 | tn(:,:,:)=tdta(:,:,:,2) |
---|
470 | sn(:,:,:)=sdta(:,:,:,2) |
---|
471 | avt(:,:,:)=avtdta(:,:,:,2) |
---|
472 | |
---|
473 | #if defined key_ldfslp |
---|
474 | CALL eos( tn, sn, rhd, rhop ) ! Time-filtered in situ density |
---|
475 | CALL bn2( tn, sn, rn2 ) ! before Brunt-Vaisala frequency |
---|
476 | CALL zdf_mxl( kt ) ! mixed layer depth |
---|
477 | CALL ldf_slp( kt, rhd, rn2 ) |
---|
478 | |
---|
479 | uslpdta(:,:,:,2)=uslp(:,:,:) |
---|
480 | vslpdta(:,:,:,2)=vslp(:,:,:) |
---|
481 | wslpidta(:,:,:,2)=wslpi(:,:,:) |
---|
482 | wslpjdta(:,:,:,2)=wslpj(:,:,:) |
---|
483 | #endif |
---|
484 | ! |
---|
485 | ! store the information of the period read |
---|
486 | ! |
---|
487 | ndyn1 = ndyn2 |
---|
488 | ndyn2 = iper |
---|
489 | ! |
---|
490 | ! we have READ another period of DATA ! |
---|
491 | IF (lwp) THEN |
---|
492 | WRITE (numout,*) ' dynamics DATA READ for the period ndyn1 =',ndyn1 |
---|
493 | WRITE (numout,*) ' and the period ndyn2 = ',ndyn2 |
---|
494 | WRITE (numout,*) ' time step is :',kt |
---|
495 | END IF |
---|
496 | |
---|
497 | END IF |
---|
498 | |
---|
499 | ! |
---|
500 | ! compute the DATA at the given time step |
---|
501 | ! |
---|
502 | IF ( nsptint == 0 ) THEN |
---|
503 | ! |
---|
504 | ! no spatial interpolation |
---|
505 | ! |
---|
506 | ! DATA are probably correct |
---|
507 | ! we have to initialize DATA IF we have changed the period |
---|
508 | ! |
---|
509 | IF (iswap.eq.1) THEN |
---|
510 | ! |
---|
511 | ! initialize now fields with the NEW DATA READ |
---|
512 | ! |
---|
513 | un(:,:,:)=udta(:,:,:,2) |
---|
514 | vn(:,:,:)=vdta(:,:,:,2) |
---|
515 | wn(:,:,:)=wdta(:,:,:,2) |
---|
516 | hdivn(:,:,:)=hdivdta(:,:,:,2) |
---|
517 | #if defined key_trc_zdfddm |
---|
518 | avs(:,:,:)=avtdta(:,:,:,2) |
---|
519 | #endif |
---|
520 | avt(:,:,:)=avtdta(:,:,:,2) |
---|
521 | tn(:,:,:)=tdta(:,:,:,2) |
---|
522 | sn(:,:,:)=sdta(:,:,:,2) |
---|
523 | #if defined key_ldfslp |
---|
524 | uslp(:,:,:)=uslpdta(:,:,:,2) |
---|
525 | vslp(:,:,:)=vslpdta(:,:,:,2) |
---|
526 | wslpi(:,:,:)=wslpidta(:,:,:,2) |
---|
527 | wslpj(:,:,:)=wslpjdta(:,:,:,2) |
---|
528 | #endif |
---|
529 | flx(:,:,:) = flxdta(:,:,:,2) |
---|
530 | hmld(:,:)=zmxldta(:,:,2) |
---|
531 | #if ! defined key_off_degrad |
---|
532 | |
---|
533 | # if defined key_traldf_c2d |
---|
534 | ahtwdta(:,:,1)= ahtwdta(:,:,2) |
---|
535 | # if defined key_trcldf_eiv |
---|
536 | aeiwdta(:,:,1)= aeiwdta(:,:,2) |
---|
537 | # endif |
---|
538 | # endif |
---|
539 | |
---|
540 | #else |
---|
541 | ahtudta(:,:,:,1) = ahtudta(:,:,:,2) |
---|
542 | ahtvdta(:,:,:,1) = ahtvdta(:,:,:,2) |
---|
543 | ahtwdta(:,:,:,1) = ahtwdta(:,:,:,2) |
---|
544 | # if defined key_trcldf_eiv |
---|
545 | aeiudta(:,:,:,1) = aeiudta(:,:,:,2) |
---|
546 | aeivdta(:,:,:,1) = aeivdta(:,:,:,2) |
---|
547 | aeiwdta(:,:,:,1) = aeiwdta(:,:,:,2) |
---|
548 | # endif |
---|
549 | |
---|
550 | #endif |
---|
551 | |
---|
552 | #if defined key_trcbbl_dif || defined key_trcbbl_adv |
---|
553 | bblx(:,:)=bblxdta(:,:,2) |
---|
554 | bbly(:,:)=bblydta(:,:,2) |
---|
555 | #endif |
---|
556 | ! |
---|
557 | ! keep needed fluxes |
---|
558 | ! |
---|
559 | wndm(:,:) = flx(:,:,jpwind) |
---|
560 | fr_i(:,:) = flx(:,:,jpice) |
---|
561 | emp(:,:) = flx(:,:,jpemp) |
---|
562 | emps(:,:) = emp(:,:) |
---|
563 | qsr(:,:) = flx(:,:,jpqsr) |
---|
564 | |
---|
565 | END IF |
---|
566 | |
---|
567 | ELSE |
---|
568 | IF ( nsptint == 1 ) THEN |
---|
569 | ! |
---|
570 | ! linear interpolation |
---|
571 | ! |
---|
572 | ! initialize now fields with a linear interpolation |
---|
573 | ! |
---|
574 | un(:,:,:) = zweighm1 * udta(:,:,:,1) + zweigh * udta(:,:,:,2) |
---|
575 | vn(:,:,:) = zweighm1 * vdta(:,:,:,1) + zweigh * vdta(:,:,:,2) |
---|
576 | wn(:,:,:) = zweighm1 * wdta(:,:,:,1) + zweigh * wdta(:,:,:,2) |
---|
577 | hdivn(:,:,:) = zweighm1 * hdivdta(:,:,:,1) + zweigh * hdivdta(:,:,:,2) |
---|
578 | #if defined key_zdfddm |
---|
579 | avs(:,:,:)= zweighm1 * avtdta(:,:,:,1) + zweigh * avtdta(:,:,:,2) |
---|
580 | #endif |
---|
581 | avt(:,:,:)= zweighm1 * avtdta(:,:,:,1) + zweigh * avtdta(:,:,:,2) |
---|
582 | tn(:,:,:) = zweighm1 * tdta(:,:,:,1) + zweigh * tdta(:,:,:,2) |
---|
583 | sn(:,:,:) = zweighm1 * sdta(:,:,:,1) + zweigh * sdta(:,:,:,2) |
---|
584 | |
---|
585 | |
---|
586 | #if defined key_ldfslp |
---|
587 | uslp(:,:,:) = zweighm1 * uslpdta(:,:,:,1) + zweigh * uslpdta(:,:,:,2) |
---|
588 | vslp(:,:,:) = zweighm1 * vslpdta(:,:,:,1) + zweigh * vslpdta(:,:,:,2) |
---|
589 | wslpi(:,:,:) = zweighm1 * wslpidta(:,:,:,1) + zweigh * wslpidta(:,:,:,2) |
---|
590 | wslpj(:,:,:) = zweighm1 * wslpjdta(:,:,:,1) + zweigh * wslpjdta(:,:,:,2) |
---|
591 | #endif |
---|
592 | flx(:,:,:) = zweighm1 * flxdta(:,:,:,1) + zweigh * flxdta(:,:,:,2) |
---|
593 | hmld(:,:) = zweighm1 * zmxldta(:,:,1) + zweigh * zmxldta(:,:,2) |
---|
594 | #if ! defined key_off_degrad |
---|
595 | |
---|
596 | # if defined key_traldf_c2d |
---|
597 | ahtw(:,:) = zweighm1 * ahtwdta(:,:,1) + zweigh * ahtwdta(:,:,2) |
---|
598 | # if defined key_trcldf_eiv |
---|
599 | aeiw(:,:) = zweighm1 * aeiwdta(:,:,1) + zweigh * aeiwdta(:,:,2) |
---|
600 | # endif |
---|
601 | # endif |
---|
602 | |
---|
603 | #else |
---|
604 | ahtu(:,:,:) = zweighm1 * ahtudta(:,:,:,1) + zweigh * ahtudta(:,:,:,2) |
---|
605 | ahtv(:,:,:) = zweighm1 * ahtvdta(:,:,:,1) + zweigh * ahtvdta(:,:,:,2) |
---|
606 | ahtw(:,:,:) = zweighm1 * ahtwdta(:,:,:,1) + zweigh * ahtwdta(:,:,:,2) |
---|
607 | # if defined key_trcldf_eiv |
---|
608 | aeiu(:,:,:) = zweighm1 * aeiudta(:,:,:,1) + zweigh * aeiudta(:,:,:,2) |
---|
609 | aeiv(:,:,:) = zweighm1 * aeivdta(:,:,:,1) + zweigh * aeivdta(:,:,:,2) |
---|
610 | aeiw(:,:,:) = zweighm1 * aeiwdta(:,:,:,1) + zweigh * aeiwdta(:,:,:,2) |
---|
611 | # endif |
---|
612 | |
---|
613 | #endif |
---|
614 | |
---|
615 | #if defined key_trcbbl_dif || defined key_trcbbl_adv |
---|
616 | bblx(:,:) = zweighm1 * bblxdta(:,:,1) + zweigh * bblxdta(:,:,2) |
---|
617 | bbly(:,:) = zweighm1 * bblydta(:,:,1) + zweigh * bblydta(:,:,2) |
---|
618 | #endif |
---|
619 | ! |
---|
620 | ! keep needed fluxes |
---|
621 | ! |
---|
622 | wndm(:,:) = flx(:,:,jpwind) |
---|
623 | fr_i(:,:) = flx(:,:,jpice) |
---|
624 | emp(:,:) = flx(:,:,jpemp) |
---|
625 | emps(:,:) = emp(:,:) |
---|
626 | qsr(:,:) = flx(:,:,jpqsr) |
---|
627 | ! |
---|
628 | ! other interpolation |
---|
629 | ! |
---|
630 | ELSE |
---|
631 | |
---|
632 | WRITE (numout,*) ' this kind of interpolation don t EXIST' |
---|
633 | WRITE (numout,*) ' at the moment. we STOP ' |
---|
634 | STOP 'dtadyn' |
---|
635 | |
---|
636 | END IF |
---|
637 | |
---|
638 | END IF |
---|
639 | ! |
---|
640 | ! lb in any case, we need rhop |
---|
641 | ! |
---|
642 | CALL eos( tn, sn, rhd, rhop ) |
---|
643 | |
---|
644 | #if ! defined key_off_degrad && defined key_traldf_c2d |
---|
645 | ! In case of 2D varying coefficients, we need aeiv and aeiu |
---|
646 | IF( lk_traldf_eiv ) CALL ldf_eiv( kt ) ! eddy induced velocity coefficient |
---|
647 | #endif |
---|
648 | |
---|
649 | END SUBROUTINE dta_dyn |
---|
650 | |
---|
651 | SUBROUTINE dynrea( kt, kenr ) |
---|
652 | !!---------------------------------------------------------------------- |
---|
653 | !! *** ROUTINE dynrea *** |
---|
654 | !! |
---|
655 | !! ** Purpose : READ dynamics fiels from OPA9 netcdf output |
---|
656 | !! |
---|
657 | !! ** Method : READ the kenr records of DATA and store in |
---|
658 | !! in udta(...,2), .... |
---|
659 | !! |
---|
660 | !! ** History : additions : M. Levy et M. Benjelloul jan 2001 |
---|
661 | !! (netcdf FORMAT) |
---|
662 | !! 05-03 (O. Aumont and A. El Moussaoui) F90 |
---|
663 | !! 06-07 : (C. Ethe) use of iom module |
---|
664 | !!---------------------------------------------------------------------- |
---|
665 | !! * Modules used |
---|
666 | USE iom |
---|
667 | |
---|
668 | !! * Arguments |
---|
669 | INTEGER, INTENT( in ) :: kt, kenr ! time index |
---|
670 | !! * Local declarations |
---|
671 | INTEGER :: ji, jj, jk, jkenr |
---|
672 | |
---|
673 | REAL(wp), DIMENSION(jpi,jpj,jpk) :: & |
---|
674 | zu, zv, zw, zt, zs, zavt , & ! 3-D dynamical fields |
---|
675 | zhdiv ! horizontal divergence |
---|
676 | |
---|
677 | REAL(wp), DIMENSION(jpi,jpj) :: & |
---|
678 | zemp, zqsr, zmld, zice, zwspd |
---|
679 | #if defined key_trcbbl_dif || defined key_trcbbl_adv |
---|
680 | REAL(wp), DIMENSION(jpi,jpj) :: & |
---|
681 | zbblx, zbbly |
---|
682 | #endif |
---|
683 | |
---|
684 | #if ! defined key_off_degrad |
---|
685 | |
---|
686 | # if defined key_traldf_c2d |
---|
687 | REAL(wp), DIMENSION(jpi,jpj) :: & |
---|
688 | zahtw |
---|
689 | # if defined key_trcldf_eiv |
---|
690 | REAL(wp), DIMENSION(jpi,jpj) :: & |
---|
691 | zaeiw |
---|
692 | # endif |
---|
693 | # endif |
---|
694 | |
---|
695 | #else |
---|
696 | |
---|
697 | REAL(wp), DIMENSION(jpi,jpj,jpk) :: & |
---|
698 | zahtu, zahtv, zahtw ! Lateral diffusivity |
---|
699 | # if defined key_trcldf_eiv |
---|
700 | REAL(wp), DIMENSION(jpi,jpj,jpk) :: & |
---|
701 | zaeiu, zaeiv, zaeiw ! G&M coefficient |
---|
702 | # endif |
---|
703 | |
---|
704 | #endif |
---|
705 | |
---|
706 | # if defined key_diaeiv |
---|
707 | !! GM Velocity : to be used latter |
---|
708 | REAL(wp), DIMENSION(jpi,jpj,jpk) :: & |
---|
709 | zeivu, zeivv, zeivw |
---|
710 | # endif |
---|
711 | |
---|
712 | CHARACTER(len=45) :: & |
---|
713 | clname_t = 'dyna_grid_T.nc', & |
---|
714 | clname_u = 'dyna_grid_U.nc', & |
---|
715 | clname_v = 'dyna_grid_V.nc', & |
---|
716 | clname_w = 'dyna_grid_W.nc' |
---|
717 | ! |
---|
718 | ! 0. Initialization |
---|
719 | ! cas d'un fichier non periodique : on utilise deux fois le premier et |
---|
720 | ! le dernier champ temporel |
---|
721 | |
---|
722 | jkenr = kenr |
---|
723 | |
---|
724 | IF(lwp) THEN |
---|
725 | WRITE(numout,*) |
---|
726 | WRITE(numout,*) 'Dynrea : reading dynamical fields, kenr = ', jkenr |
---|
727 | WRITE(numout,*) ' ~~~~~~~' |
---|
728 | #if defined key_off_degrad |
---|
729 | WRITE(numout,*) ' Degraded fields' |
---|
730 | #endif |
---|
731 | WRITE(numout,*) |
---|
732 | ENDIF |
---|
733 | |
---|
734 | |
---|
735 | IF( kt == nit000 .AND. nlecoff == 0 ) THEN |
---|
736 | |
---|
737 | nlecoff = 1 |
---|
738 | |
---|
739 | CALL iom_open ( clname_t, numfl_t ) |
---|
740 | CALL iom_open ( clname_u, numfl_u ) |
---|
741 | CALL iom_open ( clname_v, numfl_v ) |
---|
742 | CALL iom_open ( clname_w, numfl_w ) |
---|
743 | |
---|
744 | ENDIF |
---|
745 | |
---|
746 | ! file grid-T |
---|
747 | !--------------- |
---|
748 | CALL iom_get ( numfl_t, jpdom_data, 'votemper', zt (:,:,:), jkenr ) |
---|
749 | CALL iom_get ( numfl_t, jpdom_data, 'vosaline', zs (:,:,:), jkenr ) |
---|
750 | CALL iom_get ( numfl_t, jpdom_data, 'somixhgt', zmld (:,: ), jkenr ) |
---|
751 | CALL iom_get ( numfl_t, jpdom_data, 'sowaflcd', zemp (:,: ), jkenr ) |
---|
752 | CALL iom_get ( numfl_t, jpdom_data, 'soshfldo', zqsr (:,: ), jkenr ) |
---|
753 | CALL iom_get ( numfl_t, jpdom_data, 'soicecov', zice (:,: ), jkenr ) |
---|
754 | CALL iom_get ( numfl_t, jpdom_data, 'sowindsp', zwspd(:,: ), jkenr ) |
---|
755 | |
---|
756 | ! file grid-U |
---|
757 | !--------------- |
---|
758 | CALL iom_get ( numfl_u, jpdom_data, 'vozocrtx', zu (:,:,:), jkenr ) |
---|
759 | #if defined key_trcbbl_dif || defined key_trcbbl_adv |
---|
760 | CALL iom_get ( numfl_u, jpdom_data, 'sobblcox', zbblx(:,: ), jkenr ) |
---|
761 | #endif |
---|
762 | |
---|
763 | #if defined key_diaeiv |
---|
764 | !! GM Velocity : to be used latter |
---|
765 | CALL iom_get ( numfl_u, jpdom_data, 'vozoeivu', zeivu(:,:,:), jkenr ) |
---|
766 | #endif |
---|
767 | |
---|
768 | # if defined key_off_degrad |
---|
769 | CALL iom_get ( numfl_u, jpdom_data, 'vozoahtu', zahtu(:,:,:), jkenr ) |
---|
770 | # if defined key_trcldf_eiv |
---|
771 | CALL iom_get ( numfl_u, jpdom_data, 'vozoaeiu', zaeiu(:,:,:), jkenr ) |
---|
772 | # endif |
---|
773 | #endif |
---|
774 | |
---|
775 | ! file grid-V |
---|
776 | !--------------- |
---|
777 | CALL iom_get ( numfl_v, jpdom_data, 'vomecrty', zv (:,:,:), jkenr ) |
---|
778 | #if defined key_trcbbl_dif || defined key_trcbbl_adv |
---|
779 | CALL iom_get ( numfl_v, jpdom_data, 'sobblcoy', zbbly(:,: ), jkenr ) |
---|
780 | #endif |
---|
781 | |
---|
782 | #if defined key_diaeiv |
---|
783 | !! GM Velocity : to be used latter |
---|
784 | CALL iom_get ( numfl_v, jpdom_data, 'vomeeivv', zeivv(:,:,:), jkenr ) |
---|
785 | #endif |
---|
786 | |
---|
787 | #if defined key_off_degrad |
---|
788 | CALL iom_get ( numfl_v, jpdom_data, 'vomeahtv', zahtv(:,:,:), jkenr ) |
---|
789 | # if defined key_trcldf_eiv |
---|
790 | CALL iom_get ( numfl_v, jpdom_data, 'vomeaeiv', zaeiv(:,:,:), jkenr ) |
---|
791 | # endif |
---|
792 | #endif |
---|
793 | |
---|
794 | ! file grid-W |
---|
795 | !--------------- |
---|
796 | !! CALL iom_get ( numfl_w, jpdom_data, 'vovecrtz', zw (:,:,:), jkenr ) |
---|
797 | # if defined key_zdfddm |
---|
798 | CALL iom_get ( numfl_w, jpdom_data, 'voddmavs', zavt (:,:,:), jkenr ) |
---|
799 | #else |
---|
800 | CALL iom_get ( numfl_w, jpdom_data, 'votkeavt', zavt (:,:,:), jkenr ) |
---|
801 | #endif |
---|
802 | |
---|
803 | # if defined key_diaeiv |
---|
804 | !! GM Velocity : to be used latter |
---|
805 | CALL iom_get ( numfl_w, jpdom_data, 'voveeivw', zeivw(:,:,:), jkenr ) |
---|
806 | #endif |
---|
807 | |
---|
808 | #if ! defined key_off_degrad |
---|
809 | # if defined key_traldf_c2d |
---|
810 | CALL iom_get ( numfl_w, jpdom_data, 'soleahtw', zahtw (:,: ), jkenr ) |
---|
811 | # if defined key_trcldf_eiv |
---|
812 | CALL iom_get ( numfl_w, jpdom_data, 'soleaeiw', zaeiw (:,: ), jkenr ) |
---|
813 | # endif |
---|
814 | # endif |
---|
815 | #else |
---|
816 | !! degradation-integration |
---|
817 | CALL iom_get ( numfl_w, jpdom_data, 'voveahtw', zahtw(:,:,:), jkenr ) |
---|
818 | # if defined key_trcldf_eiv |
---|
819 | CALL iom_get ( numfl_w, jpdom_data, 'voveaeiw', zaeiw(:,:,:), jkenr ) |
---|
820 | # endif |
---|
821 | #endif |
---|
822 | |
---|
823 | udta(:,:,:,2) = zu(:,:,:) * umask(:,:,:) |
---|
824 | vdta(:,:,:,2) = zv(:,:,:) * vmask(:,:,:) |
---|
825 | !! wdta(:,:,:,2) = zw(:,:,:) * tmask(:,:,:) |
---|
826 | |
---|
827 | |
---|
828 | ! Computation of vertical velocity using horizontal divergence |
---|
829 | zhdiv(:,:,:) = 0. |
---|
830 | DO jk = 1, jpkm1 |
---|
831 | DO jj = 2, jpjm1 |
---|
832 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
833 | #if defined key_zco |
---|
834 | zhdiv(ji,jj,jk) = ( e2u(ji,jj) * udta(ji,jj,jk,2) - e2u(ji-1,jj ) * udta(ji-1,jj ,jk,2) & |
---|
835 | & + e1v(ji,jj) * vdta(ji,jj,jk,2) - e1v(ji ,jj-1) * vdta(ji ,jj-1,jk,2) ) & |
---|
836 | & / ( e1t(ji,jj) * e2t(ji,jj) ) |
---|
837 | #else |
---|
838 | zhdiv(ji,jj,jk) = & |
---|
839 | ( e2u(ji,jj)*fse3u(ji,jj,jk)*udta(ji,jj,jk,2) - e2u(ji-1,jj)*fse3u(ji-1,jj,jk)*udta(ji-1,jj,jk,2) & |
---|
840 | + e1v(ji,jj)*fse3v(ji,jj,jk)*vdta(ji,jj,jk,2) - e1v(ji,jj-1)*fse3v(ji,jj-1,jk)*vdta(ji,jj-1,jk,2) ) & |
---|
841 | / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) ) |
---|
842 | #endif |
---|
843 | END DO |
---|
844 | END DO |
---|
845 | ENDDO |
---|
846 | |
---|
847 | ! Lateral boundary conditions on hdiv |
---|
848 | CALL lbc_lnk( zhdiv, 'T', 1. ) |
---|
849 | |
---|
850 | |
---|
851 | ! computation of vertical velocity from the bottom |
---|
852 | zw(:,:,jpk) = 0. |
---|
853 | DO jk = jpkm1, 1, -1 |
---|
854 | zw(:,:,jk) = zw(:,:,jk+1) - fse3t(:,:,jk) * zhdiv(:,:,jk) |
---|
855 | END DO |
---|
856 | |
---|
857 | wdta(:,:,:,2) = zw(:,:,:) * tmask(:,:,:) |
---|
858 | hdivdta(:,:,:,2) = zhdiv(:,:,:) * tmask(:,:,:) |
---|
859 | |
---|
860 | tdta(:,:,:,2) = zt(:,:,:) * tmask(:,:,:) |
---|
861 | sdta(:,:,:,2) = zs(:,:,:) * tmask(:,:,:) |
---|
862 | avtdta(:,:,:,2) = zavt(:,:,:) * tmask(:,:,:) |
---|
863 | #if ! defined key_off_degrad && defined key_traldf_c2d |
---|
864 | ahtwdta(:,:,2) = zahtw(:,:) * tmask(:,:,1) |
---|
865 | #if defined key_trcldf_eiv |
---|
866 | aeiwdta(:,:,2) = zaeiw(:,:) * tmask(:,:,1) |
---|
867 | #endif |
---|
868 | #endif |
---|
869 | |
---|
870 | #if defined key_off_degrad |
---|
871 | ahtudta(:,:,:,2) = zahtu(:,:,:) * umask(:,:,:) |
---|
872 | ahtvdta(:,:,:,2) = zahtv(:,:,:) * vmask(:,:,:) |
---|
873 | ahtwdta(:,:,:,2) = zahtw(:,:,:) * tmask(:,:,:) |
---|
874 | # if defined key_trcldf_eiv |
---|
875 | aeiudta(:,:,:,2) = zaeiu(:,:,:) * umask(:,:,:) |
---|
876 | aeivdta(:,:,:,2) = zaeiv(:,:,:) * vmask(:,:,:) |
---|
877 | aeiwdta(:,:,:,2) = zaeiw(:,:,:) * tmask(:,:,:) |
---|
878 | # endif |
---|
879 | #endif |
---|
880 | |
---|
881 | ! |
---|
882 | ! flux : |
---|
883 | ! |
---|
884 | flxdta(:,:,jpwind,2) = zwspd(:,:) * tmask(:,:,1) |
---|
885 | flxdta(:,:,jpice,2) = MIN( 1., zice(:,:) ) * tmask(:,:,1) |
---|
886 | flxdta(:,:,jpemp,2) = zemp(:,:) * tmask(:,:,1) |
---|
887 | flxdta(:,:,jpqsr,2) = zqsr(:,:) * tmask(:,:,1) |
---|
888 | zmxldta(:,:,2) = zmld(:,:) * tmask(:,:,1) |
---|
889 | |
---|
890 | #if defined key_trcbbl_dif || defined key_trcbbl_adv |
---|
891 | bblxdta(:,:,2) = MAX( 0., zbblx(:,:) ) |
---|
892 | bblydta(:,:,2) = MAX( 0., zbbly(:,:) ) |
---|
893 | |
---|
894 | WHERE( bblxdta(:,:,2) > 2. ) bblxdta(:,:,2) = 0. |
---|
895 | WHERE( bblydta(:,:,2) > 2. ) bblydta(:,:,2) = 0. |
---|
896 | |
---|
897 | #endif |
---|
898 | |
---|
899 | IF( kt == nitend ) THEN |
---|
900 | CALL iom_close ( numfl_t ) |
---|
901 | CALL iom_close ( numfl_u ) |
---|
902 | CALL iom_close ( numfl_v ) |
---|
903 | CALL iom_close ( numfl_w ) |
---|
904 | ENDIF |
---|
905 | |
---|
906 | END SUBROUTINE dynrea |
---|
907 | |
---|
908 | |
---|
909 | |
---|
910 | END MODULE dtadyn |
---|