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

Last change on this file since 1359 was 1359, checked in by hoffmann, 10 years ago

new Lagrangian particle structure integrated

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