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

Last change on this file since 1783 was 1783, checked in by raasch, 8 years ago

NetCDF routines modularized; new parameter netcdf_deflate; further changes in the pmc

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