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

Last change on this file since 1414 was 1414, checked in by gronemeier, 10 years ago

Bugfix: first and last grid points as they appear in 3D volume data can now be included to masked data output

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