Re: [firedrake] Nullspace & VectorSpaceBasis
On 11 Dec 2014 19:17, "Anna Kalogirou" <a.kalogirou@leeds.ac.uk> wrote:
Thanks Andrew. About your comment 1), how can I fix this?
I think that requires a Firedrake (PyOP2) fix. Lawrence?
The two equations valid on one of the boundaries can actually be combined
into one by paying the price of a higher time derivative. This is not a problem for this simple system, so if I use
aphi = v*phi*ds(4) + g*dt**2*inner(grad(v),grad(phi))*dx Lphi = v*(2*phi1-phi0)*ds(4)
phi_problem = LinearVariationalProblem(aphi,Lphi,phi2) phi_solver = LinearVariationalSolver(phi_problem)
is equivalent to what I had before, and this time the code works and
gives a good result (compared to the exact solution). Cool...
However, I need to understand why the other case doesn't work because
that will later be generalised to become nonlinear.
Anna.
On 11/12/14 13:04, Andrew McRae wrote:
Four comments:
1) If I understand correctly, PETSc insists that the diagonal entries
of the sparse matrix are filled in. Usually, the LHS form will cause this to happen; however, that obviously isn't the case for your form!
2) It should be possible to build the nullspace you want, but I'm
clueless on the actual syntax. This wouldn't get around (1) without further fixes from our end.
3) Could the equations given not be included in the weak form for the
rest of the problem?
4) If not, there might be a hacky
"add-a-weighted-mass-matrix-to-the-LHS" fix you could use?....
Andrew
On 11 December 2014 at 12:46, Anna Kalogirou <a.kalogirou@leeds.ac.uk>
wrote:
Dear all,
I am trying to solve a simple problem (linear water wave equations), where I have the Laplace equation in a square domain with Neumann boundary conditions on the three sides and on the top boundary two equations are satisfied.
I have the following code: -------------------------------------- V = FunctionSpace(mesh,"CG",1)
eta0 = Function(V).interpolate(Expression("cos(2*pi*x[0])")) phi0 = Function(V) eta1 = Function(V) phi1 = Function(V) eta = TrialFunction(V) phi = TrialFunction(V) v = TestFunction(V)
aeta = v*eta*ds(4) Leta = v*eta0*ds(4) + dt*inner(grad(v),grad(phi0))*dx
eta_problem = LinearVariationalProblem(aeta,Leta,eta1) eta_solver = LinearVariationalSolver(eta_problem)
aphi = v*phi*ds(4) Lphi = (v*phi0 - g*dt*v*eta1)*ds(4)
phi_problem = LinearVariationalProblem(aphi,Lphi,phi1) phi_solver = LinearVariationalSolver(phi_problem) --------------------------------------
The error I get tells me that Object is in wrong state, and Matrix is missing diagonal entry 39.
I thought this had to do with the problem being singular, so I tried to build the nullspace but I don't know how to define the VectorSpaceBasis. Any help would be appreciated.
Thank you very much.
Regards, Anna.
_______________________________________________ firedrake mailing list firedrake@imperial.ac.uk https://mailman.ic.ac.uk/mailman/listinfo/firedrake
--
Dr Anna Kalogirou Research Fellow School of Mathematics University of Leeds
participants (1)
- 
                
                Andrew McRae