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

Last change on this file since 3570 was 3554, checked in by gronemeier, 5 years ago

renamed variable if to ivar; add variable description

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