Hi, I intend to implement a system of PDEs (4×4) in Firedeake package, which the function spaces are mix of tensor, vector and scalar. For this reason I prefer to understand the code of mixed formulation for Poisson equation completly. In the mixed formulation code of Poisson equation, when I change the function space from "FunctionSpace(mesh, "BDM", 1)" to "VectorFunctionSpace(mesh, "CG", 1)", there is a bug in solvig the system. I implement the following code: mesh = UnitSquareMesh(10, 10) VCG = VectorFunctionSpace(mesh, "CG", 1) DG = FunctionSpace(mesh, "DG", 0) W = VCG * DG sigma, u = TrialFunctions(W) tau, v = TestFunctions(W) x, y = SpatialCoordinate(mesh) f = Function(DG).interpolate( 10*exp(-(pow(x - 0.5, 2) + pow(y - 0.5, 2)) / 0.02)) a = (dot(sigma, tau) + div(tau)*u + div(sigma)*v)*dx L = - f*v*dx bc0 = DirichletBC(W.sub(0), as_vector([0.0, -sin(5*x)]), 1) bc1 = DirichletBC(W.sub(0), as_vector([0.0, sin(5*y)]), 2) w = Function(W) solve(a == L, w, bcs=[bc0, bc1]) sigma, u = w.split() File("poisson_mixed.pvd").write(u) try: import matplotlib.pyplot as plt except: warning("Matplotlib not imported") try: plot(u) except Exception as e: warning("Cannot plot figure. Error msg '%s'" % e) try: plt.show() except Exception as e: warning("Cannot show figure. Error msg '%s'" % e) f.interpolate(cos(x*pi*2)*cos(y*pi*2)) print(sqrt(assemble(dot(u - f, u - f) * dx))) and I get the following error: Traceback (most recent call last): File "test2.py", line 156, in <module> solve(a == L, w, bcs=[bc0, bc1]) File "/home/amireh/firedrake/src/firedrake/firedrake/solving.py", line 121, in solve _solve_varproblem(*args, **kwargs) File "/home/amireh/firedrake/src/firedrake/firedrake/solving.py", line 151, in _solve_varproblem solver.solve() File "/home/amireh/firedrake/src/firedrake/firedrake/variational_solver.py", line 242, in solve solving_utils.check_snes_convergence(self.snes) File "/home/amireh/firedrake/src/firedrake/firedrake/solving_utils.py", line 231, in check_snes_convergence %s""" % (snes.getIterationNumber(), msg)) firedrake.exceptions.ConvergenceError: Nonlinear solve failed to converge after 0 nonlinear iterations. Reason: Inner linear solve failed to converge after 10000 iterations with reason: DIVERGED_MAX_IT
Hi Amireh, First of all, changing the flux from BDM to CG1 gives you a P1-P0 Lagrange formulation for the mixed Poisson problem which is unstable and will give you an oscillatory solution. You need to either add some stabilization terms to the variational formulation or increase the polynomial order of the vector function space. Second, if you insist on using Lagrange elements for velocity, you need to apply the BCs appropriately. You would need to do something like: bcs = DirichletBC(W.sub(0).sub(0), -sin(5*x), 1) Third, the default solver/preconditioned mixed systems is GMRES/Jacobi which is a terrible solver and won’t converge. Here are some ways to solve this system of equations instead https://firedrakeproject.org/demos/saddle_point_systems.py.html Hope this helps, Justin On Sat, Jun 23, 2018 at 9:00 AM Amireh Mousavi < amireh.mousavi@math.iut.ac.ir> wrote:
Hi, I intend to implement a system of PDEs (4×4) in Firedeake package, which the function spaces are mix of tensor, vector and scalar. For this reason I prefer to understand the code of mixed formulation for Poisson equation completly. In the mixed formulation code of Poisson equation, when I change the function space from *"FunctionSpace(mesh, "BDM", 1)" *to *"VectorFunctionSpace(mesh, "CG", 1)", *there is a bug in solvig the system. I implement the following code:
mesh = UnitSquareMesh(10, 10) VCG = VectorFunctionSpace(mesh, "CG", 1) DG = FunctionSpace(mesh, "DG", 0) W = VCG * DG sigma, u = TrialFunctions(W) tau, v = TestFunctions(W) x, y = SpatialCoordinate(mesh) f = Function(DG).interpolate( 10*exp(-(pow(x - 0.5, 2) + pow(y - 0.5, 2)) / 0.02)) a = (dot(sigma, tau) + div(tau)*u + div(sigma)*v)*dx L = - f*v*dx bc0 = DirichletBC(W.sub(0), as_vector([0.0, -sin(5*x)]), 1) bc1 = DirichletBC(W.sub(0), as_vector([0.0, sin(5*y)]), 2) w = Function(W) solve(a == L, w, bcs=[bc0, bc1]) sigma, u = w.split() File("poisson_mixed.pvd").write(u) try: import matplotlib.pyplot as plt except: warning("Matplotlib not imported")
try: plot(u) except Exception as e: warning("Cannot plot figure. Error msg '%s'" % e) try: plt.show() except Exception as e: warning("Cannot show figure. Error msg '%s'" % e) f.interpolate(cos(x*pi*2)*cos(y*pi*2)) print(sqrt(assemble(dot(u - f, u - f) * dx)))
and I get the following error:
Traceback (most recent call last): File "test2.py", line 156, in <module> solve(a == L, w, bcs=[bc0, bc1]) File "/home/amireh/firedrake/src/firedrake/firedrake/solving.py", line 121, in solve _solve_varproblem(*args, **kwargs) File "/home/amireh/firedrake/src/firedrake/firedrake/solving.py", line 151, in _solve_varproblem solver.solve() File "/home/amireh/firedrake/src/firedrake/firedrake/variational_solver.py", line 242, in solve solving_utils.check_snes_convergence(self.snes) File "/home/amireh/firedrake/src/firedrake/firedrake/solving_utils.py", line 231, in check_snes_convergence %s""" % (snes.getIterationNumber(), msg)) firedrake.exceptions.ConvergenceError: Nonlinear solve failed to converge after 0 nonlinear iterations. Reason: Inner linear solve failed to converge after 10000 iterations with reason: DIVERGED_MAX_IT
_______________________________________________ firedrake mailing list firedrake@imperial.ac.uk https://mailman.ic.ac.uk/mailman/listinfo/firedrake
participants (2)
- 
                
                Amireh Mousavi
- 
                
                Justin Chang