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

Last change on this file since 1818 was 1818, checked in by maronga, 8 years ago

last commit documented / copyright update

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