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

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

preprocessor directives using machine dependent flags (lc, ibm, etc.) mostly removed from the code

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