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

Last change on this file since 1746 was 1746, checked in by gronemeier, 8 years ago

last commit documented

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