Changeset 3717 for palm/trunk/SOURCE/boundary_conds.f90
- Timestamp:
- Feb 5, 2019 5:21:16 PM (6 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
palm/trunk/SOURCE/boundary_conds.f90
r3655 r3717 25 25 ! ----------------- 26 26 ! $Id$ 27 ! Bugfix, do not set boundary conditions for potential temperature in neutral 28 ! runs. 29 ! 30 ! 3655 2019-01-07 16:51:22Z knoop 27 31 ! OpenACC port for SPEC 28 32 ! … … 235 239 humidity, ibc_pt_b, ibc_pt_t, ibc_q_b, ibc_q_t, ibc_s_b, & 236 240 ibc_s_t, ibc_uv_b, ibc_uv_t, intermediate_timestep_count, & 237 nesting_offline, n udging, ocean_mode, passive_scalar, rans_mode,&238 rans_ tke_e, tsc, salsa, use_cmax241 nesting_offline, neutral, nudging, ocean_mode, passive_scalar, & 242 rans_mode, rans_tke_e, tsc, salsa, use_cmax 239 243 240 244 USE grid_variables, & … … 328 332 !-- the sea surface temperature of the coupled ocean model. 329 333 !-- Dirichlet 330 IF ( ibc_pt_b == 0 ) THEN 331 DO l = 0, 1 332 ! 333 !-- Set kb, for upward-facing surfaces value at topography top (k-1) is set, 334 !-- for downward-facing surfaces at topography bottom (k+1). 335 kb = MERGE( -1, 1, l == 0 ) 336 !$OMP PARALLEL DO PRIVATE( i, j, k ) 337 DO m = 1, bc_h(l)%ns 338 i = bc_h(l)%i(m) 339 j = bc_h(l)%j(m) 340 k = bc_h(l)%k(m) 341 pt_p(k+kb,j,i) = pt(k+kb,j,i) 342 ENDDO 343 ENDDO 344 ! 345 !-- Neumann, zero-gradient 346 ELSEIF ( ibc_pt_b == 1 ) THEN 347 DO l = 0, 1 348 ! 349 !-- Set kb, for upward-facing surfaces value at topography top (k-1) is set, 350 !-- for downward-facing surfaces at topography bottom (k+1). 351 kb = MERGE( -1, 1, l == 0 ) 352 !$OMP PARALLEL DO PRIVATE( i, j, k ) 353 !$ACC PARALLEL LOOP PRIVATE(i, j, k) & 354 !$ACC PRESENT(bc_h, pt_p) 355 DO m = 1, bc_h(l)%ns 356 i = bc_h(l)%i(m) 357 j = bc_h(l)%j(m) 358 k = bc_h(l)%k(m) 359 pt_p(k+kb,j,i) = pt_p(k,j,i) 360 ENDDO 361 ENDDO 362 ENDIF 363 364 ! 365 !-- Temperature at top boundary 366 IF ( ibc_pt_t == 0 ) THEN 367 pt_p(nzt+1,:,:) = pt(nzt+1,:,:) 368 ! 369 !-- In case of nudging adjust top boundary to pt which is 370 !-- read in from NUDGING-DATA 371 IF ( nudging ) THEN 372 pt_p(nzt+1,:,:) = pt_init(nzt+1) 373 ENDIF 374 ELSEIF ( ibc_pt_t == 1 ) THEN 375 pt_p(nzt+1,:,:) = pt_p(nzt,:,:) 376 ELSEIF ( ibc_pt_t == 2 ) THEN 377 !$ACC KERNELS PRESENT(pt_p, dzu) 378 pt_p(nzt+1,:,:) = pt_p(nzt,:,:) + bc_pt_t_val * dzu(nzt+1) 379 !$ACC END KERNELS 334 IF ( .NOT. neutral ) THEN 335 IF ( ibc_pt_b == 0 ) THEN 336 DO l = 0, 1 337 ! 338 !-- Set kb, for upward-facing surfaces value at topography top (k-1) 339 !-- is set, for downward-facing surfaces at topography bottom (k+1). 340 kb = MERGE( -1, 1, l == 0 ) 341 !$OMP PARALLEL DO PRIVATE( i, j, k ) 342 DO m = 1, bc_h(l)%ns 343 i = bc_h(l)%i(m) 344 j = bc_h(l)%j(m) 345 k = bc_h(l)%k(m) 346 pt_p(k+kb,j,i) = pt(k+kb,j,i) 347 ENDDO 348 ENDDO 349 ! 350 !-- Neumann, zero-gradient 351 ELSEIF ( ibc_pt_b == 1 ) THEN 352 DO l = 0, 1 353 ! 354 !-- Set kb, for upward-facing surfaces value at topography top (k-1) 355 !-- is set, for downward-facing surfaces at topography bottom (k+1). 356 kb = MERGE( -1, 1, l == 0 ) 357 !$OMP PARALLEL DO PRIVATE( i, j, k ) 358 !$ACC PARALLEL LOOP PRIVATE(i, j, k) & 359 !$ACC PRESENT(bc_h, pt_p) 360 DO m = 1, bc_h(l)%ns 361 i = bc_h(l)%i(m) 362 j = bc_h(l)%j(m) 363 k = bc_h(l)%k(m) 364 pt_p(k+kb,j,i) = pt_p(k,j,i) 365 ENDDO 366 ENDDO 367 ENDIF 368 369 ! 370 !-- Temperature at top boundary 371 IF ( ibc_pt_t == 0 ) THEN 372 pt_p(nzt+1,:,:) = pt(nzt+1,:,:) 373 ! 374 !-- In case of nudging adjust top boundary to pt which is 375 !-- read in from NUDGING-DATA 376 IF ( nudging ) THEN 377 pt_p(nzt+1,:,:) = pt_init(nzt+1) 378 ENDIF 379 ELSEIF ( ibc_pt_t == 1 ) THEN 380 pt_p(nzt+1,:,:) = pt_p(nzt,:,:) 381 ELSEIF ( ibc_pt_t == 2 ) THEN 382 !$ACC KERNELS PRESENT(pt_p, dzu) 383 pt_p(nzt+1,:,:) = pt_p(nzt,:,:) + bc_pt_t_val * dzu(nzt+1) 384 !$ACC END KERNELS 385 ENDIF 380 386 ENDIF 381 387
Note: See TracChangeset
for help on using the changeset viewer.