-
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
Bug in locate_subdomain_entities
#908
Comments
RemDelaporteMathurin
added
bug
Something isn't working
fenicsx
Issue that is related to the fenicsx support
labels
Oct 31, 2024
I guess this works with the mixed domain approach since each subdomain has its own submesh.... this will not work with |
Actually running this with HydrogenTransportMaterial with several materials crashes: from mpi4py import MPI
import dolfinx
import dolfinx.fem.petsc
import numpy as np
import festim as F
def generate_mesh(n=20):
def bottom_boundary(x):
return np.isclose(x[1], 0.0)
def top_boundary(x):
return np.isclose(x[1], 1.0)
def half(x):
return x[1] <= 0.5 + 1e-14
mesh = dolfinx.mesh.create_unit_square(
MPI.COMM_WORLD, n, n, dolfinx.mesh.CellType.triangle
)
# Split domain in half and set an interface tag of 5
gdim = mesh.geometry.dim
tdim = mesh.topology.dim
fdim = tdim - 1
top_facets = dolfinx.mesh.locate_entities_boundary(mesh, fdim, top_boundary)
bottom_facets = dolfinx.mesh.locate_entities_boundary(mesh, fdim, bottom_boundary)
num_facets_local = (
mesh.topology.index_map(fdim).size_local
+ mesh.topology.index_map(fdim).num_ghosts
)
facets = np.arange(num_facets_local, dtype=np.int32)
values = np.full_like(facets, 0, dtype=np.int32)
values[top_facets] = 1
values[bottom_facets] = 2
bottom_cells = dolfinx.mesh.locate_entities(mesh, tdim, half)
num_cells_local = (
mesh.topology.index_map(tdim).size_local
+ mesh.topology.index_map(tdim).num_ghosts
)
cells = np.full(num_cells_local, 4, dtype=np.int32)
cells[bottom_cells] = 3
ct = dolfinx.mesh.meshtags(
mesh, tdim, np.arange(num_cells_local, dtype=np.int32), cells
)
all_b_facets = dolfinx.mesh.compute_incident_entities(
mesh.topology, ct.find(3), tdim, fdim
)
all_t_facets = dolfinx.mesh.compute_incident_entities(
mesh.topology, ct.find(4), tdim, fdim
)
interface = np.intersect1d(all_b_facets, all_t_facets)
values[interface] = 5
mt = dolfinx.mesh.meshtags(mesh, mesh.topology.dim - 1, facets, values)
return mesh, mt, ct
mesh, mt, ct = generate_mesh(10)
my_model = F.HydrogenTransportProblemDiscontinuousChangeVar()
my_model = F.HydrogenTransportProblem()
my_model.mesh = F.Mesh(mesh)
my_model.volume_meshtags = ct
my_model.facet_meshtags = mt
material_bottom = F.Material(D_0=1, E_D=0, K_S_0=1, E_K_S=0)
material_top = F.Material(D_0=1, E_D=0, K_S_0=1, E_K_S=0)
top_domain = F.VolumeSubdomain(4, material=material_top)
bottom_domain = F.VolumeSubdomain(3, material=material_bottom)
top_surface = F.SurfaceSubdomain(id=1)
bottom_surface = F.SurfaceSubdomain(id=2)
my_model.subdomains = [
bottom_domain,
top_domain,
top_surface,
bottom_surface,
]
H = F.SpeciesChangeVar("H", mobile=True)
my_model.species = [H]
my_model.boundary_conditions = [
F.FixedConcentrationBC(top_surface, value=1, species=H),
F.FixedConcentrationBC(bottom_surface, value=0, species=H),
]
my_model.temperature = 500.0
my_model.settings = F.Settings(atol=1e-10, rtol=1e-10, transient=False)
max_top = F.MaximumVolume(field=H, volume=top_domain)
max_bottom = F.MaximumVolume(field=H, volume=bottom_domain)
my_model.exports = [max_top, max_bottom]
my_model.initialise()
my_model.run()
assert max_top.value != max_bottom.value Produces:
|
Sign up for free
to join this conversation on GitHub.
Already have an account?
Sign in to comment
FESTIM/src/festim/subdomain/volume_subdomain.py
Lines 32 to 43 in 013937a
This will always return all the cells of the mesh, not the subdomain cells!!
We should use meshtags.find()
The text was updated successfully, but these errors were encountered: