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

Last change on this file since 1330 was 1325, checked in by suehring, 11 years ago

last commit documented

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