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

Last change on this file since 3423 was 3421, checked in by gronemeier, 5 years ago

new surface-data output; renamed output variables (pt to theta, rho_air to rho, rho_ocean to rho_sea_water)

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