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

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

last commit documented

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