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

Last change on this file since 1321 was 1321, checked in by raasch, 10 years ago

last commit documented

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