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

Last change on this file since 2001 was 2001, checked in by knoop, 8 years ago

last commit documented

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