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

Last change on this file since 3573 was 3554, checked in by gronemeier, 6 years ago

renamed variable if to ivar; add variable description

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