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

Last change on this file since 1977 was 1977, checked in by maronga, 8 years ago

last commit documented

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