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

Last change on this file since 3631 was 3589, checked in by suehring, 5 years ago

Remove erroneous UTF encoding; last commit documented

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