[1682] | 1 | !> @file data_output_3d.f90 |
---|
[2000] | 2 | !------------------------------------------------------------------------------! |
---|
[1036] | 3 | ! This file is part of PALM. |
---|
| 4 | ! |
---|
[2000] | 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. |
---|
[1036] | 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 | ! |
---|
[1818] | 17 | ! Copyright 1997-2016 Leibniz Universitaet Hannover |
---|
[2000] | 18 | !------------------------------------------------------------------------------! |
---|
[1036] | 19 | ! |
---|
[254] | 20 | ! Current revisions: |
---|
[1106] | 21 | ! ------------------ |
---|
[1961] | 22 | ! |
---|
[2012] | 23 | ! |
---|
[1552] | 24 | ! Former revisions: |
---|
| 25 | ! ----------------- |
---|
| 26 | ! $Id: data_output_3d.f90 2012 2016-09-19 17:31:38Z lvollmer $ |
---|
| 27 | ! |
---|
[2012] | 28 | ! 2011 2016-09-19 17:29:57Z kanani |
---|
| 29 | ! Flag urban_surface is now defined in module control_parameters, |
---|
| 30 | ! changed prefix for urban surface model output to "usm_", |
---|
| 31 | ! introduced control parameter varnamelength for LEN of trimvar. |
---|
| 32 | ! |
---|
[2008] | 33 | ! 2007 2016-08-24 15:47:17Z kanani |
---|
| 34 | ! Added support for new urban surface model (temporary modifications of |
---|
| 35 | ! SELECT CASE ( ) necessary, see variable trimvar) |
---|
| 36 | ! |
---|
[2001] | 37 | ! 2000 2016-08-20 18:09:15Z knoop |
---|
| 38 | ! Forced header and separation lines into 80 columns |
---|
| 39 | ! |
---|
[1981] | 40 | ! 1980 2016-07-29 15:51:57Z suehring |
---|
| 41 | ! Bugfix, in order to steer user-defined output, setting flag found explicitly |
---|
| 42 | ! to .F. |
---|
| 43 | ! |
---|
[1977] | 44 | ! 1976 2016-07-27 13:28:04Z maronga |
---|
| 45 | ! Output of radiation quantities is now done directly in the respective module |
---|
| 46 | ! |
---|
[1973] | 47 | ! 1972 2016-07-26 07:52:02Z maronga |
---|
| 48 | ! Output of land surface quantities is now done directly in the respective module. |
---|
| 49 | ! Unnecessary directive __parallel removed. |
---|
| 50 | ! |
---|
[1961] | 51 | ! 1960 2016-07-12 16:34:24Z suehring |
---|
| 52 | ! Scalar surface flux added |
---|
| 53 | ! |
---|
[1851] | 54 | ! 1849 2016-04-08 11:33:18Z hoffmann |
---|
| 55 | ! prr moved to arrays_3d |
---|
| 56 | ! |
---|
[1823] | 57 | ! 1822 2016-04-07 07:49:42Z hoffmann |
---|
| 58 | ! prr vertical dimensions set to nzb_do to nzt_do. Unused variables deleted. |
---|
| 59 | ! |
---|
[1809] | 60 | ! 1808 2016-04-05 19:44:00Z raasch |
---|
| 61 | ! test output removed |
---|
| 62 | ! |
---|
[1784] | 63 | ! 1783 2016-03-06 18:36:17Z raasch |
---|
| 64 | ! name change of netcdf routines and module + related changes |
---|
| 65 | ! |
---|
[1746] | 66 | ! 1745 2016-02-05 13:06:51Z gronemeier |
---|
| 67 | ! Bugfix: test if time axis limit exceeds moved to point after call of check_open |
---|
| 68 | ! |
---|
[1692] | 69 | ! 1691 2015-10-26 16:17:44Z maronga |
---|
| 70 | ! Added output of radiative heating rates for RRTMG |
---|
| 71 | ! |
---|
[1683] | 72 | ! 1682 2015-10-07 23:56:08Z knoop |
---|
| 73 | ! Code annotations made doxygen readable |
---|
| 74 | ! |
---|
[1586] | 75 | ! 1585 2015-04-30 07:05:52Z maronga |
---|
| 76 | ! Added support for RRTMG |
---|
| 77 | ! |
---|
[1552] | 78 | ! 1551 2015-03-03 14:18:16Z maronga |
---|
[1551] | 79 | ! Added suppport for land surface model and radiation model output. In the course |
---|
| 80 | ! of this action, the limits for vertical loops have been changed (from nzb and |
---|
| 81 | ! nzt+1 to nzb_do and nzt_do, respectively in order to allow soil model output). |
---|
| 82 | ! Moreover, a new vertical grid zs was introduced. |
---|
[1329] | 83 | ! |
---|
[1360] | 84 | ! 1359 2014-04-11 17:15:14Z hoffmann |
---|
| 85 | ! New particle structure integrated. |
---|
| 86 | ! |
---|
[1354] | 87 | ! 1353 2014-04-08 15:21:23Z heinze |
---|
| 88 | ! REAL constants provided with KIND-attribute |
---|
| 89 | ! |
---|
[1329] | 90 | ! 1327 2014-03-21 11:00:16Z raasch |
---|
| 91 | ! parts concerning avs output removed, |
---|
| 92 | ! -netcdf output queries |
---|
| 93 | ! |
---|
[1321] | 94 | ! 1320 2014-03-20 08:40:49Z raasch |
---|
[1320] | 95 | ! ONLY-attribute added to USE-statements, |
---|
| 96 | ! kind-parameters added to all INTEGER and REAL declaration statements, |
---|
| 97 | ! kinds are defined in new module kinds, |
---|
| 98 | ! old module precision_kind is removed, |
---|
| 99 | ! revision history before 2012 removed, |
---|
| 100 | ! comment fields (!:) to be used for variable explanations added to |
---|
| 101 | ! all variable declaration statements |
---|
[674] | 102 | ! |
---|
[1319] | 103 | ! 1318 2014-03-17 13:35:16Z raasch |
---|
| 104 | ! barrier argument removed from cpu_log, |
---|
| 105 | ! module interfaces removed |
---|
| 106 | ! |
---|
[1309] | 107 | ! 1308 2014-03-13 14:58:42Z fricke |
---|
| 108 | ! Check, if the limit of the time dimension is exceeded for parallel output |
---|
| 109 | ! To increase the performance for parallel output, the following is done: |
---|
| 110 | ! - Update of time axis is only done by PE0 |
---|
| 111 | ! |
---|
[1245] | 112 | ! 1244 2013-10-31 08:16:56Z raasch |
---|
| 113 | ! Bugfix for index bounds in case of 3d-parallel output |
---|
| 114 | ! |
---|
[1116] | 115 | ! 1115 2013-03-26 18:16:16Z hoffmann |
---|
| 116 | ! ql is calculated by calc_liquid_water_content |
---|
| 117 | ! |
---|
[1107] | 118 | ! 1106 2013-03-04 05:31:38Z raasch |
---|
| 119 | ! array_kind renamed precision_kind |
---|
| 120 | ! |
---|
[1077] | 121 | ! 1076 2012-12-05 08:30:18Z hoffmann |
---|
| 122 | ! Bugfix in output of ql |
---|
| 123 | ! |
---|
[1054] | 124 | ! 1053 2012-11-13 17:11:03Z hoffmann |
---|
| 125 | ! +nr, qr, prr, qc and averaged quantities |
---|
| 126 | ! |
---|
[1037] | 127 | ! 1036 2012-10-22 13:43:42Z raasch |
---|
| 128 | ! code put under GPL (PALM 3.9) |
---|
| 129 | ! |
---|
[1035] | 130 | ! 1031 2012-10-19 14:35:30Z raasch |
---|
| 131 | ! netCDF4 without parallel file support implemented |
---|
| 132 | ! |
---|
[1008] | 133 | ! 1007 2012-09-19 14:30:36Z franke |
---|
| 134 | ! Bugfix: missing calculation of ql_vp added |
---|
| 135 | ! |
---|
[1] | 136 | ! Revision 1.1 1997/09/03 06:29:36 raasch |
---|
| 137 | ! Initial revision |
---|
| 138 | ! |
---|
| 139 | ! |
---|
| 140 | ! Description: |
---|
| 141 | ! ------------ |
---|
[1682] | 142 | !> Output of the 3D-arrays in netCDF and/or AVS format. |
---|
[1] | 143 | !------------------------------------------------------------------------------! |
---|
[1682] | 144 | SUBROUTINE data_output_3d( av ) |
---|
| 145 | |
---|
[1] | 146 | |
---|
[1320] | 147 | USE arrays_3d, & |
---|
[1960] | 148 | ONLY: e, nr, p, pt, prr, q, qc, ql, ql_c, ql_v, qr, rho, s, sa, tend, & |
---|
| 149 | u, v, vpt, w |
---|
[1320] | 150 | |
---|
[1] | 151 | USE averaging |
---|
[1320] | 152 | |
---|
| 153 | USE cloud_parameters, & |
---|
[1849] | 154 | ONLY: l_d_cp, pt_d_t |
---|
[1320] | 155 | |
---|
| 156 | USE control_parameters, & |
---|
[1822] | 157 | ONLY: cloud_physics, do3d, do3d_no, do3d_time_count, io_blocks, & |
---|
| 158 | io_group, message_string, ntdim_3d, nz_do3d, psolver, & |
---|
[2011] | 159 | simulated_time, time_since_reference_point, urban_surface, & |
---|
| 160 | varnamelength |
---|
[1320] | 161 | |
---|
| 162 | USE cpulog, & |
---|
| 163 | ONLY: log_point, cpu_log |
---|
| 164 | |
---|
| 165 | USE indices, & |
---|
[1822] | 166 | ONLY: nbgp, nx, nxl, nxlg, nxr, nxrg, ny, nyn, nyng, nys, nysg, nzb |
---|
[1320] | 167 | |
---|
| 168 | USE kinds |
---|
| 169 | |
---|
[1551] | 170 | USE land_surface_model_mod, & |
---|
[1972] | 171 | ONLY: land_surface, lsm_data_output_3d, nzb_soil, nzt_soil |
---|
[1551] | 172 | |
---|
[1783] | 173 | #if defined( __netcdf ) |
---|
| 174 | USE NETCDF |
---|
| 175 | #endif |
---|
| 176 | |
---|
| 177 | USE netcdf_interface, & |
---|
| 178 | ONLY: id_set_3d, id_var_do3d, id_var_time_3d, nc_stat, & |
---|
| 179 | netcdf_data_format, netcdf_handle_error |
---|
[1320] | 180 | |
---|
| 181 | USE particle_attributes, & |
---|
[1359] | 182 | ONLY: grid_particles, number_of_particles, particles, & |
---|
| 183 | particle_advection_start, prt_count |
---|
[1320] | 184 | |
---|
[1] | 185 | USE pegrid |
---|
| 186 | |
---|
[1585] | 187 | USE radiation_model_mod, & |
---|
[1976] | 188 | ONLY: radiation, radiation_data_output_3d |
---|
[1585] | 189 | |
---|
[2007] | 190 | USE urban_surface_mod, & |
---|
[2011] | 191 | ONLY: nzub, nzut, usm_data_output_3d |
---|
[1585] | 192 | |
---|
[2007] | 193 | |
---|
[1] | 194 | IMPLICIT NONE |
---|
| 195 | |
---|
[1682] | 196 | INTEGER(iwp) :: av !< |
---|
| 197 | INTEGER(iwp) :: i !< |
---|
| 198 | INTEGER(iwp) :: if !< |
---|
| 199 | INTEGER(iwp) :: j !< |
---|
| 200 | INTEGER(iwp) :: k !< |
---|
| 201 | INTEGER(iwp) :: n !< |
---|
| 202 | INTEGER(iwp) :: nzb_do !< vertical lower limit for data output |
---|
| 203 | INTEGER(iwp) :: nzt_do !< vertical upper limit for data output |
---|
[1] | 204 | |
---|
[1682] | 205 | LOGICAL :: found !< |
---|
| 206 | LOGICAL :: resorted !< |
---|
[1] | 207 | |
---|
[1682] | 208 | REAL(wp) :: mean_r !< |
---|
| 209 | REAL(wp) :: s_r2 !< |
---|
| 210 | REAL(wp) :: s_r3 !< |
---|
[1] | 211 | |
---|
[1682] | 212 | REAL(sp), DIMENSION(:,:,:), ALLOCATABLE :: local_pf !< |
---|
[1] | 213 | |
---|
[1682] | 214 | REAL(wp), DIMENSION(:,:,:), POINTER :: to_be_resorted !< |
---|
[1] | 215 | |
---|
[2011] | 216 | CHARACTER (LEN=varnamelength) :: trimvar !< TRIM of output-variable string |
---|
[2007] | 217 | |
---|
[1] | 218 | ! |
---|
| 219 | !-- Return, if nothing to output |
---|
| 220 | IF ( do3d_no(av) == 0 ) RETURN |
---|
| 221 | |
---|
| 222 | CALL cpu_log (log_point(14),'data_output_3d','start') |
---|
| 223 | |
---|
| 224 | ! |
---|
| 225 | !-- Open output file. |
---|
| 226 | !-- Also creates coordinate and fld-file for AVS. |
---|
[1031] | 227 | !-- For classic or 64bit netCDF output or output of other (old) data formats, |
---|
[493] | 228 | !-- for a run on more than one PE, each PE opens its own file and |
---|
[1] | 229 | !-- writes the data of its subdomain in binary format (regardless of the format |
---|
| 230 | !-- the user has requested). After the run, these files are combined to one |
---|
| 231 | !-- file by combine_plot_fields in the format requested by the user (netcdf |
---|
[493] | 232 | !-- and/or avs). |
---|
[1031] | 233 | !-- For netCDF4/HDF5 output, data is written in parallel into one file. |
---|
[1327] | 234 | IF ( netcdf_data_format < 5 ) THEN |
---|
| 235 | CALL check_open( 30 ) |
---|
| 236 | IF ( myid == 0 ) CALL check_open( 106+av*10 ) |
---|
[493] | 237 | ELSE |
---|
[1327] | 238 | CALL check_open( 106+av*10 ) |
---|
[493] | 239 | ENDIF |
---|
[1] | 240 | |
---|
| 241 | ! |
---|
[1745] | 242 | !-- For parallel netcdf output the time axis must be limited. Return, if this |
---|
| 243 | !-- limit is exceeded. This could be the case, if the simulated time exceeds |
---|
| 244 | !-- the given end time by the length of the given output interval. |
---|
| 245 | IF ( netcdf_data_format > 4 ) THEN |
---|
| 246 | IF ( do3d_time_count(av) + 1 > ntdim_3d(av) ) THEN |
---|
| 247 | WRITE ( message_string, * ) 'Output of 3d data is not given at t=', & |
---|
| 248 | simulated_time, '&because the maximum ', & |
---|
| 249 | 'number of output time levels is ', & |
---|
| 250 | 'exceeded.' |
---|
| 251 | CALL message( 'data_output_3d', 'PA0387', 0, 1, 0, 6, 0 ) |
---|
| 252 | CALL cpu_log( log_point(14), 'data_output_3d', 'stop' ) |
---|
| 253 | RETURN |
---|
| 254 | ENDIF |
---|
| 255 | ENDIF |
---|
| 256 | |
---|
| 257 | ! |
---|
[1031] | 258 | !-- Update the netCDF time axis |
---|
[1308] | 259 | !-- In case of parallel output, this is only done by PE0 to increase the |
---|
| 260 | !-- performance. |
---|
[1] | 261 | #if defined( __netcdf ) |
---|
[1308] | 262 | do3d_time_count(av) = do3d_time_count(av) + 1 |
---|
| 263 | IF ( myid == 0 ) THEN |
---|
[1327] | 264 | nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_time_3d(av), & |
---|
| 265 | (/ time_since_reference_point /), & |
---|
| 266 | start = (/ do3d_time_count(av) /), & |
---|
| 267 | count = (/ 1 /) ) |
---|
[1783] | 268 | CALL netcdf_handle_error( 'data_output_3d', 376 ) |
---|
[1] | 269 | ENDIF |
---|
| 270 | #endif |
---|
| 271 | |
---|
| 272 | ! |
---|
| 273 | !-- Loop over all variables to be written. |
---|
| 274 | if = 1 |
---|
| 275 | |
---|
| 276 | DO WHILE ( do3d(av,if)(1:1) /= ' ' ) |
---|
[2007] | 277 | |
---|
[1] | 278 | ! |
---|
[2007] | 279 | !-- Temporary solution to account for data output within the new urban |
---|
| 280 | !-- surface model (urban_surface_mod.f90), see also SELECT CASE ( trimvar ). |
---|
[1] | 281 | !-- Store the array chosen on the temporary array. |
---|
[2007] | 282 | trimvar = TRIM( do3d(av,if) ) |
---|
[2011] | 283 | IF ( urban_surface .AND. trimvar(1:4) == 'usm_' ) THEN |
---|
[2007] | 284 | trimvar = 'usm_output' |
---|
| 285 | resorted = .TRUE. |
---|
| 286 | nzb_do = nzub |
---|
| 287 | nzt_do = nzut |
---|
| 288 | ELSE |
---|
| 289 | resorted = .FALSE. |
---|
| 290 | nzb_do = nzb |
---|
| 291 | nzt_do = nz_do3d |
---|
| 292 | ENDIF |
---|
[1551] | 293 | ! |
---|
[1980] | 294 | !-- Set flag to steer output of radiation, land-surface, or user-defined |
---|
| 295 | !-- quantities |
---|
| 296 | found = .FALSE. |
---|
| 297 | ! |
---|
[1551] | 298 | !-- Allocate a temporary array with the desired output dimensions. |
---|
| 299 | ALLOCATE( local_pf(nxlg:nxrg,nysg:nyng,nzb_do:nzt_do) ) |
---|
| 300 | |
---|
[2007] | 301 | SELECT CASE ( trimvar ) |
---|
[1] | 302 | |
---|
| 303 | CASE ( 'e' ) |
---|
| 304 | IF ( av == 0 ) THEN |
---|
| 305 | to_be_resorted => e |
---|
| 306 | ELSE |
---|
| 307 | to_be_resorted => e_av |
---|
| 308 | ENDIF |
---|
| 309 | |
---|
[771] | 310 | CASE ( 'lpt' ) |
---|
| 311 | IF ( av == 0 ) THEN |
---|
| 312 | to_be_resorted => pt |
---|
| 313 | ELSE |
---|
| 314 | to_be_resorted => lpt_av |
---|
| 315 | ENDIF |
---|
| 316 | |
---|
[1053] | 317 | CASE ( 'nr' ) |
---|
| 318 | IF ( av == 0 ) THEN |
---|
| 319 | to_be_resorted => nr |
---|
| 320 | ELSE |
---|
| 321 | to_be_resorted => nr_av |
---|
| 322 | ENDIF |
---|
| 323 | |
---|
[1] | 324 | CASE ( 'p' ) |
---|
| 325 | IF ( av == 0 ) THEN |
---|
[727] | 326 | IF ( psolver /= 'sor' ) CALL exchange_horiz( p, nbgp ) |
---|
[1] | 327 | to_be_resorted => p |
---|
| 328 | ELSE |
---|
[727] | 329 | IF ( psolver /= 'sor' ) CALL exchange_horiz( p_av, nbgp ) |
---|
[1] | 330 | to_be_resorted => p_av |
---|
| 331 | ENDIF |
---|
| 332 | |
---|
| 333 | CASE ( 'pc' ) ! particle concentration (requires ghostpoint exchange) |
---|
| 334 | IF ( av == 0 ) THEN |
---|
[1359] | 335 | IF ( simulated_time >= particle_advection_start ) THEN |
---|
| 336 | tend = prt_count |
---|
| 337 | CALL exchange_horiz( tend, nbgp ) |
---|
| 338 | ELSE |
---|
| 339 | tend = 0.0_wp |
---|
| 340 | ENDIF |
---|
[667] | 341 | DO i = nxlg, nxrg |
---|
| 342 | DO j = nysg, nyng |
---|
[1551] | 343 | DO k = nzb_do, nzt_do |
---|
[1] | 344 | local_pf(i,j,k) = tend(k,j,i) |
---|
| 345 | ENDDO |
---|
| 346 | ENDDO |
---|
| 347 | ENDDO |
---|
| 348 | resorted = .TRUE. |
---|
| 349 | ELSE |
---|
[667] | 350 | CALL exchange_horiz( pc_av, nbgp ) |
---|
[1] | 351 | to_be_resorted => pc_av |
---|
| 352 | ENDIF |
---|
| 353 | |
---|
[1359] | 354 | CASE ( 'pr' ) ! mean particle radius (effective radius) |
---|
[1] | 355 | IF ( av == 0 ) THEN |
---|
[1359] | 356 | IF ( simulated_time >= particle_advection_start ) THEN |
---|
| 357 | DO i = nxl, nxr |
---|
| 358 | DO j = nys, nyn |
---|
[1551] | 359 | DO k = nzb_do, nzt_do |
---|
[1359] | 360 | number_of_particles = prt_count(k,j,i) |
---|
| 361 | IF (number_of_particles <= 0) CYCLE |
---|
| 362 | particles => grid_particles(k,j,i)%particles(1:number_of_particles) |
---|
| 363 | s_r2 = 0.0_wp |
---|
| 364 | s_r3 = 0.0_wp |
---|
| 365 | DO n = 1, number_of_particles |
---|
| 366 | IF ( particles(n)%particle_mask ) THEN |
---|
| 367 | s_r2 = s_r2 + particles(n)%radius**2 * & |
---|
| 368 | particles(n)%weight_factor |
---|
| 369 | s_r3 = s_r3 + particles(n)%radius**3 * & |
---|
| 370 | particles(n)%weight_factor |
---|
| 371 | ENDIF |
---|
| 372 | ENDDO |
---|
| 373 | IF ( s_r2 > 0.0_wp ) THEN |
---|
| 374 | mean_r = s_r3 / s_r2 |
---|
| 375 | ELSE |
---|
| 376 | mean_r = 0.0_wp |
---|
| 377 | ENDIF |
---|
| 378 | tend(k,j,i) = mean_r |
---|
[1] | 379 | ENDDO |
---|
| 380 | ENDDO |
---|
| 381 | ENDDO |
---|
[1359] | 382 | CALL exchange_horiz( tend, nbgp ) |
---|
| 383 | ELSE |
---|
| 384 | tend = 0.0_wp |
---|
| 385 | ENDIF |
---|
[667] | 386 | DO i = nxlg, nxrg |
---|
| 387 | DO j = nysg, nyng |
---|
[1551] | 388 | DO k = nzb_do, nzt_do |
---|
[1] | 389 | local_pf(i,j,k) = tend(k,j,i) |
---|
| 390 | ENDDO |
---|
| 391 | ENDDO |
---|
| 392 | ENDDO |
---|
| 393 | resorted = .TRUE. |
---|
| 394 | ELSE |
---|
[667] | 395 | CALL exchange_horiz( pr_av, nbgp ) |
---|
[1] | 396 | to_be_resorted => pr_av |
---|
| 397 | ENDIF |
---|
| 398 | |
---|
[1053] | 399 | CASE ( 'prr' ) |
---|
| 400 | IF ( av == 0 ) THEN |
---|
| 401 | CALL exchange_horiz( prr, nbgp ) |
---|
| 402 | DO i = nxlg, nxrg |
---|
| 403 | DO j = nysg, nyng |
---|
[1822] | 404 | DO k = nzb_do, nzt_do |
---|
[1053] | 405 | local_pf(i,j,k) = prr(k,j,i) |
---|
| 406 | ENDDO |
---|
| 407 | ENDDO |
---|
| 408 | ENDDO |
---|
| 409 | ELSE |
---|
| 410 | CALL exchange_horiz( prr_av, nbgp ) |
---|
| 411 | DO i = nxlg, nxrg |
---|
| 412 | DO j = nysg, nyng |
---|
[1822] | 413 | DO k = nzb_do, nzt_do |
---|
[1053] | 414 | local_pf(i,j,k) = prr_av(k,j,i) |
---|
| 415 | ENDDO |
---|
| 416 | ENDDO |
---|
| 417 | ENDDO |
---|
| 418 | ENDIF |
---|
| 419 | resorted = .TRUE. |
---|
| 420 | |
---|
[1] | 421 | CASE ( 'pt' ) |
---|
| 422 | IF ( av == 0 ) THEN |
---|
| 423 | IF ( .NOT. cloud_physics ) THEN |
---|
| 424 | to_be_resorted => pt |
---|
| 425 | ELSE |
---|
[667] | 426 | DO i = nxlg, nxrg |
---|
| 427 | DO j = nysg, nyng |
---|
[1551] | 428 | DO k = nzb_do, nzt_do |
---|
[1320] | 429 | local_pf(i,j,k) = pt(k,j,i) + l_d_cp * & |
---|
| 430 | pt_d_t(k) * & |
---|
[1] | 431 | ql(k,j,i) |
---|
| 432 | ENDDO |
---|
| 433 | ENDDO |
---|
| 434 | ENDDO |
---|
| 435 | resorted = .TRUE. |
---|
| 436 | ENDIF |
---|
| 437 | ELSE |
---|
| 438 | to_be_resorted => pt_av |
---|
| 439 | ENDIF |
---|
| 440 | |
---|
| 441 | CASE ( 'q' ) |
---|
| 442 | IF ( av == 0 ) THEN |
---|
| 443 | to_be_resorted => q |
---|
| 444 | ELSE |
---|
| 445 | to_be_resorted => q_av |
---|
| 446 | ENDIF |
---|
[691] | 447 | |
---|
[1053] | 448 | CASE ( 'qc' ) |
---|
[1] | 449 | IF ( av == 0 ) THEN |
---|
[1115] | 450 | to_be_resorted => qc |
---|
[1] | 451 | ELSE |
---|
[1115] | 452 | to_be_resorted => qc_av |
---|
[1] | 453 | ENDIF |
---|
| 454 | |
---|
[1053] | 455 | CASE ( 'ql' ) |
---|
| 456 | IF ( av == 0 ) THEN |
---|
[1115] | 457 | to_be_resorted => ql |
---|
[1053] | 458 | ELSE |
---|
[1115] | 459 | to_be_resorted => ql_av |
---|
[1053] | 460 | ENDIF |
---|
| 461 | |
---|
[1] | 462 | CASE ( 'ql_c' ) |
---|
| 463 | IF ( av == 0 ) THEN |
---|
| 464 | to_be_resorted => ql_c |
---|
| 465 | ELSE |
---|
| 466 | to_be_resorted => ql_c_av |
---|
| 467 | ENDIF |
---|
| 468 | |
---|
| 469 | CASE ( 'ql_v' ) |
---|
| 470 | IF ( av == 0 ) THEN |
---|
| 471 | to_be_resorted => ql_v |
---|
| 472 | ELSE |
---|
| 473 | to_be_resorted => ql_v_av |
---|
| 474 | ENDIF |
---|
| 475 | |
---|
| 476 | CASE ( 'ql_vp' ) |
---|
| 477 | IF ( av == 0 ) THEN |
---|
[1359] | 478 | IF ( simulated_time >= particle_advection_start ) THEN |
---|
| 479 | DO i = nxl, nxr |
---|
| 480 | DO j = nys, nyn |
---|
[1551] | 481 | DO k = nzb_do, nzt_do |
---|
[1359] | 482 | number_of_particles = prt_count(k,j,i) |
---|
| 483 | IF (number_of_particles <= 0) CYCLE |
---|
| 484 | particles => grid_particles(k,j,i)%particles(1:number_of_particles) |
---|
| 485 | DO n = 1, number_of_particles |
---|
| 486 | IF ( particles(n)%particle_mask ) THEN |
---|
| 487 | tend(k,j,i) = tend(k,j,i) + & |
---|
| 488 | particles(n)%weight_factor / & |
---|
| 489 | prt_count(k,j,i) |
---|
| 490 | ENDIF |
---|
| 491 | ENDDO |
---|
[1007] | 492 | ENDDO |
---|
| 493 | ENDDO |
---|
| 494 | ENDDO |
---|
[1359] | 495 | CALL exchange_horiz( tend, nbgp ) |
---|
| 496 | ELSE |
---|
| 497 | tend = 0.0_wp |
---|
| 498 | ENDIF |
---|
[1007] | 499 | DO i = nxlg, nxrg |
---|
| 500 | DO j = nysg, nyng |
---|
[1551] | 501 | DO k = nzb_do, nzt_do |
---|
[1007] | 502 | local_pf(i,j,k) = tend(k,j,i) |
---|
| 503 | ENDDO |
---|
| 504 | ENDDO |
---|
| 505 | ENDDO |
---|
| 506 | resorted = .TRUE. |
---|
[1] | 507 | ELSE |
---|
[1007] | 508 | CALL exchange_horiz( ql_vp_av, nbgp ) |
---|
[1] | 509 | to_be_resorted => ql_vp_av |
---|
| 510 | ENDIF |
---|
| 511 | |
---|
[1053] | 512 | CASE ( 'qr' ) |
---|
| 513 | IF ( av == 0 ) THEN |
---|
| 514 | to_be_resorted => qr |
---|
| 515 | ELSE |
---|
| 516 | to_be_resorted => qr_av |
---|
| 517 | ENDIF |
---|
| 518 | |
---|
[1] | 519 | CASE ( 'qv' ) |
---|
| 520 | IF ( av == 0 ) THEN |
---|
[667] | 521 | DO i = nxlg, nxrg |
---|
| 522 | DO j = nysg, nyng |
---|
[1551] | 523 | DO k = nzb_do, nzt_do |
---|
[1] | 524 | local_pf(i,j,k) = q(k,j,i) - ql(k,j,i) |
---|
| 525 | ENDDO |
---|
| 526 | ENDDO |
---|
| 527 | ENDDO |
---|
| 528 | resorted = .TRUE. |
---|
| 529 | ELSE |
---|
| 530 | to_be_resorted => qv_av |
---|
| 531 | ENDIF |
---|
| 532 | |
---|
[96] | 533 | CASE ( 'rho' ) |
---|
| 534 | IF ( av == 0 ) THEN |
---|
| 535 | to_be_resorted => rho |
---|
| 536 | ELSE |
---|
| 537 | to_be_resorted => rho_av |
---|
| 538 | ENDIF |
---|
[691] | 539 | |
---|
[1] | 540 | CASE ( 's' ) |
---|
| 541 | IF ( av == 0 ) THEN |
---|
[1960] | 542 | to_be_resorted => s |
---|
[1] | 543 | ELSE |
---|
[355] | 544 | to_be_resorted => s_av |
---|
[1] | 545 | ENDIF |
---|
[691] | 546 | |
---|
[96] | 547 | CASE ( 'sa' ) |
---|
| 548 | IF ( av == 0 ) THEN |
---|
| 549 | to_be_resorted => sa |
---|
| 550 | ELSE |
---|
| 551 | to_be_resorted => sa_av |
---|
| 552 | ENDIF |
---|
[691] | 553 | |
---|
[1] | 554 | CASE ( 'u' ) |
---|
| 555 | IF ( av == 0 ) THEN |
---|
| 556 | to_be_resorted => u |
---|
| 557 | ELSE |
---|
| 558 | to_be_resorted => u_av |
---|
| 559 | ENDIF |
---|
| 560 | |
---|
| 561 | CASE ( 'v' ) |
---|
| 562 | IF ( av == 0 ) THEN |
---|
| 563 | to_be_resorted => v |
---|
| 564 | ELSE |
---|
| 565 | to_be_resorted => v_av |
---|
| 566 | ENDIF |
---|
| 567 | |
---|
| 568 | CASE ( 'vpt' ) |
---|
| 569 | IF ( av == 0 ) THEN |
---|
| 570 | to_be_resorted => vpt |
---|
| 571 | ELSE |
---|
| 572 | to_be_resorted => vpt_av |
---|
| 573 | ENDIF |
---|
| 574 | |
---|
| 575 | CASE ( 'w' ) |
---|
| 576 | IF ( av == 0 ) THEN |
---|
| 577 | to_be_resorted => w |
---|
| 578 | ELSE |
---|
| 579 | to_be_resorted => w_av |
---|
| 580 | ENDIF |
---|
[2007] | 581 | ! |
---|
| 582 | !-- Block of urban surface model outputs |
---|
| 583 | CASE ( 'usm_output' ) |
---|
| 584 | CALL usm_data_output_3d( av, do3d(av,if), found, local_pf, & |
---|
| 585 | nzb_do, nzt_do ) |
---|
[1] | 586 | |
---|
| 587 | CASE DEFAULT |
---|
[1972] | 588 | |
---|
[1] | 589 | ! |
---|
[1972] | 590 | !-- Land surface quantity |
---|
| 591 | IF ( land_surface ) THEN |
---|
| 592 | ! |
---|
| 593 | !-- For soil model quantities, it is required to re-allocate local_pf |
---|
| 594 | nzb_do = nzb_soil |
---|
| 595 | nzt_do = nzt_soil |
---|
| 596 | |
---|
| 597 | DEALLOCATE ( local_pf ) |
---|
| 598 | ALLOCATE( local_pf(nxlg:nxrg,nysg:nyng,nzb_do:nzt_do) ) |
---|
| 599 | |
---|
| 600 | CALL lsm_data_output_3d( av, do3d(av,if), found, local_pf ) |
---|
| 601 | resorted = .TRUE. |
---|
[1976] | 602 | |
---|
| 603 | ! |
---|
| 604 | !-- If no soil model variable was found, re-allocate local_pf |
---|
| 605 | IF ( .NOT. found ) THEN |
---|
| 606 | nzb_do = nzb |
---|
| 607 | nzt_do = nz_do3d |
---|
| 608 | |
---|
| 609 | DEALLOCATE ( local_pf ) |
---|
| 610 | ALLOCATE( local_pf(nxlg:nxrg,nysg:nyng,nzb_do:nzt_do) ) |
---|
| 611 | ENDIF |
---|
| 612 | |
---|
[1972] | 613 | ENDIF |
---|
| 614 | |
---|
| 615 | ! |
---|
[1976] | 616 | !-- Radiation quantity |
---|
| 617 | IF ( .NOT. found .AND. radiation ) THEN |
---|
| 618 | CALL radiation_data_output_3d( av, do3d(av,if), found, & |
---|
| 619 | local_pf ) |
---|
| 620 | resorted = .TRUE. |
---|
| 621 | ENDIF |
---|
| 622 | |
---|
| 623 | ! |
---|
[1] | 624 | !-- User defined quantity |
---|
[1972] | 625 | IF ( .NOT. found ) THEN |
---|
| 626 | CALL user_data_output_3d( av, do3d(av,if), found, local_pf, & |
---|
| 627 | nzb_do, nzt_do ) |
---|
| 628 | resorted = .TRUE. |
---|
| 629 | ENDIF |
---|
[1] | 630 | |
---|
[254] | 631 | IF ( .NOT. found ) THEN |
---|
[1320] | 632 | message_string = 'no output available for: ' // & |
---|
[274] | 633 | TRIM( do3d(av,if) ) |
---|
[254] | 634 | CALL message( 'data_output_3d', 'PA0182', 0, 0, 0, 6, 0 ) |
---|
[1] | 635 | ENDIF |
---|
| 636 | |
---|
| 637 | END SELECT |
---|
| 638 | |
---|
| 639 | ! |
---|
| 640 | !-- Resort the array to be output, if not done above |
---|
| 641 | IF ( .NOT. resorted ) THEN |
---|
[667] | 642 | DO i = nxlg, nxrg |
---|
| 643 | DO j = nysg, nyng |
---|
[1551] | 644 | DO k = nzb_do, nzt_do |
---|
[1] | 645 | local_pf(i,j,k) = to_be_resorted(k,j,i) |
---|
| 646 | ENDDO |
---|
| 647 | ENDDO |
---|
| 648 | ENDDO |
---|
| 649 | ENDIF |
---|
| 650 | |
---|
| 651 | ! |
---|
[1327] | 652 | !-- Output of the 3D-array |
---|
| 653 | #if defined( __parallel ) |
---|
| 654 | IF ( netcdf_data_format < 5 ) THEN |
---|
[1] | 655 | ! |
---|
[1327] | 656 | !-- Non-parallel netCDF output. Data is output in parallel in |
---|
| 657 | !-- FORTRAN binary format here, and later collected into one file by |
---|
| 658 | !-- combine_plot_fields |
---|
| 659 | IF ( myid == 0 ) THEN |
---|
| 660 | WRITE ( 30 ) time_since_reference_point, & |
---|
| 661 | do3d_time_count(av), av |
---|
[1] | 662 | ENDIF |
---|
[1327] | 663 | DO i = 0, io_blocks-1 |
---|
| 664 | IF ( i == io_group ) THEN |
---|
[1551] | 665 | WRITE ( 30 ) nxlg, nxrg, nysg, nyng, nzb_do, nzt_do |
---|
| 666 | WRITE ( 30 ) local_pf(:,:,nzb_do:nzt_do) |
---|
[1327] | 667 | ENDIF |
---|
[1972] | 668 | |
---|
[1327] | 669 | CALL MPI_BARRIER( comm2d, ierr ) |
---|
[1972] | 670 | |
---|
[1327] | 671 | ENDDO |
---|
[559] | 672 | |
---|
[1327] | 673 | ELSE |
---|
[646] | 674 | #if defined( __netcdf ) |
---|
[493] | 675 | ! |
---|
[1327] | 676 | !-- Parallel output in netCDF4/HDF5 format. |
---|
| 677 | !-- Do not output redundant ghost point data except for the |
---|
| 678 | !-- boundaries of the total domain. |
---|
| 679 | IF ( nxr == nx .AND. nyn /= ny ) THEN |
---|
| 680 | nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_do3d(av,if), & |
---|
[1551] | 681 | local_pf(nxl:nxr+1,nys:nyn,nzb_do:nzt_do), & |
---|
| 682 | start = (/ nxl+1, nys+1, nzb_do+1, do3d_time_count(av) /), & |
---|
| 683 | count = (/ nxr-nxl+2, nyn-nys+1, nzt_do-nzb_do+1, 1 /) ) |
---|
[1327] | 684 | ELSEIF ( nxr /= nx .AND. nyn == ny ) THEN |
---|
| 685 | nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_do3d(av,if), & |
---|
[1551] | 686 | local_pf(nxl:nxr,nys:nyn+1,nzb_do:nzt_do), & |
---|
| 687 | start = (/ nxl+1, nys+1, nzb_do+1, do3d_time_count(av) /), & |
---|
| 688 | count = (/ nxr-nxl+1, nyn-nys+2, nzt_do-nzb_do+1, 1 /) ) |
---|
[1327] | 689 | ELSEIF ( nxr == nx .AND. nyn == ny ) THEN |
---|
| 690 | nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_do3d(av,if), & |
---|
[1551] | 691 | local_pf(nxl:nxr+1,nys:nyn+1,nzb_do:nzt_do ), & |
---|
| 692 | start = (/ nxl+1, nys+1, nzb_do+1, do3d_time_count(av) /), & |
---|
| 693 | count = (/ nxr-nxl+2, nyn-nys+2, nzt_do-nzb_do+1, 1 /) ) |
---|
[1327] | 694 | ELSE |
---|
| 695 | nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_do3d(av,if), & |
---|
[1551] | 696 | local_pf(nxl:nxr,nys:nyn,nzb_do:nzt_do), & |
---|
| 697 | start = (/ nxl+1, nys+1, nzb_do+1, do3d_time_count(av) /), & |
---|
| 698 | count = (/ nxr-nxl+1, nyn-nys+1, nzt_do-nzb_do+1, 1 /) ) |
---|
[1327] | 699 | ENDIF |
---|
[1783] | 700 | CALL netcdf_handle_error( 'data_output_3d', 386 ) |
---|
[646] | 701 | #endif |
---|
[1327] | 702 | ENDIF |
---|
[1] | 703 | #else |
---|
| 704 | #if defined( __netcdf ) |
---|
[1327] | 705 | nc_stat = NF90_PUT_VAR( id_set_3d(av), id_var_do3d(av,if), & |
---|
[1551] | 706 | local_pf(nxl:nxr+1,nys:nyn+1,nzb_do:nzt_do), & |
---|
[1327] | 707 | start = (/ 1, 1, 1, do3d_time_count(av) /), & |
---|
[1551] | 708 | count = (/ nx+2, ny+2, nzt_do-nzb_do+1, 1 /) ) |
---|
[1783] | 709 | CALL netcdf_handle_error( 'data_output_3d', 446 ) |
---|
[1] | 710 | #endif |
---|
| 711 | #endif |
---|
| 712 | |
---|
| 713 | if = if + 1 |
---|
| 714 | |
---|
| 715 | ! |
---|
[1551] | 716 | !-- Deallocate temporary array |
---|
| 717 | DEALLOCATE ( local_pf ) |
---|
[1] | 718 | |
---|
[1551] | 719 | ENDDO |
---|
[1] | 720 | |
---|
[1318] | 721 | CALL cpu_log( log_point(14), 'data_output_3d', 'stop' ) |
---|
[1] | 722 | |
---|
| 723 | ! |
---|
| 724 | !-- Formats. |
---|
[1320] | 725 | 3300 FORMAT ('variable ',I4,' file=',A,' filetype=unformatted skip=',I12/ & |
---|
[1] | 726 | 'label = ',A,A) |
---|
| 727 | |
---|
| 728 | END SUBROUTINE data_output_3d |
---|