Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Firedrake-Netgen: Boundary IDs are not loaded correctly in parallel, when loading a curved mesh via Checkpointing #27

Open
ABaierReinio opened this issue Apr 26, 2024 · 0 comments
Labels
bug Something isn't working

Comments

@ABaierReinio
Copy link
Collaborator

When a curved firedrake-netgen mesh is created in serial, if you reload it in parallel via Checkpointing, then the boundary IDs are not labelled correctly. See the MFE below.

from firedrake import *
from firedrake.output import VTKFile

# Load the mesh and apply DirichletBCs (run this in parallel)
import os.path
if os.path.exists("out/mesh.h5"):

    with CheckpointFile("out/mesh.h5", 'r') as file:
        mesh = file.load_mesh("firedrake_default")

    V = FunctionSpace(mesh, "CG", 1)
    U = Function(V)

    DirichletBC(V, -10.0, 1).apply(U)
    DirichletBC(V, 10.0, 2).apply(U)

    VTKFile("out/func.pvd").write(U)

# Build the mesh and save it (run this in serial)
else:

    from netgen.occ import *

    box = Box(Pnt(0,0,0), Pnt(1,1,1))
    body = box.MakeFillet(box.edges, 0.1)

    body.faces.Min(X).bc("left")
    body.faces.Max(X).bc("right")

    geo = OCCGeometry(body)
    ngmesh = geo.GenerateMesh(maxh=0.15)

    print(ngmesh.GetRegionNames(dim=2))
    label_flags = {"labels" : {"left" : 1, "right" : 2, "default" : 3, "fillet" : 4}}
    mesh = Mesh(Mesh(ngmesh, netgen_flags=label_flags).curve_field(2))

    # Important: If one instead generates the mesh via the code commented below, no problems arise!
    '''
    mesh = UnitCubeMesh(10, 10, 10)

    V = VectorFunctionSpace(mesh, "CG", 4)
    newcoords = assemble(interpolate(mesh.coordinates + as_vector([0.25 * sin(np.pi*mesh.coordinates[1]), 0, 0]), V))

    mesh = Mesh(newcoords)
    '''

    with CheckpointFile("out/mesh.h5", "w") as file:
        file.save_mesh(mesh)

Running the above code twice (the first time in serial, and the second time in parallel), results in the boundary conditions being applied incorrectly, as shown below:

bad

On the other hand, if you run it the second time in serial (i.e. you load the mesh in serial), or you don't curve the mesh
(i.e. you do mesh = Mesh(ngmesh, netgen_flags=label_flags) ) then the boundary conditions are correctly applied.
Furthermore, if you instead make a curved mesh in vanilla Firedrake (see the commented code above), then the boundary conditions are again correctly applied, even when loading the mesh in parallel.

good

Note also that, if you instead make a "curved" linear mesh (i.e. you do mesh = Mesh(Mesh(ngmesh, netgen_flags=label_flags).curve_field(1)) ) then you still get the problematic behaviour with boundary conditions being applied incorrectly.

@UZerbinati UZerbinati added the bug Something isn't working label Jun 4, 2024
Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Labels
bug Something isn't working
Projects
None yet
Development

No branches or pull requests

2 participants