source: palm/trunk/SOURCE/init_masks.f90 @ 3552

Last change on this file since 3552 was 3512, checked in by gronemeier, 6 years ago

Bugfix: do not output ghost points if mask_x/y_loop is not specified (init_masks)

  • Property svn:keywords set to Id
File size: 34.4 KB
RevLine 
[1682]1!> @file init_masks.f90
[2000]2!------------------------------------------------------------------------------!
[2696]3! This file is part of the PALM model system.
[1036]4!
[2000]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.
[1036]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!
[2718]17! Copyright 1997-2018 Leibniz Universitaet Hannover
[2000]18!------------------------------------------------------------------------------!
[1036]19!
[484]20! Current revisions:
[298]21! -----------------
[1805]22!
[3049]23!
[1321]24! Former revisions:
25! -----------------
26! $Id: init_masks.f90 3512 2018-11-09 18:09:51Z suehring $
[3512]27! Bugfix: do not output ghost points if mask_x/y_loop is not specified
28!
29! 3467 2018-10-30 19:05:21Z suehring
[3467]30! Implementation of a new aerosol module salsa.
31!
32! 3435 2018-10-26 18:25:44Z gronemeier
[3435]33! Add checks for chemistry and radiation model
34! Set indices for terrain-following masked output
35!
36! 3421 2018-10-24 18:39:32Z gronemeier
[3421]37! Renamed output variables
38!
39! 3419 2018-10-24 17:27:31Z gronemeier
[3294]40! ocean renamed ocean_mode
41!
42! 3274 2018-09-24 15:42:55Z knoop
[3274]43! Modularization of all bulk cloud physics code components
44!
45! 3065 2018-06-12 07:03:02Z Giersch
[3065]46! dz_stretch_level was replaced by dz_stretch_level_start
47!
48! 3049 2018-05-29 13:52:36Z Giersch
[3045]49! Error messages revised
50!
[3049]51! 3045 2018-05-28 07:55:41Z Giersch
52! Error messages revised
53!
[3045]54! 2718 2018-01-02 08:49:38Z maronga
[2716]55! Corrected "Former revisions" section
56!
57! 2696 2017-12-14 17:12:51Z kanani
58! Change in file header (GPL part)
59!
60! 2344 2017-08-09 11:00:34Z raasch
[2344]61! explicit setting of initial values for array domask required due to a bug
62! in the Cray compiler (appears only if option -eD is used)
63!
64! 2301 2017-06-29 16:41:21Z gronemeier
[2301]65! Bugfix: set variable name length to global value
66!
67!
[2292]68! Implementation of new microphysic scheme: cloud_scheme = 'morrison'
69! includes two more prognostic equations for cloud drop concentration (nc) 
70! and cloud water content (qc).
71!
72! 2271 2017-06-09 12:34:55Z sward
[2271]73! Changed error message
74!
75! 2101 2017-01-05 16:42:31Z suehring
[1321]76!
[2032]77! 2031 2016-10-21 15:11:58Z knoop
78! renamed variable rho to rho_ocean
79!
[2001]80! 2000 2016-08-20 18:09:15Z knoop
81! Forced header and separation lines into 80 columns
82!
[1823]83! 1822 2016-04-07 07:49:42Z hoffmann
84! icloud_scheme replaced by microphysics_*
85!
[1805]86! 1804 2016-04-05 16:30:18Z maronga
87! Removed code for parameter file check (__check)
88!
[1784]89! 1783 2016-03-06 18:36:17Z raasch
90! netcdf module name changed + related changes
91!
[1683]92! 1682 2015-10-07 23:56:08Z knoop
93! Code annotations made doxygen readable
94!
[1439]95! 1438 2014-07-22 14:14:06Z heinze
96! +nr, qc, qr
97!
[1415]98! 1414 2014-05-31 11:19:48Z gronemeier
99! Bugfix: first and last grid points as they appear in 3D volume data can now
100! be included to masked data output
101!
[1354]102! 1353 2014-04-08 15:21:23Z heinze
103! REAL constants provided with KIND-attribute
104!
[1325]105! 1324 2014-03-21 09:13:16Z suehring
106! Bugfix: ONLY statement for module netcdf_control removed
107!
[1321]108! 1320 2014-03-20 08:40:49Z raasch
[1320]109! ONLY-attribute added to USE-statements,
110! kind-parameters added to all INTEGER and REAL declaration statements,
111! kinds are defined in new module kinds,
112! revision history before 2012 removed,
113! comment fields (!:) to be used for variable explanations added to
114! all variable declaration statements
[298]115!
[1187]116! 1186 2013-06-18 06:22:52Z raasch
117! bugfix: 0.0 replaced by zu(nzb) as the lowest default height level for masks,
118! because a zero value does not work in case of ocean runs
119!
[1037]120! 1036 2012-10-22 13:43:42Z raasch
121! code put under GPL (PALM 3.9)
122!
[1033]123! 1032 2012-10-21 13:03:21Z letzel
124! mask locations determined based on scalar positions
125!
[1035]126! 1031 2012-10-19 14:35:30Z raasch
127! netCDF4 without parallel file support implemented
128!
[997]129! 996 2012-09-07 10:41:47Z raasch
130! little reformatting
131!
[979]132! 978 2012-08-09 08:28:32Z fricke
133! +z0h*
134!
[810]135! 809 2012-01-30 13:32:58Z maronga
136! Bugfix: replaced .AND. and .NOT. with && and ! in the preprocessor directives
137!
[808]138! 807 2012-01-25 11:53:51Z maronga
139! New cpp directive "__check" implemented which is used by check_namelist_files
140!
[449]141! 410 2009-12-04 17:05:40Z letzel
142! Initial revision
[298]143!
[449]144!
[298]145! Description:
146! ------------
[1682]147!> Initialize masked data output
[298]148!------------------------------------------------------------------------------!
[1682]149 SUBROUTINE init_masks
[298]150
[1320]151    USE arrays_3d,                                                             &
152        ONLY:  zu, zw
153
[3435]154    USE bulk_cloud_model_mod,                                                  &
155        ONLY: bulk_cloud_model, microphysics_morrison, microphysics_seifert
156
157    USE chemistry_model_mod,                                                   &
158        ONLY:  chem_check_data_output
159
[1320]160    USE control_parameters,                                                    &
[3435]161        ONLY:  air_chemistry,                                                  &
162               constant_diffusion, cloud_droplets,                             &
[1320]163               data_output_masks, data_output_masks_user,                      &
[3065]164               doav, doav_n, domask, domask_no, dz, dz_stretch_level_start,    &
165               humidity, mask, masks, mask_scale, mask_i,                      &
[1320]166               mask_i_global, mask_j, mask_j_global, mask_k, mask_k_global,    &
[3435]167               mask_k_over_surface,                                            &
168               mask_loop, mask_size, mask_size_l, mask_start_l,                &
169               mask_surface, mask_x,                                           &
[1320]170               mask_x_loop, mask_xyz_dimension, mask_y, mask_y_loop, mask_z,   &
171               mask_z_loop, max_masks,  message_string, mid,                   &
[3294]172               passive_scalar, ocean_mode, varnamelength
[1320]173
174    USE grid_variables,                                                        &
175        ONLY:  dx, dy
176
177    USE indices,                                                               &
178        ONLY:  nx, nxl, nxr, ny, nyn, nys, nz, nzb, nzt
179
180    USE kinds
181
[1783]182    USE netcdf_interface,                                                      &
183        ONLY:  domask_unit, netcdf_data_format
[1320]184
185    USE particle_attributes,                                                   &
186        ONLY:  particle_advection
187
[3467]188    USE pegrid 
[298]189
[3435]190    USE radiation_model_mod,                                                   &
191        ONLY:  radiation, radiation_check_data_output
[3467]192       
193    USE salsa_mod,                                                             &
194        ONLY:  salsa, salsa_check_data_output 
[3435]195
[298]196    IMPLICIT NONE
197
[2301]198    CHARACTER (LEN=varnamelength) ::  var  !<
[1682]199    CHARACTER (LEN=7) ::  unit !<
[1320]200   
[2301]201    CHARACTER (LEN=varnamelength), DIMENSION(max_masks,100) ::  do_mask      !<
202    CHARACTER (LEN=varnamelength), DIMENSION(max_masks,100) ::  do_mask_user !<
[298]203
[3435]204    INTEGER(iwp) ::  count        !< counting masking indices along a dimension
[1682]205    INTEGER(iwp) ::  i            !<
206    INTEGER(iwp) ::  ilen         !<
207    INTEGER(iwp) ::  ind(6)       !<
208    INTEGER(iwp) ::  ind_array(1) !<
209    INTEGER(iwp) ::  j            !<
210    INTEGER(iwp) ::  k            !<
[3435]211    INTEGER(iwp) ::  m            !< mask index
[1682]212    INTEGER(iwp) ::  n            !<
213    INTEGER(iwp) ::  sender       !<
[1320]214   
[1682]215    INTEGER(iwp), DIMENSION(:), ALLOCATABLE ::  tmp_array !<
[298]216
[1682]217    LOGICAL ::  found !<
[2344]218
[298]219!
[2344]220!-- Initial values are explicitly set here due to a bug in the Cray compiler
221!-- in case of assignments of initial values in declaration statements for
222!-- arrays with more than 9999 elements (appears with -eD only)
223    domask = ' '
224
225!
[298]226!-- Allocation and initialization
[1414]227    ALLOCATE( tmp_array( MAX(nx,ny,nz)+2 ) )
[298]228
[1414]229    ALLOCATE( mask_i(max_masks,nxr-nxl+2), &
230              mask_j(max_masks,nyn-nys+2), &
231              mask_k(max_masks,nzt-nzb+2) )
[298]232!
233!-- internal mask arrays ("mask,dimension,selection")
234    ALLOCATE( mask(max_masks,3,mask_xyz_dimension), &
235              mask_loop(max_masks,3,3) )
[493]236   
[298]237!
[3045]238!-- Parallel mask output not yet supported. In check_parameters data format
239!-- is restricted and is switched back to non-parallel output. Therefore the
240!-- following error can not occur at the moment.
[1031]241    IF ( netcdf_data_format > 4 )  THEN
[3045]242       message_string = 'netCDF file formats '//                               &
243                        '5 and 6 (with parallel I/O support)'//                &
244                        ' are currently not supported.'
[564]245       CALL message( 'init_masks', 'PA0328', 1, 2, 0, 6, 0 )
[298]246    ENDIF
[553]247
[298]248!
249!-- Store data output parameters for masked data output in few shared arrays
[996]250    DO  mid = 1, masks
[564]251   
[553]252       do_mask     (mid,:) = data_output_masks(mid,:)
253       do_mask_user(mid,:) = data_output_masks_user(mid,:)
254       mask      (mid,1,:) = mask_x(mid,:) 
255       mask      (mid,2,:) = mask_y(mid,:)
256       mask      (mid,3,:) = mask_z(mid,:) 
[3435]257!
258!--    Flag a mask as terrain following
259       IF ( mask_k_over_surface(mid,1) /= -1_iwp )  THEN
260          mask_surface(mid) = .TRUE.
261       ENDIF
262
[1353]263       IF ( mask_x_loop(mid,1) == -1.0_wp  .AND.  mask_x_loop(mid,2) == -1.0_wp&
264            .AND.  mask_x_loop(mid,3) == -1.0_wp )  THEN
265          mask_loop(mid,1,1:2) = -1.0_wp
266          mask_loop(mid,1,3)   =  0.0_wp
[564]267       ELSE
268          mask_loop(mid,1,:) = mask_x_loop(mid,:)
269       ENDIF
[1353]270       IF ( mask_y_loop(mid,1) == -1.0_wp  .AND.  mask_y_loop(mid,2) == -1.0_wp&
271            .AND.  mask_y_loop(mid,3) == -1.0_wp )  THEN
272          mask_loop(mid,2,1:2) = -1.0_wp
273          mask_loop(mid,2,3)   =  0.0_wp
[564]274       ELSE
275          mask_loop(mid,2,:) = mask_y_loop(mid,:)
276       ENDIF
[1353]277       IF ( mask_z_loop(mid,1) == -1.0_wp  .AND.  mask_z_loop(mid,2) == -1.0_wp&
278            .AND.  mask_z_loop(mid,3) == -1.0_wp )  THEN
279          mask_loop(mid,3,1:2) = -1.0_wp
280          mask_loop(mid,3,3)   =  0.0_wp
[564]281       ELSE
282          mask_loop(mid,3,:) = mask_z_loop(mid,:)
283       ENDIF
284       
[553]285    ENDDO
286   
[298]287    mask_i = -1; mask_j = -1; mask_k = -1
[553]288   
[298]289!
290!-- Global arrays are required by define_netcdf_header.
[1031]291    IF ( myid == 0  .OR.  netcdf_data_format > 4 )  THEN
[1414]292       ALLOCATE( mask_i_global(max_masks,nx+2), &
293                 mask_j_global(max_masks,ny+2), &
294                 mask_k_global(max_masks,nz+2) )
[298]295       mask_i_global = -1; mask_j_global = -1; mask_k_global = -1
296    ENDIF
297
298!
299!-- Determine variable names for each mask
300    DO  mid = 1, masks
301!
302!--    Append user-defined data output variables to the standard data output
303       IF ( do_mask_user(mid,1) /= ' ' )  THEN
304          i = 1
305          DO  WHILE ( do_mask(mid,i) /= ' '  .AND.  i <= 100 )
306             i = i + 1
307          ENDDO
308          j = 1
309          DO  WHILE ( do_mask_user(mid,j) /= ' '  .AND.  j <= 100 )
310             IF ( i > 100 )  THEN
[3045]311                WRITE ( message_string, * ) 'number of output quantitities ',  &
312                     'given by data_output_mask and data_output_mask_user ',   &
[298]313                     'exceeds the limit of 100'
[564]314                CALL message( 'init_masks', 'PA0329', 1, 2, 0, 6, 0 )
[298]315             ENDIF
316             do_mask(mid,i) = do_mask_user(mid,j)
317             i = i + 1
318             j = j + 1
319          ENDDO
320       ENDIF
321
322!
323!--    Check and set steering parameters for mask data output and averaging
324       i   = 1
[996]325       DO WHILE ( do_mask(mid,i) /= ' '  .AND.  i <= 100 )
[298]326!
327!--       Check for data averaging
328          ilen = LEN_TRIM( do_mask(mid,i) )
329          j = 0                                              ! no data averaging
330          IF ( ilen > 3 )  THEN
331             IF ( do_mask(mid,i)(ilen-2:ilen) == '_av' )  THEN
332                j = 1                                           ! data averaging
333                do_mask(mid,i) = do_mask(mid,i)(1:ilen-3)
334             ENDIF
335          ENDIF
336          var = TRIM( do_mask(mid,i) )
337!
338!--       Check for allowed value and set units
339          SELECT CASE ( TRIM( var ) )
340
341             CASE ( 'e' )
342                IF ( constant_diffusion )  THEN
[3045]343                   WRITE ( message_string, * ) 'output of "', TRIM( var ),     &
[298]344                        '" requires constant_diffusion = .FALSE.'
[564]345                   CALL message( 'init_masks', 'PA0103', 1, 2, 0, 6, 0 )
[298]346                ENDIF
347                unit = 'm2/s2'
348
[3421]349             CASE ( 'thetal' )
[3274]350                IF ( .NOT. bulk_cloud_model )  THEN
[3045]351                   WRITE ( message_string, * ) 'output of "', TRIM( var ),     &
[3274]352                        '" requires bulk_cloud_model = .TRUE.'
[773]353                   CALL message( 'init_masks', 'PA0108', 1, 2, 0, 6, 0 )
[771]354                ENDIF
355                unit = 'K'
356
[2292]357             CASE ( 'nc' )
[3274]358                IF ( .NOT. bulk_cloud_model )  THEN
[3045]359                   WRITE ( message_string, * ) 'output of "', TRIM( var ),     &
[3274]360                        '" requires bulk_cloud_model = .TRUE.'
[2292]361                   CALL message( 'init_masks', 'PA0108', 1, 2, 0, 6, 0 )
362                 ELSEIF ( .NOT. microphysics_morrison ) THEN
[3045]363                   message_string = 'output of "' // TRIM( var ) // '" ' //    &
364                         'requires  = morrison'
[2292]365                   CALL message( 'check_parameters', 'PA0359', 1, 2, 0, 6, 0 )
366                ENDIF
367                unit = '1/m3'
368
[1438]369             CASE ( 'nr' )
[3274]370                IF ( .NOT. bulk_cloud_model )  THEN
[3045]371                   WRITE ( message_string, * ) 'output of "', TRIM( var ),     &
[3274]372                        '" requires bulk_cloud_model = .TRUE.'
[1438]373                   CALL message( 'init_masks', 'PA0108', 1, 2, 0, 6, 0 )
[1822]374                 ELSEIF ( .NOT. microphysics_seifert ) THEN
[3045]375                   message_string = 'output of "' // TRIM( var ) // '"' //     &
376                         'requires cloud_scheme = seifert_beheng'
[1438]377                   CALL message( 'check_parameters', 'PA0359', 1, 2, 0, 6, 0 )
378                ENDIF
379                unit = '1/m3'
380
[298]381             CASE ( 'pc', 'pr' )
382                IF ( .NOT. particle_advection )  THEN
[3045]383                   WRITE ( message_string, * ) 'output of "', TRIM( var ),     &
384                        '" requires a "particles_par"-NAMELIST in the ',       &
[298]385                        'parameter file (PARIN)'
[564]386                   CALL message( 'init_masks', 'PA0104', 1, 2, 0, 6, 0 )
[298]387                ENDIF
388                IF ( TRIM( var ) == 'pc' )  unit = 'number'
389                IF ( TRIM( var ) == 'pr' )  unit = 'm'
390
[3421]391             CASE ( 'q', 'thetav' )
[298]392                IF ( .NOT. humidity )  THEN
[3045]393                   WRITE ( message_string, * ) 'output of "', TRIM( var ),     &
[298]394                        '" requires humidity = .TRUE.'
[564]395                   CALL message( 'init_masks', 'PA0105', 1, 2, 0, 6, 0 )
[298]396                ENDIF
397                IF ( TRIM( var ) == 'q'   )  unit = 'kg/kg'
[3421]398                IF ( TRIM( var ) == 'thetav' )  unit = 'K'
[298]399
[1438]400             CASE ( 'qc' )
[3274]401                IF ( .NOT. bulk_cloud_model )  THEN
[3045]402                   message_string = 'output of "' // TRIM( var ) // '"' //     &
[3274]403                            'requires bulk_cloud_model = .TRUE.'
[1438]404                   CALL message( 'check_parameters', 'PA0108', 1, 2, 0, 6, 0 )
405                ENDIF
406                unit = 'kg/kg'
407
[298]408             CASE ( 'ql' )
[3274]409                IF ( .NOT. ( bulk_cloud_model  .OR.  cloud_droplets ) )  THEN
[3045]410                   WRITE ( message_string, * ) 'output of "', TRIM( var ),     &
[3274]411                        '" requires bulk_cloud_model = .TRUE. or ',            &
412                        'cloud_droplets = .TRUE.'
[3045]413                   CALL message( 'init_masks', 'PA0106', 1, 2, 0, 6, 0 )
[298]414                ENDIF
415                unit = 'kg/kg'
416
417             CASE ( 'ql_c', 'ql_v', 'ql_vp' )
418                IF ( .NOT. cloud_droplets )  THEN
[3045]419                   WRITE ( message_string, * ) 'output of "', TRIM( var ),     &
[298]420                        '" requires cloud_droplets = .TRUE.'
[564]421                   CALL message( 'init_masks', 'PA0107', 1, 2, 0, 6, 0 )
[298]422                ENDIF
423                IF ( TRIM( var ) == 'ql_c'  )  unit = 'kg/kg'
424                IF ( TRIM( var ) == 'ql_v'  )  unit = 'm3'
425                IF ( TRIM( var ) == 'ql_vp' )  unit = 'none'
426
427             CASE ( 'qv' )
[3274]428                IF ( .NOT. bulk_cloud_model )  THEN
[3045]429                   WRITE ( message_string, * ) 'output of "', TRIM( var ),     &
[3274]430                        '" requires bulk_cloud_model = .TRUE.'
[564]431                   CALL message( 'init_masks', 'PA0108', 1, 2, 0, 6, 0 )
[298]432                ENDIF
433                unit = 'kg/kg'
434
[1438]435             CASE ( 'qr' )
[3274]436                IF ( .NOT. bulk_cloud_model )  THEN
[3045]437                   message_string = 'output of "' // TRIM( var ) // '" ' //    &
[3274]438                            'requires bulk_cloud_model = .TRUE.'
[1438]439                   CALL message( 'check_parameters', 'PA0108', 1, 2, 0, 6, 0 )
[1822]440                ELSEIF ( .NOT. microphysics_seifert ) THEN
[3045]441                   message_string = 'output of "' // TRIM( var ) // '" ' //    &
442                            'requires cloud_scheme = seifert_beheng'
[1438]443                   CALL message( 'check_parameters', 'PA0359', 1, 2, 0, 6, 0 )
444                ENDIF
445                unit = 'kg/kg'
446
[3421]447             CASE ( 'rho_sea_water' )
[3294]448                IF ( .NOT. ocean_mode )  THEN
[3045]449                   WRITE ( message_string, * ) 'output of "', TRIM( var ),     &
[3294]450                        '" requires ocean mode'
[564]451                   CALL message( 'init_masks', 'PA0109', 1, 2, 0, 6, 0 )
[298]452                ENDIF
453                unit = 'kg/m3'
454
455             CASE ( 's' )
456                IF ( .NOT. passive_scalar )  THEN
[3045]457                   WRITE ( message_string, * ) 'output of "', TRIM( var ),     &
[298]458                        '" requires passive_scalar = .TRUE.'
[564]459                   CALL message( 'init_masks', 'PA0110', 1, 2, 0, 6, 0 )
[298]460                ENDIF
461                unit = 'conc'
462
463             CASE ( 'sa' )
[3294]464                IF ( .NOT. ocean_mode )  THEN
[3045]465                   WRITE ( message_string, * ) 'output of "', TRIM( var ),     &
[3294]466                        '" requires ocean mode'
[564]467                   CALL message( 'init_masks', 'PA0109', 1, 2, 0, 6, 0 )
[298]468                ENDIF
469                unit = 'psu'
470
[3421]471             CASE ( 'us*', 't*', 'lwp*', 'pra*', 'prr*', 'z0*', 'z0h*' )
[3045]472                WRITE ( message_string, * ) 'illegal value for data_',         &
473                     'output: "', TRIM( var ), '" is only allowed',            &
[298]474                     'for horizontal cross section'
[564]475                CALL message( 'init_masks', 'PA0111', 1, 2, 0, 6, 0 )
[298]476
[3421]477             CASE ( 'p', 'theta', 'u', 'v', 'w' )
[298]478                IF ( TRIM( var ) == 'p'  )  unit = 'Pa'
[3421]479                IF ( TRIM( var ) == 'theta' )  unit = 'K'
[298]480                IF ( TRIM( var ) == 'u'  )  unit = 'm/s'
481                IF ( TRIM( var ) == 'v'  )  unit = 'm/s'
482                IF ( TRIM( var ) == 'w'  )  unit = 'm/s'
483                CONTINUE
484
485             CASE DEFAULT
[3435]486             
[298]487                CALL user_check_data_output( var, unit )
[3467]488                           
489                IF ( salsa )  THEN
490                   CALL salsa_check_data_output( var, unit )
491                ENDIF               
[298]492
[3435]493                IF ( unit == 'illegal'  .AND.  air_chemistry                   &
494                     .AND.  (var(1:3) == 'kc_' .OR. var(1:3) == 'em_') )  THEN
495                   CALL chem_check_data_output( var, unit, 0, 0, 0 )
496                ENDIF
497
[298]498                IF ( unit == 'illegal' )  THEN
[3435]499                   CALL radiation_check_data_output( var, unit, 0, 0, 0 )
500                ENDIF
501
502                IF ( unit == 'illegal' )  THEN
[298]503                   IF ( do_mask_user(mid,1) /= ' ' )  THEN
[3045]504                      WRITE ( message_string, * ) 'illegal value for data_',   &
505                           'output_masks or data_output_masks_user: "',        &
[298]506                           TRIM( do_mask(mid,i) ), '"'
[3045]507                      CALL message( 'init_masks', 'PA0018', 1, 2, 0, 6, 0 )
[298]508                   ELSE
[3045]509                      WRITE ( message_string, * ) 'illegal value for data_',   &
510                           ' output_masks : "', TRIM( do_mask(mid,i) ), '"'
[564]511                      CALL message( 'init_masks', 'PA0330', 1, 2, 0, 6, 0 )
[298]512                   ENDIF
513                ENDIF
514
515          END SELECT
516!
517!--       Set the internal steering parameters appropriately
518          domask_no(mid,j)                    = domask_no(mid,j) + 1
519          domask(mid,j,domask_no(mid,j))      = do_mask(mid,i)
520          domask_unit(mid,j,domask_no(mid,j)) = unit
521
522          IF ( j == 1 )  THEN
523!
524!--          Check, if variable is already subject to averaging
525             found = .FALSE.
526             DO  k = 1, doav_n
527                IF ( TRIM( doav(k) ) == TRIM( var ) )  found = .TRUE.
528             ENDDO
529
530             IF ( .NOT. found )  THEN
531                doav_n = doav_n + 1
532                doav(doav_n) = var
533             ENDIF
534          ENDIF
535
536          i = i + 1
537
538       ENDDO   ! do_mask(mid,i)
539    ENDDO   ! mid
540
541
542!
543!-- Determine mask locations for each mask
544    DO  mid = 1, masks
545!
546!--    Set local masks for each subdomain along all three dimensions
547       CALL set_mask_locations( 1, dx, 'dx', nx, 'nx', nxl, nxr )
548       CALL set_mask_locations( 2, dy, 'dy', ny, 'ny', nys, nyn )
[3435]549       IF ( .NOT. mask_surface(mid) )  THEN
550          CALL set_mask_locations( 3, dz(1), 'dz', nz, 'nz', nzb, nzt )
551       ELSE
[298]552!
[3435]553!--       Set vertical mask locations and size in case of terrain-following
554!--       output
555          count = 0
556          DO  WHILE ( mask_k_over_surface(mid, count+1) >= 0 )
557             m = mask_k_over_surface(mid, count+1)
558             IF ( m > nz+1 )  THEN
559                WRITE ( message_string, '(I3,A,I3,A,I1,3A,I3)' )               &
560                     m,' in mask ',mid,' along dimension ', 3,                 &
561                     ' exceeds (nz+1) = ',nz+1
562                CALL message( 'init_masks', 'PA0331', 1, 2, 0, 6, 0 )
563             ENDIF
564             count = count + 1
565             mask_k(mid,count) = mask_k_over_surface(mid, count)
566             IF ( count == mask_xyz_dimension )  EXIT
567          ENDDO
568          mask_start_l(mid,3) = 1
569          mask_size(mid,3)    = count
570          mask_size_l(mid,3)  = count
571       ENDIF
572!
[298]573!--    Set global masks along all three dimensions (required by
574!--    define_netcdf_header).
[1804]575#if defined( __parallel )
[298]576!
577!--    PE0 receives partial arrays from all processors of the respective mask
578!--    and outputs them. Here a barrier has to be set, because otherwise
579!--    "-MPI- FATAL: Remote protocol queue full" may occur.
[807]580
[298]581       CALL MPI_BARRIER( comm2d, ierr )
582
583       IF ( myid == 0 )  THEN
584!
585!--       Local arrays can be relocated directly.
586          mask_i_global(mid,mask_start_l(mid,1): &
587                       mask_start_l(mid,1)+mask_size_l(mid,1)-1) = &
588                       mask_i(mid,:mask_size_l(mid,1))
589          mask_j_global(mid,mask_start_l(mid,2): &
590                       mask_start_l(mid,2)+mask_size_l(mid,2)-1) = &
591                       mask_j(mid,:mask_size_l(mid,2))
592          mask_k_global(mid,mask_start_l(mid,3): &
593                       mask_start_l(mid,3)+mask_size_l(mid,3)-1) = &
594                       mask_k(mid,:mask_size_l(mid,3))
595!
596!--       Receive data from all other PEs.
597          DO  n = 1, numprocs-1
598!
599!--          Receive index limits first, then arrays.
600!--          Index limits are received in arbitrary order from the PEs.
601             CALL MPI_RECV( ind(1), 6, MPI_INTEGER, MPI_ANY_SOURCE, 0,  &
602                  comm2d, status, ierr )
603!
604!--          Not all PEs have data for the mask.
605             IF ( ind(1) /= -9999 )  THEN
606                sender = status(MPI_SOURCE)
607                CALL MPI_RECV( tmp_array(ind(1)), ind(2)-ind(1)+1,  &
[363]608                               MPI_INTEGER, sender, 1, comm2d, status, ierr )
[298]609                mask_i_global(mid,ind(1):ind(2)) = tmp_array(ind(1):ind(2))
610                CALL MPI_RECV( tmp_array(ind(3)), ind(4)-ind(3)+1,  &
[363]611                               MPI_INTEGER, sender, 2, comm2d, status, ierr )
[298]612                mask_j_global(mid,ind(3):ind(4)) = tmp_array(ind(3):ind(4))
613                CALL MPI_RECV( tmp_array(ind(5)), ind(6)-ind(5)+1,  &
[363]614                               MPI_INTEGER, sender, 3, comm2d, status, ierr )
[298]615                mask_k_global(mid,ind(5):ind(6)) = tmp_array(ind(5):ind(6))
616             ENDIF
617          ENDDO
618
619       ELSE
620!
621!--       If at least part of the mask resides on the PE, send the index limits
622!--       for the target array, otherwise send -9999 to PE0.
[996]623          IF ( mask_size_l(mid,1) > 0  .AND.  mask_size_l(mid,2) > 0  .AND.  &
624               mask_size_l(mid,3) > 0  )  THEN
[298]625             ind(1) = mask_start_l(mid,1)
626             ind(2) = mask_start_l(mid,1) + mask_size_l(mid,1) - 1
627             ind(3) = mask_start_l(mid,2)
628             ind(4) = mask_start_l(mid,2) + mask_size_l(mid,2) - 1
629             ind(5) = mask_start_l(mid,3)
630             ind(6) = mask_start_l(mid,3) + mask_size_l(mid,3) - 1
631          ELSE
632             ind(1) = -9999; ind(2) = -9999
633             ind(3) = -9999; ind(4) = -9999
634             ind(5) = -9999; ind(6) = -9999
635          ENDIF
636          CALL MPI_SEND( ind(1), 6, MPI_INTEGER, 0, 0, comm2d, ierr )
637!
638!--       If applicable, send data to PE0.
639          IF ( ind(1) /= -9999 )  THEN
640             tmp_array(:mask_size_l(mid,1)) = mask_i(mid,:mask_size_l(mid,1))
641             CALL MPI_SEND( tmp_array(1), mask_size_l(mid,1),  &
[363]642                            MPI_INTEGER, 0, 1, comm2d, ierr )
[298]643             tmp_array(:mask_size_l(mid,2)) = mask_j(mid,:mask_size_l(mid,2))
644             CALL MPI_SEND( tmp_array(1), mask_size_l(mid,2),  &
[363]645                            MPI_INTEGER, 0, 2, comm2d, ierr )
[298]646             tmp_array(:mask_size_l(mid,3)) = mask_k(mid,:mask_size_l(mid,3))
647             CALL MPI_SEND( tmp_array(1), mask_size_l(mid,3),  &
[363]648                            MPI_INTEGER, 0, 3, comm2d, ierr )
[298]649          ENDIF
650       ENDIF
651!
652!--    A barrier has to be set, because otherwise some PEs may proceed too fast
653!--    so that PE0 may receive wrong data on tag 0.
654       CALL MPI_BARRIER( comm2d, ierr )
[409]655       
[1031]656       IF ( netcdf_data_format > 4 )  THEN
[409]657         
[1414]658          CALL MPI_BCAST( mask_i_global(mid,:), nx+2, MPI_INTEGER, 0, comm2d, &
[493]659                          ierr )
[1414]660          CALL MPI_BCAST( mask_j_global(mid,:), ny+2, MPI_INTEGER, 0, comm2d, &
[493]661                          ierr )
[1414]662          CALL MPI_BCAST( mask_k_global(mid,:), nz+2, MPI_INTEGER, 0, comm2d, &
[493]663                          ierr ) 
[1682]664     
[409]665       ENDIF
[298]666
[1804]667#else
[298]668!
669!--    Local arrays can be relocated directly.
670       mask_i_global(mid,:) = mask_i(mid,:)
671       mask_j_global(mid,:) = mask_j(mid,:)
672       mask_k_global(mid,:) = mask_k(mid,:)
673#endif
674    ENDDO   ! mid
675
676    DEALLOCATE( tmp_array )
677!
678!-- Internal mask arrays cannot be deallocated on PE 0 because they are
679!-- required for header output on PE 0.
680    IF ( myid /= 0 )  DEALLOCATE( mask, mask_loop )
681
682 CONTAINS
[493]683
[298]684!------------------------------------------------------------------------------!
685! Description:
686! ------------
[1682]687!> Set local mask for each subdomain along 'dim' direction.
[298]688!------------------------------------------------------------------------------!
[1682]689    SUBROUTINE set_mask_locations( dim, dxyz, dxyz_string, nxyz, nxyz_string, &
690                                   lb, ub )
[298]691
692       IMPLICIT NONE
693
[1682]694       CHARACTER (LEN=2) ::  dxyz_string !<
695       CHARACTER (LEN=2) ::  nxyz_string !<
[1320]696       
[1682]697       INTEGER(iwp)  ::  count       !<
698       INTEGER(iwp)  ::  count_l     !<
699       INTEGER(iwp)  ::  dim         !<
700       INTEGER(iwp)  ::  m           !<
701       INTEGER(iwp)  ::  loop_begin  !<
702       INTEGER(iwp)  ::  loop_end    !<
703       INTEGER(iwp)  ::  loop_stride !<
704       INTEGER(iwp)  ::  lb          !<
705       INTEGER(iwp)  ::  nxyz        !<
706       INTEGER(iwp)  ::  ub          !<
[1320]707       
[1682]708       REAL(wp)      ::  dxyz  !<
709       REAL(wp)      ::  ddxyz !<
710       REAL(wp)      ::  tmp1  !<
711       REAL(wp)      ::  tmp2  !<
[298]712
[1353]713       count = 0;  count_l = 0 
714       ddxyz = 1.0_wp / dxyz 
715       tmp1  = 0.0_wp
716       tmp2  = 0.0_wp
[298]717
[1353]718       IF ( mask(mid,dim,1) >= 0.0_wp )  THEN
[298]719!
720!--       use predefined mask_* array
[1353]721          DO  WHILE ( mask(mid,dim,count+1) >= 0.0_wp )
[298]722             count = count + 1
723             IF ( dim == 1 .OR. dim == 2 )  THEN
[1353]724                m = NINT( mask(mid,dim,count) * mask_scale(dim) * ddxyz - 0.5_wp )
[1414]725                IF ( m < 0 )  m = 0  ! avoid negative values
[298]726             ELSEIF ( dim == 3 )  THEN
727                ind_array =  &
[595]728                     MINLOC( ABS( mask(mid,dim,count) * mask_scale(dim) - zu ) )
[298]729                m = ind_array(1) - 1 + nzb  ! MINLOC uses lower array bound 1
730             ENDIF
[1414]731             IF ( m > (nxyz+1) )  THEN
[3045]732                WRITE ( message_string, '(I3,A,I3,A,I1,3A,I3)' )               &
733                     m,' in mask ',mid,' along dimension ',dim,                &
[1414]734                     ' exceeds (',nxyz_string,'+1) = ',nxyz+1
[564]735                CALL message( 'init_masks', 'PA0331', 1, 2, 0, 6, 0 )
[298]736             ENDIF
[1414]737             IF ( ( m >= lb .AND. m <= ub ) .OR.     &
738                  ( m == (nxyz+1) .AND. ub == nxyz )  )  THEN
[298]739                IF ( count_l == 0 )  mask_start_l(mid,dim) = count
740                count_l = count_l + 1
741                IF ( dim == 1 )  THEN
742                   mask_i(mid,count_l) = m
743                ELSEIF ( dim == 2 )  THEN
744                   mask_j(mid,count_l) = m
745                ELSEIF ( dim == 3 )  THEN
746                   mask_k(mid,count_l) = m
747                ENDIF
748             ENDIF
749             IF ( count == mask_xyz_dimension )  EXIT
750          ENDDO
751          mask_size(mid,dim)   = count
752          mask_size_l(mid,dim) = count_l
753
754       ELSE
755!
756!--       use predefined mask_loop_* array, or use the default (all grid points
757!--       along this direction)
[1353]758          IF ( mask_loop(mid,dim,1) < 0.0_wp )  THEN
[298]759             tmp1 = mask_loop(mid,dim,1)
[1186]760             mask_loop(mid,dim,1) = zw(nzb)  !   lowest level  (default)
[298]761          ENDIF
762          IF ( dim == 1 .OR. dim == 2 )  THEN
[1353]763             IF ( mask_loop(mid,dim,2) < 0.0_wp )  THEN
[298]764                tmp2 = mask_loop(mid,dim,2)
[3512]765                mask_loop(mid,dim,2) = nxyz*dxyz / mask_scale(dim)   ! (default)
[298]766             ENDIF
767             IF ( MAXVAL( mask_loop(mid,dim,1:2) )  &
[1414]768                  > (nxyz+1) * dxyz / mask_scale(dim) )  THEN
[303]769                WRITE ( message_string, '(2(A,I3,A,I1,A,F9.3),5A,I1,A,F9.3)' ) &
[3045]770                     'mask_loop(',mid,',',dim,',1)=',mask_loop(mid,dim,1),     &
[303]771                     ' and/or mask_loop(',mid,',',dim,',2)=', &
[3045]772                     mask_loop(mid,dim,2),' exceed (', &
773                     nxyz_string,'+1)*',dxyz_string,'/mask_scale(',dim,')=',   &
[1414]774                     (nxyz+1)*dxyz/mask_scale(dim)
[564]775                CALL message( 'init_masks', 'PA0332', 1, 2, 0, 6, 0 )
[298]776             ENDIF
[3045]777             loop_begin  = NINT( mask_loop(mid,dim,1) * mask_scale(dim)        &
[1353]778                  * ddxyz - 0.5_wp )
[3045]779             loop_end    = NINT( mask_loop(mid,dim,2) * mask_scale(dim)        &
[1353]780                  * ddxyz - 0.5_wp )
[3045]781             loop_stride = NINT( mask_loop(mid,dim,3) * mask_scale(dim)        &
[1032]782                  * ddxyz )
[1414]783             IF ( loop_begin == -1 )  loop_begin = 0  ! avoid negative values
[298]784          ELSEIF ( dim == 3 )  THEN
[1353]785             IF ( mask_loop(mid,dim,2) < 0.0_wp )  THEN
[298]786                tmp2 = mask_loop(mid,dim,2)
[1414]787                mask_loop(mid,dim,2) = zu(nz+1) / mask_scale(dim)   ! (default)
[298]788             ENDIF
789             IF ( MAXVAL( mask_loop(mid,dim,1:2) )  &
[1414]790                  > zu(nz+1) / mask_scale(dim) )  THEN
[3045]791                WRITE ( message_string, '(2(A,I3,A,I1,A,F9.3),A,I1,A,F9.3)' )  &
792                     'mask_loop(',mid,',',dim,',1)=',mask_loop(mid,dim,1),     &
[303]793                     ' and/or mask_loop(',mid,',',dim,',2)=', &
[3045]794                     mask_loop(mid,dim,2),' exceed zu(nz+1)/mask_scale(',dim,  &
[1414]795                     ')=',zu(nz+1)/mask_scale(dim)
[564]796                CALL message( 'init_masks', 'PA0333', 1, 2, 0, 6, 0 )
[298]797             ENDIF
[303]798             ind_array =  &
[595]799                  MINLOC( ABS( mask_loop(mid,dim,1) * mask_scale(dim) - zu ) )
[298]800             loop_begin =  &
[303]801                  ind_array(1) - 1 + nzb ! MINLOC uses lower array bound 1
802             ind_array =  &
[595]803                  MINLOC( ABS( mask_loop(mid,dim,2) * mask_scale(dim) - zu ) )
[298]804             loop_end = ind_array(1) - 1 + nzb ! MINLOC uses lower array bound 1
805!
[303]806!--          The following line assumes a constant vertical grid spacing within
807!--          the vertical mask range; it fails for vertical grid stretching.
[3065]808!--          Maybe revise later. Issue warning but continue execution. ABS(...)
809!--          within the IF statement is necessary because the default value of
810!--          dz_stretch_level_start is -9999999.9_wp.
[303]811             loop_stride = NINT( mask_loop(mid,dim,3) * mask_scale(dim) * ddxyz )
812
[3065]813             IF ( mask_loop(mid,dim,2) * mask_scale(dim) >                     &
814                  ABS( dz_stretch_level_start(1) ) )  THEN
[3045]815                WRITE ( message_string, '(A,I3,A,I1,A,F9.3,A,F8.2,3A)' )       &
816                     'mask_loop(',mid,',',dim,',2)=', mask_loop(mid,dim,2),    &
[3065]817                     ' exceeds dz_stretch_level=',dz_stretch_level_start(1),   &
[3046]818                     '.&Vertical mask locations will not ',                    &
[3045]819                     'match the desired heights within the stretching ',       &
[2271]820                     'region.'
[564]821                CALL message( 'init_masks', 'PA0334', 0, 1, 0, 6, 0 )
[303]822             ENDIF
823
[298]824          ENDIF
825!
826!--       If necessary, reset mask_loop(mid,dim,1) and mask_loop(mid,dim,2).
[1353]827          IF ( tmp1 < 0.0_wp )  mask_loop(mid,dim,1) = tmp1
828          IF ( tmp2 < 0.0_wp )  mask_loop(mid,dim,2) = tmp2
[298]829!
830!--       The default stride +/-1 (every grid point) applies if
831!--       mask_loop(mid,dim,3) is not specified (its default is zero).
832          IF ( loop_stride == 0 )  THEN
833             IF ( loop_end >= loop_begin )  THEN
834                loop_stride =  1
835             ELSE
836                loop_stride = -1
837             ENDIF
838          ENDIF
839          DO  m = loop_begin, loop_end, loop_stride
840             count = count + 1
[1414]841             IF ( ( m >= lb  .AND.  m <= ub ) .OR.   &
842                  ( m == (nxyz+1) .AND. ub == nxyz )  )  THEN
[298]843                IF ( count_l == 0 )  mask_start_l(mid,dim) = count
844                count_l = count_l + 1
845                IF ( dim == 1 )  THEN
846                   mask_i(mid,count_l) = m
847                ELSEIF ( dim == 2 )  THEN
848                   mask_j(mid,count_l) = m
849                ELSEIF ( dim == 3 )  THEN
850                   mask_k(mid,count_l) = m
851                ENDIF
852             ENDIF
853          ENDDO
854          mask_size(mid,dim)   = count
855          mask_size_l(mid,dim) = count_l
[1414]856
[298]857       ENDIF
858
859    END SUBROUTINE set_mask_locations
860
861 END SUBROUTINE init_masks
Note: See TracBrowser for help on using the repository browser.