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

Last change on this file since 1822 was 1822, checked in by hoffmann, 8 years ago

changes in LPM and bulk cloud microphysics

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