7.2
general documentation
Atmospheric examples

Local variables to be added

integer ifac, ii
integer izone
integer ilelt, nlelt
integer f_id_rough, f_id_t_rough
double precision d2s3
double precision zref, xuref
double precision ustar, rugd, rugt
double precision zent, xuent, xvent
double precision xkent, xeent
integer, allocatable, dimension(:) :: lstelt
double precision, dimension(:), pointer :: bpro_roughness
double precision, dimension(:), pointer :: bpro_roughness_t

Initialization and finalization

Initialization and finalization is similar to that of the base examples

Example 1

For boundary faces of color 11, assign an inlet boundary condition prescribed from the meteo profile with automatic choice between inlet/ outlet according to the meteo profile.

call getfbr('11',nlelt,lstelt)
! Get a new zone number (1 <= izone <= nozppm)
izone = maxval(izfppp) + 1
do ilelt = 1, nlelt
ifac = lstelt(ilelt)
! - Zone to which the face belongs
izfppp(ifac) = izone
! - Boundary conditions are prescribed from the meteo profile
iprofm(izone) = 1
! - Chemical boundary conditions are prescribed from the chemistry profile
iprofc(izone) = 1
! - boundary condition type can be set to ientre or i_convective_inlet
itypfb(ifac) = ientre
! - automatic determination of type (inlet/outlet) according to sign of
! mass flux
iautom(ifac) = 1
enddo

Example 2

For boundary faces of color 21, assign an inlet boundary condition prescribed from the meteo profile.

call getfbr('21',nlelt,lstelt)
! Get a new zone number (1 <= izone <= nozppm)
izone = maxval(izfppp) + 1
do ilelt = 1, nlelt
ifac = lstelt(ilelt)
! - Zone to which the face belongs
izfppp(ifac) = izone
! - Boundary conditions are prescribed from the meteo profile
iprofm(izone) = 1
! - Chemical boundary conditions are prescribed from the chemistry profile
iprofc(izone) = 1
! - Assign inlet boundary conditions
itypfb(ifac) = ientre
enddo

Example 3

For boundary faces of color 31, assign an inlet boundary condition prescribed from the meteo profile except for dynamical variables which are prescribed with a rough log law.

call getfbr('31',nlelt,lstelt)
! Get a new zone number (1 <= izone <= nozppm)
izone = maxval(izfppp) + 1
do ilelt = 1, nlelt
ifac = lstelt(ilelt)
! - Zone to which the face belongs
izfppp(ifac) = izone
! - Boundary conditions are prescribed from the meteo profile
iprofm(izone) = 1
! - Chemical boundary conditions are prescribed from the chemistry profile
iprofc(izone) = 1
! - Dynamical variables are prescribed with a rough log law
zent=cdgfbo(3,ifac)
ustar = xkappa*xuref/log((zref+rugd)/rugd)
xuent = ustar/xkappa*log((zent+rugd)/rugd)
xvent = 0.d0
xkent = ustar**2/sqrt(cmu)
xeent = ustar**3/xkappa/(zent+rugd)
itypfb(ifac) = ientre
rcodcl(ifac,iu,1) = xuent
rcodcl(ifac,iv,1) = xvent
rcodcl(ifac,iw,1) = 0.d0
! itytur is a flag equal to iturb/10
if (itytur.eq.2) then
rcodcl(ifac,ik,1) = xkent
rcodcl(ifac,iep,1) = xeent
elseif(itytur.eq.3) then
rcodcl(ifac,ir11,1) = d2s3*xkent
rcodcl(ifac,ir22,1) = d2s3*xkent
rcodcl(ifac,ir33,1) = d2s3*xkent
rcodcl(ifac,ir12,1) = 0.d0
rcodcl(ifac,ir13,1) = 0.d0
rcodcl(ifac,ir23,1) = 0.d0
rcodcl(ifac,iep,1) = xeent
elseif(iturb.eq.50) then
rcodcl(ifac,ik,1) = xkent
rcodcl(ifac,iep,1) = xeent
rcodcl(ifac,iphi,1) = d2s3
rcodcl(ifac,ifb,1) = 0.d0
elseif(iturb.eq.60) then
rcodcl(ifac,ik,1) = xkent
rcodcl(ifac,iomg,1) = xeent/cmu/xkent
elseif(iturb.eq.70) then
rcodcl(ifac,inusa,1) = cmu*xkent**2/xeent
endif
enddo

Example 4

Prescribe at boundary faces of color '12' an outlet.

call field_get_id_try("boundary_roughness", f_id_rough)
if (f_id_rough.ge.0) call field_get_val_s(f_id_rough, bpro_roughness)
call field_get_id_try("boundary_thermal_roughness", f_id_t_rough)
if (f_id_t_rough.ge.0) call field_get_val_s(f_id_t_rough, bpro_roughness_t)
call getfbr('15', nlelt, lstelt)
! Get a new zone number (1 <= izone <= nozppm)
izone = maxval(izfppp) + 1
do ilelt = 1, nlelt
ifac = lstelt(ilelt)
! Wall: zero flow (zero flux for pressure)
! rough friction for velocities (+ turbulent variables)
! - Zone to which the zone belongs
izfppp(ifac) = izone
itypfb(ifac) = iparug
! Roughness for velocity: rugd
bpro_roughness(ifac) = rugd
! Roughness for scalars (if required)
if (f_id_rough.ge.0) &
bpro_roughness_t(ifac) = 0.01d0
! By default zero flux for scalars
enddo

Example 5

Prescribe at boundary faces of color 4 a symmetry.

call getfbr('4', nlelt, lstelt)
! Get a new zone number (1 <= izone <= nozppm)
izone = maxval(izfppp) + 1
do ilelt = 1, nlelt
ifac = lstelt(ilelt)
! - Zone to which the zone belongs
izfppp(ifac) = izone
itypfb(ifac) = isymet
enddo