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

Last change on this file since 1973 was 1973, checked in by maronga, 8 years ago

last commit documented

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