1 | MODULE diaar5 |
---|
2 | !!====================================================================== |
---|
3 | !! *** MODULE diaar5 *** |
---|
4 | !! AR5 diagnostics |
---|
5 | !!====================================================================== |
---|
6 | !! History : 3.2 ! 2009-11 (S. Masson) Original code |
---|
7 | !! 3.3 ! 2010-10 (C. Ethe, G. Madec) reorganisation of initialisation phase + merge TRC-TRA |
---|
8 | !!---------------------------------------------------------------------- |
---|
9 | !! dia_ar5 : AR5 diagnostics |
---|
10 | !! dia_ar5_init : initialisation of AR5 diagnostics |
---|
11 | !!---------------------------------------------------------------------- |
---|
12 | USE oce ! ocean dynamics and active tracers |
---|
13 | USE dom_oce ! ocean space and time domain |
---|
14 | USE eosbn2 ! equation of state (eos_bn2 routine) |
---|
15 | USE phycst ! physical constant |
---|
16 | USE in_out_manager ! I/O manager |
---|
17 | USE zdfddm |
---|
18 | USE zdf_oce |
---|
19 | ! |
---|
20 | USE lib_mpp ! distribued memory computing library |
---|
21 | USE iom ! I/O manager library |
---|
22 | USE fldread ! type FLD_N |
---|
23 | USE timing ! preformance summary |
---|
24 | |
---|
25 | IMPLICIT NONE |
---|
26 | PRIVATE |
---|
27 | |
---|
28 | PUBLIC dia_ar5 ! routine called in step.F90 module |
---|
29 | PUBLIC dia_ar5_alloc ! routine called in nemogcm.F90 module |
---|
30 | PUBLIC dia_ar5_hst ! heat/salt transport |
---|
31 | |
---|
32 | REAL(wp) :: vol0 ! ocean volume (interior domain) |
---|
33 | REAL(wp) :: area_tot ! total ocean surface (interior domain) |
---|
34 | REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,: ) :: area ! cell surface (interior domain) |
---|
35 | REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,: ) :: thick0 ! ocean thickness (interior domain) |
---|
36 | REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: sn0 ! initial salinity |
---|
37 | |
---|
38 | LOGICAL :: l_ar5 |
---|
39 | |
---|
40 | !! * Substitutions |
---|
41 | # include "vectopt_loop_substitute.h90" |
---|
42 | !!---------------------------------------------------------------------- |
---|
43 | !! NEMO/OCE 4.0 , NEMO Consortium (2018) |
---|
44 | !! $Id$ |
---|
45 | !! Software governed by the CeCILL license (see ./LICENSE) |
---|
46 | !!---------------------------------------------------------------------- |
---|
47 | CONTAINS |
---|
48 | |
---|
49 | FUNCTION dia_ar5_alloc() |
---|
50 | !!---------------------------------------------------------------------- |
---|
51 | !! *** ROUTINE dia_ar5_alloc *** |
---|
52 | !!---------------------------------------------------------------------- |
---|
53 | INTEGER :: dia_ar5_alloc |
---|
54 | !!---------------------------------------------------------------------- |
---|
55 | ! |
---|
56 | ALLOCATE( area(jpi,jpj), thick0(jpi,jpj) , sn0(jpi,jpj,jpk) , STAT=dia_ar5_alloc ) |
---|
57 | ! |
---|
58 | CALL mpp_sum ( 'diaar5', dia_ar5_alloc ) |
---|
59 | IF( dia_ar5_alloc /= 0 ) CALL ctl_stop( 'STOP', 'dia_ar5_alloc: failed to allocate arrays' ) |
---|
60 | ! |
---|
61 | END FUNCTION dia_ar5_alloc |
---|
62 | |
---|
63 | |
---|
64 | SUBROUTINE dia_ar5( kt ) |
---|
65 | !!---------------------------------------------------------------------- |
---|
66 | !! *** ROUTINE dia_ar5 *** |
---|
67 | !! |
---|
68 | !! ** Purpose : compute and output some AR5 diagnostics |
---|
69 | !!---------------------------------------------------------------------- |
---|
70 | ! |
---|
71 | INTEGER, INTENT( in ) :: kt ! ocean time-step index |
---|
72 | ! |
---|
73 | INTEGER :: ji, jj, jk ! dummy loop arguments |
---|
74 | REAL(wp) :: zvolssh, zvol, zssh_steric, zztmp, zarho, ztemp, zsal, zmass |
---|
75 | REAL(wp) :: zaw, zbw, zrw |
---|
76 | ! |
---|
77 | REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: zarea_ssh , zbotpres ! 2D workspace |
---|
78 | REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: zpe ! 2D workspace |
---|
79 | REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) :: zrhd , zrhop ! 3D workspace |
---|
80 | REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:) :: ztsn ! 4D workspace |
---|
81 | |
---|
82 | !!-------------------------------------------------------------------- |
---|
83 | IF( ln_timing ) CALL timing_start('dia_ar5') |
---|
84 | |
---|
85 | IF( kt == nit000 ) CALL dia_ar5_init |
---|
86 | |
---|
87 | IF( l_ar5 ) THEN |
---|
88 | ALLOCATE( zarea_ssh(jpi,jpj) , zbotpres(jpi,jpj) ) |
---|
89 | ALLOCATE( zrhd(jpi,jpj,jpk) , zrhop(jpi,jpj,jpk) ) |
---|
90 | ALLOCATE( ztsn(jpi,jpj,jpk,jpts) ) |
---|
91 | zarea_ssh(:,:) = area(:,:) * sshn(:,:) |
---|
92 | ENDIF |
---|
93 | ! |
---|
94 | IF( iom_use( 'voltot' ) .OR. iom_use( 'sshtot' ) .OR. iom_use( 'sshdyn' ) ) THEN |
---|
95 | ! ! total volume of liquid seawater |
---|
96 | zvolssh = SUM( zarea_ssh(:,:) ) |
---|
97 | CALL mpp_sum( 'diaar5', zvolssh ) |
---|
98 | zvol = vol0 + zvolssh |
---|
99 | |
---|
100 | CALL iom_put( 'voltot', zvol ) |
---|
101 | CALL iom_put( 'sshtot', zvolssh / area_tot ) |
---|
102 | CALL iom_put( 'sshdyn', sshn(:,:) - (zvolssh / area_tot) ) |
---|
103 | ! |
---|
104 | ENDIF |
---|
105 | |
---|
106 | IF( iom_use( 'botpres' ) .OR. iom_use( 'sshthster' ) .OR. iom_use( 'sshsteric' ) ) THEN |
---|
107 | ! |
---|
108 | ztsn(:,:,:,jp_tem) = tsn(:,:,:,jp_tem) ! thermosteric ssh |
---|
109 | ztsn(:,:,:,jp_sal) = sn0(:,:,:) |
---|
110 | CALL eos( ztsn, zrhd, gdept_n(:,:,:) ) ! now in situ density using initial salinity |
---|
111 | ! |
---|
112 | zbotpres(:,:) = 0._wp ! no atmospheric surface pressure, levitating sea-ice |
---|
113 | DO jk = 1, jpkm1 |
---|
114 | zbotpres(:,:) = zbotpres(:,:) + e3t_n(:,:,jk) * zrhd(:,:,jk) |
---|
115 | END DO |
---|
116 | IF( ln_linssh ) THEN |
---|
117 | IF( ln_isfcav ) THEN |
---|
118 | DO ji = 1, jpi |
---|
119 | DO jj = 1, jpj |
---|
120 | zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,mikt(ji,jj)) + riceload(ji,jj) |
---|
121 | END DO |
---|
122 | END DO |
---|
123 | ELSE |
---|
124 | zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1) |
---|
125 | END IF |
---|
126 | !!gm |
---|
127 | !!gm riceload should be added in both ln_linssh=T or F, no? |
---|
128 | !!gm |
---|
129 | END IF |
---|
130 | ! |
---|
131 | zarho = SUM( area(:,:) * zbotpres(:,:) ) |
---|
132 | CALL mpp_sum( 'diaar5', zarho ) |
---|
133 | zssh_steric = - zarho / area_tot |
---|
134 | CALL iom_put( 'sshthster', zssh_steric ) |
---|
135 | |
---|
136 | ! ! steric sea surface height |
---|
137 | CALL eos( tsn, zrhd, zrhop, gdept_n(:,:,:) ) ! now in situ and potential density |
---|
138 | zrhop(:,:,jpk) = 0._wp |
---|
139 | CALL iom_put( 'rhop', zrhop ) |
---|
140 | ! |
---|
141 | zbotpres(:,:) = 0._wp ! no atmospheric surface pressure, levitating sea-ice |
---|
142 | DO jk = 1, jpkm1 |
---|
143 | zbotpres(:,:) = zbotpres(:,:) + e3t_n(:,:,jk) * zrhd(:,:,jk) |
---|
144 | END DO |
---|
145 | IF( ln_linssh ) THEN |
---|
146 | IF ( ln_isfcav ) THEN |
---|
147 | DO ji = 1,jpi |
---|
148 | DO jj = 1,jpj |
---|
149 | zbotpres(ji,jj) = zbotpres(ji,jj) + sshn(ji,jj) * zrhd(ji,jj,mikt(ji,jj)) + riceload(ji,jj) |
---|
150 | END DO |
---|
151 | END DO |
---|
152 | ELSE |
---|
153 | zbotpres(:,:) = zbotpres(:,:) + sshn(:,:) * zrhd(:,:,1) |
---|
154 | END IF |
---|
155 | END IF |
---|
156 | ! |
---|
157 | zarho = SUM( area(:,:) * zbotpres(:,:) ) |
---|
158 | CALL mpp_sum( 'diaar5', zarho ) |
---|
159 | zssh_steric = - zarho / area_tot |
---|
160 | CALL iom_put( 'sshsteric', zssh_steric ) |
---|
161 | |
---|
162 | ! ! ocean bottom pressure |
---|
163 | zztmp = rau0 * grav * 1.e-4_wp ! recover pressure from pressure anomaly and cover to dbar = 1.e4 Pa |
---|
164 | zbotpres(:,:) = zztmp * ( zbotpres(:,:) + sshn(:,:) + thick0(:,:) ) |
---|
165 | CALL iom_put( 'botpres', zbotpres ) |
---|
166 | ! |
---|
167 | ENDIF |
---|
168 | |
---|
169 | IF( iom_use( 'masstot' ) .OR. iom_use( 'temptot' ) .OR. iom_use( 'saltot' ) ) THEN |
---|
170 | ! ! Mean density anomalie, temperature and salinity |
---|
171 | ztemp = 0._wp |
---|
172 | zsal = 0._wp |
---|
173 | DO jk = 1, jpkm1 |
---|
174 | DO jj = 1, jpj |
---|
175 | DO ji = 1, jpi |
---|
176 | zztmp = area(ji,jj) * e3t_n(ji,jj,jk) |
---|
177 | ztemp = ztemp + zztmp * tsn(ji,jj,jk,jp_tem) |
---|
178 | zsal = zsal + zztmp * tsn(ji,jj,jk,jp_sal) |
---|
179 | END DO |
---|
180 | END DO |
---|
181 | END DO |
---|
182 | IF( ln_linssh ) THEN |
---|
183 | IF( ln_isfcav ) THEN |
---|
184 | DO ji = 1, jpi |
---|
185 | DO jj = 1, jpj |
---|
186 | ztemp = ztemp + zarea_ssh(ji,jj) * tsn(ji,jj,mikt(ji,jj),jp_tem) |
---|
187 | zsal = zsal + zarea_ssh(ji,jj) * tsn(ji,jj,mikt(ji,jj),jp_sal) |
---|
188 | END DO |
---|
189 | END DO |
---|
190 | ELSE |
---|
191 | ztemp = ztemp + SUM( zarea_ssh(:,:) * tsn(:,:,1,jp_tem) ) |
---|
192 | zsal = zsal + SUM( zarea_ssh(:,:) * tsn(:,:,1,jp_sal) ) |
---|
193 | END IF |
---|
194 | ENDIF |
---|
195 | IF( lk_mpp ) THEN |
---|
196 | CALL mpp_sum( 'diaar5', ztemp ) |
---|
197 | CALL mpp_sum( 'diaar5', zsal ) |
---|
198 | END IF |
---|
199 | ! |
---|
200 | zmass = rau0 * ( zarho + zvol ) ! total mass of liquid seawater |
---|
201 | ztemp = ztemp / zvol ! potential temperature in liquid seawater |
---|
202 | zsal = zsal / zvol ! Salinity of liquid seawater |
---|
203 | ! |
---|
204 | CALL iom_put( 'masstot', zmass ) |
---|
205 | CALL iom_put( 'temptot', ztemp ) |
---|
206 | CALL iom_put( 'saltot' , zsal ) |
---|
207 | ! |
---|
208 | ENDIF |
---|
209 | |
---|
210 | IF( iom_use( 'tnpeo' )) THEN |
---|
211 | ! Work done against stratification by vertical mixing |
---|
212 | ! Exclude points where rn2 is negative as convection kicks in here and |
---|
213 | ! work is not being done against stratification |
---|
214 | ALLOCATE( zpe(jpi,jpj) ) |
---|
215 | zpe(:,:) = 0._wp |
---|
216 | IF( ln_zdfddm ) THEN |
---|
217 | DO jk = 2, jpk |
---|
218 | DO jj = 1, jpj |
---|
219 | DO ji = 1, jpi |
---|
220 | IF( rn2(ji,jj,jk) > 0._wp ) THEN |
---|
221 | zrw = ( gdepw_n(ji,jj,jk ) - gdept_n(ji,jj,jk) ) & |
---|
222 | & / ( gdept_n(ji,jj,jk-1) - gdept_n(ji,jj,jk) ) |
---|
223 | !!gm this can be reduced to : (depw-dept) / e3w (NB idem dans bn2 !) |
---|
224 | ! zrw = ( gdept_n(ji,jj,jk) - gdepw_n(ji,jj,jk) ) / e3w_n(ji,jj,jk) |
---|
225 | !!gm end |
---|
226 | ! |
---|
227 | zaw = rab_n(ji,jj,jk,jp_tem) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_tem)* zrw |
---|
228 | zbw = rab_n(ji,jj,jk,jp_sal) * (1. - zrw) + rab_n(ji,jj,jk-1,jp_sal)* zrw |
---|
229 | ! |
---|
230 | zpe(ji, jj) = zpe(ji, jj) & |
---|
231 | & - grav * ( avt(ji,jj,jk) * zaw * (tsn(ji,jj,jk-1,jp_tem) - tsn(ji,jj,jk,jp_tem) ) & |
---|
232 | & - avs(ji,jj,jk) * zbw * (tsn(ji,jj,jk-1,jp_sal) - tsn(ji,jj,jk,jp_sal) ) ) |
---|
233 | ENDIF |
---|
234 | END DO |
---|
235 | END DO |
---|
236 | END DO |
---|
237 | ELSE |
---|
238 | DO jk = 1, jpk |
---|
239 | DO ji = 1, jpi |
---|
240 | DO jj = 1, jpj |
---|
241 | zpe(ji,jj) = zpe(ji,jj) + avt(ji, jj, jk) * MIN(0._wp,rn2(ji, jj, jk)) * rau0 * e3w_n(ji, jj, jk) |
---|
242 | END DO |
---|
243 | END DO |
---|
244 | END DO |
---|
245 | ENDIF |
---|
246 | !!gm useless lbc_lnk since the computation above is performed over 1:jpi & 1:jpj |
---|
247 | !!gm CALL lbc_lnk( 'diaar5', zpe, 'T', 1._wp) |
---|
248 | CALL iom_put( 'tnpeo', zpe ) |
---|
249 | DEALLOCATE( zpe ) |
---|
250 | ENDIF |
---|
251 | |
---|
252 | IF( l_ar5 ) THEN |
---|
253 | DEALLOCATE( zarea_ssh , zbotpres ) |
---|
254 | DEALLOCATE( zrhd , zrhop ) |
---|
255 | DEALLOCATE( ztsn ) |
---|
256 | ENDIF |
---|
257 | ! |
---|
258 | IF( ln_timing ) CALL timing_stop('dia_ar5') |
---|
259 | ! |
---|
260 | END SUBROUTINE dia_ar5 |
---|
261 | |
---|
262 | |
---|
263 | SUBROUTINE dia_ar5_hst( ktra, cptr, pua, pva ) |
---|
264 | !!---------------------------------------------------------------------- |
---|
265 | !! *** ROUTINE dia_ar5_htr *** |
---|
266 | !!---------------------------------------------------------------------- |
---|
267 | !! Wrapper for heat transport calculations |
---|
268 | !! Called from all advection and/or diffusion routines |
---|
269 | !!---------------------------------------------------------------------- |
---|
270 | INTEGER , INTENT(in ) :: ktra ! tracer index |
---|
271 | CHARACTER(len=3) , INTENT(in) :: cptr ! transport type 'adv'/'ldf' |
---|
272 | REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in) :: pua ! 3D input array of advection/diffusion |
---|
273 | REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in) :: pva ! 3D input array of advection/diffusion |
---|
274 | ! |
---|
275 | INTEGER :: ji, jj, jk |
---|
276 | REAL(wp), DIMENSION(jpi,jpj) :: z2d |
---|
277 | |
---|
278 | |
---|
279 | z2d(:,:) = pua(:,:,1) |
---|
280 | DO jk = 1, jpkm1 |
---|
281 | DO jj = 2, jpjm1 |
---|
282 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
283 | z2d(ji,jj) = z2d(ji,jj) + pua(ji,jj,jk) |
---|
284 | END DO |
---|
285 | END DO |
---|
286 | END DO |
---|
287 | CALL lbc_lnk( 'diaar5', z2d, 'U', -1. ) |
---|
288 | IF( cptr == 'adv' ) THEN |
---|
289 | IF( ktra == jp_tem ) CALL iom_put( "uadv_heattr" , rau0_rcp * z2d ) ! advective heat transport in i-direction |
---|
290 | IF( ktra == jp_sal ) CALL iom_put( "uadv_salttr" , rau0 * z2d ) ! advective salt transport in i-direction |
---|
291 | ENDIF |
---|
292 | IF( cptr == 'ldf' ) THEN |
---|
293 | IF( ktra == jp_tem ) CALL iom_put( "udiff_heattr" , rau0_rcp * z2d ) ! diffusive heat transport in i-direction |
---|
294 | IF( ktra == jp_sal ) CALL iom_put( "udiff_salttr" , rau0 * z2d ) ! diffusive salt transport in i-direction |
---|
295 | ENDIF |
---|
296 | ! |
---|
297 | z2d(:,:) = pva(:,:,1) |
---|
298 | DO jk = 1, jpkm1 |
---|
299 | DO jj = 2, jpjm1 |
---|
300 | DO ji = fs_2, fs_jpim1 ! vector opt. |
---|
301 | z2d(ji,jj) = z2d(ji,jj) + pva(ji,jj,jk) |
---|
302 | END DO |
---|
303 | END DO |
---|
304 | END DO |
---|
305 | CALL lbc_lnk( 'diaar5', z2d, 'V', -1. ) |
---|
306 | IF( cptr == 'adv' ) THEN |
---|
307 | IF( ktra == jp_tem ) CALL iom_put( "vadv_heattr" , rau0_rcp * z2d ) ! advective heat transport in j-direction |
---|
308 | IF( ktra == jp_sal ) CALL iom_put( "vadv_salttr" , rau0 * z2d ) ! advective salt transport in j-direction |
---|
309 | ENDIF |
---|
310 | IF( cptr == 'ldf' ) THEN |
---|
311 | IF( ktra == jp_tem ) CALL iom_put( "vdiff_heattr" , rau0_rcp * z2d ) ! diffusive heat transport in j-direction |
---|
312 | IF( ktra == jp_sal ) CALL iom_put( "vdiff_salttr" , rau0 * z2d ) ! diffusive salt transport in j-direction |
---|
313 | ENDIF |
---|
314 | |
---|
315 | END SUBROUTINE dia_ar5_hst |
---|
316 | |
---|
317 | |
---|
318 | SUBROUTINE dia_ar5_init |
---|
319 | !!---------------------------------------------------------------------- |
---|
320 | !! *** ROUTINE dia_ar5_init *** |
---|
321 | !! |
---|
322 | !! ** Purpose : initialization for AR5 diagnostic computation |
---|
323 | !!---------------------------------------------------------------------- |
---|
324 | INTEGER :: inum |
---|
325 | INTEGER :: ik |
---|
326 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
327 | REAL(wp) :: zztmp |
---|
328 | REAL(wp), ALLOCATABLE, DIMENSION(:,:,:,:) :: zsaldta ! Jan/Dec levitus salinity |
---|
329 | ! |
---|
330 | !!---------------------------------------------------------------------- |
---|
331 | ! |
---|
332 | l_ar5 = .FALSE. |
---|
333 | IF( iom_use( 'voltot' ) .OR. iom_use( 'sshtot' ) .OR. iom_use( 'sshdyn' ) .OR. & |
---|
334 | & iom_use( 'masstot' ) .OR. iom_use( 'temptot' ) .OR. iom_use( 'saltot' ) .OR. & |
---|
335 | & iom_use( 'botpres' ) .OR. iom_use( 'sshthster' ) .OR. iom_use( 'sshsteric' ) ) L_ar5 = .TRUE. |
---|
336 | |
---|
337 | IF( l_ar5 ) THEN |
---|
338 | ! |
---|
339 | ! ! allocate dia_ar5 arrays |
---|
340 | IF( dia_ar5_alloc() /= 0 ) CALL ctl_stop( 'STOP', 'dia_ar5_init : unable to allocate arrays' ) |
---|
341 | |
---|
342 | area(:,:) = e1e2t(:,:) * tmask_i(:,:) |
---|
343 | |
---|
344 | area_tot = SUM( area(:,:) ) ; CALL mpp_sum( 'diaar5', area_tot ) |
---|
345 | |
---|
346 | vol0 = 0._wp |
---|
347 | thick0(:,:) = 0._wp |
---|
348 | DO jk = 1, jpkm1 |
---|
349 | vol0 = vol0 + SUM( area (:,:) * tmask(:,:,jk) * e3t_0(:,:,jk) ) |
---|
350 | thick0(:,:) = thick0(:,:) + tmask_i(:,:) * tmask(:,:,jk) * e3t_0(:,:,jk) |
---|
351 | END DO |
---|
352 | CALL mpp_sum( 'diaar5', vol0 ) |
---|
353 | |
---|
354 | IF( iom_use( 'sshthster' ) ) THEN |
---|
355 | ALLOCATE( zsaldta(jpi,jpj,jpj,jpts) ) |
---|
356 | CALL iom_open ( 'sali_ref_clim_monthly', inum ) |
---|
357 | CALL iom_get ( inum, jpdom_data, 'vosaline' , zsaldta(:,:,:,1), 1 ) |
---|
358 | CALL iom_get ( inum, jpdom_data, 'vosaline' , zsaldta(:,:,:,2), 12 ) |
---|
359 | CALL iom_close( inum ) |
---|
360 | |
---|
361 | sn0(:,:,:) = 0.5_wp * ( zsaldta(:,:,:,1) + zsaldta(:,:,:,2) ) |
---|
362 | sn0(:,:,:) = sn0(:,:,:) * tmask(:,:,:) |
---|
363 | IF( ln_zps ) THEN ! z-coord. partial steps |
---|
364 | DO jj = 1, jpj ! interpolation of salinity at the last ocean level (i.e. the partial step) |
---|
365 | DO ji = 1, jpi |
---|
366 | ik = mbkt(ji,jj) |
---|
367 | IF( ik > 1 ) THEN |
---|
368 | zztmp = ( gdept_1d(ik) - gdept_0(ji,jj,ik) ) / ( gdept_1d(ik) - gdept_1d(ik-1) ) |
---|
369 | sn0(ji,jj,ik) = ( 1._wp - zztmp ) * sn0(ji,jj,ik) + zztmp * sn0(ji,jj,ik-1) |
---|
370 | ENDIF |
---|
371 | END DO |
---|
372 | END DO |
---|
373 | ENDIF |
---|
374 | ! |
---|
375 | DEALLOCATE( zsaldta ) |
---|
376 | ENDIF |
---|
377 | ! |
---|
378 | ENDIF |
---|
379 | ! |
---|
380 | END SUBROUTINE dia_ar5_init |
---|
381 | |
---|
382 | !!====================================================================== |
---|
383 | END MODULE diaar5 |
---|