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

Last change on this file since 1683 was 1683, checked in by knoop, 9 years ago

last commit documented

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