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

Last change on this file since 2011 was 2011, checked in by kanani, 8 years ago

changes related to steering and formating of urban surface model

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