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

Last change on this file since 1977 was 1977, checked in by maronga, 8 years ago

last commit documented

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