← Back to team overview

dolfin team mailing list archive

Re: [HG DOLFIN] First iteration of refinement of function spaces. Refinement of a function

 

On Fri, Oct 30, 2009 at 08:06:27PM +0100, Marie Rognes wrote:
> DOLFIN wrote:
>> One or more new changesets pushed to the primary dolfin repository.
>> A short summary of the last three changesets is included below.
>>
>> changeset:   7408:f32a6ec49c47
>> tag:         tip
>> user:        Anders Logg <logg@xxxxxxxxx>
>> date:        Thu Oct 29 23:28:47 2009 +0100
>> files:       dolfin/fem/DofMap.cpp dolfin/fem/DofMap.h dolfin/function/FunctionSpace.cpp dolfin/function/FunctionSpace.h dolfin/mesh/Mesh.h
>> description:
>> First iteration of refinement of function spaces. Refinement of a function
>> space should take care of refining the mesh, rebuilding the dofmap and
>> interpolating all its member functions to the new function space.
>>
>> All essential bits should be in place but needs some debugging to
>> find memory errors/leaks.
>>
>>
>
> This seems to be not working quite optimally yet.
>
> 1) The attached script fails with the error:
>
>
>    python: dolfin/function/Function.cpp:392: virtual void
>    dolfin::Function::compute_vertex_values(double*, const
>    dolfin::Mesh&) const: Assertion `&mesh == &_function_space->mesh()'
>    failed.
>    Aborted
>

I'll take a look, but probably not until after the weekend.

> 2) The original mesh does not seem to be refined. Shouldn't it be?

Yes it should. That's a bug.

> 3) Not much seems to happen if f is an Expression (since not much seems
> to happen with the mesh)
>
> 4) What happens when there are multiple function spaces associated with
> one mesh?

That requires some thought.

--
Anders



> from dolfin import *
>
> mesh = UnitSquare(5, 5)
> print "Old num cells = ", len(mesh.cells())
>
> # Mark cells for refinement
> p = Point(0.5, 0.5)
> cell_markers = MeshFunction("bool", mesh, mesh.topology().dim())
> for c in cells(mesh):
>     if c.midpoint().distance(p) < 0.2:
>         cell_markers[c] = True
>     else:
>         cell_markers[c] = False
>
> V = FunctionSpace(mesh, "CG", 1)
>
> f = Function(V)
> for i in range(len(f.vector().array())):
>     f.vector()[i] = i+1
> #f = Expression("sin(pi*x[0])", V=V)
>
> plot(f, mesh=mesh, title="Old function")
>
> V.refine(cell_markers)
> print "New num cells = ", len(mesh.cells())
>
> plot(f, mesh=mesh, title="New function")
> interactive()

> _______________________________________________
> DOLFIN-dev mailing list
> DOLFIN-dev@xxxxxxxxxx
> http://www.fenics.org/mailman/listinfo/dolfin-dev

Attachment: signature.asc
Description: Digital signature


References