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

Last change on this file since 3538 was 3467, checked in by suehring, 5 years ago

Branch salsa @3446 re-integrated into trunk

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