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

Last change on this file since 3742 was 3687, checked in by knoop, 5 years ago

Moved all user routunes that are dependencies of the PALM core only, to user_module.f90
The files that formerly contained these routines, have been deleted.
Also module_interface routines for init_mask and last_actions have been added.

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