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

Last change on this file since 3494 was 3467, checked in by suehring, 6 years ago

Branch salsa @3446 re-integrated into trunk

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