boundary_velocity Subroutine

public subroutine boundary_velocity(mesh, bc, face_id, interior_velocity, value, boundary_density)

Evaluates the velocity vector at a boundary face.

Implements standard FV boundary evaluations: - Wall/Dirichlet: Returns the specified patch velocity and therefore prescribes volumetric face flux. - mass_flux: Prescribes signed normal mass flux rhodot(U,n), converted to normal velocity using the caller-supplied boundary density. - inlet_mass_flux/outlet_mass_flux: Prescribe positive mass-flux magnitudes and apply the inward/outward sign automatically. Tangential velocity components come from patch_u/v/w. - Symmetry: Subtracts the normal component from the interior velocity. - Neumann/Periodic*: Extrapolates the interior velocity to the face.

Arguments

Type IntentOptional Attributes Name
type(mesh_t), intent(in) :: mesh

Mesh data structure.

type(bc_set_t), intent(in) :: bc

Boundary condition set.

integer, intent(in) :: face_id

ID of the boundary face.

real(kind=rk), intent(in) :: interior_velocity(3)
real(kind=rk), intent(out) :: value(3)
real(kind=rk), intent(in), optional :: boundary_density

Optional density used by mass_flux boundaries [kg/m^3].


Source Code

   subroutine boundary_velocity(mesh, bc, face_id, interior_velocity, value, boundary_density)
      type(mesh_t), intent(in) :: mesh
      type(bc_set_t), intent(in) :: bc
      integer, intent(in) :: face_id
      real(rk), intent(in) :: interior_velocity(3)
      real(rk), intent(out) :: value(3)
      real(rk), intent(in), optional :: boundary_density

      integer :: patch_id
      real(rk) :: nhat(3), un, rho_b, tangent(3)

      nhat = zero
      un = zero
      rho_b = 1.0_rk
      tangent = zero

      patch_id = mesh%faces(face_id)%patch
      if (patch_id <= 0) then
         value = interior_velocity
         return
      end if

      select case (bc%patches(patch_id)%velocity_type_id)
      case (bc_wall, bc_dirichlet)
         value = bc%patches(patch_id)%velocity
      case (bc_mass_flux, bc_inlet_mass_flux, bc_outlet_mass_flux)
         if (.not. present(boundary_density)) then
            call fatal_error('bc', 'mass-flux velocity boundary requires a boundary density')
         end if
         rho_b = max(boundary_density, tiny(1.0_rk))
         nhat = mesh%faces(face_id)%normal
         select case (bc%patches(patch_id)%velocity_type_id)
         case (bc_mass_flux)
            ! Expert/standard FV convention: signed outward-normal mass flux.
            un = bc%patches(patch_id)%mass_flux / rho_b
         case (bc_inlet_mass_flux)
            ! User-friendly inlet convention: positive magnitude into domain.
            un = -abs(bc%patches(patch_id)%mass_flux) / rho_b
         case (bc_outlet_mass_flux)
            ! User-friendly outlet convention: positive magnitude out of domain.
            un = abs(bc%patches(patch_id)%mass_flux) / rho_b
         end select
         tangent = bc%patches(patch_id)%velocity - dot_product(bc%patches(patch_id)%velocity, nhat) * nhat
         value = tangent + un * nhat
      case (bc_symmetry)
         nhat = mesh%faces(face_id)%normal
         un = dot_product(interior_velocity, nhat)
         value = interior_velocity - un * nhat
      case default
         value = interior_velocity
      end select
   end subroutine boundary_velocity