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

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

Merge of branch palm4u into trunk

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