Index
StartUpDG.MeshDataStartUpDG.RefElemDataStartUpDG.RefElemDataStartUpDG.RefElemDataStartUpDG.RefElemDataStartUpDG.RefElemDataStartUpDG.RefElemDataStartUpDG.VertexMeshPlotterStartUpDG.boundary_face_centroidsStartUpDG.ck45StartUpDG.connect_meshStartUpDG.diagE_sbp_nodesStartUpDG.estimate_hStartUpDG.hybridized_SBP_operatorsStartUpDG.inverse_trace_constantStartUpDG.readGmsh2DStartUpDG.tag_boundary_facesStartUpDG.uniform_mesh
Functions
StartUpDG.MeshData — Typestruct MeshData{Dim, Tv, Ti}MeshData: contains info for a high order piecewise polynomial discretization on an unstructured mesh.
Use @unpack to extract fields. Example:
N,K1D = 3,2
rd = RefElemData(Tri(),N)
VX,VY,EToV = uniform_mesh(Tri(),K1D)
md = MeshElemData(VX,VY,EToV,rd)
@unpack x,y = mdStartUpDG.RefElemData — Typestruct RefElemData{Dim, ElemShape <: AbstractElemShape, Nfaces, Tv}RefElemData: contains info (interpolation points, volume/face quadrature, operators) for a high order nodal basis on a given reference element.
Use @unpack to extract fields. Example:
N = 3
rd = RefElemData(Tri(),N)
@unpack r,s = rdStartUpDG.RefElemData — Methodfunction RefElemData(elem; N, kwargs...)
function RefElemData(elem, approxType; N, kwargs...)Keyword argument constructor for RefElemData (to "label" N via rd = RefElemData(Line(),N=3))
StartUpDG.RefElemData — Methodfunction RefElemData(elementType::Quad, approxType::SBP, N)
SBP reference element data for DGSEM.
StartUpDG.RefElemData — MethodRefElemData(elem::Line, N;
quad_rule_vol = quad_nodes(elem,N+1))
RefElemData(elem::Union{Tri,Quad}, N;
quad_rule_vol = quad_nodes(elem,N),
quad_rule_face = gauss_quad(0,0,N))
RefElemData(elem::Hex,N;
quad_rule_vol = quad_nodes(elem,N),
quad_rule_face = quad_nodes(Quad(),N))
RefElemData(elem; N, kwargs...) # version with keyword argConstructor for RefElemData for different element types.
StartUpDG.RefElemData — Methodfunction RefElemData(elementType::Quad, approxType::SBP, N)
SBP reference element data for DGSEM.
StartUpDG.RefElemData — Methodfunction RefElemData(elementType::Tri, approxType::SBP, N; kwargs...)kwargs = quadraturestrength=2N-1 (or 2N), quadrule_face=:Lobatto (or :Legendre)
StartUpDG.VertexMeshPlotter — TypeVertexMeshPlotter(VX,VY,EToV,fv)
VertexMeshPlotter(triout::TriangulateIO)Plot recipe to plot a quadrilateral or triangular mesh. Usage: plot(MeshPlotter(...))
StartUpDG.boundary_face_centroids — Methodfunction boundary_face_centroids(md)Returns face centroids and boundary_face_ids on the boundaries of the domain given by md::MeshData.
StartUpDG.ck45 — Methodck45()Returns coefficients rka,rkb,rkc for the 4th order 5-stage low storage Carpenter/Kennedy Runge Kutta method. Coefficients evolve the residual, solution, and local time, e.g.,
Example
res = rk4a[i]*res + dt*rhs # i = RK stage
@. u += rk4b[i]*resStartUpDG.connect_mesh — Methodconnect_mesh(EToV,fv)Initialize element connectivity matrices, element to element and element to face connectivity.
Inputs:
EToVis aKbyNvmatrix whose rows identify theNvvertices
which make up an element.
fv(an array of arrays containing unordered indices of face vertices).
Output: FToF, length(fv) by K index array containing face-to-face connectivity.
StartUpDG.diagE_sbp_nodes — MethodTriangular SBP nodes with diagonal boundary matrices. Nodes from Hicken 2019
StartUpDG.estimate_h — Methodestimate_h(rd::RefElemData,md::MeshData)Estimates the mesh size via min sizeofdomain * |J|/|sJ|, since |J| = O(hᵈ) and |sJ| = O(hᵈ⁻¹).
StartUpDG.hybridized_SBP_operators — Methodfunction hybridized_SBP_operators(rd::RefElemData{DIMS})Constructs hybridized SBP operators given a RefElemData. Returns operators Qrsth...,VhP,Ph.
StartUpDG.inverse_trace_constant — Methodfunction inverse_trace_constant(rd::RefElemData)Returns the degree-dependent constant in the inverse trace equality over the reference element (as reported in "GPU-accelerated dG methods on hybrid meshes" by Chan, Wang, Modave, Remacle, Warburton 2016).
Can be used to estimate dependence of maximum stable timestep on degree of approximation.
StartUpDG.readGmsh2D — Methodfunction readGmsh2D(filename)reads triangular GMSH 2D file format 2.2 0 8. returns VX,VY,EToV
Examples
VXY,EToV = readGmsh2D("eulerSquareCylinder2D.msh")StartUpDG.tag_boundary_faces — Methodfunction tag_boundary_faces(md,boundary_name::Symbol=:entire_boundary)
function tag_boundary_faces(md,boundary_list::Dict{Symbol,<:Function})When called without arguments, just returns Dict(:entireboundary=>boundaryfaces).
Example usage:
julia> rd = RefElemData(Tri(),N=1)
julia> md = MeshData(uniform_mesh(Tri(),2)...,rd)
julia> on_bottom_boundary(x,y,tol=1e-13) = abs(y+1) < tol
julia> on_top_boundary(x,y,tol=1e-13) = abs(y-1) < tol
julia> tag_boundary_faces(Dict(:bottom => on_bottom_boundary,
:top => on_top_boundary), md)StartUpDG.uniform_mesh — Method uniform_mesh(elem::Line,Kx)
uniform_mesh(elem::Tri,Kx,Ky)
uniform_mesh(elem::Quad,Kx,Ky)
uniform_mesh(elem::Hex,Kx,Ky,Kz)Uniform Kx (by Ky by Kz) mesh on $[-1,1]^d$, where d is the spatial dimension. Returns VX,VY,VZ,EToV. Can also use kwargs via uniform_mesh(elem; K1D=16)