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

Last change on this file since 1995 was 1981, checked in by suehring, 8 years ago

last commit documented

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