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

Last change on this file since 4180 was 4180, checked in by scharf, 5 years ago

removed comments in 'Former revisions' section that are older than 01.01.2019

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