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: Curved mesh creation sometimes fails in parallel #28

Closed
ABaierReinio opened this issue Apr 30, 2024 · 1 comment · Fixed by #29
Closed

Firedrake-Netgen: Curved mesh creation sometimes fails in parallel #28

ABaierReinio opened this issue Apr 30, 2024 · 1 comment · Fixed by #29

Comments

@ABaierReinio
Copy link
Collaborator

MFE below:

from firedrake import *
from firedrake.output import VTKFile
from netgen.occ import *

box = Box(Pnt(0, 0, 0), Pnt(1, 1, 1))
body = box.MakeFillet(box.edges, 0.1)
geo = OCCGeometry(body)
ngmesh = geo.GenerateMesh(maxh=0.1)
mesh = Mesh(Mesh(ngmesh).curve_field(2))
VTKFile("mesh.pvd").write(mesh)

Running this code in serial, the mesh is created with no issues.

But running it in parallel (e.g. mpiexec -n 2 python3 test.py) results in the following warning:

firedrake:WARNING Not able to curve element 4724, residual is: 0.18171646033208677

So it seems some of the elements are not being curved in this case.

UZerbinati pushed a commit that referenced this issue May 2, 2024
This fix #28
Requires #MR on the netgen (already opend on the TU Gitlab) waiting for approval.

Signed-off-by: Umberto Zerbinati <zerbinati@maths.ox.ac.uk>
@UZerbinati UZerbinati mentioned this issue May 2, 2024
@UZerbinati
Copy link
Collaborator

The fix is two-parted. You need to use PR #29 and the Netgen MR mirrored here.
Also, you need to change your code:

from firedrake import *
from firedrake.output import VTKFile
from netgen.occ import *
import netgen

if COMM_WORLD.rank == 0:
    box = Box(Pnt(0, 0, 0), Pnt(1, 1, 1))
    body = box.MakeFillet(box.edges, 0.1)
    geo = OCCGeometry(body)
    ngmesh = geo.GenerateMesh(maxh=0.1)
else:
    ngmesh = netgen.libngpy._meshing.Mesh(3)
mesh = Mesh(Mesh(ngmesh).curve_field(2))
VTKFile("mesh.pvd").write(mesh)

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Labels
None yet
Projects
None yet
Development

Successfully merging a pull request may close this issue.

2 participants