/usr/share/code_saturne/user_examples/cs_user_boundary_conditions-auto_inlet_profile.f90 is in code-saturne-data 4.2.0+repack-1build1.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 | !-------------------------------------------------------------------------------
! Code_Saturne version 4.2.0
! --------------------------
! This file is part of Code_Saturne, a general-purpose CFD tool.
!
! Copyright (C) 1998-2015 EDF S.A.
!
! This program is free software; you can redistribute it and/or modify it under
! the terms of the GNU General Public License as published by the Free Software
! Foundation; either version 2 of the License, or (at your option) any later
! version.
!
! This program is distributed in the hope that it will be useful, but WITHOUT
! ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
! FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
! details.
!
! You should have received a copy of the GNU General Public License along with
! this program; if not, write to the Free Software Foundation, Inc., 51 Franklin
! Street, Fifth Floor, Boston, MA 02110-1301, USA.
!-------------------------------------------------------------------------------
!===============================================================================
! Function:
! ---------
! Example of cs_user_boundary_conditions subroutine.f90 for inlet
! with automatic inlet profile.
! This example assumes the mesh is orthogonal at the inlet.
!
!-------------------------------------------------------------------------------
!-------------------------------------------------------------------------------
! Arguments
!______________________________________________________________________________.
! mode name role !
!______________________________________________________________________________!
!> \param[in] nvar total number of variables
!> \param[in] nscal total number of scalars
!> \param[out] icodcl boundary condition code:
!> - 1 Dirichlet
!> - 2 Radiative outlet
!> - 3 Neumann
!> - 4 sliding and
!> \f$ \vect{u} \cdot \vect{n} = 0 \f$
!> - 5 smooth wall and
!> \f$ \vect{u} \cdot \vect{n} = 0 \f$
!> - 6 rough wall and
!> \f$ \vect{u} \cdot \vect{n} = 0 \f$
!> - 9 free inlet/outlet
!> (input mass flux blocked to 0)
!> - 13 Dirichlet for the advection operator and
!> Neumann for the diffusion operator
!> \param[in] itrifb indirection for boundary faces ordering
!> \param[in,out] itypfb boundary face types
!> \param[out] izfppp boundary face zone number
!> \param[in] dt time step (per cell)
!> \param[in,out] rcodcl boundary condition values:
!> - rcodcl(1) value of the dirichlet
!> - rcodcl(2) value of the exterior exchange
!> coefficient (infinite if no exchange)
!> - rcodcl(3) value flux density
!> (negative if gain) in w/m2 or roughness
!> in m if icodcl=6
!> -# for the velocity \f$ (\mu+\mu_T)
!> \gradt \, \vect{u} \cdot \vect{n} \f$
!> -# for the pressure \f$ \Delta t
!> \grad P \cdot \vect{n} \f$
!> -# for a scalar \f$ cp \left( K +
!> \dfrac{K_T}{\sigma_T} \right)
!> \grad T \cdot \vect{n} \f$
!_______________________________________________________________________________
subroutine cs_user_boundary_conditions &
( nvar , nscal , &
icodcl , itrifb , itypfb , izfppp , &
dt , &
rcodcl )
!===============================================================================
!===============================================================================
! Module files
!===============================================================================
use paramx
use numvar
use optcal
use cstphy
use cstnum
use entsor
use parall
use period
use ihmpre
use ppppar
use ppthch
use coincl
use cpincl
use ppincl
use ppcpfu
use atincl
use atsoil
use ctincl
use elincl
use cs_fuel_incl
use mesh
use field
use turbomachinery
use cs_c_bindings
!===============================================================================
implicit none
! Arguments
integer nvar , nscal
integer icodcl(nfabor,nvarcl)
integer itrifb(nfabor), itypfb(nfabor)
integer izfppp(nfabor)
double precision dt(ncelet)
double precision rcodcl(nfabor,nvarcl,3)
! Local variables
!< [loc_var_dec]
integer ifac, iel, ii, ilelt, nlelt
double precision xustar2, xdh, d2s3, rhomoy
double precision acc(2), fmprsc, fmul, uref2, vnrm
integer, allocatable, dimension(:) :: lstelt, mrkcel
double precision, dimension(:), pointer :: bfpro_rom
double precision, dimension(:), pointer :: cvar_r11, cvar_r22, cvar_r33
double precision, dimension(:), pointer :: cvar_r12, cvar_r23, cvar_r13
double precision, dimension(:), pointer :: cvar_k, cvar_ep, cvar_phi
double precision, dimension(:), pointer :: cvar_omg
double precision, dimension(:), pointer :: cvar_al, cvar_fb
double precision, dimension(:), pointer :: cvar_nusa
double precision, dimension(:), pointer :: cvar_scal
double precision, dimension(:,:), pointer :: cvar_vel
!< [loc_var_dec]
!===============================================================================
! Initialization
!===============================================================================
!< [init]
allocate(lstelt(nfabor)) ! temporary array for boundary faces selection
! Map field arrays
call field_get_val_v(ivarfl(iu), cvar_vel)
call field_get_val_s(ibrom, bfpro_rom)
if (itytur.eq.2) then
call field_get_val_s(ivarfl(ik), cvar_k)
call field_get_val_s(ivarfl(iep), cvar_ep)
elseif (itytur.eq.3) then
call field_get_val_s(ivarfl(ir11), cvar_r11)
call field_get_val_s(ivarfl(ir22), cvar_r22)
call field_get_val_s(ivarfl(ir33), cvar_r33)
call field_get_val_s(ivarfl(ir12), cvar_r12)
call field_get_val_s(ivarfl(ir13), cvar_r13)
call field_get_val_s(ivarfl(ir23), cvar_r23)
call field_get_val_s(ivarfl(iep), cvar_ep)
if (iturb.eq.32) then
call field_get_val_s(ivarfl(ial), cvar_al)
endif
elseif (itytur.eq.5) then
call field_get_val_s(ivarfl(ik), cvar_k)
call field_get_val_s(ivarfl(iep), cvar_ep)
call field_get_val_s(ivarfl(iphi), cvar_phi)
if (iturb.eq.50) then
call field_get_val_s(ivarfl(ifb), cvar_fb)
elseif (iturb.eq.51) then
call field_get_val_s(ivarfl(ial), cvar_al)
endif
elseif (iturb.eq.60) then
call field_get_val_s(ivarfl(ik), cvar_k)
call field_get_val_s(ivarfl(iomg), cvar_omg)
elseif (iturb.eq.70) then
call field_get_val_s(ivarfl(inusa), cvar_nusa)
endif
d2s3 = 2.d0/3.d0
!< [init]
!===============================================================================
! Assign a pseudo-periodic channel type inlet to a set of boundary faces.
! For each subset:
! - use selection criteria to filter boundary faces of a given subset
! - loop on faces from a subset
! - set the boundary condition for each face
! A feedback loop is used so as to progressively reach a state similar
! to that of a periodic channel at the inlet.
!===============================================================================
!< [example_1]
call getfbr('INLET', nlelt, lstelt)
!==========
fmprsc = 1.d0 ! mean prescribed velocity
if (ntcabs.eq.1) then
! For the Rij-EBRSM model (and possibly V2f), we need a non-flat profile,
! so as to ensure turbulent production, and avoid laminarization;
! here, we simply divide the initial velocity by 10 for inlet
! faces adjacent to the wall.
! The loop below assumes wall conditions have been defined first
! (in the GUI, or in this file, before the current test).
if (iturb.eq.32 .or. itytur.eq.5) then
allocate(mrkcel(ncelet))
do iel = 1, ncelet
mrkcel(iel) = 0
enddo
do ifac = 1, nfabor
if (itypfb(ifac) .eq. iparoi) then
iel = ifabor(ifac)
mrkcel(iel) = 1
endif
enddo
endif
do ilelt = 1, nlelt
ifac = lstelt(ilelt)
iel = ifabor(ifac)
itypfb(ifac) = ientre
rcodcl(ifac,iu,1) = - fmprsc * surfbo(1,ifac) / surfbn(ifac)
rcodcl(ifac,iv,1) = - fmprsc * surfbo(2,ifac) / surfbn(ifac)
rcodcl(ifac,iw,1) = - fmprsc * surfbo(3,ifac) / surfbn(ifac)
if (iturb.eq.32 .or. itytur.eq.5) then
if (mrkcel(iel) .eq. 1) then
rcodcl(ifac,iu,1) = fmprsc/10.d0
endif
endif
uref2 = rcodcl(ifac,iu,1)**2 &
+ rcodcl(ifac,iv,1)**2 &
+ rcodcl(ifac,iw,1)**2
uref2 = max(uref2,1.d-12)
! Turbulence example computed using equations valid for a pipe.
! We will be careful to specify a hydraulic diameter adapted
! to the current inlet.
! We will also be careful if necessary to use a more precise
! formula for the dynamic viscosity use in the calculation of
! the Reynolds number (especially if it is variable, it may be
! useful to take the law from 'usphyv'. Here, we use by default
! the 'viscl0" value.
! Regarding the density, we have access to its value at boundary
! faces (romb) so this value is the one used here (specifically,
! it is consistent with the processing in 'usphyv', in case of
! variable density)
! Hydraulic diameter
xdh = 1.d0
! Calculation of turbulent inlet conditions using
! the turbulence intensity and standard laws for a circular pipe
! (their initialization is not needed here but is good practice)
rhomoy = bfpro_rom(ifac)
call turbulence_bc_inlet_hyd_diam(ifac, uref2, xdh, rhomoy, viscl0, &
rcodcl)
! Handle scalars
if (nscal.gt.0) then
do ii = 1, nscal
rcodcl(ifac,isca(ii),1) = 1.d0
enddo
endif
enddo
if (iturb.eq.32 .or. itytur.eq.5) then
deallocate(mrkcel)
endif
else
! Subsequent time steps
!----------------------
acc(1) = 0.d0
acc(2) = 0.d0
! Estimate multiplier
do ilelt = 1, nlelt
ifac = lstelt(ilelt)
iel = ifabor(ifac)
vnrm = sqrt(cvar_vel(1,iel)**2 + cvar_vel(2,iel)**2 + cvar_vel(3,iel)**2)
acc(1) = acc(1) + vnrm*surfbn(ifac)
acc(2) = acc(2) + surfbn(ifac)
enddo
if (irangp.ge.0) then
call parrsm(2, acc)
endif
fmul = fmprsc/(acc(1)/acc(2)) ! 1 / estimate flow multiplier
! Apply BC
do ilelt = 1, nlelt
ifac = lstelt(ilelt)
iel = ifabor(ifac)
itypfb(ifac) = ientre
vnrm = sqrt(cvar_vel(1,iel)**2 + cvar_vel(2,iel)**2 + cvar_vel(3,iel)**2)
rcodcl(ifac,iu,1) = - fmul * vnrm * surfbo(1,ifac) / surfbn(ifac)
rcodcl(ifac,iv,1) = - fmul * vnrm * surfbo(2,ifac) / surfbn(ifac)
rcodcl(ifac,iw,1) = - fmul * vnrm * surfbo(3,ifac) / surfbn(ifac)
if (itytur.eq.2) then
rcodcl(ifac,ik,1) = cvar_k(iel)
rcodcl(ifac,iep,1) = cvar_ep(iel)
elseif (itytur.eq.3) then
rcodcl(ifac,ir11,1) = cvar_r11(iel)
rcodcl(ifac,ir22,1) = cvar_r22(iel)
rcodcl(ifac,ir33,1) = cvar_r33(iel)
rcodcl(ifac,ir12,1) = cvar_r12(iel)
rcodcl(ifac,ir13,1) = cvar_r13(iel)
rcodcl(ifac,ir23,1) = cvar_r23(iel)
rcodcl(ifac,iep,1) = cvar_ep(iel)
if (iturb.eq.32) then
rcodcl(ifac,ial,1) = cvar_al(iel)
endif
elseif (itytur.eq.5) then
rcodcl(ifac,ik,1) = cvar_k(iel)
rcodcl(ifac,iep,1) = cvar_ep(iel)
rcodcl(ifac,iphi,1) = cvar_phi(iel)
if (iturb.eq.50) then
rcodcl(ifac,ifb,1) = cvar_fb(iel)
elseif (iturb.eq.51) then
rcodcl(ifac,ial,1) = cvar_al(iel)
endif
elseif (iturb.eq.60) then
rcodcl(ifac,ik,1) = cvar_k(iel)
rcodcl(ifac,iomg,1) = cvar_omg(iel)
elseif (iturb.eq.70) then
rcodcl(ifac,inusa,1) = cvar_nusa(iel)
endif
! Handle scalars (a correction similar to that of velocity is suggested
! rather than the simpler code below)
if (nscal.gt.0) then
do ii = 1, nscal
call field_get_val_s(ivarfl(isca(ii)), cvar_scal)
rcodcl(ifac,isca(ii),1) = cvar_scal(ii)
enddo
endif
enddo
endif
!< [example_1]
!--------
! Formats
!--------
!----
! End
!----
deallocate(lstelt) ! temporary array for boundary faces selection
return
end subroutine cs_user_boundary_conditions
|