On Tuesday August 9 2011 09:27:50 Anders Logg wrote: > On Tue, Aug 09, 2011 at 08:50:39AM -0700, Johan Hake wrote: > > Thanks Mikael! > > > > That worked after I scaled up the stabilizing factor for a larger mesh. > > See attached image :) > > Consider adding, and documenting, a demo for the Eikonal equation!
Will do! Consider I mean ;) Johan > -- > Anders > > > Even if this is cool, I am planing to extend the CGAL wrapper so a user > > can compute the closest distance between a vertex and the boundary. This > > might be faster, and less mesh dependent. > > > > I am also hoping I can make this Boundary marker dependent. Given a > > FacetFunction I should be able to calculate the closest distance to a > > subset of the boundary which is nice for a particular application I have > > in mind :) > > > > Cheers! > > > > Johan > > > > On Tuesday August 9 2011 01:04:11 Mikael Mortensen wrote: > > > Hi > > > > > > I've been using the signed distance function (Eikonal equation) for > > > some time > > > in various turbulence models. This equation is very stiff so in my > > > approach I > > > have found it necessary to add a stabilization term. I also solve a > > > simpler linear Poisson problem first to get a good initial guess for > > > the nonlinear Newton solver. Otherwise it's pretty straightforward. > > > You could probably use the > > > NonlinearVariational etc classes to wrap this up nicely, but I haven't > > > done that. > > > The code I use is basically (copy and pasted from > > > http://bazaar.launchpad.net/~cbc.rans/cbc.rans/mikael/view/head:/cbc/ra > > > ns/t urbsolvers/Eikonal.py): > > > > > > V = FunctionSpace(mesh, 'CG', 1) > > > v = TestFunction(V) > > > u = TrialFunction(V) > > > f = Constant(1.0) > > > y = Function(V) > > > > > > # Initialization problem to get good initial guess for nonlinear > > > problem: F1 = inner(grad(u), grad(v))*dx - f*v*dx > > > a1, L1 = lhs(F1), rhs(F1) > > > A1, b1 = assemble_system(a1, L1) > > > # Apply boundary conditions: DirichletBC = 0 on the boundary > > > for bc in bcs: bc.apply(A1, b1) > > > solve(A1, y.vector(), b1) > > > > > > # Stabilized Eikonal equation > > > eps = Constant(0.01) > > > F = sqrt(inner(grad(y), grad(y)))*v*dx - f*v*dx + eps*inner(grad(y), > > > grad(v))*dx > > > > > > J = derivative(F, y, u) > > > etc > > > Solve with Newton iterations > > > > > > This works for me for a range of geometries. Hope it is helpful. > > > > > > Best regards > > > > > > Mikael > > > > > > On 8 August 2011 22:31, Johan Hake <[email protected]> wrote: > > > > Hello! > > > > > > > > Has anyone in the FEniCS comunity used FEniCS to solve for the signed > > > > distance > > > > function? > > > > > > > > http://en.wikipedia.org/wiki/Signed_distance > > > > > > > > This functions gives the distance to the boundary for any given point > > > > inside > > > > the domain. I would like to compute this for the vertices of a mesh. > > > > It looks > > > > like it is a nontrivial problem and I wonder if any of you have put > > > > your teeth > > > > into this one. > > > > > > > > Johan > > > > > > > > _______________________________________________ > > > > Mailing list: https://launchpad.net/~dolfin > > > > Post to : [email protected] > > > > Unsubscribe : https://launchpad.net/~dolfin > > > > More help : https://help.launchpad.net/ListHelp > > > > _______________________________________________ > > Mailing list: https://launchpad.net/~dolfin > > Post to : [email protected] > > Unsubscribe : https://launchpad.net/~dolfin > > More help : https://help.launchpad.net/ListHelp _______________________________________________ Mailing list: https://launchpad.net/~dolfin Post to : [email protected] Unsubscribe : https://launchpad.net/~dolfin More help : https://help.launchpad.net/ListHelp

