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

Last change on this file since 1186 was 1186, checked in by raasch, 11 years ago

bugfix: 0.0 replaced by zu(nzb) as the lowest default height level for masks,
because a zero value does not work in case of ocean runs

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