source: palm/trunk/SOURCE/data_output_3d.f90 @ 1612

Last change on this file since 1612 was 1586, checked in by maronga, 9 years ago

last commit documented

  • Property svn:keywords set to Id
File size: 22.9 KB
Line 
1 SUBROUTINE data_output_3d( av )
2
3!--------------------------------------------------------------------------------!
4! This file is part of PALM.
5!
6! PALM is free software: you can redistribute it and/or modify it under the terms
7! of the GNU General Public License as published by the Free Software Foundation,
8! either version 3 of the License, or (at your option) any later version.
9!
10! PALM is distributed in the hope that it will be useful, but WITHOUT ANY
11! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR
12! A PARTICULAR PURPOSE.  See the GNU General Public License for more details.
13!
14! You should have received a copy of the GNU General Public License along with
15! PALM. If not, see <http://www.gnu.org/licenses/>.
16!
17! Copyright 1997-2014 Leibniz Universitaet Hannover
18!--------------------------------------------------------------------------------!
19!
20! Current revisions:
21! ------------------
22!
23!
24! Former revisions:
25! -----------------
26! $Id: data_output_3d.f90 1586 2015-04-30 07:11:33Z maronga $
27!
28! 1585 2015-04-30 07:05:52Z maronga
29! Added support for RRTMG
30!
31! 1551 2015-03-03 14:18:16Z maronga
32! Added suppport for land surface model and radiation model output. In the course
33! of this action, the limits for vertical loops have been changed (from nzb and
34! nzt+1 to nzb_do and nzt_do, respectively in order to allow soil model output).
35! Moreover, a new vertical grid zs was introduced.
36!
37! 1359 2014-04-11 17:15:14Z hoffmann
38! New particle structure integrated.
39!
40! 1353 2014-04-08 15:21:23Z heinze
41! REAL constants provided with KIND-attribute
42!
43! 1327 2014-03-21 11:00:16Z raasch
44! parts concerning avs output removed,
45! -netcdf output queries
46!
47! 1320 2014-03-20 08:40:49Z raasch
48! ONLY-attribute added to USE-statements,
49! kind-parameters added to all INTEGER and REAL declaration statements,
50! kinds are defined in new module kinds,
51! old module precision_kind is removed,
52! revision history before 2012 removed,
53! comment fields (!:) to be used for variable explanations added to
54! all variable declaration statements
55!
56! 1318 2014-03-17 13:35:16Z raasch
57! barrier argument removed from cpu_log,
58! module interfaces removed
59!
60! 1308 2014-03-13 14:58:42Z fricke
61! Check, if the limit of the time dimension is exceeded for parallel output
62! To increase the performance for parallel output, the following is done:
63! - Update of time axis is only done by PE0
64!
65! 1244 2013-10-31 08:16:56Z raasch
66! Bugfix for index bounds in case of 3d-parallel output
67!
68! 1115 2013-03-26 18:16:16Z hoffmann
69! ql is calculated by calc_liquid_water_content
70!
71! 1106 2013-03-04 05:31:38Z raasch
72! array_kind renamed precision_kind
73!
74! 1076 2012-12-05 08:30:18Z hoffmann
75! Bugfix in output of ql
76!
77! 1053 2012-11-13 17:11:03Z hoffmann
78! +nr, qr, prr, qc and averaged quantities
79!
80! 1036 2012-10-22 13:43:42Z raasch
81! code put under GPL (PALM 3.9)
82!
83! 1031 2012-10-19 14:35:30Z raasch
84! netCDF4 without parallel file support implemented
85!
86! 1007 2012-09-19 14:30:36Z franke
87! Bugfix: missing calculation of ql_vp added
88!
89! Revision 1.1  1997/09/03 06:29:36  raasch
90! Initial revision
91!
92!
93! Description:
94! ------------
95! Output of the 3D-arrays in netCDF and/or AVS format.
96!------------------------------------------------------------------------------!
97
98    USE arrays_3d,                                                             &
99        ONLY:  e, nr, p, pt, q, qc, ql, ql_c, ql_v, qr, rho, sa, tend, u, v,   &
100               vpt, w
101       
102    USE averaging
103       
104    USE cloud_parameters,                                                      &
105        ONLY:  l_d_cp, prr, pt_d_t
106       
107    USE control_parameters,                                                    &
108        ONLY:  avs_data_file, cloud_physics, do3d, do3d_avs_n,                 &
109               do3d_no, do3d_time_count, io_blocks, io_group,                  &
110               message_string, netcdf_data_format, ntdim_3d,                   &
111               nz_do3d, plot_3d_precision, psolver, simulated_time,            &
112               simulated_time_chr, skip_do_avs, time_since_reference_point
113       
114    USE cpulog,                                                                &
115        ONLY:  log_point, cpu_log
116       
117    USE indices,                                                               &
118        ONLY:  nbgp, nx, nxl, nxlg, nxr, nxrg, ny, nyn, nyng, nys, nysg, nzt,  &
119               nzb
120       
121    USE kinds
122   
123    USE land_surface_model_mod,                                                &
124        ONLY: m_soil, m_soil_av, nzb_soil, nzt_soil, t_soil, t_soil_av
125
126    USE netcdf_control
127       
128    USE particle_attributes,                                                   &
129        ONLY:  grid_particles, number_of_particles, particles,                 &
130               particle_advection_start, prt_count
131       
132    USE pegrid
133
134    USE radiation_model_mod,                                                   &
135        ONLY:  rad_lw_in, rad_lw_in_av, rad_lw_out, rad_lw_out_av,             &
136               rad_sw_in, rad_sw_in_av, rad_sw_out, rad_sw_out_av
137
138
139    IMPLICIT NONE
140
141    CHARACTER (LEN=9) ::  simulated_time_mod  !:
142
143    INTEGER(iwp) ::  av        !:
144    INTEGER(iwp) ::  i         !:
145    INTEGER(iwp) ::  if        !:
146    INTEGER(iwp) ::  j         !:
147    INTEGER(iwp) ::  k         !:
148    INTEGER(iwp) ::  n         !:
149    INTEGER(iwp) ::  nzb_do    !: vertical lower limit for data output
150    INTEGER(iwp) ::  nzt_do    !: vertical upper limit for data output
151    INTEGER(iwp) ::  pos       !:
152    INTEGER(iwp) ::  prec      !:
153    INTEGER(iwp) ::  psi       !:
154
155    LOGICAL      ::  found     !:
156    LOGICAL      ::  resorted  !:
157
158    REAL(wp)     ::  mean_r    !:
159    REAL(wp)     ::  s_r2      !:
160    REAL(wp)     ::  s_r3      !:
161
162    REAL(sp), DIMENSION(:,:,:), ALLOCATABLE ::  local_pf  !:
163
164    REAL(wp), DIMENSION(:,:,:), POINTER ::  to_be_resorted  !:
165
166!
167!-- Return, if nothing to output
168    IF ( do3d_no(av) == 0 )  RETURN
169!
170!-- For parallel netcdf output the time axis must be limited. Return, if this
171!-- limit is exceeded. This could be the case, if the simulated time exceeds
172!-- the given end time by the length of the given output interval.
173    IF ( netcdf_data_format > 4 )  THEN
174       IF ( do3d_time_count(av) + 1 > ntdim_3d(av) )  THEN
175          WRITE ( message_string, * ) 'Output of 3d data is not given at t=',  &
176                                      simulated_time, '&because the maximum ', & 
177                                      'number of output time levels is ',      &
178                                      'exceeded.'
179          CALL message( 'data_output_3d', 'PA0387', 0, 1, 0, 6, 0 )         
180          RETURN
181       ENDIF
182    ENDIF
183
184    CALL cpu_log (log_point(14),'data_output_3d','start')
185
186!
187!-- Open output file.
188!-- Also creates coordinate and fld-file for AVS.
189!-- For classic or 64bit netCDF output or output of other (old) data formats,
190!-- for a run on more than one PE, each PE opens its own file and
191!-- writes the data of its subdomain in binary format (regardless of the format
192!-- the user has requested). After the run, these files are combined to one
193!-- file by combine_plot_fields in the format requested by the user (netcdf
194!-- and/or avs).
195!-- For netCDF4/HDF5 output, data is written in parallel into one file.
196    IF ( netcdf_data_format < 5 )  THEN
197       CALL check_open( 30 )
198       IF ( myid == 0 )  CALL check_open( 106+av*10 )
199    ELSE
200       CALL check_open( 106+av*10 )
201    ENDIF
202
203!
204!-- Update the netCDF time axis
205!-- In case of parallel output, this is only done by PE0 to increase the
206!-- performance.
207#if defined( __netcdf )
208    do3d_time_count(av) = do3d_time_count(av) + 1
209    IF ( myid == 0 )  THEN
210       nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_time_3d(av),           &
211                               (/ time_since_reference_point /),            &
212                               start = (/ do3d_time_count(av) /),           &
213                               count = (/ 1 /) )
214       CALL handle_netcdf_error( 'data_output_3d', 376 )
215    ENDIF
216#endif
217
218!
219!-- Loop over all variables to be written.
220    if = 1
221
222    DO  WHILE ( do3d(av,if)(1:1) /= ' ' )
223!
224!--    Store the array chosen on the temporary array.
225       resorted = .FALSE.
226       nzb_do = nzb
227       nzt_do = nz_do3d
228
229!
230!--    Allocate a temporary array with the desired output dimensions.
231       ALLOCATE( local_pf(nxlg:nxrg,nysg:nyng,nzb_do:nzt_do) )
232
233       SELECT CASE ( TRIM( do3d(av,if) ) )
234
235          CASE ( 'e' )
236             IF ( av == 0 )  THEN
237                to_be_resorted => e
238             ELSE
239                to_be_resorted => e_av
240             ENDIF
241
242          CASE ( 'lpt' )
243             IF ( av == 0 )  THEN
244                to_be_resorted => pt
245             ELSE
246                to_be_resorted => lpt_av
247             ENDIF
248
249          CASE ( 'm_soil' )
250             nzb_do = nzb_soil
251             nzt_do = nzt_soil
252!
253!--          For soil model quantities, it is required to re-allocate local_pf
254             DEALLOCATE ( local_pf )
255             ALLOCATE( local_pf(nxlg:nxrg,nysg:nyng,nzb_do:nzt_do) )
256
257             IF ( av == 0 )  THEN
258                to_be_resorted => m_soil
259             ELSE
260                to_be_resorted => m_soil_av
261             ENDIF
262
263          CASE ( 'nr' )
264             IF ( av == 0 )  THEN
265                to_be_resorted => nr
266             ELSE
267                to_be_resorted => nr_av
268             ENDIF
269
270          CASE ( 'p' )
271             IF ( av == 0 )  THEN
272                IF ( psolver /= 'sor' )  CALL exchange_horiz( p, nbgp )
273                to_be_resorted => p
274             ELSE
275                IF ( psolver /= 'sor' )  CALL exchange_horiz( p_av, nbgp )
276                to_be_resorted => p_av
277             ENDIF
278
279          CASE ( 'pc' )  ! particle concentration (requires ghostpoint exchange)
280             IF ( av == 0 )  THEN
281                IF ( simulated_time >= particle_advection_start )  THEN
282                   tend = prt_count
283                   CALL exchange_horiz( tend, nbgp )
284                ELSE
285                   tend = 0.0_wp
286                ENDIF
287                DO  i = nxlg, nxrg
288                   DO  j = nysg, nyng
289                      DO  k = nzb_do, nzt_do
290                         local_pf(i,j,k) = tend(k,j,i)
291                      ENDDO
292                   ENDDO
293                ENDDO
294                resorted = .TRUE.
295             ELSE
296                CALL exchange_horiz( pc_av, nbgp )
297                to_be_resorted => pc_av
298             ENDIF
299
300          CASE ( 'pr' )  ! mean particle radius (effective radius)
301             IF ( av == 0 )  THEN
302                IF ( simulated_time >= particle_advection_start )  THEN
303                   DO  i = nxl, nxr
304                      DO  j = nys, nyn
305                         DO  k = nzb_do, nzt_do
306                            number_of_particles = prt_count(k,j,i)
307                            IF (number_of_particles <= 0)  CYCLE
308                            particles => grid_particles(k,j,i)%particles(1:number_of_particles)
309                            s_r2 = 0.0_wp
310                            s_r3 = 0.0_wp
311                            DO  n = 1, number_of_particles
312                               IF ( particles(n)%particle_mask )  THEN
313                                  s_r2 = s_r2 + particles(n)%radius**2 * &
314                                         particles(n)%weight_factor
315                                  s_r3 = s_r3 + particles(n)%radius**3 * &
316                                         particles(n)%weight_factor
317                               ENDIF
318                            ENDDO
319                            IF ( s_r2 > 0.0_wp )  THEN
320                               mean_r = s_r3 / s_r2
321                            ELSE
322                               mean_r = 0.0_wp
323                            ENDIF
324                            tend(k,j,i) = mean_r
325                         ENDDO
326                      ENDDO
327                   ENDDO
328                   CALL exchange_horiz( tend, nbgp )
329                ELSE
330                   tend = 0.0_wp
331                ENDIF
332                DO  i = nxlg, nxrg
333                   DO  j = nysg, nyng
334                      DO  k = nzb_do, nzt_do
335                         local_pf(i,j,k) = tend(k,j,i)
336                      ENDDO
337                   ENDDO
338                ENDDO
339                resorted = .TRUE.
340             ELSE
341                CALL exchange_horiz( pr_av, nbgp )
342                to_be_resorted => pr_av
343             ENDIF
344
345          CASE ( 'prr' )
346             IF ( av == 0 )  THEN
347                CALL exchange_horiz( prr, nbgp )
348                DO  i = nxlg, nxrg
349                   DO  j = nysg, nyng
350                      DO  k = nzb, nzt+1
351                         local_pf(i,j,k) = prr(k,j,i)
352                      ENDDO
353                   ENDDO
354                ENDDO
355             ELSE
356                CALL exchange_horiz( prr_av, nbgp )
357                DO  i = nxlg, nxrg
358                   DO  j = nysg, nyng
359                      DO  k = nzb, nzt+1
360                         local_pf(i,j,k) = prr_av(k,j,i)
361                      ENDDO
362                   ENDDO
363                ENDDO
364             ENDIF
365             resorted = .TRUE.
366
367          CASE ( 'pt' )
368             IF ( av == 0 )  THEN
369                IF ( .NOT. cloud_physics ) THEN
370                   to_be_resorted => pt
371                ELSE
372                   DO  i = nxlg, nxrg
373                      DO  j = nysg, nyng
374                         DO  k = nzb_do, nzt_do
375                            local_pf(i,j,k) = pt(k,j,i) + l_d_cp *             &
376                                                          pt_d_t(k) *          &
377                                                          ql(k,j,i)
378                         ENDDO
379                      ENDDO
380                   ENDDO
381                   resorted = .TRUE.
382                ENDIF
383             ELSE
384                to_be_resorted => pt_av
385             ENDIF
386
387          CASE ( 'q' )
388             IF ( av == 0 )  THEN
389                to_be_resorted => q
390             ELSE
391                to_be_resorted => q_av
392             ENDIF
393
394          CASE ( 'qc' )
395             IF ( av == 0 )  THEN
396                to_be_resorted => qc
397             ELSE
398                to_be_resorted => qc_av
399             ENDIF
400
401          CASE ( 'ql' )
402             IF ( av == 0 )  THEN
403                to_be_resorted => ql
404             ELSE
405                to_be_resorted => ql_av
406             ENDIF
407
408          CASE ( 'ql_c' )
409             IF ( av == 0 )  THEN
410                to_be_resorted => ql_c
411             ELSE
412                to_be_resorted => ql_c_av
413             ENDIF
414
415          CASE ( 'ql_v' )
416             IF ( av == 0 )  THEN
417                to_be_resorted => ql_v
418             ELSE
419                to_be_resorted => ql_v_av
420             ENDIF
421
422          CASE ( 'ql_vp' )
423             IF ( av == 0 )  THEN
424                IF ( simulated_time >= particle_advection_start )  THEN
425                   DO  i = nxl, nxr
426                      DO  j = nys, nyn
427                         DO  k = nzb_do, nzt_do
428                            number_of_particles = prt_count(k,j,i)
429                            IF (number_of_particles <= 0)  CYCLE
430                            particles => grid_particles(k,j,i)%particles(1:number_of_particles)
431                            DO  n = 1, number_of_particles
432                               IF ( particles(n)%particle_mask )  THEN
433                                  tend(k,j,i) =  tend(k,j,i) +                 &
434                                                 particles(n)%weight_factor /  &
435                                                 prt_count(k,j,i)
436                               ENDIF
437                            ENDDO
438                         ENDDO
439                      ENDDO
440                   ENDDO
441                   CALL exchange_horiz( tend, nbgp )
442                ELSE
443                   tend = 0.0_wp
444                ENDIF
445                DO  i = nxlg, nxrg
446                   DO  j = nysg, nyng
447                      DO  k = nzb_do, nzt_do
448                         local_pf(i,j,k) = tend(k,j,i)
449                      ENDDO
450                   ENDDO
451                ENDDO
452                resorted = .TRUE.
453             ELSE
454                CALL exchange_horiz( ql_vp_av, nbgp )
455                to_be_resorted => ql_vp_av
456             ENDIF
457
458          CASE ( 'qr' )
459             IF ( av == 0 )  THEN
460                to_be_resorted => qr
461             ELSE
462                to_be_resorted => qr_av
463             ENDIF
464
465          CASE ( 'qv' )
466             IF ( av == 0 )  THEN
467                DO  i = nxlg, nxrg
468                   DO  j = nysg, nyng
469                      DO  k = nzb_do, nzt_do
470                         local_pf(i,j,k) = q(k,j,i) - ql(k,j,i)
471                      ENDDO
472                   ENDDO
473                ENDDO
474                resorted = .TRUE.
475             ELSE
476                to_be_resorted => qv_av
477             ENDIF
478
479          CASE ( 'rad_sw_in' )
480             IF ( av == 0 )  THEN
481                to_be_resorted => rad_sw_in
482             ELSE
483                to_be_resorted => rad_sw_in_av
484             ENDIF
485
486          CASE ( 'rad_sw_out' )
487             IF ( av == 0 )  THEN
488                to_be_resorted => rad_sw_out
489             ELSE
490                to_be_resorted => rad_sw_out_av
491             ENDIF
492
493          CASE ( 'rad_lw_in' )
494             IF ( av == 0 )  THEN
495                to_be_resorted => rad_lw_in
496             ELSE
497                to_be_resorted => rad_lw_in_av
498             ENDIF
499
500          CASE ( 'rad_lw_out' )
501             IF ( av == 0 )  THEN
502                to_be_resorted => rad_lw_out
503             ELSE
504                to_be_resorted => rad_lw_out_av
505             ENDIF
506
507          CASE ( 'rho' )
508             IF ( av == 0 )  THEN
509                to_be_resorted => rho
510             ELSE
511                to_be_resorted => rho_av
512             ENDIF
513
514          CASE ( 's' )
515             IF ( av == 0 )  THEN
516                to_be_resorted => q
517             ELSE
518                to_be_resorted => s_av
519             ENDIF
520
521          CASE ( 'sa' )
522             IF ( av == 0 )  THEN
523                to_be_resorted => sa
524             ELSE
525                to_be_resorted => sa_av
526             ENDIF
527
528          CASE ( 't_soil' )
529             nzb_do = nzb_soil
530             nzt_do = nzt_soil
531!
532!--          For soil model quantities, it is required to re-allocate local_pf
533             DEALLOCATE ( local_pf )
534             ALLOCATE( local_pf(nxlg:nxrg,nysg:nyng,nzb_do:nzt_do) )
535
536             IF ( av == 0 )  THEN
537                to_be_resorted => t_soil
538             ELSE
539                to_be_resorted => t_soil_av
540             ENDIF
541
542          CASE ( 'u' )
543             IF ( av == 0 )  THEN
544                to_be_resorted => u
545             ELSE
546                to_be_resorted => u_av
547             ENDIF
548
549          CASE ( 'v' )
550             IF ( av == 0 )  THEN
551                to_be_resorted => v
552             ELSE
553                to_be_resorted => v_av
554             ENDIF
555
556          CASE ( 'vpt' )
557             IF ( av == 0 )  THEN
558                to_be_resorted => vpt
559             ELSE
560                to_be_resorted => vpt_av
561             ENDIF
562
563          CASE ( 'w' )
564             IF ( av == 0 )  THEN
565                to_be_resorted => w
566             ELSE
567                to_be_resorted => w_av
568             ENDIF
569
570          CASE DEFAULT
571!
572!--          User defined quantity
573             CALL user_data_output_3d( av, do3d(av,if), found, local_pf,       &
574                                       nzb_do, nzt_do )
575             resorted = .TRUE.
576
577             IF ( .NOT. found )  THEN
578                message_string =  'no output available for: ' //               &
579                                  TRIM( do3d(av,if) )
580                CALL message( 'data_output_3d', 'PA0182', 0, 0, 0, 6, 0 )
581             ENDIF
582
583       END SELECT
584
585!
586!--    Resort the array to be output, if not done above
587       IF ( .NOT. resorted )  THEN
588          DO  i = nxlg, nxrg
589             DO  j = nysg, nyng
590                DO  k = nzb_do, nzt_do
591                   local_pf(i,j,k) = to_be_resorted(k,j,i)
592                ENDDO
593             ENDDO
594          ENDDO
595       ENDIF
596
597!
598!--    Output of the 3D-array
599#if defined( __parallel )
600       IF ( netcdf_data_format < 5 )  THEN
601!
602!--       Non-parallel netCDF output. Data is output in parallel in
603!--       FORTRAN binary format here, and later collected into one file by
604!--       combine_plot_fields
605          IF ( myid == 0 )  THEN
606             WRITE ( 30 )  time_since_reference_point,                   &
607                           do3d_time_count(av), av
608          ENDIF
609          DO  i = 0, io_blocks-1
610             IF ( i == io_group )  THEN
611                WRITE ( 30 )  nxlg, nxrg, nysg, nyng, nzb_do, nzt_do
612                WRITE ( 30 )  local_pf(:,:,nzb_do:nzt_do)
613             ENDIF
614#if defined( __parallel )
615             CALL MPI_BARRIER( comm2d, ierr )
616#endif
617          ENDDO
618
619       ELSE
620#if defined( __netcdf )
621!
622!--       Parallel output in netCDF4/HDF5 format.
623!--       Do not output redundant ghost point data except for the
624!--       boundaries of the total domain.
625          IF ( nxr == nx  .AND.  nyn /= ny )  THEN
626             nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_do3d(av,if),  &
627                               local_pf(nxl:nxr+1,nys:nyn,nzb_do:nzt_do),    &
628                start = (/ nxl+1, nys+1, nzb_do+1, do3d_time_count(av) /),  &
629                count = (/ nxr-nxl+2, nyn-nys+1, nzt_do-nzb_do+1, 1 /) )
630          ELSEIF ( nxr /= nx  .AND.  nyn == ny )  THEN
631             nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_do3d(av,if),  &
632                               local_pf(nxl:nxr,nys:nyn+1,nzb_do:nzt_do),    &
633                start = (/ nxl+1, nys+1, nzb_do+1, do3d_time_count(av) /),  &
634                count = (/ nxr-nxl+1, nyn-nys+2, nzt_do-nzb_do+1, 1 /) )
635          ELSEIF ( nxr == nx  .AND.  nyn == ny )  THEN
636             nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_do3d(av,if),  &
637                             local_pf(nxl:nxr+1,nys:nyn+1,nzb_do:nzt_do  ),  &
638                start = (/ nxl+1, nys+1, nzb_do+1, do3d_time_count(av) /),  &
639                count = (/ nxr-nxl+2, nyn-nys+2, nzt_do-nzb_do+1, 1 /) )
640          ELSE
641             nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_do3d(av,if),  &
642                                 local_pf(nxl:nxr,nys:nyn,nzb_do:nzt_do),    &
643                start = (/ nxl+1, nys+1, nzb_do+1, do3d_time_count(av) /),  &
644                count = (/ nxr-nxl+1, nyn-nys+1, nzt_do-nzb_do+1, 1 /) )
645          ENDIF
646          CALL handle_netcdf_error( 'data_output_3d', 386 )
647#endif
648       ENDIF
649#else
650#if defined( __netcdf )
651       nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_do3d(av,if),        &
652                         local_pf(nxl:nxr+1,nys:nyn+1,nzb_do:nzt_do),        &
653                         start = (/ 1, 1, 1, do3d_time_count(av) /),     &
654                         count = (/ nx+2, ny+2, nzt_do-nzb_do+1, 1 /) )
655       CALL handle_netcdf_error( 'data_output_3d', 446 )
656#endif
657#endif
658
659       if = if + 1
660
661!
662!--    Deallocate temporary array
663       DEALLOCATE ( local_pf )
664
665    ENDDO
666
667    CALL cpu_log( log_point(14), 'data_output_3d', 'stop' )
668
669!
670!-- Formats.
6713300 FORMAT ('variable ',I4,'  file=',A,'  filetype=unformatted  skip=',I12/   &
672             'label = ',A,A)
673
674 END SUBROUTINE data_output_3d
Note: See TracBrowser for help on using the repository browser.