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

Last change on this file since 2163 was 2101, checked in by suehring, 8 years ago

last commit documented

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