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

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

further modularization of land surface model (2D/3D output and restart data). Bugfix for restart runs without land surface model

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