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

Last change on this file since 2001 was 2001, checked in by knoop, 8 years ago

last commit documented

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