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

Last change on this file since 3355 was 3298, checked in by kanani, 5 years ago

Merge chemistry branch at r3297 to trunk

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