This page was generated from unit-3.7-opsplit/opsplit.ipynb.

3.7 DG/HDG splitting

When solving unsteady problems with an operator-splitting method it might be beneficial to consider different space discretizations for different operators.

For a problem of the form

\[\partial_t u + A u + C u = 0\]

We consider the operator splitting:

  • 1.Step: Given \(u^0\), solve \(t^n \to t^{n+1}\): \(\quad \partial_t u + C u = 0 \Rightarrow u^{\frac12}\)

  • 2.Step: Given \(u^{\frac12}\), solve \(t^n \to t^{n+1}\): \(\quad \partial_t u + A u = 0 \Rightarrow u^{1}\)

This splitting is only first order accurate but allows to take different time discretizations for the substeps 1 and 2.

In this example we consider the Navier-Stokes problem again (cf. unit 3.2) and want to combine

  • an \(H(div)\) -conforming Hybrid DG method (which is a very good discretization for Stokes-type problems) and

  • a standard upwind DG method for the discretization of the convection

The weak form is: Find \((\mathbf{u},p):[0,T] \to (H_{0,D}^1)^d \times L^2\), s.t.

\begin{align} \int_{\Omega} \partial_t \mathbf{u} \cdot v + \int_{\Omega} \nu \nabla \mathbf{u} \nabla \mathbf{v} + \mathbf{u} \cdot \nabla \mathbf{u} \ \mathbf{v} - \int_{\Omega} \operatorname{div}(\mathbf{v}) p &= \int f \mathbf{v} && \forall \mathbf{v} \in (H_{0,D}^1)^d, \\ - \int_{\Omega} \operatorname{div}(\mathbf{u}) q &= 0 && \forall q \in L^2, \\ \quad \mathbf{u}(t=0) & = \mathbf{u}_0 \end{align}

Again, we consider the benchmark setup from http://www.featflow.de/en/benchmarks/cfdbenchmarking/flow/dfg_benchmark2_re100.html . The geometry:

image0

The viscosity is set to \(\nu = 10^{-3}\).

[1]:
order = 3
# import libraries, set up geometry and generate mesh
from ngsolve import *
from ngsolve.webgui import Draw
from netgen.occ import *
shape = Rectangle(2,0.41).Circle(0.2,0.2,0.05).Reverse().Face()
shape.edges.name="cyl"
shape.edges.Min(X).name="inlet"
shape.edges.Max(X).name="outlet"
shape.edges.Min(Y).name="wall"
shape.edges.Max(Y).name="wall"
mesh = Mesh(OCCGeometry(shape, dim=2).GenerateMesh(maxh=0.07)).Curve(order)

For the HDG formulation we use the product space with

  • \(BDM_k\): \(H(div)\) conforming FE space (degree k)

  • Vector facet space: facet functions of degree k (vector valued and only in tangential direction)

  • piecewise polynomials up to degree \(k-1\) for the pressure

image0

HDG spaces

[2]:
V1 = HDiv ( mesh, order = order, dirichlet = "wall|cyl|inlet" )
V2 = TangentialFacetFESpace(mesh, order = order, dirichlet = "wall|cyl|inlet" )
Q = L2( mesh, order = order-1)
V = FESpace ([V1,V2,Q])

u, uhat, p = V.TrialFunction()
v, vhat, q = V.TestFunction()

Parameter and directions (normal / tangential projection):

[3]:
nu = 0.001 # viscosity
n = specialcf.normal(mesh.dim)
h = specialcf.mesh_size
def tang(vec):
    return vec - (vec*n)*n

Stokes discretization

The bilinear form to the HDG (see unit-2.8-DG for scalar HDG) discretized Stokes problem is:

