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

Last change on this file since 1817 was 1809, checked in by raasch, 8 years ago

last commit documented

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