FEniCS / dolfinx

Next generation FEniCS problem solving environment
https://fenicsproject.org
GNU Lesser General Public License v3.0
696 stars 172 forks source link

[BUG]: cannot export several components of function with VTX #3269

Open RemDelaporteMathurin opened 3 weeks ago

RemDelaporteMathurin commented 3 weeks ago

Summarize the issue

I'm trying to export several components of a vector function to VTX. However, I notice that the components are not exported correctly when opening with Paraview.

image

Here I would expect two different fields.

How to reproduce the bug

Run the following MWE

Minimal Example (Python)

import dolfinx
from mpi4py import MPI
import basix

mesh = dolfinx.mesh.create_unit_square(MPI.COMM_WORLD, 10, 10)

element = basix.ufl.element(
    "CG",
    mesh.basix_cell(),
    1,
    basix.LagrangeVariant.equispaced,
)
W = dolfinx.fem.functionspace(mesh, basix.ufl.mixed_element([element, element]))

u = dolfinx.fem.Function(W)

u.sub(0).interpolate(lambda x: x[0] ** 2 + x[1] ** 2)
u.sub(1).interpolate(lambda x: (1 - x[0]) ** 2 + (1 - x[1]) ** 2)

writer = dolfinx.io.VTXWriter(
    MPI.COMM_WORLD,
    "results_vector.bp",
    [u.sub(0), u.sub(1)],
    "BP4",
)
writer.write(0.0)

Output (Python)

No response

Version

0.8.0

DOLFINx git commit

No response

Installation

conda install fenics-dolfinx=0.8

Additional information

No response

RemDelaporteMathurin commented 3 weeks ago

When exporting two different functions, this works correctly:

import dolfinx
from mpi4py import MPI
import basix

mesh = dolfinx.mesh.create_unit_square(MPI.COMM_WORLD, 10, 10)
V = dolfinx.fem.functionspace(mesh, ("DG", 1))
u = dolfinx.fem.Function(V)
u.name = "u"
u.interpolate(lambda x: x[0] ** 2 + x[1] ** 2)

v = dolfinx.fem.Function(V)
v.name = "v"
v.interpolate(lambda x: (1 - x[0]) ** 2 + (1 - x[1]) ** 2)

writer = dolfinx.io.VTXWriter(
    MPI.COMM_WORLD,
    "results.bp",
    [u, v],
    "BP4",
)

writer.write(0.0)

image

bleyerj commented 3 weeks ago

Did you try with collapsing the subfunctions ? u.sub(0).collapse()

RemDelaporteMathurin commented 3 weeks ago

@bleyerj the following worked, however I have to rename the functions in order to open it in Paraview:

u_sub_0 = u.sub(0).collapse()
u_sub_0.name = "u_sub_0"

u_sub_1 = u.sub(1).collapse()
u_sub_1.name = "u_sub_1"

writer = dolfinx.io.VTXWriter(
    MPI.COMM_WORLD,
    "results_vector.bp",
    [u_sub_0, u_sub_1],
    "BP4",
)
writer.write(0.0)

image

RemDelaporteMathurin commented 3 weeks ago

@bleyerj the only issue with this is that you cannot write time dependent solution. Consider:

u.sub(0).interpolate(lambda x: x[0] ** 2 + x[1] ** 2)
u.sub(1).interpolate(lambda x: (1 - x[0]) ** 2 + (1 - x[1]) ** 2)

u_sub_0 = u.sub(0).collapse()
u_sub_0.name = "u_sub_0"

u_sub_1 = u.sub(1).collapse()
u_sub_1.name = "u_sub_1"

writer = dolfinx.io.VTXWriter(
    MPI.COMM_WORLD,
    "results_vector.bp",
    [u_sub_0, u_sub_1],
    "BP4",
)
writer.write(0.0)
u.sub(0).interpolate(lambda x: 3 * x[0] ** 2 + x[1] ** 2)
writer.write(1.0)

The .sub(0) component isn't modified.

jorgensd commented 3 weeks ago

So I would do this slightly differently:

from mpi4py import MPI
import dolfinx

mesh = dolfinx.mesh.create_unit_square(
    MPI.COMM_WORLD, 10, 10, dolfinx.cpp.mesh.CellType.triangle
)
V = dolfinx.fem.functionspace(mesh, ("Lagrange", 1, (2,)))
u = dolfinx.fem.Function(V)
V0, v0_to_V = V.sub(0).collapse()
V1, v1_to_V = V.sub(1).collapse()

u.sub(0).interpolate(lambda x: x[0] ** 2 + x[1] ** 2)
u.sub(1).interpolate(lambda x: (1 - x[0]) ** 2 + (1 - x[1]) ** 2)

u_sub_0 = u.sub(0).collapse()
u_sub_0.name = "u_sub_0"

u_sub_1 = u.sub(1).collapse()
u_sub_1.name = "u_sub_1"

writer = dolfinx.io.VTXWriter(
    MPI.COMM_WORLD,
    "results_vector.bp",
    [u_sub_0, u_sub_1],
    "BP4",
)
writer.write(0.0)
u.sub(0).interpolate(lambda x: 3 * x[0] ** 2 + x[1] ** 2)
u_sub_0.x.array[:] = u.x.array[v0_to_V]
writer.write(1.0)
writer.close()

It is always beneficial to use the maps from the collapsed space to parent and its inverse.

RemDelaporteMathurin commented 3 weeks ago

It is always beneficial to use the maps from the collapsed space to parent and its inverse.

Thanks for the tip!

Is there a particular reason why giving u.sub(0) to the writer would work fine but not [u.sub(0), u.sub(1)], that looks more streamlined that having to create a collapsed function and do u_sub_0.x.array[:] = u.x.array[v0_to_V] everytime

jorgensd commented 2 weeks ago

It is always beneficial to use the maps from the collapsed space to parent and its inverse.

Thanks for the tip!

Is there a particular reason why giving u.sub(0) to the writer would work fine but not [u.sub(0), u.sub(1)], that looks more streamlined that having to create a collapsed function and do u_sub_0.x.array[:] = u.x.array[v0_to_V] everytime

The problem is that collapsed dofmaps aren’t alwaysed ordered the same (wrt ghosts). I thought this was fixed in: https://github.com/FEniCS/dolfinx/issues/2929