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

Last change on this file since 1327 was 1327, checked in by raasch, 10 years ago

Changed:


-s real64 removed (.mrun.config.hlrnIII)
-r8 removed (.mrun.config.imuk)
deleted: .mrun.config.imuk_ice2_netcdf4 .mrun.config.imuk_hlrn

REAL constants defined as wp-kind in modules

"baroclinicity" renamed "baroclinity", "ocean version" replaced by
"ocean mode"

code parts concerning old output formats "iso2d" and "avs" removed.
netCDF is the only remaining output format.

Errors:


bugfix: duplicate error message 56 removed

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