Boundary Integrals

Boundary Integrals

This page describes the functions that impose boundarya conditions. The boundary conditions are imposed weakly, using a penalty between the desired state (the boundary condition value) and the current state (the solution).

Entry Points

This function evaluates the boundary integrals in the Euler equations by calling the appropriate SBP function on eqn.bndryflux, which must be populated before calling this function. eqn.res is updated with the result

This is a mid level function

source

Functions

EulerEquationMod.calcBoundaryFlux

This function calculates the boundary flux for the portion of the boundary with a particular boundary condition. The eqn.q are converted to conservative variables if needed

Inputs: mesh : AbstractMesh sbp : AbstractSBP eqn : EulerEquation functor : a callable object that calculates the boundary flux at a node idx_range: the Range describing which Boundaries have the current BC bndry_facenums: An array with elements of type Boundary that tell which element faces have the boundary condition Outputs: bndryflux : the array to store the boundary flux, corresponds to bndry_facenums

The functor must have the signature functor( q, aux_vars, x, nrm_xy, bndryflux_i, eqn.params) where q are the conservative variables. where all arguments (except params) are vectors of values at a node.

params is the ParamType associated with the the EulerEquation object nrm = mesh.sbpface.normal[:, current_node]

This is a mid level function.

source

Staggered grid version

source

Like calcBoundaryFlux, but performs the integration and updates res rather than storing the flux.

source

EulerEquationMod.getBCFluxes

This function calls other functions to calculate the boundary fluxes, passing them pieces of the array needed. This populates eqn.bndryflux. It also calls writeBoundary() to do any requested output. If the options dictionary specifies not to precompute the boundary flux, this function will do the integration as well and update eqn.res.

This is a mid level function.

source

EulerEquationMod.getBCFunctors

This function uses the opts dictionary to populate mesh.bndry_funcs with the functors

func(params::ParamType,
     q::AbstractArray{Tsol,1},
     aux_vars::AbstractArray{Tres, 1},  coords::AbstractArray{Tmsh,1},
     nrm_xy::AbstractArray{Tmsh,1},
     bndryflux::AbstractArray{Tres, 1},
     bndry::BoundaryNode=NullBoundaryNode)

This is a high level function.

source

EulerEquationMod.interpolateBoundary

Interpolates the solution variables to the exterior boundary of the mesh and calculates any additional quantities at the boundary of the mesh. DG only

Inputs: mesh: an AbstractDGMesh sbp eqn opts q : the 3D array of solution variables for all elements, numdofpernode x numNodesPerElement x numEl

Inputs/Outputs: q_bndry: the array to be populated with the solution interpolated to the boundary, numdofpernode x numNodesPerFace x num boundary faces

eqn.aux_vars_bndry is also populated

Aliasing restrictions: none
source

EulerEquationMod.writeBoundary

This function writes information about the boundary faces and fluxes to files. It is controlled by the input argument writeboundary, of type Bool.

It generates the files: * boundaryfaces.dat : writes mesh.bndryfaces, an array with eltype Boundary to a file, one element per line * boundaryflux.dat : writes the element, local node number and boundary flux to a line in a human readable format * boundaryflux2.dat : writes the real part ofmesh.bndryflux to space delimited file

This is a high level function.

source

Boundary Condition Functors

Each functor defines a different type of boundary condition. Because the equation is solved in the weak form, the functors compute the flux associated with the penalty that imposes the boundary condition.

Maps boundary conditions names to the functor objects. Each functor should be callable with the signature

source

Null boundary node (all fields zero). This is useful as a default value for a functiona argument (if the argument is unused).

source

Type that identified a particular node on a particular face of a particular boundary. It also contains the index of the face within the array of Boundaries with the same boundary condition

Fields

  • element: the element the face is part of

  • face: the local face number

  • faceidx: the index within the array of Boundaries with this BC

  • node: the node on the face

source

EulerEquationMod.FreeStreamBC_dAlpha <: BCTypes

This functor uses the Roe solver to calculate the flux for a boundary state corresponding to the free stream velocity, using rho_free, Ma, aoa, and E_free

This is a low level functor

