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

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