source: palm/trunk/SOURCE/data_output_profiles.f90 @ 550

Last change on this file since 550 was 392, checked in by raasch, 14 years ago

New:
---

Adapted for machine lck
(mrun, mbuild, subjob)

bc_lr/bc_ns in most subroutines replaced by LOGICAL variables bc_lr_cyc,
bc_ns_cyc for speed optimization
(check_parameters, diffusion_u, diffusion_v, diffusion_w, modules)

Additional timestep criterion in case of simulations with plant canopy (timestep)

Check for illegal entries in section_xy|xz|yz that exceed nz+1|ny+1|nx+1
(check_parameters)

Clipping of dvrp output implemented. Default colourtable for particles
implemented, particle attributes (color, dvrp_size) can be set with new
parameters particle_color, particle_dvrpsize, color_interval,
dvrpsize_interval (init_dvrp, data_output_dvrp, modules, user_data_output_dvrp).
Slicer attributes (dvrp) are set with new routine set_slicer_attributes_dvrp
and are controlled with existing parameters slicer_range_limits.
(set_slicer_attributes_dvrp)

Ocean atmosphere coupling allows to use independent precursor runs in order
to account for different spin-up times. The time when coupling has to be
started is given by new inipar parameter coupling_start_time. The precursor
ocean run has to be started using new mrun option "-y" in order to add
appendix "_O" to all output files.
(check_for_restart, check_parameters, data_output_2d, data_output_3d,
data_output_profiles, data_output_ptseries, data_output_spectra,
data_output_tseries, header, init_coupling, modules, mrun,
parin, read_var_list, surface_coupler, time_integration, write_var_list)

Polygon reduction for topography and ground plate isosurface. Reduction level
for buildings can be chosen with parameter cluster_size. (init_dvrp)

External pressure gradient (check_parameters, header, init_3d_model, modules,
parin, prognostic_equations, read_var_list, write_var_list)

New topography case 'single_street_canyon' (header, init_grid, modules, parin,
read_var_list, user_check_parameters, user_header, user_init_grid, write_var_list)

Option to predefine a target bulk velocity for conserve_volume_flow
(check_parameters, header, init_3d_model, modules, parin, read_var_list,
write_var_list)

Option for user defined 2D data output in xy cross sections at z=nzb+1
(data_output_2d, user_data_output_2d)

xy cross section output of surface heatfluxes (latent, sensible)
(average_3d_data, check_parameters, data_output_2d, modules, read_3d_binary,
sum_up_3d_data, write_3d_binary)

average_3d_data, check_for_restart, check_parameters, data_output_2d, data_output_3d, data_output_dvrp, data_output_profiles, data_output_ptseries, data_output_spectra, data_output_tseries, init_coupling, init_dvrp, init_grid, init_3d_model, header, mbuild, modules, mrun, package_parin, parin, prognostic_equations, read_3d_binary, read_var_list, subjob, surface_coupler, timestep, time_integration, user_check_parameters, user_data_output_2d, user_data_output_dvrp, user_header, user_init_grid, write_3d_binary, write_var_list

New: set_particle_attributes, set_slicer_attributes_dvrp

Changed:


lcmuk changed to lc to avoid problems with Intel compiler on sgi-ice
(poisfft)

For extended NetCDF files, the updated title attribute includes an update of
time_average_text where appropriate. (netcdf)

In case of restart runs without extension, initial profiles are not written
to NetCDF-file anymore. (data_output_profiles, modules, read_var_list, write_var_list)

Small change in formatting of the message handling routine concering the output in the
job protocoll. (message)

initializing_actions='read_data_for_recycling' renamed to 'cyclic_fill', now
independent of turbulent_inflow (check_parameters, header, init_3d_model)

2 NetCDF error numbers changed. (data_output_3d)

A Link to the website appendix_a.html is printed for further information
about the possible errors. (message)

Temperature gradient criterion for estimating the boundary layer height
replaced by the gradient criterion of Sullivan et al. (1998). (flow_statistics)

NetCDF unit attribute in timeseries output in case of statistic regions added
(netcdf)

Output of NetCDF messages with aid of message handling routine.
(check_open, close_file, data_output_2d, data_output_3d,
data_output_profiles, data_output_ptseries, data_output_spectra,
data_output_tseries, netcdf, output_particles_netcdf)

