Index
StartUpDG.MeshData
StartUpDG.PIparams
StartUpDG.RefElemData
StartUpDG.RefElemData
StartUpDG.RefElemData
StartUpDG.RefElemData
StartUpDG.RefElemData
StartUpDG.RefElemData
StartUpDG.bcopy!
StartUpDG.ck45
StartUpDG.compute_adaptive_dt
StartUpDG.connect_mesh
StartUpDG.diagE_sbp_nodes
StartUpDG.dp56
StartUpDG.hybridized_SBP_operators
StartUpDG.make_periodic
StartUpDG.readGmsh2D
Functions
StartUpDG.MeshData
— Typestruct MeshData{Dim, Tv}
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 = md
StartUpDG.PIparams
— Typestruct PIparams
Struct containing PI controller parameters.
StartUpDG.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 = rd
StartUpDG.RefElemData
— Methodfunction RefElemData(elem; 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 arg
Constructor 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.bcopy!
— Methodbcopy!(x,y) = x .= y
Convenience routine for operations on tuples of arrays.
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]*res
StartUpDG.compute_adaptive_dt
— Functioncompute_adaptive_dt(Q,rhsQrk,dt,rkE,PI::PIparams,prevErrEst=nothing)
returns acceptstep (true/false), dtnew, errEst. uses PI error control method copied from Paranumal library (Warburton et al).
Inputs:
- Q: container of arrays, Q[i] = ith solution field
- rhsQrk: container whose entries are type(Q) for RK rhs evaluations
StartUpDG.connect_mesh
— Methodconnect_mesh(EToV,fv)
Initialize element connectivity matrices, element to element and element to face connectivity.
Inputs:
EToV
is aK
byNv
matrix whose rows identify theNv
vertices
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.dp56
— Methoddp56()
Dormand-prince 5th order 7 stage Runge-Kutta method (6 function evals via the FSAL property) Returns Butcher table arrays A
,c
and error evolution coefficients rkE
.
Note there is no b
array needed due to the FSAL property and because the last stage is used to compute the error estimator.
StartUpDG.hybridized_SBP_operators
— Methodfunction hybridized_SBP_operators(rd::RefElemData{DIMS})
Constructs hybridized SBP operators given a RefElemData. Returns operators Qrsth...,VhP,Ph.
StartUpDG.make_periodic
— Methodmake_periodic(md::MeshData{Dim},rd::RefElemData,is_periodic...) where {Dim}
make_periodic(md::MeshData{Dim},rd::RefElemData,is_periodic=ntuple(x->true,Dim)) where {Dim}
make_periodic(md::MeshData{1},rd::RefElemData,is_periodic=true)
Returns new MeshData such that the node maps mapP
and face maps FToF
are now periodic. Here, is_periodic
is a tuple of Bool
indicating whether or not to impose periodic BCs in the x
,y
, or z
coordinate.
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")