Arguments

  • obj : Object of type BCType used for multiple dispatch. Every new boundary condition needs to have its own type and entered in BCDict

  • q : Solution variable

  • aux_vars : Auxiliary variables

  • x : physical coordinates of the SBP node

  • nrm_xy : scaled normal vector in x-y space

  • bndryflux : Computed flux value at the boundary

source

###EulerEquationMod.FreeStreamBC_revm

Reverse mode for FreeStreamBC.

Inputs

  • obj : Type of the Boundary condition being evaluated. Its a subtype of BCType_revm

  • q : Solution variable

  • aux_vars : Auxiliary variables

  • x : Node coordinates

  • nrm_xy : scaled normal vector in x-y space

  • bndryflux_bar : Input flux value seed that is used to compute the reverse mode derivative.

  • params : equation object parameters

Output

  • nrm_bar : Derivative of bndryflux_bar w.r.t the mapping jacobian

source

EulerEquationMod.allOnesBC <: BCTypes

This functor uses the Roe solver to calculate the flux for a boundary state where all the conservative variables have a value 1.0

This is a low level functor

source

EulerEquationMod.allZerosBC <: BCTypes

This functor uses the Roe solver to calculate the flux for a boundary state where all the conservative variables have a value 0.0

This is a low level functor

source

EulerEquationMod.isentropicVortexBC <: BCTypes

This type and the associated call method define a functor to calculate the flux using the Roe Solver using the exact InsentropicVortex solution as boundary state. See calcBoundaryFlux for the arguments all functors must support.

This is a low level functor.

Arguments

  • obj : Object of type BCType used for multiple dispatch. Every new boundary condition needs to have its own type and entered in BCDict

  • q : Solution variable

  • aux_vars : Auxiliary variables

  • x : physical coordinates of the SBP node

  • nrm_xy : sclaed face normal in physical space

  • bndryflux : Computed flux value at the boundary

source

EulerEquationMod.isentropicVortexBC_physical <: BCTypes

This type and the associated call method define a functor to calculate the actual Euler flux using the exact InsentropicVortex solution as boundary state. See calcBoundaryFlux for the arguments all functors must support.

This is a low level functor.

source

###EulerEquationMod.isentropicVortexBC_revm

Reverse mode for isentropicVortexBC.

Inputs

  • obj : Type of the Boundary condition being evaluated. Its a subtype of BCType_revm

  • q : Solution variable

  • aux_vars : Auxiliary variables

  • x : Node coordinates

  • nrm : scaled normal vector in x-y space

  • bndryflux_bar : Input flux value seed that is used to compute the reverse mode derivative.

  • params : equation object parameters

Output

  • nrm_bar : Derivative of flux w.r.t the nrm

source

###EulerEquationMod.noPenetrationBC_revm

Reverse mode for noPenetrationBC.

Input

  • obj : Type of the Boundary condition being evaluated. Its a subtype of BCType_revm

  • q : Solution variable

  • aux_vars : Auxiliary variables

  • x : Node coordinates

  • dxidx : Mapping jacobian matrix for the SBP node

  • nrm : sbpface normal vector

  • bndryflux_bar : Input flux value seed that is used to compute the reverse mode derivative.

  • params : equation object parameters

source

This is a special boundary condition used for imposing a numerical solution as a boundary condition.

Fields

* bc_vals: an array of size numDofPerNode x numNodesPerFace x numFaces with
           this boundary condition on it.  This array can be accessed
           using the `faceidx` field of [`BoundaryNode`](@ref).

This BC is special because it has to store the boundary values it is imposing. Whenever this boundary condition is needed, the user should construct a new object, with the data inside it, and then register it just before constructing the EulerData object.

source

EulerEquationMod.unsteadyVortexBC <: BCTypes

This type and the associated call method define a functor to calculate the flux using the Roe Solver using the exact InsentropicVortex solution as boundary state. See calcBoundaryFlux for the arguments all functors must support.

This is a low level functor.

Arguments

  • obj : Object of type BCType used for multiple dispatch. Every new boundary condition needs to have its own type and entered in BCDict

  • q : Solution variable

  • aux_vars : Auxiliary variables

  • x : physical coordinates of the SBP node

  • dxidx : Mapping jacobian matrix for the SBP node

  • nrm_xy : SBP face normal

  • bndryflux : Computed flux value at the boundary

source