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

Last change on this file since 3259 was 3241, checked in by raasch, 6 years ago

various changes to avoid compiler warnings (mainly removal of unused variables)

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