source: palm/trunk/SOURCE/data_output_mask.f90 @ 4417

Last change on this file since 4417 was 4377, checked in by gronemeier, 4 years ago

bugfix: set fill value for output according to wall_flags_total_0 for non-terrain following output

  • Property svn:keywords set to Id
File size: 32.3 KB
Line 
1!> @file data_output_mask.f90
2!------------------------------------------------------------------------------!
3! This file is part of the PALM model system.
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-2020 Leibniz Universitaet Hannover
18!------------------------------------------------------------------------------!
19!
20! Current revisions:
21! -----------------
22!
23!
24! Former revisions:
25! -----------------
26! $Id: data_output_mask.f90 4377 2020-01-15 11:10:51Z monakurppa $
27! bugfix: set fill value for output according to wall_flags_total_0 for
28!         non-terrain following output
29!
30! 4360 2020-01-07 11:25:50Z suehring
31! Introduction of wall_flags_total_0, which currently sets bits based on static
32! topography information used in wall_flags_static_0
33!
34! 4331 2019-12-10 18:25:02Z suehring
35! Formatting adjustment
36!
37! 4329 2019-12-10 15:46:36Z motisi
38! Renamed wall_flags_0 to wall_flags_static_0
39!
40! 4246 2019-09-30 09:27:52Z pavelkrc
41! Corrected "Former revisions" section
42!
43! 4168 2019-08-16 13:50:17Z suehring
44! Remove variable grid
45!
46! 4167 2019-08-16 11:01:48Z suehring
47! Changed behaviour of masked output over surface to follow terrain and ignore
48! buildings (J.Resler, T.Gronemeier)
49!
50! 4069 2019-07-01 14:05:51Z Giersch
51! Masked output running index mid has been introduced as a local variable to
52! avoid runtime error (Loop variable has been modified) in time_integration
53!
54! 4039 2019-06-18 10:32:41Z suehring
55! Modularize diagnostic output
56!
57! 3994 2019-05-22 18:08:09Z suehring
58! output of turbulence intensity added
59!
60! 3665 2019-01-10 08:28:24Z raasch
61! unused variables removed
62!
63! 3655 2019-01-07 16:51:22Z knoop
64! Fix output time levels (use time_since_reference_point)
65!
66! 410 2009-12-04 17:05:40Z letzel
67! Initial version
68!
69! Description:
70! ------------
71!> Masked data output in netCDF format for current mask (current value of mid).
72!------------------------------------------------------------------------------!
73 SUBROUTINE data_output_mask( av, mid )
74
75
76
77#if defined( __netcdf )
78    USE arrays_3d,                                                             &
79        ONLY:  e, nc, nr, p, pt, q, qc, ql, ql_c, ql_v, qr, rho_ocean, s, sa,  &
80               tend, u, v, vpt, w, d_exner
81
82    USE averaging,                                                             &
83        ONLY:  e_av, lpt_av, nc_av, nr_av, p_av, pc_av, pr_av, pt_av, q_av,    &
84               qc_av, ql_av, ql_c_av, ql_v_av, ql_vp_av, qv_av, qr_av,         &
85               rho_ocean_av, s_av, sa_av, u_av, v_av, vpt_av, w_av
86
87    USE basic_constants_and_equations_mod,                                     &
88        ONLY:  lv_d_cp
89
90    USE chemistry_model_mod,                                                   &
91        ONLY:  chem_data_output_mask
92
93    USE control_parameters,                                                    &
94        ONLY:  air_chemistry, domask, domask_no, domask_time_count, mask_i,    &
95               mask_j, mask_k, mask_size, mask_size_l, mask_start_l,           &
96               mask_surface,                                                   &
97               max_masks, message_string, nz_do3d, salsa,                      &
98               time_since_reference_point
99
100    USE diagnostic_output_quantities_mod,                                      &
101        ONLY:  doq_output_mask
102
103    USE cpulog,                                                                &
104        ONLY:  cpu_log, log_point
105
106    USE indices,                                                               &
107        ONLY:  nbgp, nxl, nxr, nyn, nys, nzb, nzt, wall_flags_total_0
108
109    USE kinds
110
111    USE bulk_cloud_model_mod,                                                  &
112        ONLY:  bulk_cloud_model
113
114    USE NETCDF
115
116    USE netcdf_interface,                                                      &
117        ONLY:  fill_value, id_set_mask, id_var_domask, id_var_time_mask,       &
118               nc_stat, netcdf_data_format, netcdf_handle_error
119
120    USE particle_attributes,                                                   &
121        ONLY:  grid_particles, number_of_particles, particles,                 &
122               particle_advection_start, prt_count
123
124    USE pegrid
125
126    USE radiation_model_mod,                                                   &
127        ONLY:  radiation, radiation_data_output_mask
128
129    USE salsa_mod,                                                             &
130        ONLY:  salsa_data_output_mask
131
132
133    IMPLICIT NONE
134
135    INTEGER(iwp) ::  av                      !< flag for (non-)average output
136    INTEGER(iwp) ::  ngp                     !< number of grid points of an output slice
137    INTEGER(iwp) ::  flag_nr                 !< number of masking flag
138    INTEGER(iwp) ::  i                       !< loop index
139    INTEGER(iwp) ::  ivar                    !< variable index
140    INTEGER(iwp) ::  j                       !< loop index
141    INTEGER(iwp) ::  k                       !< loop index
142    INTEGER(iwp) ::  im                      !< loop index for masked variables
143    INTEGER(iwp) ::  jm                      !< loop index for masked variables
144    INTEGER(iwp) ::  kk                      !< vertical index
145    INTEGER(iwp) ::  mid                     !< masked output running index
146    INTEGER(iwp) ::  n                       !< loop index
147    INTEGER(iwp) ::  netcdf_data_format_save !< value of netcdf_data_format
148    INTEGER(iwp) ::  sender                  !< PE id of sending PE
149    INTEGER(iwp) ::  ktt                     !< k index of highest terrain surface
150    INTEGER(iwp) ::  ind(6)                  !< index limits (lower/upper bounds) of array 'local_2d'
151
152    LOGICAL ::  found      !< true if output variable was found
153    LOGICAL ::  resorted   !< true if variable is resorted
154
155    REAL(wp) ::  mean_r    !< mean particle radius
156    REAL(wp) ::  s_r2      !< sum( particle-radius**2 )
157    REAL(wp) ::  s_r3      !< sum( particle-radius**3 )
158
159    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  local_pf    !< output array
160#if defined( __parallel )
161    REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  total_pf    !< collected output array
162#endif
163    REAL(wp), DIMENSION(:,:,:), POINTER ::  to_be_resorted  !< points to array which shall be output
164
165!
166!-- Return, if nothing to output
167    IF ( domask_no(mid,av) == 0 )  RETURN
168
169    CALL cpu_log (log_point(49),'data_output_mask','start')
170
171!
172!-- Parallel netcdf output is not tested so far for masked data, hence
173!-- netcdf_data_format is switched back to non-paralell output.
174    netcdf_data_format_save = netcdf_data_format
175    IF ( netcdf_data_format == 5 ) netcdf_data_format = 3
176    IF ( netcdf_data_format == 6 ) netcdf_data_format = 4
177
178!
179!-- Open output file.
180    IF ( myid == 0  .OR.  netcdf_data_format > 4 )  THEN
181       CALL check_open( 200+mid+av*max_masks )
182    ENDIF
183
184!
185!-- Allocate total and local output arrays.
186#if defined( __parallel )
187    IF ( myid == 0 )  THEN
188       ALLOCATE( total_pf(mask_size(mid,1),mask_size(mid,2),mask_size(mid,3)) )
189    ENDIF
190#endif
191    ALLOCATE( local_pf(mask_size_l(mid,1),mask_size_l(mid,2), &
192                       mask_size_l(mid,3)) )
193
194!
195!-- Update the netCDF time axis.
196    domask_time_count(mid,av) = domask_time_count(mid,av) + 1
197    IF ( myid == 0  .OR.  netcdf_data_format > 4 )  THEN
198       nc_stat = NF90_PUT_VAR( id_set_mask(mid,av), id_var_time_mask(mid,av), &
199                               (/ time_since_reference_point /),              &
200                               start = (/ domask_time_count(mid,av) /),       &
201                               count = (/ 1 /) )
202       CALL netcdf_handle_error( 'data_output_mask', 460 )
203    ENDIF
204
205!
206!-- Loop over all variables to be written.
207    ivar = 1
208
209    DO  WHILE ( domask(mid,av,ivar)(1:1) /= ' ' )
210!
211!--    Reallocate local_pf on PE 0 since its shape changes during MPI exchange
212       IF ( netcdf_data_format < 5   .AND.  myid == 0  .AND.  ivar > 1 )  THEN
213          DEALLOCATE( local_pf )
214          ALLOCATE( local_pf(mask_size_l(mid,1),mask_size_l(mid,2), &
215                             mask_size_l(mid,3)) )
216       ENDIF
217!
218!--    Set masking flag for topography for not resorted arrays
219       flag_nr = 0
220!
221!--    Store the variable chosen.
222       resorted = .FALSE.
223       SELECT CASE ( TRIM( domask(mid,av,ivar) ) )
224
225          CASE ( 'e' )
226             IF ( av == 0 )  THEN
227                to_be_resorted => e
228             ELSE
229                to_be_resorted => e_av
230             ENDIF
231
232          CASE ( 'thetal' )
233             IF ( av == 0 )  THEN
234                to_be_resorted => pt
235             ELSE
236                to_be_resorted => lpt_av
237             ENDIF
238
239          CASE ( 'nc' )
240             IF ( av == 0 )  THEN
241                to_be_resorted => nc
242             ELSE
243                to_be_resorted => nc_av
244             ENDIF
245
246          CASE ( 'nr' )
247             IF ( av == 0 )  THEN
248                to_be_resorted => nr
249             ELSE
250                to_be_resorted => nr_av
251             ENDIF
252
253          CASE ( 'p' )
254             IF ( av == 0 )  THEN
255                to_be_resorted => p
256             ELSE
257                to_be_resorted => p_av
258             ENDIF
259
260          CASE ( 'pc' )  ! particle concentration (requires ghostpoint exchange)
261             IF ( av == 0 )  THEN
262                tend = prt_count
263                CALL exchange_horiz( tend, nbgp )
264                IF ( .NOT. mask_surface(mid) )  THEN
265                   DO  i = 1, mask_size_l(mid,1)
266                      DO  j = 1, mask_size_l(mid,2)
267                         DO  k = 1, mask_size_l(mid,3)
268                            local_pf(i,j,k) =  tend(mask_k(mid,k), &
269                                      mask_j(mid,j),mask_i(mid,i))
270                         ENDDO
271                      ENDDO
272                   ENDDO
273                ELSE
274!
275!--                Terrain-following masked output
276                   DO  i = 1, mask_size_l(mid,1)
277                      DO  j = 1, mask_size_l(mid,2)
278!--                      Get k index of the highest terraing surface
279                         im = mask_i(mid,i)
280                         jm = mask_j(mid,j)
281                         ktt = MINLOC( MERGE( 1, 0, BTEST( wall_flags_total_0(:,jm,im), 5 )),&
282                                       DIM = 1 ) - 1
283                         DO  k = 1, mask_size_l(mid,3)
284                            kk = MIN( ktt+mask_k(mid,k), nzt+1 )
285!--                         Set value if not in building
286                            IF ( BTEST( wall_flags_total_0(kk,jm,im), 6 ) )  THEN
287                               local_pf(i,j,k) = fill_value
288                            ELSE
289                               local_pf(i,j,k) =  tend(kk,jm,im)
290                            ENDIF
291                         ENDDO
292                      ENDDO
293                   ENDDO
294                ENDIF
295                resorted = .TRUE.
296             ELSE
297                CALL exchange_horiz( pc_av, nbgp )
298                to_be_resorted => pc_av
299             ENDIF
300
301          CASE ( 'pr' )  ! mean particle radius (effective radius)
302             IF ( av == 0 )  THEN
303                IF ( time_since_reference_point >= particle_advection_start )  THEN
304                   DO  i = nxl, nxr
305                      DO  j = nys, nyn
306                         DO  k = nzb, nz_do3d
307                            number_of_particles = prt_count(k,j,i)
308                            IF (number_of_particles <= 0)  CYCLE
309                            particles => grid_particles(k,j,i)%particles(1:number_of_particles)
310                            s_r2 = 0.0_wp
311                            s_r3 = 0.0_wp
312                            DO  n = 1, number_of_particles
313                               IF ( particles(n)%particle_mask )  THEN
314                                  s_r2 = s_r2 + grid_particles(k,j,i)%particles(n)%radius**2 * &
315                                         grid_particles(k,j,i)%particles(n)%weight_factor
316                                  s_r3 = s_r3 + grid_particles(k,j,i)%particles(n)%radius**3 * &
317                                         grid_particles(k,j,i)%particles(n)%weight_factor
318                               ENDIF
319                            ENDDO
320                            IF ( s_r2 > 0.0_wp )  THEN
321                               mean_r = s_r3 / s_r2
322                            ELSE
323                               mean_r = 0.0_wp
324                            ENDIF
325                            tend(k,j,i) = mean_r
326                         ENDDO
327                      ENDDO
328                   ENDDO
329                   CALL exchange_horiz( tend, nbgp )
330                ELSE
331                   tend = 0.0_wp
332                ENDIF
333                IF ( .NOT. mask_surface(mid) )  THEN
334                   DO  i = 1, mask_size_l(mid,1)
335                      DO  j = 1, mask_size_l(mid,2)
336                         DO  k = 1, mask_size_l(mid,3)
337                            local_pf(i,j,k) =  tend(mask_k(mid,k), &
338                                      mask_j(mid,j),mask_i(mid,i))
339                         ENDDO
340                      ENDDO
341                   ENDDO
342                ELSE
343!
344!--                Terrain-following masked output
345                   DO  i = 1, mask_size_l(mid,1)
346                      DO  j = 1, mask_size_l(mid,2)
347!--                      Get k index of the highest terraing surface
348                         im = mask_i(mid,i)
349                         jm = mask_j(mid,j)
350                         ktt = MINLOC( MERGE( 1, 0, BTEST( wall_flags_total_0(:,jm,im), 5 )), &
351                                         DIM = 1 ) - 1
352                         DO  k = 1, mask_size_l(mid,3)
353                            kk = MIN( ktt+mask_k(mid,k), nzt+1 )
354!--                         Set value if not in building
355                            IF ( BTEST( wall_flags_total_0(kk,jm,im), 6 ) )  THEN
356                               local_pf(i,j,k) = fill_value
357                            ELSE
358                               local_pf(i,j,k) =  tend(kk,jm,im)
359                            ENDIF
360                         ENDDO
361                      ENDDO
362                   ENDDO
363                ENDIF
364                resorted = .TRUE.
365             ELSE
366                CALL exchange_horiz( pr_av, nbgp )
367                to_be_resorted => pr_av
368             ENDIF
369
370          CASE ( 'theta' )
371             IF ( av == 0 )  THEN
372                IF ( .NOT. bulk_cloud_model ) THEN
373                   to_be_resorted => pt
374                ELSE
375                   IF ( .NOT. mask_surface(mid) )  THEN
376                      DO  i = 1, mask_size_l(mid,1)
377                         DO  j = 1, mask_size_l(mid,2)
378                            DO  k = 1, mask_size_l(mid,3)
379                               local_pf(i,j,k) =  &
380                                  pt(mask_k(mid,k),mask_j(mid,j),mask_i(mid,i)) &
381                                  + lv_d_cp * d_exner(mask_k(mid,k)) *          &
382                                    ql(mask_k(mid,k),mask_j(mid,j),mask_i(mid,i))
383                            ENDDO
384                         ENDDO
385                      ENDDO
386                   ELSE
387!
388!--                   Terrain-following masked output
389                      DO  i = 1, mask_size_l(mid,1)
390                         DO  j = 1, mask_size_l(mid,2)
391!--                         Get k index of the highest terraing surface
392                            im = mask_i(mid,i)
393                            jm = mask_j(mid,j)
394                            ktt = MINLOC( MERGE( 1, 0, BTEST( wall_flags_total_0(:,jm,im), 5 )), &
395                                             DIM = 1 ) - 1
396                            DO  k = 1, mask_size_l(mid,3)
397                               kk = MIN( ktt+mask_k(mid,k), nzt+1 )
398!--                            Set value if not in building
399                               IF ( BTEST( wall_flags_total_0(kk,jm,im), 6 ) )  THEN
400                                  local_pf(i,j,k) = fill_value
401                               ELSE
402                                  local_pf(i,j,k) = pt(kk,jm,im) + lv_d_cp * d_exner(kk) * ql(kk,jm,im)
403                               ENDIF
404                            ENDDO
405                         ENDDO
406                      ENDDO
407                   ENDIF
408                   resorted = .TRUE.
409                ENDIF
410             ELSE
411                to_be_resorted => pt_av
412             ENDIF
413
414          CASE ( 'q' )
415             IF ( av == 0 )  THEN
416                to_be_resorted => q
417             ELSE
418                to_be_resorted => q_av
419             ENDIF
420
421          CASE ( 'qc' )
422             IF ( av == 0 )  THEN
423                to_be_resorted => qc
424             ELSE
425                to_be_resorted => qc_av
426             ENDIF
427
428          CASE ( 'ql' )
429             IF ( av == 0 )  THEN
430                to_be_resorted => ql
431             ELSE
432                to_be_resorted => ql_av
433             ENDIF
434
435          CASE ( 'ql_c' )
436             IF ( av == 0 )  THEN
437                to_be_resorted => ql_c
438             ELSE
439                to_be_resorted => ql_c_av
440             ENDIF
441
442          CASE ( 'ql_v' )
443             IF ( av == 0 )  THEN
444                to_be_resorted => ql_v
445             ELSE
446                to_be_resorted => ql_v_av
447             ENDIF
448
449          CASE ( 'ql_vp' )
450             IF ( av == 0 )  THEN
451                IF ( time_since_reference_point >= particle_advection_start )  THEN
452                   DO  i = nxl, nxr
453                      DO  j = nys, nyn
454                         DO  k = nzb, nz_do3d
455                            number_of_particles = prt_count(k,j,i)
456                            IF (number_of_particles <= 0)  CYCLE
457                            particles => grid_particles(k,j,i)%particles(1:number_of_particles)
458                            DO  n = 1, number_of_particles
459                               IF ( particles(n)%particle_mask )  THEN
460                                  tend(k,j,i) = tend(k,j,i) + &
461                                          particles(n)%weight_factor / &
462                                          prt_count(k,j,i)
463                               ENDIF
464                            ENDDO
465                         ENDDO
466                      ENDDO
467                   ENDDO
468                   CALL exchange_horiz( tend, nbgp )
469                ELSE
470                   tend = 0.0_wp
471                ENDIF
472                IF ( .NOT. mask_surface(mid) )  THEN
473                   DO  i = 1, mask_size_l(mid,1)
474                      DO  j = 1, mask_size_l(mid,2)
475                         DO  k = 1, mask_size_l(mid,3)
476                            local_pf(i,j,k) =  tend(mask_k(mid,k), &
477                                      mask_j(mid,j),mask_i(mid,i))
478                         ENDDO
479                      ENDDO
480                   ENDDO
481                ELSE
482!
483!--                Terrain-following masked output
484                   DO  i = 1, mask_size_l(mid,1)
485                      DO  j = 1, mask_size_l(mid,2)
486!--                      Get k index of the highest terraing surface
487                         im = mask_i(mid,i)
488                         jm = mask_j(mid,j)
489                         ktt = MINLOC( MERGE( 1, 0, BTEST( wall_flags_total_0(:,jm,im), 5 )), &
490                                          DIM = 1 ) - 1
491                         DO  k = 1, mask_size_l(mid,3)
492                            kk = MIN( ktt+mask_k(mid,k), nzt+1 )
493!--                         Set value if not in building
494                            IF ( BTEST( wall_flags_total_0(kk,jm,im), 6 ) )  THEN
495                               local_pf(i,j,k) = fill_value
496                            ELSE
497                               local_pf(i,j,k) = tend(kk,jm,im)
498                            ENDIF
499                         ENDDO
500                      ENDDO
501                   ENDDO
502                ENDIF
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 ( 'qv' )
510             IF ( av == 0 )  THEN
511                IF ( .NOT. mask_surface(mid) )  THEN
512                   DO  i = 1, mask_size_l(mid,1)
513                      DO  j = 1, mask_size_l(mid,2)
514                         DO  k = 1, mask_size_l(mid,3)
515                            local_pf(i,j,k) =  &
516                                 q(mask_k(mid,k),mask_j(mid,j),mask_i(mid,i)) -  &
517                                 ql(mask_k(mid,k),mask_j(mid,j),mask_i(mid,i))
518                         ENDDO
519                      ENDDO
520                   ENDDO
521                ELSE
522!
523!--                Terrain-following masked output
524                   DO  i = 1, mask_size_l(mid,1)
525                      DO  j = 1, mask_size_l(mid,2)
526!--                      Get k index of the highest terraing surface
527                         im = mask_i(mid,i)
528                         jm = mask_j(mid,j)
529                         ktt = MINLOC( MERGE( 1, 0, BTEST( wall_flags_total_0(:,jm,im), 5 )), &
530                                          DIM = 1 ) - 1
531                         DO  k = 1, mask_size_l(mid,3)
532                            kk = MIN( ktt+mask_k(mid,k), nzt+1 )
533!--                         Set value if not in building
534                            IF ( BTEST( wall_flags_total_0(kk,jm,im), 6 ) )  THEN
535                               local_pf(i,j,k) = fill_value
536                            ELSE
537                               local_pf(i,j,k) = q(kk,jm,im) - ql(kk,jm,im)
538                            ENDIF
539                         ENDDO
540                      ENDDO
541                   ENDDO
542                ENDIF
543                resorted = .TRUE.
544             ELSE
545                to_be_resorted => qv_av
546             ENDIF
547
548          CASE ( 'qr' )
549             IF ( av == 0 )  THEN
550                to_be_resorted => qr
551             ELSE
552                to_be_resorted => qr_av
553             ENDIF
554
555          CASE ( 'rho_sea_water' )
556             IF ( av == 0 )  THEN
557                to_be_resorted => rho_ocean
558             ELSE
559                to_be_resorted => rho_ocean_av
560             ENDIF
561
562          CASE ( 's' )
563             IF ( av == 0 )  THEN
564                to_be_resorted => s
565             ELSE
566                to_be_resorted => s_av
567             ENDIF
568
569          CASE ( 'sa' )
570             IF ( av == 0 )  THEN
571                to_be_resorted => sa
572             ELSE
573                to_be_resorted => sa_av
574             ENDIF
575
576          CASE ( 'u' )
577             flag_nr = 1
578             IF ( av == 0 )  THEN
579                to_be_resorted => u
580             ELSE
581                to_be_resorted => u_av
582             ENDIF
583
584          CASE ( 'v' )
585             flag_nr = 2
586             IF ( av == 0 )  THEN
587                to_be_resorted => v
588             ELSE
589                to_be_resorted => v_av
590             ENDIF
591
592          CASE ( 'thetav' )
593             IF ( av == 0 )  THEN
594                to_be_resorted => vpt
595             ELSE
596                to_be_resorted => vpt_av
597             ENDIF
598
599          CASE ( 'w' )
600             flag_nr = 3
601             IF ( av == 0 )  THEN
602                to_be_resorted => w
603             ELSE
604                to_be_resorted => w_av
605             ENDIF
606
607          CASE DEFAULT
608!
609!--          Set flag to steer output of radiation, land-surface, or user-defined
610!--          quantities
611             found = .FALSE.
612!
613!--          Radiation quantity
614             IF ( .NOT. found  .AND. radiation )  THEN
615                CALL radiation_data_output_mask(av, domask(mid,av,ivar), found,&
616                                                local_pf, mid )
617             ENDIF
618
619             IF ( .NOT. found  .AND. air_chemistry )  THEN
620                CALL chem_data_output_mask(av, domask(mid,av,ivar), found,     &
621                                           local_pf, mid )
622             ENDIF
623!
624!--          Check for diagnostic quantities
625             IF ( .NOT. found )  THEN
626                CALL doq_output_mask( av, domask(mid,av,ivar), found, local_pf,   &
627                                      mid)
628             ENDIF
629!
630!--          SALSA quantities
631             IF ( .NOT. found .AND. salsa )  THEN
632                CALL salsa_data_output_mask( av, domask(mid,av,ivar), found,   &
633                                             local_pf, mid )
634             ENDIF
635!
636!--          User defined quantity
637             IF ( .NOT. found )  THEN
638                CALL user_data_output_mask(av, domask(mid,av,ivar), found,     &
639                                           local_pf, mid )
640             ENDIF
641
642             resorted = .TRUE.
643
644             IF ( .NOT. found )  THEN
645                WRITE ( message_string, * ) 'no masked output available for: ',&
646                                            TRIM( domask(mid,av,ivar) )
647                CALL message( 'data_output_mask', 'PA0327', 0, 0, 0, 6, 0 )
648             ENDIF
649
650       END SELECT
651
652!
653!--    Resort the array to be output, if not done above
654       IF ( .NOT. resorted )  THEN
655          IF ( .NOT. mask_surface(mid) )  THEN
656!
657!--          Default masked output
658             DO  i = 1, mask_size_l(mid,1)
659                DO  j = 1, mask_size_l(mid,2)
660                   DO  k = 1, mask_size_l(mid,3)
661                      local_pf(i,j,k) = MERGE( to_be_resorted(mask_k(mid,k),  &
662                                                              mask_j(mid,j),  &
663                                                              mask_i(mid,i)), &
664                                               REAL( fill_value, KIND = wp ), &
665                                               BTEST( wall_flags_total_0(     &
666                                                              mask_k(mid,k),  &
667                                                              mask_j(mid,j),  &
668                                                              mask_i(mid,i)), &
669                                                      flag_nr ) )
670                   ENDDO
671                ENDDO
672             ENDDO
673
674          ELSE
675!
676!--          Terrain-following masked output
677             DO  i = 1, mask_size_l(mid,1)
678                DO  j = 1, mask_size_l(mid,2)
679!--                Get k index of the highest terraing surface
680                   im = mask_i(mid,i)
681                   jm = mask_j(mid,j)
682                   ktt = MINLOC( MERGE( 1, 0, BTEST( wall_flags_total_0(:,jm,im), 5 )), &
683                                    DIM = 1 ) - 1
684                   DO  k = 1, mask_size_l(mid,3)
685                      kk = MIN( ktt+mask_k(mid,k), nzt+1 )
686!--                   Set value if not in building
687                      IF ( BTEST( wall_flags_total_0(kk,jm,im), 6 ) )  THEN
688                         local_pf(i,j,k) = fill_value
689                      ELSE
690                         local_pf(i,j,k) = to_be_resorted(kk,jm,im)
691                      ENDIF
692                   ENDDO
693                ENDDO
694             ENDDO
695
696          ENDIF
697       ENDIF
698
699!
700!--    I/O block. I/O methods are implemented
701!--    (1) for parallel execution
702!--     a. with netCDF 4 parallel I/O-enabled library
703!--     b. with netCDF 3 library
704!--    (2) for serial execution.
705!--    The choice of method depends on the correct setting of preprocessor
706!--    directives __parallel and __netcdf4_parallel as well as on the parameter
707!--    netcdf_data_format.
708#if defined( __parallel )
709#if defined( __netcdf4_parallel )
710       IF ( netcdf_data_format > 4 )  THEN
711!
712!--       (1) a. Parallel I/O using netCDF 4 (not yet tested)
713          nc_stat = NF90_PUT_VAR( id_set_mask(mid,av),                         &
714               id_var_domask(mid,av,ivar), local_pf,                           &
715               start = (/ mask_start_l(mid,1), mask_start_l(mid,2),            &
716                          mask_start_l(mid,3), domask_time_count(mid,av) /),   &
717               count = (/ mask_size_l(mid,1), mask_size_l(mid,2),              &
718                          mask_size_l(mid,3), 1 /) )
719          CALL netcdf_handle_error( 'data_output_mask', 461 )
720       ELSE
721#endif
722!
723!--       (1) b. Conventional I/O only through PE0
724!--       PE0 receives partial arrays from all processors of the respective mask
725!--       and outputs them. Here a barrier has to be set, because otherwise
726!--       "-MPI- FATAL: Remote protocol queue full" may occur.
727          CALL MPI_BARRIER( comm2d, ierr )
728
729          ngp = mask_size_l(mid,1) * mask_size_l(mid,2) * mask_size_l(mid,3)
730          IF ( myid == 0 )  THEN
731!
732!--          Local array can be relocated directly.
733             total_pf( &
734               mask_start_l(mid,1):mask_start_l(mid,1)+mask_size_l(mid,1)-1, &
735               mask_start_l(mid,2):mask_start_l(mid,2)+mask_size_l(mid,2)-1, &
736               mask_start_l(mid,3):mask_start_l(mid,3)+mask_size_l(mid,3)-1 ) &
737               = local_pf
738!
739!--          Receive data from all other PEs.
740             DO  n = 1, numprocs-1
741!
742!--             Receive index limits first, then array.
743!--             Index limits are received in arbitrary order from the PEs.
744                CALL MPI_RECV( ind(1), 6, MPI_INTEGER, MPI_ANY_SOURCE, 0,  &
745                     comm2d, status, ierr )
746!
747!--             Not all PEs have data for the mask
748                IF ( ind(1) /= -9999 )  THEN
749                   ngp = ( ind(2)-ind(1)+1 ) * (ind(4)-ind(3)+1 ) *  &
750                         ( ind(6)-ind(5)+1 )
751                   sender = status(MPI_SOURCE)
752                   DEALLOCATE( local_pf )
753                   ALLOCATE(local_pf(ind(1):ind(2),ind(3):ind(4),ind(5):ind(6)))
754                   CALL MPI_RECV( local_pf(ind(1),ind(3),ind(5)), ngp,  &
755                        MPI_REAL, sender, 1, comm2d, status, ierr )
756                   total_pf(ind(1):ind(2),ind(3):ind(4),ind(5):ind(6)) &
757                        = local_pf
758                ENDIF
759             ENDDO
760
761             nc_stat = NF90_PUT_VAR( id_set_mask(mid,av),                      &
762                  id_var_domask(mid,av,ivar), total_pf,                        &
763                  start = (/ 1, 1, 1, domask_time_count(mid,av) /),            &
764                  count = (/ mask_size(mid,1), mask_size(mid,2),               &
765                             mask_size(mid,3), 1 /) )
766             CALL netcdf_handle_error( 'data_output_mask', 462 )
767
768          ELSE
769!
770!--          If at least part of the mask resides on the PE, send the index
771!--          limits for the target array, otherwise send -9999 to PE0.
772             IF ( mask_size_l(mid,1) > 0 .AND.  mask_size_l(mid,2) > 0 .AND. &
773                  mask_size_l(mid,3) > 0  ) &
774                  THEN
775                ind(1) = mask_start_l(mid,1)
776                ind(2) = mask_start_l(mid,1) + mask_size_l(mid,1) - 1
777                ind(3) = mask_start_l(mid,2)
778                ind(4) = mask_start_l(mid,2) + mask_size_l(mid,2) - 1
779                ind(5) = mask_start_l(mid,3)
780                ind(6) = mask_start_l(mid,3) + mask_size_l(mid,3) - 1
781             ELSE
782                ind(1) = -9999; ind(2) = -9999
783                ind(3) = -9999; ind(4) = -9999
784                ind(5) = -9999; ind(6) = -9999
785             ENDIF
786             CALL MPI_SEND( ind(1), 6, MPI_INTEGER, 0, 0, comm2d, ierr )
787!
788!--          If applicable, send data to PE0.
789             IF ( ind(1) /= -9999 )  THEN
790                CALL MPI_SEND( local_pf(1,1,1), ngp, MPI_REAL, 0, 1, comm2d, &
791                     ierr )
792             ENDIF
793          ENDIF
794!
795!--       A barrier has to be set, because otherwise some PEs may proceed too
796!--       fast so that PE0 may receive wrong data on tag 0.
797          CALL MPI_BARRIER( comm2d, ierr )
798#if defined( __netcdf4_parallel )
799       ENDIF
800#endif
801#else
802!
803!--    (2) For serial execution of PALM, the single processor (PE0) holds all
804!--    data and writes them directly to file.
805       nc_stat = NF90_PUT_VAR( id_set_mask(mid,av),                            &
806                               id_var_domask(mid,av,ivar), local_pf,           &
807                             start = (/ 1, 1, 1, domask_time_count(mid,av) /), &
808                             count = (/ mask_size_l(mid,1), mask_size_l(mid,2),&
809                               mask_size_l(mid,3), 1 /) )
810       CALL netcdf_handle_error( 'data_output_mask', 463 )
811#endif
812
813       ivar = ivar + 1
814
815    ENDDO
816
817!
818!-- Deallocate temporary arrays.
819    DEALLOCATE( local_pf )
820#if defined( __parallel )
821    IF ( myid == 0 )  THEN
822       DEALLOCATE( total_pf )
823    ENDIF
824#endif
825
826!
827!-- Switch back to original format given by user (see beginning of this routine)
828    netcdf_data_format = netcdf_data_format_save
829
830    CALL cpu_log( log_point(49), 'data_output_mask', 'stop' )
831#endif
832
833
834 END SUBROUTINE data_output_mask
Note: See TracBrowser for help on using the repository browser.