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

Last change on this file since 2716 was 2716, checked in by kanani, 6 years ago

Correction of "Former revisions" section

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