Output of messages replaced by message handling routine.
(advec_particles, advec_s_bc, buoyancy, calc_spectra, check_for_restart,
check_open, coriolis, cpu_log, data_output_2d, data_output_3d, data_output_dvrp,
data_output_profiles, data_output_spectra, fft_xy, flow_statistics, header,
init_1d_model, init_3d_model, init_dvrp, init_grid, init_particles, init_pegrid,
netcdf, parin, plant_canopy_model, poisfft_hybrid, poismg, read_3d_binary,
read_var_list, surface_coupler, temperton_fft, timestep, user_actions,
user_data_output_dvrp, user_dvrp_coltab, user_init_grid, user_init_plant_canopy,
user_parin, user_read_restart_data, user_spectra )

Maximum number of tails is calculated from maximum number of particles and
skip_particles_for_tail (init_particles)

Value of vertical_particle_advection may differ for each particle group
(advec_particles, header, modules)

First constant in array den also defined as type double. (eqn_state_seawater)

Parameter dvrp_psize moved from particles_par to dvrp_graphics_par. (package_parin)

topography_grid_convention moved from userpar to inipar (check_parameters,
header, parin, read_var_list, user_check_parameters, user_header,
user_init_grid, user_parin, write_var_list)

Default value of grid_matching changed to strict.

Adjustments for runs on lcxt4 (necessary due to an software update on CRAY) and
for coupled runs on ibmy (mrun, subjob)

advec_particles, advec_s_bc, buoyancy, calc_spectra, check_for_restart, check_open, check_parameters, close_file, coriolis, cpu_log, data_output_2d, data_output_3d, data_output_dvrp, data_output_profiles, data_output_ptseries, data_output_spectra, data_output_tseries, eqn_state_seawater, fft_xy, flow_statistics, header, init_1d_model, init_3d_model, init_dvrp, init_grid, init_particles, init_pegrid, message, mrun, netcdf, output_particles_netcdf, package_parin, parin, plant_canopy_model, poisfft, poisfft_hybrid, poismg, read_3d_binary, read_var_list, sort_particles, subjob, user_check_parameters, user_header, user_init_grid, user_parin, surface_coupler, temperton_fft, timestep, user_actions, user_data_output_dvrp, user_dvrp_coltab, user_init_grid, user_init_plant_canopy, user_parin, user_read_restart_data, user_spectra, write_var_list

Errors:


Bugfix: Initial hydrostatic pressure profile in case of ocean runs is now
calculated in 5 iteration steps. (init_ocean)

Bugfix: wrong sign in buoyancy production of ocean part in case of not using
the reference density (only in 3D routine production_e) (production_e)

Bugfix: output of averaged 2d/3d quantities requires that an avaraging
interval has been set, respective error message is included (check_parameters)

Bugfix: Output on unit 14 only if requested by write_binary.
(user_last_actions)

Bugfix to avoid zero division by km_neutral (production_e)

Bugfix for extended NetCDF files: In order to avoid 'data mode' errors if
updated attributes are larger than their original size, NF90_PUT_ATT is called
in 'define mode' enclosed by NF90_REDEF and NF90_ENDDEF calls. This implies a
possible performance loss; an alternative strategy would be to ensure equal
attribute size in a job chain. (netcdf)

Bugfix: correction of initial volume flow for non-flat topography (init_3d_model)
Bugfix: zero initialization of arrays within buildings for 'cyclic_fill' (init_3d_model)

Bugfix: to_be_resorted => s_av for time-averaged scalars (data_output_2d, data_output_3d)

Bugfix: error in formatting the output (message)

Bugfix: avoid that ngp_2dh_s_inner becomes zero (init_3_model)

Typographical error: unit of wpt in dots_unit (modules)

Bugfix: error in check, if particles moved further than one subdomain length.
This check must not be applied for newly released particles. (advec_particles)

Bugfix: several tail counters are initialized, particle_tail_coordinates is
only written to file if its third index is > 0, arrays for tails are allocated
with a minimum size of 10 tails if there is no tail initially (init_particles,
advec_particles)

Bugfix: pressure included for profile output (check_parameters)

Bugfix: Type of count and count_rate changed to default INTEGER on NEC machines
(cpu_log)

Bugfix: output if particle time series only if particle advection is switched
on. (time_integration)

Bugfix: qsws was calculated in case of constant heatflux = .FALSE. (prandtl_fluxes)

Bugfix: averaging along z is not allowed for 2d quantities (e.g. u* and z0) (data_output_2d)

Typographical errors (netcdf)

If the inversion height calculated by the prerun is zero, inflow_damping_height
must be explicitly specified (init_3d_model)

Small bugfix concerning 3d 64bit netcdf output format (header)

