A Matlab DG Method For Hyperbolic Function
A Matlab DG Method For Hyperbolic Function
Abstract
This is the first in a series of papers on implementing a discontinuous Galerkin method as a MATLAB / GNU Oct-
ave toolbox. The main goal is the development of techniques that deliver optimized computational performance
combined with a compact, user-friendly interface. Our implementation relies on fully vectorized matrix / vector op-
erations and is carefully documented; in addition, a direct mapping between discretization terms and code routines
is maintained throughout. The present work focuses on a two-dimensional time-dependent diffusion equation with
space / time-varying coefficients. The spatial discretization is based on the local discontinuous Galerkin formulation
and is locally mass conservative. Approximations of orders zero through four based on orthogonal polynomials have
been implemented; more spaces of arbitrary type and order can be easily accommodated by the code structure. Time
discretization is performed using an implicit Euler method.
Keywords: MATLAB, GNU Octave, local discontinuous Galerkin method
1. Introduction
The discontinuous Galerkin (DG) methods first introduced in [1] for a hyperbolic equation started gaining in
popularity with the appearance of techniques to deal with second order terms such as Laplace operators. Two different
approaches to the discretization of second order terms are known in the literature. The first originates from the interior
penalty (IP) methods introduced in the late 1970s and early 1980s for elliptic and parabolic equations (cf. [2] for an
overview). The IP methods discretize the second order operators directly, similarly to the classical finite element
method. To produce a stable scheme, however, they need stabilization terms in the discrete formulation.
In our MATLAB [3] / GNU Octave [4] implementation FESTUNG (F inite E lement S imulation T oolbox for
UN structured G rids) available at [5], we use the second type of DG method called the local discontinuous Galer-
kin (LDG) method first proposed in [6] and further developed in [7, 8]. This scheme relies on a mixed formulation
that replaces each second order equation with two first order equations introducing in the process an auxiliary flux
variable. As opposed to methods from the IP family the LDG method is also consistent for piecewise constant ap-
proximation spaces.
In developing this toolbox we pursue a number of goals:
1. Design a general-purpose software package using the discontinuous Galerkin method—the diffusion operator
described here is only the first step—for a range of standard CFD applications and provide this toolbox as
a research tool in the open source format (cf. [5]).
2. Investigate performance optimized software development strategies in MATLAB / GNU Octave for discontinu-
ous finite element methods on unstructured grids.
∗ Correspondingauthor
Email addresses: [email protected] (Florian Frank), [email protected] (Balthasar Reuter), [email protected]
(Vadym Aizinger), [email protected] (Peter Knabner)
respectively. For a boundary edge E ∈ E∂Ω , only the definition on the left is meaningful. The one-sided values of
a vector-valued quantity v are defined analogously. The average and the jump of w on E are then given by
z = −∇c in J × Ω , (2a)
∂t c + ∇ · (d z) = f in J × Ω , (2b)
c = cD on J × ∂ΩD , (2c)
z · ν = gN on J × ∂ΩN , (2d)
0
c = c on {0} × Ω . (2e)
Placing the test functions left to the solutions has an advantage when assembling the discrete system and will be
discussed at a later time.
denote the broken polynomial space on the triangulation Th . For the semi-discrete formulation, we assume the fol-
lowing approximation for coefficients (for t ∈ J fixed): dh (t), fh (t), c0h ∈ P p (Th ). A specific way to compute these
approximations will be given in Sec. 3.4; here we only state that it is done using the L2 -projection into P p (T ), there-
fore the accuracy improves with increasing polynomial order p.
3
Problem 1 (Semi-discrete problem). Seek (zh (t), ch (t)) ∈ [P p (Th )]2 × P p (Th ) such that the following holds for t ∈ J
and ∀T − ∈ Th , ∀vh ∈ [P p (Th )]2 , ∀wh ∈ P p (Th ) :
Z Z Z
{|ch (t)|} on EΩ
v−h · νT −
vh · zh (t) − ∇ · vh ch (t) + cD (t) on ED = 0,
T − T − ∂T − c− (t) on E
h N
Z Z Z
{|dh (t) zh (t)|} · νT − + h η− ~ch (t) · νT − on EΩ
Z
T
− − − η −
wh ∂t ch (t) − ∇wh · dh (t) zh (t) + wh dh (t) zh (t) · νT − + h − ch (t) − cD (t) on ED = wh fh (t)
T
T− T− ∂T − T−
dh− (t) gN (t)
on EN
where η > 0 is a penalty coefficient.
In the second equation of Problem 1, the presence of the penalty terms is required to ensure a full rank of the system
in the absence of the time derivative [9, Lem. 2.15]. For analysis purposes, the above equations are usually summed
over all triangles T ∈ Th . In the implementation that follows, however, it is sufficient to work with local equations.
Thus far, we used an algebraic indexing style. In the remainder we switch to a mixture of algebraic and numerical
style: for instance, Ekn ∈ ∂T k ∩EΩ means all possible element indices k ∈ {1, . . . , K} and local edge indices n ∈ {1, 2, 3}
such that Ekn lies in ∂T k ∩ EΩ . This implicitly fixes the numerical indices which accordingly can be used to index
matrices or arrays.
We use a bracket notation followed by a subscript to index matrices and multidimensional arrays. Thus, for an
n-dimensional array X, the symbol [X]i1 ,...,in stands for the component of X with index il in the lth dimension. As in
MATLAB / GNU Octave, a colon is used to abbreviate all indices within a single dimension. For example, [X]:,:,i3,...,in
is a two-dimensional array / matrix.
4
We condense the coefficients associated with unknowns into two-dimensional arrays C(t), Z1 (t), Z2 (t), such that Ck j (t) ≔
[C(t)]k, j etc. The vectors [C]k,: and [Zm ]k,: , m ∈ {1, 2}, are called local representation vectors with respect to basis
functions ϕki i∈{1,...,N} for ch and for components of zh , correspondingly. In a similar way, we express the coefficient
functions as linear combinations of the basis functions: On T k , we use the local representation vectors [C0 ]k,: for c0h ,
[D]k,: for dh , and [F]k,: for fh .
Since the basis functions ϕki , i ∈ {1, . . . , N} are supported only on T k , M has a block-diagonal structure
MT 1 Z ϕk1 ϕk1 · · · ϕk1 ϕkN
..
..
.. .. ,
M = .
with MT k ≔
. . . (5)
Tk
MT K ϕkN ϕk1 · · · ϕkN ϕkN
i. e., it consists of K local mass matrices MT k ∈ RN×N . Henceforth we write M = diag MT 1 , . . . , MT K .
m KN×KN
The block matrices H ∈ R , m ∈ {1, 2} from term II are given by
Z
[Hm ](k−1)N+i,(k−1)N+ j ≔ − ∂ xm ϕki ϕk j .
Tk
Hence follows the reason for placing the test function left to the solution: otherwise, we would be assembling the
m
transpose of Hm instead. Similarly to M, matrices Hm are block-diagonal Hm = diag Hm T 1 , . . . , HT K with local matrices
given by
Z ∂ xm ϕk1 ϕk1 · · · ∂ xm ϕk1 ϕkN
Hm
.. .. ..
Tk ≔ − . . . .
Tk
∂ xm ϕkN ϕk1 · · · ∂ xm ϕkN ϕkN
In fact, all block matrices for volume integrals have block-diagonal structure due to the local support of the integrands.
The block matrices Gm ∈ RKN×KN , m ∈ {1, 2} from term V
N
X Z
m
[G ](k−1)N+i,(k−1)N+ j ≔ − Dkl (t) ∂ xm ϕki ϕkl ϕk j
l=1 Tk
are similar except that we have a non-stationary and a stationary factor. The block-diagonal matrices read Gm =
m
diag GmT 1 , . . . , GT K with local matrices
N Z ∂ xm ϕk1 ϕkl ϕk1 · · · ∂ xm ϕk1 ϕkl ϕkN
X
Gm
.. .. ..
Tk ≔ − Dkl (t) . (6)
. . .
l=1 Tk
∂ xm ϕkN ϕkl ϕk1 · · · ∂ xm ϕkN ϕkl ϕkN
Vector L(t) resulting from VII is obtained by multiplication of the representation vector of fh (t) to the global mass
matrix:
h iT
L(t) = M F11 (t) · · · F1N (t) · · · · · · F K1 (t) · · · F KN (t) .
6
Ek− n−
T k−
T k+
νk− n−
Figure 1: Two triangles adjacent to edge Ek− n− . It holds Ek− n− = Ek+ n+ and νk− n− = −νk+ n+ .
N×N
where local matrix SEkn ∈ R corresponds to interior edge Ekn of T k , n ∈ {1, 2, 3}:
Z ϕk1 ϕk1 · · · ϕk1 ϕkN
. .. .. .
SEkn = .. . . (7b)
Ekn
ϕkN ϕk1 · · · ϕkN ϕkN
Entries in off-diagonal blocks in Qm are only non-zero for pairs of triangles T k− , T k+ with ∂T k− ∩ ∂T k+ , ∅. They
consist of the mixed terms containing basis functions from both adjacent triangles and are given as
Z
1
[Qm ](k− −1)N+i,(k+ −1)N+ j ≔ νm − − ϕk− i ϕk+ j . (7c)
2 k n E k − n−
Note that the local edge index n− is given implicitly since ∂T k− ∩ ∂T k+ , ∅ consist of exactly one edge Ek− n− = Ek+ n+ .
Next, consider term VI in (3b) containing average and jump terms that produce contributions to multiple block
matrices for ϕk− i ,
N N Z N N Z
1 m X X 1 X X
νk− n− Dk− l (t) Zkm− j (t) ϕk− i ϕk− l ϕk− j + νm
k− n− Dk+ l (t) Zkm+ j (t) ϕk− i ϕk+ l ϕk+ j
2 l=1 j=1 E k − n− 2 l=1 j=1 E k − n−
and
N Z N Z
η X η X
Ck− j (t) ϕk− i ϕk− j − Ck+ j (t) ϕk− i ϕk+ j .
hT k− j=1 E k − n− hT k− j=1 E k − n−
The first integrals are responsible for entries in the diagonal and off-diagonal blocks of a block matrix Rm ∈ RKN×KN
that end up in the last row of system (15). Entries in diagonal blocks are given component-wise by
N Z
m 1 X m
X 1 X
[R ](k−1)N+i,(k−1)N+ j ≔ ν Dkl (t) ϕki ϕkl ϕk j = νm
kn [REkn ]i, j (8a)
2 E ∈∂T ∩E kn l=1 Ekn 2 Ekn ∈∂T k ∩EΩ
kn k Ω
with
N Z ϕk1 ϕkl ϕk1 · · · ϕk1 ϕkl ϕkN
X
REkn = Dkl (t)
.. .. .. . (8b)
. . .
l=1 Ekn
ϕkN ϕkl ϕk1 · · · ϕkN ϕkl ϕkN
Once again, entries in non-zero off-diagonal blocks consist of the mixed terms:
N Z
1 m X
[Rm ](k− −1)N+i,(k+ −1)N+ j ≔ ν−− Dk+ l (t) ϕk− i ϕk+ l ϕk+ j . (8c)
2 k n l=1 E k − n−
7
All off-diagonal blocks corresponding to pairs of triangles not sharing an edge are zero.
The second integral from term VI results in a block matrix S ∈ RKN×KN similar to Qm . Its entries differ only in the
coefficient and the lack of the normal. hT k from the definition of the penalty term in Problem 1 is replaced here with
the local edge length |Ekn | to ensure the uniqueness of the flux over the edge that is necessary to ensure the local mass
conservation. In the diagonal blocks we can reuse the previously defined SEkn and have
Z
X 1 X 1
[S](k−1)N+i,(k−1)N+ j ≔ η ϕki ϕk j = η [SEkn ]i, j , (9a)
E ∈∂T ∩E
|Ekn | Ekn E ∈∂T ∩E
|Ekn |
kn k Ω kn k Ω
Term VI of (3b) contains dependencies on zh (t), ch (t), and the prescribed data cD (t); thus, it produces three contribu-
tions to system (4): on the left-hand side the blocks Rm D , SD ∈ R
KN×KN
, m ∈ {1, 2} (cf. (8a), (9a))
X
[ Rm
D ](k−1)N+i,(k−1)N+ j ≔ νm
kn [REkn ]i, j , (11)
Ekn ∈∂T k ∩ED
X 1
[SD ](k−1)N+i,(k−1)N+ j ≔ η [SEkn ]i, j , (12)
Ekn ∈∂T k ∩ED
|Ekn |
Neumann Edges EN . Consider the Neumann boundary ∂ΩN . Term III of (3a) replaces the average of the primary
variable over the edge by the interior value resulting in the block-diagonal matrix Qm
N ∈R
KN×KN
(cf. (7a), (7b)) with
X
[QmN ](k−1)N+i,(k−1)N+ j ≔ νm
kn [SEkn ]i, j . (14)
Ekn ∈∂T k ∩EN
Term VI contributes to the right-hand side of system (4) since it contains given data only. The corresponding vec-
tor KN ∈ RKN reads
X X N Z
[KN ](k−1)N+i ≔ − Dkl (t) ϕki ϕkl gN (t) .
Ekn ∈∂T k ∩EN l=1 Ekn
3. Implementation
We obey the following implementation conventions:
1. Compute every piece of information only once. In particular, this means that stationary parts of the linear system
to be solved in a time step should be kept in the memory and not repeatedly assembled and that the evaluation
of functions at quadrature points should be carried out only once.
2. Avoid long for loops. With “long” loops we mean loops that scale with the mesh size, e. g., loops over the
triangles T k ∈ Th or edges Ek ∈ EΩ . Use vectorization instead.
3. Avoid changing the nonzero pattern of sparse matrices. Assemble global block matrices with the command
sparse( , , , , ), kron, or comparable commands.
Furthermore, we try to name variables as close to the theory as possible. Whenever we mention a non-built-in MAT-
LAB / GNU Octave routine they are to be found in Sec. 4.
Table 1: Arrays generated by the routine generateGridData storing topological (top, middle) and geometric (bottom) grid descriptions.
5 3 6 3 3 3
5 5 13 1 11 5 2 4 12 8 22 12 32 16
62 9 3 1 2 3 1 2 3 1 2
6 2 63 7 2
1 1 35 1 5 11 18220 3 7 3 2 3 3 3
1 9 2 53 1
21 3 237 1 6 12 15 18
6 3 2 5 391 36
101 2 11
3 2161 215
1 3 3 3317
232 2 1 32 3 14 4 32 2
110
313
2
2
120
323
2
2
130
3332
2 15 41 238 1 8
10 3 72 136 2 3 1 23212 3 1 3 6 9
3 33 1
4 2 2
6 83 8 3 2 1 3 103 7 34 4 1 1 1
3 91 17
2 3 3 1 32 1 2 11 421 24 1 3 2 1 3 2 1 3 2
2 5 22 6 20
2 2 3 10 322 340
2 12 1 1 34 6 3
3
9
1
7
2 3
19
1
11
2 3
29
1 2
15
1 1 2 3 1132 114 3 2 33 19 3 3 3
12
3 39
3 4 1
331
2 135 334 1 11 3 2 18 2 1
2
1 12 1 3 1
2 33 2 4 1 2 1 9 1 141 25 2 312 323 2
11
311 2
14
321 2
17
3312
82 11 222 21 3 2 1 1
1 1 27 114 4 22 7 17 127
1 213 2 3 15 16 245 1 2 1 2
1 21 1 220 2 3
33152
3
3 130
28 3 29 33 243 3 2 5 8
2 33 1 1 8 2 222
2 212 1 1 1
3 171
8 10
218
3 13 1
46 87 1 3 2 1 3 2 1 3 2
3 113 14 53
2 13 2 1
19 2 2 1
3 193 31344 1 2 6 6 16 10 26 14
3 3 1 2 3 1 2 3 1 2
21 3262 2 24
1 5 3 1 3 2 3 3 3
6 73 2 2 3 29 2 3 27 3
12 1 3 21 225 1 2111
1
18
2
19 1 512 3243 154 2216 3 493 150 1
10 13 16
2
2
1 281 2 30 1 23 22 2 203
326 4 12 5
1
2 38
2
2
15
1
318
2
2
125
3282
3 3
2 1 812 52 17 1 316 23 223
3 48 2
97 1 4 7
1 2 3
1
1 3
1
2
2 2 13
13 2 2 321 2110
47 1
1 3
1
2 1 3
1
2 1 3 2
1
1 11 111 7 3 1 4 5 14 9 24 13
8 5 7 1 1 1
Figure 2: Three examples for commands to build a grid. Each is visualized using visualizeGrid(g) (cf. Sec. 3.10).
Fk : T̂ ∋ x̂ 7→ x ∈ T k .
10
x̂2 vk3
b
Thus any function w : T k → R implies ŵ : T̂ → R by ŵ = w ◦ Fk , i. e., w(x) = ŵ( x̂) . The transformation of the
gradient is obtained by the chain rule:
" #
ˆ ∂ x1 w(x) ∂ x̂1 Fk1 ( x̂) + ∂ x2 w(x) ∂ x̂1 Fk2 ( x̂) ˆ k ( x̂)T ∇w(x) = ∇Fk T ∇w(x) ,
∇ŵ( x̂) = ∇w ◦ Fk ( x̂) = = ∇F
∂ x1 w(x) ∂ x̂2 Fk1 ( x̂) + ∂ x2 w(x) ∂ x̂2 Fk2 ( x̂)
on T k . Since T̂ was explicitly defined, the affine mapping can be expressed explicitly in terms of the vertices vk1 , vk2 , vk3
of T k by
h i
Fk : T̂ ∋ x̂ 7→ Bk x̂ + vk1 ∈ T k with R2×2 ∋ Bk ≔ vk2 − vk1 vk3 − vk1 . (17a)
F−1
k : T k ∋ x 7→ B−1
k (x − vk1 ) ∈ T̂ . (17b)
Since the orientations of all physical triangles are the same as of the reference triangle (cf. Fig. 3), 0 < det Bk = 2|T k |
holds. For a function w : Ω → R, we use transformation formula
Z Z
w(x) dx = 2|T k | w ◦ Fk ( x̂) d x̂ . (18a)
Tk T̂
The transformation rule for an integral over the edge Ekn ⊂ T k reads
Z Z
|Ekn |
w(x) dx = w ◦ Fk ( x̂) d x̂ . (18b)
Ekn |Ên | Ên
The rule (18b) is derived as follows: Denote by γkn : [0, 1] ∋ s 7→ γkn (s) ∈ Ekn a parametrization of the edge Ekn with
derivative γ′kn . For instance γk2 (s) ≔ (1 − s) vk3 + s vk1 . Let γ̂n : [0, 1] ∋ s 7→ γ̂n (s) ∈ Ên be defined analogously. From
Z Z 1 Z 1
w(x) dx = w ◦ γkn (s) |γ′kn (s)| ds = w ◦ γkn (s) |Ekn | ds
Ekn 0 0
and
Z Z 1 Z 1
w ◦ Fk ( x̂) d x̂ = w ◦ Fk ◦ γ̂n (s) |γ̂′n (s)| ds = w ◦ γkn (s) |Ên | ds
Ên 0 0
11
3.3. Numerical integration
As an alternative to the symbolic integration functions provided by MATLAB we implemented a quadrature inte-
gration functionality for triangle and edge integrals. In addition, this functionality is required to produce L2 -projections
(cf. Sec. 3.4) of all nonlinear functions (initial conditions, right-hand side, etc.) used in the system.
Since we transform all integrals on T k ∈ Th to the reference triangle T̂ (cf. Sec. 3.2), it is sufficient to define the
quadrature rules on T̂ (which, of course, can be rewritten to apply for every physical triangle T = FT (T̂ )):
Z XR
ĝ( q̂) ≈ ωr ĝ( q̂r ) (19)
T̂ r=1
with R quadrature points q̂r ∈ T̂ and quadrature weights ωr ∈ R. The order of a quadrature rule is the largest integer s
such that (19) is exact for polynomials g ∈ P s (T̂ ). Note that we exclusively rely on quadrature rules with positive
weights and quadrature points located strictly in the interior of T̂ and not on ∂T̂ . The rules used in the implementation
are found in the routine quadRule2D. The positions of q̂r for some quadrature formulas are illustrated in Fig. 4.
An overview of quadrature rules on triangles is found in the “Encyclopaedia of Cubature Formulas” [20]. For edge
integration, we rely on standard Gauss quadrature rules of required order.
The integrals in (3) contain integrands that are polynomials of maximum order 3p−1 on triangles and of maximum
order 3p on edges. Using quadrature integration one could choose rules that integrate all such terms exactly; however,
sufficient accuracy can be achieved with quadrature rules that are exact for polynomials of order 2p on triangles and
2p + 1 on edges (see [21]).
b b b b
bc bc
bc bc
bc bc bc bc
bc
bc
bc bc bc
bc bc bc bc bc
bc bc
b b b b b b b b
Figure 4: Positions of the quadrature points q̂r on the reference triangle T̂ as used in the routine quadRule2D for quadrature rules of order 2 to 5.
such that dh (t) is an adequate approximation of d(t). A simple way (also used in this work) to produce dh is the
L2 -projection defined locally for T k ∈ Th by
Z Z
∀wh ∈ Pd (T ) , wh dh (t) = wh d(t) .
Tk Tk
where the factor of 2|T k | canceled out. Written in matrix form, this is equivalent to
Dk1 Z ϕ̂1 ( x̂) d t, Fk ( x̂)
. ..
M̂ .. = . d x̂
T̂
DkN ϕ̂N ( x̂) d t, Fk ( x̂)
12
with local mass matrix on the reference triangle M̂ ∈ RN×N defined as in (21). This N × N system of equations
can be solved locally for every k ∈ {1, . . . , K}. Approximating the right-hand side by numerical quadrature (19) and
transposing the equation yields
XR d t, F1 ( q̂r ) d t, F1 ( q̂1 ) . . . d t, F1 ( q̂R ) ω1 ϕ̂1 ( q̂1 ) . . . ω1 ϕ̂N ( q̂1 )
.. h i .. .. .. ..
..
..
D(t) M̂ = ωr
. ϕ̂1 ( q̂r ), . . . , ϕ̂N ( q̂r ) = . . . . . . .
r=1
d t, FK ( q̂r ) d t, FK ( q̂1 ) . . . d t, FK ( q̂R ) ωR ϕ̂1 ( q̂R ) . . . ωR ϕ̂N ( q̂R )
This is the global matrix-valued (transposed) system of equations with unknown D(t) ∈ RK×N and a right-hand side
of dimension K × N. The corresponding routine is projectAlg2DG.
X Z 2 X N
Z X 2
kch − ck2L2 (Ω) = ch (x) − c(x) dx = 2 |T k | Ckl ϕ̂l ( x̂) − c ◦ Fk ( x̂) d x̂
T k ∈Th Tk T k ∈Th T̂ l=1
2
X XR N |T 1 | C11 . . . C1N ϕ̂1 ( q̂1 ) . . . ϕ̂1 ( q̂R ) ω1
X 2 . . . . .
.. − c X , X ... ,
1 2
≈2 |T k | ωr Ckl ϕ̂l ( q̂r ) − c ◦ Fk ( q̂r ) = 2 .. · .. .. ..
T k ∈Th r=1 l=1
|T K | C K1 . . . C KN ϕ̂N ( q̂1 ) . . . ϕ̂N ( q̂R ) ωR
can be assembled using a Kronecker product. The abuse of notation c(X1 , X2 ) clearly stands for the K × R matrix
with the entry c([X1 ]k,r , [X2 ]k,r ) in the kth row and rth column. The above procedure is implemented in the rou-
tine computeL2Error .
3.6. Assembly
The aim of this section is to transform the terms required to build the block matrices in (4) to the reference
triangle T̂ and then to compute those either via numerical quadrature or analytically. The assembly of the block
matrices from the local contributions is then performed in vectorized operations.
For the implementation, we need the explicit form for the components of the mappings Fk : T̂ → T k and their
inverses F−1 k : T k → T̂ as defined in (17). Recalling that 0 < det Bk = 2|T k | (cf. Sec. 3.2), we obtain
" 11 1 # " #
Bk x̂ + B12 2 1
k x̂ + v1k −1 1 B22 1 1 12 2 2
k (x − v1k ) − Bk (x − v1k ) .
Fk ( x̂) = and F (x) =
B21 1 22 2 2
k x̂ + Bk x̂ + v1k
k
2 |T k | B11 2 2 21 1 1
k (x − v1k ) − Bk (x − v1k )
In the following, we present the necessary transformation for all blocks of system (4) and name the corresponding
MATLAB / GNU Octave routines that can be found in Sec. 4.
13
3.6.1. Assembly of M
Using the transformation rule (18a) we see that for the local mass matrix MT k as defined in (5) holds
Z ϕ̂1 ϕ̂1 · · · ϕ̂1 ϕ̂N
.
MT k = 2|T k | M̂ with M̂ ≔ .. .. .. , (21)
. .
T̂
ϕ̂N ϕ̂1 · · · ϕ̂N ϕ̂N
where M̂ ∈ RN×N is the representation of the local mass matrix on the reference triangle T̂ . With (5) we see that the
global mass matrix M can be expressed as a Kronecker product of a matrix containing the areas |T k | and the local
matrix M̂:
MT 1 |T 1 |
M =
..
= 2
..
⊗ M̂ .
.
.
MT K |T K |
In the corresponding assembly routine assembleGlobM , the sparse block-diagonal matrix is generated using the com-
mand spdiags with the list g.areaT (cf. Tab. 1).
3.6.2. Assembly of Hm
The transformation rules (18a) and (20) yield
H1T k = −B22 21
k [Ĥ]:,:,1 + Bk [Ĥ]:,:,2 and H2T k = B12 11
k [Ĥ]:,:,1 − Bk [Ĥ]:,:,2
with
Z ∂ x̂m ϕ̂1 ϕ̂1 · · · ∂ x̂m ϕ̂1 ϕ̂N
[Ĥ]:,:,m ≔
.. .. .. ∈ RN×N (22)
. . .
T̂
∂ x̂m ϕ̂N ϕ̂1 · · · ∂ x̂m ϕ̂N ϕ̂N
for m ∈ {1, 2}. Similar to M, the global matrices Hm are assembled by Kronecker products in the routine assembleGlobH .
3.6.3. Assembly of Gm
Application of the product rule, (18a), and (20) give us
Z Z
∂ x1 ϕki ϕk j = B22
k [ Ĝ] i, j,l,1 − B 21
k [ Ĝ] i, j,l,2 , ∂ x2 ϕki ϕk j = −B12 11
k [Ĝ]i, j,l,1 + Bk [Ĝ]i, j,l,2 ,
Tk Tk
with a multidimensional array Ĝ ∈ RN×N×N×2 representing the transformed integral on the reference triangle T̂ :
Z
[Ĝ]i, j,l,m ≔ ∂ x̂m ϕ̂i ϕ̂ j ϕ̂l (23)
T̂
for m ∈ {1, 2}. Now we can express the local matrix G1T k from (6) as
N Z ∂ x̂1 ϕ̂1 ϕ̂1 ϕ̂l · · · ∂ x̂1 ϕ̂1 ϕ̂N ϕ̂l Z ∂ x̂2 ϕ̂1 ϕ̂1 ϕ̂l · · · ∂ x̂2 ϕ̂1 ϕ̂N ϕ̂l
X
G1T k = −
22
Dkl (t) Bk
.. .. .. − B21 .. .. ..
. . . k
. . .
l=1 T̂
T̂
∂ x̂1 ϕ̂N ϕ̂1 ϕ̂l · · · ∂ x̂1 ϕ̂N ϕ̂N ϕ̂l ∂ x̂2 ϕ̂N ϕ̂1 ϕ̂l · · · ∂ x̂2 ϕ̂N ϕ̂N ϕ̂l
N
X
=− Dkl (t) B22 21
k [Ĝ]:,:,l,1 − Bk [Ĝ]:,:,l,2
l=1
Z Z Z 1 Z 1
|Ekn | |Ekn |
ϕki ϕk j = ϕ̂i ( x̂) ϕ̂ j ( x̂) d x̂ = ϕ̂i ◦ γ̂n (s) ϕ̂ j ◦ γ̂n (s) |γ̂′n (s)| ds = |Ekn | ϕ̂i ◦ γ̂n (s) ϕ̂ j ◦ γ̂n (s) ds ,
Ekn |Ên | Ên |Ên | 0
| 0
{z }
≕[Ŝdiag ]i, j,n
(24)
where we used transformation rule (18b) and |γ̂′n (s)| = |Ên |. The explicit forms of the mappings γ̂n : [0, 1] → Ên can
be easily derived:
" # " # " #
1−s 0 s
γ̂1 (s) ≔ , γ̂2 (s) ≔ , γ̂3 (s) ≔ . (25)
s 1−s 0
Thus, we have SEkn = |Ekn |[Ŝdiag ]:,:,n allowing to define the diagonal blocks of the global matrix Sdiag using the Kron-
ecker product:
X δE1n ∈EΩ
3
Sdiag ≔ η
..
⊗ [Ŝdiag ]:,:,n ,
.
n=1
δEKn ∈EΩ
Z Z =I Z
z }| {
|Ek− n− | |Ek− n− |
ϕk− i ϕk+ j = ϕk− i ◦ Fk− ( x̂) ϕk+ j ◦ Fk+ ◦ F−1
k + ◦Fk − ( x̂) d x̂ = ϕ̂i ( x̂) ϕ̂ j ◦ F−1
k+ ◦ Fk− ( x̂) d x̂
E k − n− |Ên− | Ên− |Ên− | Ên−
Z 1
= |Ek− n− | ϕ̂i ◦ γ̂n− (s) ϕ̂ j ◦ F−1
k+ ◦ Fk− ◦ γ̂n− (s) ds .
0
Here one has to note that all maps ϑ̂n− n+ reverse the edge orientation because the affine mappings Fk maintain
the counter-clockwise vertex orientation for each T k ∈ Th . We define Ŝoffdiag ∈ RN×N×3×3 by
Z 1
[Ŝoffdiag ]i, j,n− ,n+ ≔ ϕ̂i ◦ γ̂n− (s) ϕ̂ j ◦ ϑ̂n− n+ ◦ γ̂n− (s) ds , (27)
0
15
and thus we arrive at
0 δE1n− =E2n+ ... ... δE1n− =EKn+
.. ..
3 X 3
δE − =E + 0 . .
X 2n . 1n ..
S offdiag
≔ −η .. .. .. .. offdiag
. . . . ⊗ [Ŝ ]:,:,n− ,n+ .
n =1 n =1
− + . ..
.. . 0 δE(K−1)n− =EKn+
δEKn− =EKn+ ... . . . δEKn− =E(K−1)n+ 0
The sparsity structure for off-diagonal blocks depends on the numbering of mesh entities and is given for each combi-
nation of n− and n+ by the list markE0TE0T (cf. Tab. 1). The routine assembleGlobS assembles the matrices Sdiag and
Soffdiag directly into S with a code very similar to the formulation above.
3.6.5. Assembly of Qm
The assembly of Qm from equations (7a), (7c) is analogous to S since both are constructed from the same terms
only differing in constant coefficients. Consequently, we can choose the same approach as described in 3.6.4. Again,
we split the matrix into diagonal and off-diagonal blocks Qm = Qm,diag + Qm,offdiag and assemble each separately
exploiting transformation rule (18b). This allows to write the diagonal blocks as follows:
m
3 δE1n ∈EΩ ν1n |E1n |
1 X
. .
Q m,diag
≔
..
◦
.. ⊗ [Ŝdiag ]:,:,n ,
2 n=1
δEKn ∈EΩ νmKn |E Kn |
where “◦” is the operator for the Hadamard product.
The off-diagonal blocks are assembled as before, using the mapping ϑ̂n− n+ from (26), leading to a similar repre-
sentation as for Soffdiag ,
0 δE1n− =E2n+ ... ... δE1n− =EKn+
..
..
3 3
δE − =E +
2n 1n 0 . .
m,offdiag 1 X X .. .. .. .. ..
Q ≔ . . . . .
2 n− =1 n+ =1
.. ..
. . 0 δE(K−1)n− =EKn+
δEKn− =EKn+ ... . . . δEKn− =E(K−1)n+ 0
m
ν1n− |E1n− | . . . νm 1n− |E 1n |
−
◦
.. .. offdiag
. . ⊗ [Ŝ ]:,:,n− ,n+ .
m m
νKn− |E Kn− | . . . νKn− |E Kn− |
Once again, we can use a code close to the mathematical formulation to assemble the matrices Qm,diag and Qm,offdiag .
This is realized in the routine assembleGlobQ . In the implementation the Hadamard product is replaced by a call to
the built-in function bsxfun which applies a certain element-by-element operation (here: @times) to arrays. Since all
columns in the second matrix of the product are the same this makes superfluous explicitly creating this matrix and
permits the use of a single list of all required values instead.
3.6.6. Assembly of Rm
Just as before we split the block matrices Rm from (8a), (8c) into diagonal and off-diagonal parts as Rm = Rm,diag +
m,offdiag
R . Here, integrals consist of three basis functions due to the diffusion coefficient but still can be transformed in
the same way. In diagonal blocks, this takes the form
Z Z 1
ϕki ϕkl ϕk j = |Ekn | ϕ̂i ◦ γ̂n (s) ϕ̂l ◦ γ̂n (s) ϕ̂ j ◦ γ̂n (s) ds , (28)
Ekn
|0 {z }
≕[R̂diag ]i, j,l,n
16
which can be used to define a common multidimensional array R̂diag ∈ RN×N×N×3 . This allows to rewrite the local
block matrix from (8b) as
N
X
REkn = Dkl (t)|Ekn |[R̂diag ]:,:,l,n .
l=1
The off-diagonal entries consist of integrals over triples of basis functions two of which belong to the adjacent
triangle T k+ ; thus making it necessary to apply the mapping ϑ̂n− n+ from (26). Once again, this can be written as
Z Z 1
ϕk− i ϕk+ l ϕk+ j = |Ek− n− | ϕ̂i ◦ γ̂n− (s) ϕ̂l ◦ ϑ̂n− n+ ◦ γ̂n− (s) ϕ̂ j ◦ ϑ̂n− n+ ◦ γ̂n− (s) ds , (29)
E k − n−
|0 {z }
≕ [R̂offdiag ]i, j,l,n− ,n+
with a multidimensional array R̂offdiag ∈ RN×N×N×3×3 whose help allows us to denote the assembly of Rm,offdiag
(component-wise given in (8c)) by
0 δE1n− =E2n+ · · · ··· δE1n− =EKn+
..
..
3 X 3 X N
δE − =E +
2n 1n
0 . .
1 X
.. .. .. .. ..
Rm,offdiag ≔ . . . . .
2 n− =1 n+ =1 l=1
.. ..
. . 0 δE(K−1)n− =EKn+
δEKn− =EKn+ ··· · · · δEKn− =E(K−1)n+ 0
m
ν1n− |E1n− | · · · νm 1n− |E 1n |
−
D1l (t) · · · DKl (t)
◦
.. .. .
◦ .. .. ⊗ [R̂offdiag ] − + .
. . . :,:,l,n ,n
m m
νKn− |E Kn | · · · νKn− |E Kn |
− − D1l (t) · · · DKl (t)
The corresponding code is found in assembleGlobR . Again, we make use of the function bsxfun to carry out the
Hadamard product; it is used twice, first to apply the row vector as before and then to apply the column vector of the
diffusion coefficient.
3.6.7. Assembly of Rm D
Since the entries of Rm D in (11) are computed in precisely the same way as for R
m,diag
(cf. (8a)), the corresponding
assembly routine assembleGlobRD consists only of the first part of assembleGlobR which is responsible for the assem-
bly of Rm,diag . It only differs in a factor and the list of edges, namely all edges in the set ED for which non-zero entries
(here given by markE0Tbdr) are generated.
3.6.8. Assembly of SD
diag
The same as for RmD holds for SD in (12) which, in fact, has the same entries in the diagonal blocks as S
(cf. (9a)). Consequently, the corresponding routine assembleGlobSD is again a subset of assembleGlobS .
17
3.6.10. Assembly of J Dm
The entries of J Dm in (10) are transformed using transformation rule (18b)
Z Z
m |E kn |
X X
m m
[J D ](k−1)N+i = − νkn ϕki cD (t) = − νkn ϕki ◦ Fk ( x̂) cD t, Fk ( x̂) d x̂
Ekn ∈∂T k ∩ED Ekn Ekn ∈∂T k ∩ED |Ên | Ên
X Z 1
= − νm
kn |E kn | ϕ̂i ◦ γ̂ n (s) c D t, F k ◦ γ̂ n (s) ds .
Ekn ∈∂T k ∩ED 0
This integral is then approximated using a 1D quadrature rule (19) on the interval (0, 1)
X R
X
[J Dm ](k−1)N+i ≈ − νm
kn |E kn | ωr ϕ̂i ◦ γ̂n (qr ) cD (t, Fk ◦ γ̂n (qr )) ,
Ekn ∈∂T k ∩ED r=1
allowing to vectorize the computation over all triangles and resulting in the routine assembleGlobJD .
3.6.11. Assembly of KD
The computation of KD is done similarly to the assembly of J Dm . The component-wise integrals from (13) are,
once again, transformed to the interval [0, 1] using (18b):
Z X Z 1
X 1
[KD ](k−1)N+i = η ϕki cD (t) = η ϕ̂i ◦ γ̂n (s) cD (t, Fk ◦ γ̂n (s)) ds ,
E ∈∂T ∩E
|Ekn | Ekn E ∈∂T ∩E 0
kn k D kn k D
effectively canceling out the edge length. Using a quadrature rule and vectorization, KD is assembled in the rou-
tine assembleGlobKD .
3.6.12. Assembly of KN
In the integral terms of KN an additional basis function from the diffusion coefficient appears. As before, the
integrals are transformed using transformation rules (18b), (25) and a 1D quadrature rule (19):
X N
X Z X N
X Z 1
[KN ](k−1)N+i = − Dkl (t) ϕki ϕkl gN (t) = − |Ekn | Dkl (t) ϕ̂i ◦ γ̂n (s) ϕ̂l ◦ γ̂n (s) gN (t, Fk ◦ γ̂n (s)) ds
Ekn ∈∂T k ∩EN l=1 Ekn Ekn ∈∂T k ∩EN l=1 0
X N
X R
X
≈ − |Ekn | Dkl (t) ωr ϕ̂i ◦ γ̂n (qr ) ϕ̂l ◦ γ̂n (qr ) gN (t, Fk ◦ γ̂n (qr )) .
Ekn ∈∂T k ∩EN l=1 r=1
Once again, using vectorization over all triangles the assembly routine is assembleGlobKN .
Grid data structures. The size of the grid data structures depends on the number of mesh entities. To give a rough
estimate, we will assume certain simplifications: Each triangle has three incident edges and each edge (disregarding
boundary edges) has two incident triangles, hence it holds 3 #T ≈ 2 #E. Additionally, the number of vertices is usually
less than the number of triangles, i. e., #V . #T . Using those assumptions, the memory requirement of the grid data
structures (cf. Tab. 1 plus additional lists not shown there) amounts to ≈ 89 #T · 8 Bytes.
Degrees of freedom. The memory requirements for system (15) largely depend on the sparsity structure of the matrix
which varies with the numbering of the mesh entities, the number of boundary edges, etc. In MATLAB / GNU Octave,
the memory requirements for a sparse matrix A ∈ Rn×n on a 64-bit machine can be approximated by 16·nnz(A)+8·n+8
Bytes, with nnz(A) being the number of non-zero entries in A. Coefficients, like Ck j , Dkl , Fkl , k ∈ {1, . . . , K}, j, l ∈
{1, . . . , N}, and the right-hand side entries J Dm , KD , KN , L are stored in full vectors, each of which requires KN ·8 Bytes.
The blocks of the system matrices A(t) and W can be divided into two groups: (i) blocks built from element-wise
integrations and (ii) blocks built from edge contributions. For the first case we showed in Sec. 2.4.3 that these have
a block diagonal structure due to the local support of the basis functions. Consequently, each contains K blocks of
size N × N with nonzero entries. Blocks from edge integrals also have block diagonal entries but additionally for each
of the three edges of an element two nonzero blocks exist. We neglect the blocks for boundary edges here (i. e., QN ,
RD , SD ), since these hold only entries for edges that are not contained in the interior edge blocks.
Before solving for the next time step, these blocks are assembled into the system matrices and the right-hand side
vector, effectively doubling the memory requirement. Combining these estimates, this sums to a memory requirement
of ≈ (27KN + 82KN 2 ) · 8 Bytes. Compared to (18(KN)2 + 9KN) · 8 Bytes alone for the assembled system, when using
full matrices, this is still a reasonable number and emphasizes once more the need for sparse data structures.
Total memory usage. Note that all these values are highly dependent on the shape of the mesh and additional overhead
introduced by MATLAB / GNU Octave (e. g., for GUI, interpreter, cell-data structures, temporary storage of built-in
routines, etc.). Other blocks, e. g., blocks on the reference element, like M̂, counters, helper variables, lookup tables
for the basis functions, etc., don’t scale with the mesh size and are left out in these estimates. Hence, the numbers
given here should be understood as a lower bound. Combining the partial results for the memory usage, we obtain
the total amount of ≈ (90 + 43N + 68N 2 )K · 8 Bytes. This means, computing with quadratic basis functions on a grid
of 10 000 triangles requires at least 214 MBytes of memory and should therefore be possible on any current machine.
However, computing with a grid of half a million elements and polynomials of order 4 requires more than 60 GBytes
of memory, requiring a well-equipped workstation.
19
other part takes up less than 5 % of the total computation time. Although our code is not parallelized, MATLAB’s
built-in routines (in particular, mldivide and bsxfun) make extensive use of multithreading. Hence, these results are
machine dependent and, especially on machines with more cores, the runtime distribution might be different.
For large enough a number of time steps one can disregard the execution time of the initial computations (gen-
eration of grid data, computation of basis function lookup tables and reference element blocks, etc.); then the total
runtime scales linearly with the number of time steps. Therefore, we investigate the runtime behavior of the station-
ary test case described in Sec. 3.9 and plot the computation times with respect to the grid size and number of basis
functions in Fig. 5. This shows that the computation time increases with the number of basis functions but primarily
depends on the mesh size. We also observe that doubling the number of elements increases the computation time by
more than a factor of two which is related to the fact that some steps of the algorithm (e. g., the linear solver) have a
complexity that does not depend linearly on the number of degrees of freedom.
Computation time for different grid sizes Computation time for different polynomial orders
6
10
104 p=0 p=1 p=2 K = 36 K = 144 K = 576
p=3 p=4 105 K = 2304 K = 9216 K = 36864
103 K = 147456
104
102 103
time [s]
time [s]
102
101
101
100
100
10−1 10−1
Figure 5: Computation times for varying grid size (left) and approximation order (right). Polynomial orders p = 0, 1, 2, 3, 4 correspond to
N = 1, 3, 6, 10, 15 basis functions, respectively. Measured on 4× Intel Xeon E7-4870 CPUs (each 10 cores, 20 threads) with 500 GBytes RAM and
MATLAB R2014a (8.3.0.532).
20
p 0 1 2 3 4
j kch − ck α kch − ck α kch − ck α kch − ck α kch − ck α
0 2.37E–1 — 7.70E–2 — 1.45E–2 — 1.97E–3 — 4.18E–4 —
1 7.67E–2 1.63 2.46E–2 1.65 1.26E–3 3.52 1.53E–4 3.69 1.26E–5 5.05
2 8.94E–2 – 0.22 6.67E–3 1.88 1.11E–4 3.51 1.02E–5 3.90 3.78E–7 5.06
3 9.39E–2 – 0.07 1.71E–3 1.96 1.10E–5 3.33 6.42E–7 3.99 1.16E–8 5.03
4 9.48E–2 – 0.01 4.32E–4 1.99 1.22E–6 3.18 3.94E–8 4.03 3.59E–10 5.01
5 9.49E–2 0.00 1.08E–4 2.00 1.44E–7 3.09 2.43E–9 4.02 1.11E–11 5.01
6 9.49E–2 0.00 2.71E–5 2.00 1.75E–8 3.04 1.51E–10 4.01 3.64E–13 4.94
1
Table 2: Discretization errors measured in L2 (Ω) and estimated orders of convergences using the penalty η = 1. We have h j = 3·2 j
and K = 36 · 4 j
triangles in the jth refinement level.
3.10. Visualization
In order to get a deeper insight into the data associated with a grid Th or with a discrete variable from P p (Th ), we
provide the routines visualizeGrid and visualizeData , respectively (see Sec. 4 for documentation). Since our code
accepts arbitrary basis functions, in particular the modal basis functions of Sec. 2.4.1, we have to sample those at the
Lagrangian points on each triangle (i. e. the barycenter of T for P0 (T ), the vertices of T for P1 (T ) and the vertices and
edge barycenters for P2 (T )). This mapping from the DG to the Lagrangian basis is realized in projectDG2Lagrange .
The representation of a discrete quantity in the latter basis which is as the DG representation a K × N matrix, is then
used to write a VTK-file [22]. These can be visualized and post-processed, e. g., by Paraview [23]. Unfortunately, the
current version 4.0.1 does not visualize quadratic functions as such but instead uses piecewise linear approximations
consisting of four pieces per element.
As an example, the following code generates a grid with two triangles and visualizes it using visualizeGrid ;
then a quadratic, discontinuous function is projected into the DG space for p ∈ {0, 1, 2} and written to VTK-files
using visualizeData . The resulting outputs are shown in Fig. 6 and Fig. 7.
g = g e n e r a t e G r i d D a t a([0 , -1; sqrt (3) ,0; 0 ,1; - sqrt (3) ,0] , [4 ,1 ,3; 1 ,2 ,3]);
g . idE = ( abs ( g . nuE (: ,2)) >0).*(( g . nuE (: ,1) >0)+(g . nuE (: ,2) >0)*2+1);
v i s u a l i z e G r i d( g )
fAlg = @ ( X1 , X2 ) ( X1 <0).*( X1 .^2 - X2 .^2 -1) + ( X1 >=0).*( - X1 .^2 - X2 .^2+1);
for N = [1 , 3 , 6]
p = ( sqrt (8* N +1) -3)/2; quadOrd = max (2*p , 1); c o m p u t e B a s e s O n Q u a d( N );
fDG = p r o j e c t A l g 2 D G(g , fAlg , quadOrd , c o m p u t e H a t M( N ));
fLagr = p r o j e c t D G 2 L a g r a n g e( fDG );
v i s u a l i z e D a t a(g , fLagr , ' funname ' , [ ' fDOF ' , int2str ( N )] , 1)
end % for
3 3 3 4
5 3
2 1
4 1 1 1 2 2 2 2 2
3 3
4 1
1 2 1 2
√ T √ T
Figure 6: Visualization of a grid with vertices v11 = [− 3, 0] , v12 = v21 = [0, −1]T , v13 = v23 = [0, 1]T , v22 = [ 3, 0] (cf. Sec. 3.10) with
triangle numbers, global and local vertex numbers, global and local edge numbers, and (scaled) edge normals. Global numbers are printed on, local
numbers next to the respective mesh entity. The boundary ID’s which are used to associate parts of the boundary with specific boundary conditions,
are printed on the exterior of each boundary edge.
21
range: [−2/3, 1/3] range: [−8/5, 6/5] range: [−2, 2]
codefDOF1.1.vtu fDOF3.1.vtu fDOF6.1.vtu
Figure 7: Visualization of a discontinuous function, represented with constant, linear, and quadratic basis functions (left to right) using Par-
aview (cf. Sec. 3.10). The underlying mesh is drawn in black.
4. Register of Routines
We list here all routines of our implementation in alphabetic order. For the reason of compactness, we waive
the check for correct arguments, e. g., by means of routines as assert. However, it is strongly recommended to
catch exceptions if the code is to be extended. The argument g is always a struct representing the triangulation Th
(cf. Sec. 3.1), the argument N is always the number of local basis functions N. A script that demonstrates the applica-
tion of the presented routines is given in main.m.
globG = assembleGlobG(g, hatG, dDG) assembles the matrices Gm , m ∈ {1, 2} according to Sec. 3.6.3 which are returned in the 2 × 1
cell variable globG. The input argument hatG stores the local matrices Ĝ (multidimensional array) as defined in (23) and can be computed
by computeHatG . The coefficients of the projection of the algebraic diffusion coefficient d into the broken polynomial space are stored in the
input argument dDG as explained in Sec. 3.4 and can be computed by projectAlg2DG .
function globG = a s s e m b l e G l o b G(g , hatG , dDG )
[K , N ] = size ( dDG );
globG = cell (2 , 1);
globG {1} = sparse ( K *N , K * N ); globG {2} = sparse ( K *N , K * N );
for l = 1 : N
globG {1} = globG {1} + kron ( spdiags ( - dDG (: , l ) .* g . B (: ,2 ,2) , 0 ,K , K ) , hatG (: ,: ,l ,1)) ...
- kron ( spdiags ( - dDG (: , l ) .* g . B (: ,2 ,1) , 0 ,K , K ) , hatG (: ,: ,l ,2));
globG {2} = globG {2} - kron ( spdiags ( - dDG (: , l ) .* g . B (: ,1 ,2) , 0 ,K , K ) , hatG (: ,: ,l ,1)) ...
+ kron ( spdiags ( - dDG (: , l ) .* g . B (: ,1 ,1) , 0 ,K , K ) , hatG (: ,: ,l ,2));
end % for
end % f u n c t i o n
globH = assembleGlobH(g, hatH) assembles the matrices Hm , m ∈ {1, 2} according to Sec. 3.6.2 which are returned in the 2×1 cell vari-
able globH. The input argument hatH stores the local matrices Ĥ as defined in (22) and can be computed by computeHatH.
function globH = a s s e m b l e G l o b H(g , hatH )
K = g . numT ; N = size ( hatH , 1);
globH = cell (2 , 1);
globH {1} = sparse ( K *N , K * N ); globH {2} = sparse ( K *N , K * N );
globH {1} = - kron ( spdiags ( g . B (: ,2 ,2) , 0 ,K , K ) , hatH (: ,: ,1)) ...
+ kron ( spdiags ( g . B (: ,2 ,1) , 0 ,K , K ) , hatH (: ,: ,2));
globH {2} = kron ( spdiags ( g . B (: ,1 ,2) , 0 ,K , K ) , hatH (: ,: ,1)) ...
- kron ( spdiags ( g . B (: ,1 ,1) , 0 ,K , K ) , hatH (: ,: ,2));
end % f u n c t i o n
globJD = assembleGlobJD(g, markE0Tbdr, cDAlg, N) assembles the contributions of Dirichlet boundaries J Dm , m ∈ {1, 2} to the
right-hand side of Equation (3a) according to Sec. 3.6.10 which are returned in the 2×1 cell variable globJD. The input argument markE0Tbdr
is a K × 3 logical array marking for each triangle the local edges on the Dirichlet boundary. Assuming the Dirichlet boundary has the ID 1 then
it can be computed as markE0Tbdr = g.idE0T==1. cDAlg is a function_handle for the algebraic representation of cD (cf. main.m).
Here, Q2X1 and Q2X2 carry out mappings Fk from the reference triangle T̂ to physical triangles T k for all k ∈ {1, . . . , K} at once. This allows to
evaluate cD (t) at all quadrature points in the physical domain at the same time and, thus, vectorize the quadrature with respect to the triangle index.
function globJD = a s s e m b l e G l o b J D(g , markE0Tbdr , cDAlg , N )
22
global gPhi1D
K = g . numT ; p = ( sqrt (8* N +1) -3)/2;
qOrd = 2* p +1; [Q , W ] = q u a d R u l e 1 D( qOrd );
Q2X1 = @ ( X1 , X2 ) g . B (: ,1 ,1)* X1 + g . B (: ,1 ,2)* X2 + g . coordV0T (: ,1 ,1)* ones ( size ( X1 ));
Q2X2 = @ ( X1 , X2 ) g . B (: ,2 ,1)* X1 + g . B (: ,2 ,2)* X2 + g . coordV0T (: ,1 ,2)* ones ( size ( X1 ));
globJD = cell (2 , 1);
globJD {1} = zeros ( K , N ); globJD {2} = zeros (K , N );
for n = 1 : 3
[ Q1 , Q2 ] = gammaMap (n , Q );
cDn = cDAlg ( Q2X1 ( Q1 , Q2 ) , Q2X2 ( Q1 , Q2 ));
Jkn = - m a r k E 0 T b d r(: , n ) .* g . sigE0T (: , n ) .* g . areaE0T (: , n );
for i = 1 : N
integral = cDn * ( W .* gPhi1D { qOrd }(: , i , n ) ' ) ' ;
globJD {1}(: ,i ) = globJD {1}(: ,i ) + Jkn .* g . nuE0T (: ,n ,1) .* integral ;
globJD {2}(: ,i ) = globJD {2}(: ,i ) + Jkn .* g . nuE0T (: ,n ,2) .* integral ;
end % for
end % for
globJD {1} = reshape ( globJD {1} ' ,K *N ,1); globJD {2} = reshape ( globJD {2} ' ,K *N ,1);
end % f u n c t i o n
globKD = assembleGlobKD(g, markE0Tbdr, cDAlg, eta, N) assembles the contributions of Dirichlet boundaries KD to the right-
hand side of Equation (3b) according to Sec. 3.6.11 which are returned in variable globKD. Input arguments markE0Tbdr and cDAlg are as
described in assembleGlobJD . The penalty coefficient η is given in eta.
function globKD = a s s e m b l e G l o b K D(g , markE0Tbdr , cDAlg , eta , N )
global gPhi1D
p = ( sqrt (8* N +1) -3)/2;
qOrd = 2* p +1; [Q , W ] = q u a d R u l e 1 D( qOrd );
Q2X1 = @ ( X1 , X2 ) g . B (: ,1 ,1)* X1 + g . B (: ,1 ,2)* X2 + g . coordV0T (: ,1 ,1)* ones ( size ( X1 ));
Q2X2 = @ ( X1 , X2 ) g . B (: ,2 ,1)* X1 + g . B (: ,2 ,2)* X2 + g . coordV0T (: ,1 ,2)* ones ( size ( X1 ));
globKD = zeros ( g . numT , N );
for n = 1 : 3
[ Q1 , Q2 ] = gammaMap (n , Q );
cDn = cDAlg ( Q2X1 ( Q1 , Q2 ) , Q2X2 ( Q1 , Q2 ));
for i = 1 : N
globKD (: , i ) = globKD (: , i ) + eta * m a r k E 0 T b d r(: , n ) .* ( cDn * (W ' .* gPhi1D { qOrd }(: , i , n )) );
end % for
end % for
globKD = reshape ( globKD ' , g . numT *N , 1);
end % f u n c t i o n
globKN = assembleGlobKN(g, markE0Tbdr, dDG, gNAlg) assembles the contributions of Neumann boundaries KN to the right-
hand side of Equation (3b) according to Sec. 3.6.12 which are returned in variable globKN. The input argument markE0Tbdr is described
in assembleGlobJD , dDG is the representation of the diffusion coefficient in the polynomial space (cf. projectAlg2DG ), and gNAlg is
a function_handle for the algebraic representation of gN (cf. main.m).
function globKN = a s s e m b l e G l o b K N(g , markE0Tbdr , dDG , gNAlg )
global gPhi1D
[K , N ] = size ( dDG ); p = ( sqrt (8* N +1) -3)/2;
qOrd = 2* p +1; [Q , W ] = q u a d R u l e 1 D( qOrd );
Q2X1 = @ ( X1 , X2 ) g . B (: ,1 ,1)* X1 + g . B (: ,1 ,2)* X2 + g . coordV0T (: ,1 ,1)* ones ( size ( X1 ));
Q2X2 = @ ( X1 , X2 ) g . B (: ,2 ,1)* X1 + g . B (: ,2 ,2)* X2 + g . coordV0T (: ,1 ,2)* ones ( size ( X1 ));
globKN = zeros (K , N );
for n = 1 : 3
[ Q1 , Q2 ] = gammaMap (n , Q );
gNAlgn = gNAlg ( Q2X1 ( Q1 , Q2 ) , Q2X2 ( Q1 , Q2 ));
Kkn = - m a r k E 0 T b d r(: , n ) .* g . sigE0T (: , n ) .* g . areaE0T (: , n );
for i = 1 : N
for l = 1 : N
integral = gNAlgn * ( W .* gPhi1D { qOrd }(: , i , n ) ' .* gPhi1D { qOrd }(: , l , n ) ' ) ' ;
globKN (: , i ) = globKN (: , i ) + Kkn .* dDG (: , l ) .* integral ;
end % for
end % for
end % for
globKN = reshape ( globKN ' , K *N ,1);
end % f u n c t i o n
globM = assembleGlobM(g, hatM) assembles the global mass matrix M according to Sec. 3.6.1 which is returned in variable globM. The
input argument hatM stores the local matrices M̂ as defined in (21) and can be computed by computeHatM .
function globM = a s s e m b l e G l o b M(g , hatM )
K = g . numT ;
globM = 2* kron ( spdiags ( g . areaT , 0 , K , K ) , hatM );
end % f u n c t i o n
23
globQ = assembleGlobQ(g, markE0Tint, hatQdiag, hatQoffdiag) assembles the matrices Qm , m ∈ {1, 2} according to
Sec. 3.6.5 which are returned in the 2 × 1 cell variable globQ. The input arguments hatQdiag and hatQoffdiag store the local ma-
trices Q̂m,diag = Ŝdiag and Q̂m,offdiag = Ŝoffdiag as given in (24) and (27), respectively. They can be computed by computeHatSdiag
and computeHatSoffdiag . Similarly to assembleGlobJD , the argument markE0Tint is a K × 3 logical array that marks each tri-
angle’s interior edges. Note the use of bsxfun(@times,...) to carry out the Hadamard product without building the full coefficient matrix.
function globQ = a s s e m b l e G l o b Q(g , markE0Tint , hatQdiag , h a t Q o f f d i a g)
K = g . numT ; N = size ( hatQdiag , 1);
globQ = cell (2 , 1);
globQ {1} = sparse ( K *N , K * N ); globQ {2} = sparse ( K *N , K * N );
for n = 1 : 3
Qkn = 0.5 * m a r k E 0 T i n t(: , n ) .* g . areaE0T (: , n ) .* g . sigE0T (: , n );
globQ {1} = globQ {1} + kron ( spdiags ( Qkn .* g . nuE0T (: ,n ,1) , 0 ,K , K ) , hatQdiag (: ,: , n ));
globQ {2} = globQ {2} + kron ( spdiags ( Qkn .* g . nuE0T (: ,n ,2) , 0 ,K , K ) , hatQdiag (: ,: , n ));
end % for
for nn = 1 : 3
Qkn = 0.5 * g . areaE0T (: , nn ) .* g . sigE0T (: , nn );
for np = 1 : 3
m a r k T i m e s Q k n = bsxfun ( @times , g . m a r k E 0 T E 0 T{ nn , np } , Qkn .* g . nuE0T (: , nn ,1));
globQ {1} = globQ {1} + kron ( markTimesQkn , h a t Q o f f d i a g(: ,: , nn , np ));
m a r k T i m e s Q k n = bsxfun ( @times , g . m a r k E 0 T E 0 T{ nn , np } , Qkn .* g . nuE0T (: , nn ,2));
globQ {2} = globQ {2} + kron ( markTimesQkn , h a t Q o f f d i a g(: ,: , nn , np ));
end % for
end % for
end % f u n c t i o n
globQN = assembleGlobQN(g, markE0Tbdr, hatQdiag) assembles the matrices QmN , m ∈ {1, 2} according to Sec. 3.6.9. It is es-
sentially the same routine as the diagonal part of assembleGlobQ , only with markE0Tbdr marking the Neumann boundary edges instead of
interior edges.
function globQN = a s s e m b l e G l o b Q N(g , markE0Tbdr , hatQdiag )
K = g . numT ; N = size ( hatQdiag , 1);
globQN = cell (2 , 1);
globQN {1} = sparse ( K *N , K * N ); globQN {2} = sparse ( K *N , K * N );
for n = 1 : 3
QNkn = m a r k E 0 T b d r(: , n ) .* g . areaE0T (: , n ) .* g . sigE0T (: , n );
globQN {1} = globQN {1} + kron ( spdiags ( QNkn .* g . nuE0T (: ,n ,1) , 0 ,K , K ) , hatQdiag (: ,: , n ));
globQN {2} = globQN {2} + kron ( spdiags ( QNkn .* g . nuE0T (: ,n ,2) , 0 ,K , K ) , hatQdiag (: ,: , n ));
end % for
end % f u n c t i o n
globR = assembleGlobR(g, markE0Tint, hatRdiag, hatRoffdiag, dDG) assembles the matrices Rm , m ∈ {1, 2} according to
Sec. 3.6.6 which are returned in the 2 × 1 cell variable globR. The input arguments hatRdiag and hatRoffdiag store the local matri-
ces R̂m,diag and R̂m,offdiag as defined in (28) and (29), respectively. They can be computed by computeHatRdiag and computeHatRoffdiag .
Each triangle’s interior edges are marked in markE0Tint as in assembleGlobQ . A representation of the diffusion coefficient in the polynomial
space is stored in dDG and can be computed by projectAlg2DG . The Hadamard product is carried out by bsxfun(@times,...). Note the
transposed application of Dkl (t) in the second part of the routine as a result of the diffusion coefficient being taken from the neighboring triangle T k+ .
function globR = a s s e m b l e G l o b R(g , markE0Tint , hatRdiag , hatRoffdiag , dDG )
[K , N ] = size ( dDG ); globR = cell (2 , 1);
globR {1} = sparse ( K *N , K * N ); globR {2} = sparse ( K *N , K * N );
for n = 1 : 3
Rkn = 0.5 * m a r k E 0 T i n t(: , n ) .* g . areaE0T (: , n ) .* g . sigE0T (: , n );
for l = 1 : N
globR {1} = globR {1} + kron ( spdiags ( Rkn .* g . nuE0T (: ,n ,1).* dDG (: , l ) ,0 ,K , K ) , hatRdiag (: ,: ,l , n ));
globR {2} = globR {2} + kron ( spdiags ( Rkn .* g . nuE0T (: ,n ,2).* dDG (: , l ) ,0 ,K , K ) , hatRdiag (: ,: ,l , n ));
end % for
end % for
for nn = 1 : 3
Rkn = 0.5 * g . areaE0T (: , nn ) .* g . sigE0T (: , nn );
for np = 1 : 3
m a r k E 0 T E 0 T t i m e s R k n 1 = bsxfun ( @times , g . m a r k E 0 T E 0 T{ nn , np } , Rkn .* g . nuE0T (: , nn ,1));
m a r k E 0 T E 0 T t i m e s R k n 2 = bsxfun ( @times , g . m a r k E 0 T E 0 T{ nn , np } , Rkn .* g . nuE0T (: , nn ,2));
for l = 1 : N
m a r k E 0 T E 0 T t i m e s D t i m e s R k n = bsxfun ( @times , m a r k E 0 T E 0 T t im esR kn 1 , dDG (: , l ). ' );
globR {1} = globR {1} + kron ( m a r k E 0 T E 0 T t i m e s Dt im es R kn , h a t R o f f d i a g(: ,: ,l , nn , np ));
m a r k E 0 T E 0 T t i m e s D t i m e s R k n = bsxfun ( @times , m a r k E 0 T E 0 T t im esR kn 2 , dDG (: , l ). ' );
globR {2} = globR {2} + kron ( m a r k E 0 T E 0 T t i m e s Dt im es R kn , h a t R o f f d i a g(: ,: ,l , nn , np ));
end % for
end % for
end % for
end % f u n c t i o n
24
globRD = assembleGlobRD(g, markE0Tbdr, hatRdiag, dDG) assembles the matrices RmD , m ∈ {1, 2} according to Sec. 3.6.7. It is
essentially the same routine as the first part of assembleGlobR but carried out only for the Dirichlet boundary edges marked by markE0Tbdr.
function globRD = a s s e m b l e G l o b R D(g , markE0Tbdr , hatRdiag , dDG )
[K , N ] = size ( dDG ); globRD = cell (2 , 1);
globRD {1} = sparse ( K *N , K * N ); globRD {2} = sparse ( K *N , K * N );
for n = 1 : 3
RDkn = m a r k E 0 T b d r(: , n ) .* g . areaE0T (: , n ) .* g . sigE0T (: , n );
for l = 1 : N
globRD {1} = globRD {1} + kron ( spdiags ( RDkn .* g . nuE0T (: ,n ,1).* dDG (: , l ) ,0 ,K , K ) , hatRdiag (: ,: ,l , n ));
globRD {2} = globRD {2} + kron ( spdiags ( RDkn .* g . nuE0T (: ,n ,2).* dDG (: , l ) ,0 ,K , K ) , hatRdiag (: ,: ,l , n ));
end
end % for
end % f u n c t i o n
globS = assembleGlobS(g, markE0Tint, hatSdiag, hatSoffdiag, eta) assembles the matrix S according to Sec. 3.6.4
which is returned in the variable globS. The arguments are the same as for assembleGlobQ plus the penalty coefficient η given in eta.
function globS = a s s e m b l e G l o b S(g , markE0Tint , hatSdiag , hatSoffdiag , eta )
K = g . numT ; N = size ( hatSdiag , 1);
globS = sparse ( K *N , K * N );
for n = 1 : 3
globS = globS + eta * kron ( spdiags ( m a r k E 0 T i n t(: , n ) ,0 ,K , K ) , hatSdiag (: ,: ,n ));
end % for
for nn = 1 : 3
for np = 1 : 3
globS = globS - eta * kron ( g . m a r k E 0 T E 0 T{ nn , np } , h a t S o f f d i a g(: ,: ,nn , np ));
end % for
end % for
end % f u n c t i o n
globSD = assembleGlobSD(g, markE0Tbdr, hatSdiag, eta) assembles the matrix SD according to Sec. 3.6.8. It is similar to the
first part of assembleGlobS but carried out for Dirichlet boundary edges, which are marked in markE0Tbdr .
function globSD = a s s e m b l e G l o b S D(g , markE0Tbdr , hatSdiag , eta )
K = g . numT ; N = size ( hatSdiag , 1);
globSD = sparse ( K *N , K * N );
for n = 1 : 3
globSD = globSD + eta * kron ( spdiags ( m a r k E 0 T b d r(: , n ) ,0 ,K , K ) , hatSdiag (: ,: ,n ));
end % for
end % f u n c t i o n
computeBasesOnQuad(N) evaluates the basis functions and their gradients in the quadrature points on the reference triangle T̂ according to
Sec. 3.8 for all required orders as described in Sec. 3.3. It computes the values of ϕ̂i ( x̂), ∇ϕ̂i ( x̂), ϕ̂i ◦ γ̂(s), and ϕ̂i ◦ ϑ̂n− n+ ◦ γ̂(s) and stores them
in global arrays gPhi2D, gGradPhi2D , gPhi1D, and gThetaPhi1D, respectively, which are 2p + 1 × 1 cell variables.
function c o m p u t e B a s e s O n Q u a d( N )
global gPhi2D g G r a d P h i 2 D gPhi1D g T h e t a P h i 1 D
p = ( sqrt (8* N +1) -3)/2;
if p > 0 , r e q u i r e d O r d e r s = [2*p , 2* p +1]; else r e q u i r e d O r d e r s = 1; end
gPhi2D = cell ( max ( r e q u i r e d O r d e r s ) ,1); g G r a d P h i 2 D = cell ( max ( r e q u i r e d O r d e r s ) ,1);
gPhi1D = cell ( max ( r e q u i r e d O r d e r s ) ,1); g T h e t a P h i 1 D = cell ( max ( r e q u i r e d O r d e r s ) ,1);
for it = 1 : length ( r e q u i r e d O r d e r s)
ord = r e q u i r e d O r d e r s( it );
[ Q1 , Q2 , ~] = q u a d R u l e 2 D( ord );
gPhi2D { ord } = zeros ( length ( Q1 ) , N );
for i = 1 : N
gPhi2D { ord }(: , i ) = phi (i , Q1 , Q2 );
end % for
g G r a d P h i 2 D{ ord } = zeros ( length ( Q1 ) , N , 2);
for m = 1 : 2
for i = 1 : N
g G r a d P h i 2 D{ ord }(: , i , m ) = gradPhi (i , m , Q1 , Q2 );
end % for
end % for
[Q , ~] = q u a d R u l e 1 D( ord );
gPhi1D { ord } = zeros ( length ( Q ) , N , 3);
for nn = 1 : 3
[ Q1 , Q2 ] = gammaMap ( nn , Q );
for i = 1 : N
gPhi1D { ord }(: , i , nn ) = phi (i , Q1 , Q2 );
end
for np = 1 : 3
[ QP1 , QP2 ] = theta ( nn , np , Q1 , Q2 );
25
for i = 1 : N
g T h e t a P h i 1 D{ ord }(: , i , nn , np ) = phi (i , QP1 , QP2 );
end % for
end % for
end % for
end % for
end % f u n c t i o n
hatG = computeHatG(N) computes the local matrix Ĝ (multidimensional array) on the reference triangle T̂ as given in (23) using quadrature.
function hatG = c o m p u t e H a t G( N )
global gPhi2D g G r a d P h i 2 D
hatG = zeros (N , N , N , 2); % [ N x N x N x 2]
if N > 1 % p > 0
p = ( sqrt (8* N +1) -3)/2; qOrd = max (2* p , 1); [~ ,~ ,W ] = q u a d R u l e 2 D( qOrd );
for i = 1 : N
for j = 1 : N
for l = 1 : N
for m = 1 : 2
hatG (i ,j ,l , m ) = sum ( W ' .* g G r a d P h i 2 D{ qOrd }(: , i , m ) .* gPhi2D { qOrd }(: , j ) .* gPhi2D { qOrd }(: , l ) );
end % for
end % for
end % for
end % for
end % f u n c t i o n
hatH = computeHatH(N) computes the local matrix Ĥ (multidimensional array) on the reference triangle T̂ as given in (22) using quadrature.
function hatH = c o m p u t e H a t H( N )
global gPhi2D g G r a d P h i 2 D
hatH = zeros (N , N , 2); % [ N x N x 2]
if N > 1 % p > 0
p = ( sqrt (8* N +1) -3)/2; qOrd = max (2* p , 1); [~ ,~ ,W ] = q u a d R u l e 2 D( qOrd );
for i = 1 : N
for j = 1 : N
for m = 1 : 2
hatH (i , j , m ) = sum ( W ' .* g G r a d P h i 2 D{ qOrd }(: , i , m ) .* gPhi2D { qOrd }(: , j ) );
end % for
end % for
end % for
end % f u n c t i o n
hatM = computeHatM(N) computes the local mass matrix M̂ on the reference triangle T̂ as given in (21) using quadrature.
function hatM = c o m p u t e H a t M( N )
global gPhi2D
p = ( sqrt (8* N +1) -3)/2; qOrd = max (2*p , 1); [~ ,~ ,W ] = q u a d R u l e 2 D( qOrd );
hatM = zeros ( N ); % [ N x N ]
for i = 1 : N
for j = 1 : N
hatM (i , j ) = sum ( W ' .* gPhi2D { qOrd }(: , i ) .* gPhi2D { qOrd }(: , j ) );
end % for
end % for
end % f u n c t i o n
hatRdiag = computeHatRdiag(N) computes the local matrix R̂diag (multidimensional array) on the reference triangle T̂ as given in (28)
using quadrature.
function hatRdiag = c o m p u t e H a t R d i a g( N )
global gPhi1D
p = ( sqrt (8* N +1) -3)/2; qOrd = max (2* p +1 ,1); [~ , W ] = q u a d R u l e 1 D( qOrd );
hatRdiag = zeros (N , N , N , 3); % [ N x N x N x 3]
for n = 1 : 3 % 3 edges
for l = 1 : N % N b a s i s f c t s for D ( t )
for i = 1 : N
for j = 1 : N
hatRdiag (i ,j ,l , n ) = sum (W ' .* gPhi1D { qOrd }(: , i , n ) .* gPhi1D { qOrd }(: , l , n ) .* gPhi1D { qOrd }(: , j , n ));
end % for
end % for
end % for
end % for
end % f u n c t i o n
26
hatRoffdiag = computeHatRoffdiag(N) computes the local matrix R̂offdiag (multidimensional array) on the reference triangle T̂ as given
in (29) using quadrature.
function h a t R o f f d i a g = c o m p u t e H a t R o f f d i a g( N )
global gPhi1D g T h e t a P h i 1 D
p = ( sqrt (8* N +1) -3)/2; qOrd = 2* p +1; [~ , W ] = q u a d R u l e 1 D( qOrd );
h a t R o f f d i a g = zeros (N ,N ,N ,3 ,3); % [ N x N x N x 3 x 3]
for nn = 1 : 3 % 3 edges
for np = 1 : 3
for l = 1 : N
for i = 1 : N
for j = 1 : N
h a t R o f f d i a g(i , j , l , nn , np ) = sum ( W ' .* gPhi1D { qOrd }(: , i , nn ) .* ...
g T h e t a P h i 1 D{ qOrd }(: , l , nn , np ) .* g T h e t a P h i 1 D{ qOrd }(: , j , nn , np ) );
end % for
end % for
end % for
end % for
end % for
end
hatSdiag = computeHatSdiag(N) computes the local matrix Ŝdiag (multidimensional array) on the reference triangle T̂ as given in (24)
using quadrature.
function hatSdiag = c o m p u t e H a t S d i a g( N )
global gPhi1D
p = ( sqrt (8* N +1) -3)/2; qOrd = 2* p +1; [~ , W ] = q u a d R u l e 1 D( qOrd );
hatSdiag = zeros (N , N , 3); % [ N x N x 3]
for n = 1 : 3 % 3 edges
for i = 1 : N
for j = 1 : N
hatSdiag (i , j , n ) = sum ( W ' .* gPhi1D { qOrd }(: , i , n ) .* gPhi1D { qOrd }(: , j , n ) );
end % for
end % for
end % for
end % f u n c t i o n
hatSoffdiag = computeHatSoffdiag(N) computes the local matrix Ŝoffdiag (multidimensional array) on the reference triangle T̂ as given
in (27) using quadrature.
function h a t S o f f d i a g = c o m p u t e H a t S o f f d i a g( N )
global gPhi1D g T h e t a P h i 1 D
p = ( sqrt (8* N +1) -3)/2; qOrd = 2* p +1; [~ , W ] = q u a d R u l e 1 D( qOrd );
h a t S o f f d i a g = zeros (N , N , 3 , 3); % [ N x N x 3 x 3]
for nn = 1 : 3
for np = 1 : 3
for i = 1 : N
for j = 1 : N
h a t S o f f d i a g(i , j , nn , np ) = sum ( W ' .* gPhi1D { qOrd }(: , i , nn ) .* g T h e t a P h i 1 D{ qOrd }(: , j , nn , np ) );
end % for
end % for
end % for
end % for
end % f u n c t i o n
err = computeL2Error(g, cDG, cAlg, qOrd) computes the L2 -error according to Sec. 3.5 using a quadrature of order qOrd, where
g is a struct generated by generateGridData , cDG is the K × N representation matrix of the approximate solution, cAlg is the algebraic
formulation of the exact solution (see example below).
Example: N = 6; h = 1/2; qOrd = 5;
g = domainSquare(h); computeBasesOnQuad(N); hatM = computeHatM(N);
cAlg = @(X1, X2) sin(X1).*sin(X2); cDG = projectAlg2DG(g,cAlg,N,qOrd,hatM);
computeL2Error(g, cDG, cAlg, qOrd)
function err = c o m p u t e L 2 E r r o r(g , cDG , cAlg , qOrd )
global gPhi2D
N = size ( cDG , 2); qOrd = max ( qOrd ,1);
[ Q1 , Q2 , W ] = q u a d R u l e 2 D( qOrd );
R = length ( W );
X1 = kron ( g . B (: ,1 ,1) , Q1 )+ kron ( g . B (: ,1 ,2) , Q2 )+ kron ( g . coordV0T (: ,1 ,1) , ones (1 , R ));
X2 = kron ( g . B (: ,2 ,1) , Q1 )+ kron ( g . B (: ,2 ,2) , Q2 )+ kron ( g . coordV0T (: ,1 ,2) , ones (1 , R ));
c E x O n Q u a d P t s = cAlg ( X1 , X2 ); % [ K x R ]
c A p p r x O n Q u a d P t s = cDG * gPhi2D { qOrd } ' ; % [ K x R ] = [ K x N ] * [ N x R ]
err = sqrt (2* dot (( c A p p r x O n Q u a d P t s - c E x O n Q u a d P t s ).^2 * W . ' , g . areaT ));
27
end % f u n c t i o n
g = domainCircle(h) This method uses the grid generator Gmsh [19]. A system call to Gmsh generates the ASCII
file domainCircle.mesh based on geometry information of the domain Ω stored in domainCircle.geo. The basic grid data is extracted
from domainCircle.mesh to call the routine generateGridData (cf. Sec. 3.1) and to set the boundary IDs (from 1 to 4).
function g = d o m a i n C i r c l e( h )
% % g e n e r a t i o n of d o m a i n C i r c l e. mesh using d o m a i n C i r c l e. geo
cmd = sprintf ( ' gmsh -2 - format mesh - clscale % f -o " d o m a i n C i r c l e. mesh " " d o m a i n C i r c l e. geo " ' , h )
system ( cmd );
% % e x t r a c t data from d o m a i n C i r c l e. mesh
fid = fopen ( ' d o m a i n C i r c l e. mesh ' , ' r ' );
tline = fgets ( fid );
while ischar ( tline )
if strfind ( tline , ' Vertices ' )
numV = fscanf ( fid , ' % d ' , [1 , 1]);
coordV = reshape ( fscanf ( fid , ' % f ' ) , 4 , numV ) ' ; coordV (: , 3:4) = [];
end % if
if strfind ( tline , ' Edges ' )
numEbdry = fscanf ( fid , ' % d ' , [1 , 1]);
tmp = reshape ( fscanf ( fid , ' % f ' ) , 3 , numEbdry ) ' ;
V0Ebdry = tmp (: , 1:2); idEbdry = tmp (: , 3);
end % if
if strfind ( tline , ' T r i a n g l e s ' )
numT = fscanf ( fid , ' % d ' , [1 , 1]);
V0T = reshape ( fscanf ( fid , ' % f ' ) , 4 , numT ) ' ; V0T (: , 4) = [];
end % if
tline = fgets ( fid );
end % while
fclose ( fid );
% % g e n e r a t e lists and set b o u n d a r y IDs
g = g e n e r a t e G r i d D a t a( coordV , V0T );
g . idE = zeros ( g . numE , 1);
g . idE ( g . V2E ( sub2ind ([ g . numV , g . numV ] , V0Ebdry (: ,1) , V0Ebdry (: ,2)))) = idEbdry ;
g . idE0T = g . idE ( g . E0T ); % local edge IDs
end % f u n c t i o n
domainCircle.geo Geometry description of a circle with center (0, 0) and radius 0.5 serving as input for the grid generator Gmsh [19] which
is called by the routine domainCircle .
cx = 0.0; cy = 0.0; r = 0.5;
Point (0) = { cx , cy , 0.0 , 1.0}; Point (1) = { cx -r , cy , 0.0 , 1.0}; Point (2) = { cx , cy +r , 0.0 , 1.0};
Point (3) = { cx +r , cy , 0.0 , 1.0}; Point (4) = { cx , cy -r , 0.0 , 1.0};
Circle (5) = {2 , 0 , 1}; Circle (6) = {3 , 0 , 2}; Circle (7) = {4 , 0 , 3}; Circle (8) = {1 , 0 , 4};
Line Loop (9) = {5 , 6 , 7 , 8}; // boundary edge IDs will be 5 to 8
Plane Surface (10) = {9};
g = domainSquare(h) Friedrichs–Keller triangulation on the unit square. The input argument h specifies the upper bound for the heights
of the triangles (not for the diameters). The output variable g representing the triangulation Th is of type struct according to Sec. 3.1. The
boundary IDs are set from 1 to 4.
function g = d o m a i n S q u a r e( h )
dim = ceil (1/ h ); % n u m b e r of edges per side of the unit s q u a r e
h = 1/ dim ;
% coordV
[X , Y ] = meshgrid (0: h :1);
Xlist = reshape (X , length ( X )^2 , 1); Ylist = reshape (Y , length ( X )^2 , 1);
coordV = [ Xlist , Ylist ];
% V0T
pat1 = [1 , dim +2 ,2]; % p a t t e r n of " lower - left " t r i a n g l e s
V0T1 = repmat ( pat1 , dim *( dim +1) , 1) + repmat ((0: dim *( dim +1) -1) ' , 1 , 3);
V0T1 ( dim +1 : dim +1 : dim *( dim +1) , :) = [];
pat2 = [ dim +2 , dim +3 ,2];
V0T2 = repmat ( pat2 , dim *( dim +1) , 1) + repmat ((0: dim *( dim +1) -1) ' , 1 , 3);
V0T2 ( dim +1 : dim +1 : dim *( dim +1) , :) = [];
% % g e n e r a t e grid data and b o u n d a r y IDs
g = g e n e r a t e G r i d D a t a( coordV , [ V0T1 ; V0T2 ]);
g . idE = zeros ( g . numE , 1);
g . idE ( g . baryE (: , 2) == 0) = 1; % south
g . idE ( g . baryE (: , 1) == 1) = 2; % east
g . idE ( g . baryE (: , 2) == 1) = 3; % north
g . idE ( g . baryE (: , 1) == 0) = 4; % west
g . idE0T = g . idE ( g . E0T ); % local edge IDs
end % f u n c t i o n
28
g = domainPolygon(X1, X2, h) triangulates a polygonally bounded domain Ω, where X1 and X2 are lists of x1 and x2 coordinates of the
corners. The output variable g representing the triangulation Th is of type struct according to Sec. 3.1. The boundary IDs are set to the values
produced by initmesh. This method uses the MATLAB grid generator initmesh.
function g = d o m a i n P o l y g o n( X1 , X2 , h )
gd = [2; length ( X1 (:)); X1 (:); X2 (:)]; % g e o m e t r y d e s c r i p t i o n
sf = ' polygon ' ; % set f o r m u l a
ns = double ( ' polygon ' ) ' ; % name space
[p , e , t ] = initmesh ( decsg ( gd , sf , ns ) , ' Hmax ' , h );
g = g e n e r a t e G r i d D a t a(p ' , t (1:3 , :) ' );
g . idE = zeros ( g . numE , 1);
g . idE ( g . V2E ( sub2ind ( size ( g . V2E ) , e (1 ,:) , e (2 ,:)))) = e (5 ,:);
g . idE0T = g . idE ( g . E0T ); % local edge IDs
end % f u n c t i o n
[X1, X2] = gammaMap(n, S) evaluates the parametrization of the nth edge Ên of the reference triangle T̂ at parameter values specified by
a list of parameters S, cf. (25).
function [ X1 , X2 ] = gammaMap (n , S )
S = S (:) ' ;
switch n
case 1 , X1 = 1 - S ; X2 = S ;
case 2 , X1 = zeros ( size ( S )); X2 = S ;
case 3 , X1 = S ; X2 = zeros ( size ( S ));
end % s w i t c h
end % f u n c t i o n
g = generateGridData(coordV, V0T) assembles lists describing the geometric and topological properties of a triangulation Th according
to Sec. 3.1 and stores them in the output variable g of type struct. The input arguments are the array coordV of dimension #V × 2 that contains
the x1 and x2 coordinates of the grid vertices (using a global index) and the array V0T of dimension #Th × 3 storing the global vertex indices for
each triangle with a counter-clockwise ordering. A usage example is found in Sec. 3.10. Note that the lists g.idE and g.idE0T (cf. Tab. 1)
storing the global and local edge indices are not generated and have to be defined manually after calling generateGridData .
function g = g e n e r a t e G r i d D a t a( coordV , V0T )
g . coordV = coordV ;
g . V0T = V0T ;
g . numT = size ( g . V0T , 1);
g . numV = size ( g . coordV , 1);
% the f o l l o w i n g i m p l i c i t e l y d e f i n e s the signs of the edges
g . V2T = sparse ( g . V0T (: , [1 2 3 1 2 3 1 2 3]) , g . V0T (: , [2 3 1 2 3 1 2 3 1]) , ...
[(1: g . numT ) ' , zeros ( g . numT ,3) ,(1:g . numT ) ' , zeros ( g . numT ,3) ,(1: g . numT ) ' ] , g . numV , g . numV );
% the f o l l o w i n g i m p l i c i t e l y d e f i n e s the edge n u m b e r s
[r , c ] = find ( triu ( g . V2T + g . V2T ' ));
g . V2E = sparse (r , c , 1 : size (r , 1) , g . numV , g . numV );
g . V2E = g . V2E + g . V2E ' ;
idxE = full ( g . V2E ( sub2ind ([ g . numV , g . numV ] , g . V0T ( end : -1:1 ,[1 ,2 ,3]) , g . V0T ( end : -1:1 ,[2 ,3 ,1])))) ' ;
g . V0E ( idxE (:) , 1) = reshape ( g . V0T ( end : -1:1 , [1 ,2 ,3]) ' , 3* g . numT , 1);
g . V0E ( idxE (:) , 2) = reshape ( g . V0T ( end : -1:1 , [2 ,3 ,1]) ' , 3* g . numT , 1);
g . T0E ( idxE (:) , 1) = reshape ( full ( g . V2T ( sub2ind ([ g . numV , g . numV ] , ...
g . V0T ( end : -1:1 ,[1 ,2 ,3]) , g . V0T ( end : -1:1 ,[2 ,3 ,1])))) ' , 3* g . numT , 1);
g . T0E ( idxE (:) , 2) = reshape ( full ( g . V2T ( sub2ind ([ g . numV , g . numV ] , ...
g . V0T ( end : -1:1 ,[2 ,3 ,1]) , g . V0T ( end : -1:1 ,[1 ,2 ,3])))) ' , 3* g . numT , 1);
g . numE = size ( g . V0E , 1);
g . vecE = g . coordV ( g . V0E (: , 2) , :) - g . coordV ( g . V0E (: , 1) , :);
g . areaE = ( g . vecE (: , 1).^2 + g . vecE (: , 2 ) . ^ 2 ) . ^ ( 1 / 2 ) ;
g . nuE = g . vecE * [0 , -1; 1 ,0] ./ g . areaE (: , [1 , 1]);
g . areaT = ...
( g . coordV ( g . V0T (: ,1) ,1).*g . coordV ( g . V0T (: ,2) ,2) + g . coordV ( g . V0T (: ,2) ,1).*g . coordV ( g . V0T (: ,3) ,2) ...
+ g . coordV ( g . V0T (: ,3) ,1).*g . coordV ( g . V0T (: ,1) ,2) - g . coordV ( g . V0T (: ,1) ,1).*g . coordV ( g . V0T (: ,3) ,2) ...
- g . coordV ( g . V0T (: ,2) ,1).*g . coordV ( g . V0T (: ,1) ,2) - g . coordV ( g . V0T (: ,3) ,1).*g . coordV ( g . V0T (: ,2) ,2) )/2;
g . baryT = ( g . coordV ( g . V0T (: ,1) ,:)+g . coordV ( g . V0T (: ,2) ,:)+g . coordV ( g . V0T (: ,3) ,:))/3;
g . E0T = full ( g . V2E ( sub2ind ([ g . numV , g . numV ] , g . V0T (: ,[2 ,3 ,1]) , g . V0T (: ,[3 ,1 ,2]))));
g . areaE0T = g . areaE ( g . E0T );
g . sigE0T = 1 -2*( bsxfun ( @eq , reshape ( g . T0E ( g . E0T ,2) , g . numT ,3) ,(1: g . numT ) ' ));
g . baryE = 0.5 * ( g . coordV ( g . V0E (: , 1) , :) + g . coordV ( g . V0E (: , 2) , :));
for n = 1 : 3
g . coordV0T (: , n , 1) = g . coordV ( g . V0T (: , n ) , 1) ' ; g . coordV0T (: , n , 2) = g . coordV ( g . V0T (: , n ) , 2) ' ;
g . baryE0T (: , n , 1) = g . baryE ( g . E0T (: , n ) , 1) ' ; g . baryE0T (: , n , 2) = g . baryE ( g . E0T (: , n ) , 2) ' ;
g . nuE0T (: , n , 1) = g . nuE ( g . E0T (: , n ) , 1) ' ; g . nuE0T (: , n , 2) = g . nuE ( g . E0T (: , n ) , 2) ' ;
Tn = g . sigE0T (: , n ) == 1; Tp = ~ Tn ;
g . E0E ( g . E0T ( Tn , n ) , 1) = n ; g . E0E ( g . E0T ( Tp , n ) , 2) = n ;
end % for
for m = 1 : 2
g . B (: , m , 1) = g . coordV0T (: , 2 , m ) - g . coordV0T (: , 1 , m );
g . B (: , m , 2) = g . coordV0T (: , 3 , m ) - g . coordV0T (: , 1 , m );
29
end % for
markEint = g . E0E (: , 2) ~= 0; % mark i n t e r i o r edges
g . m a r k E 0 T E 0 T = cell (3 , 3);
for nn = 1 : 3
for np = 1 : 3
g . m a r k E 0 T E 0 T{ nn , np } = sparse ( g . numT , g . numT );
markEn = g . E0E (: , 1) == nn ; markEp = g . E0E (: , 2) == np ;
idx = markEn & markEp & markEint ;
g . m a r k E 0 T E 0 T{ nn , np }( sub2ind ([ g . numT , g . numT ] , g . T0E ( idx , 1) , g . T0E ( idx , 2))) = 1;
markEn = g . E0E (: , 2) == nn ; markEp = g . E0E (: , 1) == np ;
idx = markEn & markEp & markEint ;
g . m a r k E 0 T E 0 T{ nn , np }( sub2ind ([ g . numT , g . numT ] , g . T0E ( idx , 2) , g . T0E ( idx , 1))) = 1;
end % for
end % for
end % f u n c t i o n
ret = gradPhi(i, m, X1, X2) evaluates the mth component of the gradient of the ith basis function ϕ̂i on the reference triangle T̂
(cf. Sec. 2.4.1) at points specified by a list of x̂1 coordinates X1 and x̂2 coordinates X2.
function ret = gradPhi (i , m , X1 , X2 )
switch m
case 1
switch i
case 1 , ret =
zeros ( size ( X1 ));
case 2 , ret =
-6* ones ( size ( X1 ));
case 3 , ret =
-2* sqrt (3)* ones ( size ( X1 ));
case 4 , ret =
sqrt (6)*(20* X1 - 8);
case 5 , ret =
sqrt (3)*(10* X1 - 4);
case 6 , ret =
6* sqrt (5)*(3* X1 + 4* X2 - 2);
case 7 , ret =
2* sqrt ( 2 ) * ( 1 5 + ( - 9 0 + 1 0 5 *X1 ).* X1 );
case 8 , ret =
2* sqrt (6)*(13+( -66+63* X1 ).* X1 +( -24+84* X1 ).* X2 );
case 9 , ret =
2* sqrt (10)*(9+( -30+21* X1 ).* X1 +( -48+84* X1 +42* X2 ).* X2 );
case 10 , ret =
2* sqrt (14)*(3+( -6+3* X1 ).* X1 +( -24+24* X1 +30* X2 ).* X2 );
case 11 , ret =
sqrt ( 1 0 ) * ( - 2 4 + ( 2 5 2 + ( - 6 7 2 + 5 0 4 *X1 ).* X1 ).* X1 );
case 12 , ret =
sqrt ( 3 0 ) * ( - 2 2 + ( 2 1 0 + ( - 5 0 4 + 3 3 6 *X1 ).* X1 ).* X1 +(42+( -336+504* X1 ).* X1 ).* X2 );
case 13 , ret =
5* sqrt (2)*( -18+(138+( -264+144* X1 ).* X1 ).* X1 ...
+(102+( -624+648* X1 ).* X1 +( -96+432* X1 ).* X2 ).* X2 );
case 14 , ret = sqrt (70)*( -12+(60+( -84+36* X1 ).* X1 ).* X1 ...
+(132+( -456+324* X1 ).* X1 +( -300+540* X1 +180* X2 ).* X2 ).* X2 );
case 15 , ret = 3* sqrt (10)*( -4+(12+( -12+4* X1 ).* X1 ).* X1 ...
+(60+( -120+60* X1 ).* X1 +( -180+180* X1 +140* X2 ).* X2 ).* X2 );
end
case 2
switch i
case 1 , ret = zeros ( size ( X1 ));
case 2 , ret = zeros ( size ( X1 ));
case 3 , ret = -4* sqrt (3)* ones ( size ( X1 ));
case 4 , ret = zeros ( size ( X1 ));
case 5 , ret = 2* sqrt (3)*( -15* X2 + 6);
case 6 , ret = 6* sqrt (5)*(4* X1 + 3* X2 - 2);
case 7 , ret = zeros ( size ( X1 ));
case 8 , ret = 2* sqrt (6)*(2+( -24+42* X1 ).* X1 );
case 9 , ret = 2* sqrt (10)*(6+( -48+42* X1 ).* X1 +( -12+84* X1 ).* X2 );
case 10 , ret = 2* sqrt ( 1 4 ) * ( 1 2 + (- 2 4 + 1 2 * X1 ).* X1 +( -60+60* X1 +60* X2 ).* X2 );
case 11 , ret = zeros ( size ( X1 ));
case 12 , ret = sqrt (30)*( -2+(42+( -168+168* X1 ).* X1 ).* X1 );
case 13 , ret = 5* sqrt (2)*( -6+(102+( -312+216* X1 ).* X1 ).* X1 ...
+(12+( -192+432* X1 ).* X1 ).* X2 );
case 14 , ret = sqrt ( 7 0 ) * ( - 1 2 + ( 1 3 2 + ( - 2 2 8 + 1 0 8 *X1 ).* X1 ).* X1 ...
+(60+( -600+540* X1 ).* X1 +( -60+540* X1 ).* X2 ).* X2 );
case 15 , ret = 3* sqrt (10)*( -20+(60+( -60+20* X1 ).* X1 ).* X1 ...
+(180+( -360+180* X1 ).* X1 +( -420+420* X1 +280* X2 ).* X2 ).* X2 );
end % s w i t c h
end % s w i t c h
end % f u n c t i o n
main.m This is the main script to solve Equation (2) which can be used as a template for further modifications.
function main ()
%% parameters
hmax = 1/8; % m a x i m u m edge l e n g t h of t r i a n g l e
p = 2; % local p o l y n o m i a l d e g r e e
tEnd = pi /2; % end time
numSteps = 10; % n u m b e r of time steps
isVisGrid = false ; % v i s u a l i z a t i o n of grid
30
isVisSol = true ; % v i s u a l i z a t i o n of s o l u t i o n
eta = 1; % p e n a l t y p a r a m e t e r ( eta >0)
% % c o e f f i c i e n t s and b o u n d a r y data
c0 = @ ( x1 , x2 ) sin ( x1 ).* cos ( x2 );
dAlg = @ (t , x1 , x2 ) ( x1 <3/4& x1 >1/4& x2 <3/4& x2 >1/4) + 0.01;
fAlg = @ (t , x1 , x2 ) 0.1* t *( x1 == x1 );
cDAlg = @ (t , x1 , x2 ) sin (2* pi * x2 + t );
gNAlg = @ (t , x1 , x2 ) x2 ;
%% triangulation
g = d o m a i n S q u a r e( hmax );
if isVisGrid , v i s u a l i z e G r i d( g ); end
%% globally constant parameters
K = g . numT ; % n u m b e r of t r i a n g l e s
N = nchoosek ( p + 2 , p ); % n u m b e r of local DOFs
tau = tEnd / numSteps ; % time step size
m a r k E 0 T i n t = g . idE0T == 0; % [ K x 3] mark local edges that are i n t e r i o r
m a r k E 0 T b d r N = g . idE0T == 1 | g . idE0T == 3; % [ K x 3] mark local edges on the N e u m a n n b o u n d a r y
m a r k E 0 T b d r D = ~( m a r k E 0 T i n t | m a r k E 0 T b d r N); % [ K x 3] mark local edges on the D i r i c h l e t b o u n d a r y
% % l o o k u p table for basis f u n c t i o n
c o m p u t e B a s e s O n Q u a d( N );
% % c o m p u t a t i o n of m a t r i c e s on the r e f e r e n c e t r i a n g l e
hatM = c o m p u t e H a t M( N );
hatG = c o m p u t e H a t G( N );
hatH = c o m p u t e H a t H( N );
hatRdiag = c o m p u t e H a t R d i a g( N );
h a t R o f f d i a g = c o m p u t e H a t R o f f d i a g( N );
hatSdiag = c o m p u t e H a t S d i a g( N );
h a t S o f f d i a g = c o m p u t e H a t S o f f d i a g( N );
% % a s s e m b l y of time - i n d e p e n d e n t g l o b a l m a t r i c e s
globM = a s s e m b l e G l o b M(g , hatM );
globH = a s s e m b l e G l o b H(g , hatH );
globQ = a s s e m b l e G l o b Q(g , markE0Tint , hatSdiag , h a t S o f f d i a g);
globQN = a s s e m b l e G l o b Q N(g , markE0TbdrN , hatSdiag );
globS = a s s e m b l e G l o b S(g , markE0Tint , hatSdiag , hatSoffdiag , eta );
globSD = a s s e m b l e G l o b S D(g , markE0TbdrD , hatSdiag , eta );
sysW = [ sparse (2* K *N ,3* K * N ) ; sparse ( K *N ,2* K * N ) , globM ];
% % i n i t i a l data
cDG = p r o j e c t A l g 2 D G(g , c0 , 2* p , hatM );
sysY = [ zeros (2* K *N ,1) ; reshape ( cDG ' , K *N , 1) ];
% % time s t e p p i n g
for nStep = 1 : numSteps
t = nStep * tau ;
% % L2 p r o j e c t i o n s of a l g e b r a i c c o e f f i c i e n t s
dDG = p r o j e c t A l g 2 D G(g , @ ( x1 , x2 ) dAlg (t , x1 , x2 ) , 2* p , hatM );
fDG = p r o j e c t A l g 2 D G(g , @ ( x1 , x2 ) fAlg (t , x1 , x2 ) , 2* p , hatM );
% % a s s e m b l y of time - d e p e n d e n t g l o b a l m a t r i c e s
globG = a s s e m b l e G l o b G(g , hatG , dDG );
globR = a s s e m b l e G l o b R(g , markE0Tint , hatRdiag , hatRoffdiag , dDG );
% % a s s e m b l y of D i r i c h l e t b o u n d a r y c o n t r i b u t i o n s
globRD = a s s e m b l e G l o b R D(g , markE0TbdrD , hatRdiag , dDG );
globJD = a s s e m b l e G l o b J D(g , markE0TbdrD , @ ( x1 , x2 ) cDAlg (t , x1 , x2 ) , N );
globKD = a s s e m b l e G l o b K D(g , markE0TbdrD , @ ( x1 , x2 ) cDAlg (t , x1 , x2 ) , eta , N );
% % a s s e m b l y of N e u m a n n b o u n d a r y c o n t r i b u t i o n s
globKN = a s s e m b l e G l o b K N(g , markE0TbdrN , dDG , @ ( x1 , x2 ) gNAlg (t , x1 , x2 ));
% % a s s e m b l y of the s o u r c e c o n t r i b u t i o n
globL = globM * reshape ( fDG ' , K *N , 1);
% % b u i l d i n g and s o l v i n g the s y s t e m
sysA = [ globM , sparse ( K *N , K * N ) , globH {1}+ globQ {1}+ globQN {1};
sparse ( K *N , K * N ) , globM , globH {2}+ globQ {2}+ globQN {2};
globG {1}+ globR {1}+ globRD {1} , globG {2}+ globR {2}+ globRD {2} , globS + globSD ];
sysV = [ globJD {1}; globJD {2}; globKD + globKN + globL ];
sysY = ( sysW + tau * sysA ) \ ( sysW * sysY + tau * sysV );
%% visualization
if isVisSol && p <= 2
cDG = reshape ( sysY (2* K * N +1 : 3* K * N ) , N , K ) ' ;
c L a g r a n g e = p r o j e c t D G 2 L a g r a n g e( cDG );
v i s u a l i z e D a t a(g , cLagrange , ' c_h ' , ' solution ' , nStep )
end % if
end % for
end % f u n c t i o n
ret = phi(i, X1, X2) evaluates the ith basis function ϕ̂i on the reference triangle T̂ (cf. Sec. 2.4.1) at points specified by a list of x̂1 coor-
dinates X1 and x̂2 coordinates X2.
function ret = phi (i , X1 , X2 )
switch i
31
case 1, ret =
sqrt (2)* ones ( size ( X1 ));
case 2, ret =
2 - 6* X1 ;
case 3, ret =
2* sqrt (3)*(1 - X1 - 2* X2 );
case 4, ret =
sqrt ( 6 ) * ( ( 1 0 *X1 - 8).* X1 + 1);
case 5, ret =
sqrt (3)*((5* X1 - 4).* X1 + ( -15* X2 + 12).* X2 - 1);
case 6, ret =
3* sqrt (5)*((3* X1 + 8* X2 - 4).* X1 + (3* X2 - 4).* X2 + 1);
case 7, ret =
2* sqrt (2)*( -1+(15+( -45+35* X1 ).* X1 ).* X1 );
case 8, ret =
2* sqrt (6)*( -1+(13+( -33+21* X1 ).* X1 ).* X1 +(2+( -24+42* X1 ).* X1 ).* X2 );
case 9, ret =
2* sqrt (10)*( -1+(9+( -15+7* X1 ).* X1 ).* X1 +(6+( -48+42* X1 ).* X1 +( -6+42* X1 ).* X2 ).* X2 );
case 10 , ret =
2* sqrt (14)*( -1+(3+( -3+ X1 ).* X1 ).* X1 +(12+( -24+12* X1 ).* X1 +( -30+30* X1 +20* X2 ).* X2 ).* X2 );
case 11 , ret =
sqrt ( 1 0 ) * ( 1 + ( - 2 4 + ( 1 2 6 + ( - 2 2 4 + 1 2 6 *X1 ).* X1 ).* X1 ).* X1 );
case 12 , ret =
sqrt ( 3 0 ) * ( 1 + ( - 2 2 + ( 1 0 5 + ( - 1 6 8 + 8 4 *X1 ).* X1 ).* X1 ).* X1 +( -2+(42+( -168+168* X1 ).* X1 ).* X1 ).* X2 );
case 13 , ret =
5* sqrt (2)*(1+( -18+(69+( -88+36* X1 ).* X1 ).* X1 ).* X1 +( -6+(102+( -312+216* X1 ).* X1 ).* X1 ...
+(6+( -96+216* X1 ).* X1 ).* X2 ).* X2 );
case 14 , ret = sqrt (70)*(1+( -12+(30+( -28+9* X1 ).* X1 ).* X1 ).* X1 +( -12+(132+( -228+108* X1 ).* X1 ).* X1 ...
+(30+( -300+270* X1 ).* X1 +( -20+180* X1 ).* X2 ).* X2 ).* X2 );
case 15 , ret = 3* sqrt (10)*(1+( -4+(6+( -4+ X1 ).* X1 ).* X1 ).* X1 +( -20+(60+( -60+20* X1 ).* X1 ).* X1 ...
+(90+( -180+90* X1 ).* X1 +( -140+140* X1 +70* X2 ).* X2 ).* X2 ).* X2 );
end % s w i t c h
end % f u n c t i o n
cDG = projectAlg2DG(g, cAlg, R, hatM) computes the representation matrix cDG of an algebraic function in the DG / modal basis by
performing the L2 -projection described in Sec. 3.4. The non-obvious input arguments are as follows: cAlg is of type function_handle taking
rows of coordinates and hatM the matrix M̂ computed by computeHatM.
function cDG = p r o j e c t A l g 2 D G(g , cAlg , ord , hatM )
global gPhi2D
ord = max ( ord ,1); [ Q1 , Q2 , W ] = q u a d R u l e 2 D( ord );
K = g . numT ; N = size ( hatM , 1);
F1 = @ ( X1 , X2 ) g . B (: ,1 ,1)* X1 + g . B (: ,1 ,2)* X2 + g . coordV0T (: ,1 ,1)* ones ( size ( X1 ));
F2 = @ ( X1 , X2 ) g . B (: ,2 ,1)* X1 + g . B (: ,2 ,2)* X2 + g . coordV0T (: ,1 ,2)* ones ( size ( X1 ));
rhs = cAlg ( F1 ( Q1 , Q2 ) , F2 ( Q1 , Q2 )) * ( repmat ( W . ' , 1 , N ) .* gPhi2D { ord });
cDG = rhs / hatM ;
end % f u n c t i o n
cLagr = projectDG2Lagrange(cDG) converts the representation matrix in the DG / modal basis to the respective representation matrix in
a Lagrange / nodal basis, both of dimension K × N. In this routine, the local basis functions ϕ̂i are sampled at the Lagrange nodes on the reference
triangle T̂ (cf. visualizeData ), whose x̂1 and x̂2 coordinates are stored in the variables L1 and L2. Note that this routine is defined only for
polynomial degrees p ∈ {0, 1, 2}.
function cLagr = p r o j e c t D G 2 L a g r a n g e( cDG )
[K , N ] = size ( cDG );
switch N
case 1 , L1 = 1/3; L2 = 1/3; % locally constant
case 3 , L1 = [0 , 1 , 0]; L2 = [0 , 0 , 1]; % locally linear
case 6 , L1 = [0 , 1 , 0 , 1/2 , 0 , 1/2]; L2 = [0 , 0 , 1 , 1/2 , 1/2 , 0]; % l o c a l l y q u a d r a t i c
end % s w i t c h
cLagr = zeros (K , N );
for i = 1 : N
cLagr = cLagr + cDG (: , i ) * phi (i , L1 , L2 );
end % for
end % f u n c t i o n
[Q, W] = quadRule1D(qOrd) returns a list of Gauss–Legendre quadrature points Q within the interval (0, 1) and associated (positive)
weights W. The quadrature rule is exact for polynomials of order qOrd. The length of the interval [0, 1] is incorporated in the weights. A rule with
n points is exact for polynomials up to order 2n − 1.
Example: f = @(s) s.^2; [Q, W] = quadRule1D(2); dot(W, f(Q))
function [Q , W ] = q u a d R u l e 1 D( qOrd )
switch qOrd
case {0 , 1} % R = 1 , n u m b e r of q u a d r a t u r e p o i n t s
Q = 0;
W = 2;
case {2 , 3} % R = 2
Q = sqrt (1/3)*[ -1 , 1];
W = [1 , 1];
case {4 , 5} % R = 3
Q = sqrt (3/5)*[ -1 , 0 , 1];
W = 1/9*[5 , 8 , 5];
case {6 , 7} % R = 4
Q = [ -1 , -1 ,1 ,1].* sqrt (3/7+[1 , -1 , -1 ,1]*2/7*sqrt (6/5));
W = 1/36*(18 + sqrt (30)*[ -1 ,1 ,1 , -1]);
32
case {8 , 9} % R = 5
Q = [ -1 , -1 ,0 ,1 ,1].* sqrt (5+[2 , -2 ,0 , -2 ,2]* sqrt ( 1 0 / 7 ) ) / 3 ;
W = 1 / 9 0 0 * ( 3 2 2 + 1 3 * sqrt (70)*[ -1 ,1 ,0 ,1 , -1]+[0 ,0 ,190 ,0 ,0]);
case {10 , 11} % R = 6
Q = [ 0 . 6 6 1 2 0 9 3 8 6 46 62 64 5, -0.6612093864662645 , -0.2386191860831969 , ...
0 . 2 3 8 6 1 9 1 8 6 08 31 96 9, -0.9324695142031521 , 0 . 9 3 2 4 6 9 5 1 4 2 0 3 1 5 2 1 ] ;
W = [ 0 . 3 6 0 7 6 1 5 7 3 04 81 38 6, 0 . 3 6 0 7 6 1 5 7 3 04 81 38 6, 0 . 4 6 7 9 1 3 9 3 4 57 26 91 0, ...
0 . 4 6 7 9 1 3 9 3 4 57 26 91 0, 0 . 1 7 1 3 2 4 4 9 2 37 91 70 4, 0 . 1 7 1 3 2 4 4 9 2 3 7 9 1 7 0 ] ;
case {12 , 13} % R = 7
Q = [ 0 . 0 0 0 0 0 0 0 0 0 00 00 00 0, 0 . 4 0 5 8 4 5 1 5 1 37 73 97 2, -0.4058451513773972 , ...
-0.7415311855993945 , 0 . 7 4 1 5 3 1 1 8 5 59 93 94 5, -0.9491079123427585 , ...
0.9491079123427585];
W = [ 0 . 4 1 7 9 5 9 1 8 3 67 34 69 4, 0 . 3 8 1 8 3 0 0 5 0 50 51 18 9, 0 . 3 8 1 8 3 0 0 5 0 50 51 18 9, ...
0 . 2 7 9 7 0 5 3 9 1 48 92 76 6, 0 . 2 7 9 7 0 5 3 9 1 48 92 76 6, 0 . 1 2 9 4 8 4 9 6 6 16 88 69 7, ...
0.1294849661688697];
case {14 , 15} % R = 8
Q = [ -0.1834346424956498 , 0 . 1 8 3 4 3 4 6 4 2 49 56 49 8, -0.5255324099163290 , ...
0 . 5 2 5 5 3 2 4 0 9 91 63 29 0, -0.7966664774136267 , 0 . 7 9 6 6 6 6 4 7 7 41 36 26 7, ...
-0.9602898564975363 , 0 . 9 6 0 2 8 9 8 5 6 4 9 7 5 3 6 3 ] ;
W = [ 0 . 3 6 2 6 8 3 7 8 3 37 83 62 0, 0 . 3 6 2 6 8 3 7 8 3 37 83 62 0, 0 . 3 1 3 7 0 6 6 4 5 87 78 87 3, ...
0 . 3 1 3 7 0 6 6 4 5 87 78 87 3, 0 . 2 2 2 3 8 1 0 3 4 45 33 74 5, 0 . 2 2 2 3 8 1 0 3 4 45 33 74 5, ...
0 . 1 0 1 2 2 8 5 3 6 29 03 76 3, 0 . 1 0 1 2 2 8 5 3 6 2 9 0 3 7 6 3 ] ;
case {16 , 17} % R = 9
Q = [ 0 . 0 0 0 0 0 0 0 0 0 00 00 00 0, -0.8360311073266358 , 0 . 8 3 6 0 3 1 1 0 7 32 66 35 8, ...
-0.9681602395076261 , 0 . 9 6 8 1 6 0 2 3 9 50 76 26 1, -0.3242534234038089 , ...
0 . 3 2 4 2 5 3 4 2 3 40 38 08 9, -0.6133714327005904 , 0 . 6 1 3 3 7 1 4 3 2 7 0 0 5 9 0 4 ] ;
W = [ 0 . 3 3 0 2 3 9 3 5 5 00 12 59 8, 0 . 1 8 0 6 4 8 1 6 0 69 48 57 4, 0 . 1 8 0 6 4 8 1 6 0 69 48 57 4, ...
0 . 0 8 1 2 7 4 3 8 8 36 15 74 4, 0 . 0 8 1 2 7 4 3 8 8 36 15 74 4, 0 . 3 1 2 3 4 7 0 7 7 04 00 02 9, ...
0 . 3 1 2 3 4 7 0 7 7 04 00 02 9, 0 . 2 6 0 6 1 0 6 9 6 40 29 35 4, 0 . 2 6 0 6 1 0 6 9 6 4 0 2 9 3 5 4 ] ;
end % s w i t c h
Q = ( Q + 1)/2; W = W /2; % t r a n s f o r m a t i o n [ -1; 1] -> [0 , 1]
end % f u n c t i o n
T
[Q1, Q2, W] = quadRule2D(qOrd) returns quadrature points q̂r = [q̂1r , q̂2r ] within the reference triangle T̂ in lists of x̂1 and x̂2 coordi-
nates Q1 and Q2, respectively, and the associated weights Q (cf. Sec. 3.3). The quadrature rule is exact for polynomials of order qOrd (cf. [24] for
orders 1, 2, 5, [25] for order 3, and [26] for order 4, 6). If qOrd is greater than 6, we call the third party function triquad [27] that uses Gaussian
quadrature points on a square which is collapsed to a triangle. The area 1/2 of the reference triangle T̂ is incorporated in the weights such that the
integral over one is 1/2.
Example: [Q1, Q2, W] = quadRule2D(2); N = 3; M = eye(N);
for i=1:N,for j=1:N,M(i,j)=sum(W.*phi(i,Q1,Q2).*phi(j,Q1,Q2));end,end
function [ Q1 , Q2 , W ] = q u a d R u l e 2 D( qOrd )
switch qOrd
case {0 , 1} % R = 1
Q1 = 1/3; Q2 = 1/3; W = 1/2;
case 2 % R = 3
Q1 = [1/6 , 2/3 , 1/6]; Q2 = [1/6 , 1/6 , 2/3]; W = [1/6 , 1/6 , 1/6];
case 3 % R = 4
Q1 = [0.666390246 , 0.178558728 , 0.280019915 , 0 . 0 7 5 0 3 1 1 0 9 ] ;
Q2 = [0.178558728 , 0.666390246 , 0.075031109 , 0 . 2 8 0 0 1 9 9 1 5 ] ;
W = [0.159020691 , 0.159020691 , 0.090979309 , 0 . 0 9 0 9 7 9 3 0 9 ] ;
case 4 % R = 6
Q1 = [ 0 . 4 4 5 9 4 8 4 9 09 15 965 , 0 . 1 0 8 1 0 3 0 1 81 68 07 0, 0 . 4 4 5 9 4 8 4 9 09 15 965, ...
0 . 0 9 1 5 7 6 2 1 35 09 77 1, 0 . 8 1 6 8 4 7 5 7 29 80 45 8, 0 . 0 9 1 5 7 6 2 1 3 5 0 9 7 7 1 ] ;
Q2 = [ 0 . 1 0 8 1 0 3 0 1 81 68 070 , 0 . 4 4 5 9 4 8 4 9 09 15 96 5, 0 . 4 4 5 9 4 8 4 9 09 15 965, ...
0 . 8 1 6 8 4 7 5 7 29 80 45 8, 0 . 0 9 1 5 7 6 2 1 35 09 77 1, 0 . 0 9 1 5 7 6 2 1 3 5 0 9 7 7 1 ] ;
W = [ 0 . 1 1 1 6 9 0 7 9 48 39 005 , 0 . 1 1 1 6 9 0 7 9 48 39 00 5, 0 . 1 1 1 6 9 0 7 9 48 39 005, ...
0 . 0 5 4 9 7 5 8 7 18 27 66 1, 0 . 0 5 4 9 7 5 8 7 18 27 66 1, 0 . 0 5 4 9 7 5 8 7 1 8 2 7 6 6 1 ] ;
case 5 % R = 7
a1 = (6 - sqrt ( 1 5 ) ) / 2 1 ; a2 = (6+ sqrt ( 1 5 ) ) / 2 1 ;
w1 = (155 - sqrt ( 1 5 ) ) / 2 4 0 0 ; w2 = (155+ sqrt ( 1 5 ) ) / 2 4 0 0 ;
Q1 = [1/3 , a1 , 1 -2* a1 , a1 , a2 , 1 -2* a2 , a2 ];
Q2 = [1/3 , 1 -2* a1 , a1 , a1 , 1 -2* a2 , a2 , a2 ];
W = [9/80 , w1 , w1 , w1 , w2 , w2 , w2 ];
case 6 % R = 12
Q1 = [ 0 . 0 6 3 0 8 9 0 1 44 91 502 , 0 . 8 7 3 8 2 1 9 7 10 16 99 6, 0 . 0 6 3 0 8 9 0 1 44 91 502, ...
0 . 2 4 9 2 8 6 7 4 51 70 91 0, 0 . 5 0 1 4 2 6 5 0 96 58 17 9, 0 . 2 4 9 2 8 6 7 4 51 70 910, ...
0 . 3 1 0 3 5 2 4 5 10 33 78 5, 0 . 0 5 3 1 4 5 0 4 98 44 81 6, 0 . 6 3 6 5 0 2 4 9 91 21 399, ...
0 . 0 5 3 1 4 5 0 4 98 44 81 6, 0 . 6 3 6 5 0 2 4 9 91 21 39 9, 0 . 3 1 0 3 5 2 4 5 1 0 3 3 7 8 5 ] ;
Q2 = [ 0 . 0 6 3 0 8 9 0 1 44 91 502 , 0 . 0 6 3 0 8 9 0 1 44 91 50 2, 0 . 8 7 3 8 2 1 9 7 10 16 996, ...
0 . 2 4 9 2 8 6 7 4 51 70 91 0, 0 . 2 4 9 2 8 6 7 4 51 70 91 0, 0 . 5 0 1 4 2 6 5 0 96 58 179, ...
0 . 0 5 3 1 4 5 0 4 98 44 81 6, 0 . 3 1 0 3 5 2 4 5 10 33 78 5, 0 . 0 5 3 1 4 5 0 4 98 44 816, ...
0 . 6 3 6 5 0 2 4 9 91 21 39 9, 0 . 3 1 0 3 5 2 4 5 10 33 78 5, 0 . 6 3 6 5 0 2 4 9 9 1 2 1 3 9 9 ] ;
W = [ 0 . 0 2 5 4 2 2 4 5 31 85 103 , 0 . 0 2 5 4 2 2 4 5 31 85 10 3, 0 . 0 2 5 4 2 2 4 5 31 85 103, ...
0 . 0 5 8 3 9 3 1 3 78 63 18 9, 0 . 0 5 8 3 9 3 1 3 78 63 18 9, 0 . 0 5 8 3 9 3 1 3 78 63 189, ...
33
0 . 0 4 1 4 2 5 5 3 78 09 18 7, 0 . 0 4 1 4 2 5 5 3 78 09 18 7, 0 . 0 4 1 4 2 5 5 3 78 09 187, ...
0 . 0 4 1 4 2 5 5 3 78 09 18 7, 0 . 0 4 1 4 2 5 5 3 78 09 18 7, 0 . 0 4 1 4 2 5 5 3 7 8 0 9 1 8 7 ] ;
o t h e r w i s e % use Gauss q u a d r a t u r e p o i n t s on s q u a r e
[X ,Y , Wx , Wy ] = triquad ( qOrd , [0 0; 0 1; 1 0]); % third party function , see r e f e r e n c e s
Q1 = X (:). ' ;
Q2 = Y (:). ' ;
W = Wx * Wy . ' ; W = W (:). ' ;
end % s w i t c h
end % f u n c t i o n
[XP1, XP2] = theta(nn, np, X1, X2) returns the mapped points from n− th edge to the n+ th edge of the reference triangle T̂ , cf. (26).
function [ XP1 , XP2 ] = theta ( nn , np , X1 , X2 )
switch nn
case 1
switch np
case 1 , XP1 = 1 - X1 ; XP2 = 1 - X2 ;
case 2 , XP1 = zeros ( size ( X1 )); XP2 = X2 ;
case 3 , XP1 = X1 ; XP2 = zeros ( size ( X1 ));
end % s w i t c h
case 2
switch np
case 1 , XP1 = 1 - X2 ; XP2 = X2 ;
case 2 , XP1 = zeros ( size ( X1 )); XP2 = 1 - X2 ;
case 3 , XP1 = X2 ; XP2 = zeros ( size ( X1 ));
end % s w i t c h
case 3
switch np
case 1 , XP1 = X1 ; XP2 = 1 - X1 ;
case 2 , XP1 = zeros ( size ( X1 )); XP2 = X1 ;
case 3 , XP1 = 1 - X1 ; XP2 = zeros ( size ( X1 ));
end % s w i t c h
end % s w i t c h
end % f u n c t i o n
visualizeGrid(g) visualizes the triangulation Th along with global and local indices and edge normals. The input argument g is the output
of the routine generateGridData .
function v i s u a l i z e G r i d( g )
figure ( ' Color ' , [1 , 1 , 1]); % white b a c k g r o u n d
hold ( ' on ' ) , axis ( ' off ' )
daspect ([1 , 1 , 1]) % a d j u s t a s p e c t ration , r e q u i r e s O c t a v e >= 3.8
t e x t a r r a y = @ ( x1 , x2 , s ) arrayfun ( @ (a ,b , c ) text (a ,b , int2str ( c )) , x1 , x2 , s );
%% triangle boundaries
trisurf ( g . V0T , g . coordV (: ,1) , g . coordV (: ,2) , zeros ( g . numV ,1) , ' f a c e c o l o r ' , ' none ' );
% % edge n o r m a l s ( s c a l e d by the f a c t o r of | E |/4 to fit the s c a l e s)
quiver ( g . baryE (: ,1) , g . baryE (: ,2) , ...
g . nuE (: ,1).*g . areaE (:)/4 , g . nuE (: , 2).* g . areaE (:)/4 , 0)
% % local edge n u m b e r s
w = [1/12 , 11/24 , 11/24; 11/24 , 1/12 , 11/24; 11/24 , 11/24 , 1/12];
for kE = 1 : 3
t e x t a r r a y( reshape ( g . coordV ( g . V0T ,1) , g . numT ,3)* w (: , kE ) , ...
reshape ( g . coordV ( g . V0T ,2) , g . numT ,3)* w (: , kE ) , kE * ones ( g . numT , 1))
end % for
%% global vertex numbers
t e x t a r r a y( g . coordV (: ,1) , g . coordV (: ,2) , (1: g . numV ) ' );
% % local v e r t e x n u m b e r s
w = [5/6 , 1/12 , 1/12; 1/12 , 5/6 , 1/12; 1/12 , 1/12 , 5/6];
for kV = 1 : 3
t e x t a r r a y( reshape ( g . coordV ( g . V0T ,1) , g . numT ,3)* w (: , kV ) , ...
reshape ( g . coordV ( g . V0T ,2) , g . numT ,3)* w (: , kV ) , kV * ones ( g . numT , 1))
end % for
% % g l o b a l edge n u m b e r s
t e x t a r r a y( g . baryE (: ,1) , g . baryE (: ,2) , (1: g . numE ) ' );
%% triangle numbers
t e x t a r r a y( g . baryT (: ,1) , g . baryT (: ,2) , (1: g . numT ) ' );
% % edge IDs
markEext = g . idE ~= 0; % mark b o u n d a r y edges
t e x t a r r a y( g . baryE ( markEext ,1) + g . nuE ( markEext ,1).* g . areaE ( markEext )/8 , ...
g . baryE ( markEext ,2) + g . nuE ( markEext ,2).* g . areaE ( markEext )/8 , g . idE ( markEext ))
end % f u n c t i o n
34
visualizeData(g, repLagr, varName, fileName, tLvl) writes a .vtu file for the visualization of a discrete quantity
in P p (Th ), p ∈ {0, 1, 2} defined on the triangulation g according to generateGridData . The name of the generated file is
fileName.tLvl.vtu, where tLvl stands for time level. The name of the quantity within the file is specified by varName. The argu-
ment repLagr should be a list of dimension K × N containing the Lagragian representation of the quantity. The kth row of repLagr has to hold
the value on T k for P0 (Th ), the values on the vertices of T k for P1 (Th ), and the values on the vertices and on the edge barycenters of T k for P2 (Th ).
Note that we treat functions of P0 (Th ) as if they were in P1 (Th ) as we assign the constant value on T k to the verices vk1 , vk2 , vk3 . An alternative
was using CellData instead of PointData . A usage example is found in Sec. 3.10.
function v i s u a l i z e D a t a(g , repLagr , varName , fileName , tLvl )
[K , N ] = size ( repLagr ) ;
% % open file
fileName = [ fileName , ' . ' , num2str ( tLvl ) , ' . vtu ' ];
file = fopen ( fileName , ' wt ' ) ; % if this file exists , then o v e r w r i t e
%% header
fprintf ( file , ' <? xml version ="1.0"? >\ n ' ) ;
fprintf ( file , ' < VTKFile type =" U n s t r u c t u r e d G r i d" version ="0.1" b y t e _ o r d e r=" L i t t l e E n d i a n" c o m p r e s s o r="
֒→v t k Z L i b D a t a C o m p r e s s o r" >\ n ' ) ;
fprintf ( file , ' < UnstructuredGrid >\ n ' ) ;
% % p o i n t s & cells
switch N
case {1 , 3}
P1 = reshape ( g . coordV0T (: , : , 1) ' , 3* K , 1) ;
P2 = reshape ( g . coordV0T (: , : , 2) ' , 3* K , 1) ;
numP = 3; % n u m b e r of local p o i n t s
id = 5; % vtk ID for l i n e a r p o l y n o m i a l s
case 6
P1 = reshape ([ g . coordV0T (: ,: ,1) , g . baryE0T (: ,[3 ,1 ,2] ,1) ] ' ,6*K ,1) ;
P2 = reshape ([ g . coordV0T (: ,: ,2) , g . baryE0T (: ,[3 ,1 ,2] ,2) ] ' ,6*K ,1) ;
numP = 6; % n u m b e r of local p o i n t s
id = 22; % vtk ID for q u a d r a t i c p o l y n o m i a l s
end % s w i t c h
fprintf ( file , ' < Piece N u m b e r O f P o i n t s ="% d " N u m b e r O f C e l l s ="% d " >\ n ' ,K * numP , K ) ;
fprintf ( file , ' < Points >\ n ' ) ;
fprintf ( file , ' < D a t a A r r a y type =" Float32 " N u m b e r O f C o m p o n e n t s ="3" format =" ascii " >\ n ' ) ;
fprintf ( file , ' %.3 e %.3 e %.3 e \ n ' , [ P1 , P2 , zeros ( numP *K , 1) ] ' ) ;
fprintf ( file , ' </ DataArray >\ n ' ) ;
fprintf ( file , ' </ Points >\ n ' ) ;
fprintf ( file , ' < Cells >\ n ' ) ;
fprintf ( file , ' < D a t a A r r a y type =" Int32 " Name =" c o n n e c t i v i t y" format =" ascii " >\ n ' ) ;
fprintf ( file , ' ' ) ; fprintf ( file , ' % d ' , 0: K * numP -1) ;
fprintf ( file , ' \ n </ DataArray >\ n ' ) ;
fprintf ( file , ' < D a t a A r r a y type =" Int32 " Name =" offsets " format =" ascii " >\ n ' ) ;
fprintf ( file , ' % d \ n ' , numP : numP : numP * K ) ;
fprintf ( file , ' </ DataArray >\ n ' ) ;
fprintf ( file , ' < D a t a A r r a y type =" UInt8 " Name =" types " format =" ascii " >\ n ' ) ;
fprintf ( file , ' % d \ n ' , id * ones (K , 1) ) ;
fprintf ( file , ' </ DataArray >\ n ' ) ;
fprintf ( file , ' </ Cells >\ n ' ) ;
% % data
switch N
case 1 % l o c a l l y c o n s t a n t
dataLagr = kron ( repLagr , [1;1;1]) ' ;
case 3 % l o c a l l y q u a d r a t i c
dataLagr = reshape ( repLagr ' , 1 , K * N ) ;
case 6 % l o c a l l y q u a d r a t i c ( p e r m u t a t i o n of local edge i n d i c e s due to vtk f o r m a t)
dataLagr = reshape ( repLagr (: , [1 ,2 ,3 ,6 ,4 ,5]) ' , 1 , K * N ) ;
end % s w i t c h
fprintf ( file , ' < P o i n t D a t a Scalars ="% s " >\ n ' , varName ) ;
fprintf ( file , ' < D a t a A r r a y type =" Float32 " Name ="% s " N u m b e r O f C o m p o n e n t s="1" format =" ascii " >\ n ' ,
֒→varName ) ;
fprintf ( file , ' %.3 e \ n ' , dataLagr ) ;
fprintf ( file , ' </ DataArray >\ n ' ) ;
fprintf ( file , ' </ PointData >\ n ' ) ;
%% footer
fprintf ( file , ' </ Piece >\ n ' ) ;
fprintf ( file , ' </ UnstructuredGrid >\ n ' ) ;
fprintf ( file , ' </ VTKFile >\ n ' ) ;
% % close file
fclose ( file ) ;
disp ([ ' Data written to ' fileName ])
end % f u n c t i o n
35
5. Conclusion and Outlook
The MATLAB / GNU Octave toolbox described in this work represents the first step of a multi-purpose package
that will include performance optimized discretizations for a range of standard problems, first, from the CFD, and then,
conceivably, from other application areas. Several important features will be added in the upcoming parts of this paper
and in the new releases of the toolbox: slope limiters based on the nodal slope limiting procedure proposed in [28, 29],
convection terms, nonlinear advection operators, and higher order time solvers. Furthermore, the object orientation
capabilites provided by MATLAB’s classdef, which are expected to be supported in future GNU Octave versions,
will be exploited to provide a more powerful and comfortable user interface.
Acknowledgments
The work of B. Reuter was supported by the German Research Foundation (DFG) under grant AI 117/1-1.
Index of notation
Symbol Definition
{| · |}, ~· average
" and
# jump on an edge
A
diag(A, B) ≔ , block-diagonal matrix with blocks A, B
B
#M cardinality of a set M
P
a·b ≔ 2m=1 am bm , Euclidean scalar product in R2
◦ composition of functions or Hadamard product
⊗ Kronecker product
c concentration (scalar-valued unknown)
d diffusion coefficient
η penalty parameter
δm∈M ≔ {1 if m ∈ M, 0 if m < M}, Kronecker delta
Ekn , Ên nth edge of the physical triangle T k , nth edge of the reference triangle T̂
V, E, T sets of vertices, edges, and triangles
ED , EN set of boundary edges, E∂Ω = ED ∪ EN
EΩ , E∂Ω set of interior edges, set of boundary edges
Fk affine mapping from T̂ to T k
h mesh fineness
hT ≔ diam(h), diameter of triangle T ∈ Th
J ≔ (0, tend ), open time interval
K ≔ #Th , number of triangles
νT unit normal on ∂T pointing outward of T
νk ≔ νT k
N ≔ (p + 1)(p + 2)/2, number of local degrees of freedom
ωr quadrature weight associated with q̂r
Ω, ∂Ω spatial domain in two dimensions, boundary of Ω
∂ΩD , ∂ΩN Dirichlet
√ and Neumann boundaries, ∂Ω = ∂ΩD ∪ ∂ΩN
p = ( 8N + 1 − 3)/2, polynomial degree
ϕik , ϕ̂i ith hierarchical basis function on T k , ith hierarchical basis function on T̂
P p (T ) space of polynomials of degree at most p
P p (Th ) ≔ {wh : Ω → R ; ∀T ∈ Th , wh |T ∈ P p (T )}
q̂r rth quadrature point in T̂
R number of quadrature points
R+ , R+0 set of (strictly) positive real numbers, set of nonnegative real numbers
t time variable
tn nth time level
tend end time
ϑ̂n− n+ mapping from Ên− to Ên+
∆tn ≔ tn+1 − tn , time step size
T k , ∂T k kth physical triangle, boundary of T k
T̂ bi-unit reference triangle
vki v̂i ith vertex of the physical triangle T k , ith vertex of the reference triangle T̂
T
x = [x1 , x2 ] , space variable in the physical domain Ω
T
x̂ = [ x̂1 , x̂2 ] , space variable in the reference triangle T̂
z diffusion mass flux (vector-valued unknown)
36
[1] H. Reed, T. R. Hill, Triangular mesh methods for the neutron transport equation, Tech. Rep. LA-UR-73-479, Los Alamos Scientific Labora-
tory, NM (1973).
[2] D. N. Arnold, F. Brezzi, B. Cockburn, L. D. Marini, Analysis of discontinuous Galerkin methods for elliptic problems, SIAM J. Num. Anal.
39 (2002) 1749–1779.
[3] The Mathworks, Inc., Natick, MA, USA, MATLAB and Statistics Toolbox Release 2014a.
[4] Octave community, GNU Octave 3.8.1 (2014).
URL http://www.gnu.org/software/octave/
[5] F. Frank, B. Reuter, V. Aizinger (2014). [link].
URL http://www.math.fau.de/FESTUNG
[6] B. Cockburn, C. Shu, The local discontinuous Galerkin method for time-dependent convection–diffusion systems, SIAM Journal on Numer-
ical Analysis 35 (6) (1998) 2440–2463. doi:10.1137/S0036142997316712.
[7] V. Aizinger, C. Dawson, B. Cockburn, P. Castillo, The local discontinuous Galerkin method for contaminant transport, Advances in Water
Resources 24 (1) (2000) 73–87. doi:10.1016/S0309-1708(00)00022-1.
[8] V. Aizinger, C. Dawson, A discontinuous Galerkin method for two-dimensional flow and transport in shallow water, Advances in Water
Resources 25 (1) (2002) 67–84. doi:10.1016/S0309-1708(01)00019-7.
[9] B. Rivière, Discontinuous Galerkin Methods for Solving Elliptic and Parabolic Equations: Theory and Implementation, Society for Industrial
and Applied Mathematics, 2008.
[10] M. G. Larson, F. Bengzon, The Finite Element Method: Theory, Implementation, and Applications, Vol. 10 of Texts in Computational Science
and Engineering, Springer, 2013.
[11] J. S. Hesthaven, T. Warburton, Nodal Discontinuous Galerkin Methods: Algorithms, Analysis, and Applications, Texts in Applied Mathe-
matics, Springer, 2008.
[12] T. Warburton (2014). [link].
URL http://www.nudg.org/
[13] Z. Fu, L. F. Gatica, F.-J. Sayas, Matlab tools for HDG in three dimensions, ArXiv e-prints arXiv:1305.1489.
[14] L. Chen, iFEM: an integrated finite element methods package in MATLAB, Tech. rep., Technical Report, University of California at Irvine
(2009).
[15] S. Funken, D. Praetorius, P. Wissgott, Efficient implementation of adaptive p1-fem in matlab, Comput. Methods Appl. Math. 11 (4) (2011)
460–490.
[16] T. Rahman, J. Valdman, Fast MATLAB assembly of FEM matrices in 2d and 3d: nodal elements, Appl. Math. Comput. 219 (13) (2013)
7151–7158. doi:10.1016/j.amc.2011.08.043.
[17] B. Cockburn, C.-W. Shu, TVB Runge–Kutta local projection discontinuous Galerkin finite element method for conservation laws ii: general
framework, Mathematics of Computation 52 (186) (1989) 411–435.
[18] C. Bahriawati, C. Carstensen, Three matlab implementations of the lowest-order Raviart–Thomas MFEM with a posteriori error control,
Computational Methods in Applied Mathematics 5 (4) (2005) 333–361.
[19] C. Geuzaine, J.-F. Remacle, Gmsh: A 3-D finite element mesh generator with built-in pre-and post-processing facilities, International Journal
for Numerical Methods in Engineering 79 (11) (2009) 1309–1331.
[20] R. Cools, An encyclopaedia of cubature formulas, Journal of Complexity 19 (3) (2003) 445–453. doi:10.1016/S0885-064X(03)00011-6.
[21] B. Cockburn, C.-W. Shu, The Runge–Kutta discontinuous Galerkin method for conservation laws V: multidimensional systems, J. Comput.
Phys. 141 (2) (1998) 199–224. doi:10.1006/jcph.1998.5892.
[22] L. S. Avila, Kitware, Inc, The VTK User’s Guide, Kitware, Incorporated, 2010.
[23] A. H. Squillacote, The ParaView Guide: A Parallel Visualization Application, Kitware, 2007.
[24] A. H. Stroud, Approximate Calculation of Multiple Integrals, Prentice-Hall Series in Automatic Computation, Prentice-Hall, 1971.
[25] P. Hillion, Numerical integration on a triangle, International Journal for Numerical Methods in Engineering 11 (5) (1977) 797–815.
doi:10.1002/nme.1620110504.
[26] G. R. Cowper, Gaussian quadrature formulas for triangles, International Journal for Numerical Methods in Engineering 7 (3) (1973) 405–408.
doi:10.1002/nme.1620070316.
[27] G. von Winckel, Gaussian Quadrature for Triangles. MATLAB Central File Exchange. Retrieved July 30, 2014. (2005). [link].
URL http://www.mathworks.com/matlabcentral/fileexchange/9230-gaussian-quadrature-for-triangles
[28] V. Aizinger, A geometry independent slope limiter for the discontinuous Galerkin method, in: E. Krause, Y. Shokin, M. Resch, D. Krner,
N. Shokina (Eds.), Computational Science and High Performance Computing IV, Vol. 115 of Notes on Numerical Fluid Mechanics and
Multidisciplinary Design, Springer Berlin Heidelberg, 2011, pp. 207–217. doi:10.1007/978-3-642-17770-5_16.
[29] D. Kuzmin, A vertex-based hierarchical slope limiter for adaptive discontinuous Galerkin methods, Journal of Computational
and Applied Mathematics 233 (12) (2010) 3077–3085, Finite Element Methods in Engineering and Science (FEMTEC 2009).
doi:10.1016/j.cam.2009.05.028.
37