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

Last change on this file since 4356 was 4346, checked in by motisi, 4 years ago

Introduction of wall_flags_total_0, which currently sets bits based on static topography information used in wall_flags_static_0

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