Bugfix: dt_fixed removed from the restart file, because otherwise, no change
from a fixed to a variable timestep would be possible in restart runs.
(read_var_list, write_var_list)

Bugfix: initial setting of time_coupling in coupled restart runs (time_integration)

advec_particles, check_parameters, cpu_log, data_output_2d, data_output_3d, header, init_3d_model, init_particles, init_ocean, modules, netcdf, prandtl_fluxes, production_e, read_var_list, time_integration, user_last_actions, write_var_list

  • Property svn:keywords set to Id
File size: 25.7 KB
Line 
1 SUBROUTINE data_output_profiles
2
3!------------------------------------------------------------------------------!
4! Current revisions:
5! -----------------
6!
7!
8! Former revisions:
9! -----------------
10! $Id: data_output_profiles.f90 392 2009-09-24 10:39:14Z maronga $
11!
12! 345 2009-07-01 14:37:56Z heinze
13! In case of restart runs without extension, initial profiles are not written
14! to NetCDF-file anymore.
15! simulated_time in NetCDF output replaced by time_since_reference_point.
16! Output of NetCDF messages with aid of message handling routine.
17! Output of messages replaced by message handling routine.
18!
19! 197 2008-09-16 15:29:03Z raasch
20! Time coordinate t=0 stored on netcdf-file only if an output is required for
21! this time for at least one of the profiles
22!
23! February 2007
24! RCS Log replace by Id keyword, revision history cleaned up
25!
26! 87 2007-05-22 15:46:47Z raasch
27! var_hom renamed pr_palm
28!
29! Revision 1.18  2006/08/16 14:27:04  raasch
30! PRINT* statements for testing removed
31!
32! Revision 1.1  1997/09/12 06:28:48  raasch
33! Initial revision
34!
35!
36! Description:
37! ------------
38! Plot output of 1D-profiles for PROFIL
39!------------------------------------------------------------------------------!
40
41    USE control_parameters
42    USE cpulog
43    USE indices
44    USE interfaces
45    USE netcdf_control
46    USE pegrid
47    USE profil_parameter
48    USE statistics
49
50    IMPLICIT NONE
51
52
53    INTEGER ::  i, id, ilc, ils, j, k, sr
54    REAL    ::  uxma, uxmi
55
56
57!
58!-- If required, compute statistics
59    IF ( .NOT. flow_statistics_called )  CALL flow_statistics
60
61!
62!-- Flow_statistics has its own CPU time measurement
63    CALL cpu_log( log_point(15), 'data_output_profiles', 'start' )
64
65!
66!-- If required, compute temporal average
67    IF ( averaging_interval_pr == 0.0 )  THEN
68       hom_sum(:,:,:) = hom(:,1,:,:)
69    ELSE
70       IF ( average_count_pr > 0 )  THEN
71          hom_sum = hom_sum / REAL( average_count_pr )
72       ELSE
73!
74!--       This case may happen if dt_dopr is changed in the d3par-list of
75!--       a restart run
76          RETURN
77       ENDIF
78    ENDIF
79
80   
81    IF ( myid == 0 )  THEN
82
83!
84!--    Plot-output for each (sub-)region
85
86!
87!--    Open file for profile output in NetCDF format
88       IF ( netcdf_output )  THEN
89          CALL check_open( 104 )
90       ENDIF
91
92!
93!--    Open PROFIL-output files for each (sub-)region
94       IF ( profil_output )  THEN
95          DO  sr = 0, statistic_regions
96             CALL check_open( 40 + sr )
97          ENDDO
98       ENDIF
99
100!
101!--    Increment the counter for number of output times
102       dopr_time_count = dopr_time_count + 1
103
104!
105!--    Re-set to zero the counter for the number of profiles already written
106!--    at the current output time into the respective crosses
107       cross_pnc_local = 0
108
109!
110!--    Output of initial profiles
111       IF ( dopr_time_count == 1 )  THEN
112       
113          IF ( .NOT. output_for_t0 ) THEN
114
115             IF ( netcdf_output )  THEN
116#if defined( __netcdf )         
117!
118!--             Store initial time (t=0) to time axis, but only if an output
119!--             is required for at least one of the profiles
120                DO  i = 1, dopr_n
121                IF ( dopr_initial_index(i) /= 0 )  THEN
122                   nc_stat = NF90_PUT_VAR( id_set_pr, id_var_time_pr,  &
123                                              (/ 0.0 /), start = (/ 1 /), &
124                                              count = (/ 1 /) )
125                      CALL handle_netcdf_error( 'data_output_profiles', 329 )
126                      output_for_t0 = .TRUE.
127                      EXIT
128                   ENDIF
129                ENDDO
130
131!
132!--             Store normalization factors
133                nc_stat = NF90_PUT_VAR( id_set_pr, id_var_norm_dopr(1), & ! wpt0
134                                     (/ hom_sum(nzb,18,normalizing_region) /), &
135                                        start = (/ 1 /), count = (/ 1 /) )
136                CALL handle_netcdf_error( 'data_output_profiles', 330 )
137
138                nc_stat = NF90_PUT_VAR( id_set_pr, id_var_norm_dopr(2), & ! ws2
139                           (/ hom_sum(nzb+8,pr_palm,normalizing_region)**2 /), &
140                                        start = (/ 1 /), count = (/ 1 /) )
141                CALL handle_netcdf_error( 'data_output_profiles', 331 )
142                nc_stat = NF90_PUT_VAR( id_set_pr, id_var_norm_dopr(3), & ! tsw2
143                           (/ hom_sum(nzb+3,pr_palm,normalizing_region)**2 /), &
144                                     start = (/ 1 /), count = (/ 1 /) )
145                CALL handle_netcdf_error( 'data_output_profiles', 332 )
146                nc_stat = NF90_PUT_VAR( id_set_pr, id_var_norm_dopr(4), & ! ws3
147                           (/ hom_sum(nzb+8,pr_palm,normalizing_region)**3 /), &
148                                        start = (/ 1 /), count = (/ 1 /) )
149                CALL handle_netcdf_error( 'data_output_profiles', 333 )
150
151                nc_stat = NF90_PUT_VAR( id_set_pr, id_var_norm_dopr(5), &!ws2tsw
152                           (/ hom_sum(nzb+8,pr_palm,normalizing_region)**3 *   &
153                              hom_sum(nzb+3,pr_palm,normalizing_region)    /), &
154                                        start = (/ 1 /), count = (/ 1 /) )
155                CALL handle_netcdf_error( 'data_output_profiles', 334 )
156
157                nc_stat = NF90_PUT_VAR( id_set_pr, id_var_norm_dopr(6), &!wstsw2
158                           (/ hom_sum(nzb+8,pr_palm,normalizing_region) *      &
159                              hom_sum(nzb+3,pr_palm,normalizing_region)**2 /), &
160                                        start = (/ 1 /), count = (/ 1 /) )
161                CALL handle_netcdf_error( 'data_output_profiles', 335 )
162
163                nc_stat = NF90_PUT_VAR( id_set_pr, id_var_norm_dopr(7), & ! z_i
164                              (/ hom_sum(nzb+6,pr_palm,normalizing_region) /), &
165                                        start = (/ 1 /), count = (/ 1 /) )
166                CALL handle_netcdf_error( 'data_output_profiles', 336 )
167             
168#endif
169             ENDIF
170!
171!--          Loop over all 1D variables
172             DO  i = 1, dopr_n
173
174                IF ( dopr_initial_index(i) /= 0 )  THEN
175
176!
177!--                Output for the individual (sub-)regions
178                   DO  sr = 0, statistic_regions
179
180                      IF ( profil_output )  THEN
181                         id = 40 + sr
182!
183!--                      Write Label-Header
184                         WRITE ( id, 100 )  TRIM( data_output_pr(i) ), '(t=0)'
185!
186!--                      Write total profile
187                         DO  k = nzb, nzt+1
188                            WRITE ( id, 101 )  hom(k,2,dopr_initial_index(i),sr), &
189                                               hom(k,1,dopr_initial_index(i),sr)
190                         ENDDO
191!
192!--                      Write separation label
193                         WRITE ( id, 102 )
194                      ENDIF
195
196                      IF ( netcdf_output )  THEN
197#if defined( __netcdf )
198!
199!--                      Write data to netcdf file
200                         nc_stat = NF90_PUT_VAR( id_set_pr, id_var_dopr(i,sr),    &
201                                       hom(nzb:nzt+1,1,dopr_initial_index(i),sr), &
202                                                 start = (/ 1, 1 /),              &
203                                                 count = (/ nzt-nzb+2, 1 /) )
204                         CALL handle_netcdf_error( 'data_output_profiles', 337 )
205#endif
206                      ENDIF
207
208                   ENDDO
209
210                   IF ( profil_output )  THEN
211!
212!--                   Determine indices for later NAMELIST-output (s. below)
213                      profile_number = profile_number + 1
214                      j = dopr_crossindex(i)
215                      IF ( j /= 0 )  THEN
216                         cross_profile_number_count(j) = &
217                                                  cross_profile_number_count(j) + 1
218                         k = cross_profile_number_count(j)
219                         cross_profile_numbers(k,j) = profile_number
220!
221!--                      Initial profiles are always drawn as solid lines in
222!--                      anti-background colour.
223                         cross_linecolors(k,j) = 1
224                         cross_linestyles(k,j) = 0
225!
226!--                      If required, extend x-value range of the respective
227!--                      cross, provided it has not been specified in &
228!--                      check_parameters. Determination over all (sub-)regions.
229                         IF ( cross_uxmin(j) == 0.0  .AND. &
230                              cross_uxmax(j) == 0.0 )  THEN
231
232                            DO  sr = 0, statistic_regions
233
234                               uxmi = &
235                               MINVAL( hom(:nz_do1d,1,dopr_initial_index(i),sr) )
236 
237                               uxma = &
238                               MAXVAL( hom(:nz_do1d,1,dopr_initial_index(i),sr) )
239!
240!--                            When the value range of the first line in the
241!--                            corresponding cross is determined, its value range
242!--                            is simply adopted.
243                               IF ( cross_uxmin_computed(j) > &
244                                    cross_uxmax_computed(j) )  THEN
245                                  cross_uxmin_computed(j) = uxmi
246                                  cross_uxmax_computed(j) = uxma
247                               ELSE
248                                  cross_uxmin_computed(j) = &
249                                               MIN( cross_uxmin_computed(j), uxmi )
250                                  cross_uxmax_computed(j) = &
251                                               MAX( cross_uxmax_computed(j), uxma )
252                               ENDIF
253
254                            ENDDO
255
256                         ENDIF
257!
258!--                      If required, determine and note normalizing factors
259                         SELECT CASE ( cross_normalized_x(j) )
260   
261                             CASE ( 'ts2' )
262                               cross_normx_factor(k,j) = &
263                                ( hom_sum(nzb+3,pr_palm,normalizing_region) )**2
264                            CASE ( 'wpt0' )
265                               cross_normx_factor(k,j) = &
266                                hom_sum(nzb,18,normalizing_region)
267                            CASE ( 'wsts2' )
268                               cross_normx_factor(k,j) = &
269                                hom_sum(nzb+8,pr_palm,normalizing_region)  &
270                              * ( hom_sum(nzb+3,pr_palm,normalizing_region) )**2
271                            CASE ( 'ws2' )
272                               cross_normx_factor(k,j) = &
273                                ( hom_sum(nzb+8,pr_palm,normalizing_region) )**2
274                            CASE ( 'ws2ts' )
275                               cross_normx_factor(k,j) = &
276                              ( hom_sum(nzb+8,pr_palm,normalizing_region) )**2 &
277                              * hom_sum(nzb+3,pr_palm,normalizing_region)
278                             CASE ( 'ws3' )
279                               cross_normx_factor(k,j) = &
280                                ( hom_sum(nzb+8,pr_palm,normalizing_region) )**3
281
282                         END SELECT
283
284                         SELECT CASE ( cross_normalized_y(j) )
285 
286                            CASE ( 'z_i' )
287                                cross_normy_factor(k,j) = &
288                                        hom_sum(nzb+6,pr_palm,normalizing_region)
289
290                         END SELECT
291
292!
293!--                      Check the normalizing factors for zeros and deactivate
294!--                      the normalization, if required.
295                         IF ( cross_normx_factor(k,j) == 0.0  .OR. &
296                              cross_normy_factor(k,j) == 0.0 )  THEN
297                            WRITE( message_string, * ) 'data_output_profiles: ',  &
298                                'normalizing cross ',j, ' is not possible ',     &
299                                'since one of the & normalizing factors ',       &
300                                 'is zero! & cross_normx_factor(',k,',',j,') = ', &
301                                                         cross_normx_factor(k,j), &
302                                 ' & cross_normy_factor(',k,',',j,') = ',         &
303                                                         cross_normy_factor(k,j)
304                            CALL message( 'data_output_profiles', 'PA0185',&
305                                                                    0, 1, 0, 6, 0 )
306                            cross_normx_factor(k,j) = 1.0
307                            cross_normy_factor(k,j) = 1.0
308                            cross_normalized_x(j) = ' '
309                            cross_normalized_y(j) = ' '
310                         ENDIF
311
312!
313!--                      If required, extend normalized x-value range of the
314!--                      respective cross, provided it has not been specified in
315!--                      check_parameters. Determination over all (sub-)regions.
316                         IF ( cross_uxmin_normalized(j) == 0.0  .AND. &
317                              cross_uxmax_normalized(j) == 0.0 )  THEN
318
319                            DO  sr = 0, statistic_regions
320
321                               uxmi = MINVAL( hom(:nz_do1d,1,             &
322                                              dopr_initial_index(i),sr) ) / &
323                                      cross_normx_factor(k,j)
324                               uxma = MAXVAL( hom(:nz_do1d,1,             &
325                                              dopr_initial_index(i),sr) ) / &
326                                       cross_normx_factor(k,j)
327!
328!--                            When the value range of the first line in the
329!--                            corresponding cross is determined, its value range
330!--                            is simply adopted.
331                               IF ( cross_uxmin_normalized_computed(j) > &
332                                    cross_uxmax_normalized_computed(j) )  THEN
333                                  cross_uxmin_normalized_computed(j) = uxmi
334                                  cross_uxmax_normalized_computed(j) = uxma
335                               ELSE
336                                  cross_uxmin_normalized_computed(j) = &
337                                    MIN( cross_uxmin_normalized_computed(j), uxmi )
338                                  cross_uxmax_normalized_computed(j) = &
339                                    MAX( cross_uxmax_normalized_computed(j), uxma )
340                               ENDIF
341
342                            ENDDO
343
344                         ENDIF
345
346                      ENDIF   ! Index determination
347
348                   ENDIF   ! profil output
349
350                ENDIF   ! Initial profile available
351
352             ENDDO   ! Loop over dopr_n for initial profiles
353
354             IF ( netcdf_output  .AND.  output_for_t0 )  THEN
355                dopr_time_count = dopr_time_count + 1
356             ENDIF
357
358          END IF
359       ENDIF   ! Initial profiles
360
361       IF ( netcdf_output )  THEN
362#if defined( __netcdf )
363
364!
365!--       Store time to time axis         
366          nc_stat = NF90_PUT_VAR( id_set_pr, id_var_time_pr,        &
367                                  (/ time_since_reference_point /), &
368                                  start = (/ dopr_time_count /),    &
369                                  count = (/ 1 /) )
370          CALL handle_netcdf_error( 'data_output_profiles', 338 )
371
372!
373!--       Store normalization factors
374          nc_stat = NF90_PUT_VAR( id_set_pr, id_var_norm_dopr(1), &  ! wpt0
375                                  (/ hom_sum(nzb,18,normalizing_region) /), &
376                                  start = (/ dopr_time_count /),               &
377                                  count = (/ 1 /) )
378          CALL handle_netcdf_error( 'data_output_profiles', 339 )
379
380          nc_stat = NF90_PUT_VAR( id_set_pr, id_var_norm_dopr(2), &  ! ws2
381                        (/ hom_sum(nzb+8,pr_palm,normalizing_region)**2 /), &
382                                  start = (/ dopr_time_count /),               &
383                                  count = (/ 1 /) )
384          CALL handle_netcdf_error( 'data_output_profiles', 340 )
385
386          nc_stat = NF90_PUT_VAR( id_set_pr, id_var_norm_dopr(3), &  ! tsw2
387                        (/ hom_sum(nzb+3,pr_palm,normalizing_region)**2 /), &
388                                  start = (/ dopr_time_count /),               &
389                                  count = (/ 1 /) )
390          CALL handle_netcdf_error( 'data_output_profiles', 341 )
391
392          nc_stat = NF90_PUT_VAR( id_set_pr, id_var_norm_dopr(4), &  ! ws3
393                        (/ hom_sum(nzb+8,pr_palm,normalizing_region)**3 /), &
394                                  start = (/ dopr_time_count /),               &
395                                  count = (/ 1 /) )
396          CALL handle_netcdf_error( 'data_output_profiles', 342 )
397
398          nc_stat = NF90_PUT_VAR( id_set_pr, id_var_norm_dopr(5), &  ! ws2tsw
399                        (/ hom_sum(nzb+8,pr_palm,normalizing_region)**3 *   &
400                           hom_sum(nzb+3,pr_palm,normalizing_region)    /), &
401                                  start = (/ dopr_time_count /),               &
402                                  count = (/ 1 /) )
403          CALL handle_netcdf_error( 'data_output_profiles', 343 )
404         
405          nc_stat = NF90_PUT_VAR( id_set_pr, id_var_norm_dopr(6), &  ! wstsw2
406                        (/ hom_sum(nzb+8,pr_palm,normalizing_region) *      &
407                           hom_sum(nzb+3,pr_palm,normalizing_region)**2 /), &
408                                  start = (/ dopr_time_count /),               &
409                                  count = (/ 1 /) )
410          CALL handle_netcdf_error( 'data_output_profiles', 344 )
411
412          nc_stat = NF90_PUT_VAR( id_set_pr, id_var_norm_dopr(7), &  ! z_i
413                           (/ hom_sum(nzb+6,pr_palm,normalizing_region) /), &
414                                  start = (/ dopr_time_count /),               &
415                                  count = (/ 1 /) )
416          CALL handle_netcdf_error( 'data_output_profiles', 345 )
417#endif
418       ENDIF
419
420!
421!--    Output of the individual (non-initial) profiles
422       DO  i = 1, dopr_n
423
424!
425!--       Output for the individual (sub-)domains
426          DO  sr = 0, statistic_regions
427
428             IF ( profil_output )  THEN
429                id = 40 + sr
430!
431!--             Write Label-Header
432                WRITE ( id, 100 )  TRIM( dopr_label(i) ), simulated_time_chr
433!
434!--             Output of total profile
435                DO  k = nzb, nzt+1
436                   WRITE ( id, 101 )  hom(k,2,dopr_index(i),sr), &
437                                      hom_sum(k,dopr_index(i),sr)
438                ENDDO
439!
440!--             Write separation label
441                WRITE ( id, 102 )
442             ENDIF
443
444             IF ( netcdf_output )  THEN
445#if defined( __netcdf )
446!
447!--             Write data to netcdf file
448                nc_stat = NF90_PUT_VAR( id_set_pr, id_var_dopr(i,sr),          &
449                                        hom_sum(nzb:nzt+1,dopr_index(i),sr),&
450                                        start = (/ 1, dopr_time_count /),      &
451                                        count = (/ nzt-nzb+2, 1 /) )
452                CALL handle_netcdf_error( 'data_output_profiles', 346 )
453#endif
454             ENDIF
455
456          ENDDO
457
458          IF ( profil_output )  THEN
459!
460!--          Determine profile number on file and note the data for later
461!--          NAMELIST output, if the respective profile is to be drawn by
462!--          PROFIL (if it shall not be drawn, the variable dopr_crossindex has
463!--          the value 0, otherwise the number of the coordinate cross)
464             profile_number = profile_number + 1
465             j = dopr_crossindex(i)
466
467             IF ( j /= 0 )  THEN
468                cross_profile_number_count(j) = cross_profile_number_count(j) +1
469                k = cross_profile_number_count(j)
470                cross_pnc_local(j)            = cross_pnc_local(j)            +1
471                cross_profile_numbers(k,j) = profile_number
472                ilc = MOD( dopr_time_count, 10 )
473                IF ( ilc == 0 )  ilc = 10
474                cross_linecolors(k,j) = linecolors(ilc)
475                ils = MOD( cross_pnc_local(j), 11 )
476                IF ( ils == 0 )  ils = 11
477                cross_linestyles(k,j) = linestyles(ils)
478!
479!--             If required, extend x-value range of the respective coordinate
480!--             cross, provided it has not been specified in check_parameters.
481!--             Determination over all (sub-)regions.
482                IF ( cross_uxmin(j) == 0.0  .AND.  cross_uxmax(j) == 0.0 )  THEN
483
484                   DO  sr = 0, statistic_regions
485
486                      uxmi = MINVAL( hom_sum(:nz_do1d,dopr_index(i),sr) )
487                      uxma = MAXVAL( hom_sum(:nz_do1d,dopr_index(i),sr) )
488!
489!--                   When the value range of the first line in the
490!--                   corresponding cross is determined, its value range is
491!--                   simply adopted.
492                      IF ( cross_uxmin_computed(j) > cross_uxmax_computed(j) ) &
493                      THEN
494                         cross_uxmin_computed(j) = uxmi
495                         cross_uxmax_computed(j) = uxma
496                      ELSE
497                         cross_uxmin_computed(j) = &
498                                           MIN( cross_uxmin_computed(j), uxmi )
499                         cross_uxmax_computed(j) = &
500                                           MAX( cross_uxmax_computed(j), uxma )
501                      ENDIF
502
503                   ENDDO
504
505                ENDIF
506!
507!--             If required, store the normalizing factors
508                SELECT CASE ( cross_normalized_x(j) )
509
510                   CASE ( 'tsw2' )
511                      cross_normx_factor(k,j) = &
512                            ( hom_sum(nzb+11,pr_palm,normalizing_region) )**2
513                   CASE ( 'wpt0' )
514                      cross_normx_factor(k,j) = &
515                              hom_sum(nzb,18,normalizing_region)
516                   CASE ( 'wstsw2' )
517                      cross_normx_factor(k,j) = &
518                              hom_sum(nzb+8,pr_palm,normalizing_region)  &
519                          * ( hom_sum(nzb+11,pr_palm,normalizing_region) )**2
520                   CASE ( 'ws2' )
521                      cross_normx_factor(k,j) = &
522                            ( hom_sum(nzb+8,pr_palm,normalizing_region) )**2
523                   CASE ( 'ws2tsw' )
524                      cross_normx_factor(k,j) = &
525                            ( hom_sum(nzb+8,pr_palm,normalizing_region) )**2&
526                            * hom_sum(nzb+11,pr_palm,normalizing_region)
527                   CASE ( 'ws3' )
528                      cross_normx_factor(k,j) = &
529                            ( hom_sum(nzb+8,pr_palm,normalizing_region) )**3
530
531                END SELECT
532                SELECT CASE ( cross_normalized_y(j) )
533
534                   CASE ( 'z_i' )
535                      cross_normy_factor(k,j) = &
536                                   hom_sum(nzb+6,pr_palm,normalizing_region)
537
538                END SELECT
539
540!
541!--             Check the normalizing factors for zeros and deactivate
542!--             the normalization, if required.
543                IF ( cross_normx_factor(k,j) == 0.0  .OR. &
544                     cross_normy_factor(k,j) == 0.0 )  THEN
545                   WRITE( message_string, * ) 'data_output_profiles: ',  &
546                              'normalizing cross ',j, ' is not possible ',     &
547                              'since one of the & normalizing factors ',       &
548                              'is zero! & cross_normx_factor(',k,',',j,') = ', &
549                                                      cross_normx_factor(k,j), &
550                              ' & cross_normy_factor(',k,',',j,') = ',         &
551                                                      cross_normy_factor(k,j)
552                         CALL message( 'data_output_profiles', 'PA0185',&
553                                                                0, 1, 0, 6, 0 )
554                    cross_normx_factor(k,j) = 1.0
555                    cross_normy_factor(k,j) = 1.0
556                    cross_normalized_x(j) = ' '
557                    cross_normalized_y(j) = ' '
558                ENDIF
559
560!
561!--             If required, extend normalized x-value range of the respective 
562!--             cross, provided it has not been specified in check_parameters.
563!--             Determination over all (sub-)regions.
564                IF ( cross_uxmin_normalized(j) == 0.0  .AND. &
565                     cross_uxmax_normalized(j) == 0.0 )  THEN
566
567                   DO  sr = 0, statistic_regions
568
569                      uxmi = MINVAL( hom(:nz_do1d,1,dopr_index(i),sr) ) / &
570                             cross_normx_factor(k,j)
571                      uxma = MAXVAL( hom(:nz_do1d,1,dopr_index(i),sr) ) / &
572                             cross_normx_factor(k,j)
573!
574!--                   When the value range of the first line in the
575!--                   corresponding cross is determined, its value range is
576!--                   simply adopted.
577                      IF ( cross_uxmin_normalized_computed(j) > &
578                           cross_uxmax_normalized_computed(j) )  THEN
579                         cross_uxmin_normalized_computed(j) = uxmi
580                         cross_uxmax_normalized_computed(j) = uxma
581                      ELSE
582                         cross_uxmin_normalized_computed(j) = &
583                                MIN( cross_uxmin_normalized_computed(j), uxmi )
584                         cross_uxmax_normalized_computed(j) = &
585                                MAX( cross_uxmax_normalized_computed(j), uxma )
586                      ENDIF
587
588                   ENDDO
589
590                ENDIF
591
592             ENDIF   ! Index determination
593
594          ENDIF   ! profil output
595
596       ENDDO   ! Loop over dopr_n
597
598    ENDIF  ! Output on PE0
599
600!
601!-- If averaging has been done above, the summation counter must be re-set.
602    IF ( averaging_interval_pr /= 0.0 )  THEN
603       average_count_pr = 0
604    ENDIF
605
606    CALL cpu_log( log_point(15), 'data_output_profiles','stop', 'nobarrier' )
607
608!
609!-- Formats
610100 FORMAT ('#1 ',A,1X,A)
611101 FORMAT (E15.7,1X,E15.7)
612102 FORMAT ('NEXT')
613
614 END SUBROUTINE data_output_profiles
Note: See TracBrowser for help on using the repository browser.