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

Last change on this file since 1054 was 1054, checked in by hoffmann, 11 years ago

last commit documented

  • Property svn:keywords set to Id
File size: 20.7 KB
Line 
1 SUBROUTINE data_output_3d( av )
2
3!--------------------------------------------------------------------------------!
4! This file is part of PALM.
5!
6! PALM is free software: you can redistribute it and/or modify it under the terms
7! of the GNU General Public License as published by the Free Software Foundation,
8! either version 3 of the License, or (at your option) any later 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-2012  Leibniz University Hannover
18!--------------------------------------------------------------------------------!
19!
20! Current revisions:
21! -----------------
22!
23! Former revisions:
24! -----------------
25! $Id: data_output_3d.f90 1054 2012-11-13 17:30:09Z hoffmann $
26!
27! 1053 2012-11-13 17:11:03Z hoffmann
28! +nr, qr, prr, qc and averaged quantities
29!
30! 1036 2012-10-22 13:43:42Z raasch
31! code put under GPL (PALM 3.9)
32!
33! 1031 2012-10-19 14:35:30Z raasch
34! netCDF4 without parallel file support implemented
35!
36! 1007 2012-09-19 14:30:36Z franke
37! Bugfix: missing calculation of ql_vp added
38!
39! 790 2011-11-29 03:11:20Z raasch
40! bugfix: calculation of 'pr' must depend on the particle weighting factor,
41! nzt+1 replaced by nz_do3d for 'pr'
42!
43! 771 2011-10-27 10:56:21Z heinze
44! +lpt
45!
46! 759 2011-09-15 13:58:31Z raasch
47! Splitting of parallel I/O
48!
49! 727 2011-04-20 20:05:25Z suehring
50! Exchange ghost layers also for p_av.
51!
52! 725 2011-04-11 09:37:01Z suehring
53! Exchange ghost layers for p regardless of used pressure solver (except SOR).
54!
55! 691 2011-03-04 08:45:30Z maronga
56! Replaced simulated_time by time_since_reference_point
57!
58! 673 2011-01-18 16:19:48Z suehring
59! When using Multigrid or SOR solver an additional CALL exchange_horiz is
60! is needed for pressure output.
61!
62! 667 2010-12-23 12:06:00Z suehring/gryschka
63! nxl-1, nxr+1, nys-1, nyn+1 replaced by nxlg, nxrg, nysg, nyng in loops and
64! allocation of arrays.  Calls of exchange_horiz are modified.
65! Skip-value skip_do_avs changed to a dynamic adaption of ghost points.
66!
67! 646 2010-12-15 13:03:52Z raasch
68! bugfix: missing define statements for netcdf added
69!
70! 493 2010-03-01 08:30:24Z raasch
71! netCDF4 support (parallel output)
72!
73! 355 2009-07-17 01:03:01Z letzel
74! simulated_time in netCDF output replaced by time_since_reference_point.
75! Output of netCDF messages with aid of message handling routine.
76! Output of messages replaced by message handling routine.
77! Bugfix: to_be_resorted => s_av for time-averaged scalars
78!
79! 96 2007-06-04 08:07:41Z raasch
80! Output of density and salinity
81!
82! 75 2007-03-22 09:54:05Z raasch
83! 2nd+3rd argument removed from exchange horiz
84!
85! RCS Log replace by Id keyword, revision history cleaned up
86!
87! Revision 1.3  2006/06/02 15:18:59  raasch
88! +argument "found", -argument grid in call of routine user_data_output_3d
89!
90! Revision 1.2  2006/02/23 10:23:07  raasch
91! Former subroutine plot_3d renamed data_output_3d, pl.. renamed do..,
92! .._anz renamed .._n,
93! output extended to (almost) all quantities, output of user-defined quantities
94!
95! Revision 1.1  1997/09/03 06:29:36  raasch
96! Initial revision
97!
98!
99! Description:
100! ------------
101! Output of the 3D-arrays in netCDF and/or AVS format.
102!------------------------------------------------------------------------------!
103
104    USE array_kind
105    USE arrays_3d
106    USE averaging
107    USE cloud_parameters
108    USE control_parameters
109    USE cpulog
110    USE indices
111    USE interfaces
112    USE netcdf_control
113    USE particle_attributes
114    USE pegrid
115
116    IMPLICIT NONE
117
118    CHARACTER (LEN=9) ::  simulated_time_mod
119
120    INTEGER           ::  av, i, if, j, k, n, pos, prec, psi
121
122    LOGICAL           ::  found, resorted
123
124    REAL              ::  mean_r, s_r3, s_r4
125
126    REAL(spk), DIMENSION(:,:,:), ALLOCATABLE  ::  local_pf
127
128    REAL, DIMENSION(:,:,:), POINTER ::  to_be_resorted
129
130!
131!-- Return, if nothing to output
132    IF ( do3d_no(av) == 0 )  RETURN
133
134    CALL cpu_log (log_point(14),'data_output_3d','start')
135
136!
137!-- Open output file.
138!-- Also creates coordinate and fld-file for AVS.
139!-- For classic or 64bit netCDF output or output of other (old) data formats,
140!-- for a run on more than one PE, each PE opens its own file and
141!-- writes the data of its subdomain in binary format (regardless of the format
142!-- the user has requested). After the run, these files are combined to one
143!-- file by combine_plot_fields in the format requested by the user (netcdf
144!-- and/or avs).
145!-- For netCDF4/HDF5 output, data is written in parallel into one file.
146    IF ( netcdf_output )  THEN
147       IF ( netcdf_data_format < 5 )  THEN
148          CALL check_open( 30 )
149          IF ( myid == 0 )  CALL check_open( 106+av*10 )
150       ELSE
151          CALL check_open( 106+av*10 )
152       ENDIF
153    ELSE
154       IF ( avs_output  .OR.  ( numprocs > 1 ) )  CALL check_open( 30 )
155    ENDIF
156
157!
158!-- Allocate a temporary array with the desired output dimensions.
159    ALLOCATE( local_pf(nxlg:nxrg,nysg:nyng,nzb:nz_do3d) )
160
161!
162!-- Update the netCDF time axis
163#if defined( __netcdf )
164    IF ( myid == 0  .OR.  netcdf_data_format > 4 )  THEN
165       do3d_time_count(av) = do3d_time_count(av) + 1
166       IF ( netcdf_output )  THEN
167          nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_time_3d(av), &
168                                  (/ time_since_reference_point /),  &
169                                  start = (/ do3d_time_count(av) /), &
170                                  count = (/ 1 /) )
171          CALL handle_netcdf_error( 'data_output_3d', 376 )
172       ENDIF
173    ENDIF
174#endif
175
176!
177!-- Loop over all variables to be written.
178    if = 1
179
180    DO  WHILE ( do3d(av,if)(1:1) /= ' ' )
181!
182!--    Set the precision for data compression.
183       IF ( do3d_compress )  THEN
184          DO  i = 1, 100
185             IF ( plot_3d_precision(i)%variable == do3d(av,if) )  THEN
186                prec = plot_3d_precision(i)%precision
187                EXIT
188             ENDIF
189          ENDDO
190       ENDIF
191
192!
193!--    Store the array chosen on the temporary array.
194       resorted = .FALSE.
195       SELECT CASE ( TRIM( do3d(av,if) ) )
196
197          CASE ( 'e' )
198             IF ( av == 0 )  THEN
199                to_be_resorted => e
200             ELSE
201                to_be_resorted => e_av
202             ENDIF
203
204          CASE ( 'lpt' )
205             IF ( av == 0 )  THEN
206                to_be_resorted => pt
207             ELSE
208                to_be_resorted => lpt_av
209             ENDIF
210
211          CASE ( 'nr' )
212             IF ( av == 0 )  THEN
213                to_be_resorted => nr
214             ELSE
215                to_be_resorted => nr_av
216             ENDIF
217
218          CASE ( 'p' )
219             IF ( av == 0 )  THEN
220                IF ( psolver /= 'sor' )  CALL exchange_horiz( p, nbgp )
221                to_be_resorted => p
222             ELSE
223                IF ( psolver /= 'sor' )  CALL exchange_horiz( p_av, nbgp )
224                to_be_resorted => p_av
225             ENDIF
226
227          CASE ( 'pc' )  ! particle concentration (requires ghostpoint exchange)
228             IF ( av == 0 )  THEN
229                tend = prt_count
230                CALL exchange_horiz( tend, nbgp )
231                DO  i = nxlg, nxrg
232                   DO  j = nysg, nyng
233                      DO  k = nzb, nz_do3d
234                         local_pf(i,j,k) = tend(k,j,i)
235                      ENDDO
236                   ENDDO
237                ENDDO
238                resorted = .TRUE.
239             ELSE
240                CALL exchange_horiz( pc_av, nbgp )
241                to_be_resorted => pc_av
242             ENDIF
243
244          CASE ( 'pr' )  ! mean particle radius
245             IF ( av == 0 )  THEN
246                DO  i = nxl, nxr
247                   DO  j = nys, nyn
248                      DO  k = nzb, nz_do3d
249                         psi = prt_start_index(k,j,i)
250                         s_r3 = 0.0
251                         s_r4 = 0.0
252                         DO  n = psi, psi+prt_count(k,j,i)-1
253                         s_r3 = s_r3 + particles(n)%radius**3 * &
254                                       particles(n)%weight_factor
255                         s_r4 = s_r4 + particles(n)%radius**4 * &
256                                       particles(n)%weight_factor
257                         ENDDO
258                         IF ( s_r3 /= 0.0 )  THEN
259                            mean_r = s_r4 / s_r3
260                         ELSE
261                            mean_r = 0.0
262                         ENDIF
263                         tend(k,j,i) = mean_r
264                      ENDDO
265                   ENDDO
266                ENDDO
267                CALL exchange_horiz( tend, nbgp )
268                DO  i = nxlg, nxrg
269                   DO  j = nysg, nyng
270                      DO  k = nzb, nz_do3d
271                         local_pf(i,j,k) = tend(k,j,i)
272                      ENDDO
273                   ENDDO
274                ENDDO
275                resorted = .TRUE.
276             ELSE
277                CALL exchange_horiz( pr_av, nbgp )
278                to_be_resorted => pr_av
279             ENDIF
280
281          CASE ( 'prr' )
282             IF ( av == 0 )  THEN
283                CALL exchange_horiz( prr, nbgp )
284                DO  i = nxlg, nxrg
285                   DO  j = nysg, nyng
286                      DO  k = nzb, nzt+1
287                         local_pf(i,j,k) = prr(k,j,i)
288                      ENDDO
289                   ENDDO
290                ENDDO
291             ELSE
292                CALL exchange_horiz( prr_av, nbgp )
293                DO  i = nxlg, nxrg
294                   DO  j = nysg, nyng
295                      DO  k = nzb, nzt+1
296                         local_pf(i,j,k) = prr_av(k,j,i)
297                      ENDDO
298                   ENDDO
299                ENDDO
300             ENDIF
301             resorted = .TRUE.
302
303          CASE ( 'pt' )
304             IF ( av == 0 )  THEN
305                IF ( .NOT. cloud_physics ) THEN
306                   to_be_resorted => pt
307                ELSE
308                   DO  i = nxlg, nxrg
309                      DO  j = nysg, nyng
310                         DO  k = nzb, nz_do3d
311                            local_pf(i,j,k) = pt(k,j,i) + l_d_cp *    &
312                                                          pt_d_t(k) * &
313                                                          ql(k,j,i)
314                         ENDDO
315                      ENDDO
316                   ENDDO
317                   resorted = .TRUE.
318                ENDIF
319             ELSE
320                to_be_resorted => pt_av
321             ENDIF
322
323          CASE ( 'q' )
324             IF ( av == 0 )  THEN
325                to_be_resorted => q
326             ELSE
327                to_be_resorted => q_av
328             ENDIF
329
330          CASE ( 'qc' )
331             IF ( av == 0 )  THEN
332                to_be_resorted => ql
333             ELSE
334                to_be_resorted => ql_av
335             ENDIF
336
337          CASE ( 'ql' )
338             IF ( av == 0 )  THEN
339                IF ( icloud_scheme == 0 )  THEN
340                   DO  i = nxlg, nxrg
341                      DO  j = nysg, nyng
342                         DO  k = nzb, nz_do3d
343                            local_pf(i,j,k) = ql(k,j,i) + qr(k,j,i)
344                         ENDDO
345                      ENDDO
346                   ENDDO
347                   resorted = .TRUE.
348                ELSE
349                   to_be_resorted => ql
350                ENDIF
351             ELSE
352                IF ( icloud_scheme == 0 )  THEN
353                   DO  i = nxlg, nxrg
354                      DO  j = nysg, nyng
355                         DO  k = nzb, nz_do3d
356                            local_pf(i,j,k) = ql_av(k,j,i) + qr_av(k,j,i)
357                         ENDDO
358                      ENDDO
359                   ENDDO
360                   resorted = .TRUE.
361                ELSE
362                   to_be_resorted => ql_av
363                ENDIF
364             ENDIF
365
366          CASE ( 'ql_c' )
367             IF ( av == 0 )  THEN
368                to_be_resorted => ql_c
369             ELSE
370                to_be_resorted => ql_c_av
371             ENDIF
372
373          CASE ( 'ql_v' )
374             IF ( av == 0 )  THEN
375                to_be_resorted => ql_v
376             ELSE
377                to_be_resorted => ql_v_av
378             ENDIF
379
380          CASE ( 'ql_vp' )
381             IF ( av == 0 )  THEN
382                DO  i = nxl, nxr
383                   DO  j = nys, nyn
384                      DO  k = nzb, nz_do3d
385                         psi = prt_start_index(k,j,i)
386                         DO  n = psi, psi+prt_count(k,j,i)-1
387                            tend(k,j,i) = tend(k,j,i) + &
388                                          particles(n)%weight_factor / &
389                                          prt_count(k,j,i)
390                         ENDDO
391                      ENDDO
392                   ENDDO
393                ENDDO
394                CALL exchange_horiz( tend, nbgp )
395                DO  i = nxlg, nxrg
396                   DO  j = nysg, nyng
397                      DO  k = nzb, nz_do3d
398                         local_pf(i,j,k) = tend(k,j,i)
399                      ENDDO
400                   ENDDO
401                ENDDO
402                resorted = .TRUE.
403             ELSE
404                CALL exchange_horiz( ql_vp_av, nbgp )
405                to_be_resorted => ql_vp_av
406             ENDIF
407
408          CASE ( 'qr' )
409             IF ( av == 0 )  THEN
410                to_be_resorted => qr
411             ELSE
412                to_be_resorted => qr_av
413             ENDIF
414
415          CASE ( 'qv' )
416             IF ( av == 0 )  THEN
417                DO  i = nxlg, nxrg
418                   DO  j = nysg, nyng
419                      DO  k = nzb, nz_do3d
420                         local_pf(i,j,k) = q(k,j,i) - ql(k,j,i)
421                      ENDDO
422                   ENDDO
423                ENDDO
424                resorted = .TRUE.
425             ELSE
426                to_be_resorted => qv_av
427             ENDIF
428
429          CASE ( 'rho' )
430             IF ( av == 0 )  THEN
431                to_be_resorted => rho
432             ELSE
433                to_be_resorted => rho_av
434             ENDIF
435
436          CASE ( 's' )
437             IF ( av == 0 )  THEN
438                to_be_resorted => q
439             ELSE
440                to_be_resorted => s_av
441             ENDIF
442
443          CASE ( 'sa' )
444             IF ( av == 0 )  THEN
445                to_be_resorted => sa
446             ELSE
447                to_be_resorted => sa_av
448             ENDIF
449
450          CASE ( 'u' )
451             IF ( av == 0 )  THEN
452                to_be_resorted => u
453             ELSE
454                to_be_resorted => u_av
455             ENDIF
456
457          CASE ( 'v' )
458             IF ( av == 0 )  THEN
459                to_be_resorted => v
460             ELSE
461                to_be_resorted => v_av
462             ENDIF
463
464          CASE ( 'vpt' )
465             IF ( av == 0 )  THEN
466                to_be_resorted => vpt
467             ELSE
468                to_be_resorted => vpt_av
469             ENDIF
470
471          CASE ( 'w' )
472             IF ( av == 0 )  THEN
473                to_be_resorted => w
474             ELSE
475                to_be_resorted => w_av
476             ENDIF
477
478          CASE DEFAULT
479!
480!--          User defined quantity
481             CALL user_data_output_3d( av, do3d(av,if), found, local_pf, &
482                                       nz_do3d )
483             resorted = .TRUE.
484
485             IF ( .NOT. found )  THEN
486                message_string =  'no output available for: ' //   &
487                                  TRIM( do3d(av,if) )
488                CALL message( 'data_output_3d', 'PA0182', 0, 0, 0, 6, 0 )
489             ENDIF
490
491       END SELECT
492
493!
494!--    Resort the array to be output, if not done above
495       IF ( .NOT. resorted )  THEN
496          DO  i = nxlg, nxrg
497             DO  j = nysg, nyng
498                DO  k = nzb, nz_do3d
499                   local_pf(i,j,k) = to_be_resorted(k,j,i)
500                ENDDO
501             ENDDO
502          ENDDO
503       ENDIF
504
505!
506!--    Output of the volume data information for the AVS-FLD-file.
507       do3d_avs_n = do3d_avs_n + 1
508       IF ( myid == 0  .AND.  avs_output )  THEN
509!
510!--       AVS-labels must not contain any colons. Hence they must be removed
511!--       from the time character string.
512          simulated_time_mod = simulated_time_chr
513          DO  WHILE ( SCAN( simulated_time_mod, ':' ) /= 0 )
514             pos = SCAN( simulated_time_mod, ':' )
515             simulated_time_mod(pos:pos) = '/'
516          ENDDO
517
518          IF ( av == 0 )  THEN
519             WRITE ( 33, 3300 )  do3d_avs_n, TRIM( avs_data_file ), &
520                                 skip_do_avs, TRIM( do3d(av,if) ), &
521                                 TRIM( simulated_time_mod )
522          ELSE
523             WRITE ( 33, 3300 )  do3d_avs_n, TRIM( avs_data_file ), &
524                                 skip_do_avs, TRIM( do3d(av,if) ) // &
525                                 ' averaged', TRIM( simulated_time_mod )
526          ENDIF
527!
528!--       Determine the Skip-value for the next array. Record end and start
529!--       require 4 byte each.
530          skip_do_avs = skip_do_avs + ( ( ( nx+2*nbgp ) * ( ny+2*nbgp ) * &
531                                          ( nz_do3d+1 ) ) * 4 + 8 )
532       ENDIF
533
534!
535!--    Output of the 3D-array. (compressed/uncompressed)
536       IF ( do3d_compress )  THEN
537!
538!--       Compression, output of compression information on FLD-file and output
539!--       of compressed data.
540          CALL write_compressed( local_pf, 30, 33, myid, nxl, nxr, nyn, nys, &
541                                 nzb, nz_do3d, prec, nbgp )
542       ELSE
543!
544!--       Uncompressed output.
545#if defined( __parallel )
546          IF ( netcdf_output )  THEN
547             IF ( netcdf_data_format < 5 )  THEN
548!
549!--             Non-parallel netCDF output. Data is output in parallel in
550!--             FORTRAN binary format here, and later collected into one file by
551!--             combine_plot_fields
552                IF ( myid == 0 )  THEN
553                   WRITE ( 30 )  time_since_reference_point,                   &
554                                 do3d_time_count(av), av
555                ENDIF
556                DO  i = 0, io_blocks-1
557                   IF ( i == io_group )  THEN
558                      WRITE ( 30 )  nxlg, nxrg, nysg, nyng, nzb, nz_do3d
559                      WRITE ( 30 )  local_pf
560                   ENDIF
561#if defined( __parallel )
562                   CALL MPI_BARRIER( comm2d, ierr )
563#endif
564                ENDDO
565
566             ELSE
567#if defined( __netcdf )
568!
569!--             Parallel output in netCDF4/HDF5 format.
570!--             Do not output redundant ghost point data except for the
571!--             boundaries of the total domain.
572                IF ( nxr == nx  .AND.  nyn /= ny )  THEN
573                   nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_do3d(av,if), &
574                                  local_pf(nxl:nxrg,nys:nyn,nzb:nz_do3d),    &
575                      start = (/ nxl+1, nys+1, nzb+1, do3d_time_count(av) /), &
576                      count = (/ nxr-nxl+1+nbgp, nyn-nys+1, nz_do3d-nzb+1, 1 /) )
577                ELSEIF ( nxr /= nx  .AND.  nyn == ny )  THEN
578                   nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_do3d(av,if), &
579                                  local_pf(nxl:nxr,nys:nyng,nzb:nz_do3d),    &
580                      start = (/ nxl+1, nys+1, nzb+1, do3d_time_count(av) /), &
581                      count = (/ nxr-nxl+1, nyn-nys+1+nbgp, nz_do3d-nzb+1, 1 /) )
582                ELSEIF ( nxr == nx  .AND.  nyn == ny )  THEN
583                   nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_do3d(av,if), &
584                                  local_pf(nxl:nxrg,nys:nyng,nzb:nz_do3d),  &
585                      start = (/ nxl+1, nys+1, nzb+1, do3d_time_count(av) /), &
586                      count = (/ nxr-nxl+1+nbgp, nyn-nys+1+nbgp, nz_do3d-nzb+1, 1 /) )
587                ELSE
588                   nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_do3d(av,if), &
589                                  local_pf(nxl:nxr,nys:nyn,nzb:nz_do3d),      &
590                      start = (/ nxl+1, nys+1, nzb+1, do3d_time_count(av) /), &
591                      count = (/ nxr-nxl+1, nyn-nys+1, nz_do3d-nzb+1, 1 /) )
592                ENDIF
593                CALL handle_netcdf_error( 'data_output_3d', 386 )
594#endif
595             ENDIF
596          ENDIF
597#else
598          IF ( avs_output )  THEN
599             WRITE ( 30 )  local_pf(nxl:nxr+1,nys:nyn+1,:)
600          ENDIF
601#if defined( __netcdf )
602          IF ( netcdf_output )  THEN
603
604             nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_do3d(av,if),    &
605                               local_pf(nxl:nxr+1,nys:nyn+1,nzb:nz_do3d),  &
606                               start = (/ 1, 1, 1, do3d_time_count(av) /), &
607                               count = (/ nx+2, ny+2, nz_do3d-nzb+1, 1 /) )
608             CALL handle_netcdf_error( 'data_output_3d', 446 )
609
610          ENDIF
611#endif
612#endif
613       ENDIF
614
615       if = if + 1
616
617    ENDDO
618
619!
620!-- Deallocate temporary array.
621    DEALLOCATE( local_pf )
622
623
624    CALL cpu_log (log_point(14),'data_output_3d','stop','nobarrier')
625
626!
627!-- Formats.
6283300 FORMAT ('variable ',I4,'  file=',A,'  filetype=unformatted  skip=',I12/ &
629             'label = ',A,A)
630
631 END SUBROUTINE data_output_3d
Note: See TracBrowser for help on using the repository browser.