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

Last change on this file since 3130 was 3083, checked in by gronemeier, 6 years ago

merge with branch rans

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