\begin{align*} K_h((\mathbf{u}_T,\mathbf{u}_F,p),(\mathbf{v}_T,\mathbf{v}_F,q) := & \displaystyle \sum_{T\in\mathcal{T}} \left\{ \int_{T} \nu {\nabla} {\mathbf{u}_T} \! : \! {\nabla} {\mathbf{v}_T} \ d {x} - \int_{\Omega} \operatorname{div}(\mathbf{u}) q \ d {x} - \int_{\Omega} \operatorname{div}(\mathbf{v}) p \ d {x} \right. \\ & \qquad \left. - \int_{\partial T} \nu \frac{\partial {\mathbf{u}_T}}{\partial {n} } [ \mathbf{v} ]_t \ d {s} - \int_{\partial T} \nu \frac{\partial {\mathbf{v}_T}}{\partial {n} } [ \mathbf{u} ]_t \ d {s} + \int_{\partial T} \nu \frac{\alpha k^2}{h} [\mathbf{u}]_t [\mathbf{v}]_t \ d {s} \right\} \end{align*}

where \([\cdot]_t\) is the tangential projection of the jump \((\cdot)_T - (\cdot)_F\).

[4]:
alpha = 4  # SIP parameter
dS = dx(element_boundary=True)
a = BilinearForm ( V, symmetric=True)
a += nu * InnerProduct ( Grad(u), Grad(v) ) * dx
a += nu * InnerProduct ( Grad(u)*n, tang(vhat-v) ) * dS
a += nu * InnerProduct ( Grad(v)*n, tang(uhat-u) ) * dS
a += nu * alpha*order*order/h * InnerProduct(tang(vhat-v), tang(uhat-u)) * dS
a += (-div(u)*q - div(v)*p) * dx
a.Assemble()
[4]:
<ngsolve.comp.BilinearForm at 0x7fcd683782f0>

The mass matrix is simply

\begin{align*} M_h((\mathbf{u}_T,\mathbf{u}_F,p),(\mathbf{v}_T,\mathbf{v}_F,q) := \int_{\Omega} \mathbf{u}_T \cdot \mathbf{v}_T dx \end{align*}

[5]:
m = BilinearForm(u * v * dx, symmetric=True).Assemble()

Initial conditions

As initial condition we solve the Stokes problem:

[6]:
gfu = GridFunction(V)

U0 = 1.5
uin = CF( (U0*4*y*(0.41-y)/(0.41*0.41),0) )
gfu.components[0].Set(uin, definedon=mesh.Boundaries("inlet"))

invstokes = a.mat.Inverse(V.FreeDofs(), inverse="umfpack")
gfu.vec.data += invstokes @ -a.mat * gfu.vec
[7]:
scenes = [ \
    Draw( gfu.components[0], mesh, "velocity"), \
    Draw( gfu.components[2], mesh, "pressure")];

Application of convection

In the operator splitted approach we want to apply only operator applications for the convection part. Further, we want to do this in a usual DG setting. As a model problem we use the following procedure:

  • Given \((\mathbf{u},p)\) in HDG space: project into \(\hat{\mathbf{u}}\) in usual DG space

  • Solve \(\partial_t \hat{\mathbf{u}} = C \hat{\mathbf{u}}\) by explicit scheme (involves convection evaluations and mass matrix operations only)

  • Solve Unsteady Stokes step with r.h.s. from convection sub-problem. To this end evaluate mixed mass matrix \(\int \hat{\mathbf{u}} \cdot \mathbf{v}\) to obtain a functional on the HDG space

For the projection steps we use mixed mass matrices:

  • \(M_m\) : \(HDG \times DG \to \mathbb{R}\)

  • \(M_m^T\) : \(DG \times HDG \to \mathbb{R}\)

  • \(M_{DG}\) : \(DG \times DG \to \mathbb{R}\) (block diagonal)

To set up mixed mass matrices we use a bilinear form with two different FESpaces.

We do not assemble the matrices as we will only need the matrix-vector applications of \(M_m\), \(M_m^T\) and \(M_{DG}^{-1}\).

There is a more elegant version (ConvertL2Operator, similar to TraceOperator) that we discuss below.

mixed mass matrices

[8]:
VL2 = VectorL2(mesh, order=order, piola=True)
uL2, vL2 = VL2.TnT()
bfmixed = BilinearForm ( vL2*u*dx, nonassemble=True )
bfmixedT = BilinearForm ( uL2*v*dx, nonassemble=True)

convection operator

  • convection operation with standard Upwinding

  • No set up of the matrix (only interested in operator applications)

  • for the advection velocity we use the \(H(div)\)-conforming velocity (which is pointwise divergence free).

[9]:
vel = gfu.components[0]
convL2 = BilinearForm( (-InnerProduct(Grad(vL2) * vel, uL2)) * dx, nonassemble=True )
un = InnerProduct(vel,n)
upwindL2 = IfPos(un, un*uL2, un*uL2.Other(bnd=uin))

dskel_inner  = dx(skeleton=True)
dskel_bound  = ds(skeleton=True)

convL2 += InnerProduct (upwindL2, vL2-vL2.Other()) * dskel_inner
convL2 += InnerProduct (upwindL2, vL2) * dskel_bound

Solution of convection steps as an operator (BaseMatrix)

We now define the solution of the convection problem for an initial data \(u\) in the HDG space. The return value (y) is \(M_m \hat{u}\) where \(\hat{u}\) is the solution of several explicit Euler steps of the convection problem

\[\partial_t \hat{u} = - M_{DG}^{-1} C \hat{u}\]
[10]:
class PropagateConvection(BaseMatrix):
    def __init__(self,tau,steps):
        super(PropagateConvection, self).__init__()
        self.tau = tau; self.steps = steps
        self.h = V.ndof; self.w = V.ndof # operator domain and range
        self.mL2 = VL2.Mass(Id(mesh.dim)); self.invmL2 = self.mL2.Inverse()
        self.vecL2 = bfmixed.mat.CreateColVector() # temp vector
    def Mult(self, x, y):
        self.vecL2.data = self.invmL2 @ bfmixed.mat * x # <- project from Hdiv to L2
        for i in range(self.steps):
            self.vecL2.data -= self.tau/self.steps * self.invmL2 @ convL2.mat * self.vecL2
        y.data = bfmixedT.mat * self.vecL2
    def CreateColVector(self):
        return CreateVVector(self.h)

The operator splitted method (Yanenko splitting)

Setup of \(M^*\):

[11]:
t = 0; tend = 0
tau = 0.01; substeps = 10

mstar = m.mat.CreateMatrix()
mstar.AsVector().data = m.mat.AsVector() + tau * a.mat.AsVector()
inv = mstar.Inverse(V.FreeDofs(), inverse="umfpack")

The time loop:

[12]:
for s in scenes : s.Draw()
[13]:
%%time
tend += 1
res = gfu.vec.CreateVector()
convpropagator = PropagateConvection(tau,substeps)
while t < tend:
    gfu.vec.data += inv @ (convpropagator - mstar) * gfu.vec
    t += tau
    print ("\r  t =", t, end="")
    for s in scenes : s.Redraw()
more-to-come:

t = 0.01 t = 0.02

</pre>

t = 0.01 t = 0.02

end{sphinxVerbatim}

t = 0.01 t = 0.02

more-to-come:

t = 0.03 t = 0.04

</pre>

t = 0.03 t = 0.04

end{sphinxVerbatim}

t = 0.03 t = 0.04

more-to-come:

t = 0.05 t = 0.060000000000000005

</pre>

t = 0.05 t = 0.060000000000000005

end{sphinxVerbatim}

t = 0.05 t = 0.060000000000000005

more-to-come:

t = 0.07 t = 0.08

</pre>

t = 0.07 t = 0.08

end{sphinxVerbatim}

t = 0.07 t = 0.08

more-to-come:

t = 0.09 t = 0.09999999999999999

</pre>

t = 0.09 t = 0.09999999999999999

end{sphinxVerbatim}

t = 0.09 t = 0.09999999999999999

more-to-come:

t = 0.10999999999999999 t = 0.11999999999999998

</pre>

t = 0.10999999999999999 t = 0.11999999999999998

end{sphinxVerbatim}

t = 0.10999999999999999 t = 0.11999999999999998

more-to-come:

t = 0.12999999999999998 t = 0.13999999999999999

</pre>

t = 0.12999999999999998 t = 0.13999999999999999

end{sphinxVerbatim}

t = 0.12999999999999998 t = 0.13999999999999999

more-to-come:

t = 0.15 t = 0.16

</pre>

t = 0.15 t = 0.16

end{sphinxVerbatim}

t = 0.15 t = 0.16

more-to-come:

t = 0.17 t = 0.18000000000000002

</pre>

t = 0.17 t = 0.18000000000000002

end{sphinxVerbatim}

t = 0.17 t = 0.18000000000000002

more-to-come:

t = 0.19000000000000003 t = 0.20000000000000004

</pre>

t = 0.19000000000000003 t = 0.20000000000000004

end{sphinxVerbatim}

t = 0.19000000000000003 t = 0.20000000000000004

more-to-come:

t = 0.21000000000000005 t = 0.22000000000000006

</pre>

t = 0.21000000000000005 t = 0.22000000000000006

end{sphinxVerbatim}

t = 0.21000000000000005 t = 0.22000000000000006

more-to-come:

t = 0.23000000000000007 t = 0.24000000000000007

</pre>

t = 0.23000000000000007 t = 0.24000000000000007

end{sphinxVerbatim}

t = 0.23000000000000007 t = 0.24000000000000007

more-to-come:

t = 0.25000000000000006 t = 0.26000000000000006

</pre>

t = 0.25000000000000006 t = 0.26000000000000006

end{sphinxVerbatim}

t = 0.25000000000000006 t = 0.26000000000000006

more-to-come:

t = 0.2700000000000001 t = 0.2800000000000001

</pre>

t = 0.2700000000000001 t = 0.2800000000000001

end{sphinxVerbatim}

t = 0.2700000000000001 t = 0.2800000000000001

more-to-come:

t = 0.2900000000000001 t = 0.3000000000000001

</pre>

t = 0.2900000000000001 t = 0.3000000000000001

end{sphinxVerbatim}

t = 0.2900000000000001 t = 0.3000000000000001

more-to-come:

t = 0.3100000000000001 t = 0.3200000000000001

</pre>

t = 0.3100000000000001 t = 0.3200000000000001

end{sphinxVerbatim}

t = 0.3100000000000001 t = 0.3200000000000001

more-to-come:

t = 0.3300000000000001 t = 0.34000000000000014

</pre>

t = 0.3300000000000001 t = 0.34000000000000014

end{sphinxVerbatim}

t = 0.3300000000000001 t = 0.34000000000000014

more-to-come:

t = 0.35000000000000014 t = 0.36000000000000015

</pre>

t = 0.35000000000000014 t = 0.36000000000000015

end{sphinxVerbatim}

t = 0.35000000000000014 t = 0.36000000000000015

more-to-come:

t = 0.37000000000000016 t = 0.38000000000000017

</pre>

t = 0.37000000000000016 t = 0.38000000000000017

end{sphinxVerbatim}

t = 0.37000000000000016 t = 0.38000000000000017

more-to-come:

t = 0.3900000000000002 t = 0.4000000000000002

</pre>

t = 0.3900000000000002 t = 0.4000000000000002

end{sphinxVerbatim}

t = 0.3900000000000002 t = 0.4000000000000002

more-to-come:

t = 0.4100000000000002 t = 0.4200000000000002

</pre>

t = 0.4100000000000002 t = 0.4200000000000002

end{sphinxVerbatim}

t = 0.4100000000000002 t = 0.4200000000000002

more-to-come:

t = 0.4300000000000002 t = 0.4400000000000002

</pre>

t = 0.4300000000000002 t = 0.4400000000000002

end{sphinxVerbatim}

t = 0.4300000000000002 t = 0.4400000000000002

more-to-come:

t = 0.45000000000000023 t = 0.46000000000000024

</pre>

t = 0.45000000000000023 t = 0.46000000000000024

end{sphinxVerbatim}

t = 0.45000000000000023 t = 0.46000000000000024

more-to-come:

t = 0.47000000000000025 t = 0.48000000000000026

</pre>

t = 0.47000000000000025 t = 0.48000000000000026

end{sphinxVerbatim}

t = 0.47000000000000025 t = 0.48000000000000026

more-to-come:

t = 0.49000000000000027 t = 0.5000000000000002

</pre>

t = 0.49000000000000027 t = 0.5000000000000002

end{sphinxVerbatim}

t = 0.49000000000000027 t = 0.5000000000000002

more-to-come:

t = 0.5100000000000002 t = 0.5200000000000002

</pre>

t = 0.5100000000000002 t = 0.5200000000000002

end{sphinxVerbatim}

t = 0.5100000000000002 t = 0.5200000000000002

more-to-come:

t = 0.5300000000000002 t = 0.5400000000000003

</pre>

t = 0.5300000000000002 t = 0.5400000000000003

end{sphinxVerbatim}

t = 0.5300000000000002 t = 0.5400000000000003

more-to-come:

t = 0.5500000000000003 t = 0.5600000000000003

</pre>

t = 0.5500000000000003 t = 0.5600000000000003

end{sphinxVerbatim}

t = 0.5500000000000003 t = 0.5600000000000003

more-to-come:

t = 0.5700000000000003 t = 0.5800000000000003

</pre>

t = 0.5700000000000003 t = 0.5800000000000003

end{sphinxVerbatim}

t = 0.5700000000000003 t = 0.5800000000000003

more-to-come:

t = 0.5900000000000003 t = 0.6000000000000003

</pre>

t = 0.5900000000000003 t = 0.6000000000000003

end{sphinxVerbatim}

t = 0.5900000000000003 t = 0.6000000000000003

more-to-come:

t = 0.6100000000000003 t = 0.6200000000000003

</pre>

t = 0.6100000000000003 t = 0.6200000000000003

end{sphinxVerbatim}

t = 0.6100000000000003 t = 0.6200000000000003

more-to-come:

t = 0.6300000000000003 t = 0.6400000000000003

</pre>

t = 0.6300000000000003 t = 0.6400000000000003

end{sphinxVerbatim}

t = 0.6300000000000003 t = 0.6400000000000003

more-to-come:

t = 0.6500000000000004 t = 0.6600000000000004

</pre>

t = 0.6500000000000004 t = 0.6600000000000004

end{sphinxVerbatim}

t = 0.6500000000000004 t = 0.6600000000000004

more-to-come:

t = 0.6700000000000004 t = 0.6800000000000004

</pre>

t = 0.6700000000000004 t = 0.6800000000000004

end{sphinxVerbatim}

t = 0.6700000000000004 t = 0.6800000000000004

more-to-come:

t = 0.6900000000000004 t = 0.7000000000000004

</pre>

t = 0.6900000000000004 t = 0.7000000000000004

end{sphinxVerbatim}

t = 0.6900000000000004 t = 0.7000000000000004

more-to-come:

t = 0.7100000000000004 t = 0.7200000000000004

</pre>

t = 0.7100000000000004 t = 0.7200000000000004

end{sphinxVerbatim}

t = 0.7100000000000004 t = 0.7200000000000004

more-to-come:

t = 0.7300000000000004 t = 0.7400000000000004

</pre>

t = 0.7300000000000004 t = 0.7400000000000004

end{sphinxVerbatim}

t = 0.7300000000000004 t = 0.7400000000000004

more-to-come:

t = 0.7500000000000004 t = 0.7600000000000005

</pre>

t = 0.7500000000000004 t = 0.7600000000000005

end{sphinxVerbatim}

t = 0.7500000000000004 t = 0.7600000000000005

more-to-come:

t = 0.7700000000000005 t = 0.7800000000000005

</pre>

t = 0.7700000000000005 t = 0.7800000000000005

end{sphinxVerbatim}

t = 0.7700000000000005 t = 0.7800000000000005

more-to-come:

t = 0.7900000000000005 t = 0.8000000000000005

</pre>

t = 0.7900000000000005 t = 0.8000000000000005

end{sphinxVerbatim}

t = 0.7900000000000005 t = 0.8000000000000005

more-to-come:

t = 0.8100000000000005 t = 0.8200000000000005

</pre>

t = 0.8100000000000005 t = 0.8200000000000005

end{sphinxVerbatim}

t = 0.8100000000000005 t = 0.8200000000000005

more-to-come:

t = 0.8300000000000005 t = 0.8400000000000005

</pre>

t = 0.8300000000000005 t = 0.8400000000000005

end{sphinxVerbatim}

t = 0.8300000000000005 t = 0.8400000000000005

more-to-come:

t = 0.8500000000000005 t = 0.8600000000000005

</pre>

t = 0.8500000000000005 t = 0.8600000000000005

end{sphinxVerbatim}

t = 0.8500000000000005 t = 0.8600000000000005

more-to-come:

t = 0.8700000000000006 t = 0.8800000000000006

</pre>

t = 0.8700000000000006 t = 0.8800000000000006

end{sphinxVerbatim}

t = 0.8700000000000006 t = 0.8800000000000006

more-to-come:

t = 0.8900000000000006 t = 0.9000000000000006

</pre>

t = 0.8900000000000006 t = 0.9000000000000006

end{sphinxVerbatim}

t = 0.8900000000000006 t = 0.9000000000000006

more-to-come:

t = 0.9100000000000006 t = 0.9200000000000006 t = 0.9300000000000006

</pre>

t = 0.9100000000000006 t = 0.9200000000000006 t = 0.9300000000000006

end{sphinxVerbatim}

t = 0.9100000000000006 t = 0.9200000000000006 t = 0.9300000000000006

more-to-come:

t = 0.9400000000000006 t = 0.9500000000000006

</pre>

t = 0.9400000000000006 t = 0.9500000000000006

end{sphinxVerbatim}

t = 0.9400000000000006 t = 0.9500000000000006

more-to-come:

t = 0.9600000000000006 t = 0.9700000000000006

</pre>

t = 0.9600000000000006 t = 0.9700000000000006

end{sphinxVerbatim}

t = 0.9600000000000006 t = 0.9700000000000006

more-to-come:

t = 0.9800000000000006 t = 0.9900000000000007

</pre>

t = 0.9800000000000006 t = 0.9900000000000007

end{sphinxVerbatim}

t = 0.9800000000000006 t = 0.9900000000000007

t = 1.0000000000000007CPU times: user 11.6 s, sys: 47.1 ms, total: 11.6 s

Wall time: 11.6 s </pre>

t = 1.0000000000000007CPU times: user 11.6 s, sys: 47.1 ms, total: 11.6 s

Wall time: 11.6 s end{sphinxVerbatim}

t = 1.0000000000000007CPU times: user 11.6 s, sys: 47.1 ms, total: 11.6 s

Wall time: 11.6 s

ConvertL2Operator

Instead of doing the conversion between the spaces "by hand", we can use the convenience operator ConvertL2Operator that realizes $ M_{DG}^{-1} :nbsphinx-math:`cdot `M_m $:

[14]:
HDG_to_L2 = V1.ConvertL2Operator(VL2) @ Embedding(V.ndof, V.Range(0)).T

The transpose of \(M_{DG}^{-1} \cdot M_m\) is \(M_m^T \cdot M_{DG}^{-1}\) * V1.ConvertL2Operator(VL2): V1 \(\to\) VL2 * V1.ConvertL2Operator(VL2).T: VL2' \(\to\) V1'

Overwrite the application and use the ConvertL2Operator instead ot the mixed mass matrices:

[15]:
def NewMult(self, x, y):
    self.vecL2.data = HDG_to_L2 * x # <- project from Hdiv to L2
    for i in range(self.steps):
        self.vecL2.data -= tau/self.steps*self.invmL2@convL2.mat*self.vecL2
    y.data = HDG_to_L2.T @ self.mL2 * self.vecL2

PropagateConvection.Mult = NewMult

The time loop:

[16]:
for s in scenes : s.Draw()
[17]:
%%time
tend += 1
convpropagator = PropagateConvection(tau,substeps)
while t < tend:
    gfu.vec.data += inv @ (convpropagator - mstar) * gfu.vec
    t += tau
    print ("\r  t =", t, end="")
    for s in scenes : s.Redraw()
more-to-come:

t = 1.0100000000000007 t = 1.0200000000000007

</pre>

t = 1.0100000000000007 t = 1.0200000000000007

end{sphinxVerbatim}

t = 1.0100000000000007 t = 1.0200000000000007

more-to-come:

t = 1.0300000000000007 t = 1.0400000000000007 t = 1.0500000000000007

</pre>

t = 1.0300000000000007 t = 1.0400000000000007 t = 1.0500000000000007

end{sphinxVerbatim}

t = 1.0300000000000007 t = 1.0400000000000007 t = 1.0500000000000007

more-to-come:

t = 1.0600000000000007 t = 1.0700000000000007

</pre>

t = 1.0600000000000007 t = 1.0700000000000007

end{sphinxVerbatim}

t = 1.0600000000000007 t = 1.0700000000000007

more-to-come:

t = 1.0800000000000007 t = 1.0900000000000007 t = 1.1000000000000008

</pre>

t = 1.0800000000000007 t = 1.0900000000000007 t = 1.1000000000000008

end{sphinxVerbatim}

t = 1.0800000000000007 t = 1.0900000000000007 t = 1.1000000000000008

more-to-come:

t = 1.1100000000000008 t = 1.1200000000000008

</pre>

t = 1.1100000000000008 t = 1.1200000000000008

end{sphinxVerbatim}

t = 1.1100000000000008 t = 1.1200000000000008

more-to-come:

t = 1.1300000000000008 t = 1.1400000000000008 t = 1.1500000000000008

</pre>

t = 1.1300000000000008 t = 1.1400000000000008 t = 1.1500000000000008

end{sphinxVerbatim}

t = 1.1300000000000008 t = 1.1400000000000008 t = 1.1500000000000008

more-to-come:

t = 1.1600000000000008 t = 1.1700000000000008

</pre>

t = 1.1600000000000008 t = 1.1700000000000008

end{sphinxVerbatim}

t = 1.1600000000000008 t = 1.1700000000000008

more-to-come:

t = 1.1800000000000008 t = 1.1900000000000008

</pre>

t = 1.1800000000000008 t = 1.1900000000000008

end{sphinxVerbatim}

t = 1.1800000000000008 t = 1.1900000000000008

more-to-come:

t = 1.2000000000000008 t = 1.2100000000000009 t = 1.2200000000000009

</pre>

t = 1.2000000000000008 t = 1.2100000000000009 t = 1.2200000000000009

end{sphinxVerbatim}

t = 1.2000000000000008 t = 1.2100000000000009 t = 1.2200000000000009

more-to-come:

t = 1.2300000000000009 t = 1.2400000000000009

</pre>

t = 1.2300000000000009 t = 1.2400000000000009

end{sphinxVerbatim}

t = 1.2300000000000009 t = 1.2400000000000009

more-to-come:

t = 1.2500000000000009 t = 1.260000000000001 t = 1.270000000000001

</pre>

t = 1.2500000000000009 t = 1.260000000000001 t = 1.270000000000001

end{sphinxVerbatim}

t = 1.2500000000000009 t = 1.260000000000001 t = 1.270000000000001

more-to-come:

t = 1.280000000000001 t = 1.290000000000001

</pre>

t = 1.280000000000001 t = 1.290000000000001

end{sphinxVerbatim}

t = 1.280000000000001 t = 1.290000000000001

more-to-come:

t = 1.300000000000001 t = 1.310000000000001

</pre>

t = 1.300000000000001 t = 1.310000000000001

end{sphinxVerbatim}

t = 1.300000000000001 t = 1.310000000000001

more-to-come:

t = 1.320000000000001 t = 1.330000000000001

</pre>

t = 1.320000000000001 t = 1.330000000000001

end{sphinxVerbatim}

t = 1.320000000000001 t = 1.330000000000001

more-to-come:

t = 1.340000000000001 t = 1.350000000000001 t = 1.360000000000001

</pre>

t = 1.340000000000001 t = 1.350000000000001 t = 1.360000000000001

end{sphinxVerbatim}

t = 1.340000000000001 t = 1.350000000000001 t = 1.360000000000001

more-to-come:

t = 1.370000000000001 t = 1.380000000000001

</pre>

t = 1.370000000000001 t = 1.380000000000001

end{sphinxVerbatim}

t = 1.370000000000001 t = 1.380000000000001

more-to-come:

t = 1.390000000000001 t = 1.400000000000001 t = 1.410000000000001

</pre>

t = 1.390000000000001 t = 1.400000000000001 t = 1.410000000000001

end{sphinxVerbatim}

t = 1.390000000000001 t = 1.400000000000001 t = 1.410000000000001

more-to-come:

t = 1.420000000000001 t = 1.430000000000001

</pre>

t = 1.420000000000001 t = 1.430000000000001

end{sphinxVerbatim}

t = 1.420000000000001 t = 1.430000000000001

more-to-come:

t = 1.440000000000001 t = 1.450000000000001 t = 1.460000000000001

</pre>

t = 1.440000000000001 t = 1.450000000000001 t = 1.460000000000001

end{sphinxVerbatim}

t = 1.440000000000001 t = 1.450000000000001 t = 1.460000000000001

more-to-come:

t = 1.470000000000001 t = 1.480000000000001

</pre>

t = 1.470000000000001 t = 1.480000000000001

end{sphinxVerbatim}

t = 1.470000000000001 t = 1.480000000000001

more-to-come:

t = 1.490000000000001 t = 1.500000000000001

</pre>

t = 1.490000000000001 t = 1.500000000000001

end{sphinxVerbatim}

t = 1.490000000000001 t = 1.500000000000001

more-to-come:

t = 1.5100000000000011 t = 1.5200000000000011 t = 1.5300000000000011

</pre>

t = 1.5100000000000011 t = 1.5200000000000011 t = 1.5300000000000011

end{sphinxVerbatim}

t = 1.5100000000000011 t = 1.5200000000000011 t = 1.5300000000000011

more-to-come:

t = 1.5400000000000011 t = 1.5500000000000012

</pre>

t = 1.5400000000000011 t = 1.5500000000000012

end{sphinxVerbatim}

t = 1.5400000000000011 t = 1.5500000000000012

more-to-come:

t = 1.5600000000000012 t = 1.5700000000000012 t = 1.5800000000000012

</pre>

t = 1.5600000000000012 t = 1.5700000000000012 t = 1.5800000000000012

end{sphinxVerbatim}

t = 1.5600000000000012 t = 1.5700000000000012 t = 1.5800000000000012

more-to-come:

t = 1.5900000000000012 t = 1.6000000000000012

</pre>

t = 1.5900000000000012 t = 1.6000000000000012

end{sphinxVerbatim}

t = 1.5900000000000012 t = 1.6000000000000012

more-to-come:

t = 1.6100000000000012 t = 1.6200000000000012 t = 1.6300000000000012

</pre>

t = 1.6100000000000012 t = 1.6200000000000012 t = 1.6300000000000012

end{sphinxVerbatim}

t = 1.6100000000000012 t = 1.6200000000000012 t = 1.6300000000000012

more-to-come:

t = 1.6400000000000012 t = 1.6500000000000012

</pre>

t = 1.6400000000000012 t = 1.6500000000000012

end{sphinxVerbatim}

t = 1.6400000000000012 t = 1.6500000000000012

more-to-come:

t = 1.6600000000000013 t = 1.6700000000000013 t = 1.6800000000000013

</pre>

t = 1.6600000000000013 t = 1.6700000000000013 t = 1.6800000000000013

end{sphinxVerbatim}

t = 1.6600000000000013 t = 1.6700000000000013 t = 1.6800000000000013

more-to-come:

t = 1.6900000000000013 t = 1.7000000000000013

</pre>

t = 1.6900000000000013 t = 1.7000000000000013

end{sphinxVerbatim}

t = 1.6900000000000013 t = 1.7000000000000013

more-to-come:

t = 1.7100000000000013 t = 1.7200000000000013

</pre>

t = 1.7100000000000013 t = 1.7200000000000013

end{sphinxVerbatim}

t = 1.7100000000000013 t = 1.7200000000000013

more-to-come:

t = 1.7300000000000013 t = 1.7400000000000013 t = 1.7500000000000013

</pre>

t = 1.7300000000000013 t = 1.7400000000000013 t = 1.7500000000000013

end{sphinxVerbatim}

t = 1.7300000000000013 t = 1.7400000000000013 t = 1.7500000000000013

more-to-come:

t = 1.7600000000000013 t = 1.7700000000000014

</pre>

t = 1.7600000000000013 t = 1.7700000000000014

end{sphinxVerbatim}

t = 1.7600000000000013 t = 1.7700000000000014

more-to-come:

t = 1.7800000000000014 t = 1.7900000000000014

</pre>

t = 1.7800000000000014 t = 1.7900000000000014

end{sphinxVerbatim}

t = 1.7800000000000014 t = 1.7900000000000014

more-to-come:

t = 1.8000000000000014 t = 1.8100000000000014

</pre>

t = 1.8000000000000014 t = 1.8100000000000014

end{sphinxVerbatim}

t = 1.8000000000000014 t = 1.8100000000000014

more-to-come:

t = 1.8200000000000014 t = 1.8300000000000014 t = 1.8400000000000014

</pre>

t = 1.8200000000000014 t = 1.8300000000000014 t = 1.8400000000000014

end{sphinxVerbatim}

t = 1.8200000000000014 t = 1.8300000000000014 t = 1.8400000000000014

more-to-come:

t = 1.8500000000000014 t = 1.8600000000000014

</pre>

t = 1.8500000000000014 t = 1.8600000000000014

end{sphinxVerbatim}

t = 1.8500000000000014 t = 1.8600000000000014

more-to-come:

t = 1.8700000000000014 t = 1.8800000000000014 t = 1.8900000000000015

</pre>

t = 1.8700000000000014 t = 1.8800000000000014 t = 1.8900000000000015

end{sphinxVerbatim}

t = 1.8700000000000014 t = 1.8800000000000014 t = 1.8900000000000015

more-to-come:

t = 1.9000000000000015 t = 1.9100000000000015

</pre>

t = 1.9000000000000015 t = 1.9100000000000015

end{sphinxVerbatim}

t = 1.9000000000000015 t = 1.9100000000000015

more-to-come:

t = 1.9200000000000015 t = 1.9300000000000015

</pre>

t = 1.9200000000000015 t = 1.9300000000000015

end{sphinxVerbatim}

t = 1.9200000000000015 t = 1.9300000000000015

more-to-come:

t = 1.9400000000000015 t = 1.9500000000000015

</pre>

t = 1.9400000000000015 t = 1.9500000000000015

end{sphinxVerbatim}

t = 1.9400000000000015 t = 1.9500000000000015

more-to-come:

t = 1.9600000000000015 t = 1.9700000000000015

</pre>

t = 1.9600000000000015 t = 1.9700000000000015

end{sphinxVerbatim}

t = 1.9600000000000015 t = 1.9700000000000015

more-to-come:

t = 1.9800000000000015 t = 1.9900000000000015

</pre>

t = 1.9800000000000015 t = 1.9900000000000015

end{sphinxVerbatim}

t = 1.9800000000000015 t = 1.9900000000000015

t = 2.0000000000000013CPU times: user 11.2 s, sys: 28.2 ms, total: 11.2 s

Wall time: 11.2 s </pre>

t = 2.0000000000000013CPU times: user 11.2 s, sys: 28.2 ms, total: 11.2 s

Wall time: 11.2 s end{sphinxVerbatim}

t = 2.0000000000000013CPU times: user 11.2 s, sys: 28.2 ms, total: 11.2 s

Wall time: 11.2 s

Tasks:

  • Implement higher order Splitting (e.g. Strang) schemes

  • Implement IMEX schemes and compare