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

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