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

Last change on this file since 3046 was 3045, checked in by Giersch, 6 years ago

Code adjusted according to coding standards, renamed namelists, error messages revised until PA0347, output CASE 108 disabled

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