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

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

further modularization of land surface model (2D/3D output and restart data)

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