spectralDNS / shenfun

High performance computational platform in Python for the spectral Galerkin method

Home Page:http://shenfun.readthedocs.org

Geek Repo:Geek Repo

Github PK Tool:Github PK Tool

Two questions in the demo MixedPoisson.py

yoczhang opened this issue · comments

Dear Mikael,
In the MixedPoisson.py, I have two questions:

1 . How to modify it to periodic boundary conditions in both x- and y-directions. I tried to change SD = FunctionSpace(N, family, bc =(0, 0)) to SD = FunctionSpace(N, family), but an error occurred.

2 . For the mixed formulation
$$
g - \nabla(u) &= 0 \
\nabla \cdot g &= f \
$$
g is a vector denoted by g=(g0, g1), g and u are trial functions, let p=(p0, p1) and q are test functions.

If I don't want to use the VT = VectorSpace(TT) to get the space g belongs to, and then M = BlockMatrix(A00+A01+A10) to generate the final matrix, I want to get the final matrix (by myself) like this:

[ A0   O   B0' ]
[  O   A1  B1' ]
[ B0   B1   O  ]

where A0=(p0, g0), A1=(p1, g1), O is zero matrix, B0=(q, g1_x), B1=(q, g2_y).

Is this possible, and if so, can a complete code or demo be provided?

Thank you very much!

Hi
Very good questions! I have not used a purely Fourier composite space before, but there is no reason why it should not work! Unfortunately, I need to add some code to the BlockMatrix and the BlockMatrixSolver for this to work. I will add it to master and push as soon as possible.

Thanks!

It should work now. Following code runs fine for me:

import os
import sys
import numpy as np
from sympy import symbols, sin, cos
from shenfun import *

x, y = symbols("x,y", real=True)

# Create a manufactured solution for verification
ue = sin(2*x)*cos(1*y)
dux = ue.diff(x, 1)
duy = ue.diff(y, 1)
fe = ue.diff(x, 2) + ue.diff(y, 2)

N = 24
K0 = FunctionSpace(N, 'Fourier', dtype='D')
SD = FunctionSpace(N, 'Fourier', dtype='d')
TD = TensorProductSpace(comm, (K0, SD))
VT = VectorSpace(TD)
Q = CompositeSpace([VT, TD])
gu = TrialFunction(Q)
pq = TestFunction(Q)
g, u = gu
p, q = pq
A00 = inner(p, g)
A01 = inner(div(p), u)
A10 = inner(q, div(g))
# Get f and g on quad points
vfj = Array(Q, buffer=(0, 0, fe))
vj, fj = vfj
vf_hat = Function(Q)
v_hat, f_hat = vf_hat
f_hat = inner(q, fj, output_array=f_hat)
M = BlockMatrix(A00+A01+A10)
gu_hat = M.solve(vf_hat, constraints=((2, 0, 0),))
gu = gu_hat.backward()
g_, u_ = gu
uj = Array(TD, buffer=ue)
duxj = Array(TD, buffer=dux)
duyj = Array(TD, buffer=duy)
error = [np.sqrt(inner(1, (uj-u_)**2)),
         np.sqrt(inner(1, (duxj-g_[0])**2)),
         np.sqrt(inner(1, (duyj-g_[1])**2))]
import matplotlib.pyplot as plt
plt.figure()
X = TD.local_mesh(True)
plt.contourf(X[0], X[1], u_)
plt.figure()
plt.quiver(X[1], X[0], g_[1], g_[0])
plt.figure()
plt.spy(M.diags(0, format='csr').toarray()) 
plt.show()

Yes, it works for me.

@mikaem Hi Mikael and other question about all periodic bcs, I don't quite understand the choice of solver, if I want to solve the biharmonic2D problem with all periodic bcs in (x, y) directions, I modified part of the code according to the above, but got the error prompt again.

So may other examples (such as Poisson and biharmonic equations with all periodic bcs) be provided here?

Thanks!

Hi
The biharmonic2D problem is scalar and the assembled Fourier matrix will be diagonal. You should not try to modify biharmonic2D, but rather fourier_poisson2D, where you simply need to change the equation:

import os
from sympy import symbols, cos, sin, lambdify
import numpy as np
from shenfun import inner, div, grad, TestFunction, TrialFunction, Array, FunctionSpace, \
    TensorProductSpace, Function, comm, la

# Use sympy to compute a rhs, given an analytical solution
x, y = symbols("x,y", real=True)
ue = cos(4*x) + sin(8*y)

# Change right hand side to reflect biharmonic equation
fe = ue.diff(x, 4) + ue.diff(y, 4) + 2*ue.diff(x, 2, y, 2)

# Size of discretization
N = (64, 64)

K0 = FunctionSpace(N[0], family='F', dtype='D', domain=(-2*np.pi, 2*np.pi))
K1 = FunctionSpace(N[1], family='F', dtype='d', domain=(-2*np.pi, 2*np.pi))
T = TensorProductSpace(comm, (K0, K1), axes=(0, 1))
u = TrialFunction(T)
v = TestFunction(T)

# Get f on quad points
fj = Array(T, buffer=fe)

# Compute right hand side
f_hat = Function(T)
f_hat = inner(v, fj, output_array=f_hat)

# Solve biharmonic equation
u_hat = Function(T)
A = inner(div(grad(v)), div(grad(u)))
#A = inner(v, div(grad(div(grad(u)))))
sol = la.SolverDiagonal(A)
u_hat = sol(f_hat, u_hat, constraints=((0, 0),))

uq = Array(T)
uq = T.backward(u_hat, uq)
uj = Array(T, buffer=ue)
error = np.sqrt(inner(1, (uj-uq)**2))
assert abs(error) < 1e-6
if comm.Get_rank() == 0:
    print(f"fourier_biharmonic2D L2 error = {abs(error):2.6e}")

OK, thank you!