-
Notifications
You must be signed in to change notification settings - Fork 25
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
Multi-material issue on fenicsx
#898
Comments
This is the mesh I use |
If we comment out the |
The 1D equivalent works fine.... import numpy as np
import festim as F
import dolfinx
dolfinx.log.set_log_level(dolfinx.log.LogLevel.INFO)
my_model = F.HydrogenTransportProblem()
interface_1 = 0.2
interface_2 = 0.8
vertices = np.concatenate(
[
np.linspace(0, interface_1, num=100),
np.linspace(interface_1, interface_2, num=100),
np.linspace(interface_2, 1, num=100),
]
)
my_model.mesh = F.Mesh1D(vertices)
mat = F.Material(D_0=1, E_D=0)
vol1 = F.VolumeSubdomain1D(id=1, material=mat, borders=[0, interface_1])
vol2 = F.VolumeSubdomain1D(id=2, material=mat, borders=[interface_1, interface_2])
vol3 = F.VolumeSubdomain1D(id=3, material=mat, borders=[interface_2, 1])
surface1 = F.SurfaceSubdomain1D(id=4, x=vertices[0])
surface2 = F.SurfaceSubdomain1D(id=5, x=vertices[-1])
my_model.subdomains = [vol1, vol2, vol3, surface1, surface2]
mobile = F.Species(name="H", mobile=True)
trapped_1 = F.Species(name="H_trapped_W1", mobile=False)
empty_1a = F.ImplicitSpecies(n=0.5, others=[trapped_1])
my_model.species = [mobile, trapped_1]
my_model.reactions = [
F.Reaction(
reactant=[mobile, empty_1a],
product=[trapped_1],
k_0=1,
E_k=mat.E_D,
p_0=0.1,
E_p=0.87,
volume=vol,
)
for vol in my_model.volume_subdomains
]
my_model.temperature = 600
my_model.boundary_conditions = [
F.FixedConcentrationBC(subdomain=surface1, species=mobile, value=2),
F.FixedConcentrationBC(subdomain=surface2, species=mobile, value=0),
]
my_model.settings = F.Settings(
atol=1e-6,
rtol=1e-6,
transient=False,
)
my_model.initialise()
my_model.run()
|
I tried to make a 2D MWE with a dolfinx mesh but couldn't reproduce the issue with this: import dolfinx
import festim as F
import numpy as np
import mpi4py.MPI as MPI
def make_mesh(N=10):
mesh = dolfinx.mesh.create_unit_square(MPI.COMM_WORLD, N, N)
def top_half(x):
return x[1] > 0.5 - 1e-10
def bottom_half(x):
return x[1] < 0.5 + 1e-10
vdim = mesh.topology.dim
fdim = mesh.topology.dim - 1
num_cells = mesh.topology.index_map(vdim).size_local
mesh_cell_indices = np.arange(num_cells, dtype=np.int32)
tags_volumes = np.full(num_cells, 0, dtype=np.int32)
entities_top = dolfinx.mesh.locate_entities(mesh, mesh.topology.dim, top_half)
entities_bot = dolfinx.mesh.locate_entities(mesh, mesh.topology.dim, bottom_half)
tags_volumes[entities_bot] = 1
tags_volumes[entities_top] = 2
ct = dolfinx.mesh.meshtags(
mesh, dim=mesh.topology.dim, entities=mesh_cell_indices, values=tags_volumes
)
# find all facets in domain and mark them as 0
def top_surface(x):
return np.isclose(x[1], 1)
def bottom_surface(x):
return np.isclose(x[1], 0)
entities_top = dolfinx.mesh.locate_entities(mesh, fdim, top_surface)
entities_bot = dolfinx.mesh.locate_entities(mesh, fdim, bottom_surface)
# concatenate entities in one array and mark them as 1 and 2
entities = np.concatenate([entities_bot, entities_top])
tags_facets = np.concatenate(
[np.full_like(entities_bot, 4), np.full_like(entities_top, 5)]
)
ft = dolfinx.mesh.meshtags(mesh, dim=fdim, entities=entities, values=tags_facets)
return mesh, ct, ft
dolfinx.log.set_log_level(dolfinx.log.LogLevel.INFO)
my_model = F.HydrogenTransportProblem()
mesh, ct, ft = make_mesh(N=10)
# with dolfinx.io.XDMFFile(MPI.COMM_WORLD, "mt.xdmf", "w") as xdmf:
# xdmf.write_mesh(mesh)
# xdmf.write_meshtags(ct, mesh.geometry)
# with dolfinx.io.XDMFFile(MPI.COMM_WORLD, "ft.xdmf", "w") as xdmf:
# xdmf.write_mesh(mesh)
# xdmf.write_meshtags(ft, mesh.geometry)
my_model.mesh = F.Mesh(mesh=mesh)
my_model.volume_meshtags = ct
my_model.facet_meshtags = ft
mat = F.Material(D_0=1, E_D=0.1)
vol1 = F.VolumeSubdomain(id=1, material=mat)
vol2 = F.VolumeSubdomain(id=2, material=mat)
surface1 = F.SurfaceSubdomain(id=4)
surface2 = F.SurfaceSubdomain(id=5)
my_model.subdomains = [vol1, vol2, surface1, surface2]
mobile = F.Species(name="H", mobile=True)
trapped = F.Species(name="H_trapped", mobile=False)
empty = F.ImplicitSpecies(n=0.5, others=[trapped])
my_model.species = [mobile, trapped]
my_model.reactions = [
F.Reaction(
reactant=[mobile, empty],
product=[trapped],
k_0=1,
E_k=mat.E_D,
p_0=0.1,
E_p=0.87,
volume=vol,
)
for vol in my_model.volume_subdomains
]
my_model.temperature = 600
my_model.boundary_conditions = [
F.FixedConcentrationBC(subdomain=surface1, species=mobile, value=1),
F.FixedConcentrationBC(subdomain=surface2, species=mobile, value=0),
]
my_model.settings = F.Settings(
atol=1e-6,
rtol=1e-6,
transient=False,
)
my_model.initialise()
my_model.run() |
Setting the following in ksp = self.solver.krylov_solver
ksp.setType("preonly")
ksp.getPC().setType("lu")
ksp.getPC().setFactorSolverType("mumps")
ksp.setErrorIfNotConverged(True) |
Sign up for free
to join this conversation on GitHub.
Already have an account?
Sign in to comment
Run on the
fenicsx
branchThe text was updated successfully, but these errors were encountered: